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