xref: /petsc/src/snes/impls/nasm/nasm.c (revision 79f8d099bee1bd3e2e41b7725692d5d7f212fcc1)
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   PetscBool  same_local_solves;   /* flag to determine if the solvers have been individually modified */
18 
19   /* logging events */
20   PetscLogEvent eventrestrictinterp;
21   PetscLogEvent eventsubsolve;
22 } SNES_NASM;
23 
24 const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0};
25 
26 #undef __FUNCT__
27 #define __FUNCT__ "SNESReset_NASM"
28 PetscErrorCode SNESReset_NASM(SNES snes)
29 {
30   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
31   PetscErrorCode ierr;
32   PetscInt       i;
33 
34   PetscFunctionBegin;
35   for (i=0; i<nasm->n; i++) {
36     if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); }
37     if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); }
38     if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); }
39     if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); }
40 
41     if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); }
42     if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); }
43     if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); }
44     if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); }
45   }
46 
47   if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);}
48   if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);}
49   if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);}
50   if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);}
51 
52   if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);}
53   if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);}
54   if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);}
55   if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);}
56 
57   nasm->eventrestrictinterp = 0;
58   nasm->eventsubsolve = 0;
59   PetscFunctionReturn(0);
60 }
61 
62 #undef __FUNCT__
63 #define __FUNCT__ "SNESDestroy_NASM"
64 PetscErrorCode SNESDestroy_NASM(SNES snes)
65 {
66   PetscErrorCode ierr;
67 
68   PetscFunctionBegin;
69   ierr = SNESReset_NASM(snes);CHKERRQ(ierr);
70   ierr = PetscFree(snes->data);CHKERRQ(ierr);
71   PetscFunctionReturn(0);
72 }
73 
74 #undef __FUNCT__
75 #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private"
76 PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx)
77 {
78   PetscErrorCode ierr;
79   Vec            bcs = (Vec)ctx;
80 
81   PetscFunctionBegin;
82   ierr = VecCopy(bcs,l);CHKERRQ(ierr);
83   PetscFunctionReturn(0);
84 }
85 
86 #undef __FUNCT__
87 #define __FUNCT__ "SNESSetUp_NASM"
88 PetscErrorCode SNESSetUp_NASM(SNES snes)
89 {
90   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
91   PetscErrorCode ierr;
92   DM             dm,subdm;
93   DM             *subdms;
94   PetscInt       i;
95   const char     *optionsprefix;
96   Vec            F;
97   PetscMPIInt    size;
98   KSP            ksp;
99   PC             pc;
100 
101   PetscFunctionBegin;
102   if (!nasm->subsnes) {
103     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
104     if (dm) {
105       nasm->usesdm = PETSC_TRUE;
106       ierr         = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr);
107       if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
108       ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr);
109 
110       ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
111       ierr = PetscMalloc(nasm->n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr);
112       for (i=0; i<nasm->n; i++) {
113         ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr);
114         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr);
115         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr);
116         ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr);
117         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRQ(ierr);
118         if (size == 1) {
119           ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr);
120           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
121           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
122           ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
123         }
124         ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr);
125         ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr);
126       }
127       ierr = PetscFree(subdms);CHKERRQ(ierr);
128     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!");
129   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!");
130   /* allocate the global vectors */
131   if (!nasm->x) {
132     ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->x);CHKERRQ(ierr);
133     ierr = PetscMemzero(nasm->x,nasm->n*sizeof(Vec));CHKERRQ(ierr);
134   }
135   if (!nasm->xl) {
136     ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->xl);CHKERRQ(ierr);
137     ierr = PetscMemzero(nasm->xl,nasm->n*sizeof(Vec));CHKERRQ(ierr);
138   }
139   if (!nasm->y) {
140     ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->y);CHKERRQ(ierr);
141     ierr = PetscMemzero(nasm->y,nasm->n*sizeof(Vec));CHKERRQ(ierr);
142   }
143   if (!nasm->b) {
144     ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->b);CHKERRQ(ierr);
145     ierr = PetscMemzero(nasm->b,nasm->n*sizeof(Vec));CHKERRQ(ierr);
146   }
147 
148   for (i=0; i<nasm->n; i++) {
149     ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr);
150     if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);}
151     if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);}
152     if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);}
153     if (!nasm->xl[i]) {
154       ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr);
155       ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr);
156     }
157     ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr);
158   }
159   if (nasm->finaljacobian) {ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);}
160   PetscFunctionReturn(0);
161 }
162 
163 #undef __FUNCT__
164 #define __FUNCT__ "SNESSetFromOptions_NASM"
165 PetscErrorCode SNESSetFromOptions_NASM(SNES snes)
166 {
167   PetscErrorCode    ierr;
168   PCASMType         asmtype;
169   PetscBool         flg,monflg,subviewflg;
170   SNES_NASM         *nasm = (SNES_NASM*)snes->data;
171 
172   PetscFunctionBegin;
173   ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr);
174   ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
175   if (flg) nasm->type = asmtype;
176   flg    = PETSC_FALSE;
177   monflg = PETSC_TRUE;
178   subviewflg = PETSC_FALSE;
179   ierr   = PetscOptionsBool("-snes_nasm_sub_view","Print detailed information for every processor when using -snes_view","",subviewflg,&subviewflg,&flg);CHKERRQ(ierr);
180   if (flg) {
181     nasm->same_local_solves = PETSC_FALSE;
182     if (!subviewflg) {
183       nasm->same_local_solves = PETSC_TRUE;
184     }
185   }
186   ierr   = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr);
187   ierr   = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr);
188   if (flg) {
189     ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr);
190     ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr);
191   }
192   ierr = PetscOptionsTail();CHKERRQ(ierr);
193   PetscFunctionReturn(0);
194 }
195 
196 #undef __FUNCT__
197 #define __FUNCT__ "SNESView_NASM"
198 PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
199 {
200   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
201   PetscErrorCode ierr;
202   PetscMPIInt    rank,size;
203   PetscInt       i,j,N,bsz;
204   PetscBool      iascii,isstring;
205   PetscViewer    sviewer;
206   MPI_Comm       comm;
207 
208   PetscFunctionBegin;
209   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
210   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
211   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
212   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
213   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
214   ierr = MPI_Allreduce(&nasm->n,&N,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
215   if (iascii) {
216     ierr = PetscViewerASCIIPrintf(viewer, "  Nonlinear Additive Schwarz: total subdomain blocks = %D\n",N);CHKERRQ(ierr);
217     if (nasm->same_local_solves) {
218       if (nasm->subsnes) {
219         ierr = PetscViewerASCIIPrintf(viewer,"  Local solve is the same for all blocks:\n");CHKERRQ(ierr);
220         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
221         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
222         if (!rank) {
223           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
224           ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);
225           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
226         }
227         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
228         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
229       }
230     } else {
231       /* print the solver on each block */
232       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
233       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr);
234       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
235       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
236       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following SNES objects:\n");CHKERRQ(ierr);
237       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
238       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
239       for (j=0; j<size; j++) {
240         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
241         if (rank == j) {
242           for (i=0; i<nasm->n; i++) {
243             ierr = VecGetLocalSize(nasm->x[i],&bsz);CHKERRQ(ierr);
244             ierr = PetscViewerASCIIPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr);
245             ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr);
246             ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
247           }
248         }
249         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
250         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
251       }
252       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
253     }
254   } else if (isstring) {
255     ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr);
256     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
257     if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);}
258     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
259   }
260   PetscFunctionReturn(0);
261 }
262 
263 #undef __FUNCT__
264 #define __FUNCT__ "SNESNASMSetSubdomains"
265 /*@
266    SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems.
267 
268    Not Collective
269 
270    Input Parameters:
271 +  SNES - the SNES context
272 .  n - the number of local subdomains
273 .  subsnes - solvers defined on the local subdomains
274 .  iscatter - scatters into the nonoverlapping portions of the local subdomains
275 .  oscatter - scatters into the overlapping portions of the local subdomains
276 -  gscatter - scatters into the (ghosted) local vector of the local subdomain
277 
278    Level: intermediate
279 
280 .keywords: SNES, NASM
281 
282 .seealso: SNESNASM, SNESNASMGetSubdomains()
283 @*/
284 PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
285 {
286   PetscErrorCode ierr;
287   PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*);
288 
289   PetscFunctionBegin;
290   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f);CHKERRQ(ierr);
291   ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);
292   PetscFunctionReturn(0);
293 }
294 
295 #undef __FUNCT__
296 #define __FUNCT__ "SNESNASMSetSubdomains_NASM"
297 PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
298 {
299   PetscInt       i;
300   PetscErrorCode ierr;
301   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
302 
303   PetscFunctionBegin;
304   if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
305 
306   /* tear down the previously set things */
307   ierr = SNESReset(snes);CHKERRQ(ierr);
308 
309   nasm->n = n;
310   if (oscatter) {
311     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);}
312   }
313   if (iscatter) {
314     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);}
315   }
316   if (gscatter) {
317     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);}
318   }
319   if (oscatter) {
320     ierr = PetscMalloc(n*sizeof(IS),&nasm->oscatter);CHKERRQ(ierr);
321     for (i=0; i<n; i++) {
322       nasm->oscatter[i] = oscatter[i];
323     }
324   }
325   if (iscatter) {
326     ierr = PetscMalloc(n*sizeof(IS),&nasm->iscatter);CHKERRQ(ierr);
327     for (i=0; i<n; i++) {
328       nasm->iscatter[i] = iscatter[i];
329     }
330   }
331   if (gscatter) {
332     ierr = PetscMalloc(n*sizeof(IS),&nasm->gscatter);CHKERRQ(ierr);
333     for (i=0; i<n; i++) {
334       nasm->gscatter[i] = gscatter[i];
335     }
336   }
337 
338   if (subsnes) {
339     ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr);
340     for (i=0; i<n; i++) {
341       nasm->subsnes[i] = subsnes[i];
342     }
343     nasm->same_local_solves = PETSC_FALSE;
344   }
345   PetscFunctionReturn(0);
346 }
347 
348 #undef __FUNCT__
349 #define __FUNCT__ "SNESNASMGetSubdomains"
350 /*@
351    SNESNASMGetSubdomains - Get the local subdomain context.
352 
353    Not Collective
354 
355    Input Parameters:
356 .  SNES - the SNES context
357 
358    Output Parameters:
359 +  n - the number of local subdomains
360 .  subsnes - solvers defined on the local subdomains
361 .  iscatter - scatters into the nonoverlapping portions of the local subdomains
362 .  oscatter - scatters into the overlapping portions of the local subdomains
363 -  gscatter - scatters into the (ghosted) local vector of the local subdomain
364 
365    Level: intermediate
366 
367 .keywords: SNES, NASM
368 
369 .seealso: SNESNASM, SNESNASMSetSubdomains()
370 @*/
371 PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
372 {
373   PetscErrorCode ierr;
374   PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**);
375 
376   PetscFunctionBegin;
377   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f);CHKERRQ(ierr);
378   ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);
379   PetscFunctionReturn(0);
380 }
381 
382 #undef __FUNCT__
383 #define __FUNCT__ "SNESNASMGetSubdomains_NASM"
384 PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
385 {
386   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
387 
388   PetscFunctionBegin;
389   if (n) *n = nasm->n;
390   if (oscatter) *oscatter = nasm->oscatter;
391   if (iscatter) *iscatter = nasm->iscatter;
392   if (gscatter) *gscatter = nasm->gscatter;
393   if (subsnes)  {
394     *subsnes  = nasm->subsnes;
395     nasm->same_local_solves = PETSC_FALSE;
396   }
397   PetscFunctionReturn(0);
398 }
399 
400 #undef __FUNCT__
401 #define __FUNCT__ "SNESNASMGetSubdomainVecs"
402 /*@
403    SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors
404 
405    Not Collective
406 
407    Input Parameters:
408 .  SNES - the SNES context
409 
410    Output Parameters:
411 +  n - the number of local subdomains
412 .  x - The subdomain solution vector
413 .  y - The subdomain step vector
414 .  b - The subdomain RHS vector
415 -  xl - The subdomain local vectors (ghosted)
416 
417    Level: developer
418 
419 .keywords: SNES, NASM
420 
421 .seealso: SNESNASM, SNESNASMGetSubdomains()
422 @*/
423 PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl)
424 {
425   PetscErrorCode ierr;
426   PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**);
427 
428   PetscFunctionBegin;
429   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f);CHKERRQ(ierr);
430   ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);
431   PetscFunctionReturn(0);
432 }
433 
434 #undef __FUNCT__
435 #define __FUNCT__ "SNESNASMGetSubdomainVecs_NASM"
436 PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl)
437 {
438   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
439 
440   PetscFunctionBegin;
441   if (n)  *n  = nasm->n;
442   if (x)  *x  = nasm->x;
443   if (y)  *y  = nasm->y;
444   if (b)  *b  = nasm->b;
445   if (xl) *xl = nasm->xl;
446   PetscFunctionReturn(0);
447 }
448 
449 #undef __FUNCT__
450 #define __FUNCT__ "SNESNASMSetComputeFinalJacobian"
451 /*@
452    SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence
453 
454    Collective on SNES
455 
456    Input Parameters:
457 +  SNES - the SNES context
458 -  flg - indication of whether to compute the jacobians or not
459 
460    Level: developer
461 
462    Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian
463    is needed at each linear iteration.
464 
465 .keywords: SNES, NASM, ASPIN
466 
467 .seealso: SNESNASM, SNESNASMGetSubdomains()
468 @*/
469 PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg)
470 {
471   PetscErrorCode (*f)(SNES,PetscBool);
472   PetscErrorCode ierr;
473 
474   PetscFunctionBegin;
475   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f);CHKERRQ(ierr);
476   ierr = (f)(snes,flg);CHKERRQ(ierr);
477   PetscFunctionReturn(0);
478 }
479 
480 #undef __FUNCT__
481 #define __FUNCT__ "SNESNASMSetComputeFinalJacobian_NASM"
482 PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg)
483 {
484   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
485 
486   PetscFunctionBegin;
487   nasm->finaljacobian = flg;
488   PetscFunctionReturn(0);
489 }
490 
491 #undef __FUNCT__
492 #define __FUNCT__ "SNESNASMSolveLocal_Private"
493 PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X)
494 {
495   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
496   SNES           subsnes;
497   PetscInt       i;
498   PetscErrorCode ierr;
499   Vec            Xlloc,Xl,Bl,Yl;
500   VecScatter     iscat,oscat,gscat;
501   DM             dm,subdm;
502 
503   PetscFunctionBegin;
504   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
505   ierr = VecSet(Y,0);CHKERRQ(ierr);
506   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
507   for (i=0; i<nasm->n; i++) {
508     /* scatter the solution to the local solution */
509     Xlloc = nasm->xl[i];
510     gscat   = nasm->gscatter[i];
511     oscat   = nasm->oscatter[i];
512     ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
513     if (B) {
514       /* scatter the RHS to the local RHS */
515       Bl   = nasm->b[i];
516       ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
517     }
518   }
519   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
520 
521 
522   if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
523   for (i=0; i<nasm->n; i++) {
524     Xl    = nasm->x[i];
525     Xlloc = nasm->xl[i];
526     Yl    = nasm->y[i];
527     subsnes = nasm->subsnes[i];
528     ierr    = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
529     iscat   = nasm->iscatter[i];
530     oscat   = nasm->oscatter[i];
531     gscat   = nasm->gscatter[i];
532     ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
533     if (B) {
534       Bl   = nasm->b[i];
535       ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
536     } else Bl = NULL;
537     ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
538     ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
539     ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
540     ierr = VecCopy(Xl,Yl);CHKERRQ(ierr);
541     ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr);
542     ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr);
543     if (nasm->type == PC_ASM_BASIC) {
544       ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
545     } else if (nasm->type == PC_ASM_RESTRICT) {
546       ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
547     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
548   }
549   if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
550   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
551   for (i=0; i<nasm->n; i++) {
552     Yl    = nasm->y[i];
553     iscat   = nasm->iscatter[i];
554     oscat   = nasm->oscatter[i];
555     if (nasm->type == PC_ASM_BASIC) {
556       ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
557     } else if (nasm->type == PC_ASM_RESTRICT) {
558       ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
559     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
560   }
561   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
562   ierr = VecAXPY(X,1.0,Y);CHKERRQ(ierr);
563   PetscFunctionReturn(0);
564 }
565 
566 #undef __FUNCT__
567 #define __FUNCT__ "SNESNASMComputeFinalJacobian_Private"
568 PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec X)
569 {
570   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
571   SNES           subsnes;
572   PetscInt       i,lag = 1;
573   PetscErrorCode ierr;
574   Vec            Xlloc,Xl,Fl,F;
575   VecScatter     oscat,gscat;
576   DM             dm,subdm;
577   MatStructure   flg = DIFFERENT_NONZERO_PATTERN;
578   PetscFunctionBegin;
579   F = snes->vec_func;
580   if (snes->normtype == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);}
581   ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
582   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
583   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
584   for (i=0; i<nasm->n; i++) {
585     Xlloc = nasm->xl[i];
586     Fl    = nasm->subsnes[i]->vec_func;
587     gscat = nasm->gscatter[i];
588     oscat = nasm->oscatter[i];
589     ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
590     ierr = VecScatterBegin(oscat,F,Fl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
591   }
592   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
593   for (i=0; i<nasm->n; i++) {
594     Fl      = nasm->subsnes[i]->vec_func;
595     Xl      = nasm->x[i];
596     Xlloc   = nasm->xl[i];
597     subsnes = nasm->subsnes[i];
598     oscat   = nasm->oscatter[i];
599     gscat   = nasm->gscatter[i];
600     ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
601     ierr = VecScatterEnd(oscat,F,Fl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
602     ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
603     ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
604     ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
605     ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
606 
607     if (subsnes->lagjacobian == -1)    subsnes->lagjacobian = -2;
608     else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian;
609     ierr = SNESComputeJacobian(subsnes,Xl,&subsnes->jacobian,&subsnes->jacobian_pre,&flg);CHKERRQ(ierr);
610     if (lag > 1) subsnes->lagjacobian = lag;
611     ierr = KSPSetOperators(subsnes->ksp,subsnes->jacobian,subsnes->jacobian_pre,flg);CHKERRQ(ierr);
612   }
613   PetscFunctionReturn(0);
614 }
615 
616 #undef __FUNCT__
617 #define __FUNCT__ "SNESSolve_NASM"
618 PetscErrorCode SNESSolve_NASM(SNES snes)
619 {
620   Vec            F;
621   Vec            X;
622   Vec            B;
623   Vec            Y;
624   PetscInt       i;
625   PetscReal      fnorm = 0.0;
626   PetscErrorCode ierr;
627   SNESNormType   normtype;
628   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
629 
630   PetscFunctionBegin;
631   X = snes->vec_sol;
632   Y = snes->vec_sol_update;
633   F = snes->vec_func;
634   B = snes->vec_rhs;
635 
636   ierr         = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
637   snes->iter   = 0;
638   snes->norm   = 0.;
639   ierr         = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
640   snes->reason = SNES_CONVERGED_ITERATING;
641   ierr         = SNESGetNormType(snes, &normtype);CHKERRQ(ierr);
642   if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
643     /* compute the initial function and preconditioned update delX */
644     if (!snes->vec_func_init_set) {
645       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
646       if (snes->domainerror) {
647         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
648         PetscFunctionReturn(0);
649       }
650     } else snes->vec_func_init_set = PETSC_FALSE;
651 
652     /* convergence test */
653     if (!snes->norm_init_set) {
654       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
655       if (PetscIsInfOrNanReal(fnorm)) {
656         snes->reason = SNES_DIVERGED_FNORM_NAN;
657         PetscFunctionReturn(0);
658       }
659     } else {
660       fnorm               = snes->norm_init;
661       snes->norm_init_set = PETSC_FALSE;
662     }
663     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
664     snes->iter = 0;
665     snes->norm = fnorm;
666     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
667     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
668     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
669 
670     /* set parameter for default relative tolerance convergence test */
671     snes->ttol = fnorm*snes->rtol;
672 
673     /* test convergence */
674     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
675     if (snes->reason) PetscFunctionReturn(0);
676   } else {
677     ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
678     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
679     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
680   }
681 
682   /* Call general purpose update function */
683   if (snes->ops->update) {
684     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
685   }
686 
687   for (i = 0; i < snes->max_its; i++) {
688     ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr);
689     if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) {
690       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
691       if (snes->domainerror) {
692         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
693         break;
694       }
695       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
696       if (PetscIsInfOrNanReal(fnorm)) {
697         snes->reason = SNES_DIVERGED_FNORM_NAN;
698         break;
699       }
700     }
701     /* Monitor convergence */
702     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
703     snes->iter = i+1;
704     snes->norm = fnorm;
705     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
706     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
707     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
708     /* Test for convergence */
709     if (normtype == SNES_NORM_FUNCTION) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
710     if (snes->reason) break;
711     /* Call general purpose update function */
712     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
713   }
714   if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);}
715   if (normtype == SNES_NORM_FUNCTION) {
716     if (i == snes->max_its) {
717       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
718       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
719     }
720   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */
721   PetscFunctionReturn(0);
722 }
723 
724 /*MC
725   SNESNASM - Nonlinear Additive Schwartz
726 
727    Options Database:
728 +  -snes_nasm_log - enable logging events for the communication and solve stages
729 .  -snes_nasm_type <basic,restrict> - type of subdomain update used
730 .  -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate
731 .  -sub_snes_ - options prefix of the subdomain nonlinear solves
732 .  -sub_ksp_ - options prefix of the subdomain Krylov solver
733 -  -sub_pc_ - options prefix of the subdomain preconditioner
734 
735    Level: advanced
736 
737 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
738 M*/
739 
740 #undef __FUNCT__
741 #define __FUNCT__ "SNESCreate_NASM"
742 PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes)
743 {
744   SNES_NASM      *nasm;
745   PetscErrorCode ierr;
746 
747   PetscFunctionBegin;
748   ierr       = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr);
749   snes->data = (void*)nasm;
750 
751   nasm->n        = PETSC_DECIDE;
752   nasm->subsnes  = 0;
753   nasm->x        = 0;
754   nasm->xl       = 0;
755   nasm->y        = 0;
756   nasm->b        = 0;
757   nasm->oscatter = 0;
758   nasm->iscatter = 0;
759   nasm->gscatter = 0;
760 
761   nasm->type = PC_ASM_BASIC;
762   nasm->finaljacobian = PETSC_FALSE;
763   nasm->same_local_solves = PETSC_TRUE;
764 
765   snes->ops->destroy        = SNESDestroy_NASM;
766   snes->ops->setup          = SNESSetUp_NASM;
767   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
768   snes->ops->view           = SNESView_NASM;
769   snes->ops->solve          = SNESSolve_NASM;
770   snes->ops->reset          = SNESReset_NASM;
771 
772   snes->usesksp = PETSC_FALSE;
773   snes->usespc  = PETSC_FALSE;
774 
775   nasm->eventrestrictinterp = 0;
776   nasm->eventsubsolve       = 0;
777 
778   if (!snes->tolerancesset) {
779     snes->max_its   = 10000;
780     snes->max_funcs = 10000;
781   }
782 
783   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr);
784   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr);
785   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr);
786   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr);
787   PetscFunctionReturn(0);
788 }
789 
790