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