xref: /petsc/src/snes/impls/nasm/nasm.c (revision 2891b2c670da7fd7c7fda60af7b64c4e470daea4)
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 
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 
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 
88d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm, Vec g, InsertMode mode, Vec l, void *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 
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;
105ed3c11a9SPeter Brune   PetscMPIInt size;
106ed3c11a9SPeter Brune   KSP         ksp;
107ed3c11a9SPeter Brune   PC          pc;
108eaedb033SPeter Brune 
109eaedb033SPeter Brune   PetscFunctionBegin;
110eaedb033SPeter Brune   if (!nasm->subsnes) {
1119566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(snes, &dm));
1120a696f66SPeter Brune     if (dm) {
113eaedb033SPeter Brune       nasm->usesdm = PETSC_TRUE;
1149566063dSJacob Faibussowitsch       PetscCall(DMCreateDomainDecomposition(dm, &nasm->n, NULL, NULL, NULL, &subdms));
11528b400f6SJacob Faibussowitsch       PetscCheck(subdms, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
1169566063dSJacob Faibussowitsch       PetscCall(DMCreateDomainDecompositionScatters(dm, nasm->n, subdms, &nasm->iscatter, &nasm->oscatter, &nasm->gscatter));
1179566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nasm->n, &nasm->oscatter_copy));
11848a46eb9SPierre Jolivet       for (i = 0; i < nasm->n; i++) PetscCall(VecScatterCopy(nasm->oscatter[i], &nasm->oscatter_copy[i]));
119eaedb033SPeter Brune 
1209566063dSJacob Faibussowitsch       PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
1219566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nasm->n, &nasm->subsnes));
122111ade9eSPeter Brune       for (i = 0; i < nasm->n; i++) {
1239566063dSJacob Faibussowitsch         PetscCall(SNESCreate(PETSC_COMM_SELF, &nasm->subsnes[i]));
1249566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)nasm->subsnes[i], (PetscObject)snes, 1));
1259566063dSJacob Faibussowitsch         PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i], optionsprefix));
1269566063dSJacob Faibussowitsch         PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i], "sub_"));
1279566063dSJacob Faibussowitsch         PetscCall(SNESSetDM(nasm->subsnes[i], subdms[i]));
1289566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]), &size));
129ed3c11a9SPeter Brune         if (size == 1) {
1309566063dSJacob Faibussowitsch           PetscCall(SNESGetKSP(nasm->subsnes[i], &ksp));
1319566063dSJacob Faibussowitsch           PetscCall(KSPGetPC(ksp, &pc));
1329566063dSJacob Faibussowitsch           PetscCall(KSPSetType(ksp, KSPPREONLY));
1339566063dSJacob Faibussowitsch           PetscCall(PCSetType(pc, PCLU));
134ed3c11a9SPeter Brune         }
135*2891b2c6SStefano Zampini         if (snes->ops->usercompute) {
136*2891b2c6SStefano Zampini           PetscCall(SNESSetComputeApplicationContext(nasm->subsnes[i], snes->ops->usercompute, snes->ops->userdestroy));
137*2891b2c6SStefano Zampini         } else {
138*2891b2c6SStefano Zampini           void *ctx;
139*2891b2c6SStefano Zampini 
140*2891b2c6SStefano Zampini           PetscCall(SNESGetApplicationContext(snes, &ctx));
141*2891b2c6SStefano Zampini           PetscCall(SNESSetApplicationContext(nasm->subsnes[i], ctx));
142*2891b2c6SStefano Zampini         }
1439566063dSJacob Faibussowitsch         PetscCall(SNESSetFromOptions(nasm->subsnes[i]));
1449566063dSJacob Faibussowitsch         PetscCall(DMDestroy(&subdms[i]));
145111ade9eSPeter Brune       }
1469566063dSJacob Faibussowitsch       PetscCall(PetscFree(subdms));
147ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Cannot construct local problems automatically without a DM!");
148ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must set subproblems manually if there is no DM!");
149111ade9eSPeter Brune   /* allocate the global vectors */
15048a46eb9SPierre Jolivet   if (!nasm->x) PetscCall(PetscCalloc1(nasm->n, &nasm->x));
15148a46eb9SPierre Jolivet   if (!nasm->xl) PetscCall(PetscCalloc1(nasm->n, &nasm->xl));
15248a46eb9SPierre Jolivet   if (!nasm->y) PetscCall(PetscCalloc1(nasm->n, &nasm->y));
15348a46eb9SPierre Jolivet   if (!nasm->b) PetscCall(PetscCalloc1(nasm->n, &nasm->b));
154111ade9eSPeter Brune 
155111ade9eSPeter Brune   for (i = 0; i < nasm->n; i++) {
1569566063dSJacob Faibussowitsch     PetscCall(SNESGetFunction(nasm->subsnes[i], &F, NULL, NULL));
1579566063dSJacob Faibussowitsch     if (!nasm->x[i]) PetscCall(VecDuplicate(F, &nasm->x[i]));
1589566063dSJacob Faibussowitsch     if (!nasm->y[i]) PetscCall(VecDuplicate(F, &nasm->y[i]));
1599566063dSJacob Faibussowitsch     if (!nasm->b[i]) PetscCall(VecDuplicate(F, &nasm->b[i]));
16076857b2aSPeter Brune     if (!nasm->xl[i]) {
1619566063dSJacob Faibussowitsch       PetscCall(SNESGetDM(nasm->subsnes[i], &subdm));
1629566063dSJacob Faibussowitsch       PetscCall(DMCreateLocalVector(subdm, &nasm->xl[i]));
1639566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalHookAdd(subdm, DMGlobalToLocalSubDomainDirichletHook_Private, NULL, nasm->xl[i]));
164111ade9eSPeter Brune     }
16561ba4676SBarry Smith   }
166602bec5dSPeter Brune   if (nasm->finaljacobian) {
1679566063dSJacob Faibussowitsch     PetscCall(SNESSetUpMatrices(snes));
16848a46eb9SPierre Jolivet     if (nasm->fjtype == 2) PetscCall(VecDuplicate(snes->vec_sol, &nasm->xinit));
16948a46eb9SPierre Jolivet     for (i = 0; i < nasm->n; i++) PetscCall(SNESSetUpMatrices(nasm->subsnes[i]));
170602bec5dSPeter Brune   }
1713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
172eaedb033SPeter Brune }
173eaedb033SPeter Brune 
174d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NASM(SNES snes, PetscOptionItems *PetscOptionsObject)
175d71ae5a4SJacob Faibussowitsch {
176111ade9eSPeter Brune   PCASMType  asmtype;
17783dc3634SPierre Jolivet   PetscBool  flg, monflg;
178111ade9eSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
1796e111a19SKarl Rupp 
180eaedb033SPeter Brune   PetscFunctionBegin;
181d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Nonlinear Additive Schwarz options");
1829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_nasm_type", "Type of restriction/extension", "", SNESNASMTypes, (PetscEnum)nasm->type, (PetscEnum *)&asmtype, &flg));
1839566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESNASMSetType(snes, asmtype));
184b20c023fSPeter Brune   flg    = PETSC_FALSE;
185b20c023fSPeter Brune   monflg = PETSC_TRUE;
1869566063dSJacob 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));
1879566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESNASMSetDamping(snes, nasm->damping));
1889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDeprecated("-snes_nasm_sub_view", NULL, "3.15", "Use -snes_view ::ascii_info_detail"));
1899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_nasm_finaljacobian", "Compute the global jacobian of the final iterate (for ASPIN)", "", nasm->finaljacobian, &nasm->finaljacobian, NULL));
1909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-snes_nasm_finaljacobian_type", "The type of the final jacobian computed.", "", SNESNASMFJTypes, 3, SNESNASMFJTypes[0], &nasm->fjtype, NULL));
1919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_nasm_log", "Log times for subSNES solves and restriction", "", monflg, &monflg, &flg));
192b20c023fSPeter Brune   if (flg) {
1939566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister("SNESNASMSubSolve", ((PetscObject)snes)->classid, &nasm->eventsubsolve));
1949566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister("SNESNASMRestrict", ((PetscObject)snes)->classid, &nasm->eventrestrictinterp));
195b20c023fSPeter Brune   }
196d0609cedSBarry Smith   PetscOptionsHeadEnd();
1973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
198eaedb033SPeter Brune }
199eaedb033SPeter Brune 
200d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
201d71ae5a4SJacob Faibussowitsch {
202b20c023fSPeter Brune   SNES_NASM        *nasm = (SNES_NASM *)snes->data;
203a4f17876SPeter Brune   PetscMPIInt       rank, size;
204dd2fa690SBarry Smith   PetscInt          i, N, bsz;
205b20c023fSPeter Brune   PetscBool         iascii, isstring;
206b20c023fSPeter Brune   PetscViewer       sviewer;
207ce94432eSBarry Smith   MPI_Comm          comm;
20883dc3634SPierre Jolivet   PetscViewerFormat format;
20983dc3634SPierre Jolivet   const char       *prefix;
210b20c023fSPeter Brune 
211b20c023fSPeter Brune   PetscFunctionBegin;
2129566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
2139566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2149566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
2159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2171c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&nasm->n, &N, 1, MPIU_INT, MPI_SUM, comm));
218b20c023fSPeter Brune   if (iascii) {
21963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  total subdomain blocks = %" PetscInt_FMT "\n", N));
2209566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
22183dc3634SPierre Jolivet     if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
222a4f17876SPeter Brune       if (nasm->subsnes) {
2239566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for first block on rank 0:\n"));
2249566063dSJacob Faibussowitsch         PetscCall(SNESGetOptionsPrefix(snes, &prefix));
2259566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Use -%ssnes_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : ""));
2269566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
2279566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
228dd400576SPatrick Sanan         if (rank == 0) {
2299566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
2309566063dSJacob Faibussowitsch           PetscCall(SNESView(nasm->subsnes[0], sviewer));
2319566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
232a4f17876SPeter Brune         }
2339566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2349566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
235a4f17876SPeter Brune       }
236a4f17876SPeter Brune     } else {
237a4f17876SPeter Brune       /* print the solver on each block */
2389566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
23963a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] number of local blocks = %" PetscInt_FMT "\n", (int)rank, nasm->n));
2409566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2419566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2429566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for each block is in the following SNES objects:\n"));
2439566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
2449566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n"));
2459566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
246a4f17876SPeter Brune       for (i = 0; i < nasm->n; i++) {
2479566063dSJacob Faibussowitsch         PetscCall(VecGetLocalSize(nasm->x[i], &bsz));
24863a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n", (int)rank, i, bsz));
2499566063dSJacob Faibussowitsch         PetscCall(SNESView(nasm->subsnes[i], sviewer));
2509566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(sviewer, "- - - - - - - - - - - - - - - - - -\n"));
251b20c023fSPeter Brune       }
2529566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2539566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
254a4f17876SPeter Brune     }
255b20c023fSPeter Brune   } else if (isstring) {
25663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " blocks=%" PetscInt_FMT ",type=%s", N, SNESNASMTypes[nasm->type]));
2579566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2589566063dSJacob Faibussowitsch     if (nasm->subsnes && rank == 0) PetscCall(SNESView(nasm->subsnes[0], sviewer));
2599566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
260b20c023fSPeter Brune   }
2613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
262eaedb033SPeter Brune }
263eaedb033SPeter Brune 
264e0331734SPeter Brune /*@
265420bcc1bSBarry Smith   SNESNASMSetType - Set the type of subdomain update used for the nonlinear additive Schwarz solver `SNESNASM`
266e0331734SPeter Brune 
267c3339decSBarry Smith   Logically Collective
268e0331734SPeter Brune 
269e0331734SPeter Brune   Input Parameters:
270f6dfbefdSBarry Smith + snes - the `SNES` context
271f6dfbefdSBarry Smith - type - the type of update, `PC_ASM_BASIC` or `PC_ASM_RESTRICT`
272e0331734SPeter Brune 
273420bcc1bSBarry Smith   Options Database Key:
274420bcc1bSBarry Smith . -snes_nasm_type <basic,restrict> - type of subdomain update used
275420bcc1bSBarry Smith 
276e0331734SPeter Brune   Level: intermediate
277e0331734SPeter Brune 
278420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetType()`, `PCASMSetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType`
279e0331734SPeter Brune @*/
280d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetType(SNES snes, PCASMType type)
281d71ae5a4SJacob Faibussowitsch {
282e0331734SPeter Brune   PetscErrorCode (*f)(SNES, PCASMType);
283e0331734SPeter Brune 
284e0331734SPeter Brune   PetscFunctionBegin;
2859566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetType_C", &f));
2869566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, type));
2873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
288e0331734SPeter Brune }
289e0331734SPeter Brune 
290d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetType_NASM(SNES snes, PCASMType type)
291d71ae5a4SJacob Faibussowitsch {
292e0331734SPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
293e0331734SPeter Brune 
294e0331734SPeter Brune   PetscFunctionBegin;
2950b121fc5SBarry Smith   PetscCheck(type == PC_ASM_BASIC || type == PC_ASM_RESTRICT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESNASM only supports basic and restrict types");
296e0331734SPeter Brune   nasm->type = type;
2973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
298e0331734SPeter Brune }
299e0331734SPeter Brune 
300e0331734SPeter Brune /*@
301420bcc1bSBarry Smith   SNESNASMGetType - Get the type of subdomain update used for the nonlinear additive Schwarz solver `SNESNASM`
302e0331734SPeter Brune 
303c3339decSBarry Smith   Logically Collective
304e0331734SPeter Brune 
3052fe279fdSBarry Smith   Input Parameter:
306f6dfbefdSBarry Smith . snes - the `SNES` context
307e0331734SPeter Brune 
3082fe279fdSBarry Smith   Output Parameter:
309e0331734SPeter Brune . type - the type of update
310e0331734SPeter Brune 
311e0331734SPeter Brune   Level: intermediate
312e0331734SPeter Brune 
313420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetType()`, `PCASMGetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType`
314e0331734SPeter Brune @*/
315d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetType(SNES snes, PCASMType *type)
316d71ae5a4SJacob Faibussowitsch {
317e0331734SPeter Brune   PetscFunctionBegin;
318cac4c232SBarry Smith   PetscUseMethod(snes, "SNESNASMGetType_C", (SNES, PCASMType *), (snes, type));
3193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
320e0331734SPeter Brune }
321e0331734SPeter Brune 
322d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetType_NASM(SNES snes, PCASMType *type)
323d71ae5a4SJacob Faibussowitsch {
324e0331734SPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
325e0331734SPeter Brune 
326e0331734SPeter Brune   PetscFunctionBegin;
327e0331734SPeter Brune   *type = nasm->type;
3283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
329e0331734SPeter Brune }
330e0331734SPeter Brune 
33176857b2aSPeter Brune /*@
332f6dfbefdSBarry Smith   SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems in the nonlinear additive Schwarz solver
33376857b2aSPeter Brune 
334f6dfbefdSBarry Smith   Logically Collective
33576857b2aSPeter Brune 
33676857b2aSPeter Brune   Input Parameters:
337420bcc1bSBarry Smith + snes     - the `SNES` context
33876857b2aSPeter Brune . n        - the number of local subdomains
33976857b2aSPeter Brune . subsnes  - solvers defined on the local subdomains
34076857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains
34176857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains
34276857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain
34376857b2aSPeter Brune 
34476857b2aSPeter Brune   Level: intermediate
34576857b2aSPeter Brune 
346420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()`
34776857b2aSPeter Brune @*/
348d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetSubdomains(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[])
349d71ae5a4SJacob Faibussowitsch {
350111ade9eSPeter Brune   PetscErrorCode (*f)(SNES, PetscInt, SNES *, VecScatter *, VecScatter *, VecScatter *);
3516e111a19SKarl Rupp 
352eaedb033SPeter Brune   PetscFunctionBegin;
3539566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", &f));
3549566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter));
3553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
356eaedb033SPeter Brune }
357eaedb033SPeter Brune 
358d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[])
359d71ae5a4SJacob Faibussowitsch {
360eaedb033SPeter Brune   PetscInt   i;
361eaedb033SPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
3626e111a19SKarl Rupp 
363eaedb033SPeter Brune   PetscFunctionBegin;
36428b400f6SJacob Faibussowitsch   PetscCheck(!snes->setupcalled, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
365eaedb033SPeter Brune 
366111ade9eSPeter Brune   /* tear down the previously set things */
3679566063dSJacob Faibussowitsch   PetscCall(SNESReset(snes));
368111ade9eSPeter Brune 
369eaedb033SPeter Brune   nasm->n = n;
370111ade9eSPeter Brune   if (oscatter) {
3719566063dSJacob Faibussowitsch     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)oscatter[i]));
372eaedb033SPeter Brune   }
373111ade9eSPeter Brune   if (iscatter) {
3749566063dSJacob Faibussowitsch     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)iscatter[i]));
375eaedb033SPeter Brune   }
376111ade9eSPeter Brune   if (gscatter) {
3779566063dSJacob Faibussowitsch     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)gscatter[i]));
378111ade9eSPeter Brune   }
379111ade9eSPeter Brune   if (oscatter) {
3809566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->oscatter));
3819566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->oscatter_copy));
382eaedb033SPeter Brune     for (i = 0; i < n; i++) {
383111ade9eSPeter Brune       nasm->oscatter[i] = oscatter[i];
3849566063dSJacob Faibussowitsch       PetscCall(VecScatterCopy(oscatter[i], &nasm->oscatter_copy[i]));
385eaedb033SPeter Brune     }
386111ade9eSPeter Brune   }
387111ade9eSPeter Brune   if (iscatter) {
3889566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->iscatter));
389ad540459SPierre Jolivet     for (i = 0; i < n; i++) nasm->iscatter[i] = iscatter[i];
390eaedb033SPeter Brune   }
391111ade9eSPeter Brune   if (gscatter) {
3929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->gscatter));
393ad540459SPierre Jolivet     for (i = 0; i < n; i++) nasm->gscatter[i] = gscatter[i];
394eaedb033SPeter Brune   }
395111ade9eSPeter Brune 
396eaedb033SPeter Brune   if (subsnes) {
3979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->subsnes));
398ad540459SPierre Jolivet     for (i = 0; i < n; i++) nasm->subsnes[i] = subsnes[i];
399eaedb033SPeter Brune   }
4003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
401eaedb033SPeter Brune }
402eaedb033SPeter Brune 
40376857b2aSPeter Brune /*@
404f6dfbefdSBarry Smith   SNESNASMGetSubdomains - Get the local subdomain contexts for the nonlinear additive Schwarz solver
40576857b2aSPeter Brune 
406f6dfbefdSBarry Smith   Not Collective but some of the objects returned will be parallel
40776857b2aSPeter Brune 
408f899ff85SJose E. Roman   Input Parameter:
409f6dfbefdSBarry Smith . snes - the `SNES` context
41076857b2aSPeter Brune 
41176857b2aSPeter Brune   Output Parameters:
41276857b2aSPeter Brune + n        - the number of local subdomains
41376857b2aSPeter Brune . subsnes  - solvers defined on the local subdomains
41476857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains
41576857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains
41676857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain
41776857b2aSPeter Brune 
41876857b2aSPeter Brune   Level: intermediate
41976857b2aSPeter Brune 
420420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetSubdomains()`
42176857b2aSPeter Brune @*/
422d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetSubdomains(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[])
423d71ae5a4SJacob Faibussowitsch {
42476857b2aSPeter Brune   PetscErrorCode (*f)(SNES, PetscInt *, SNES **, VecScatter **, VecScatter **, VecScatter **);
42576857b2aSPeter Brune 
42676857b2aSPeter Brune   PetscFunctionBegin;
4279566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", &f));
4289566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter));
4293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43076857b2aSPeter Brune }
43176857b2aSPeter Brune 
432d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[])
433d71ae5a4SJacob Faibussowitsch {
43476857b2aSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
43576857b2aSPeter Brune 
43676857b2aSPeter Brune   PetscFunctionBegin;
43776857b2aSPeter Brune   if (n) *n = nasm->n;
43876857b2aSPeter Brune   if (oscatter) *oscatter = nasm->oscatter;
43976857b2aSPeter Brune   if (iscatter) *iscatter = nasm->iscatter;
44076857b2aSPeter Brune   if (gscatter) *gscatter = nasm->gscatter;
44183dc3634SPierre Jolivet   if (subsnes) *subsnes = nasm->subsnes;
4423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44376857b2aSPeter Brune }
44476857b2aSPeter Brune 
44576857b2aSPeter Brune /*@
446f6dfbefdSBarry Smith   SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors for the nonlinear additive Schwarz solver
44776857b2aSPeter Brune 
44876857b2aSPeter Brune   Not Collective
44976857b2aSPeter Brune 
450f899ff85SJose E. Roman   Input Parameter:
451f6dfbefdSBarry Smith . snes - the `SNES` context
45276857b2aSPeter Brune 
45376857b2aSPeter Brune   Output Parameters:
45476857b2aSPeter Brune + n  - the number of local subdomains
45576857b2aSPeter Brune . x  - The subdomain solution vector
45676857b2aSPeter Brune . y  - The subdomain step vector
45776857b2aSPeter Brune . b  - The subdomain RHS vector
45876857b2aSPeter Brune - xl - The subdomain local vectors (ghosted)
45976857b2aSPeter Brune 
46076857b2aSPeter Brune   Level: developer
46176857b2aSPeter Brune 
462420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()`
46376857b2aSPeter Brune @*/
464d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl)
465d71ae5a4SJacob Faibussowitsch {
46676857b2aSPeter Brune   PetscErrorCode (*f)(SNES, PetscInt *, Vec **, Vec **, Vec **, Vec **);
46776857b2aSPeter Brune 
46876857b2aSPeter Brune   PetscFunctionBegin;
4699566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", &f));
4709566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, n, x, y, b, xl));
4713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47276857b2aSPeter Brune }
47376857b2aSPeter Brune 
474d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl)
475d71ae5a4SJacob Faibussowitsch {
47676857b2aSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
47776857b2aSPeter Brune 
47876857b2aSPeter Brune   PetscFunctionBegin;
47976857b2aSPeter Brune   if (n) *n = nasm->n;
48076857b2aSPeter Brune   if (x) *x = nasm->x;
48176857b2aSPeter Brune   if (y) *y = nasm->y;
48276857b2aSPeter Brune   if (b) *b = nasm->b;
48376857b2aSPeter Brune   if (xl) *xl = nasm->xl;
4843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
48576857b2aSPeter Brune }
48676857b2aSPeter Brune 
487d728fb7dSPeter Brune /*@
488f6dfbefdSBarry Smith   SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain Jacobians upon convergence for the
489f6dfbefdSBarry Smith   nonlinear additive Schwarz solver
490d728fb7dSPeter Brune 
491c3339decSBarry Smith   Collective
492d728fb7dSPeter Brune 
493d728fb7dSPeter Brune   Input Parameters:
494f6dfbefdSBarry Smith + snes - the SNES context
495f6dfbefdSBarry Smith - flg  - `PETSC_TRUE` to compute the Jacobians
496d728fb7dSPeter Brune 
497d728fb7dSPeter Brune   Level: developer
498d728fb7dSPeter Brune 
49995452b02SPatrick Sanan   Notes:
500f6dfbefdSBarry Smith   This is used almost exclusively in the implementation of `SNESASPIN`, where the converged subdomain and global Jacobian
501d728fb7dSPeter Brune   is needed at each linear iteration.
502d728fb7dSPeter Brune 
503420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()`
504d728fb7dSPeter Brune @*/
505d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes, PetscBool flg)
506d71ae5a4SJacob Faibussowitsch {
507d728fb7dSPeter Brune   PetscErrorCode (*f)(SNES, PetscBool);
508d728fb7dSPeter Brune 
509d728fb7dSPeter Brune   PetscFunctionBegin;
5109566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", &f));
5119566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, flg));
5123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
513d728fb7dSPeter Brune }
514d728fb7dSPeter Brune 
515d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes, PetscBool flg)
516d71ae5a4SJacob Faibussowitsch {
517d728fb7dSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
518d728fb7dSPeter Brune 
519d728fb7dSPeter Brune   PetscFunctionBegin;
520d728fb7dSPeter Brune   nasm->finaljacobian = flg;
5213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
522d728fb7dSPeter Brune }
52376857b2aSPeter Brune 
524610116beSPeter Brune /*@
525f6dfbefdSBarry Smith   SNESNASMSetDamping - Sets the update damping for `SNESNASM` the nonlinear additive Schwarz solver
526610116beSPeter Brune 
52720f4b53cSBarry Smith   Logically Collective
528610116beSPeter Brune 
529610116beSPeter Brune   Input Parameters:
530f6dfbefdSBarry Smith + snes - the `SNES` context
531610116beSPeter Brune - dmp  - damping
532610116beSPeter Brune 
533420bcc1bSBarry Smith   Options Database Key:
534420bcc1bSBarry Smith . -snes_nasm_damping <dmp> - the new solution is obtained as old solution plus `dmp` times (sum of the solutions on the subdomains)
535420bcc1bSBarry Smith 
536610116beSPeter Brune   Level: intermediate
537610116beSPeter Brune 
538420bcc1bSBarry Smith   Note:
53995452b02SPatrick Sanan   The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
5405dfa0f3bSBarry Smith 
541420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetDamping()`
542610116beSPeter Brune @*/
543d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetDamping(SNES snes, PetscReal dmp)
544d71ae5a4SJacob Faibussowitsch {
545610116beSPeter Brune   PetscErrorCode (*f)(SNES, PetscReal);
546610116beSPeter Brune 
547610116beSPeter Brune   PetscFunctionBegin;
5489566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetDamping_C", (void (**)(void)) & f));
5499566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, dmp));
5503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
551610116beSPeter Brune }
552610116beSPeter Brune 
553d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes, PetscReal dmp)
554d71ae5a4SJacob Faibussowitsch {
555610116beSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
556610116beSPeter Brune 
557610116beSPeter Brune   PetscFunctionBegin;
558610116beSPeter Brune   nasm->damping = dmp;
5593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
560610116beSPeter Brune }
561610116beSPeter Brune 
562610116beSPeter Brune /*@
563f6dfbefdSBarry Smith   SNESNASMGetDamping - Gets the update damping for `SNESNASM` the nonlinear additive Schwarz solver
564610116beSPeter Brune 
565610116beSPeter Brune   Not Collective
566610116beSPeter Brune 
567f6dfbefdSBarry Smith   Input Parameter:
568420bcc1bSBarry Smith . snes - the `SNES` context
569f6dfbefdSBarry Smith 
570f6dfbefdSBarry Smith   Output Parameter:
571f6dfbefdSBarry Smith . dmp - damping
572610116beSPeter Brune 
573610116beSPeter Brune   Level: intermediate
574610116beSPeter Brune 
575420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetDamping()`
576610116beSPeter Brune @*/
577d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetDamping(SNES snes, PetscReal *dmp)
578d71ae5a4SJacob Faibussowitsch {
579610116beSPeter Brune   PetscFunctionBegin;
580cac4c232SBarry Smith   PetscUseMethod(snes, "SNESNASMGetDamping_C", (SNES, PetscReal *), (snes, dmp));
5813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
582610116beSPeter Brune }
583610116beSPeter Brune 
584d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes, PetscReal *dmp)
585d71ae5a4SJacob Faibussowitsch {
586610116beSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
587610116beSPeter Brune 
588610116beSPeter Brune   PetscFunctionBegin;
589610116beSPeter Brune   *dmp = nasm->damping;
5903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
591610116beSPeter Brune }
592610116beSPeter Brune 
59314eb1c5cSMatthew G. Knepley /*
59414eb1c5cSMatthew G. Knepley   Input Parameters:
59514eb1c5cSMatthew G. Knepley + snes - The solver
59614eb1c5cSMatthew G. Knepley . B - The RHS vector
59714eb1c5cSMatthew G. Knepley - X - The initial guess
59814eb1c5cSMatthew G. Knepley 
5992fe279fdSBarry Smith   Output Parameter:
60014eb1c5cSMatthew G. Knepley . Y - The solution update
60114eb1c5cSMatthew G. Knepley 
60214eb1c5cSMatthew G. Knepley   TODO: All scatters should be packed into one
60314eb1c5cSMatthew G. Knepley */
60466976f2fSJacob Faibussowitsch static PetscErrorCode SNESNASMSolveLocal_Private(SNES snes, Vec B, Vec Y, Vec X)
605d71ae5a4SJacob Faibussowitsch {
606eaedb033SPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
607258e1594SPeter Brune   SNES       subsnes;
608eaedb033SPeter Brune   PetscInt   i;
609610116beSPeter Brune   PetscReal  dmp;
610f10b3e88SPatrick Farrell   Vec        Xl, Bl, Yl, Xlloc;
611f10b3e88SPatrick Farrell   VecScatter iscat, oscat, gscat, oscat_copy;
612111ade9eSPeter Brune   DM         dm, subdm;
613e0331734SPeter Brune   PCASMType  type;
6140adebc6cSBarry Smith 
615eaedb033SPeter Brune   PetscFunctionBegin;
6169566063dSJacob Faibussowitsch   PetscCall(SNESNASMGetType(snes, &type));
6179566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
6189566063dSJacob Faibussowitsch   PetscCall(VecSet(Y, 0));
6199566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
620eaedb033SPeter Brune   for (i = 0; i < nasm->n; i++) {
621f10b3e88SPatrick Farrell     /* scatter the solution to the global solution and the local solution */
622f10b3e88SPatrick Farrell     Xl         = nasm->x[i];
62370c78f05SPeter Brune     Xlloc      = nasm->xl[i];
62470c78f05SPeter Brune     oscat      = nasm->oscatter[i];
625f10b3e88SPatrick Farrell     oscat_copy = nasm->oscatter_copy[i];
626f10b3e88SPatrick Farrell     gscat      = nasm->gscatter[i];
6279566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD));
6289566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
62970c78f05SPeter Brune     if (B) {
63070c78f05SPeter Brune       /* scatter the RHS to the local RHS */
63170c78f05SPeter Brune       Bl = nasm->b[i];
6329566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD));
63370c78f05SPeter Brune     }
63470c78f05SPeter Brune   }
6359566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
636b20c023fSPeter Brune 
6379566063dSJacob Faibussowitsch   if (nasm->eventsubsolve) PetscCall(PetscLogEventBegin(nasm->eventsubsolve, snes, 0, 0, 0));
63870c78f05SPeter Brune   for (i = 0; i < nasm->n; i++) {
63970c78f05SPeter Brune     Xl      = nasm->x[i];
64070c78f05SPeter Brune     Xlloc   = nasm->xl[i];
64170c78f05SPeter Brune     Yl      = nasm->y[i];
642258e1594SPeter Brune     subsnes = nasm->subsnes[i];
6439566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(subsnes, &subdm));
644111ade9eSPeter Brune     iscat      = nasm->iscatter[i];
645111ade9eSPeter Brune     oscat      = nasm->oscatter[i];
646f10b3e88SPatrick Farrell     oscat_copy = nasm->oscatter_copy[i];
647111ade9eSPeter Brune     gscat      = nasm->gscatter[i];
6489566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD));
6499566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
65024b7f281SPeter Brune     if (B) {
65124b7f281SPeter Brune       Bl = nasm->b[i];
6529566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD));
653ed3c11a9SPeter Brune     } else Bl = NULL;
654f10b3e88SPatrick Farrell 
6559566063dSJacob Faibussowitsch     PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm));
6569566063dSJacob Faibussowitsch     PetscCall(VecCopy(Xl, Yl));
6579566063dSJacob Faibussowitsch     PetscCall(SNESSolve(subsnes, Bl, Xl));
6589566063dSJacob Faibussowitsch     PetscCall(VecAYPX(Yl, -1.0, Xl));
6599566063dSJacob Faibussowitsch     PetscCall(VecScale(Yl, nasm->damping));
660e0331734SPeter Brune     if (type == PC_ASM_BASIC) {
6619566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
6629566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
663e0331734SPeter Brune     } else if (type == PC_ASM_RESTRICT) {
6649566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
6659566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
666ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Only basic and restrict types are supported for SNESNASM");
667eaedb033SPeter Brune   }
6689566063dSJacob Faibussowitsch   if (nasm->eventsubsolve) PetscCall(PetscLogEventEnd(nasm->eventsubsolve, snes, 0, 0, 0));
6699566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
6701baa6e33SBarry Smith   if (nasm->weight_set) PetscCall(VecPointwiseMult(Y, Y, nasm->weight));
6719566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
6729566063dSJacob Faibussowitsch   PetscCall(SNESNASMGetDamping(snes, &dmp));
6739566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, dmp, Y));
6743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
675eaedb033SPeter Brune }
676eaedb033SPeter Brune 
677d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal)
678d71ae5a4SJacob Faibussowitsch {
679602bec5dSPeter Brune   Vec        X    = Xfinal;
680d728fb7dSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
681d728fb7dSPeter Brune   SNES       subsnes;
682ca44f815SPeter Brune   PetscInt   i, lag = 1;
683e59f0a30SPeter Brune   Vec        Xlloc, Xl, Fl, F;
684d728fb7dSPeter Brune   VecScatter oscat, gscat;
685d728fb7dSPeter Brune   DM         dm, subdm;
686d1e9a80fSBarry Smith 
687d728fb7dSPeter Brune   PetscFunctionBegin;
688602bec5dSPeter Brune   if (nasm->fjtype == 2) X = nasm->xinit;
689e59f0a30SPeter Brune   F = snes->vec_func;
6909566063dSJacob Faibussowitsch   if (snes->normschedule == SNES_NORM_NONE) PetscCall(SNESComputeFunction(snes, X, F));
6919566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
6929566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
6939566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
694602bec5dSPeter Brune   if (nasm->fjtype != 1) {
695d728fb7dSPeter Brune     for (i = 0; i < nasm->n; i++) {
696d728fb7dSPeter Brune       Xlloc = nasm->xl[i];
697d728fb7dSPeter Brune       gscat = nasm->gscatter[i];
6989566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
699602bec5dSPeter Brune     }
700d728fb7dSPeter Brune   }
7019566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
702d728fb7dSPeter Brune   for (i = 0; i < nasm->n; i++) {
703e59f0a30SPeter Brune     Fl      = nasm->subsnes[i]->vec_func;
704d728fb7dSPeter Brune     Xl      = nasm->x[i];
705d728fb7dSPeter Brune     Xlloc   = nasm->xl[i];
706d728fb7dSPeter Brune     subsnes = nasm->subsnes[i];
707d728fb7dSPeter Brune     oscat   = nasm->oscatter[i];
708d728fb7dSPeter Brune     gscat   = nasm->gscatter[i];
7099566063dSJacob Faibussowitsch     if (nasm->fjtype != 1) PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
7109566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(subsnes, &subdm));
7119566063dSJacob Faibussowitsch     PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm));
712602bec5dSPeter Brune     if (nasm->fjtype != 1) {
7139566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(subdm, Xlloc, INSERT_VALUES, Xl));
7149566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(subdm, Xlloc, INSERT_VALUES, Xl));
715602bec5dSPeter Brune     }
716ca44f815SPeter Brune     if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2;
717ca44f815SPeter Brune     else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian;
7189566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(subsnes, Xl, Fl));
7199566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(subsnes, Xl, subsnes->jacobian, subsnes->jacobian_pre));
720ca44f815SPeter Brune     if (lag > 1) subsnes->lagjacobian = lag;
721d728fb7dSPeter Brune   }
7223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
723d728fb7dSPeter Brune }
724d728fb7dSPeter Brune 
725d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NASM(SNES snes)
726d71ae5a4SJacob Faibussowitsch {
727eaedb033SPeter Brune   Vec              F;
728eaedb033SPeter Brune   Vec              X;
729eaedb033SPeter Brune   Vec              B;
730111ade9eSPeter Brune   Vec              Y;
731eaedb033SPeter Brune   PetscInt         i;
732ed3c11a9SPeter Brune   PetscReal        fnorm = 0.0;
733365a6726SPeter Brune   SNESNormSchedule normschedule;
734d728fb7dSPeter Brune   SNES_NASM       *nasm = (SNES_NASM *)snes->data;
735eaedb033SPeter Brune 
736eaedb033SPeter Brune   PetscFunctionBegin;
737c579b300SPatrick Farrell 
7380b121fc5SBarry 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);
739c579b300SPatrick Farrell 
7409566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
741eaedb033SPeter Brune   X = snes->vec_sol;
742111ade9eSPeter Brune   Y = snes->vec_sol_update;
743eaedb033SPeter Brune   F = snes->vec_func;
744eaedb033SPeter Brune   B = snes->vec_rhs;
745eaedb033SPeter Brune 
7469566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
747eaedb033SPeter Brune   snes->iter = 0;
748eaedb033SPeter Brune   snes->norm = 0.;
7499566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
750eaedb033SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
7519566063dSJacob Faibussowitsch   PetscCall(SNESGetNormSchedule(snes, &normschedule));
752365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
753eaedb033SPeter Brune     /* compute the initial function and preconditioned update delX */
754eaedb033SPeter Brune     if (!snes->vec_func_init_set) {
7559566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
7561aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
757eaedb033SPeter Brune 
7589566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
759422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
7609566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
761eaedb033SPeter Brune     snes->iter = 0;
762eaedb033SPeter Brune     snes->norm = fnorm;
7639566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7649566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
765eaedb033SPeter Brune 
766eaedb033SPeter Brune     /* test convergence */
7672d157150SStefano Zampini     PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
7682d157150SStefano Zampini     PetscCall(SNESMonitor(snes, 0, snes->norm));
7693ba16761SJacob Faibussowitsch     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
770eaedb033SPeter Brune   } else {
7719566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7729566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
7732d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
774eaedb033SPeter Brune   }
775eaedb033SPeter Brune 
776eaedb033SPeter Brune   /* Call general purpose update function */
777dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
778602bec5dSPeter Brune   /* copy the initial solution over for later */
7799566063dSJacob Faibussowitsch   if (nasm->fjtype == 2) PetscCall(VecCopy(X, nasm->xinit));
780eaedb033SPeter Brune 
781eaedb033SPeter Brune   for (i = 0; i < snes->max_its; i++) {
7829566063dSJacob Faibussowitsch     PetscCall(SNESNASMSolveLocal_Private(snes, B, Y, X));
783365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
7849566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
7859566063dSJacob Faibussowitsch       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
786422a814eSBarry Smith       SNESCheckFunctionNorm(snes, fnorm);
787eaedb033SPeter Brune     }
788eaedb033SPeter Brune     /* Monitor convergence */
7899566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
790eaedb033SPeter Brune     snes->iter = i + 1;
791eaedb033SPeter Brune     snes->norm = fnorm;
7929566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7939566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
794eaedb033SPeter Brune     /* Test for convergence */
7952d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm));
7962d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
797d728fb7dSPeter Brune     if (snes->reason) break;
798eaedb033SPeter Brune     /* Call general purpose update function */
799dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
800eaedb033SPeter Brune   }
80107b62357SFande Kong   if (nasm->finaljacobian) {
8029566063dSJacob Faibussowitsch     PetscCall(SNESNASMComputeFinalJacobian_Private(snes, X));
80307b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
80407b62357SFande Kong   }
8053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
806eaedb033SPeter Brune }
807eaedb033SPeter Brune 
808eaedb033SPeter Brune /*MC
809f6dfbefdSBarry Smith   SNESNASM - Nonlinear Additive Schwarz solver
810eaedb033SPeter Brune 
811f6dfbefdSBarry Smith    Options Database Keys:
81270c78603SPeter Brune +  -snes_nasm_log - enable logging events for the communication and solve stages
81370c78603SPeter Brune .  -snes_nasm_type <basic,restrict> - type of subdomain update used
814420bcc1bSBarry Smith .  -snes_nasm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
81570c78603SPeter Brune .  -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate
816150d49b7SHong Zhang .  -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the jacobian is calculated at
81770c78603SPeter Brune .  -sub_snes_ - options prefix of the subdomain nonlinear solves
81870c78603SPeter Brune .  -sub_ksp_ - options prefix of the subdomain Krylov solver
81970c78603SPeter Brune -  -sub_pc_ - options prefix of the subdomain preconditioner
82070c78603SPeter Brune 
821eaedb033SPeter Brune    Level: advanced
822eaedb033SPeter Brune 
823f6dfbefdSBarry Smith    Note:
824f6dfbefdSBarry Smith    This is not often used directly as a solver, it converges too slowly. However it works well as a nonlinear preconditioner for
825f6dfbefdSBarry Smith    the `SNESASPIN` solver
826f6dfbefdSBarry Smith 
827f6dfbefdSBarry Smith    Developer Note:
828f6dfbefdSBarry Smith    This is a non-Newton based nonlinear solver that does not directly require a Jacobian; hence the flag snes->usesksp is set to
829f6dfbefdSBarry Smith        false and `SNESView()` and -snes_view do not display a `KSP` object. However, if the flag nasm->finaljacobian is set (for example, if
830f6dfbefdSBarry Smith        `SNESNASM` is used as a nonlinear preconditioner for  `SNESASPIN`) then `SNESSetUpMatrices()` is called to generate the
831f6dfbefdSBarry Smith        Jacobian (needed by `SNESASPIN`)
832f6dfbefdSBarry 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
833f6dfbefdSBarry Smith        used by `SNESASPIN` they share the same Jacobian matrices because `SNESSetUp()` (called on the outer `SNESASPIN`) causes the inner `SNES`
834f6dfbefdSBarry Smith        object (in this case `SNESNASM`) to inherit the outer Jacobian matrices.
835951fe5abSBarry Smith 
8364f02bc6aSBarry Smith    References:
837606c0280SSatish Balay .  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
8384f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
8394f02bc6aSBarry Smith 
840420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESNASMSetType()`, `SNESNASMGetType()`, `SNESNASMSetSubdomains()`, `SNESNASMGetSubdomains()`,
841420bcc1bSBarry Smith           `SNESNASMGetSubdomainVecs()`, `SNESNASMSetComputeFinalJacobian()`, `SNESNASMSetDamping()`, `SNESNASMGetDamping()`, `SNESNASMSetWeight()`,
842420bcc1bSBarry Smith           `SNESNASMGetSNES()`, `SNESNASMGetNumber()`
843eaedb033SPeter Brune M*/
844eaedb033SPeter Brune 
845d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes)
846d71ae5a4SJacob Faibussowitsch {
847eaedb033SPeter Brune   SNES_NASM *nasm;
848eaedb033SPeter Brune 
849eaedb033SPeter Brune   PetscFunctionBegin;
8504dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&nasm));
851eaedb033SPeter Brune   snes->data = (void *)nasm;
852eaedb033SPeter Brune 
853eaedb033SPeter Brune   nasm->n             = PETSC_DECIDE;
8549e5d0892SLisandro Dalcin   nasm->subsnes       = NULL;
8559e5d0892SLisandro Dalcin   nasm->x             = NULL;
8569e5d0892SLisandro Dalcin   nasm->xl            = NULL;
8579e5d0892SLisandro Dalcin   nasm->y             = NULL;
8589e5d0892SLisandro Dalcin   nasm->b             = NULL;
8599e5d0892SLisandro Dalcin   nasm->oscatter      = NULL;
8609e5d0892SLisandro Dalcin   nasm->oscatter_copy = NULL;
8619e5d0892SLisandro Dalcin   nasm->iscatter      = NULL;
8629e5d0892SLisandro Dalcin   nasm->gscatter      = NULL;
863610116beSPeter Brune   nasm->damping       = 1.;
864111ade9eSPeter Brune 
865111ade9eSPeter Brune   nasm->type          = PC_ASM_BASIC;
866d728fb7dSPeter Brune   nasm->finaljacobian = PETSC_FALSE;
867f10b3e88SPatrick Farrell   nasm->weight_set    = PETSC_FALSE;
868eaedb033SPeter Brune 
869eaedb033SPeter Brune   snes->ops->destroy        = SNESDestroy_NASM;
870eaedb033SPeter Brune   snes->ops->setup          = SNESSetUp_NASM;
871eaedb033SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
872eaedb033SPeter Brune   snes->ops->view           = SNESView_NASM;
873eaedb033SPeter Brune   snes->ops->solve          = SNESSolve_NASM;
874eaedb033SPeter Brune   snes->ops->reset          = SNESReset_NASM;
875eaedb033SPeter Brune 
876eaedb033SPeter Brune   snes->usesksp = PETSC_FALSE;
877efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
878eaedb033SPeter Brune 
8794fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
8804fc747eaSLawrence Mitchell 
881602bec5dSPeter Brune   nasm->fjtype              = 0;
882602bec5dSPeter Brune   nasm->xinit               = NULL;
8830298fd71SBarry Smith   nasm->eventrestrictinterp = 0;
8840298fd71SBarry Smith   nasm->eventsubsolve       = 0;
885b20c023fSPeter Brune 
886eaedb033SPeter Brune   if (!snes->tolerancesset) {
887eaedb033SPeter Brune     snes->max_its   = 10000;
888eaedb033SPeter Brune     snes->max_funcs = 10000;
889eaedb033SPeter Brune   }
890eaedb033SPeter Brune 
8919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", SNESNASMSetType_NASM));
8929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", SNESNASMGetType_NASM));
8939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", SNESNASMSetSubdomains_NASM));
8949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", SNESNASMGetSubdomains_NASM));
8959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", SNESNASMSetDamping_NASM));
8969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", SNESNASMGetDamping_NASM));
8979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", SNESNASMGetSubdomainVecs_NASM));
8989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", SNESNASMSetComputeFinalJacobian_NASM));
8993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
900eaedb033SPeter Brune }
90199e0435eSBarry Smith 
902448b6425SPatrick Farrell /*@
903448b6425SPatrick Farrell   SNESNASMGetSNES - Gets a subsolver
904448b6425SPatrick Farrell 
90520f4b53cSBarry Smith   Not Collective
906448b6425SPatrick Farrell 
907448b6425SPatrick Farrell   Input Parameters:
90820f4b53cSBarry Smith + snes - the `SNES` context
909448b6425SPatrick Farrell - i    - the number of the subsnes to get
910448b6425SPatrick Farrell 
9112fe279fdSBarry Smith   Output Parameter:
912448b6425SPatrick Farrell . subsnes - the subsolver context
913448b6425SPatrick Farrell 
914448b6425SPatrick Farrell   Level: intermediate
915448b6425SPatrick Farrell 
916420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNASM`, `SNESNASMGetNumber()`
917448b6425SPatrick Farrell @*/
918d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetSNES(SNES snes, PetscInt i, SNES *subsnes)
919d71ae5a4SJacob Faibussowitsch {
920448b6425SPatrick Farrell   SNES_NASM *nasm = (SNES_NASM *)snes->data;
921448b6425SPatrick Farrell 
922448b6425SPatrick Farrell   PetscFunctionBegin;
9230b121fc5SBarry Smith   PetscCheck(i >= 0 && i < nasm->n, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "No such subsolver");
924448b6425SPatrick Farrell   *subsnes = nasm->subsnes[i];
9253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
926448b6425SPatrick Farrell }
927448b6425SPatrick Farrell 
928448b6425SPatrick Farrell /*@
929448b6425SPatrick Farrell   SNESNASMGetNumber - Gets number of subsolvers
930448b6425SPatrick Farrell 
93120f4b53cSBarry Smith   Not Collective
932448b6425SPatrick Farrell 
9332fe279fdSBarry Smith   Input Parameter:
93420f4b53cSBarry Smith . snes - the `SNES` context
935448b6425SPatrick Farrell 
9362fe279fdSBarry Smith   Output Parameter:
937448b6425SPatrick Farrell . n - the number of subsolvers
938448b6425SPatrick Farrell 
939448b6425SPatrick Farrell   Level: intermediate
940448b6425SPatrick Farrell 
941420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNASM`, `SNESNASMGetSNES()`
942448b6425SPatrick Farrell @*/
943d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetNumber(SNES snes, PetscInt *n)
944d71ae5a4SJacob Faibussowitsch {
945448b6425SPatrick Farrell   SNES_NASM *nasm = (SNES_NASM *)snes->data;
946448b6425SPatrick Farrell 
947448b6425SPatrick Farrell   PetscFunctionBegin;
948448b6425SPatrick Farrell   *n = nasm->n;
9493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
950448b6425SPatrick Farrell }
951448b6425SPatrick Farrell 
952f10b3e88SPatrick Farrell /*@
953f10b3e88SPatrick Farrell   SNESNASMSetWeight - Sets weight to use when adding overlapping updates
954f10b3e88SPatrick Farrell 
955f10b3e88SPatrick Farrell   Collective
956f10b3e88SPatrick Farrell 
957f10b3e88SPatrick Farrell   Input Parameters:
95820f4b53cSBarry Smith + snes   - the `SNES` context
959f10b3e88SPatrick Farrell - weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in)
960f10b3e88SPatrick Farrell 
961f10b3e88SPatrick Farrell   Level: intermediate
962f10b3e88SPatrick Farrell 
963420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNASM`
964f10b3e88SPatrick Farrell @*/
965d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetWeight(SNES snes, Vec weight)
966d71ae5a4SJacob Faibussowitsch {
967f10b3e88SPatrick Farrell   SNES_NASM *nasm = (SNES_NASM *)snes->data;
968f10b3e88SPatrick Farrell 
969f10b3e88SPatrick Farrell   PetscFunctionBegin;
970f10b3e88SPatrick Farrell 
9719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&nasm->weight));
972f10b3e88SPatrick Farrell   nasm->weight_set = PETSC_TRUE;
973f10b3e88SPatrick Farrell   nasm->weight     = weight;
9749566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)nasm->weight));
975f10b3e88SPatrick Farrell 
9763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
977f10b3e88SPatrick Farrell }
978