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