xref: /petsc/src/snes/impls/nasm/nasm.c (revision 2ff79c18c26c94ed8cb599682f680f231dca6444)
1 #include <petsc/private/snesimpl.h> /*I   "petscsnes.h"   I*/
2 #include <petscdm.h>
3 
4 typedef struct {
5   PetscInt    n;             /* local subdomains */
6   SNES       *subsnes;       /* nonlinear solvers for each subdomain */
7   Vec        *x;             /* solution vectors */
8   Vec        *xl;            /* solution local vectors */
9   Vec        *y;             /* step vectors */
10   Vec        *b;             /* rhs vectors */
11   Vec         weight;        /* weighting for adding updates on overlaps, in global space */
12   VecScatter *oscatter;      /* scatter from global space to the subdomain global space */
13   VecScatter *oscatter_copy; /* copy of the above */
14   VecScatter *iscatter;      /* scatter from global space to the nonoverlapping subdomain space */
15   VecScatter *gscatter;      /* scatter from global space to the subdomain local space */
16   PCASMType   type;          /* ASM type */
17   PetscBool   usesdm;        /* use the DM for setting up the subproblems */
18   PetscBool   finaljacobian; /* compute the jacobian of the converged solution */
19   PetscReal   damping;       /* damping parameter for updates from the blocks */
20   PetscBool   weight_set;    /* use a weight in the overlap updates */
21 
22   /* logging events */
23   PetscLogEvent eventrestrictinterp;
24   PetscLogEvent eventsubsolve;
25 
26   PetscInt fjtype; /* type of computed jacobian */
27   Vec      xinit;  /* initial solution in case the final jacobian type is computed as first */
28 } SNES_NASM;
29 
30 const char *const SNESNASMTypes[]   = {"NONE", "RESTRICT", "INTERPOLATE", "BASIC", "PCASMType", "PC_ASM_", NULL};
31 const char *const SNESNASMFJTypes[] = {"FINALOUTER", "FINALINNER", "INITIAL"};
32 
33 static PetscErrorCode SNESReset_NASM(SNES snes)
34 {
35   SNES_NASM *nasm = (SNES_NASM *)snes->data;
36   PetscInt   i;
37 
38   PetscFunctionBegin;
39   for (i = 0; i < nasm->n; i++) {
40     if (nasm->xl) PetscCall(VecDestroy(&nasm->xl[i]));
41     if (nasm->x) PetscCall(VecDestroy(&nasm->x[i]));
42     if (nasm->y) PetscCall(VecDestroy(&nasm->y[i]));
43     if (nasm->b) PetscCall(VecDestroy(&nasm->b[i]));
44 
45     if (nasm->subsnes) PetscCall(SNESDestroy(&nasm->subsnes[i]));
46     if (nasm->oscatter) PetscCall(VecScatterDestroy(&nasm->oscatter[i]));
47     if (nasm->oscatter_copy) PetscCall(VecScatterDestroy(&nasm->oscatter_copy[i]));
48     if (nasm->iscatter) PetscCall(VecScatterDestroy(&nasm->iscatter[i]));
49     if (nasm->gscatter) PetscCall(VecScatterDestroy(&nasm->gscatter[i]));
50   }
51 
52   PetscCall(PetscFree(nasm->x));
53   PetscCall(PetscFree(nasm->xl));
54   PetscCall(PetscFree(nasm->y));
55   PetscCall(PetscFree(nasm->b));
56 
57   if (nasm->xinit) PetscCall(VecDestroy(&nasm->xinit));
58 
59   PetscCall(PetscFree(nasm->subsnes));
60   PetscCall(PetscFree(nasm->oscatter));
61   PetscCall(PetscFree(nasm->oscatter_copy));
62   PetscCall(PetscFree(nasm->iscatter));
63   PetscCall(PetscFree(nasm->gscatter));
64 
65   if (nasm->weight_set) PetscCall(VecDestroy(&nasm->weight));
66 
67   nasm->eventrestrictinterp = 0;
68   nasm->eventsubsolve       = 0;
69   PetscFunctionReturn(PETSC_SUCCESS);
70 }
71 
72 static PetscErrorCode SNESDestroy_NASM(SNES snes)
73 {
74   PetscFunctionBegin;
75   PetscCall(SNESReset_NASM(snes));
76   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", NULL));
77   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", NULL));
78   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", NULL));
79   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", NULL));
80   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", NULL));
81   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", NULL));
82   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", NULL));
83   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", NULL));
84   PetscCall(PetscFree(snes->data));
85   PetscFunctionReturn(PETSC_SUCCESS);
86 }
87 
88 static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm, Vec g, InsertMode mode, Vec l, void *ctx)
89 {
90   Vec bcs = (Vec)ctx;
91 
92   PetscFunctionBegin;
93   PetscCall(VecCopy(bcs, l));
94   PetscFunctionReturn(PETSC_SUCCESS);
95 }
96 
97 static PetscErrorCode SNESSetUp_NASM(SNES snes)
98 {
99   SNES_NASM  *nasm = (SNES_NASM *)snes->data;
100   DM          dm, subdm;
101   DM         *subdms;
102   PetscInt    i;
103   const char *optionsprefix;
104   Vec         F;
105 
106   PetscFunctionBegin;
107   if (!nasm->subsnes) {
108     PetscCall(SNESGetDM(snes, &dm));
109     if (dm) {
110       nasm->usesdm = PETSC_TRUE;
111       PetscCall(DMCreateDomainDecomposition(dm, &nasm->n, NULL, NULL, NULL, &subdms));
112       PetscCheck(subdms, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
113       PetscCall(DMCreateDomainDecompositionScatters(dm, nasm->n, subdms, &nasm->iscatter, &nasm->oscatter, &nasm->gscatter));
114       PetscCall(PetscMalloc1(nasm->n, &nasm->oscatter_copy));
115       for (i = 0; i < nasm->n; i++) PetscCall(VecScatterCopy(nasm->oscatter[i], &nasm->oscatter_copy[i]));
116 
117       PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
118       PetscCall(PetscMalloc1(nasm->n, &nasm->subsnes));
119       for (i = 0; i < nasm->n; i++) {
120         PetscCall(SNESCreate(PetscObjectComm((PetscObject)subdms[i]), &nasm->subsnes[i]));
121         PetscCall(PetscObjectIncrementTabLevel((PetscObject)nasm->subsnes[i], (PetscObject)snes, 1));
122         PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i], optionsprefix));
123         PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i], "sub_"));
124         PetscCall(SNESSetDM(nasm->subsnes[i], subdms[i]));
125         if (snes->ops->usercompute) {
126           PetscCall(SNESSetComputeApplicationContext(nasm->subsnes[i], snes->ops->usercompute, snes->ops->ctxdestroy));
127         } else {
128           void *ctx;
129 
130           PetscCall(SNESGetApplicationContext(snes, &ctx));
131           PetscCall(SNESSetApplicationContext(nasm->subsnes[i], ctx));
132         }
133         PetscCall(SNESSetFromOptions(nasm->subsnes[i]));
134         PetscCall(DMDestroy(&subdms[i]));
135       }
136       PetscCall(PetscFree(subdms));
137     } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Cannot construct local problems automatically without a DM!");
138   } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must set subproblems manually if there is no DM!");
139   /* allocate the global vectors */
140   if (!nasm->x) PetscCall(PetscCalloc1(nasm->n, &nasm->x));
141   if (!nasm->xl) PetscCall(PetscCalloc1(nasm->n, &nasm->xl));
142   if (!nasm->y) PetscCall(PetscCalloc1(nasm->n, &nasm->y));
143   if (!nasm->b) PetscCall(PetscCalloc1(nasm->n, &nasm->b));
144 
145   for (i = 0; i < nasm->n; i++) {
146     PetscCall(SNESGetFunction(nasm->subsnes[i], &F, NULL, NULL));
147     if (!nasm->x[i]) PetscCall(VecDuplicate(F, &nasm->x[i]));
148     if (!nasm->y[i]) PetscCall(VecDuplicate(F, &nasm->y[i]));
149     if (!nasm->b[i]) PetscCall(VecDuplicate(F, &nasm->b[i]));
150     if (!nasm->xl[i]) {
151       PetscCall(SNESGetDM(nasm->subsnes[i], &subdm));
152       PetscCall(DMCreateLocalVector(subdm, &nasm->xl[i]));
153       PetscCall(DMGlobalToLocalHookAdd(subdm, DMGlobalToLocalSubDomainDirichletHook_Private, NULL, nasm->xl[i]));
154     }
155   }
156   if (nasm->finaljacobian) {
157     PetscCall(SNESSetUpMatrices(snes));
158     if (nasm->fjtype == 2) PetscCall(VecDuplicate(snes->vec_sol, &nasm->xinit));
159     for (i = 0; i < nasm->n; i++) PetscCall(SNESSetUpMatrices(nasm->subsnes[i]));
160   }
161   PetscFunctionReturn(PETSC_SUCCESS);
162 }
163 
164 static PetscErrorCode SNESSetFromOptions_NASM(SNES snes, PetscOptionItems PetscOptionsObject)
165 {
166   PCASMType  asmtype;
167   PetscBool  flg, monflg;
168   SNES_NASM *nasm = (SNES_NASM *)snes->data;
169 
170   PetscFunctionBegin;
171   PetscOptionsHeadBegin(PetscOptionsObject, "Nonlinear Additive Schwarz options");
172   PetscCall(PetscOptionsEnum("-snes_nasm_type", "Type of restriction/extension", "", SNESNASMTypes, (PetscEnum)nasm->type, (PetscEnum *)&asmtype, &flg));
173   if (flg) PetscCall(SNESNASMSetType(snes, asmtype));
174   flg    = PETSC_FALSE;
175   monflg = PETSC_TRUE;
176   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));
177   if (flg) PetscCall(SNESNASMSetDamping(snes, nasm->damping));
178   PetscCall(PetscOptionsDeprecated("-snes_nasm_sub_view", NULL, "3.15", "Use -snes_view ::ascii_info_detail"));
179   PetscCall(PetscOptionsBool("-snes_nasm_finaljacobian", "Compute the global jacobian of the final iterate (for ASPIN)", "", nasm->finaljacobian, &nasm->finaljacobian, NULL));
180   PetscCall(PetscOptionsEList("-snes_nasm_finaljacobian_type", "The type of the final jacobian computed.", "", SNESNASMFJTypes, 3, SNESNASMFJTypes[0], &nasm->fjtype, NULL));
181   PetscCall(PetscOptionsBool("-snes_nasm_log", "Log times for subSNES solves and restriction", "", monflg, &monflg, &flg));
182   if (flg) {
183     PetscCall(PetscLogEventRegister("SNESNASMSubSolve", ((PetscObject)snes)->classid, &nasm->eventsubsolve));
184     PetscCall(PetscLogEventRegister("SNESNASMRestrict", ((PetscObject)snes)->classid, &nasm->eventrestrictinterp));
185   }
186   PetscOptionsHeadEnd();
187   PetscFunctionReturn(PETSC_SUCCESS);
188 }
189 
190 static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
191 {
192   SNES_NASM        *nasm = (SNES_NASM *)snes->data;
193   PetscMPIInt       rank, size;
194   PetscInt          i, N, bsz;
195   PetscBool         isascii, isstring;
196   PetscViewer       sviewer;
197   MPI_Comm          comm;
198   PetscViewerFormat format;
199   const char       *prefix;
200 
201   PetscFunctionBegin;
202   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
203   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
204   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
205   PetscCallMPI(MPI_Comm_rank(comm, &rank));
206   PetscCallMPI(MPI_Comm_size(comm, &size));
207   PetscCallMPI(MPIU_Allreduce(&nasm->n, &N, 1, MPIU_INT, MPI_SUM, comm));
208   if (isascii) {
209     PetscCall(PetscViewerASCIIPrintf(viewer, "  total subdomain blocks = %" PetscInt_FMT "\n", N));
210     PetscCall(PetscViewerGetFormat(viewer, &format));
211     if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
212       if (nasm->subsnes) {
213         PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for first block on rank 0:\n"));
214         PetscCall(SNESGetOptionsPrefix(snes, &prefix));
215         PetscCall(PetscViewerASCIIPrintf(viewer, "  Use -%ssnes_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : ""));
216         PetscCall(PetscViewerASCIIPushTab(viewer));
217         PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
218         if (rank == 0) {
219           PetscCall(PetscViewerASCIIPushTab(sviewer));
220           PetscCall(SNESView(nasm->subsnes[0], sviewer));
221           PetscCall(PetscViewerASCIIPopTab(sviewer));
222         }
223         PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
224         PetscCall(PetscViewerASCIIPopTab(viewer));
225       }
226     } else {
227       /* print the solver on each block */
228       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
229       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] number of local blocks = %" PetscInt_FMT "\n", rank, nasm->n));
230       PetscCall(PetscViewerFlush(viewer));
231       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
232       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for each block is in the following SNES objects:\n"));
233       PetscCall(PetscViewerASCIIPushTab(viewer));
234       PetscCall(PetscViewerASCIIPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n"));
235       PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
236       for (i = 0; i < nasm->n; i++) {
237         PetscCall(VecGetLocalSize(nasm->x[i], &bsz));
238         PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n", rank, i, bsz));
239         PetscCall(SNESView(nasm->subsnes[i], sviewer));
240         PetscCall(PetscViewerASCIIPrintf(sviewer, "- - - - - - - - - - - - - - - - - -\n"));
241       }
242       PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
243       PetscCall(PetscViewerASCIIPopTab(viewer));
244     }
245   } else if (isstring) {
246     PetscCall(PetscViewerStringSPrintf(viewer, " blocks=%" PetscInt_FMT ",type=%s", N, SNESNASMTypes[nasm->type]));
247     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
248     if (nasm->subsnes && rank == 0) PetscCall(SNESView(nasm->subsnes[0], sviewer));
249     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
250   }
251   PetscFunctionReturn(PETSC_SUCCESS);
252 }
253 
254 /*@
255   SNESNASMSetType - Set the type of subdomain update used for the nonlinear additive Schwarz solver `SNESNASM`
256 
257   Logically Collective
258 
259   Input Parameters:
260 + snes - the `SNES` context
261 - type - the type of update, `PC_ASM_BASIC` or `PC_ASM_RESTRICT`
262 
263   Options Database Key:
264 . -snes_nasm_type <basic,restrict> - type of subdomain update used
265 
266   Level: intermediate
267 
268 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetType()`, `PCASMSetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType`
269 @*/
270 PetscErrorCode SNESNASMSetType(SNES snes, PCASMType type)
271 {
272   PetscFunctionBegin;
273   PetscTryMethod(snes, "SNESNASMSetType_C", (SNES, PCASMType), (snes, type));
274   PetscFunctionReturn(PETSC_SUCCESS);
275 }
276 
277 static PetscErrorCode SNESNASMSetType_NASM(SNES snes, PCASMType type)
278 {
279   SNES_NASM *nasm = (SNES_NASM *)snes->data;
280 
281   PetscFunctionBegin;
282   PetscCheck(type == PC_ASM_BASIC || type == PC_ASM_RESTRICT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESNASM only supports basic and restrict types");
283   nasm->type = type;
284   PetscFunctionReturn(PETSC_SUCCESS);
285 }
286 
287 /*@
288   SNESNASMGetType - Get the type of subdomain update used for the nonlinear additive Schwarz solver `SNESNASM`
289 
290   Logically Collective
291 
292   Input Parameter:
293 . snes - the `SNES` context
294 
295   Output Parameter:
296 . type - the type of update
297 
298   Level: intermediate
299 
300 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetType()`, `PCASMGetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType`
301 @*/
302 PetscErrorCode SNESNASMGetType(SNES snes, PCASMType *type)
303 {
304   PetscFunctionBegin;
305   PetscUseMethod(snes, "SNESNASMGetType_C", (SNES, PCASMType *), (snes, type));
306   PetscFunctionReturn(PETSC_SUCCESS);
307 }
308 
309 static PetscErrorCode SNESNASMGetType_NASM(SNES snes, PCASMType *type)
310 {
311   SNES_NASM *nasm = (SNES_NASM *)snes->data;
312 
313   PetscFunctionBegin;
314   *type = nasm->type;
315   PetscFunctionReturn(PETSC_SUCCESS);
316 }
317 
318 /*@
319   SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems in the nonlinear additive Schwarz solver
320 
321   Logically Collective
322 
323   Input Parameters:
324 + snes     - the `SNES` context
325 . n        - the number of local subdomains
326 . subsnes  - solvers defined on the local subdomains
327 . iscatter - scatters into the nonoverlapping portions of the local subdomains
328 . oscatter - scatters into the overlapping portions of the local subdomains
329 - gscatter - scatters into the (ghosted) local vector of the local subdomain
330 
331   Level: intermediate
332 
333 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()`
334 @*/
335 PetscErrorCode SNESNASMSetSubdomains(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[])
336 {
337   PetscErrorCode (*f)(SNES, PetscInt, SNES *, VecScatter *, VecScatter *, VecScatter *);
338 
339   PetscFunctionBegin;
340   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", &f));
341   if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter));
342   PetscFunctionReturn(PETSC_SUCCESS);
343 }
344 
345 static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[])
346 {
347   PetscInt   i;
348   SNES_NASM *nasm = (SNES_NASM *)snes->data;
349 
350   PetscFunctionBegin;
351   PetscCheck(!snes->setupcalled, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
352 
353   /* tear down the previously set things */
354   PetscCall(SNESReset(snes));
355 
356   nasm->n = n;
357   if (oscatter) {
358     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)oscatter[i]));
359   }
360   if (iscatter) {
361     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)iscatter[i]));
362   }
363   if (gscatter) {
364     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)gscatter[i]));
365   }
366   if (oscatter) {
367     PetscCall(PetscMalloc1(n, &nasm->oscatter));
368     PetscCall(PetscMalloc1(n, &nasm->oscatter_copy));
369     for (i = 0; i < n; i++) {
370       nasm->oscatter[i] = oscatter[i];
371       PetscCall(VecScatterCopy(oscatter[i], &nasm->oscatter_copy[i]));
372     }
373   }
374   if (iscatter) {
375     PetscCall(PetscMalloc1(n, &nasm->iscatter));
376     for (i = 0; i < n; i++) nasm->iscatter[i] = iscatter[i];
377   }
378   if (gscatter) {
379     PetscCall(PetscMalloc1(n, &nasm->gscatter));
380     for (i = 0; i < n; i++) nasm->gscatter[i] = gscatter[i];
381   }
382 
383   if (subsnes) {
384     PetscCall(PetscMalloc1(n, &nasm->subsnes));
385     for (i = 0; i < n; i++) nasm->subsnes[i] = subsnes[i];
386   }
387   PetscFunctionReturn(PETSC_SUCCESS);
388 }
389 
390 /*@
391   SNESNASMGetSubdomains - Get the local subdomain contexts for the nonlinear additive Schwarz solver
392 
393   Not Collective but some of the objects returned will be parallel
394 
395   Input Parameter:
396 . snes - the `SNES` context
397 
398   Output Parameters:
399 + n        - the number of local subdomains
400 . subsnes  - solvers defined on the local subdomains
401 . iscatter - scatters into the nonoverlapping portions of the local subdomains
402 . oscatter - scatters into the overlapping portions of the local subdomains
403 - gscatter - scatters into the (ghosted) local vector of the local subdomain
404 
405   Level: intermediate
406 
407 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetSubdomains()`
408 @*/
409 PetscErrorCode SNESNASMGetSubdomains(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[])
410 {
411   PetscErrorCode (*f)(SNES, PetscInt *, SNES **, VecScatter **, VecScatter **, VecScatter **);
412 
413   PetscFunctionBegin;
414   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", &f));
415   if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter));
416   PetscFunctionReturn(PETSC_SUCCESS);
417 }
418 
419 static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[])
420 {
421   SNES_NASM *nasm = (SNES_NASM *)snes->data;
422 
423   PetscFunctionBegin;
424   if (n) *n = nasm->n;
425   if (oscatter) *oscatter = nasm->oscatter;
426   if (iscatter) *iscatter = nasm->iscatter;
427   if (gscatter) *gscatter = nasm->gscatter;
428   if (subsnes) *subsnes = nasm->subsnes;
429   PetscFunctionReturn(PETSC_SUCCESS);
430 }
431 
432 /*@
433   SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors for the nonlinear additive Schwarz solver
434 
435   Not Collective
436 
437   Input Parameter:
438 . snes - the `SNES` context
439 
440   Output Parameters:
441 + n  - the number of local subdomains
442 . x  - The subdomain solution vector
443 . y  - The subdomain step vector
444 . b  - The subdomain RHS vector
445 - xl - The subdomain local vectors (ghosted)
446 
447   Level: developer
448 
449 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()`
450 @*/
451 PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes, PetscInt *n, Vec *x[], Vec *y[], Vec *b[], Vec *xl[])
452 {
453   PetscErrorCode (*f)(SNES, PetscInt *, Vec **, Vec **, Vec **, Vec **);
454 
455   PetscFunctionBegin;
456   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", &f));
457   if (f) PetscCall((f)(snes, n, x, y, b, xl));
458   PetscFunctionReturn(PETSC_SUCCESS);
459 }
460 
461 static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl)
462 {
463   SNES_NASM *nasm = (SNES_NASM *)snes->data;
464 
465   PetscFunctionBegin;
466   if (n) *n = nasm->n;
467   if (x) *x = nasm->x;
468   if (y) *y = nasm->y;
469   if (b) *b = nasm->b;
470   if (xl) *xl = nasm->xl;
471   PetscFunctionReturn(PETSC_SUCCESS);
472 }
473 
474 /*@
475   SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain Jacobians upon convergence for the
476   nonlinear additive Schwarz solver
477 
478   Collective
479 
480   Input Parameters:
481 + snes - the SNES context
482 - flg  - `PETSC_TRUE` to compute the Jacobians
483 
484   Level: developer
485 
486   Notes:
487   This is used almost exclusively in the implementation of `SNESASPIN`, where the converged subdomain and global Jacobian
488   is needed at each linear iteration.
489 
490 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()`
491 @*/
492 PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes, PetscBool flg)
493 {
494   PetscErrorCode (*f)(SNES, PetscBool);
495 
496   PetscFunctionBegin;
497   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", &f));
498   if (f) PetscCall((f)(snes, flg));
499   PetscFunctionReturn(PETSC_SUCCESS);
500 }
501 
502 static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes, PetscBool flg)
503 {
504   SNES_NASM *nasm = (SNES_NASM *)snes->data;
505 
506   PetscFunctionBegin;
507   nasm->finaljacobian = flg;
508   PetscFunctionReturn(PETSC_SUCCESS);
509 }
510 
511 /*@
512   SNESNASMSetDamping - Sets the update damping for `SNESNASM` the nonlinear additive Schwarz solver
513 
514   Logically Collective
515 
516   Input Parameters:
517 + snes - the `SNES` context
518 - dmp  - damping
519 
520   Options Database Key:
521 . -snes_nasm_damping <dmp> - the new solution is obtained as old solution plus `dmp` times (sum of the solutions on the subdomains)
522 
523   Level: intermediate
524 
525   Note:
526   The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
527 
528 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetDamping()`
529 @*/
530 PetscErrorCode SNESNASMSetDamping(SNES snes, PetscReal dmp)
531 {
532   PetscFunctionBegin;
533   PetscTryMethod(snes, "SNESNASMSetDamping_C", (SNES, PetscReal), (snes, dmp));
534   PetscFunctionReturn(PETSC_SUCCESS);
535 }
536 
537 static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes, PetscReal dmp)
538 {
539   SNES_NASM *nasm = (SNES_NASM *)snes->data;
540 
541   PetscFunctionBegin;
542   nasm->damping = dmp;
543   PetscFunctionReturn(PETSC_SUCCESS);
544 }
545 
546 /*@
547   SNESNASMGetDamping - Gets the update damping for `SNESNASM` the nonlinear additive Schwarz solver
548 
549   Not Collective
550 
551   Input Parameter:
552 . snes - the `SNES` context
553 
554   Output Parameter:
555 . dmp - damping
556 
557   Level: intermediate
558 
559 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetDamping()`
560 @*/
561 PetscErrorCode SNESNASMGetDamping(SNES snes, PetscReal *dmp)
562 {
563   PetscFunctionBegin;
564   PetscUseMethod(snes, "SNESNASMGetDamping_C", (SNES, PetscReal *), (snes, dmp));
565   PetscFunctionReturn(PETSC_SUCCESS);
566 }
567 
568 static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes, PetscReal *dmp)
569 {
570   SNES_NASM *nasm = (SNES_NASM *)snes->data;
571 
572   PetscFunctionBegin;
573   *dmp = nasm->damping;
574   PetscFunctionReturn(PETSC_SUCCESS);
575 }
576 
577 /*
578   Input Parameters:
579 + snes - The solver
580 . B - The RHS vector
581 - X - The initial guess
582 
583   Output Parameter:
584 . Y - The solution update
585 
586   TODO: All scatters should be packed into one
587 */
588 static PetscErrorCode SNESNASMSolveLocal_Private(SNES snes, Vec B, Vec Y, Vec X)
589 {
590   SNES_NASM *nasm = (SNES_NASM *)snes->data;
591   SNES       subsnes;
592   PetscInt   i;
593   PetscReal  dmp;
594   Vec        Xl, Bl, Yl, Xlloc;
595   VecScatter iscat, oscat, gscat, oscat_copy;
596   DM         dm, subdm;
597   PCASMType  type;
598 
599   PetscFunctionBegin;
600   PetscCall(SNESNASMGetType(snes, &type));
601   PetscCall(SNESGetDM(snes, &dm));
602   PetscCall(VecSet(Y, 0));
603   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
604   for (i = 0; i < nasm->n; i++) { /* scatter the global solution to the overlap solution and the local solution */
605     Xl         = nasm->x[i];
606     Xlloc      = nasm->xl[i];
607     oscat      = nasm->oscatter[i];
608     oscat_copy = nasm->oscatter_copy[i];
609     gscat      = nasm->gscatter[i];
610     PetscCall(VecScatterBegin(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD));
611     PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
612 
613     if (B) {
614       /* scatter the RHS to the local RHS */
615       Bl = nasm->b[i];
616       PetscCall(VecScatterBegin(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD));
617     }
618   }
619   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
620 
621   if (nasm->eventsubsolve) PetscCall(PetscLogEventBegin(nasm->eventsubsolve, snes, 0, 0, 0));
622   for (i = 0; i < nasm->n; i++) {
623     PetscErrorCode (*bl)(DM, Vec, void *);
624     void *bctx;
625 
626     Xl      = nasm->x[i];
627     Xlloc   = nasm->xl[i];
628     Yl      = nasm->y[i];
629     subsnes = nasm->subsnes[i];
630     PetscCall(SNESGetDM(subsnes, &subdm));
631     iscat      = nasm->iscatter[i];
632     oscat      = nasm->oscatter[i];
633     oscat_copy = nasm->oscatter_copy[i];
634     gscat      = nasm->gscatter[i];
635     PetscCall(VecScatterEnd(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD));
636     PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
637     if (B) {
638       Bl = nasm->b[i];
639       PetscCall(VecScatterEnd(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD));
640     } else Bl = NULL;
641 
642     PetscCall(SNESGetDM(subsnes, &subdm));
643     PetscCall(DMSNESGetBoundaryLocal(subdm, &bl, &bctx));
644     if (bl) PetscCall((*bl)(subdm, Xlloc, bctx));
645 
646     PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm));
647     PetscCall(VecCopy(Xl, Yl));
648     PetscCall(SNESSolve(subsnes, Bl, Xl));
649     PetscCall(VecAYPX(Yl, -1.0, Xl));
650     PetscCall(VecScale(Yl, nasm->damping));
651     if (type == PC_ASM_BASIC) {
652       PetscCall(VecScatterBegin(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
653       PetscCall(VecScatterEnd(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
654     } else if (type == PC_ASM_RESTRICT) {
655       PetscCall(VecScatterBegin(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
656       PetscCall(VecScatterEnd(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
657     } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Only basic and restrict types are supported for SNESNASM");
658   }
659   if (nasm->eventsubsolve) PetscCall(PetscLogEventEnd(nasm->eventsubsolve, snes, 0, 0, 0));
660   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
661   if (nasm->weight_set) PetscCall(VecPointwiseMult(Y, Y, nasm->weight));
662   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
663   PetscCall(SNESNASMGetDamping(snes, &dmp));
664   PetscCall(VecAXPY(X, dmp, Y));
665   PetscFunctionReturn(PETSC_SUCCESS);
666 }
667 
668 static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal)
669 {
670   Vec        X    = Xfinal;
671   SNES_NASM *nasm = (SNES_NASM *)snes->data;
672   SNES       subsnes;
673   PetscInt   i, lag = 1;
674   Vec        Xlloc, Xl, Fl, F;
675   VecScatter oscat, gscat;
676   DM         dm, subdm;
677 
678   PetscFunctionBegin;
679   if (nasm->fjtype == 2) X = nasm->xinit;
680   F = snes->vec_func;
681   if (snes->normschedule == SNES_NORM_NONE) PetscCall(SNESComputeFunction(snes, X, F));
682   PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
683   PetscCall(SNESGetDM(snes, &dm));
684   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
685   if (nasm->fjtype != 1) {
686     for (i = 0; i < nasm->n; i++) {
687       Xlloc = nasm->xl[i];
688       gscat = nasm->gscatter[i];
689       PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
690     }
691   }
692   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
693   for (i = 0; i < nasm->n; i++) {
694     Fl      = nasm->subsnes[i]->vec_func;
695     Xl      = nasm->x[i];
696     Xlloc   = nasm->xl[i];
697     subsnes = nasm->subsnes[i];
698     oscat   = nasm->oscatter[i];
699     gscat   = nasm->gscatter[i];
700     if (nasm->fjtype != 1) PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
701     PetscCall(SNESGetDM(subsnes, &subdm));
702     PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm));
703     if (nasm->fjtype != 1) {
704       PetscCall(DMLocalToGlobalBegin(subdm, Xlloc, INSERT_VALUES, Xl));
705       PetscCall(DMLocalToGlobalEnd(subdm, Xlloc, INSERT_VALUES, Xl));
706     }
707     if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2;
708     else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian;
709     PetscCall(SNESComputeFunction(subsnes, Xl, Fl));
710     PetscCall(SNESComputeJacobian(subsnes, Xl, subsnes->jacobian, subsnes->jacobian_pre));
711     if (lag > 1) subsnes->lagjacobian = lag;
712   }
713   PetscFunctionReturn(PETSC_SUCCESS);
714 }
715 
716 static PetscErrorCode SNESSolve_NASM(SNES snes)
717 {
718   Vec              F;
719   Vec              X;
720   Vec              B;
721   Vec              Y;
722   PetscInt         i;
723   PetscReal        fnorm = 0.0;
724   SNESNormSchedule normschedule;
725   SNES_NASM       *nasm = (SNES_NASM *)snes->data;
726 
727   PetscFunctionBegin;
728   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);
729 
730   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
731   X = snes->vec_sol;
732   Y = snes->vec_sol_update;
733   F = snes->vec_func;
734   B = snes->vec_rhs;
735 
736   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
737   snes->iter = 0;
738   snes->norm = 0.;
739   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
740   snes->reason = SNES_CONVERGED_ITERATING;
741   PetscCall(SNESGetNormSchedule(snes, &normschedule));
742   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY || !snes->max_its) {
743     /* compute the initial function and preconditioned update delX */
744     if (!snes->vec_func_init_set) PetscCall(SNESComputeFunction(snes, X, F));
745     else snes->vec_func_init_set = PETSC_FALSE;
746 
747     PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
748     SNESCheckFunctionNorm(snes, fnorm);
749     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
750     snes->iter = 0;
751     snes->norm = fnorm;
752     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
753     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
754 
755     /* test convergence */
756     PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
757     PetscCall(SNESMonitor(snes, 0, snes->norm));
758     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
759   } else {
760     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
761     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
762     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
763   }
764 
765   /* Call general purpose update function */
766   PetscTryTypeMethod(snes, update, snes->iter);
767   /* copy the initial solution over for later */
768   if (nasm->fjtype == 2) PetscCall(VecCopy(X, nasm->xinit));
769 
770   for (i = 0; i < snes->max_its; i++) {
771     PetscCall(SNESNASMSolveLocal_Private(snes, B, Y, X));
772     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
773       PetscCall(SNESComputeFunction(snes, X, F));
774       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
775       SNESCheckFunctionNorm(snes, fnorm);
776     }
777     /* Monitor convergence */
778     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
779     snes->iter = i + 1;
780     snes->norm = fnorm;
781     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
782     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
783     /* Test for convergence */
784     PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm));
785     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
786     if (snes->reason) break;
787     /* Call general purpose update function */
788     PetscTryTypeMethod(snes, update, snes->iter);
789   }
790   if (nasm->finaljacobian) {
791     PetscCall(SNESNASMComputeFinalJacobian_Private(snes, X));
792     SNESCheckJacobianDomainerror(snes);
793   }
794   PetscFunctionReturn(PETSC_SUCCESS);
795 }
796 
797 /*MC
798   SNESNASM - Nonlinear Additive Schwarz solver {cite}`ck02`, {cite}`bruneknepleysmithtu15`
799 
800    Options Database Keys:
801 +  -snes_nasm_log                                                - enable logging events for the communication and solve stages
802 .  -snes_nasm_type <basic,restrict>                              - type of subdomain update used
803 .  -snes_nasm_damping <dmp>                                      - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
804 .  -snes_nasm_finaljacobian                                      - compute the local and global Jacobians of the final iterate
805 .  -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the Jacobian is calculated at
806 .  -sub_snes_                                                    - options prefix of the subdomain nonlinear solves
807 .  -sub_ksp_                                                     - options prefix of the subdomain Krylov solver
808 -  -sub_pc_                                                      - options prefix of the subdomain preconditioner
809 
810    Level: advanced
811 
812    Note:
813    This is not often used directly as a solver, it converges too slowly. However it works well as a nonlinear preconditioner for
814    the `SNESASPIN` solver
815 
816    Developer Note:
817    This is a non-Newton based nonlinear solver that does not directly require a Jacobian; hence the flag snes->usesksp is set to
818    false and `SNESView()` and -snes_view do not display a `KSP` object. However, if the flag nasm->finaljacobian is set (for example, if
819    `SNESNASM` is used as a nonlinear preconditioner for  `SNESASPIN`) then `SNESSetUpMatrices()` is called to generate the
820    Jacobian (needed by `SNESASPIN`)
821    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
822    used by `SNESASPIN` they share the same Jacobian matrices because `SNESSetUp()` (called on the outer `SNESASPIN`) causes the inner `SNES`
823    object (in this case `SNESNASM`) to inherit the outer Jacobian matrices.
824 
825 .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESNASMSetType()`, `SNESNASMGetType()`, `SNESNASMSetSubdomains()`, `SNESNASMGetSubdomains()`,
826           `SNESNASMGetSubdomainVecs()`, `SNESNASMSetComputeFinalJacobian()`, `SNESNASMSetDamping()`, `SNESNASMGetDamping()`, `SNESNASMSetWeight()`,
827           `SNESNASMGetSNES()`, `SNESNASMGetNumber()`
828 M*/
829 
830 PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes)
831 {
832   SNES_NASM *nasm;
833 
834   PetscFunctionBegin;
835   PetscCall(PetscNew(&nasm));
836   snes->data = (void *)nasm;
837 
838   nasm->n             = PETSC_DECIDE;
839   nasm->subsnes       = NULL;
840   nasm->x             = NULL;
841   nasm->xl            = NULL;
842   nasm->y             = NULL;
843   nasm->b             = NULL;
844   nasm->oscatter      = NULL;
845   nasm->oscatter_copy = NULL;
846   nasm->iscatter      = NULL;
847   nasm->gscatter      = NULL;
848   nasm->damping       = 1.;
849 
850   nasm->type          = PC_ASM_BASIC;
851   nasm->finaljacobian = PETSC_FALSE;
852   nasm->weight_set    = PETSC_FALSE;
853 
854   snes->ops->destroy        = SNESDestroy_NASM;
855   snes->ops->setup          = SNESSetUp_NASM;
856   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
857   snes->ops->view           = SNESView_NASM;
858   snes->ops->solve          = SNESSolve_NASM;
859   snes->ops->reset          = SNESReset_NASM;
860 
861   snes->usesksp = PETSC_FALSE;
862   snes->usesnpc = PETSC_FALSE;
863 
864   snes->alwayscomputesfinalresidual = PETSC_FALSE;
865 
866   nasm->fjtype              = 0;
867   nasm->xinit               = NULL;
868   nasm->eventrestrictinterp = 0;
869   nasm->eventsubsolve       = 0;
870 
871   PetscCall(SNESParametersInitialize(snes));
872   PetscObjectParameterSetDefault(snes, max_funcs, 10000);
873   PetscObjectParameterSetDefault(snes, max_its, 10000);
874 
875   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", SNESNASMSetType_NASM));
876   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", SNESNASMGetType_NASM));
877   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", SNESNASMSetSubdomains_NASM));
878   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", SNESNASMGetSubdomains_NASM));
879   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", SNESNASMSetDamping_NASM));
880   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", SNESNASMGetDamping_NASM));
881   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", SNESNASMGetSubdomainVecs_NASM));
882   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", SNESNASMSetComputeFinalJacobian_NASM));
883   PetscFunctionReturn(PETSC_SUCCESS);
884 }
885 
886 /*@
887   SNESNASMGetSNES - Gets a subsolver
888 
889   Not Collective
890 
891   Input Parameters:
892 + snes - the `SNES` context
893 - i    - the number of the subsnes to get
894 
895   Output Parameter:
896 . subsnes - the subsolver context
897 
898   Level: intermediate
899 
900 .seealso: [](ch_snes), `SNESNASM`, `SNESNASMGetNumber()`
901 @*/
902 PetscErrorCode SNESNASMGetSNES(SNES snes, PetscInt i, SNES *subsnes)
903 {
904   SNES_NASM *nasm = (SNES_NASM *)snes->data;
905 
906   PetscFunctionBegin;
907   PetscCheck(i >= 0 && i < nasm->n, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "No such subsolver");
908   *subsnes = nasm->subsnes[i];
909   PetscFunctionReturn(PETSC_SUCCESS);
910 }
911 
912 /*@
913   SNESNASMGetNumber - Gets number of subsolvers
914 
915   Not Collective
916 
917   Input Parameter:
918 . snes - the `SNES` context
919 
920   Output Parameter:
921 . n - the number of subsolvers
922 
923   Level: intermediate
924 
925 .seealso: [](ch_snes), `SNESNASM`, `SNESNASMGetSNES()`
926 @*/
927 PetscErrorCode SNESNASMGetNumber(SNES snes, PetscInt *n)
928 {
929   SNES_NASM *nasm = (SNES_NASM *)snes->data;
930 
931   PetscFunctionBegin;
932   *n = nasm->n;
933   PetscFunctionReturn(PETSC_SUCCESS);
934 }
935 
936 /*@
937   SNESNASMSetWeight - Sets weight to use when adding overlapping updates
938 
939   Collective
940 
941   Input Parameters:
942 + snes   - the `SNES` context
943 - weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in)
944 
945   Level: intermediate
946 
947 .seealso: [](ch_snes), `SNESNASM`
948 @*/
949 PetscErrorCode SNESNASMSetWeight(SNES snes, Vec weight)
950 {
951   SNES_NASM *nasm = (SNES_NASM *)snes->data;
952 
953   PetscFunctionBegin;
954   PetscCall(VecDestroy(&nasm->weight));
955   nasm->weight_set = PETSC_TRUE;
956   nasm->weight     = weight;
957   PetscCall(PetscObjectReference((PetscObject)nasm->weight));
958   PetscFunctionReturn(PETSC_SUCCESS);
959 }
960