xref: /petsc/src/snes/impls/nasm/nasm.c (revision 5dc64a9727a0cc1147601b21df2b9dcd374fa7c9)
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   PetscFunctionReturn(0);
559 }
560 
561 /*@
562    SNESNASMSetDamping - Sets the update damping for NASM
563 
564    Logically collective on SNES
565 
566    Input Parameters:
567 +  SNES - the SNES context
568 -  dmp - damping
569 
570    Level: intermediate
571 
572    Notes:
573     The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
574 
575 .keywords: SNES, NASM, damping
576 
577 .seealso: SNESNASM, SNESNASMGetDamping()
578 @*/
579 PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp)
580 {
581   PetscErrorCode (*f)(SNES,PetscReal);
582   PetscErrorCode ierr;
583 
584   PetscFunctionBegin;
585   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr);
586   if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);}
587   PetscFunctionReturn(0);
588 }
589 
590 static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp)
591 {
592   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
593 
594   PetscFunctionBegin;
595   nasm->damping = dmp;
596   PetscFunctionReturn(0);
597 }
598 
599 /*@
600    SNESNASMGetDamping - Gets the update damping for NASM
601 
602    Not Collective
603 
604    Input Parameters:
605 +  SNES - the SNES context
606 -  dmp - damping
607 
608    Level: intermediate
609 
610 .keywords: SNES, NASM, damping
611 
612 .seealso: SNESNASM, SNESNASMSetDamping()
613 @*/
614 PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp)
615 {
616   PetscErrorCode ierr;
617 
618   PetscFunctionBegin;
619   ierr = PetscUseMethod(snes,"SNESNASMGetDamping_C",(SNES,PetscReal*),(snes,dmp));CHKERRQ(ierr);
620   PetscFunctionReturn(0);
621 }
622 
623 static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp)
624 {
625   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
626 
627   PetscFunctionBegin;
628   *dmp = nasm->damping;
629   PetscFunctionReturn(0);
630 }
631 
632 
633 /*
634   Input Parameters:
635 + snes - The solver
636 . B - The RHS vector
637 - X - The initial guess
638 
639   Output Parameters:
640 . Y - The solution update
641 
642   TODO: All scatters should be packed into one
643 */
644 PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X)
645 {
646   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
647   SNES           subsnes;
648   PetscInt       i;
649   PetscReal      dmp;
650   PetscErrorCode ierr;
651   Vec            Xl,Bl,Yl,Xlloc;
652   VecScatter     iscat,oscat,gscat,oscat_copy;
653   DM             dm,subdm;
654   PCASMType      type;
655 
656   PetscFunctionBegin;
657   ierr = SNESNASMGetType(snes,&type);CHKERRQ(ierr);
658   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
659   ierr = VecSet(Y,0);CHKERRQ(ierr);
660   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
661   for (i=0; i<nasm->n; i++) {
662     /* scatter the solution to the global solution and the local solution */
663     Xl      = nasm->x[i];
664     Xlloc   = nasm->xl[i];
665     oscat   = nasm->oscatter[i];
666     oscat_copy = nasm->oscatter_copy[i];
667     gscat   = nasm->gscatter[i];
668     ierr = VecScatterBegin(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
669     ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
670     if (B) {
671       /* scatter the RHS to the local RHS */
672       Bl   = nasm->b[i];
673       ierr = VecScatterBegin(oscat_copy,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
674     }
675   }
676   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
677 
678 
679   if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
680   for (i=0; i<nasm->n; i++) {
681     Xl    = nasm->x[i];
682     Xlloc = nasm->xl[i];
683     Yl    = nasm->y[i];
684     subsnes = nasm->subsnes[i];
685     ierr    = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
686     iscat   = nasm->iscatter[i];
687     oscat   = nasm->oscatter[i];
688     oscat_copy = nasm->oscatter_copy[i];
689     gscat   = nasm->gscatter[i];
690     ierr = VecScatterEnd(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
691     ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
692     if (B) {
693       Bl   = nasm->b[i];
694       ierr = VecScatterEnd(oscat_copy,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
695     } else Bl = NULL;
696 
697     ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
698     ierr = VecCopy(Xl,Yl);CHKERRQ(ierr);
699     ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr);
700     ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr);
701     ierr = VecScale(Yl, nasm->damping);CHKERRQ(ierr);
702     if (type == PC_ASM_BASIC) {
703       ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
704     } else if (type == PC_ASM_RESTRICT) {
705       ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
706     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
707   }
708   if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
709   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
710   for (i=0; i<nasm->n; i++) {
711     Yl    = nasm->y[i];
712     iscat   = nasm->iscatter[i];
713     oscat   = nasm->oscatter[i];
714     if (type == PC_ASM_BASIC) {
715       ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
716     } else if (type == PC_ASM_RESTRICT) {
717       ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
718     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
719   }
720   if (nasm->weight_set) {
721     ierr = VecPointwiseMult(Y,Y,nasm->weight);CHKERRQ(ierr);
722   }
723   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
724   ierr = SNESNASMGetDamping(snes,&dmp);CHKERRQ(ierr);
725   ierr = VecAXPY(X,dmp,Y);CHKERRQ(ierr);
726   PetscFunctionReturn(0);
727 }
728 
729 static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal)
730 {
731   Vec            X = Xfinal;
732   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
733   SNES           subsnes;
734   PetscInt       i,lag = 1;
735   PetscErrorCode ierr;
736   Vec            Xlloc,Xl,Fl,F;
737   VecScatter     oscat,gscat;
738   DM             dm,subdm;
739 
740   PetscFunctionBegin;
741   if (nasm->fjtype == 2) X = nasm->xinit;
742   F = snes->vec_func;
743   if (snes->normschedule == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);}
744   ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
745   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
746   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
747   if (nasm->fjtype != 1) {
748     for (i=0; i<nasm->n; i++) {
749       Xlloc = nasm->xl[i];
750       gscat = nasm->gscatter[i];
751       ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
752     }
753   }
754   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
755   for (i=0; i<nasm->n; i++) {
756     Fl      = nasm->subsnes[i]->vec_func;
757     Xl      = nasm->x[i];
758     Xlloc   = nasm->xl[i];
759     subsnes = nasm->subsnes[i];
760     oscat   = nasm->oscatter[i];
761     gscat   = nasm->gscatter[i];
762     if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);}
763     ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
764     ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
765     if (nasm->fjtype != 1) {
766       ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
767       ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
768     }
769     if (subsnes->lagjacobian == -1)    subsnes->lagjacobian = -2;
770     else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian;
771     ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr);
772     ierr = SNESComputeJacobian(subsnes,Xl,subsnes->jacobian,subsnes->jacobian_pre);CHKERRQ(ierr);
773     if (lag > 1) subsnes->lagjacobian = lag;
774   }
775   PetscFunctionReturn(0);
776 }
777 
778 static PetscErrorCode SNESSolve_NASM(SNES snes)
779 {
780   Vec              F;
781   Vec              X;
782   Vec              B;
783   Vec              Y;
784   PetscInt         i;
785   PetscReal        fnorm = 0.0;
786   PetscErrorCode   ierr;
787   SNESNormSchedule normschedule;
788   SNES_NASM        *nasm = (SNES_NASM*)snes->data;
789 
790   PetscFunctionBegin;
791 
792   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);
793 
794   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
795   X = snes->vec_sol;
796   Y = snes->vec_sol_update;
797   F = snes->vec_func;
798   B = snes->vec_rhs;
799 
800   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
801   snes->iter   = 0;
802   snes->norm   = 0.;
803   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
804   snes->reason = SNES_CONVERGED_ITERATING;
805   ierr         = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr);
806   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
807     /* compute the initial function and preconditioned update delX */
808     if (!snes->vec_func_init_set) {
809       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
810     } else snes->vec_func_init_set = PETSC_FALSE;
811 
812     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
813     SNESCheckFunctionNorm(snes,fnorm);
814     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
815     snes->iter = 0;
816     snes->norm = fnorm;
817     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
818     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
819     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
820 
821     /* test convergence */
822     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
823     if (snes->reason) PetscFunctionReturn(0);
824   } else {
825     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
826     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
827     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
828   }
829 
830   /* Call general purpose update function */
831   if (snes->ops->update) {
832     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
833   }
834   /* copy the initial solution over for later */
835   if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);}
836 
837   for (i=0; i < snes->max_its; i++) {
838     ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr);
839     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
840       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
841       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
842       SNESCheckFunctionNorm(snes,fnorm);
843     }
844     /* Monitor convergence */
845     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
846     snes->iter = i+1;
847     snes->norm = fnorm;
848     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
849     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
850     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
851     /* Test for convergence */
852     if (normschedule == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
853     if (snes->reason) break;
854     /* Call general purpose update function */
855     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
856   }
857   if (nasm->finaljacobian) {
858     ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);
859     SNESCheckJacobianDomainerror(snes);
860   }
861   if (normschedule == SNES_NORM_ALWAYS) {
862     if (i == snes->max_its) {
863       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
864       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
865     }
866   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */
867   PetscFunctionReturn(0);
868 }
869 
870 /*MC
871   SNESNASM - Nonlinear Additive Schwartz
872 
873    Options Database:
874 +  -snes_nasm_log - enable logging events for the communication and solve stages
875 .  -snes_nasm_type <basic,restrict> - type of subdomain update used
876 .  -snes_asm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
877 .  -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate
878 .  -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the jacobian is calculated at
879 .  -sub_snes_ - options prefix of the subdomain nonlinear solves
880 .  -sub_ksp_ - options prefix of the subdomain Krylov solver
881 -  -sub_pc_ - options prefix of the subdomain preconditioner
882 
883    Level: advanced
884 
885    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
886        false and SNESView() and -snes_view do not display a KSP object. However the flag nasm->finaljacobian is set (for example if
887        NASM is used as a nonlinear preconditioner for  KSPASPIN) then SNESSetUpMatrices() is called to generate the Jacobian (needed by KSPASPIN)
888        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
889        used by SNESASPIN they share the same Jacobian matrices because SNESSetUp() (called on the outer SNES KSPASPIN) causes the inner SNES
890        object (in this case SNESNASM) to inherit the outer Jacobian matrices.
891 
892    References:
893 .  1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
894    SIAM Review, 57(4), 2015
895 
896 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types), SNESNASMSetType(), SNESNASMGetType(), SNESNASMSetSubdomains(), SNESNASMGetSubdomains(), SNESNASMGetSubdomainVecs(), SNESNASMSetComputeFinalJacobian(), SNESNASMSetDamping(), SNESNASMGetDamping()
897 M*/
898 
899 PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes)
900 {
901   SNES_NASM      *nasm;
902   PetscErrorCode ierr;
903 
904   PetscFunctionBegin;
905   ierr       = PetscNewLog(snes,&nasm);CHKERRQ(ierr);
906   snes->data = (void*)nasm;
907 
908   nasm->n        = PETSC_DECIDE;
909   nasm->subsnes  = 0;
910   nasm->x        = 0;
911   nasm->xl       = 0;
912   nasm->y        = 0;
913   nasm->b        = 0;
914   nasm->oscatter = 0;
915   nasm->oscatter_copy = 0;
916   nasm->iscatter = 0;
917   nasm->gscatter = 0;
918   nasm->damping  = 1.;
919 
920   nasm->type              = PC_ASM_BASIC;
921   nasm->finaljacobian     = PETSC_FALSE;
922   nasm->same_local_solves = PETSC_TRUE;
923   nasm->weight_set        = PETSC_FALSE;
924 
925   snes->ops->destroy        = SNESDestroy_NASM;
926   snes->ops->setup          = SNESSetUp_NASM;
927   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
928   snes->ops->view           = SNESView_NASM;
929   snes->ops->solve          = SNESSolve_NASM;
930   snes->ops->reset          = SNESReset_NASM;
931 
932   snes->usesksp = PETSC_FALSE;
933   snes->usesnpc = PETSC_FALSE;
934 
935   snes->alwayscomputesfinalresidual = PETSC_FALSE;
936 
937   nasm->fjtype              = 0;
938   nasm->xinit               = NULL;
939   nasm->eventrestrictinterp = 0;
940   nasm->eventsubsolve       = 0;
941 
942   if (!snes->tolerancesset) {
943     snes->max_its   = 10000;
944     snes->max_funcs = 10000;
945   }
946 
947   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetType_C",SNESNASMSetType_NASM);CHKERRQ(ierr);
948   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetType_C",SNESNASMGetType_NASM);CHKERRQ(ierr);
949   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr);
950   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr);
951   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",SNESNASMSetDamping_NASM);CHKERRQ(ierr);
952   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",SNESNASMGetDamping_NASM);CHKERRQ(ierr);
953   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr);
954   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr);
955   PetscFunctionReturn(0);
956 }
957 
958 /*@
959    SNESNASMGetSNES - Gets a subsolver
960 
961    Not collective
962 
963    Input Parameters:
964 +  snes - the SNES context
965 -  i - the number of the subsnes to get
966 
967    Output Parameters:
968 .  subsnes - the subsolver context
969 
970    Level: intermediate
971 
972 .keywords: SNES, NASM
973 
974 .seealso: SNESNASM, SNESNASMGetNumber()
975 @*/
976 PetscErrorCode SNESNASMGetSNES(SNES snes,PetscInt i,SNES *subsnes)
977 {
978   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
979 
980   PetscFunctionBegin;
981   if (i < 0 || i >= nasm->n) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"No such subsolver");
982   *subsnes = nasm->subsnes[i];
983   PetscFunctionReturn(0);
984 }
985 
986 /*@
987    SNESNASMGetNumber - Gets number of subsolvers
988 
989    Not collective
990 
991    Input Parameters:
992 .  snes - the SNES context
993 
994    Output Parameters:
995 .  n - the number of subsolvers
996 
997    Level: intermediate
998 
999 .keywords: SNES, NASM
1000 
1001 .seealso: SNESNASM, SNESNASMGetSNES()
1002 @*/
1003 PetscErrorCode SNESNASMGetNumber(SNES snes,PetscInt *n)
1004 {
1005   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
1006 
1007   PetscFunctionBegin;
1008   *n = nasm->n;
1009   PetscFunctionReturn(0);
1010 }
1011 
1012 /*@
1013    SNESNASMSetWeight - Sets weight to use when adding overlapping updates
1014 
1015    Collective
1016 
1017    Input Parameters:
1018 +  snes - the SNES context
1019 -  weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in)
1020 
1021    Level: intermediate
1022 
1023 .keywords: SNES, NASM
1024 
1025 .seealso: SNESNASM
1026 @*/
1027 PetscErrorCode SNESNASMSetWeight(SNES snes,Vec weight)
1028 {
1029   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
1030   PetscErrorCode ierr;
1031 
1032   PetscFunctionBegin;
1033 
1034   ierr = VecDestroy(&nasm->weight);CHKERRQ(ierr);
1035   nasm->weight_set = PETSC_TRUE;
1036   nasm->weight     = weight;
1037   ierr = PetscObjectReference((PetscObject)nasm->weight);CHKERRQ(ierr);
1038 
1039   PetscFunctionReturn(0);
1040 }
1041