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