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