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