xref: /petsc/src/snes/impls/nasm/nasm.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 
SNESReset_NASM(SNES snes)33d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_NASM(SNES snes)
34d71ae5a4SJacob Faibussowitsch {
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 
6548a46eb9SPierre Jolivet   if (nasm->weight_set) PetscCall(VecDestroy(&nasm->weight));
66f10b3e88SPatrick Farrell 
67b20c023fSPeter Brune   nasm->eventrestrictinterp = 0;
68b20c023fSPeter Brune   nasm->eventsubsolve       = 0;
693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70eaedb033SPeter Brune }
71eaedb033SPeter Brune 
SNESDestroy_NASM(SNES snes)72d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NASM(SNES snes)
73d71ae5a4SJacob Faibussowitsch {
74eaedb033SPeter Brune   PetscFunctionBegin;
759566063dSJacob Faibussowitsch   PetscCall(SNESReset_NASM(snes));
762e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", NULL));
772e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", NULL));
782e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", NULL));
792e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", NULL));
802e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", NULL));
812e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", NULL));
822e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", NULL));
832e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", NULL));
849566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86eaedb033SPeter Brune }
87eaedb033SPeter Brune 
DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,PetscCtx ctx)88*2a8381b2SBarry Smith static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm, Vec g, InsertMode mode, Vec l, PetscCtx ctx)
89d71ae5a4SJacob Faibussowitsch {
90111ade9eSPeter Brune   Vec bcs = (Vec)ctx;
916e111a19SKarl Rupp 
92111ade9eSPeter Brune   PetscFunctionBegin;
939566063dSJacob Faibussowitsch   PetscCall(VecCopy(bcs, l));
943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
95111ade9eSPeter Brune }
96111ade9eSPeter Brune 
SNESSetUp_NASM(SNES snes)97d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NASM(SNES snes)
98d71ae5a4SJacob Faibussowitsch {
99eaedb033SPeter Brune   SNES_NASM  *nasm = (SNES_NASM *)snes->data;
10076857b2aSPeter Brune   DM          dm, subdm;
101111ade9eSPeter Brune   DM         *subdms;
102111ade9eSPeter Brune   PetscInt    i;
103eaedb033SPeter Brune   const char *optionsprefix;
104111ade9eSPeter Brune   Vec         F;
105eaedb033SPeter Brune 
106eaedb033SPeter Brune   PetscFunctionBegin;
107eaedb033SPeter Brune   if (!nasm->subsnes) {
1089566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(snes, &dm));
1090a696f66SPeter Brune     if (dm) {
110eaedb033SPeter Brune       nasm->usesdm = PETSC_TRUE;
1119566063dSJacob Faibussowitsch       PetscCall(DMCreateDomainDecomposition(dm, &nasm->n, NULL, NULL, NULL, &subdms));
11228b400f6SJacob Faibussowitsch       PetscCheck(subdms, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
1139566063dSJacob Faibussowitsch       PetscCall(DMCreateDomainDecompositionScatters(dm, nasm->n, subdms, &nasm->iscatter, &nasm->oscatter, &nasm->gscatter));
1149566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nasm->n, &nasm->oscatter_copy));
11548a46eb9SPierre Jolivet       for (i = 0; i < nasm->n; i++) PetscCall(VecScatterCopy(nasm->oscatter[i], &nasm->oscatter_copy[i]));
116eaedb033SPeter Brune 
1179566063dSJacob Faibussowitsch       PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
1189566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nasm->n, &nasm->subsnes));
119111ade9eSPeter Brune       for (i = 0; i < nasm->n; i++) {
120de54d9edSStefano Zampini         PetscCall(SNESCreate(PetscObjectComm((PetscObject)subdms[i]), &nasm->subsnes[i]));
1219566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)nasm->subsnes[i], (PetscObject)snes, 1));
1229566063dSJacob Faibussowitsch         PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i], optionsprefix));
1239566063dSJacob Faibussowitsch         PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i], "sub_"));
1249566063dSJacob Faibussowitsch         PetscCall(SNESSetDM(nasm->subsnes[i], subdms[i]));
125*2a8381b2SBarry Smith         if (snes->ops->ctxcompute) {
126*2a8381b2SBarry Smith           PetscCall(SNESSetComputeApplicationContext(nasm->subsnes[i], snes->ops->ctxcompute, snes->ops->ctxdestroy));
1272891b2c6SStefano Zampini         } else {
128*2a8381b2SBarry Smith           PetscCtx ctx;
1292891b2c6SStefano Zampini 
1302891b2c6SStefano Zampini           PetscCall(SNESGetApplicationContext(snes, &ctx));
1312891b2c6SStefano Zampini           PetscCall(SNESSetApplicationContext(nasm->subsnes[i], ctx));
1322891b2c6SStefano Zampini         }
1339566063dSJacob Faibussowitsch         PetscCall(SNESSetFromOptions(nasm->subsnes[i]));
1349566063dSJacob Faibussowitsch         PetscCall(DMDestroy(&subdms[i]));
135111ade9eSPeter Brune       }
1369566063dSJacob Faibussowitsch       PetscCall(PetscFree(subdms));
137ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Cannot construct local problems automatically without a DM!");
138ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must set subproblems manually if there is no DM!");
139111ade9eSPeter Brune   /* allocate the global vectors */
14048a46eb9SPierre Jolivet   if (!nasm->x) PetscCall(PetscCalloc1(nasm->n, &nasm->x));
14148a46eb9SPierre Jolivet   if (!nasm->xl) PetscCall(PetscCalloc1(nasm->n, &nasm->xl));
14248a46eb9SPierre Jolivet   if (!nasm->y) PetscCall(PetscCalloc1(nasm->n, &nasm->y));
14348a46eb9SPierre Jolivet   if (!nasm->b) PetscCall(PetscCalloc1(nasm->n, &nasm->b));
144111ade9eSPeter Brune 
145111ade9eSPeter Brune   for (i = 0; i < nasm->n; i++) {
1469566063dSJacob Faibussowitsch     PetscCall(SNESGetFunction(nasm->subsnes[i], &F, NULL, NULL));
1479566063dSJacob Faibussowitsch     if (!nasm->x[i]) PetscCall(VecDuplicate(F, &nasm->x[i]));
1489566063dSJacob Faibussowitsch     if (!nasm->y[i]) PetscCall(VecDuplicate(F, &nasm->y[i]));
1499566063dSJacob Faibussowitsch     if (!nasm->b[i]) PetscCall(VecDuplicate(F, &nasm->b[i]));
15076857b2aSPeter Brune     if (!nasm->xl[i]) {
1519566063dSJacob Faibussowitsch       PetscCall(SNESGetDM(nasm->subsnes[i], &subdm));
1529566063dSJacob Faibussowitsch       PetscCall(DMCreateLocalVector(subdm, &nasm->xl[i]));
1539566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalHookAdd(subdm, DMGlobalToLocalSubDomainDirichletHook_Private, NULL, nasm->xl[i]));
154111ade9eSPeter Brune     }
15561ba4676SBarry Smith   }
156602bec5dSPeter Brune   if (nasm->finaljacobian) {
1579566063dSJacob Faibussowitsch     PetscCall(SNESSetUpMatrices(snes));
15848a46eb9SPierre Jolivet     if (nasm->fjtype == 2) PetscCall(VecDuplicate(snes->vec_sol, &nasm->xinit));
15948a46eb9SPierre Jolivet     for (i = 0; i < nasm->n; i++) PetscCall(SNESSetUpMatrices(nasm->subsnes[i]));
160602bec5dSPeter Brune   }
1613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
162eaedb033SPeter Brune }
163eaedb033SPeter Brune 
SNESSetFromOptions_NASM(SNES snes,PetscOptionItems PetscOptionsObject)164ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_NASM(SNES snes, PetscOptionItems PetscOptionsObject)
165d71ae5a4SJacob Faibussowitsch {
166111ade9eSPeter Brune   PCASMType  asmtype;
16783dc3634SPierre Jolivet   PetscBool  flg, monflg;
168111ade9eSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
1696e111a19SKarl Rupp 
170eaedb033SPeter Brune   PetscFunctionBegin;
171d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Nonlinear Additive Schwarz options");
1729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_nasm_type", "Type of restriction/extension", "", SNESNASMTypes, (PetscEnum)nasm->type, (PetscEnum *)&asmtype, &flg));
1739566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESNASMSetType(snes, asmtype));
174b20c023fSPeter Brune   flg    = PETSC_FALSE;
175b20c023fSPeter Brune   monflg = PETSC_TRUE;
1769566063dSJacob 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));
1779566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESNASMSetDamping(snes, nasm->damping));
1789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDeprecated("-snes_nasm_sub_view", NULL, "3.15", "Use -snes_view ::ascii_info_detail"));
1799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_nasm_finaljacobian", "Compute the global jacobian of the final iterate (for ASPIN)", "", nasm->finaljacobian, &nasm->finaljacobian, NULL));
1809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-snes_nasm_finaljacobian_type", "The type of the final jacobian computed.", "", SNESNASMFJTypes, 3, SNESNASMFJTypes[0], &nasm->fjtype, NULL));
1819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_nasm_log", "Log times for subSNES solves and restriction", "", monflg, &monflg, &flg));
182b20c023fSPeter Brune   if (flg) {
1839566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister("SNESNASMSubSolve", ((PetscObject)snes)->classid, &nasm->eventsubsolve));
1849566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister("SNESNASMRestrict", ((PetscObject)snes)->classid, &nasm->eventrestrictinterp));
185b20c023fSPeter Brune   }
186d0609cedSBarry Smith   PetscOptionsHeadEnd();
1873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
188eaedb033SPeter Brune }
189eaedb033SPeter Brune 
SNESView_NASM(SNES snes,PetscViewer viewer)190d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
191d71ae5a4SJacob Faibussowitsch {
192b20c023fSPeter Brune   SNES_NASM        *nasm = (SNES_NASM *)snes->data;
193a4f17876SPeter Brune   PetscMPIInt       rank, size;
194dd2fa690SBarry Smith   PetscInt          i, N, bsz;
1959f196a02SMartin Diehl   PetscBool         isascii, isstring;
196b20c023fSPeter Brune   PetscViewer       sviewer;
197ce94432eSBarry Smith   MPI_Comm          comm;
19883dc3634SPierre Jolivet   PetscViewerFormat format;
19983dc3634SPierre Jolivet   const char       *prefix;
200b20c023fSPeter Brune 
201b20c023fSPeter Brune   PetscFunctionBegin;
2029566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
2039f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2049566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
2059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2069566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
207462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&nasm->n, &N, 1, MPIU_INT, MPI_SUM, comm));
2089f196a02SMartin Diehl   if (isascii) {
20963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  total subdomain blocks = %" PetscInt_FMT "\n", N));
2109566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
21183dc3634SPierre Jolivet     if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
212a4f17876SPeter Brune       if (nasm->subsnes) {
2139566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for first block on rank 0:\n"));
2149566063dSJacob Faibussowitsch         PetscCall(SNESGetOptionsPrefix(snes, &prefix));
2159566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Use -%ssnes_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : ""));
2169566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
2179566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
218dd400576SPatrick Sanan         if (rank == 0) {
219b4025f61SBarry Smith           PetscCall(PetscViewerASCIIPushTab(sviewer));
2209566063dSJacob Faibussowitsch           PetscCall(SNESView(nasm->subsnes[0], sviewer));
221b4025f61SBarry Smith           PetscCall(PetscViewerASCIIPopTab(sviewer));
222a4f17876SPeter Brune         }
2239566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2249566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
225a4f17876SPeter Brune       }
226a4f17876SPeter Brune     } else {
227a4f17876SPeter Brune       /* print the solver on each block */
2289566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
229835f2295SStefano Zampini       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] number of local blocks = %" PetscInt_FMT "\n", rank, nasm->n));
2309566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2319566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2329566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for each block is in the following SNES objects:\n"));
2339566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
2349566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n"));
2359566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
236a4f17876SPeter Brune       for (i = 0; i < nasm->n; i++) {
2379566063dSJacob Faibussowitsch         PetscCall(VecGetLocalSize(nasm->x[i], &bsz));
238835f2295SStefano Zampini         PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n", rank, i, bsz));
2399566063dSJacob Faibussowitsch         PetscCall(SNESView(nasm->subsnes[i], sviewer));
2409566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(sviewer, "- - - - - - - - - - - - - - - - - -\n"));
241b20c023fSPeter Brune       }
2429566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2439566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
244a4f17876SPeter Brune     }
245b20c023fSPeter Brune   } else if (isstring) {
24663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " blocks=%" PetscInt_FMT ",type=%s", N, SNESNASMTypes[nasm->type]));
2479566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2489566063dSJacob Faibussowitsch     if (nasm->subsnes && rank == 0) PetscCall(SNESView(nasm->subsnes[0], sviewer));
2499566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
250b20c023fSPeter Brune   }
2513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
252eaedb033SPeter Brune }
253eaedb033SPeter Brune 
254e0331734SPeter Brune /*@
255420bcc1bSBarry Smith   SNESNASMSetType - Set the type of subdomain update used for the nonlinear additive Schwarz solver `SNESNASM`
256e0331734SPeter Brune 
257c3339decSBarry Smith   Logically Collective
258e0331734SPeter Brune 
259e0331734SPeter Brune   Input Parameters:
260f6dfbefdSBarry Smith + snes - the `SNES` context
261f6dfbefdSBarry Smith - type - the type of update, `PC_ASM_BASIC` or `PC_ASM_RESTRICT`
262e0331734SPeter Brune 
263420bcc1bSBarry Smith   Options Database Key:
264420bcc1bSBarry Smith . -snes_nasm_type <basic,restrict> - type of subdomain update used
265420bcc1bSBarry Smith 
266e0331734SPeter Brune   Level: intermediate
267e0331734SPeter Brune 
268420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetType()`, `PCASMSetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType`
269e0331734SPeter Brune @*/
SNESNASMSetType(SNES snes,PCASMType type)270d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetType(SNES snes, PCASMType type)
271d71ae5a4SJacob Faibussowitsch {
272e0331734SPeter Brune   PetscFunctionBegin;
273835f2295SStefano Zampini   PetscTryMethod(snes, "SNESNASMSetType_C", (SNES, PCASMType), (snes, type));
2743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
275e0331734SPeter Brune }
276e0331734SPeter Brune 
SNESNASMSetType_NASM(SNES snes,PCASMType type)277d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetType_NASM(SNES snes, PCASMType type)
278d71ae5a4SJacob Faibussowitsch {
279e0331734SPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
280e0331734SPeter Brune 
281e0331734SPeter Brune   PetscFunctionBegin;
2820b121fc5SBarry Smith   PetscCheck(type == PC_ASM_BASIC || type == PC_ASM_RESTRICT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESNASM only supports basic and restrict types");
283e0331734SPeter Brune   nasm->type = type;
2843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
285e0331734SPeter Brune }
286e0331734SPeter Brune 
287e0331734SPeter Brune /*@
288420bcc1bSBarry Smith   SNESNASMGetType - Get the type of subdomain update used for the nonlinear additive Schwarz solver `SNESNASM`
289e0331734SPeter Brune 
290c3339decSBarry Smith   Logically Collective
291e0331734SPeter Brune 
2922fe279fdSBarry Smith   Input Parameter:
293f6dfbefdSBarry Smith . snes - the `SNES` context
294e0331734SPeter Brune 
2952fe279fdSBarry Smith   Output Parameter:
296e0331734SPeter Brune . type - the type of update
297e0331734SPeter Brune 
298e0331734SPeter Brune   Level: intermediate
299e0331734SPeter Brune 
300420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetType()`, `PCASMGetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType`
301e0331734SPeter Brune @*/
SNESNASMGetType(SNES snes,PCASMType * type)302d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetType(SNES snes, PCASMType *type)
303d71ae5a4SJacob Faibussowitsch {
304e0331734SPeter Brune   PetscFunctionBegin;
305cac4c232SBarry Smith   PetscUseMethod(snes, "SNESNASMGetType_C", (SNES, PCASMType *), (snes, type));
3063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
307e0331734SPeter Brune }
308e0331734SPeter Brune 
SNESNASMGetType_NASM(SNES snes,PCASMType * type)309d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetType_NASM(SNES snes, PCASMType *type)
310d71ae5a4SJacob Faibussowitsch {
311e0331734SPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
312e0331734SPeter Brune 
313e0331734SPeter Brune   PetscFunctionBegin;
314e0331734SPeter Brune   *type = nasm->type;
3153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
316e0331734SPeter Brune }
317e0331734SPeter Brune 
31876857b2aSPeter Brune /*@
319f6dfbefdSBarry Smith   SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems in the nonlinear additive Schwarz solver
32076857b2aSPeter Brune 
321f6dfbefdSBarry Smith   Logically Collective
32276857b2aSPeter Brune 
32376857b2aSPeter Brune   Input Parameters:
324420bcc1bSBarry Smith + snes     - the `SNES` context
32576857b2aSPeter Brune . n        - the number of local subdomains
32676857b2aSPeter Brune . subsnes  - solvers defined on the local subdomains
32776857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains
32876857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains
32976857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain
33076857b2aSPeter Brune 
33176857b2aSPeter Brune   Level: intermediate
33276857b2aSPeter Brune 
333420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()`
33476857b2aSPeter Brune @*/
SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])335d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetSubdomains(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[])
336d71ae5a4SJacob Faibussowitsch {
337111ade9eSPeter Brune   PetscErrorCode (*f)(SNES, PetscInt, SNES *, VecScatter *, VecScatter *, VecScatter *);
3386e111a19SKarl Rupp 
339eaedb033SPeter Brune   PetscFunctionBegin;
3409566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", &f));
3419566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter));
3423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
343eaedb033SPeter Brune }
344eaedb033SPeter Brune 
SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])345d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[])
346d71ae5a4SJacob Faibussowitsch {
347eaedb033SPeter Brune   PetscInt   i;
348eaedb033SPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
3496e111a19SKarl Rupp 
350eaedb033SPeter Brune   PetscFunctionBegin;
35128b400f6SJacob Faibussowitsch   PetscCheck(!snes->setupcalled, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
352eaedb033SPeter Brune 
353111ade9eSPeter Brune   /* tear down the previously set things */
3549566063dSJacob Faibussowitsch   PetscCall(SNESReset(snes));
355111ade9eSPeter Brune 
356eaedb033SPeter Brune   nasm->n = n;
357111ade9eSPeter Brune   if (oscatter) {
3589566063dSJacob Faibussowitsch     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)oscatter[i]));
359eaedb033SPeter Brune   }
360111ade9eSPeter Brune   if (iscatter) {
3619566063dSJacob Faibussowitsch     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)iscatter[i]));
362eaedb033SPeter Brune   }
363111ade9eSPeter Brune   if (gscatter) {
3649566063dSJacob Faibussowitsch     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)gscatter[i]));
365111ade9eSPeter Brune   }
366111ade9eSPeter Brune   if (oscatter) {
3679566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->oscatter));
3689566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->oscatter_copy));
369eaedb033SPeter Brune     for (i = 0; i < n; i++) {
370111ade9eSPeter Brune       nasm->oscatter[i] = oscatter[i];
3719566063dSJacob Faibussowitsch       PetscCall(VecScatterCopy(oscatter[i], &nasm->oscatter_copy[i]));
372eaedb033SPeter Brune     }
373111ade9eSPeter Brune   }
374111ade9eSPeter Brune   if (iscatter) {
3759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->iscatter));
376ad540459SPierre Jolivet     for (i = 0; i < n; i++) nasm->iscatter[i] = iscatter[i];
377eaedb033SPeter Brune   }
378111ade9eSPeter Brune   if (gscatter) {
3799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->gscatter));
380ad540459SPierre Jolivet     for (i = 0; i < n; i++) nasm->gscatter[i] = gscatter[i];
381eaedb033SPeter Brune   }
382111ade9eSPeter Brune 
383eaedb033SPeter Brune   if (subsnes) {
3849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->subsnes));
385ad540459SPierre Jolivet     for (i = 0; i < n; i++) nasm->subsnes[i] = subsnes[i];
386eaedb033SPeter Brune   }
3873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
388eaedb033SPeter Brune }
389eaedb033SPeter Brune 
39076857b2aSPeter Brune /*@
391f6dfbefdSBarry Smith   SNESNASMGetSubdomains - Get the local subdomain contexts for the nonlinear additive Schwarz solver
39276857b2aSPeter Brune 
393f6dfbefdSBarry Smith   Not Collective but some of the objects returned will be parallel
39476857b2aSPeter Brune 
395f899ff85SJose E. Roman   Input Parameter:
396f6dfbefdSBarry Smith . snes - the `SNES` context
39776857b2aSPeter Brune 
39876857b2aSPeter Brune   Output Parameters:
39976857b2aSPeter Brune + n        - the number of local subdomains
40076857b2aSPeter Brune . subsnes  - solvers defined on the local subdomains
40176857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains
40276857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains
40376857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain
40476857b2aSPeter Brune 
40576857b2aSPeter Brune   Level: intermediate
40676857b2aSPeter Brune 
407420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetSubdomains()`
40876857b2aSPeter Brune @*/
SNESNASMGetSubdomains(SNES snes,PetscInt * n,SNES * subsnes[],VecScatter * iscatter[],VecScatter * oscatter[],VecScatter * gscatter[])409d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetSubdomains(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[])
410d71ae5a4SJacob Faibussowitsch {
41176857b2aSPeter Brune   PetscErrorCode (*f)(SNES, PetscInt *, SNES **, VecScatter **, VecScatter **, VecScatter **);
41276857b2aSPeter Brune 
41376857b2aSPeter Brune   PetscFunctionBegin;
4149566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", &f));
4159566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter));
4163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41776857b2aSPeter Brune }
41876857b2aSPeter Brune 
SNESNASMGetSubdomains_NASM(SNES snes,PetscInt * n,SNES * subsnes[],VecScatter * iscatter[],VecScatter * oscatter[],VecScatter * gscatter[])419d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[])
420d71ae5a4SJacob Faibussowitsch {
42176857b2aSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
42276857b2aSPeter Brune 
42376857b2aSPeter Brune   PetscFunctionBegin;
42476857b2aSPeter Brune   if (n) *n = nasm->n;
42576857b2aSPeter Brune   if (oscatter) *oscatter = nasm->oscatter;
42676857b2aSPeter Brune   if (iscatter) *iscatter = nasm->iscatter;
42776857b2aSPeter Brune   if (gscatter) *gscatter = nasm->gscatter;
42883dc3634SPierre Jolivet   if (subsnes) *subsnes = nasm->subsnes;
4293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43076857b2aSPeter Brune }
43176857b2aSPeter Brune 
43276857b2aSPeter Brune /*@
433f6dfbefdSBarry Smith   SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors for the nonlinear additive Schwarz solver
43476857b2aSPeter Brune 
43576857b2aSPeter Brune   Not Collective
43676857b2aSPeter Brune 
437f899ff85SJose E. Roman   Input Parameter:
438f6dfbefdSBarry Smith . snes - the `SNES` context
43976857b2aSPeter Brune 
44076857b2aSPeter Brune   Output Parameters:
44176857b2aSPeter Brune + n  - the number of local subdomains
44276857b2aSPeter Brune . x  - The subdomain solution vector
44376857b2aSPeter Brune . y  - The subdomain step vector
44476857b2aSPeter Brune . b  - The subdomain RHS vector
44576857b2aSPeter Brune - xl - The subdomain local vectors (ghosted)
44676857b2aSPeter Brune 
44776857b2aSPeter Brune   Level: developer
44876857b2aSPeter Brune 
449420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()`
45076857b2aSPeter Brune @*/
SNESNASMGetSubdomainVecs(SNES snes,PetscInt * n,Vec * x[],Vec * y[],Vec * b[],Vec * xl[])451ce78bad3SBarry Smith PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes, PetscInt *n, Vec *x[], Vec *y[], Vec *b[], Vec *xl[])
452d71ae5a4SJacob Faibussowitsch {
45376857b2aSPeter Brune   PetscErrorCode (*f)(SNES, PetscInt *, Vec **, Vec **, Vec **, Vec **);
45476857b2aSPeter Brune 
45576857b2aSPeter Brune   PetscFunctionBegin;
4569566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", &f));
4579566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, n, x, y, b, xl));
4583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45976857b2aSPeter Brune }
46076857b2aSPeter Brune 
SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt * n,Vec ** x,Vec ** y,Vec ** b,Vec ** xl)461d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl)
462d71ae5a4SJacob Faibussowitsch {
46376857b2aSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
46476857b2aSPeter Brune 
46576857b2aSPeter Brune   PetscFunctionBegin;
46676857b2aSPeter Brune   if (n) *n = nasm->n;
46776857b2aSPeter Brune   if (x) *x = nasm->x;
46876857b2aSPeter Brune   if (y) *y = nasm->y;
46976857b2aSPeter Brune   if (b) *b = nasm->b;
47076857b2aSPeter Brune   if (xl) *xl = nasm->xl;
4713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47276857b2aSPeter Brune }
47376857b2aSPeter Brune 
474d728fb7dSPeter Brune /*@
475f6dfbefdSBarry Smith   SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain Jacobians upon convergence for the
476f6dfbefdSBarry Smith   nonlinear additive Schwarz solver
477d728fb7dSPeter Brune 
478c3339decSBarry Smith   Collective
479d728fb7dSPeter Brune 
480d728fb7dSPeter Brune   Input Parameters:
481f6dfbefdSBarry Smith + snes - the SNES context
482f6dfbefdSBarry Smith - flg  - `PETSC_TRUE` to compute the Jacobians
483d728fb7dSPeter Brune 
484d728fb7dSPeter Brune   Level: developer
485d728fb7dSPeter Brune 
48695452b02SPatrick Sanan   Notes:
487f6dfbefdSBarry Smith   This is used almost exclusively in the implementation of `SNESASPIN`, where the converged subdomain and global Jacobian
488d728fb7dSPeter Brune   is needed at each linear iteration.
489d728fb7dSPeter Brune 
490420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()`
491d728fb7dSPeter Brune @*/
SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg)492d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes, PetscBool flg)
493d71ae5a4SJacob Faibussowitsch {
494d728fb7dSPeter Brune   PetscErrorCode (*f)(SNES, PetscBool);
495d728fb7dSPeter Brune 
496d728fb7dSPeter Brune   PetscFunctionBegin;
4979566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", &f));
4989566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, flg));
4993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
500d728fb7dSPeter Brune }
501d728fb7dSPeter Brune 
SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg)502d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes, PetscBool flg)
503d71ae5a4SJacob Faibussowitsch {
504d728fb7dSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
505d728fb7dSPeter Brune 
506d728fb7dSPeter Brune   PetscFunctionBegin;
507d728fb7dSPeter Brune   nasm->finaljacobian = flg;
5083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
509d728fb7dSPeter Brune }
51076857b2aSPeter Brune 
511610116beSPeter Brune /*@
512f6dfbefdSBarry Smith   SNESNASMSetDamping - Sets the update damping for `SNESNASM` the nonlinear additive Schwarz solver
513610116beSPeter Brune 
51420f4b53cSBarry Smith   Logically Collective
515610116beSPeter Brune 
516610116beSPeter Brune   Input Parameters:
517f6dfbefdSBarry Smith + snes - the `SNES` context
518610116beSPeter Brune - dmp  - damping
519610116beSPeter Brune 
520420bcc1bSBarry Smith   Options Database Key:
521420bcc1bSBarry Smith . -snes_nasm_damping <dmp> - the new solution is obtained as old solution plus `dmp` times (sum of the solutions on the subdomains)
522420bcc1bSBarry Smith 
523610116beSPeter Brune   Level: intermediate
524610116beSPeter Brune 
525420bcc1bSBarry Smith   Note:
52695452b02SPatrick Sanan   The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
5275dfa0f3bSBarry Smith 
528420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetDamping()`
529610116beSPeter Brune @*/
SNESNASMSetDamping(SNES snes,PetscReal dmp)530d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetDamping(SNES snes, PetscReal dmp)
531d71ae5a4SJacob Faibussowitsch {
532610116beSPeter Brune   PetscFunctionBegin;
533835f2295SStefano Zampini   PetscTryMethod(snes, "SNESNASMSetDamping_C", (SNES, PetscReal), (snes, dmp));
5343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
535610116beSPeter Brune }
536610116beSPeter Brune 
SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp)537d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes, PetscReal dmp)
538d71ae5a4SJacob Faibussowitsch {
539610116beSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
540610116beSPeter Brune 
541610116beSPeter Brune   PetscFunctionBegin;
542610116beSPeter Brune   nasm->damping = dmp;
5433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
544610116beSPeter Brune }
545610116beSPeter Brune 
546610116beSPeter Brune /*@
547f6dfbefdSBarry Smith   SNESNASMGetDamping - Gets the update damping for `SNESNASM` the nonlinear additive Schwarz solver
548610116beSPeter Brune 
549610116beSPeter Brune   Not Collective
550610116beSPeter Brune 
551f6dfbefdSBarry Smith   Input Parameter:
552420bcc1bSBarry Smith . snes - the `SNES` context
553f6dfbefdSBarry Smith 
554f6dfbefdSBarry Smith   Output Parameter:
555f6dfbefdSBarry Smith . dmp - damping
556610116beSPeter Brune 
557610116beSPeter Brune   Level: intermediate
558610116beSPeter Brune 
559420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetDamping()`
560610116beSPeter Brune @*/
SNESNASMGetDamping(SNES snes,PetscReal * dmp)561d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetDamping(SNES snes, PetscReal *dmp)
562d71ae5a4SJacob Faibussowitsch {
563610116beSPeter Brune   PetscFunctionBegin;
564cac4c232SBarry Smith   PetscUseMethod(snes, "SNESNASMGetDamping_C", (SNES, PetscReal *), (snes, dmp));
5653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
566610116beSPeter Brune }
567610116beSPeter Brune 
SNESNASMGetDamping_NASM(SNES snes,PetscReal * dmp)568d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes, PetscReal *dmp)
569d71ae5a4SJacob Faibussowitsch {
570610116beSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
571610116beSPeter Brune 
572610116beSPeter Brune   PetscFunctionBegin;
573610116beSPeter Brune   *dmp = nasm->damping;
5743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
575610116beSPeter Brune }
576610116beSPeter Brune 
57714eb1c5cSMatthew G. Knepley /*
57814eb1c5cSMatthew G. Knepley   Input Parameters:
57914eb1c5cSMatthew G. Knepley + snes - The solver
58014eb1c5cSMatthew G. Knepley . B - The RHS vector
58114eb1c5cSMatthew G. Knepley - X - The initial guess
58214eb1c5cSMatthew G. Knepley 
5832fe279fdSBarry Smith   Output Parameter:
58414eb1c5cSMatthew G. Knepley . Y - The solution update
58514eb1c5cSMatthew G. Knepley 
58614eb1c5cSMatthew G. Knepley   TODO: All scatters should be packed into one
58714eb1c5cSMatthew G. Knepley */
SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X)58866976f2fSJacob Faibussowitsch static PetscErrorCode SNESNASMSolveLocal_Private(SNES snes, Vec B, Vec Y, Vec X)
589d71ae5a4SJacob Faibussowitsch {
590eaedb033SPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
591258e1594SPeter Brune   SNES       subsnes;
592eaedb033SPeter Brune   PetscInt   i;
593610116beSPeter Brune   PetscReal  dmp;
594f10b3e88SPatrick Farrell   Vec        Xl, Bl, Yl, Xlloc;
595f10b3e88SPatrick Farrell   VecScatter iscat, oscat, gscat, oscat_copy;
596111ade9eSPeter Brune   DM         dm, subdm;
597e0331734SPeter Brune   PCASMType  type;
5980adebc6cSBarry Smith 
599eaedb033SPeter Brune   PetscFunctionBegin;
6009566063dSJacob Faibussowitsch   PetscCall(SNESNASMGetType(snes, &type));
6019566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
6029566063dSJacob Faibussowitsch   PetscCall(VecSet(Y, 0));
6039566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
604de54d9edSStefano Zampini   for (i = 0; i < nasm->n; i++) { /* scatter the global solution to the overlap solution and the local solution */
605f10b3e88SPatrick Farrell     Xl         = nasm->x[i];
60670c78f05SPeter Brune     Xlloc      = nasm->xl[i];
60770c78f05SPeter Brune     oscat      = nasm->oscatter[i];
608f10b3e88SPatrick Farrell     oscat_copy = nasm->oscatter_copy[i];
609f10b3e88SPatrick Farrell     gscat      = nasm->gscatter[i];
6109566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD));
6119566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
612de54d9edSStefano Zampini 
61370c78f05SPeter Brune     if (B) {
61470c78f05SPeter Brune       /* scatter the RHS to the local RHS */
61570c78f05SPeter Brune       Bl = nasm->b[i];
6169566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD));
61770c78f05SPeter Brune     }
61870c78f05SPeter Brune   }
6199566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
620b20c023fSPeter Brune 
6219566063dSJacob Faibussowitsch   if (nasm->eventsubsolve) PetscCall(PetscLogEventBegin(nasm->eventsubsolve, snes, 0, 0, 0));
62270c78f05SPeter Brune   for (i = 0; i < nasm->n; i++) {
623de54d9edSStefano Zampini     PetscErrorCode (*bl)(DM, Vec, void *);
624de54d9edSStefano Zampini     void *bctx;
625de54d9edSStefano Zampini 
62670c78f05SPeter Brune     Xl      = nasm->x[i];
62770c78f05SPeter Brune     Xlloc   = nasm->xl[i];
62870c78f05SPeter Brune     Yl      = nasm->y[i];
629258e1594SPeter Brune     subsnes = nasm->subsnes[i];
6309566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(subsnes, &subdm));
631111ade9eSPeter Brune     iscat      = nasm->iscatter[i];
632111ade9eSPeter Brune     oscat      = nasm->oscatter[i];
633f10b3e88SPatrick Farrell     oscat_copy = nasm->oscatter_copy[i];
634111ade9eSPeter Brune     gscat      = nasm->gscatter[i];
6359566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD));
6369566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
63724b7f281SPeter Brune     if (B) {
63824b7f281SPeter Brune       Bl = nasm->b[i];
6399566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD));
640ed3c11a9SPeter Brune     } else Bl = NULL;
641f10b3e88SPatrick Farrell 
642de54d9edSStefano Zampini     PetscCall(SNESGetDM(subsnes, &subdm));
643de54d9edSStefano Zampini     PetscCall(DMSNESGetBoundaryLocal(subdm, &bl, &bctx));
644de54d9edSStefano Zampini     if (bl) PetscCall((*bl)(subdm, Xlloc, bctx));
645de54d9edSStefano Zampini 
6469566063dSJacob Faibussowitsch     PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm));
6479566063dSJacob Faibussowitsch     PetscCall(VecCopy(Xl, Yl));
6489566063dSJacob Faibussowitsch     PetscCall(SNESSolve(subsnes, Bl, Xl));
6499566063dSJacob Faibussowitsch     PetscCall(VecAYPX(Yl, -1.0, Xl));
6509566063dSJacob Faibussowitsch     PetscCall(VecScale(Yl, nasm->damping));
651e0331734SPeter Brune     if (type == PC_ASM_BASIC) {
6529566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
6539566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
654e0331734SPeter Brune     } else if (type == PC_ASM_RESTRICT) {
6559566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
6569566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
657ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Only basic and restrict types are supported for SNESNASM");
658eaedb033SPeter Brune   }
6599566063dSJacob Faibussowitsch   if (nasm->eventsubsolve) PetscCall(PetscLogEventEnd(nasm->eventsubsolve, snes, 0, 0, 0));
6609566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
6611baa6e33SBarry Smith   if (nasm->weight_set) PetscCall(VecPointwiseMult(Y, Y, nasm->weight));
6629566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
6639566063dSJacob Faibussowitsch   PetscCall(SNESNASMGetDamping(snes, &dmp));
6649566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, dmp, Y));
6653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
666eaedb033SPeter Brune }
667eaedb033SPeter Brune 
SNESNASMComputeFinalJacobian_Private(SNES snes,Vec Xfinal)668d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal)
669d71ae5a4SJacob Faibussowitsch {
670602bec5dSPeter Brune   Vec        X    = Xfinal;
671d728fb7dSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
672d728fb7dSPeter Brune   SNES       subsnes;
673ca44f815SPeter Brune   PetscInt   i, lag = 1;
674e59f0a30SPeter Brune   Vec        Xlloc, Xl, Fl, F;
675d728fb7dSPeter Brune   VecScatter oscat, gscat;
676d728fb7dSPeter Brune   DM         dm, subdm;
677d1e9a80fSBarry Smith 
678d728fb7dSPeter Brune   PetscFunctionBegin;
679602bec5dSPeter Brune   if (nasm->fjtype == 2) X = nasm->xinit;
680e59f0a30SPeter Brune   F = snes->vec_func;
6819566063dSJacob Faibussowitsch   if (snes->normschedule == SNES_NORM_NONE) PetscCall(SNESComputeFunction(snes, X, F));
6829566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
68376c63389SBarry Smith   SNESCheckJacobianDomainError(snes);
6849566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
6859566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
686602bec5dSPeter Brune   if (nasm->fjtype != 1) {
687d728fb7dSPeter Brune     for (i = 0; i < nasm->n; i++) {
688d728fb7dSPeter Brune       Xlloc = nasm->xl[i];
689d728fb7dSPeter Brune       gscat = nasm->gscatter[i];
6909566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
691602bec5dSPeter Brune     }
692d728fb7dSPeter Brune   }
6939566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
694d728fb7dSPeter Brune   for (i = 0; i < nasm->n; i++) {
695e59f0a30SPeter Brune     Fl      = nasm->subsnes[i]->vec_func;
696d728fb7dSPeter Brune     Xl      = nasm->x[i];
697d728fb7dSPeter Brune     Xlloc   = nasm->xl[i];
698d728fb7dSPeter Brune     subsnes = nasm->subsnes[i];
699d728fb7dSPeter Brune     oscat   = nasm->oscatter[i];
700d728fb7dSPeter Brune     gscat   = nasm->gscatter[i];
7019566063dSJacob Faibussowitsch     if (nasm->fjtype != 1) PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
7029566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(subsnes, &subdm));
7039566063dSJacob Faibussowitsch     PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm));
704602bec5dSPeter Brune     if (nasm->fjtype != 1) {
7059566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(subdm, Xlloc, INSERT_VALUES, Xl));
7069566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(subdm, Xlloc, INSERT_VALUES, Xl));
707602bec5dSPeter Brune     }
708ca44f815SPeter Brune     if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2;
709ca44f815SPeter Brune     else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian;
7109566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(subsnes, Xl, Fl));
7119566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(subsnes, Xl, subsnes->jacobian, subsnes->jacobian_pre));
712ca44f815SPeter Brune     if (lag > 1) subsnes->lagjacobian = lag;
713d728fb7dSPeter Brune   }
7143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
715d728fb7dSPeter Brune }
716d728fb7dSPeter Brune 
SNESSolve_NASM(SNES snes)717d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NASM(SNES snes)
718d71ae5a4SJacob Faibussowitsch {
719eaedb033SPeter Brune   Vec              F;
720eaedb033SPeter Brune   Vec              X;
721eaedb033SPeter Brune   Vec              B;
722111ade9eSPeter Brune   Vec              Y;
723eaedb033SPeter Brune   PetscInt         i;
724ed3c11a9SPeter Brune   PetscReal        fnorm = 0.0;
725365a6726SPeter Brune   SNESNormSchedule normschedule;
726d728fb7dSPeter Brune   SNES_NASM       *nasm = (SNES_NASM *)snes->data;
727eaedb033SPeter Brune 
728eaedb033SPeter Brune   PetscFunctionBegin;
7290b121fc5SBarry 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);
730c579b300SPatrick Farrell 
7319566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
732eaedb033SPeter Brune   X = snes->vec_sol;
733111ade9eSPeter Brune   Y = snes->vec_sol_update;
734eaedb033SPeter Brune   F = snes->vec_func;
735eaedb033SPeter Brune   B = snes->vec_rhs;
736eaedb033SPeter Brune 
7379566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
738eaedb033SPeter Brune   snes->iter = 0;
739eaedb033SPeter Brune   snes->norm = 0.;
7409566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
741eaedb033SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
7429566063dSJacob Faibussowitsch   PetscCall(SNESGetNormSchedule(snes, &normschedule));
743de54d9edSStefano Zampini   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY || !snes->max_its) {
744eaedb033SPeter Brune     /* compute the initial function and preconditioned update delX */
745ac530a7eSPierre Jolivet     if (!snes->vec_func_init_set) PetscCall(SNESComputeFunction(snes, X, F));
746ac530a7eSPierre Jolivet     else snes->vec_func_init_set = PETSC_FALSE;
747eaedb033SPeter Brune 
7489566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
74976c63389SBarry Smith     SNESCheckFunctionDomainError(snes, fnorm);
7509566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
751eaedb033SPeter Brune     snes->iter = 0;
752eaedb033SPeter Brune     snes->norm = fnorm;
7539566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7549566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
755eaedb033SPeter Brune 
756eaedb033SPeter Brune     /* test convergence */
7572d157150SStefano Zampini     PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
7582d157150SStefano Zampini     PetscCall(SNESMonitor(snes, 0, snes->norm));
7593ba16761SJacob Faibussowitsch     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
760eaedb033SPeter Brune   } else {
7619566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7629566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
7632d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
764eaedb033SPeter Brune   }
765eaedb033SPeter Brune 
766eaedb033SPeter Brune   /* Call general purpose update function */
767dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
768602bec5dSPeter Brune   /* copy the initial solution over for later */
7699566063dSJacob Faibussowitsch   if (nasm->fjtype == 2) PetscCall(VecCopy(X, nasm->xinit));
770eaedb033SPeter Brune 
771eaedb033SPeter Brune   for (i = 0; i < snes->max_its; i++) {
7729566063dSJacob Faibussowitsch     PetscCall(SNESNASMSolveLocal_Private(snes, B, Y, X));
773365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
7749566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
7759566063dSJacob Faibussowitsch       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
77676c63389SBarry Smith       SNESCheckFunctionDomainError(snes, fnorm);
777eaedb033SPeter Brune     }
778eaedb033SPeter Brune     /* Monitor convergence */
7799566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
780eaedb033SPeter Brune     snes->iter = i + 1;
781eaedb033SPeter Brune     snes->norm = fnorm;
7829566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7839566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
784eaedb033SPeter Brune     /* Test for convergence */
7852d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm));
7862d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
787d728fb7dSPeter Brune     if (snes->reason) break;
788eaedb033SPeter Brune     /* Call general purpose update function */
789dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
790eaedb033SPeter Brune   }
79107b62357SFande Kong   if (nasm->finaljacobian) {
7929566063dSJacob Faibussowitsch     PetscCall(SNESNASMComputeFinalJacobian_Private(snes, X));
79376c63389SBarry Smith     SNESCheckJacobianDomainError(snes);
79407b62357SFande Kong   }
7953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
796eaedb033SPeter Brune }
797eaedb033SPeter Brune 
798eaedb033SPeter Brune /*MC
7991d27aa22SBarry Smith   SNESNASM - Nonlinear Additive Schwarz solver {cite}`ck02`, {cite}`bruneknepleysmithtu15`
800eaedb033SPeter Brune 
801f6dfbefdSBarry Smith    Options Database Keys:
80270c78603SPeter Brune +  -snes_nasm_log                                                - enable logging events for the communication and solve stages
80370c78603SPeter Brune .  -snes_nasm_type <basic,restrict>                              - type of subdomain update used
804420bcc1bSBarry Smith .  -snes_nasm_damping <dmp>                                      - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
8051d27aa22SBarry Smith .  -snes_nasm_finaljacobian                                      - compute the local and global Jacobians of the final iterate
8061d27aa22SBarry Smith .  -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the Jacobian is calculated at
80770c78603SPeter Brune .  -sub_snes_                                                    - options prefix of the subdomain nonlinear solves
80870c78603SPeter Brune .  -sub_ksp_                                                     - options prefix of the subdomain Krylov solver
80970c78603SPeter Brune -  -sub_pc_                                                      - options prefix of the subdomain preconditioner
81070c78603SPeter Brune 
811eaedb033SPeter Brune    Level: advanced
812eaedb033SPeter Brune 
813f6dfbefdSBarry Smith    Note:
814f6dfbefdSBarry Smith    This is not often used directly as a solver, it converges too slowly. However it works well as a nonlinear preconditioner for
815f6dfbefdSBarry Smith    the `SNESASPIN` solver
816f6dfbefdSBarry Smith 
817f6dfbefdSBarry Smith    Developer Note:
818f6dfbefdSBarry Smith    This is a non-Newton based nonlinear solver that does not directly require a Jacobian; hence the flag snes->usesksp is set to
819f6dfbefdSBarry Smith    false and `SNESView()` and -snes_view do not display a `KSP` object. However, if the flag nasm->finaljacobian is set (for example, if
820f6dfbefdSBarry Smith    `SNESNASM` is used as a nonlinear preconditioner for  `SNESASPIN`) then `SNESSetUpMatrices()` is called to generate the
821f6dfbefdSBarry Smith    Jacobian (needed by `SNESASPIN`)
822f6dfbefdSBarry Smith    and this utilizes the inner `KSP` object for storing the matrices, but the `KSP` is never used for solving a linear system. When `SNESNASM` is
823f6dfbefdSBarry Smith    used by `SNESASPIN` they share the same Jacobian matrices because `SNESSetUp()` (called on the outer `SNESASPIN`) causes the inner `SNES`
824f6dfbefdSBarry Smith    object (in this case `SNESNASM`) to inherit the outer Jacobian matrices.
825951fe5abSBarry Smith 
826420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESNASMSetType()`, `SNESNASMGetType()`, `SNESNASMSetSubdomains()`, `SNESNASMGetSubdomains()`,
827420bcc1bSBarry Smith           `SNESNASMGetSubdomainVecs()`, `SNESNASMSetComputeFinalJacobian()`, `SNESNASMSetDamping()`, `SNESNASMGetDamping()`, `SNESNASMSetWeight()`,
828420bcc1bSBarry Smith           `SNESNASMGetSNES()`, `SNESNASMGetNumber()`
829eaedb033SPeter Brune M*/
830eaedb033SPeter Brune 
SNESCreate_NASM(SNES snes)831d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes)
832d71ae5a4SJacob Faibussowitsch {
833eaedb033SPeter Brune   SNES_NASM *nasm;
834eaedb033SPeter Brune 
835eaedb033SPeter Brune   PetscFunctionBegin;
8364dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&nasm));
837eaedb033SPeter Brune   snes->data = (void *)nasm;
838eaedb033SPeter Brune 
839eaedb033SPeter Brune   nasm->n             = PETSC_DECIDE;
8409e5d0892SLisandro Dalcin   nasm->subsnes       = NULL;
8419e5d0892SLisandro Dalcin   nasm->x             = NULL;
8429e5d0892SLisandro Dalcin   nasm->xl            = NULL;
8439e5d0892SLisandro Dalcin   nasm->y             = NULL;
8449e5d0892SLisandro Dalcin   nasm->b             = NULL;
8459e5d0892SLisandro Dalcin   nasm->oscatter      = NULL;
8469e5d0892SLisandro Dalcin   nasm->oscatter_copy = NULL;
8479e5d0892SLisandro Dalcin   nasm->iscatter      = NULL;
8489e5d0892SLisandro Dalcin   nasm->gscatter      = NULL;
849610116beSPeter Brune   nasm->damping       = 1.;
850111ade9eSPeter Brune 
851111ade9eSPeter Brune   nasm->type          = PC_ASM_BASIC;
852d728fb7dSPeter Brune   nasm->finaljacobian = PETSC_FALSE;
853f10b3e88SPatrick Farrell   nasm->weight_set    = PETSC_FALSE;
854eaedb033SPeter Brune 
855eaedb033SPeter Brune   snes->ops->destroy        = SNESDestroy_NASM;
856eaedb033SPeter Brune   snes->ops->setup          = SNESSetUp_NASM;
857eaedb033SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
858eaedb033SPeter Brune   snes->ops->view           = SNESView_NASM;
859eaedb033SPeter Brune   snes->ops->solve          = SNESSolve_NASM;
860eaedb033SPeter Brune   snes->ops->reset          = SNESReset_NASM;
861eaedb033SPeter Brune 
862eaedb033SPeter Brune   snes->usesksp = PETSC_FALSE;
863efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
864eaedb033SPeter Brune 
8654fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
8664fc747eaSLawrence Mitchell 
867602bec5dSPeter Brune   nasm->fjtype              = 0;
868602bec5dSPeter Brune   nasm->xinit               = NULL;
8690298fd71SBarry Smith   nasm->eventrestrictinterp = 0;
8700298fd71SBarry Smith   nasm->eventsubsolve       = 0;
871b20c023fSPeter Brune 
87277e5a1f9SBarry Smith   PetscCall(SNESParametersInitialize(snes));
87377e5a1f9SBarry Smith   PetscObjectParameterSetDefault(snes, max_funcs, 10000);
87477e5a1f9SBarry Smith   PetscObjectParameterSetDefault(snes, max_its, 10000);
875eaedb033SPeter Brune 
8769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", SNESNASMSetType_NASM));
8779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", SNESNASMGetType_NASM));
8789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", SNESNASMSetSubdomains_NASM));
8799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", SNESNASMGetSubdomains_NASM));
8809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", SNESNASMSetDamping_NASM));
8819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", SNESNASMGetDamping_NASM));
8829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", SNESNASMGetSubdomainVecs_NASM));
8839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", SNESNASMSetComputeFinalJacobian_NASM));
8843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
885eaedb033SPeter Brune }
88699e0435eSBarry Smith 
887448b6425SPatrick Farrell /*@
888448b6425SPatrick Farrell   SNESNASMGetSNES - Gets a subsolver
889448b6425SPatrick Farrell 
89020f4b53cSBarry Smith   Not Collective
891448b6425SPatrick Farrell 
892448b6425SPatrick Farrell   Input Parameters:
89320f4b53cSBarry Smith + snes - the `SNES` context
894448b6425SPatrick Farrell - i    - the number of the subsnes to get
895448b6425SPatrick Farrell 
8962fe279fdSBarry Smith   Output Parameter:
897448b6425SPatrick Farrell . subsnes - the subsolver context
898448b6425SPatrick Farrell 
899448b6425SPatrick Farrell   Level: intermediate
900448b6425SPatrick Farrell 
901420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNASM`, `SNESNASMGetNumber()`
902448b6425SPatrick Farrell @*/
SNESNASMGetSNES(SNES snes,PetscInt i,SNES * subsnes)903d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetSNES(SNES snes, PetscInt i, SNES *subsnes)
904d71ae5a4SJacob Faibussowitsch {
905448b6425SPatrick Farrell   SNES_NASM *nasm = (SNES_NASM *)snes->data;
906448b6425SPatrick Farrell 
907448b6425SPatrick Farrell   PetscFunctionBegin;
9080b121fc5SBarry Smith   PetscCheck(i >= 0 && i < nasm->n, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "No such subsolver");
909448b6425SPatrick Farrell   *subsnes = nasm->subsnes[i];
9103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
911448b6425SPatrick Farrell }
912448b6425SPatrick Farrell 
913448b6425SPatrick Farrell /*@
914448b6425SPatrick Farrell   SNESNASMGetNumber - Gets number of subsolvers
915448b6425SPatrick Farrell 
91620f4b53cSBarry Smith   Not Collective
917448b6425SPatrick Farrell 
9182fe279fdSBarry Smith   Input Parameter:
91920f4b53cSBarry Smith . snes - the `SNES` context
920448b6425SPatrick Farrell 
9212fe279fdSBarry Smith   Output Parameter:
922448b6425SPatrick Farrell . n - the number of subsolvers
923448b6425SPatrick Farrell 
924448b6425SPatrick Farrell   Level: intermediate
925448b6425SPatrick Farrell 
926420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNASM`, `SNESNASMGetSNES()`
927448b6425SPatrick Farrell @*/
SNESNASMGetNumber(SNES snes,PetscInt * n)928d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetNumber(SNES snes, PetscInt *n)
929d71ae5a4SJacob Faibussowitsch {
930448b6425SPatrick Farrell   SNES_NASM *nasm = (SNES_NASM *)snes->data;
931448b6425SPatrick Farrell 
932448b6425SPatrick Farrell   PetscFunctionBegin;
933448b6425SPatrick Farrell   *n = nasm->n;
9343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
935448b6425SPatrick Farrell }
936448b6425SPatrick Farrell 
937f10b3e88SPatrick Farrell /*@
938f10b3e88SPatrick Farrell   SNESNASMSetWeight - Sets weight to use when adding overlapping updates
939f10b3e88SPatrick Farrell 
940f10b3e88SPatrick Farrell   Collective
941f10b3e88SPatrick Farrell 
942f10b3e88SPatrick Farrell   Input Parameters:
94320f4b53cSBarry Smith + snes   - the `SNES` context
944f10b3e88SPatrick Farrell - weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in)
945f10b3e88SPatrick Farrell 
946f10b3e88SPatrick Farrell   Level: intermediate
947f10b3e88SPatrick Farrell 
948420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNASM`
949f10b3e88SPatrick Farrell @*/
SNESNASMSetWeight(SNES snes,Vec weight)950d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetWeight(SNES snes, Vec weight)
951d71ae5a4SJacob Faibussowitsch {
952f10b3e88SPatrick Farrell   SNES_NASM *nasm = (SNES_NASM *)snes->data;
953f10b3e88SPatrick Farrell 
954f10b3e88SPatrick Farrell   PetscFunctionBegin;
9559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&nasm->weight));
956f10b3e88SPatrick Farrell   nasm->weight_set = PETSC_TRUE;
957f10b3e88SPatrick Farrell   nasm->weight     = weight;
9589566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)nasm->weight));
9593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
960f10b3e88SPatrick Farrell }
961