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