1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 2111ade9eSPeter Brune #include <petscdm.h> 3eaedb033SPeter Brune 4eaedb033SPeter Brune typedef struct { 5eaedb033SPeter Brune PetscInt n; /* local subdomains */ 6eaedb033SPeter Brune SNES *subsnes; /* nonlinear solvers for each subdomain */ 7eaedb033SPeter Brune Vec *x; /* solution vectors */ 8111ade9eSPeter Brune Vec *xl; /* solution local vectors */ 9111ade9eSPeter Brune Vec *y; /* step vectors */ 10eaedb033SPeter Brune Vec *b; /* rhs vectors */ 11f10b3e88SPatrick Farrell Vec weight; /* weighting for adding updates on overlaps, in global space */ 12111ade9eSPeter Brune VecScatter *oscatter; /* scatter from global space to the subdomain global space */ 13f10b3e88SPatrick Farrell VecScatter *oscatter_copy; /* copy of the above */ 14111ade9eSPeter Brune VecScatter *iscatter; /* scatter from global space to the nonoverlapping subdomain space */ 15111ade9eSPeter Brune VecScatter *gscatter; /* scatter from global space to the subdomain local space */ 16111ade9eSPeter Brune PCASMType type; /* ASM type */ 17111ade9eSPeter Brune PetscBool usesdm; /* use the DM for setting up the subproblems */ 18d728fb7dSPeter Brune PetscBool finaljacobian; /* compute the jacobian of the converged solution */ 19610116beSPeter Brune PetscReal damping; /* damping parameter for updates from the blocks */ 20f10b3e88SPatrick Farrell PetscBool weight_set; /* use a weight in the overlap updates */ 21b20c023fSPeter Brune 22b20c023fSPeter Brune /* logging events */ 23b20c023fSPeter Brune PetscLogEvent eventrestrictinterp; 24b20c023fSPeter Brune PetscLogEvent eventsubsolve; 25602bec5dSPeter Brune 26602bec5dSPeter Brune PetscInt fjtype; /* type of computed jacobian */ 27602bec5dSPeter Brune Vec xinit; /* initial solution in case the final jacobian type is computed as first */ 28eaedb033SPeter Brune } SNES_NASM; 29eaedb033SPeter Brune 309e5d0892SLisandro Dalcin const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",NULL}; 31602bec5dSPeter Brune const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"}; 32b20c023fSPeter Brune 3340244768SBarry Smith static PetscErrorCode SNESReset_NASM(SNES snes) 34eaedb033SPeter Brune { 35eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 36eaedb033SPeter Brune PetscInt i; 376e111a19SKarl Rupp 38eaedb033SPeter Brune PetscFunctionBegin; 39eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 409566063dSJacob Faibussowitsch if (nasm->xl) PetscCall(VecDestroy(&nasm->xl[i])); 419566063dSJacob Faibussowitsch if (nasm->x) PetscCall(VecDestroy(&nasm->x[i])); 429566063dSJacob Faibussowitsch if (nasm->y) PetscCall(VecDestroy(&nasm->y[i])); 439566063dSJacob Faibussowitsch if (nasm->b) PetscCall(VecDestroy(&nasm->b[i])); 44eaedb033SPeter Brune 459566063dSJacob Faibussowitsch if (nasm->subsnes) PetscCall(SNESDestroy(&nasm->subsnes[i])); 469566063dSJacob Faibussowitsch if (nasm->oscatter) PetscCall(VecScatterDestroy(&nasm->oscatter[i])); 479566063dSJacob Faibussowitsch if (nasm->oscatter_copy) PetscCall(VecScatterDestroy(&nasm->oscatter_copy[i])); 489566063dSJacob Faibussowitsch if (nasm->iscatter) PetscCall(VecScatterDestroy(&nasm->iscatter[i])); 499566063dSJacob Faibussowitsch if (nasm->gscatter) PetscCall(VecScatterDestroy(&nasm->gscatter[i])); 50eaedb033SPeter Brune } 51111ade9eSPeter Brune 529566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->x)); 539566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->xl)); 549566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->y)); 559566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->b)); 56111ade9eSPeter Brune 579566063dSJacob Faibussowitsch if (nasm->xinit) PetscCall(VecDestroy(&nasm->xinit)); 58602bec5dSPeter Brune 599566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->subsnes)); 609566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->oscatter)); 619566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->oscatter_copy)); 629566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->iscatter)); 639566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->gscatter)); 64b20c023fSPeter Brune 65f10b3e88SPatrick Farrell if (nasm->weight_set) { 669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&nasm->weight)); 67f10b3e88SPatrick Farrell } 68f10b3e88SPatrick Farrell 69b20c023fSPeter Brune nasm->eventrestrictinterp = 0; 70b20c023fSPeter Brune nasm->eventsubsolve = 0; 71eaedb033SPeter Brune PetscFunctionReturn(0); 72eaedb033SPeter Brune } 73eaedb033SPeter Brune 7440244768SBarry Smith static PetscErrorCode SNESDestroy_NASM(SNES snes) 75eaedb033SPeter Brune { 76eaedb033SPeter Brune PetscFunctionBegin; 779566063dSJacob Faibussowitsch PetscCall(SNESReset_NASM(snes)); 78*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetType_C",NULL)); 79*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetType_C",NULL)); 80*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",NULL)); 81*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",NULL)); 82*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",NULL)); 83*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",NULL)); 84*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",NULL)); 85*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",NULL)); 869566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 87eaedb033SPeter Brune PetscFunctionReturn(0); 88eaedb033SPeter Brune } 89eaedb033SPeter Brune 9040244768SBarry Smith static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx) 910adebc6cSBarry Smith { 92111ade9eSPeter Brune Vec bcs = (Vec)ctx; 936e111a19SKarl Rupp 94111ade9eSPeter Brune PetscFunctionBegin; 959566063dSJacob Faibussowitsch PetscCall(VecCopy(bcs,l)); 96111ade9eSPeter Brune PetscFunctionReturn(0); 97111ade9eSPeter Brune } 98111ade9eSPeter Brune 9940244768SBarry Smith static PetscErrorCode SNESSetUp_NASM(SNES snes) 100eaedb033SPeter Brune { 101eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 10276857b2aSPeter Brune DM dm,subdm; 103111ade9eSPeter Brune DM *subdms; 104111ade9eSPeter Brune PetscInt i; 105eaedb033SPeter Brune const char *optionsprefix; 106111ade9eSPeter Brune Vec F; 107ed3c11a9SPeter Brune PetscMPIInt size; 108ed3c11a9SPeter Brune KSP ksp; 109ed3c11a9SPeter Brune PC pc; 110eaedb033SPeter Brune 111eaedb033SPeter Brune PetscFunctionBegin; 112eaedb033SPeter Brune if (!nasm->subsnes) { 1139566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 1140a696f66SPeter Brune if (dm) { 115eaedb033SPeter Brune nasm->usesdm = PETSC_TRUE; 1169566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms)); 11728b400f6SJacob Faibussowitsch PetscCheck(subdms,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 1189566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter)); 1199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nasm->n, &nasm->oscatter_copy)); 120f10b3e88SPatrick Farrell for (i=0; i<nasm->n; i++) { 1219566063dSJacob Faibussowitsch PetscCall(VecScatterCopy(nasm->oscatter[i], &nasm->oscatter_copy[i])); 122f10b3e88SPatrick Farrell } 123eaedb033SPeter Brune 1249566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 1259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nasm->n,&nasm->subsnes)); 126111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 1279566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i])); 1289566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)nasm->subsnes[i], (PetscObject)snes, 1)); 1299566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix)); 1309566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_")); 1319566063dSJacob Faibussowitsch PetscCall(SNESSetDM(nasm->subsnes[i],subdms[i])); 1329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size)); 133ed3c11a9SPeter Brune if (size == 1) { 1349566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(nasm->subsnes[i],&ksp)); 1359566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 1369566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp,KSPPREONLY)); 1379566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCLU)); 138ed3c11a9SPeter Brune } 1399566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(nasm->subsnes[i])); 1409566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdms[i])); 141111ade9eSPeter Brune } 1429566063dSJacob Faibussowitsch PetscCall(PetscFree(subdms)); 143ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!"); 144ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!"); 145111ade9eSPeter Brune /* allocate the global vectors */ 146e245e0aaSPeter Brune if (!nasm->x) { 1479566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nasm->n,&nasm->x)); 148e245e0aaSPeter Brune } 149e245e0aaSPeter Brune if (!nasm->xl) { 1509566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nasm->n,&nasm->xl)); 151e245e0aaSPeter Brune } 152e245e0aaSPeter Brune if (!nasm->y) { 1539566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nasm->n,&nasm->y)); 154e245e0aaSPeter Brune } 155e245e0aaSPeter Brune if (!nasm->b) { 1569566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nasm->n,&nasm->b)); 157e245e0aaSPeter Brune } 158111ade9eSPeter Brune 159111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 1609566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL)); 1619566063dSJacob Faibussowitsch if (!nasm->x[i]) PetscCall(VecDuplicate(F,&nasm->x[i])); 1629566063dSJacob Faibussowitsch if (!nasm->y[i]) PetscCall(VecDuplicate(F,&nasm->y[i])); 1639566063dSJacob Faibussowitsch if (!nasm->b[i]) PetscCall(VecDuplicate(F,&nasm->b[i])); 16476857b2aSPeter Brune if (!nasm->xl[i]) { 1659566063dSJacob Faibussowitsch PetscCall(SNESGetDM(nasm->subsnes[i],&subdm)); 1669566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(subdm,&nasm->xl[i])); 1679566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i])); 168111ade9eSPeter Brune } 16961ba4676SBarry Smith } 170602bec5dSPeter Brune if (nasm->finaljacobian) { 1719566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 172602bec5dSPeter Brune if (nasm->fjtype == 2) { 1739566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol,&nasm->xinit)); 174602bec5dSPeter Brune } 175602bec5dSPeter Brune for (i=0; i<nasm->n;i++) { 1769566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(nasm->subsnes[i])); 177602bec5dSPeter Brune } 178602bec5dSPeter Brune } 179eaedb033SPeter Brune PetscFunctionReturn(0); 180eaedb033SPeter Brune } 181eaedb033SPeter Brune 18240244768SBarry Smith static PetscErrorCode SNESSetFromOptions_NASM(PetscOptionItems *PetscOptionsObject,SNES snes) 183eaedb033SPeter Brune { 184111ade9eSPeter Brune PCASMType asmtype; 18583dc3634SPierre Jolivet PetscBool flg,monflg; 186111ade9eSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 1876e111a19SKarl Rupp 188eaedb033SPeter Brune PetscFunctionBegin; 189d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Nonlinear Additive Schwarz options"); 1909566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg)); 1919566063dSJacob Faibussowitsch if (flg) PetscCall(SNESNASMSetType(snes,asmtype)); 192b20c023fSPeter Brune flg = PETSC_FALSE; 193b20c023fSPeter Brune monflg = PETSC_TRUE; 1949566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_nasm_damping","The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)","SNESNASMSetDamping",nasm->damping,&nasm->damping,&flg)); 1959566063dSJacob Faibussowitsch if (flg) PetscCall(SNESNASMSetDamping(snes,nasm->damping)); 1969566063dSJacob Faibussowitsch PetscCall(PetscOptionsDeprecated("-snes_nasm_sub_view",NULL,"3.15","Use -snes_view ::ascii_info_detail")); 1979566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL)); 1989566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL)); 1999566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg)); 200b20c023fSPeter Brune if (flg) { 2019566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve)); 2029566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp)); 203b20c023fSPeter Brune } 204d0609cedSBarry Smith PetscOptionsHeadEnd(); 205eaedb033SPeter Brune PetscFunctionReturn(0); 206eaedb033SPeter Brune } 207eaedb033SPeter Brune 20840244768SBarry Smith static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 209eaedb033SPeter Brune { 210b20c023fSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 211a4f17876SPeter Brune PetscMPIInt rank,size; 212dd2fa690SBarry Smith PetscInt i,N,bsz; 213b20c023fSPeter Brune PetscBool iascii,isstring; 214b20c023fSPeter Brune PetscViewer sviewer; 215ce94432eSBarry Smith MPI_Comm comm; 21683dc3634SPierre Jolivet PetscViewerFormat format; 21783dc3634SPierre Jolivet const char *prefix; 218b20c023fSPeter Brune 219b20c023fSPeter Brune PetscFunctionBegin; 2209566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes,&comm)); 2219566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 2229566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring)); 2239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 2249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 2251c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&nasm->n,&N,1,MPIU_INT,MPI_SUM,comm)); 226b20c023fSPeter Brune if (iascii) { 22763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " total subdomain blocks = %" PetscInt_FMT "\n",N)); 2289566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer,&format)); 22983dc3634SPierre Jolivet if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 230a4f17876SPeter Brune if (nasm->subsnes) { 2319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Local solver information for first block on rank 0:\n")); 2329566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes,&prefix)); 2339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Use -%ssnes_view ::ascii_info_detail to display information for all blocks\n",prefix?prefix:"")); 2349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2359566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 236dd400576SPatrick Sanan if (rank == 0) { 2379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2389566063dSJacob Faibussowitsch PetscCall(SNESView(nasm->subsnes[0],sviewer)); 2399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 240a4f17876SPeter Brune } 2419566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 243a4f17876SPeter Brune } 244a4f17876SPeter Brune } else { 245a4f17876SPeter Brune /* print the solver on each block */ 2469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 24763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %" PetscInt_FMT "\n",(int)rank,nasm->n)); 2489566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 2509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Local solver information for each block is in the following SNES objects:\n")); 2519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n")); 2539566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 254a4f17876SPeter Brune for (i=0; i<nasm->n; i++) { 2559566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(nasm->x[i],&bsz)); 25663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(sviewer,"[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n",(int)rank,i,bsz)); 2579566063dSJacob Faibussowitsch PetscCall(SNESView(nasm->subsnes[i],sviewer)); 2589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n")); 259b20c023fSPeter Brune } 2609566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 2619566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 263a4f17876SPeter Brune } 264b20c023fSPeter Brune } else if (isstring) { 26563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer," blocks=%" PetscInt_FMT ",type=%s",N,SNESNASMTypes[nasm->type])); 2669566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 2679566063dSJacob Faibussowitsch if (nasm->subsnes && rank == 0) PetscCall(SNESView(nasm->subsnes[0],sviewer)); 2689566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 269b20c023fSPeter Brune } 270eaedb033SPeter Brune PetscFunctionReturn(0); 271eaedb033SPeter Brune } 272eaedb033SPeter Brune 273e0331734SPeter Brune /*@ 274e0331734SPeter Brune SNESNASMSetType - Set the type of subdomain update used 275e0331734SPeter Brune 276e0331734SPeter Brune Logically Collective on SNES 277e0331734SPeter Brune 278e0331734SPeter Brune Input Parameters: 279e0331734SPeter Brune + SNES - the SNES context 280e0331734SPeter Brune - type - the type of update, PC_ASM_BASIC or PC_ASM_RESTRICT 281e0331734SPeter Brune 282e0331734SPeter Brune Level: intermediate 283e0331734SPeter Brune 284db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetType()`, `PCASMSetType()` 285e0331734SPeter Brune @*/ 286e0331734SPeter Brune PetscErrorCode SNESNASMSetType(SNES snes,PCASMType type) 287e0331734SPeter Brune { 288e0331734SPeter Brune PetscErrorCode (*f)(SNES,PCASMType); 289e0331734SPeter Brune 290e0331734SPeter Brune PetscFunctionBegin; 2919566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetType_C",&f)); 2929566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes,type)); 293e0331734SPeter Brune PetscFunctionReturn(0); 294e0331734SPeter Brune } 295e0331734SPeter Brune 29640244768SBarry Smith static PetscErrorCode SNESNASMSetType_NASM(SNES snes,PCASMType type) 297e0331734SPeter Brune { 298e0331734SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 299e0331734SPeter Brune 300e0331734SPeter Brune PetscFunctionBegin; 3010b121fc5SBarry Smith PetscCheck(type == PC_ASM_BASIC || type == PC_ASM_RESTRICT,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESNASM only supports basic and restrict types"); 302e0331734SPeter Brune nasm->type = type; 303e0331734SPeter Brune PetscFunctionReturn(0); 304e0331734SPeter Brune } 305e0331734SPeter Brune 306e0331734SPeter Brune /*@ 307e0331734SPeter Brune SNESNASMGetType - Get the type of subdomain update used 308e0331734SPeter Brune 309e0331734SPeter Brune Logically Collective on SNES 310e0331734SPeter Brune 311e0331734SPeter Brune Input Parameters: 312e0331734SPeter Brune . SNES - the SNES context 313e0331734SPeter Brune 314e0331734SPeter Brune Output Parameters: 315e0331734SPeter Brune . type - the type of update 316e0331734SPeter Brune 317e0331734SPeter Brune Level: intermediate 318e0331734SPeter Brune 319db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMSetType()`, `PCASMGetType()` 320e0331734SPeter Brune @*/ 321e0331734SPeter Brune PetscErrorCode SNESNASMGetType(SNES snes,PCASMType *type) 322e0331734SPeter Brune { 323e0331734SPeter Brune PetscFunctionBegin; 324cac4c232SBarry Smith PetscUseMethod(snes,"SNESNASMGetType_C",(SNES,PCASMType*),(snes,type)); 325e0331734SPeter Brune PetscFunctionReturn(0); 326e0331734SPeter Brune } 327e0331734SPeter Brune 32840244768SBarry Smith static PetscErrorCode SNESNASMGetType_NASM(SNES snes,PCASMType *type) 329e0331734SPeter Brune { 330e0331734SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 331e0331734SPeter Brune 332e0331734SPeter Brune PetscFunctionBegin; 333e0331734SPeter Brune *type = nasm->type; 334e0331734SPeter Brune PetscFunctionReturn(0); 335e0331734SPeter Brune } 336e0331734SPeter Brune 33776857b2aSPeter Brune /*@ 33876857b2aSPeter Brune SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems. 33976857b2aSPeter Brune 34076857b2aSPeter Brune Not Collective 34176857b2aSPeter Brune 34276857b2aSPeter Brune Input Parameters: 34376857b2aSPeter Brune + SNES - the SNES context 34476857b2aSPeter Brune . n - the number of local subdomains 34576857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 34676857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 34776857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 34876857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 34976857b2aSPeter Brune 35076857b2aSPeter Brune Level: intermediate 35176857b2aSPeter Brune 352db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetSubdomains()` 35376857b2aSPeter Brune @*/ 354a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 355a6dfd86eSKarl Rupp { 356111ade9eSPeter Brune PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 3576e111a19SKarl Rupp 358eaedb033SPeter Brune PetscFunctionBegin; 3599566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f)); 3609566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes,n,subsnes,iscatter,oscatter,gscatter)); 361eaedb033SPeter Brune PetscFunctionReturn(0); 362eaedb033SPeter Brune } 363eaedb033SPeter Brune 36440244768SBarry Smith static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 365a6dfd86eSKarl Rupp { 366eaedb033SPeter Brune PetscInt i; 367eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 3686e111a19SKarl Rupp 369eaedb033SPeter Brune PetscFunctionBegin; 37028b400f6SJacob Faibussowitsch PetscCheck(!snes->setupcalled,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 371eaedb033SPeter Brune 372111ade9eSPeter Brune /* tear down the previously set things */ 3739566063dSJacob Faibussowitsch PetscCall(SNESReset(snes)); 374111ade9eSPeter Brune 375eaedb033SPeter Brune nasm->n = n; 376111ade9eSPeter Brune if (oscatter) { 3779566063dSJacob Faibussowitsch for (i=0; i<n; i++) PetscCall(PetscObjectReference((PetscObject)oscatter[i])); 378eaedb033SPeter Brune } 379111ade9eSPeter Brune if (iscatter) { 3809566063dSJacob Faibussowitsch for (i=0; i<n; i++) PetscCall(PetscObjectReference((PetscObject)iscatter[i])); 381eaedb033SPeter Brune } 382111ade9eSPeter Brune if (gscatter) { 3839566063dSJacob Faibussowitsch for (i=0; i<n; i++) PetscCall(PetscObjectReference((PetscObject)gscatter[i])); 384111ade9eSPeter Brune } 385111ade9eSPeter Brune if (oscatter) { 3869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&nasm->oscatter)); 3879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&nasm->oscatter_copy)); 388eaedb033SPeter Brune for (i=0; i<n; i++) { 389111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 3909566063dSJacob Faibussowitsch PetscCall(VecScatterCopy(oscatter[i], &nasm->oscatter_copy[i])); 391eaedb033SPeter Brune } 392111ade9eSPeter Brune } 393111ade9eSPeter Brune if (iscatter) { 3949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&nasm->iscatter)); 395eaedb033SPeter Brune for (i=0; i<n; i++) { 396111ade9eSPeter Brune nasm->iscatter[i] = iscatter[i]; 397eaedb033SPeter Brune } 398eaedb033SPeter Brune } 399111ade9eSPeter Brune if (gscatter) { 4009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&nasm->gscatter)); 401eaedb033SPeter Brune for (i=0; i<n; i++) { 402111ade9eSPeter Brune nasm->gscatter[i] = gscatter[i]; 403eaedb033SPeter Brune } 404eaedb033SPeter Brune } 405111ade9eSPeter Brune 406eaedb033SPeter Brune if (subsnes) { 4079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&nasm->subsnes)); 408eaedb033SPeter Brune for (i=0; i<n; i++) { 409eaedb033SPeter Brune nasm->subsnes[i] = subsnes[i]; 410eaedb033SPeter Brune } 411eaedb033SPeter Brune } 412eaedb033SPeter Brune PetscFunctionReturn(0); 413eaedb033SPeter Brune } 414eaedb033SPeter Brune 41576857b2aSPeter Brune /*@ 41676857b2aSPeter Brune SNESNASMGetSubdomains - Get the local subdomain context. 41776857b2aSPeter Brune 41876857b2aSPeter Brune Not Collective 41976857b2aSPeter Brune 420f899ff85SJose E. Roman Input Parameter: 42176857b2aSPeter Brune . SNES - the SNES context 42276857b2aSPeter Brune 42376857b2aSPeter Brune Output Parameters: 42476857b2aSPeter Brune + n - the number of local subdomains 42576857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 42676857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 42776857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 42876857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 42976857b2aSPeter Brune 43076857b2aSPeter Brune Level: intermediate 43176857b2aSPeter Brune 432db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMSetSubdomains()` 43376857b2aSPeter Brune @*/ 43476857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 43576857b2aSPeter Brune { 43676857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**); 43776857b2aSPeter Brune 43876857b2aSPeter Brune PetscFunctionBegin; 4399566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f)); 4409566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes,n,subsnes,iscatter,oscatter,gscatter)); 44176857b2aSPeter Brune PetscFunctionReturn(0); 44276857b2aSPeter Brune } 44376857b2aSPeter Brune 44440244768SBarry Smith static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 44576857b2aSPeter Brune { 44676857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 44776857b2aSPeter Brune 44876857b2aSPeter Brune PetscFunctionBegin; 44976857b2aSPeter Brune if (n) *n = nasm->n; 45076857b2aSPeter Brune if (oscatter) *oscatter = nasm->oscatter; 45176857b2aSPeter Brune if (iscatter) *iscatter = nasm->iscatter; 45276857b2aSPeter Brune if (gscatter) *gscatter = nasm->gscatter; 45383dc3634SPierre Jolivet if (subsnes) *subsnes = nasm->subsnes; 45476857b2aSPeter Brune PetscFunctionReturn(0); 45576857b2aSPeter Brune } 45676857b2aSPeter Brune 45776857b2aSPeter Brune /*@ 45876857b2aSPeter Brune SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors 45976857b2aSPeter Brune 46076857b2aSPeter Brune Not Collective 46176857b2aSPeter Brune 462f899ff85SJose E. Roman Input Parameter: 46376857b2aSPeter Brune . SNES - the SNES context 46476857b2aSPeter Brune 46576857b2aSPeter Brune Output Parameters: 46676857b2aSPeter Brune + n - the number of local subdomains 46776857b2aSPeter Brune . x - The subdomain solution vector 46876857b2aSPeter Brune . y - The subdomain step vector 46976857b2aSPeter Brune . b - The subdomain RHS vector 47076857b2aSPeter Brune - xl - The subdomain local vectors (ghosted) 47176857b2aSPeter Brune 47276857b2aSPeter Brune Level: developer 47376857b2aSPeter Brune 474db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetSubdomains()` 47576857b2aSPeter Brune @*/ 47676857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl) 47776857b2aSPeter Brune { 47876857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**); 47976857b2aSPeter Brune 48076857b2aSPeter Brune PetscFunctionBegin; 4819566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f)); 4829566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes,n,x,y,b,xl)); 48376857b2aSPeter Brune PetscFunctionReturn(0); 48476857b2aSPeter Brune } 48576857b2aSPeter Brune 48640244768SBarry Smith static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl) 48776857b2aSPeter Brune { 48876857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 48976857b2aSPeter Brune 49076857b2aSPeter Brune PetscFunctionBegin; 49176857b2aSPeter Brune if (n) *n = nasm->n; 49276857b2aSPeter Brune if (x) *x = nasm->x; 49376857b2aSPeter Brune if (y) *y = nasm->y; 49476857b2aSPeter Brune if (b) *b = nasm->b; 49576857b2aSPeter Brune if (xl) *xl = nasm->xl; 49676857b2aSPeter Brune PetscFunctionReturn(0); 49776857b2aSPeter Brune } 49876857b2aSPeter Brune 499d728fb7dSPeter Brune /*@ 5008bf45196SRichard Tran Mills SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain Jacobians upon convergence 501d728fb7dSPeter Brune 502d728fb7dSPeter Brune Collective on SNES 503d728fb7dSPeter Brune 504d728fb7dSPeter Brune Input Parameters: 505d728fb7dSPeter Brune + SNES - the SNES context 5068bf45196SRichard Tran Mills - flg - indication of whether to compute the Jacobians or not 507d728fb7dSPeter Brune 508d728fb7dSPeter Brune Level: developer 509d728fb7dSPeter Brune 51095452b02SPatrick Sanan Notes: 5118bf45196SRichard Tran Mills This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global Jacobian 512d728fb7dSPeter Brune is needed at each linear iteration. 513d728fb7dSPeter Brune 514db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetSubdomains()` 515d728fb7dSPeter Brune @*/ 516d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg) 517d728fb7dSPeter Brune { 518d728fb7dSPeter Brune PetscErrorCode (*f)(SNES,PetscBool); 519d728fb7dSPeter Brune 520d728fb7dSPeter Brune PetscFunctionBegin; 5219566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f)); 5229566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes,flg)); 523d728fb7dSPeter Brune PetscFunctionReturn(0); 524d728fb7dSPeter Brune } 525d728fb7dSPeter Brune 52640244768SBarry Smith static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg) 527d728fb7dSPeter Brune { 528d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 529d728fb7dSPeter Brune 530d728fb7dSPeter Brune PetscFunctionBegin; 531d728fb7dSPeter Brune nasm->finaljacobian = flg; 532d728fb7dSPeter Brune PetscFunctionReturn(0); 533d728fb7dSPeter Brune } 53476857b2aSPeter Brune 535610116beSPeter Brune /*@ 536610116beSPeter Brune SNESNASMSetDamping - Sets the update damping for NASM 537610116beSPeter Brune 538610116beSPeter Brune Logically collective on SNES 539610116beSPeter Brune 540610116beSPeter Brune Input Parameters: 541610116beSPeter Brune + SNES - the SNES context 542610116beSPeter Brune - dmp - damping 543610116beSPeter Brune 544610116beSPeter Brune Level: intermediate 545610116beSPeter Brune 54695452b02SPatrick Sanan Notes: 54795452b02SPatrick Sanan The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 5485dfa0f3bSBarry Smith 549db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetDamping()` 550610116beSPeter Brune @*/ 551610116beSPeter Brune PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp) 552610116beSPeter Brune { 553610116beSPeter Brune PetscErrorCode (*f)(SNES,PetscReal); 554610116beSPeter Brune 555610116beSPeter Brune PetscFunctionBegin; 5569566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f)); 5579566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes,dmp)); 558610116beSPeter Brune PetscFunctionReturn(0); 559610116beSPeter Brune } 560610116beSPeter Brune 56140244768SBarry Smith static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp) 562610116beSPeter Brune { 563610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 564610116beSPeter Brune 565610116beSPeter Brune PetscFunctionBegin; 566610116beSPeter Brune nasm->damping = dmp; 567610116beSPeter Brune PetscFunctionReturn(0); 568610116beSPeter Brune } 569610116beSPeter Brune 570610116beSPeter Brune /*@ 571610116beSPeter Brune SNESNASMGetDamping - Gets the update damping for NASM 572610116beSPeter Brune 573610116beSPeter Brune Not Collective 574610116beSPeter Brune 575610116beSPeter Brune Input Parameters: 576610116beSPeter Brune + SNES - the SNES context 577610116beSPeter Brune - dmp - damping 578610116beSPeter Brune 579610116beSPeter Brune Level: intermediate 580610116beSPeter Brune 581db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMSetDamping()` 582610116beSPeter Brune @*/ 583610116beSPeter Brune PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp) 584610116beSPeter Brune { 585610116beSPeter Brune PetscFunctionBegin; 586cac4c232SBarry Smith PetscUseMethod(snes,"SNESNASMGetDamping_C",(SNES,PetscReal*),(snes,dmp)); 587610116beSPeter Brune PetscFunctionReturn(0); 588610116beSPeter Brune } 589610116beSPeter Brune 59040244768SBarry Smith static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp) 591610116beSPeter Brune { 592610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 593610116beSPeter Brune 594610116beSPeter Brune PetscFunctionBegin; 595610116beSPeter Brune *dmp = nasm->damping; 596610116beSPeter Brune PetscFunctionReturn(0); 597610116beSPeter Brune } 598610116beSPeter Brune 59914eb1c5cSMatthew G. Knepley /* 60014eb1c5cSMatthew G. Knepley Input Parameters: 60114eb1c5cSMatthew G. Knepley + snes - The solver 60214eb1c5cSMatthew G. Knepley . B - The RHS vector 60314eb1c5cSMatthew G. Knepley - X - The initial guess 60414eb1c5cSMatthew G. Knepley 60514eb1c5cSMatthew G. Knepley Output Parameters: 60614eb1c5cSMatthew G. Knepley . Y - The solution update 60714eb1c5cSMatthew G. Knepley 60814eb1c5cSMatthew G. Knepley TODO: All scatters should be packed into one 60914eb1c5cSMatthew G. Knepley */ 6100adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) 6110adebc6cSBarry Smith { 612eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 613258e1594SPeter Brune SNES subsnes; 614eaedb033SPeter Brune PetscInt i; 615610116beSPeter Brune PetscReal dmp; 616f10b3e88SPatrick Farrell Vec Xl,Bl,Yl,Xlloc; 617f10b3e88SPatrick Farrell VecScatter iscat,oscat,gscat,oscat_copy; 618111ade9eSPeter Brune DM dm,subdm; 619e0331734SPeter Brune PCASMType type; 6200adebc6cSBarry Smith 621eaedb033SPeter Brune PetscFunctionBegin; 6229566063dSJacob Faibussowitsch PetscCall(SNESNASMGetType(snes,&type)); 6239566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 6249566063dSJacob Faibussowitsch PetscCall(VecSet(Y,0)); 6259566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0)); 626eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 627f10b3e88SPatrick Farrell /* scatter the solution to the global solution and the local solution */ 628f10b3e88SPatrick Farrell Xl = nasm->x[i]; 62970c78f05SPeter Brune Xlloc = nasm->xl[i]; 63070c78f05SPeter Brune oscat = nasm->oscatter[i]; 631f10b3e88SPatrick Farrell oscat_copy = nasm->oscatter_copy[i]; 632f10b3e88SPatrick Farrell gscat = nasm->gscatter[i]; 6339566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD)); 6349566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD)); 63570c78f05SPeter Brune if (B) { 63670c78f05SPeter Brune /* scatter the RHS to the local RHS */ 63770c78f05SPeter Brune Bl = nasm->b[i]; 6389566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscat_copy,B,Bl,INSERT_VALUES,SCATTER_FORWARD)); 63970c78f05SPeter Brune } 64070c78f05SPeter Brune } 6419566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0)); 642b20c023fSPeter Brune 6439566063dSJacob Faibussowitsch if (nasm->eventsubsolve) PetscCall(PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0)); 64470c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 64570c78f05SPeter Brune Xl = nasm->x[i]; 64670c78f05SPeter Brune Xlloc = nasm->xl[i]; 64770c78f05SPeter Brune Yl = nasm->y[i]; 648258e1594SPeter Brune subsnes = nasm->subsnes[i]; 6499566063dSJacob Faibussowitsch PetscCall(SNESGetDM(subsnes,&subdm)); 650111ade9eSPeter Brune iscat = nasm->iscatter[i]; 651111ade9eSPeter Brune oscat = nasm->oscatter[i]; 652f10b3e88SPatrick Farrell oscat_copy = nasm->oscatter_copy[i]; 653111ade9eSPeter Brune gscat = nasm->gscatter[i]; 6549566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD)); 6559566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD)); 65624b7f281SPeter Brune if (B) { 65724b7f281SPeter Brune Bl = nasm->b[i]; 6589566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscat_copy,B,Bl,INSERT_VALUES,SCATTER_FORWARD)); 659ed3c11a9SPeter Brune } else Bl = NULL; 660f10b3e88SPatrick Farrell 6619566063dSJacob Faibussowitsch PetscCall(DMSubDomainRestrict(dm,oscat,gscat,subdm)); 6629566063dSJacob Faibussowitsch PetscCall(VecCopy(Xl,Yl)); 6639566063dSJacob Faibussowitsch PetscCall(SNESSolve(subsnes,Bl,Xl)); 6649566063dSJacob Faibussowitsch PetscCall(VecAYPX(Yl,-1.0,Xl)); 6659566063dSJacob Faibussowitsch PetscCall(VecScale(Yl, nasm->damping)); 666e0331734SPeter Brune if (type == PC_ASM_BASIC) { 6679566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE)); 6689566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE)); 669e0331734SPeter Brune } else if (type == PC_ASM_RESTRICT) { 6709566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE)); 6719566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE)); 672ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 673eaedb033SPeter Brune } 6749566063dSJacob Faibussowitsch if (nasm->eventsubsolve) PetscCall(PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0)); 6759566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0)); 676f10b3e88SPatrick Farrell if (nasm->weight_set) { 6779566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Y,Y,nasm->weight)); 678f10b3e88SPatrick Farrell } 6799566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0)); 6809566063dSJacob Faibussowitsch PetscCall(SNESNASMGetDamping(snes,&dmp)); 6819566063dSJacob Faibussowitsch PetscCall(VecAXPY(X,dmp,Y)); 682eaedb033SPeter Brune PetscFunctionReturn(0); 683eaedb033SPeter Brune } 684eaedb033SPeter Brune 68540244768SBarry Smith static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 686d728fb7dSPeter Brune { 687602bec5dSPeter Brune Vec X = Xfinal; 688d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 689d728fb7dSPeter Brune SNES subsnes; 690ca44f815SPeter Brune PetscInt i,lag = 1; 691e59f0a30SPeter Brune Vec Xlloc,Xl,Fl,F; 692d728fb7dSPeter Brune VecScatter oscat,gscat; 693d728fb7dSPeter Brune DM dm,subdm; 694d1e9a80fSBarry Smith 695d728fb7dSPeter Brune PetscFunctionBegin; 696602bec5dSPeter Brune if (nasm->fjtype == 2) X = nasm->xinit; 697e59f0a30SPeter Brune F = snes->vec_func; 6989566063dSJacob Faibussowitsch if (snes->normschedule == SNES_NORM_NONE) PetscCall(SNESComputeFunction(snes,X,F)); 6999566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre)); 7009566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 7019566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0)); 702602bec5dSPeter Brune if (nasm->fjtype != 1) { 703d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 704d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 705d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 7069566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD)); 707602bec5dSPeter Brune } 708d728fb7dSPeter Brune } 7099566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0)); 710d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 711e59f0a30SPeter Brune Fl = nasm->subsnes[i]->vec_func; 712d728fb7dSPeter Brune Xl = nasm->x[i]; 713d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 714d728fb7dSPeter Brune subsnes = nasm->subsnes[i]; 715d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 716d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 7179566063dSJacob Faibussowitsch if (nasm->fjtype != 1) PetscCall(VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD)); 7189566063dSJacob Faibussowitsch PetscCall(SNESGetDM(subsnes,&subdm)); 7199566063dSJacob Faibussowitsch PetscCall(DMSubDomainRestrict(dm,oscat,gscat,subdm)); 720602bec5dSPeter Brune if (nasm->fjtype != 1) { 7219566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl)); 7229566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl)); 723602bec5dSPeter Brune } 724ca44f815SPeter Brune if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 725ca44f815SPeter Brune else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 7269566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(subsnes,Xl,Fl)); 7279566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(subsnes,Xl,subsnes->jacobian,subsnes->jacobian_pre)); 728ca44f815SPeter Brune if (lag > 1) subsnes->lagjacobian = lag; 729d728fb7dSPeter Brune } 730d728fb7dSPeter Brune PetscFunctionReturn(0); 731d728fb7dSPeter Brune } 732d728fb7dSPeter Brune 73340244768SBarry Smith static PetscErrorCode SNESSolve_NASM(SNES snes) 734eaedb033SPeter Brune { 735eaedb033SPeter Brune Vec F; 736eaedb033SPeter Brune Vec X; 737eaedb033SPeter Brune Vec B; 738111ade9eSPeter Brune Vec Y; 739eaedb033SPeter Brune PetscInt i; 740ed3c11a9SPeter Brune PetscReal fnorm = 0.0; 741365a6726SPeter Brune SNESNormSchedule normschedule; 742d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 743eaedb033SPeter Brune 744eaedb033SPeter Brune PetscFunctionBegin; 745c579b300SPatrick Farrell 7460b121fc5SBarry Smith PetscCheck(!snes->xl & !snes->xu && !snes->ops->computevariablebounds,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 747c579b300SPatrick Farrell 7489566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite)); 749eaedb033SPeter Brune X = snes->vec_sol; 750111ade9eSPeter Brune Y = snes->vec_sol_update; 751eaedb033SPeter Brune F = snes->vec_func; 752eaedb033SPeter Brune B = snes->vec_rhs; 753eaedb033SPeter Brune 7549566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 755eaedb033SPeter Brune snes->iter = 0; 756eaedb033SPeter Brune snes->norm = 0.; 7579566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 758eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 7599566063dSJacob Faibussowitsch PetscCall(SNESGetNormSchedule(snes, &normschedule)); 760365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 761eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 762eaedb033SPeter Brune if (!snes->vec_func_init_set) { 7639566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,X,F)); 7641aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 765eaedb033SPeter Brune 7669566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 767422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 7689566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 769eaedb033SPeter Brune snes->iter = 0; 770eaedb033SPeter Brune snes->norm = fnorm; 7719566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7729566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0)); 7739566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,0,snes->norm)); 774eaedb033SPeter Brune 775eaedb033SPeter Brune /* test convergence */ 7769566063dSJacob Faibussowitsch PetscCall((*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP)); 777eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 778eaedb033SPeter Brune } else { 7799566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7809566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0)); 7819566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,0,snes->norm)); 782eaedb033SPeter Brune } 783eaedb033SPeter Brune 784eaedb033SPeter Brune /* Call general purpose update function */ 785eaedb033SPeter Brune if (snes->ops->update) { 7869566063dSJacob Faibussowitsch PetscCall((*snes->ops->update)(snes, snes->iter)); 787eaedb033SPeter Brune } 788602bec5dSPeter Brune /* copy the initial solution over for later */ 7899566063dSJacob Faibussowitsch if (nasm->fjtype == 2) PetscCall(VecCopy(X,nasm->xinit)); 790eaedb033SPeter Brune 791eaedb033SPeter Brune for (i=0; i < snes->max_its; i++) { 7929566063dSJacob Faibussowitsch PetscCall(SNESNASMSolveLocal_Private(snes,B,Y,X)); 793365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 7949566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,X,F)); 7959566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 796422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 797eaedb033SPeter Brune } 798eaedb033SPeter Brune /* Monitor convergence */ 7999566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 800eaedb033SPeter Brune snes->iter = i+1; 801eaedb033SPeter Brune snes->norm = fnorm; 8029566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 8039566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0)); 8049566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,snes->iter,snes->norm)); 805eaedb033SPeter Brune /* Test for convergence */ 8069566063dSJacob Faibussowitsch if (normschedule == SNES_NORM_ALWAYS) PetscCall((*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP)); 807d728fb7dSPeter Brune if (snes->reason) break; 808eaedb033SPeter Brune /* Call general purpose update function */ 8099566063dSJacob Faibussowitsch if (snes->ops->update) PetscCall((*snes->ops->update)(snes, snes->iter)); 810eaedb033SPeter Brune } 81107b62357SFande Kong if (nasm->finaljacobian) { 8129566063dSJacob Faibussowitsch PetscCall(SNESNASMComputeFinalJacobian_Private(snes,X)); 81307b62357SFande Kong SNESCheckJacobianDomainerror(snes); 81407b62357SFande Kong } 815365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) { 816eaedb033SPeter Brune if (i == snes->max_its) { 81763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes,"Maximum number of iterations has been reached: %" PetscInt_FMT "\n",snes->max_its)); 818eaedb033SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 819eaedb033SPeter Brune } 8201aa26658SKarl Rupp } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 821eaedb033SPeter Brune PetscFunctionReturn(0); 822eaedb033SPeter Brune } 823eaedb033SPeter Brune 824eaedb033SPeter Brune /*MC 825f3f228e0SPierre Jolivet SNESNASM - Nonlinear Additive Schwarz 826eaedb033SPeter Brune 82770c78603SPeter Brune Options Database: 82870c78603SPeter Brune + -snes_nasm_log - enable logging events for the communication and solve stages 82970c78603SPeter Brune . -snes_nasm_type <basic,restrict> - type of subdomain update used 8305dfa0f3bSBarry Smith . -snes_asm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 83170c78603SPeter Brune . -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate 832150d49b7SHong Zhang . -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the jacobian is calculated at 83370c78603SPeter Brune . -sub_snes_ - options prefix of the subdomain nonlinear solves 83470c78603SPeter Brune . -sub_ksp_ - options prefix of the subdomain Krylov solver 83570c78603SPeter Brune - -sub_pc_ - options prefix of the subdomain preconditioner 83670c78603SPeter Brune 837eaedb033SPeter Brune Level: advanced 838eaedb033SPeter Brune 839951fe5abSBarry Smith Developer Note: This is a non-Newton based nonlinear solver that does not directly require a Jacobian; hence the flag snes->usesksp is set to 8407f851d37SRichard Tran Mills false and SNESView() and -snes_view do not display a KSP object. However, if the flag nasm->finaljacobian is set (for example, if 841951fe5abSBarry Smith NASM is used as a nonlinear preconditioner for KSPASPIN) then SNESSetUpMatrices() is called to generate the Jacobian (needed by KSPASPIN) 842951fe5abSBarry Smith and this utilizes the KSP for storing the matrices, but the KSP is never used for solving a linear system. Note that when SNESNASM is 843951fe5abSBarry Smith used by SNESASPIN they share the same Jacobian matrices because SNESSetUp() (called on the outer SNES KSPASPIN) causes the inner SNES 844951fe5abSBarry Smith object (in this case SNESNASM) to inherit the outer Jacobian matrices. 845951fe5abSBarry Smith 8464f02bc6aSBarry Smith References: 847606c0280SSatish Balay . * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 8484f02bc6aSBarry Smith SIAM Review, 57(4), 2015 8494f02bc6aSBarry Smith 850db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESNASMSetType()`, `SNESNASMGetType()`, `SNESNASMSetSubdomains()`, `SNESNASMGetSubdomains()`, `SNESNASMGetSubdomainVecs()`, `SNESNASMSetComputeFinalJacobian()`, `SNESNASMSetDamping()`, `SNESNASMGetDamping()` 851eaedb033SPeter Brune M*/ 852eaedb033SPeter Brune 8538cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) 854eaedb033SPeter Brune { 855eaedb033SPeter Brune SNES_NASM *nasm; 856eaedb033SPeter Brune 857eaedb033SPeter Brune PetscFunctionBegin; 8589566063dSJacob Faibussowitsch PetscCall(PetscNewLog(snes,&nasm)); 859eaedb033SPeter Brune snes->data = (void*)nasm; 860eaedb033SPeter Brune 861eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 8629e5d0892SLisandro Dalcin nasm->subsnes = NULL; 8639e5d0892SLisandro Dalcin nasm->x = NULL; 8649e5d0892SLisandro Dalcin nasm->xl = NULL; 8659e5d0892SLisandro Dalcin nasm->y = NULL; 8669e5d0892SLisandro Dalcin nasm->b = NULL; 8679e5d0892SLisandro Dalcin nasm->oscatter = NULL; 8689e5d0892SLisandro Dalcin nasm->oscatter_copy = NULL; 8699e5d0892SLisandro Dalcin nasm->iscatter = NULL; 8709e5d0892SLisandro Dalcin nasm->gscatter = NULL; 871610116beSPeter Brune nasm->damping = 1.; 872111ade9eSPeter Brune 873111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 874d728fb7dSPeter Brune nasm->finaljacobian = PETSC_FALSE; 875f10b3e88SPatrick Farrell nasm->weight_set = PETSC_FALSE; 876eaedb033SPeter Brune 877eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 878eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 879eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 880eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 881eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 882eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 883eaedb033SPeter Brune 884eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 885efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 886eaedb033SPeter Brune 8874fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 8884fc747eaSLawrence Mitchell 889602bec5dSPeter Brune nasm->fjtype = 0; 890602bec5dSPeter Brune nasm->xinit = NULL; 8910298fd71SBarry Smith nasm->eventrestrictinterp = 0; 8920298fd71SBarry Smith nasm->eventsubsolve = 0; 893b20c023fSPeter Brune 894eaedb033SPeter Brune if (!snes->tolerancesset) { 895eaedb033SPeter Brune snes->max_its = 10000; 896eaedb033SPeter Brune snes->max_funcs = 10000; 897eaedb033SPeter Brune } 898eaedb033SPeter Brune 8999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetType_C",SNESNASMSetType_NASM)); 9009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetType_C",SNESNASMGetType_NASM)); 9019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM)); 9029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM)); 9039566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",SNESNASMSetDamping_NASM)); 9049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",SNESNASMGetDamping_NASM)); 9059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM)); 9069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM)); 907eaedb033SPeter Brune PetscFunctionReturn(0); 908eaedb033SPeter Brune } 90999e0435eSBarry Smith 910448b6425SPatrick Farrell /*@ 911448b6425SPatrick Farrell SNESNASMGetSNES - Gets a subsolver 912448b6425SPatrick Farrell 913448b6425SPatrick Farrell Not collective 914448b6425SPatrick Farrell 915448b6425SPatrick Farrell Input Parameters: 916448b6425SPatrick Farrell + snes - the SNES context 917448b6425SPatrick Farrell - i - the number of the subsnes to get 918448b6425SPatrick Farrell 919448b6425SPatrick Farrell Output Parameters: 920448b6425SPatrick Farrell . subsnes - the subsolver context 921448b6425SPatrick Farrell 922448b6425SPatrick Farrell Level: intermediate 923448b6425SPatrick Farrell 924db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetNumber()` 925448b6425SPatrick Farrell @*/ 926448b6425SPatrick Farrell PetscErrorCode SNESNASMGetSNES(SNES snes,PetscInt i,SNES *subsnes) 927448b6425SPatrick Farrell { 928448b6425SPatrick Farrell SNES_NASM *nasm = (SNES_NASM*)snes->data; 929448b6425SPatrick Farrell 930448b6425SPatrick Farrell PetscFunctionBegin; 9310b121fc5SBarry Smith PetscCheck(i >= 0 && i < nasm->n,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"No such subsolver"); 932448b6425SPatrick Farrell *subsnes = nasm->subsnes[i]; 933448b6425SPatrick Farrell PetscFunctionReturn(0); 934448b6425SPatrick Farrell } 935448b6425SPatrick Farrell 936448b6425SPatrick Farrell /*@ 937448b6425SPatrick Farrell SNESNASMGetNumber - Gets number of subsolvers 938448b6425SPatrick Farrell 939448b6425SPatrick Farrell Not collective 940448b6425SPatrick Farrell 941448b6425SPatrick Farrell Input Parameters: 942448b6425SPatrick Farrell . snes - the SNES context 943448b6425SPatrick Farrell 944448b6425SPatrick Farrell Output Parameters: 945448b6425SPatrick Farrell . n - the number of subsolvers 946448b6425SPatrick Farrell 947448b6425SPatrick Farrell Level: intermediate 948448b6425SPatrick Farrell 949db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetSNES()` 950448b6425SPatrick Farrell @*/ 951448b6425SPatrick Farrell PetscErrorCode SNESNASMGetNumber(SNES snes,PetscInt *n) 952448b6425SPatrick Farrell { 953448b6425SPatrick Farrell SNES_NASM *nasm = (SNES_NASM*)snes->data; 954448b6425SPatrick Farrell 955448b6425SPatrick Farrell PetscFunctionBegin; 956448b6425SPatrick Farrell *n = nasm->n; 957448b6425SPatrick Farrell PetscFunctionReturn(0); 958448b6425SPatrick Farrell } 959448b6425SPatrick Farrell 960f10b3e88SPatrick Farrell /*@ 961f10b3e88SPatrick Farrell SNESNASMSetWeight - Sets weight to use when adding overlapping updates 962f10b3e88SPatrick Farrell 963f10b3e88SPatrick Farrell Collective 964f10b3e88SPatrick Farrell 965f10b3e88SPatrick Farrell Input Parameters: 966f10b3e88SPatrick Farrell + snes - the SNES context 967f10b3e88SPatrick Farrell - weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in) 968f10b3e88SPatrick Farrell 969f10b3e88SPatrick Farrell Level: intermediate 970f10b3e88SPatrick Farrell 971db781477SPatrick Sanan .seealso: `SNESNASM` 972f10b3e88SPatrick Farrell @*/ 973f10b3e88SPatrick Farrell PetscErrorCode SNESNASMSetWeight(SNES snes,Vec weight) 974f10b3e88SPatrick Farrell { 975f10b3e88SPatrick Farrell SNES_NASM *nasm = (SNES_NASM*)snes->data; 976f10b3e88SPatrick Farrell 977f10b3e88SPatrick Farrell PetscFunctionBegin; 978f10b3e88SPatrick Farrell 9799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&nasm->weight)); 980f10b3e88SPatrick Farrell nasm->weight_set = PETSC_TRUE; 981f10b3e88SPatrick Farrell nasm->weight = weight; 9829566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)nasm->weight)); 983f10b3e88SPatrick Farrell 984f10b3e88SPatrick Farrell PetscFunctionReturn(0); 985f10b3e88SPatrick Farrell } 986