xref: /petsc/src/snes/impls/nasm/nasm.c (revision 1bd63e3efbce1229b47e2ec5bb1e7d27e9e9b2dd)
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 `SNESNASM`
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   Options Database Key:
266 . -snes_nasm_type <basic,restrict> - type of subdomain update used
267 
268   Level: intermediate
269 
270 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetType()`, `PCASMSetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType`
271 @*/
272 PetscErrorCode SNESNASMSetType(SNES snes, PCASMType type)
273 {
274   PetscErrorCode (*f)(SNES, PCASMType);
275 
276   PetscFunctionBegin;
277   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetType_C", &f));
278   if (f) PetscCall((f)(snes, type));
279   PetscFunctionReturn(PETSC_SUCCESS);
280 }
281 
282 static PetscErrorCode SNESNASMSetType_NASM(SNES snes, PCASMType type)
283 {
284   SNES_NASM *nasm = (SNES_NASM *)snes->data;
285 
286   PetscFunctionBegin;
287   PetscCheck(type == PC_ASM_BASIC || type == PC_ASM_RESTRICT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESNASM only supports basic and restrict types");
288   nasm->type = type;
289   PetscFunctionReturn(PETSC_SUCCESS);
290 }
291 
292 /*@
293   SNESNASMGetType - Get the type of subdomain update used for the nonlinear additive Schwarz solver `SNESNASM`
294 
295   Logically Collective
296 
297   Input Parameter:
298 . snes - the `SNES` context
299 
300   Output Parameter:
301 . type - the type of update
302 
303   Level: intermediate
304 
305 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetType()`, `PCASMGetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType`
306 @*/
307 PetscErrorCode SNESNASMGetType(SNES snes, PCASMType *type)
308 {
309   PetscFunctionBegin;
310   PetscUseMethod(snes, "SNESNASMGetType_C", (SNES, PCASMType *), (snes, type));
311   PetscFunctionReturn(PETSC_SUCCESS);
312 }
313 
314 static PetscErrorCode SNESNASMGetType_NASM(SNES snes, PCASMType *type)
315 {
316   SNES_NASM *nasm = (SNES_NASM *)snes->data;
317 
318   PetscFunctionBegin;
319   *type = nasm->type;
320   PetscFunctionReturn(PETSC_SUCCESS);
321 }
322 
323 /*@
324   SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems in the nonlinear additive Schwarz solver
325 
326   Logically Collective
327 
328   Input Parameters:
329 + snes     - the `SNES` context
330 . n        - the number of local subdomains
331 . subsnes  - solvers defined on the local subdomains
332 . iscatter - scatters into the nonoverlapping portions of the local subdomains
333 . oscatter - scatters into the overlapping portions of the local subdomains
334 - gscatter - scatters into the (ghosted) local vector of the local subdomain
335 
336   Level: intermediate
337 
338 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()`
339 @*/
340 PetscErrorCode SNESNASMSetSubdomains(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[])
341 {
342   PetscErrorCode (*f)(SNES, PetscInt, SNES *, VecScatter *, VecScatter *, VecScatter *);
343 
344   PetscFunctionBegin;
345   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", &f));
346   if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter));
347   PetscFunctionReturn(PETSC_SUCCESS);
348 }
349 
350 static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[])
351 {
352   PetscInt   i;
353   SNES_NASM *nasm = (SNES_NASM *)snes->data;
354 
355   PetscFunctionBegin;
356   PetscCheck(!snes->setupcalled, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
357 
358   /* tear down the previously set things */
359   PetscCall(SNESReset(snes));
360 
361   nasm->n = n;
362   if (oscatter) {
363     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)oscatter[i]));
364   }
365   if (iscatter) {
366     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)iscatter[i]));
367   }
368   if (gscatter) {
369     for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)gscatter[i]));
370   }
371   if (oscatter) {
372     PetscCall(PetscMalloc1(n, &nasm->oscatter));
373     PetscCall(PetscMalloc1(n, &nasm->oscatter_copy));
374     for (i = 0; i < n; i++) {
375       nasm->oscatter[i] = oscatter[i];
376       PetscCall(VecScatterCopy(oscatter[i], &nasm->oscatter_copy[i]));
377     }
378   }
379   if (iscatter) {
380     PetscCall(PetscMalloc1(n, &nasm->iscatter));
381     for (i = 0; i < n; i++) nasm->iscatter[i] = iscatter[i];
382   }
383   if (gscatter) {
384     PetscCall(PetscMalloc1(n, &nasm->gscatter));
385     for (i = 0; i < n; i++) nasm->gscatter[i] = gscatter[i];
386   }
387 
388   if (subsnes) {
389     PetscCall(PetscMalloc1(n, &nasm->subsnes));
390     for (i = 0; i < n; i++) nasm->subsnes[i] = subsnes[i];
391   }
392   PetscFunctionReturn(PETSC_SUCCESS);
393 }
394 
395 /*@
396   SNESNASMGetSubdomains - Get the local subdomain contexts for the nonlinear additive Schwarz solver
397 
398   Not Collective but some of the objects returned will be parallel
399 
400   Input Parameter:
401 . snes - the `SNES` context
402 
403   Output Parameters:
404 + n        - the number of local subdomains
405 . subsnes  - solvers defined on the local subdomains
406 . iscatter - scatters into the nonoverlapping portions of the local subdomains
407 . oscatter - scatters into the overlapping portions of the local subdomains
408 - gscatter - scatters into the (ghosted) local vector of the local subdomain
409 
410   Level: intermediate
411 
412 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetSubdomains()`
413 @*/
414 PetscErrorCode SNESNASMGetSubdomains(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[])
415 {
416   PetscErrorCode (*f)(SNES, PetscInt *, SNES **, VecScatter **, VecScatter **, VecScatter **);
417 
418   PetscFunctionBegin;
419   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", &f));
420   if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter));
421   PetscFunctionReturn(PETSC_SUCCESS);
422 }
423 
424 static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[])
425 {
426   SNES_NASM *nasm = (SNES_NASM *)snes->data;
427 
428   PetscFunctionBegin;
429   if (n) *n = nasm->n;
430   if (oscatter) *oscatter = nasm->oscatter;
431   if (iscatter) *iscatter = nasm->iscatter;
432   if (gscatter) *gscatter = nasm->gscatter;
433   if (subsnes) *subsnes = nasm->subsnes;
434   PetscFunctionReturn(PETSC_SUCCESS);
435 }
436 
437 /*@
438   SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors for the nonlinear additive Schwarz solver
439 
440   Not Collective
441 
442   Input Parameter:
443 . snes - the `SNES` context
444 
445   Output Parameters:
446 + n  - the number of local subdomains
447 . x  - The subdomain solution vector
448 . y  - The subdomain step vector
449 . b  - The subdomain RHS vector
450 - xl - The subdomain local vectors (ghosted)
451 
452   Level: developer
453 
454 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()`
455 @*/
456 PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl)
457 {
458   PetscErrorCode (*f)(SNES, PetscInt *, Vec **, Vec **, Vec **, Vec **);
459 
460   PetscFunctionBegin;
461   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", &f));
462   if (f) PetscCall((f)(snes, n, x, y, b, xl));
463   PetscFunctionReturn(PETSC_SUCCESS);
464 }
465 
466 static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl)
467 {
468   SNES_NASM *nasm = (SNES_NASM *)snes->data;
469 
470   PetscFunctionBegin;
471   if (n) *n = nasm->n;
472   if (x) *x = nasm->x;
473   if (y) *y = nasm->y;
474   if (b) *b = nasm->b;
475   if (xl) *xl = nasm->xl;
476   PetscFunctionReturn(PETSC_SUCCESS);
477 }
478 
479 /*@
480   SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain Jacobians upon convergence for the
481   nonlinear additive Schwarz solver
482 
483   Collective
484 
485   Input Parameters:
486 + snes - the SNES context
487 - flg  - `PETSC_TRUE` to compute the Jacobians
488 
489   Level: developer
490 
491   Notes:
492   This is used almost exclusively in the implementation of `SNESASPIN`, where the converged subdomain and global Jacobian
493   is needed at each linear iteration.
494 
495 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()`
496 @*/
497 PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes, PetscBool flg)
498 {
499   PetscErrorCode (*f)(SNES, PetscBool);
500 
501   PetscFunctionBegin;
502   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", &f));
503   if (f) PetscCall((f)(snes, flg));
504   PetscFunctionReturn(PETSC_SUCCESS);
505 }
506 
507 static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes, PetscBool flg)
508 {
509   SNES_NASM *nasm = (SNES_NASM *)snes->data;
510 
511   PetscFunctionBegin;
512   nasm->finaljacobian = flg;
513   PetscFunctionReturn(PETSC_SUCCESS);
514 }
515 
516 /*@
517   SNESNASMSetDamping - Sets the update damping for `SNESNASM` the nonlinear additive Schwarz solver
518 
519   Logically Collective
520 
521   Input Parameters:
522 + snes - the `SNES` context
523 - dmp  - damping
524 
525   Options Database Key:
526 . -snes_nasm_damping <dmp> - the new solution is obtained as old solution plus `dmp` times (sum of the solutions on the subdomains)
527 
528   Level: intermediate
529 
530   Note:
531   The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
532 
533 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetDamping()`
534 @*/
535 PetscErrorCode SNESNASMSetDamping(SNES snes, PetscReal dmp)
536 {
537   PetscErrorCode (*f)(SNES, PetscReal);
538 
539   PetscFunctionBegin;
540   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetDamping_C", (void (**)(void)) & f));
541   if (f) PetscCall((f)(snes, dmp));
542   PetscFunctionReturn(PETSC_SUCCESS);
543 }
544 
545 static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes, PetscReal dmp)
546 {
547   SNES_NASM *nasm = (SNES_NASM *)snes->data;
548 
549   PetscFunctionBegin;
550   nasm->damping = dmp;
551   PetscFunctionReturn(PETSC_SUCCESS);
552 }
553 
554 /*@
555   SNESNASMGetDamping - Gets the update damping for `SNESNASM` the nonlinear additive Schwarz solver
556 
557   Not Collective
558 
559   Input Parameter:
560 . snes - the `SNES` context
561 
562   Output Parameter:
563 . dmp - damping
564 
565   Level: intermediate
566 
567 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetDamping()`
568 @*/
569 PetscErrorCode SNESNASMGetDamping(SNES snes, PetscReal *dmp)
570 {
571   PetscFunctionBegin;
572   PetscUseMethod(snes, "SNESNASMGetDamping_C", (SNES, PetscReal *), (snes, dmp));
573   PetscFunctionReturn(PETSC_SUCCESS);
574 }
575 
576 static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes, PetscReal *dmp)
577 {
578   SNES_NASM *nasm = (SNES_NASM *)snes->data;
579 
580   PetscFunctionBegin;
581   *dmp = nasm->damping;
582   PetscFunctionReturn(PETSC_SUCCESS);
583 }
584 
585 /*
586   Input Parameters:
587 + snes - The solver
588 . B - The RHS vector
589 - X - The initial guess
590 
591   Output Parameter:
592 . Y - The solution update
593 
594   TODO: All scatters should be packed into one
595 */
596 static PetscErrorCode SNESNASMSolveLocal_Private(SNES snes, Vec B, Vec Y, Vec X)
597 {
598   SNES_NASM *nasm = (SNES_NASM *)snes->data;
599   SNES       subsnes;
600   PetscInt   i;
601   PetscReal  dmp;
602   Vec        Xl, Bl, Yl, Xlloc;
603   VecScatter iscat, oscat, gscat, oscat_copy;
604   DM         dm, subdm;
605   PCASMType  type;
606 
607   PetscFunctionBegin;
608   PetscCall(SNESNASMGetType(snes, &type));
609   PetscCall(SNESGetDM(snes, &dm));
610   PetscCall(VecSet(Y, 0));
611   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
612   for (i = 0; i < nasm->n; i++) {
613     /* scatter the solution to the global solution and the local solution */
614     Xl         = nasm->x[i];
615     Xlloc      = nasm->xl[i];
616     oscat      = nasm->oscatter[i];
617     oscat_copy = nasm->oscatter_copy[i];
618     gscat      = nasm->gscatter[i];
619     PetscCall(VecScatterBegin(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD));
620     PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
621     if (B) {
622       /* scatter the RHS to the local RHS */
623       Bl = nasm->b[i];
624       PetscCall(VecScatterBegin(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD));
625     }
626   }
627   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
628 
629   if (nasm->eventsubsolve) PetscCall(PetscLogEventBegin(nasm->eventsubsolve, snes, 0, 0, 0));
630   for (i = 0; i < nasm->n; i++) {
631     Xl      = nasm->x[i];
632     Xlloc   = nasm->xl[i];
633     Yl      = nasm->y[i];
634     subsnes = nasm->subsnes[i];
635     PetscCall(SNESGetDM(subsnes, &subdm));
636     iscat      = nasm->iscatter[i];
637     oscat      = nasm->oscatter[i];
638     oscat_copy = nasm->oscatter_copy[i];
639     gscat      = nasm->gscatter[i];
640     PetscCall(VecScatterEnd(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD));
641     PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
642     if (B) {
643       Bl = nasm->b[i];
644       PetscCall(VecScatterEnd(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD));
645     } else Bl = NULL;
646 
647     PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm));
648     PetscCall(VecCopy(Xl, Yl));
649     PetscCall(SNESSolve(subsnes, Bl, Xl));
650     PetscCall(VecAYPX(Yl, -1.0, Xl));
651     PetscCall(VecScale(Yl, nasm->damping));
652     if (type == PC_ASM_BASIC) {
653       PetscCall(VecScatterBegin(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
654       PetscCall(VecScatterEnd(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
655     } else if (type == PC_ASM_RESTRICT) {
656       PetscCall(VecScatterBegin(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
657       PetscCall(VecScatterEnd(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE));
658     } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Only basic and restrict types are supported for SNESNASM");
659   }
660   if (nasm->eventsubsolve) PetscCall(PetscLogEventEnd(nasm->eventsubsolve, snes, 0, 0, 0));
661   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
662   if (nasm->weight_set) PetscCall(VecPointwiseMult(Y, Y, nasm->weight));
663   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
664   PetscCall(SNESNASMGetDamping(snes, &dmp));
665   PetscCall(VecAXPY(X, dmp, Y));
666   PetscFunctionReturn(PETSC_SUCCESS);
667 }
668 
669 static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal)
670 {
671   Vec        X    = Xfinal;
672   SNES_NASM *nasm = (SNES_NASM *)snes->data;
673   SNES       subsnes;
674   PetscInt   i, lag = 1;
675   Vec        Xlloc, Xl, Fl, F;
676   VecScatter oscat, gscat;
677   DM         dm, subdm;
678 
679   PetscFunctionBegin;
680   if (nasm->fjtype == 2) X = nasm->xinit;
681   F = snes->vec_func;
682   if (snes->normschedule == SNES_NORM_NONE) PetscCall(SNESComputeFunction(snes, X, F));
683   PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
684   PetscCall(SNESGetDM(snes, &dm));
685   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0));
686   if (nasm->fjtype != 1) {
687     for (i = 0; i < nasm->n; i++) {
688       Xlloc = nasm->xl[i];
689       gscat = nasm->gscatter[i];
690       PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
691     }
692   }
693   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0));
694   for (i = 0; i < nasm->n; i++) {
695     Fl      = nasm->subsnes[i]->vec_func;
696     Xl      = nasm->x[i];
697     Xlloc   = nasm->xl[i];
698     subsnes = nasm->subsnes[i];
699     oscat   = nasm->oscatter[i];
700     gscat   = nasm->gscatter[i];
701     if (nasm->fjtype != 1) PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD));
702     PetscCall(SNESGetDM(subsnes, &subdm));
703     PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm));
704     if (nasm->fjtype != 1) {
705       PetscCall(DMLocalToGlobalBegin(subdm, Xlloc, INSERT_VALUES, Xl));
706       PetscCall(DMLocalToGlobalEnd(subdm, Xlloc, INSERT_VALUES, Xl));
707     }
708     if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2;
709     else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian;
710     PetscCall(SNESComputeFunction(subsnes, Xl, Fl));
711     PetscCall(SNESComputeJacobian(subsnes, Xl, subsnes->jacobian, subsnes->jacobian_pre));
712     if (lag > 1) subsnes->lagjacobian = lag;
713   }
714   PetscFunctionReturn(PETSC_SUCCESS);
715 }
716 
717 static PetscErrorCode SNESSolve_NASM(SNES snes)
718 {
719   Vec              F;
720   Vec              X;
721   Vec              B;
722   Vec              Y;
723   PetscInt         i;
724   PetscReal        fnorm = 0.0;
725   SNESNormSchedule normschedule;
726   SNES_NASM       *nasm = (SNES_NASM *)snes->data;
727 
728   PetscFunctionBegin;
729 
730   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);
731 
732   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
733   X = snes->vec_sol;
734   Y = snes->vec_sol_update;
735   F = snes->vec_func;
736   B = snes->vec_rhs;
737 
738   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
739   snes->iter = 0;
740   snes->norm = 0.;
741   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
742   snes->reason = SNES_CONVERGED_ITERATING;
743   PetscCall(SNESGetNormSchedule(snes, &normschedule));
744   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
745     /* compute the initial function and preconditioned update delX */
746     if (!snes->vec_func_init_set) {
747       PetscCall(SNESComputeFunction(snes, X, F));
748     } else snes->vec_func_init_set = PETSC_FALSE;
749 
750     PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
751     SNESCheckFunctionNorm(snes, fnorm);
752     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
753     snes->iter = 0;
754     snes->norm = fnorm;
755     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
756     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
757 
758     /* test convergence */
759     PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
760     PetscCall(SNESMonitor(snes, 0, snes->norm));
761     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
762   } else {
763     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
764     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
765     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
766   }
767 
768   /* Call general purpose update function */
769   PetscTryTypeMethod(snes, update, snes->iter);
770   /* copy the initial solution over for later */
771   if (nasm->fjtype == 2) PetscCall(VecCopy(X, nasm->xinit));
772 
773   for (i = 0; i < snes->max_its; i++) {
774     PetscCall(SNESNASMSolveLocal_Private(snes, B, Y, X));
775     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
776       PetscCall(SNESComputeFunction(snes, X, F));
777       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
778       SNESCheckFunctionNorm(snes, fnorm);
779     }
780     /* Monitor convergence */
781     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
782     snes->iter = i + 1;
783     snes->norm = fnorm;
784     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
785     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
786     /* Test for convergence */
787     PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm));
788     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
789     if (snes->reason) break;
790     /* Call general purpose update function */
791     PetscTryTypeMethod(snes, update, snes->iter);
792   }
793   if (nasm->finaljacobian) {
794     PetscCall(SNESNASMComputeFinalJacobian_Private(snes, X));
795     SNESCheckJacobianDomainerror(snes);
796   }
797   PetscFunctionReturn(PETSC_SUCCESS);
798 }
799 
800 /*MC
801   SNESNASM - Nonlinear Additive Schwarz solver
802 
803    Options Database Keys:
804 +  -snes_nasm_log - enable logging events for the communication and solve stages
805 .  -snes_nasm_type <basic,restrict> - type of subdomain update used
806 .  -snes_nasm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
807 .  -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate
808 .  -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the jacobian is calculated at
809 .  -sub_snes_ - options prefix of the subdomain nonlinear solves
810 .  -sub_ksp_ - options prefix of the subdomain Krylov solver
811 -  -sub_pc_ - options prefix of the subdomain preconditioner
812 
813    Level: advanced
814 
815    Note:
816    This is not often used directly as a solver, it converges too slowly. However it works well as a nonlinear preconditioner for
817    the `SNESASPIN` solver
818 
819    Developer Note:
820    This is a non-Newton based nonlinear solver that does not directly require a Jacobian; hence the flag snes->usesksp is set to
821        false and `SNESView()` and -snes_view do not display a `KSP` object. However, if the flag nasm->finaljacobian is set (for example, if
822        `SNESNASM` is used as a nonlinear preconditioner for  `SNESASPIN`) then `SNESSetUpMatrices()` is called to generate the
823        Jacobian (needed by `SNESASPIN`)
824        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
825        used by `SNESASPIN` they share the same Jacobian matrices because `SNESSetUp()` (called on the outer `SNESASPIN`) causes the inner `SNES`
826        object (in this case `SNESNASM`) to inherit the outer Jacobian matrices.
827 
828    References:
829 .  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
830    SIAM Review, 57(4), 2015
831 
832 .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESNASMSetType()`, `SNESNASMGetType()`, `SNESNASMSetSubdomains()`, `SNESNASMGetSubdomains()`,
833           `SNESNASMGetSubdomainVecs()`, `SNESNASMSetComputeFinalJacobian()`, `SNESNASMSetDamping()`, `SNESNASMGetDamping()`, `SNESNASMSetWeight()`,
834           `SNESNASMGetSNES()`, `SNESNASMGetNumber()`
835 M*/
836 
837 PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes)
838 {
839   SNES_NASM *nasm;
840 
841   PetscFunctionBegin;
842   PetscCall(PetscNew(&nasm));
843   snes->data = (void *)nasm;
844 
845   nasm->n             = PETSC_DECIDE;
846   nasm->subsnes       = NULL;
847   nasm->x             = NULL;
848   nasm->xl            = NULL;
849   nasm->y             = NULL;
850   nasm->b             = NULL;
851   nasm->oscatter      = NULL;
852   nasm->oscatter_copy = NULL;
853   nasm->iscatter      = NULL;
854   nasm->gscatter      = NULL;
855   nasm->damping       = 1.;
856 
857   nasm->type          = PC_ASM_BASIC;
858   nasm->finaljacobian = PETSC_FALSE;
859   nasm->weight_set    = PETSC_FALSE;
860 
861   snes->ops->destroy        = SNESDestroy_NASM;
862   snes->ops->setup          = SNESSetUp_NASM;
863   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
864   snes->ops->view           = SNESView_NASM;
865   snes->ops->solve          = SNESSolve_NASM;
866   snes->ops->reset          = SNESReset_NASM;
867 
868   snes->usesksp = PETSC_FALSE;
869   snes->usesnpc = PETSC_FALSE;
870 
871   snes->alwayscomputesfinalresidual = PETSC_FALSE;
872 
873   nasm->fjtype              = 0;
874   nasm->xinit               = NULL;
875   nasm->eventrestrictinterp = 0;
876   nasm->eventsubsolve       = 0;
877 
878   if (!snes->tolerancesset) {
879     snes->max_its   = 10000;
880     snes->max_funcs = 10000;
881   }
882 
883   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", SNESNASMSetType_NASM));
884   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", SNESNASMGetType_NASM));
885   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", SNESNASMSetSubdomains_NASM));
886   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", SNESNASMGetSubdomains_NASM));
887   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", SNESNASMSetDamping_NASM));
888   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", SNESNASMGetDamping_NASM));
889   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", SNESNASMGetSubdomainVecs_NASM));
890   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", SNESNASMSetComputeFinalJacobian_NASM));
891   PetscFunctionReturn(PETSC_SUCCESS);
892 }
893 
894 /*@
895   SNESNASMGetSNES - Gets a subsolver
896 
897   Not Collective
898 
899   Input Parameters:
900 + snes - the `SNES` context
901 - i    - the number of the subsnes to get
902 
903   Output Parameter:
904 . subsnes - the subsolver context
905 
906   Level: intermediate
907 
908 .seealso: [](ch_snes), `SNESNASM`, `SNESNASMGetNumber()`
909 @*/
910 PetscErrorCode SNESNASMGetSNES(SNES snes, PetscInt i, SNES *subsnes)
911 {
912   SNES_NASM *nasm = (SNES_NASM *)snes->data;
913 
914   PetscFunctionBegin;
915   PetscCheck(i >= 0 && i < nasm->n, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "No such subsolver");
916   *subsnes = nasm->subsnes[i];
917   PetscFunctionReturn(PETSC_SUCCESS);
918 }
919 
920 /*@
921   SNESNASMGetNumber - Gets number of subsolvers
922 
923   Not Collective
924 
925   Input Parameter:
926 . snes - the `SNES` context
927 
928   Output Parameter:
929 . n - the number of subsolvers
930 
931   Level: intermediate
932 
933 .seealso: [](ch_snes), `SNESNASM`, `SNESNASMGetSNES()`
934 @*/
935 PetscErrorCode SNESNASMGetNumber(SNES snes, PetscInt *n)
936 {
937   SNES_NASM *nasm = (SNES_NASM *)snes->data;
938 
939   PetscFunctionBegin;
940   *n = nasm->n;
941   PetscFunctionReturn(PETSC_SUCCESS);
942 }
943 
944 /*@
945   SNESNASMSetWeight - Sets weight to use when adding overlapping updates
946 
947   Collective
948 
949   Input Parameters:
950 + snes   - the `SNES` context
951 - weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in)
952 
953   Level: intermediate
954 
955 .seealso: [](ch_snes), `SNESNASM`
956 @*/
957 PetscErrorCode SNESNASMSetWeight(SNES snes, Vec weight)
958 {
959   SNES_NASM *nasm = (SNES_NASM *)snes->data;
960 
961   PetscFunctionBegin;
962 
963   PetscCall(VecDestroy(&nasm->weight));
964   nasm->weight_set = PETSC_TRUE;
965   nasm->weight     = weight;
966   PetscCall(PetscObjectReference((PetscObject)nasm->weight));
967 
968   PetscFunctionReturn(PETSC_SUCCESS);
969 }
970