xref: /petsc/src/snes/impls/nasm/nasm.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
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 
339371c9d4SSatish Balay static PetscErrorCode SNESReset_NASM(SNES snes) {
34eaedb033SPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
35eaedb033SPeter Brune   PetscInt   i;
366e111a19SKarl Rupp 
37eaedb033SPeter Brune   PetscFunctionBegin;
38eaedb033SPeter Brune   for (i = 0; i < nasm->n; i++) {
399566063dSJacob Faibussowitsch     if (nasm->xl) PetscCall(VecDestroy(&nasm->xl[i]));
409566063dSJacob Faibussowitsch     if (nasm->x) PetscCall(VecDestroy(&nasm->x[i]));
419566063dSJacob Faibussowitsch     if (nasm->y) PetscCall(VecDestroy(&nasm->y[i]));
429566063dSJacob Faibussowitsch     if (nasm->b) PetscCall(VecDestroy(&nasm->b[i]));
43eaedb033SPeter Brune 
449566063dSJacob Faibussowitsch     if (nasm->subsnes) PetscCall(SNESDestroy(&nasm->subsnes[i]));
459566063dSJacob Faibussowitsch     if (nasm->oscatter) PetscCall(VecScatterDestroy(&nasm->oscatter[i]));
469566063dSJacob Faibussowitsch     if (nasm->oscatter_copy) PetscCall(VecScatterDestroy(&nasm->oscatter_copy[i]));
479566063dSJacob Faibussowitsch     if (nasm->iscatter) PetscCall(VecScatterDestroy(&nasm->iscatter[i]));
489566063dSJacob Faibussowitsch     if (nasm->gscatter) PetscCall(VecScatterDestroy(&nasm->gscatter[i]));
49eaedb033SPeter Brune   }
50111ade9eSPeter Brune 
519566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->x));
529566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->xl));
539566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->y));
549566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->b));
55111ade9eSPeter Brune 
569566063dSJacob Faibussowitsch   if (nasm->xinit) PetscCall(VecDestroy(&nasm->xinit));
57602bec5dSPeter Brune 
589566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->subsnes));
599566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->oscatter));
609566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->oscatter_copy));
619566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->iscatter));
629566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->gscatter));
63b20c023fSPeter Brune 
6448a46eb9SPierre Jolivet   if (nasm->weight_set) PetscCall(VecDestroy(&nasm->weight));
65f10b3e88SPatrick Farrell 
66b20c023fSPeter Brune   nasm->eventrestrictinterp = 0;
67b20c023fSPeter Brune   nasm->eventsubsolve       = 0;
68eaedb033SPeter Brune   PetscFunctionReturn(0);
69eaedb033SPeter Brune }
70eaedb033SPeter Brune 
719371c9d4SSatish Balay static PetscErrorCode SNESDestroy_NASM(SNES snes) {
72eaedb033SPeter Brune   PetscFunctionBegin;
739566063dSJacob Faibussowitsch   PetscCall(SNESReset_NASM(snes));
742e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", NULL));
752e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", NULL));
762e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", NULL));
772e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", NULL));
782e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", NULL));
792e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", NULL));
802e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", NULL));
812e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", NULL));
829566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
83eaedb033SPeter Brune   PetscFunctionReturn(0);
84eaedb033SPeter Brune }
85eaedb033SPeter Brune 
869371c9d4SSatish Balay static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm, Vec g, InsertMode mode, Vec l, void *ctx) {
87111ade9eSPeter Brune   Vec bcs = (Vec)ctx;
886e111a19SKarl Rupp 
89111ade9eSPeter Brune   PetscFunctionBegin;
909566063dSJacob Faibussowitsch   PetscCall(VecCopy(bcs, l));
91111ade9eSPeter Brune   PetscFunctionReturn(0);
92111ade9eSPeter Brune }
93111ade9eSPeter Brune 
949371c9d4SSatish Balay static PetscErrorCode SNESSetUp_NASM(SNES snes) {
95eaedb033SPeter Brune   SNES_NASM  *nasm = (SNES_NASM *)snes->data;
9676857b2aSPeter Brune   DM          dm, subdm;
97111ade9eSPeter Brune   DM         *subdms;
98111ade9eSPeter Brune   PetscInt    i;
99eaedb033SPeter Brune   const char *optionsprefix;
100111ade9eSPeter Brune   Vec         F;
101ed3c11a9SPeter Brune   PetscMPIInt size;
102ed3c11a9SPeter Brune   KSP         ksp;
103ed3c11a9SPeter Brune   PC          pc;
104eaedb033SPeter Brune 
105eaedb033SPeter Brune   PetscFunctionBegin;
106eaedb033SPeter Brune   if (!nasm->subsnes) {
1079566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(snes, &dm));
1080a696f66SPeter Brune     if (dm) {
109eaedb033SPeter Brune       nasm->usesdm = PETSC_TRUE;
1109566063dSJacob Faibussowitsch       PetscCall(DMCreateDomainDecomposition(dm, &nasm->n, NULL, NULL, NULL, &subdms));
11128b400f6SJacob Faibussowitsch       PetscCheck(subdms, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
1129566063dSJacob Faibussowitsch       PetscCall(DMCreateDomainDecompositionScatters(dm, nasm->n, subdms, &nasm->iscatter, &nasm->oscatter, &nasm->gscatter));
1139566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nasm->n, &nasm->oscatter_copy));
11448a46eb9SPierre Jolivet       for (i = 0; i < nasm->n; i++) PetscCall(VecScatterCopy(nasm->oscatter[i], &nasm->oscatter_copy[i]));
115eaedb033SPeter Brune 
1169566063dSJacob Faibussowitsch       PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
1179566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nasm->n, &nasm->subsnes));
118111ade9eSPeter Brune       for (i = 0; i < nasm->n; i++) {
1199566063dSJacob Faibussowitsch         PetscCall(SNESCreate(PETSC_COMM_SELF, &nasm->subsnes[i]));
1209566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)nasm->subsnes[i], (PetscObject)snes, 1));
1219566063dSJacob Faibussowitsch         PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i], optionsprefix));
1229566063dSJacob Faibussowitsch         PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i], "sub_"));
1239566063dSJacob Faibussowitsch         PetscCall(SNESSetDM(nasm->subsnes[i], subdms[i]));
1249566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]), &size));
125ed3c11a9SPeter Brune         if (size == 1) {
1269566063dSJacob Faibussowitsch           PetscCall(SNESGetKSP(nasm->subsnes[i], &ksp));
1279566063dSJacob Faibussowitsch           PetscCall(KSPGetPC(ksp, &pc));
1289566063dSJacob Faibussowitsch           PetscCall(KSPSetType(ksp, KSPPREONLY));
1299566063dSJacob Faibussowitsch           PetscCall(PCSetType(pc, PCLU));
130ed3c11a9SPeter Brune         }
1319566063dSJacob Faibussowitsch         PetscCall(SNESSetFromOptions(nasm->subsnes[i]));
1329566063dSJacob Faibussowitsch         PetscCall(DMDestroy(&subdms[i]));
133111ade9eSPeter Brune       }
1349566063dSJacob Faibussowitsch       PetscCall(PetscFree(subdms));
135ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Cannot construct local problems automatically without a DM!");
136ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must set subproblems manually if there is no DM!");
137111ade9eSPeter Brune   /* allocate the global vectors */
13848a46eb9SPierre Jolivet   if (!nasm->x) PetscCall(PetscCalloc1(nasm->n, &nasm->x));
13948a46eb9SPierre Jolivet   if (!nasm->xl) PetscCall(PetscCalloc1(nasm->n, &nasm->xl));
14048a46eb9SPierre Jolivet   if (!nasm->y) PetscCall(PetscCalloc1(nasm->n, &nasm->y));
14148a46eb9SPierre Jolivet   if (!nasm->b) PetscCall(PetscCalloc1(nasm->n, &nasm->b));
142111ade9eSPeter Brune 
143111ade9eSPeter Brune   for (i = 0; i < nasm->n; i++) {
1449566063dSJacob Faibussowitsch     PetscCall(SNESGetFunction(nasm->subsnes[i], &F, NULL, NULL));
1459566063dSJacob Faibussowitsch     if (!nasm->x[i]) PetscCall(VecDuplicate(F, &nasm->x[i]));
1469566063dSJacob Faibussowitsch     if (!nasm->y[i]) PetscCall(VecDuplicate(F, &nasm->y[i]));
1479566063dSJacob Faibussowitsch     if (!nasm->b[i]) PetscCall(VecDuplicate(F, &nasm->b[i]));
14876857b2aSPeter Brune     if (!nasm->xl[i]) {
1499566063dSJacob Faibussowitsch       PetscCall(SNESGetDM(nasm->subsnes[i], &subdm));
1509566063dSJacob Faibussowitsch       PetscCall(DMCreateLocalVector(subdm, &nasm->xl[i]));
1519566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalHookAdd(subdm, DMGlobalToLocalSubDomainDirichletHook_Private, NULL, nasm->xl[i]));
152111ade9eSPeter Brune     }
15361ba4676SBarry Smith   }
154602bec5dSPeter Brune   if (nasm->finaljacobian) {
1559566063dSJacob Faibussowitsch     PetscCall(SNESSetUpMatrices(snes));
15648a46eb9SPierre Jolivet     if (nasm->fjtype == 2) PetscCall(VecDuplicate(snes->vec_sol, &nasm->xinit));
15748a46eb9SPierre Jolivet     for (i = 0; i < nasm->n; i++) PetscCall(SNESSetUpMatrices(nasm->subsnes[i]));
158602bec5dSPeter Brune   }
159eaedb033SPeter Brune   PetscFunctionReturn(0);
160eaedb033SPeter Brune }
161eaedb033SPeter Brune 
1629371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_NASM(SNES snes, PetscOptionItems *PetscOptionsObject) {
163111ade9eSPeter Brune   PCASMType  asmtype;
16483dc3634SPierre Jolivet   PetscBool  flg, monflg;
165111ade9eSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
1666e111a19SKarl Rupp 
167eaedb033SPeter Brune   PetscFunctionBegin;
168d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Nonlinear Additive Schwarz options");
1699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_nasm_type", "Type of restriction/extension", "", SNESNASMTypes, (PetscEnum)nasm->type, (PetscEnum *)&asmtype, &flg));
1709566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESNASMSetType(snes, asmtype));
171b20c023fSPeter Brune   flg    = PETSC_FALSE;
172b20c023fSPeter Brune   monflg = PETSC_TRUE;
1739566063dSJacob 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));
1749566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESNASMSetDamping(snes, nasm->damping));
1759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDeprecated("-snes_nasm_sub_view", NULL, "3.15", "Use -snes_view ::ascii_info_detail"));
1769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_nasm_finaljacobian", "Compute the global jacobian of the final iterate (for ASPIN)", "", nasm->finaljacobian, &nasm->finaljacobian, NULL));
1779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-snes_nasm_finaljacobian_type", "The type of the final jacobian computed.", "", SNESNASMFJTypes, 3, SNESNASMFJTypes[0], &nasm->fjtype, NULL));
1789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_nasm_log", "Log times for subSNES solves and restriction", "", monflg, &monflg, &flg));
179b20c023fSPeter Brune   if (flg) {
1809566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister("SNESNASMSubSolve", ((PetscObject)snes)->classid, &nasm->eventsubsolve));
1819566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister("SNESNASMRestrict", ((PetscObject)snes)->classid, &nasm->eventrestrictinterp));
182b20c023fSPeter Brune   }
183d0609cedSBarry Smith   PetscOptionsHeadEnd();
184eaedb033SPeter Brune   PetscFunctionReturn(0);
185eaedb033SPeter Brune }
186eaedb033SPeter Brune 
1879371c9d4SSatish Balay static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) {
188b20c023fSPeter Brune   SNES_NASM        *nasm = (SNES_NASM *)snes->data;
189a4f17876SPeter Brune   PetscMPIInt       rank, size;
190dd2fa690SBarry Smith   PetscInt          i, N, bsz;
191b20c023fSPeter Brune   PetscBool         iascii, isstring;
192b20c023fSPeter Brune   PetscViewer       sviewer;
193ce94432eSBarry Smith   MPI_Comm          comm;
19483dc3634SPierre Jolivet   PetscViewerFormat format;
19583dc3634SPierre Jolivet   const char       *prefix;
196b20c023fSPeter Brune 
197b20c023fSPeter Brune   PetscFunctionBegin;
1989566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
1999566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2009566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
2019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2031c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&nasm->n, &N, 1, MPIU_INT, MPI_SUM, comm));
204b20c023fSPeter Brune   if (iascii) {
20563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  total subdomain blocks = %" PetscInt_FMT "\n", N));
2069566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
20783dc3634SPierre Jolivet     if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
208a4f17876SPeter Brune       if (nasm->subsnes) {
2099566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for first block on rank 0:\n"));
2109566063dSJacob Faibussowitsch         PetscCall(SNESGetOptionsPrefix(snes, &prefix));
2119566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Use -%ssnes_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : ""));
2129566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
2139566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
214dd400576SPatrick Sanan         if (rank == 0) {
2159566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
2169566063dSJacob Faibussowitsch           PetscCall(SNESView(nasm->subsnes[0], sviewer));
2179566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
218a4f17876SPeter Brune         }
2199566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2209566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
221a4f17876SPeter Brune       }
222a4f17876SPeter Brune     } else {
223a4f17876SPeter Brune       /* print the solver on each block */
2249566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
22563a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] number of local blocks = %" PetscInt_FMT "\n", (int)rank, nasm->n));
2269566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2279566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2289566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for each block is in the following SNES objects:\n"));
2299566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
2309566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n"));
2319566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
232a4f17876SPeter Brune       for (i = 0; i < nasm->n; i++) {
2339566063dSJacob Faibussowitsch         PetscCall(VecGetLocalSize(nasm->x[i], &bsz));
23463a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n", (int)rank, i, bsz));
2359566063dSJacob Faibussowitsch         PetscCall(SNESView(nasm->subsnes[i], sviewer));
2369566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(sviewer, "- - - - - - - - - - - - - - - - - -\n"));
237b20c023fSPeter Brune       }
2389566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2399566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2409566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
241a4f17876SPeter Brune     }
242b20c023fSPeter Brune   } else if (isstring) {
24363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " blocks=%" PetscInt_FMT ",type=%s", N, SNESNASMTypes[nasm->type]));
2449566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2459566063dSJacob Faibussowitsch     if (nasm->subsnes && rank == 0) PetscCall(SNESView(nasm->subsnes[0], sviewer));
2469566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
247b20c023fSPeter Brune   }
248eaedb033SPeter Brune   PetscFunctionReturn(0);
249eaedb033SPeter Brune }
250eaedb033SPeter Brune 
251e0331734SPeter Brune /*@
252*f6dfbefdSBarry Smith    SNESNASMSetType - Set the type of subdomain update used for the nonlinear additive Schwarz solver
253e0331734SPeter Brune 
254*f6dfbefdSBarry Smith    Logically Collective on snes
255e0331734SPeter Brune 
256e0331734SPeter Brune    Input Parameters:
257*f6dfbefdSBarry Smith +  snes - the `SNES` context
258*f6dfbefdSBarry Smith -  type - the type of update, `PC_ASM_BASIC` or `PC_ASM_RESTRICT`
259e0331734SPeter Brune 
260e0331734SPeter Brune    Level: intermediate
261e0331734SPeter Brune 
262*f6dfbefdSBarry Smith .seealso: `SNESNASM`, `SNESNASMGetType()`, `PCASMSetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType`
263e0331734SPeter Brune @*/
2649371c9d4SSatish Balay PetscErrorCode SNESNASMSetType(SNES snes, PCASMType type) {
265e0331734SPeter Brune   PetscErrorCode (*f)(SNES, PCASMType);
266e0331734SPeter Brune 
267e0331734SPeter Brune   PetscFunctionBegin;
2689566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetType_C", &f));
2699566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, type));
270e0331734SPeter Brune   PetscFunctionReturn(0);
271e0331734SPeter Brune }
272e0331734SPeter Brune 
2739371c9d4SSatish Balay static PetscErrorCode SNESNASMSetType_NASM(SNES snes, PCASMType type) {
274e0331734SPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
275e0331734SPeter Brune 
276e0331734SPeter Brune   PetscFunctionBegin;
2770b121fc5SBarry Smith   PetscCheck(type == PC_ASM_BASIC || type == PC_ASM_RESTRICT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESNASM only supports basic and restrict types");
278e0331734SPeter Brune   nasm->type = type;
279e0331734SPeter Brune   PetscFunctionReturn(0);
280e0331734SPeter Brune }
281e0331734SPeter Brune 
282e0331734SPeter Brune /*@
283*f6dfbefdSBarry Smith    SNESNASMGetType - Get the type of subdomain update used for the nonlinear additive Schwarz solver
284e0331734SPeter Brune 
285*f6dfbefdSBarry Smith    Logically Collective on snes
286e0331734SPeter Brune 
287e0331734SPeter Brune    Input Parameters:
288*f6dfbefdSBarry Smith .  snes - the `SNES` context
289e0331734SPeter Brune 
290e0331734SPeter Brune    Output Parameters:
291e0331734SPeter Brune .  type - the type of update
292e0331734SPeter Brune 
293e0331734SPeter Brune    Level: intermediate
294e0331734SPeter Brune 
295*f6dfbefdSBarry Smith .seealso: `SNESNASM`, `SNESNASMSetType()`, `PCASMGetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType`
296e0331734SPeter Brune @*/
2979371c9d4SSatish Balay PetscErrorCode SNESNASMGetType(SNES snes, PCASMType *type) {
298e0331734SPeter Brune   PetscFunctionBegin;
299cac4c232SBarry Smith   PetscUseMethod(snes, "SNESNASMGetType_C", (SNES, PCASMType *), (snes, type));
300e0331734SPeter Brune   PetscFunctionReturn(0);
301e0331734SPeter Brune }
302e0331734SPeter Brune 
3039371c9d4SSatish Balay static PetscErrorCode SNESNASMGetType_NASM(SNES snes, PCASMType *type) {
304e0331734SPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
305e0331734SPeter Brune 
306e0331734SPeter Brune   PetscFunctionBegin;
307e0331734SPeter Brune   *type = nasm->type;
308e0331734SPeter Brune   PetscFunctionReturn(0);
309e0331734SPeter Brune }
310e0331734SPeter Brune 
31176857b2aSPeter Brune /*@
312*f6dfbefdSBarry Smith    SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems in the nonlinear additive Schwarz solver
31376857b2aSPeter Brune 
314*f6dfbefdSBarry Smith    Logically Collective
31576857b2aSPeter Brune 
31676857b2aSPeter Brune    Input Parameters:
317*f6dfbefdSBarry Smith +  snes - the SNES context
31876857b2aSPeter Brune .  n - the number of local subdomains
31976857b2aSPeter Brune .  subsnes - solvers defined on the local subdomains
32076857b2aSPeter Brune .  iscatter - scatters into the nonoverlapping portions of the local subdomains
32176857b2aSPeter Brune .  oscatter - scatters into the overlapping portions of the local subdomains
32276857b2aSPeter Brune -  gscatter - scatters into the (ghosted) local vector of the local subdomain
32376857b2aSPeter Brune 
32476857b2aSPeter Brune    Level: intermediate
32576857b2aSPeter Brune 
326db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetSubdomains()`
32776857b2aSPeter Brune @*/
3289371c9d4SSatish Balay PetscErrorCode SNESNASMSetSubdomains(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[]) {
329111ade9eSPeter Brune   PetscErrorCode (*f)(SNES, PetscInt, SNES *, VecScatter *, VecScatter *, VecScatter *);
3306e111a19SKarl Rupp 
331eaedb033SPeter Brune   PetscFunctionBegin;
3329566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", &f));
3339566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter));
334eaedb033SPeter Brune   PetscFunctionReturn(0);
335eaedb033SPeter Brune }
336eaedb033SPeter Brune 
3379371c9d4SSatish Balay static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[]) {
338eaedb033SPeter Brune   PetscInt   i;
339eaedb033SPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
3406e111a19SKarl Rupp 
341eaedb033SPeter Brune   PetscFunctionBegin;
34228b400f6SJacob Faibussowitsch   PetscCheck(!snes->setupcalled, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
343eaedb033SPeter Brune 
344111ade9eSPeter Brune   /* tear down the previously set things */
3459566063dSJacob Faibussowitsch   PetscCall(SNESReset(snes));
346111ade9eSPeter Brune 
347eaedb033SPeter Brune   nasm->n = n;
348111ade9eSPeter Brune   if (oscatter) {
3499566063dSJacob Faibussowitsch     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)oscatter[i]));
350eaedb033SPeter Brune   }
351111ade9eSPeter Brune   if (iscatter) {
3529566063dSJacob Faibussowitsch     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)iscatter[i]));
353eaedb033SPeter Brune   }
354111ade9eSPeter Brune   if (gscatter) {
3559566063dSJacob Faibussowitsch     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)gscatter[i]));
356111ade9eSPeter Brune   }
357111ade9eSPeter Brune   if (oscatter) {
3589566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->oscatter));
3599566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->oscatter_copy));
360eaedb033SPeter Brune     for (i = 0; i < n; i++) {
361111ade9eSPeter Brune       nasm->oscatter[i] = oscatter[i];
3629566063dSJacob Faibussowitsch       PetscCall(VecScatterCopy(oscatter[i], &nasm->oscatter_copy[i]));
363eaedb033SPeter Brune     }
364111ade9eSPeter Brune   }
365111ade9eSPeter Brune   if (iscatter) {
3669566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->iscatter));
367ad540459SPierre Jolivet     for (i = 0; i < n; i++) nasm->iscatter[i] = iscatter[i];
368eaedb033SPeter Brune   }
369111ade9eSPeter Brune   if (gscatter) {
3709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->gscatter));
371ad540459SPierre Jolivet     for (i = 0; i < n; i++) nasm->gscatter[i] = gscatter[i];
372eaedb033SPeter Brune   }
373111ade9eSPeter Brune 
374eaedb033SPeter Brune   if (subsnes) {
3759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &nasm->subsnes));
376ad540459SPierre Jolivet     for (i = 0; i < n; i++) nasm->subsnes[i] = subsnes[i];
377eaedb033SPeter Brune   }
378eaedb033SPeter Brune   PetscFunctionReturn(0);
379eaedb033SPeter Brune }
380eaedb033SPeter Brune 
38176857b2aSPeter Brune /*@
382*f6dfbefdSBarry Smith    SNESNASMGetSubdomains - Get the local subdomain contexts for the nonlinear additive Schwarz solver
38376857b2aSPeter Brune 
384*f6dfbefdSBarry Smith    Not Collective but some of the objects returned will be parallel
38576857b2aSPeter Brune 
386f899ff85SJose E. Roman    Input Parameter:
387*f6dfbefdSBarry Smith .  snes - the `SNES` context
38876857b2aSPeter Brune 
38976857b2aSPeter Brune    Output Parameters:
39076857b2aSPeter Brune +  n - the number of local subdomains
39176857b2aSPeter Brune .  subsnes - solvers defined on the local subdomains
39276857b2aSPeter Brune .  iscatter - scatters into the nonoverlapping portions of the local subdomains
39376857b2aSPeter Brune .  oscatter - scatters into the overlapping portions of the local subdomains
39476857b2aSPeter Brune -  gscatter - scatters into the (ghosted) local vector of the local subdomain
39576857b2aSPeter Brune 
39676857b2aSPeter Brune    Level: intermediate
39776857b2aSPeter Brune 
398db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMSetSubdomains()`
39976857b2aSPeter Brune @*/
4009371c9d4SSatish Balay PetscErrorCode SNESNASMGetSubdomains(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[]) {
40176857b2aSPeter Brune   PetscErrorCode (*f)(SNES, PetscInt *, SNES **, VecScatter **, VecScatter **, VecScatter **);
40276857b2aSPeter Brune 
40376857b2aSPeter Brune   PetscFunctionBegin;
4049566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", &f));
4059566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter));
40676857b2aSPeter Brune   PetscFunctionReturn(0);
40776857b2aSPeter Brune }
40876857b2aSPeter Brune 
4099371c9d4SSatish Balay static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[]) {
41076857b2aSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
41176857b2aSPeter Brune 
41276857b2aSPeter Brune   PetscFunctionBegin;
41376857b2aSPeter Brune   if (n) *n = nasm->n;
41476857b2aSPeter Brune   if (oscatter) *oscatter = nasm->oscatter;
41576857b2aSPeter Brune   if (iscatter) *iscatter = nasm->iscatter;
41676857b2aSPeter Brune   if (gscatter) *gscatter = nasm->gscatter;
41783dc3634SPierre Jolivet   if (subsnes) *subsnes = nasm->subsnes;
41876857b2aSPeter Brune   PetscFunctionReturn(0);
41976857b2aSPeter Brune }
42076857b2aSPeter Brune 
42176857b2aSPeter Brune /*@
422*f6dfbefdSBarry Smith    SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors for the nonlinear additive Schwarz solver
42376857b2aSPeter Brune 
42476857b2aSPeter Brune    Not Collective
42576857b2aSPeter Brune 
426f899ff85SJose E. Roman    Input Parameter:
427*f6dfbefdSBarry Smith .  snes - the `SNES` context
42876857b2aSPeter Brune 
42976857b2aSPeter Brune    Output Parameters:
43076857b2aSPeter Brune +  n - the number of local subdomains
43176857b2aSPeter Brune .  x - The subdomain solution vector
43276857b2aSPeter Brune .  y - The subdomain step vector
43376857b2aSPeter Brune .  b - The subdomain RHS vector
43476857b2aSPeter Brune -  xl - The subdomain local vectors (ghosted)
43576857b2aSPeter Brune 
43676857b2aSPeter Brune    Level: developer
43776857b2aSPeter Brune 
438db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetSubdomains()`
43976857b2aSPeter Brune @*/
4409371c9d4SSatish Balay PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl) {
44176857b2aSPeter Brune   PetscErrorCode (*f)(SNES, PetscInt *, Vec **, Vec **, Vec **, Vec **);
44276857b2aSPeter Brune 
44376857b2aSPeter Brune   PetscFunctionBegin;
4449566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", &f));
4459566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, n, x, y, b, xl));
44676857b2aSPeter Brune   PetscFunctionReturn(0);
44776857b2aSPeter Brune }
44876857b2aSPeter Brune 
4499371c9d4SSatish Balay static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl) {
45076857b2aSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
45176857b2aSPeter Brune 
45276857b2aSPeter Brune   PetscFunctionBegin;
45376857b2aSPeter Brune   if (n) *n = nasm->n;
45476857b2aSPeter Brune   if (x) *x = nasm->x;
45576857b2aSPeter Brune   if (y) *y = nasm->y;
45676857b2aSPeter Brune   if (b) *b = nasm->b;
45776857b2aSPeter Brune   if (xl) *xl = nasm->xl;
45876857b2aSPeter Brune   PetscFunctionReturn(0);
45976857b2aSPeter Brune }
46076857b2aSPeter Brune 
461d728fb7dSPeter Brune /*@
462*f6dfbefdSBarry Smith    SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain Jacobians upon convergence for the
463*f6dfbefdSBarry Smith    nonlinear additive Schwarz solver
464d728fb7dSPeter Brune 
465*f6dfbefdSBarry Smith    Collective on snes
466d728fb7dSPeter Brune 
467d728fb7dSPeter Brune    Input Parameters:
468*f6dfbefdSBarry Smith +  snes - the SNES context
469*f6dfbefdSBarry Smith -  flg - `PETSC_TRUE` to compute the Jacobians
470d728fb7dSPeter Brune 
471d728fb7dSPeter Brune    Level: developer
472d728fb7dSPeter Brune 
47395452b02SPatrick Sanan    Notes:
474*f6dfbefdSBarry Smith    This is used almost exclusively in the implementation of `SNESASPIN`, where the converged subdomain and global Jacobian
475d728fb7dSPeter Brune    is needed at each linear iteration.
476d728fb7dSPeter Brune 
477db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetSubdomains()`
478d728fb7dSPeter Brune @*/
4799371c9d4SSatish Balay PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes, PetscBool flg) {
480d728fb7dSPeter Brune   PetscErrorCode (*f)(SNES, PetscBool);
481d728fb7dSPeter Brune 
482d728fb7dSPeter Brune   PetscFunctionBegin;
4839566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", &f));
4849566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, flg));
485d728fb7dSPeter Brune   PetscFunctionReturn(0);
486d728fb7dSPeter Brune }
487d728fb7dSPeter Brune 
4889371c9d4SSatish Balay static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes, PetscBool flg) {
489d728fb7dSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
490d728fb7dSPeter Brune 
491d728fb7dSPeter Brune   PetscFunctionBegin;
492d728fb7dSPeter Brune   nasm->finaljacobian = flg;
493d728fb7dSPeter Brune   PetscFunctionReturn(0);
494d728fb7dSPeter Brune }
49576857b2aSPeter Brune 
496610116beSPeter Brune /*@
497*f6dfbefdSBarry Smith    SNESNASMSetDamping - Sets the update damping for `SNESNASM` the nonlinear additive Schwarz solver
498610116beSPeter Brune 
499*f6dfbefdSBarry Smith    Logically collective on snes
500610116beSPeter Brune 
501610116beSPeter Brune    Input Parameters:
502*f6dfbefdSBarry Smith +  snes - the `SNES` context
503610116beSPeter Brune -  dmp - damping
504610116beSPeter Brune 
505610116beSPeter Brune    Level: intermediate
506610116beSPeter Brune 
50795452b02SPatrick Sanan    Notes:
50895452b02SPatrick Sanan     The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
5095dfa0f3bSBarry Smith 
510db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetDamping()`
511610116beSPeter Brune @*/
5129371c9d4SSatish Balay PetscErrorCode SNESNASMSetDamping(SNES snes, PetscReal dmp) {
513610116beSPeter Brune   PetscErrorCode (*f)(SNES, PetscReal);
514610116beSPeter Brune 
515610116beSPeter Brune   PetscFunctionBegin;
5169566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetDamping_C", (void (**)(void)) & f));
5179566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, dmp));
518610116beSPeter Brune   PetscFunctionReturn(0);
519610116beSPeter Brune }
520610116beSPeter Brune 
5219371c9d4SSatish Balay static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes, PetscReal dmp) {
522610116beSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
523610116beSPeter Brune 
524610116beSPeter Brune   PetscFunctionBegin;
525610116beSPeter Brune   nasm->damping = dmp;
526610116beSPeter Brune   PetscFunctionReturn(0);
527610116beSPeter Brune }
528610116beSPeter Brune 
529610116beSPeter Brune /*@
530*f6dfbefdSBarry Smith    SNESNASMGetDamping - Gets the update damping for `SNESNASM` the nonlinear additive Schwarz solver
531610116beSPeter Brune 
532610116beSPeter Brune    Not Collective
533610116beSPeter Brune 
534*f6dfbefdSBarry Smith    Input Parameter:
535*f6dfbefdSBarry Smith .  snes - the SNES context
536*f6dfbefdSBarry Smith 
537*f6dfbefdSBarry Smith    Output Parameter:
538*f6dfbefdSBarry Smith .  dmp - damping
539610116beSPeter Brune 
540610116beSPeter Brune    Level: intermediate
541610116beSPeter Brune 
542db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMSetDamping()`
543610116beSPeter Brune @*/
5449371c9d4SSatish Balay PetscErrorCode SNESNASMGetDamping(SNES snes, PetscReal *dmp) {
545610116beSPeter Brune   PetscFunctionBegin;
546cac4c232SBarry Smith   PetscUseMethod(snes, "SNESNASMGetDamping_C", (SNES, PetscReal *), (snes, dmp));
547610116beSPeter Brune   PetscFunctionReturn(0);
548610116beSPeter Brune }
549610116beSPeter Brune 
5509371c9d4SSatish Balay static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes, PetscReal *dmp) {
551610116beSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
552610116beSPeter Brune 
553610116beSPeter Brune   PetscFunctionBegin;
554610116beSPeter Brune   *dmp = nasm->damping;
555610116beSPeter Brune   PetscFunctionReturn(0);
556610116beSPeter Brune }
557610116beSPeter Brune 
55814eb1c5cSMatthew G. Knepley /*
55914eb1c5cSMatthew G. Knepley   Input Parameters:
56014eb1c5cSMatthew G. Knepley + snes - The solver
56114eb1c5cSMatthew G. Knepley . B - The RHS vector
56214eb1c5cSMatthew G. Knepley - X - The initial guess
56314eb1c5cSMatthew G. Knepley 
56414eb1c5cSMatthew G. Knepley   Output Parameters:
56514eb1c5cSMatthew G. Knepley . Y - The solution update
56614eb1c5cSMatthew G. Knepley 
56714eb1c5cSMatthew G. Knepley   TODO: All scatters should be packed into one
56814eb1c5cSMatthew G. Knepley */
5699371c9d4SSatish Balay PetscErrorCode SNESNASMSolveLocal_Private(SNES snes, Vec B, Vec Y, Vec X) {
570eaedb033SPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
571258e1594SPeter Brune   SNES       subsnes;
572eaedb033SPeter Brune   PetscInt   i;
573610116beSPeter Brune   PetscReal  dmp;
574f10b3e88SPatrick Farrell   Vec        Xl, Bl, Yl, Xlloc;
575f10b3e88SPatrick Farrell   VecScatter iscat, oscat, gscat, oscat_copy;
576111ade9eSPeter Brune   DM         dm, subdm;
577e0331734SPeter Brune   PCASMType  type;
5780adebc6cSBarry Smith 
579eaedb033SPeter Brune   PetscFunctionBegin;
5809566063dSJacob Faibussowitsch   PetscCall(SNESNASMGetType(snes, &type));
5819566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
5829566063dSJacob Faibussowitsch   PetscCall(VecSet(Y, 0));
5839566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
584eaedb033SPeter Brune   for (i = 0; i < nasm->n; i++) {
585f10b3e88SPatrick Farrell     /* scatter the solution to the global solution and the local solution */
586f10b3e88SPatrick Farrell     Xl         = nasm->x[i];
58770c78f05SPeter Brune     Xlloc      = nasm->xl[i];
58870c78f05SPeter Brune     oscat      = nasm->oscatter[i];
589f10b3e88SPatrick Farrell     oscat_copy = nasm->oscatter_copy[i];
590f10b3e88SPatrick Farrell     gscat      = nasm->gscatter[i];
5919566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD));
5929566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
59370c78f05SPeter Brune     if (B) {
59470c78f05SPeter Brune       /* scatter the RHS to the local RHS */
59570c78f05SPeter Brune       Bl = nasm->b[i];
5969566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD));
59770c78f05SPeter Brune     }
59870c78f05SPeter Brune   }
5999566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
600b20c023fSPeter Brune 
6019566063dSJacob Faibussowitsch   if (nasm->eventsubsolve) PetscCall(PetscLogEventBegin(nasm->eventsubsolve, snes, 0, 0, 0));
60270c78f05SPeter Brune   for (i = 0; i < nasm->n; i++) {
60370c78f05SPeter Brune     Xl      = nasm->x[i];
60470c78f05SPeter Brune     Xlloc   = nasm->xl[i];
60570c78f05SPeter Brune     Yl      = nasm->y[i];
606258e1594SPeter Brune     subsnes = nasm->subsnes[i];
6079566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(subsnes, &subdm));
608111ade9eSPeter Brune     iscat      = nasm->iscatter[i];
609111ade9eSPeter Brune     oscat      = nasm->oscatter[i];
610f10b3e88SPatrick Farrell     oscat_copy = nasm->oscatter_copy[i];
611111ade9eSPeter Brune     gscat      = nasm->gscatter[i];
6129566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD));
6139566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
61424b7f281SPeter Brune     if (B) {
61524b7f281SPeter Brune       Bl = nasm->b[i];
6169566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD));
617ed3c11a9SPeter Brune     } else Bl = NULL;
618f10b3e88SPatrick Farrell 
6199566063dSJacob Faibussowitsch     PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm));
6209566063dSJacob Faibussowitsch     PetscCall(VecCopy(Xl, Yl));
6219566063dSJacob Faibussowitsch     PetscCall(SNESSolve(subsnes, Bl, Xl));
6229566063dSJacob Faibussowitsch     PetscCall(VecAYPX(Yl, -1.0, Xl));
6239566063dSJacob Faibussowitsch     PetscCall(VecScale(Yl, nasm->damping));
624e0331734SPeter Brune     if (type == PC_ASM_BASIC) {
6259566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
6269566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
627e0331734SPeter Brune     } else if (type == PC_ASM_RESTRICT) {
6289566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
6299566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
630ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Only basic and restrict types are supported for SNESNASM");
631eaedb033SPeter Brune   }
6329566063dSJacob Faibussowitsch   if (nasm->eventsubsolve) PetscCall(PetscLogEventEnd(nasm->eventsubsolve, snes, 0, 0, 0));
6339566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
6341baa6e33SBarry Smith   if (nasm->weight_set) PetscCall(VecPointwiseMult(Y, Y, nasm->weight));
6359566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
6369566063dSJacob Faibussowitsch   PetscCall(SNESNASMGetDamping(snes, &dmp));
6379566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, dmp, Y));
638eaedb033SPeter Brune   PetscFunctionReturn(0);
639eaedb033SPeter Brune }
640eaedb033SPeter Brune 
6419371c9d4SSatish Balay static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) {
642602bec5dSPeter Brune   Vec        X    = Xfinal;
643d728fb7dSPeter Brune   SNES_NASM *nasm = (SNES_NASM *)snes->data;
644d728fb7dSPeter Brune   SNES       subsnes;
645ca44f815SPeter Brune   PetscInt   i, lag = 1;
646e59f0a30SPeter Brune   Vec        Xlloc, Xl, Fl, F;
647d728fb7dSPeter Brune   VecScatter oscat, gscat;
648d728fb7dSPeter Brune   DM         dm, subdm;
649d1e9a80fSBarry Smith 
650d728fb7dSPeter Brune   PetscFunctionBegin;
651602bec5dSPeter Brune   if (nasm->fjtype == 2) X = nasm->xinit;
652e59f0a30SPeter Brune   F = snes->vec_func;
6539566063dSJacob Faibussowitsch   if (snes->normschedule == SNES_NORM_NONE) PetscCall(SNESComputeFunction(snes, X, F));
6549566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
6559566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
6569566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
657602bec5dSPeter Brune   if (nasm->fjtype != 1) {
658d728fb7dSPeter Brune     for (i = 0; i < nasm->n; i++) {
659d728fb7dSPeter Brune       Xlloc = nasm->xl[i];
660d728fb7dSPeter Brune       gscat = nasm->gscatter[i];
6619566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
662602bec5dSPeter Brune     }
663d728fb7dSPeter Brune   }
6649566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
665d728fb7dSPeter Brune   for (i = 0; i < nasm->n; i++) {
666e59f0a30SPeter Brune     Fl      = nasm->subsnes[i]->vec_func;
667d728fb7dSPeter Brune     Xl      = nasm->x[i];
668d728fb7dSPeter Brune     Xlloc   = nasm->xl[i];
669d728fb7dSPeter Brune     subsnes = nasm->subsnes[i];
670d728fb7dSPeter Brune     oscat   = nasm->oscatter[i];
671d728fb7dSPeter Brune     gscat   = nasm->gscatter[i];
6729566063dSJacob Faibussowitsch     if (nasm->fjtype != 1) PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
6739566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(subsnes, &subdm));
6749566063dSJacob Faibussowitsch     PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm));
675602bec5dSPeter Brune     if (nasm->fjtype != 1) {
6769566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(subdm, Xlloc, INSERT_VALUES, Xl));
6779566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(subdm, Xlloc, INSERT_VALUES, Xl));
678602bec5dSPeter Brune     }
679ca44f815SPeter Brune     if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2;
680ca44f815SPeter Brune     else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian;
6819566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(subsnes, Xl, Fl));
6829566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(subsnes, Xl, subsnes->jacobian, subsnes->jacobian_pre));
683ca44f815SPeter Brune     if (lag > 1) subsnes->lagjacobian = lag;
684d728fb7dSPeter Brune   }
685d728fb7dSPeter Brune   PetscFunctionReturn(0);
686d728fb7dSPeter Brune }
687d728fb7dSPeter Brune 
6889371c9d4SSatish Balay static PetscErrorCode SNESSolve_NASM(SNES snes) {
689eaedb033SPeter Brune   Vec              F;
690eaedb033SPeter Brune   Vec              X;
691eaedb033SPeter Brune   Vec              B;
692111ade9eSPeter Brune   Vec              Y;
693eaedb033SPeter Brune   PetscInt         i;
694ed3c11a9SPeter Brune   PetscReal        fnorm = 0.0;
695365a6726SPeter Brune   SNESNormSchedule normschedule;
696d728fb7dSPeter Brune   SNES_NASM       *nasm = (SNES_NASM *)snes->data;
697eaedb033SPeter Brune 
698eaedb033SPeter Brune   PetscFunctionBegin;
699c579b300SPatrick Farrell 
7000b121fc5SBarry 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);
701c579b300SPatrick Farrell 
7029566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
703eaedb033SPeter Brune   X = snes->vec_sol;
704111ade9eSPeter Brune   Y = snes->vec_sol_update;
705eaedb033SPeter Brune   F = snes->vec_func;
706eaedb033SPeter Brune   B = snes->vec_rhs;
707eaedb033SPeter Brune 
7089566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
709eaedb033SPeter Brune   snes->iter = 0;
710eaedb033SPeter Brune   snes->norm = 0.;
7119566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
712eaedb033SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
7139566063dSJacob Faibussowitsch   PetscCall(SNESGetNormSchedule(snes, &normschedule));
714365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
715eaedb033SPeter Brune     /* compute the initial function and preconditioned update delX */
716eaedb033SPeter Brune     if (!snes->vec_func_init_set) {
7179566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
7181aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
719eaedb033SPeter Brune 
7209566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
721422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
7229566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
723eaedb033SPeter Brune     snes->iter = 0;
724eaedb033SPeter Brune     snes->norm = fnorm;
7259566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7269566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
7279566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, 0, snes->norm));
728eaedb033SPeter Brune 
729eaedb033SPeter Brune     /* test convergence */
730dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
731eaedb033SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
732eaedb033SPeter Brune   } else {
7339566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7349566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
7359566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, 0, snes->norm));
736eaedb033SPeter Brune   }
737eaedb033SPeter Brune 
738eaedb033SPeter Brune   /* Call general purpose update function */
739dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
740602bec5dSPeter Brune   /* copy the initial solution over for later */
7419566063dSJacob Faibussowitsch   if (nasm->fjtype == 2) PetscCall(VecCopy(X, nasm->xinit));
742eaedb033SPeter Brune 
743eaedb033SPeter Brune   for (i = 0; i < snes->max_its; i++) {
7449566063dSJacob Faibussowitsch     PetscCall(SNESNASMSolveLocal_Private(snes, B, Y, X));
745365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
7469566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
7479566063dSJacob Faibussowitsch       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
748422a814eSBarry Smith       SNESCheckFunctionNorm(snes, fnorm);
749eaedb033SPeter Brune     }
750eaedb033SPeter Brune     /* Monitor convergence */
7519566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
752eaedb033SPeter Brune     snes->iter = i + 1;
753eaedb033SPeter Brune     snes->norm = fnorm;
7549566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7559566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
7569566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
757eaedb033SPeter Brune     /* Test for convergence */
758dbbe0bcdSBarry Smith     if (normschedule == SNES_NORM_ALWAYS) PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
759d728fb7dSPeter Brune     if (snes->reason) break;
760eaedb033SPeter Brune     /* Call general purpose update function */
761dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
762eaedb033SPeter Brune   }
76307b62357SFande Kong   if (nasm->finaljacobian) {
7649566063dSJacob Faibussowitsch     PetscCall(SNESNASMComputeFinalJacobian_Private(snes, X));
76507b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
76607b62357SFande Kong   }
767365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS) {
768eaedb033SPeter Brune     if (i == snes->max_its) {
76963a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", snes->max_its));
770eaedb033SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
771eaedb033SPeter Brune     }
772*f6dfbefdSBarry Smith   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a nonlinear preconditioner */
773eaedb033SPeter Brune   PetscFunctionReturn(0);
774eaedb033SPeter Brune }
775eaedb033SPeter Brune 
776eaedb033SPeter Brune /*MC
777*f6dfbefdSBarry Smith   SNESNASM - Nonlinear Additive Schwarz solver
778eaedb033SPeter Brune 
779*f6dfbefdSBarry Smith    Options Database Keys:
78070c78603SPeter Brune +  -snes_nasm_log - enable logging events for the communication and solve stages
78170c78603SPeter Brune .  -snes_nasm_type <basic,restrict> - type of subdomain update used
7825dfa0f3bSBarry Smith .  -snes_asm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
78370c78603SPeter Brune .  -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate
784150d49b7SHong Zhang .  -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the jacobian is calculated at
78570c78603SPeter Brune .  -sub_snes_ - options prefix of the subdomain nonlinear solves
78670c78603SPeter Brune .  -sub_ksp_ - options prefix of the subdomain Krylov solver
78770c78603SPeter Brune -  -sub_pc_ - options prefix of the subdomain preconditioner
78870c78603SPeter Brune 
789eaedb033SPeter Brune    Level: advanced
790eaedb033SPeter Brune 
791*f6dfbefdSBarry Smith    Note:
792*f6dfbefdSBarry Smith    This is not often used directly as a solver, it converges too slowly. However it works well as a nonlinear preconditioner for
793*f6dfbefdSBarry Smith    the `SNESASPIN` solver
794*f6dfbefdSBarry Smith 
795*f6dfbefdSBarry Smith    Developer Note:
796*f6dfbefdSBarry Smith    This is a non-Newton based nonlinear solver that does not directly require a Jacobian; hence the flag snes->usesksp is set to
797*f6dfbefdSBarry Smith        false and `SNESView()` and -snes_view do not display a `KSP` object. However, if the flag nasm->finaljacobian is set (for example, if
798*f6dfbefdSBarry Smith        `SNESNASM` is used as a nonlinear preconditioner for  `SNESASPIN`) then `SNESSetUpMatrices()` is called to generate the
799*f6dfbefdSBarry Smith        Jacobian (needed by `SNESASPIN`)
800*f6dfbefdSBarry 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
801*f6dfbefdSBarry Smith        used by `SNESASPIN` they share the same Jacobian matrices because `SNESSetUp()` (called on the outer `SNESASPIN`) causes the inner `SNES`
802*f6dfbefdSBarry Smith        object (in this case `SNESNASM`) to inherit the outer Jacobian matrices.
803951fe5abSBarry Smith 
8044f02bc6aSBarry Smith    References:
805606c0280SSatish Balay .  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
8064f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
8074f02bc6aSBarry Smith 
808db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESNASMSetType()`, `SNESNASMGetType()`, `SNESNASMSetSubdomains()`, `SNESNASMGetSubdomains()`, `SNESNASMGetSubdomainVecs()`, `SNESNASMSetComputeFinalJacobian()`, `SNESNASMSetDamping()`, `SNESNASMGetDamping()`
809eaedb033SPeter Brune M*/
810eaedb033SPeter Brune 
8119371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) {
812eaedb033SPeter Brune   SNES_NASM *nasm;
813eaedb033SPeter Brune 
814eaedb033SPeter Brune   PetscFunctionBegin;
8159566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes, &nasm));
816eaedb033SPeter Brune   snes->data = (void *)nasm;
817eaedb033SPeter Brune 
818eaedb033SPeter Brune   nasm->n             = PETSC_DECIDE;
8199e5d0892SLisandro Dalcin   nasm->subsnes       = NULL;
8209e5d0892SLisandro Dalcin   nasm->x             = NULL;
8219e5d0892SLisandro Dalcin   nasm->xl            = NULL;
8229e5d0892SLisandro Dalcin   nasm->y             = NULL;
8239e5d0892SLisandro Dalcin   nasm->b             = NULL;
8249e5d0892SLisandro Dalcin   nasm->oscatter      = NULL;
8259e5d0892SLisandro Dalcin   nasm->oscatter_copy = NULL;
8269e5d0892SLisandro Dalcin   nasm->iscatter      = NULL;
8279e5d0892SLisandro Dalcin   nasm->gscatter      = NULL;
828610116beSPeter Brune   nasm->damping       = 1.;
829111ade9eSPeter Brune 
830111ade9eSPeter Brune   nasm->type          = PC_ASM_BASIC;
831d728fb7dSPeter Brune   nasm->finaljacobian = PETSC_FALSE;
832f10b3e88SPatrick Farrell   nasm->weight_set    = PETSC_FALSE;
833eaedb033SPeter Brune 
834eaedb033SPeter Brune   snes->ops->destroy        = SNESDestroy_NASM;
835eaedb033SPeter Brune   snes->ops->setup          = SNESSetUp_NASM;
836eaedb033SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
837eaedb033SPeter Brune   snes->ops->view           = SNESView_NASM;
838eaedb033SPeter Brune   snes->ops->solve          = SNESSolve_NASM;
839eaedb033SPeter Brune   snes->ops->reset          = SNESReset_NASM;
840eaedb033SPeter Brune 
841eaedb033SPeter Brune   snes->usesksp = PETSC_FALSE;
842efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
843eaedb033SPeter Brune 
8444fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
8454fc747eaSLawrence Mitchell 
846602bec5dSPeter Brune   nasm->fjtype              = 0;
847602bec5dSPeter Brune   nasm->xinit               = NULL;
8480298fd71SBarry Smith   nasm->eventrestrictinterp = 0;
8490298fd71SBarry Smith   nasm->eventsubsolve       = 0;
850b20c023fSPeter Brune 
851eaedb033SPeter Brune   if (!snes->tolerancesset) {
852eaedb033SPeter Brune     snes->max_its   = 10000;
853eaedb033SPeter Brune     snes->max_funcs = 10000;
854eaedb033SPeter Brune   }
855eaedb033SPeter Brune 
8569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", SNESNASMSetType_NASM));
8579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", SNESNASMGetType_NASM));
8589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", SNESNASMSetSubdomains_NASM));
8599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", SNESNASMGetSubdomains_NASM));
8609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", SNESNASMSetDamping_NASM));
8619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", SNESNASMGetDamping_NASM));
8629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", SNESNASMGetSubdomainVecs_NASM));
8639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", SNESNASMSetComputeFinalJacobian_NASM));
864eaedb033SPeter Brune   PetscFunctionReturn(0);
865eaedb033SPeter Brune }
86699e0435eSBarry Smith 
867448b6425SPatrick Farrell /*@
868448b6425SPatrick Farrell    SNESNASMGetSNES - Gets a subsolver
869448b6425SPatrick Farrell 
870448b6425SPatrick Farrell    Not collective
871448b6425SPatrick Farrell 
872448b6425SPatrick Farrell    Input Parameters:
873448b6425SPatrick Farrell +  snes - the SNES context
874448b6425SPatrick Farrell -  i - the number of the subsnes to get
875448b6425SPatrick Farrell 
876448b6425SPatrick Farrell    Output Parameters:
877448b6425SPatrick Farrell .  subsnes - the subsolver context
878448b6425SPatrick Farrell 
879448b6425SPatrick Farrell    Level: intermediate
880448b6425SPatrick Farrell 
881db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetNumber()`
882448b6425SPatrick Farrell @*/
8839371c9d4SSatish Balay PetscErrorCode SNESNASMGetSNES(SNES snes, PetscInt i, SNES *subsnes) {
884448b6425SPatrick Farrell   SNES_NASM *nasm = (SNES_NASM *)snes->data;
885448b6425SPatrick Farrell 
886448b6425SPatrick Farrell   PetscFunctionBegin;
8870b121fc5SBarry Smith   PetscCheck(i >= 0 && i < nasm->n, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "No such subsolver");
888448b6425SPatrick Farrell   *subsnes = nasm->subsnes[i];
889448b6425SPatrick Farrell   PetscFunctionReturn(0);
890448b6425SPatrick Farrell }
891448b6425SPatrick Farrell 
892448b6425SPatrick Farrell /*@
893448b6425SPatrick Farrell    SNESNASMGetNumber - Gets number of subsolvers
894448b6425SPatrick Farrell 
895448b6425SPatrick Farrell    Not collective
896448b6425SPatrick Farrell 
897448b6425SPatrick Farrell    Input Parameters:
898448b6425SPatrick Farrell .  snes - the SNES context
899448b6425SPatrick Farrell 
900448b6425SPatrick Farrell    Output Parameters:
901448b6425SPatrick Farrell .  n - the number of subsolvers
902448b6425SPatrick Farrell 
903448b6425SPatrick Farrell    Level: intermediate
904448b6425SPatrick Farrell 
905db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetSNES()`
906448b6425SPatrick Farrell @*/
9079371c9d4SSatish Balay PetscErrorCode SNESNASMGetNumber(SNES snes, PetscInt *n) {
908448b6425SPatrick Farrell   SNES_NASM *nasm = (SNES_NASM *)snes->data;
909448b6425SPatrick Farrell 
910448b6425SPatrick Farrell   PetscFunctionBegin;
911448b6425SPatrick Farrell   *n = nasm->n;
912448b6425SPatrick Farrell   PetscFunctionReturn(0);
913448b6425SPatrick Farrell }
914448b6425SPatrick Farrell 
915f10b3e88SPatrick Farrell /*@
916f10b3e88SPatrick Farrell    SNESNASMSetWeight - Sets weight to use when adding overlapping updates
917f10b3e88SPatrick Farrell 
918f10b3e88SPatrick Farrell    Collective
919f10b3e88SPatrick Farrell 
920f10b3e88SPatrick Farrell    Input Parameters:
921f10b3e88SPatrick Farrell +  snes - the SNES context
922f10b3e88SPatrick Farrell -  weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in)
923f10b3e88SPatrick Farrell 
924f10b3e88SPatrick Farrell    Level: intermediate
925f10b3e88SPatrick Farrell 
926db781477SPatrick Sanan .seealso: `SNESNASM`
927f10b3e88SPatrick Farrell @*/
9289371c9d4SSatish Balay PetscErrorCode SNESNASMSetWeight(SNES snes, Vec weight) {
929f10b3e88SPatrick Farrell   SNES_NASM *nasm = (SNES_NASM *)snes->data;
930f10b3e88SPatrick Farrell 
931f10b3e88SPatrick Farrell   PetscFunctionBegin;
932f10b3e88SPatrick Farrell 
9339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&nasm->weight));
934f10b3e88SPatrick Farrell   nasm->weight_set = PETSC_TRUE;
935f10b3e88SPatrick Farrell   nasm->weight     = weight;
9369566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)nasm->weight));
937f10b3e88SPatrick Farrell 
938f10b3e88SPatrick Farrell   PetscFunctionReturn(0);
939f10b3e88SPatrick Farrell }
940