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