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