xref: /petsc/src/snes/impls/nasm/nasm.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
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 
17   /* logging events */
18   PetscLogEvent eventrestrictinterp;
19   PetscLogEvent eventsubsolve;
20 } SNES_NASM;
21 
22 const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0};
23 
24 #undef __FUNCT__
25 #define __FUNCT__ "SNESReset_NASM"
26 PetscErrorCode SNESReset_NASM(SNES snes)
27 {
28   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
29   PetscErrorCode ierr;
30   PetscInt       i;
31 
32   PetscFunctionBegin;
33   for (i=0; i<nasm->n; i++) {
34     if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); }
35     if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); }
36     if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); }
37     if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); }
38 
39     if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); }
40     if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); }
41     if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); }
42     if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); }
43   }
44 
45   if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);}
46   if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);}
47   if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);}
48   if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);}
49 
50   if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);}
51   if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);}
52   if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);}
53   if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);}
54 
55   nasm->eventrestrictinterp = 0;
56   nasm->eventsubsolve = 0;
57   PetscFunctionReturn(0);
58 }
59 
60 #undef __FUNCT__
61 #define __FUNCT__ "SNESDestroy_NASM"
62 PetscErrorCode SNESDestroy_NASM(SNES snes)
63 {
64   PetscErrorCode ierr;
65 
66   PetscFunctionBegin;
67   ierr = SNESReset_NASM(snes);CHKERRQ(ierr);
68   ierr = PetscFree(snes->data);CHKERRQ(ierr);
69   PetscFunctionReturn(0);
70 }
71 
72 #undef __FUNCT__
73 #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private"
74 PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx)
75 {
76   PetscErrorCode ierr;
77   Vec            bcs = (Vec)ctx;
78 
79   PetscFunctionBegin;
80   ierr = VecCopy(bcs,l);CHKERRQ(ierr);
81   PetscFunctionReturn(0);
82 }
83 
84 #undef __FUNCT__
85 #define __FUNCT__ "SNESSetUp_NASM"
86 PetscErrorCode SNESSetUp_NASM(SNES snes)
87 {
88   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
89   PetscErrorCode ierr;
90   DM             dm;
91   DM             *subdms;
92   PetscInt       i;
93   const char     *optionsprefix;
94   Vec            F;
95 
96   PetscFunctionBegin;
97   if (!nasm->subsnes) {
98     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
99     if (dm) {
100       nasm->usesdm = PETSC_TRUE;
101       ierr         = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr);
102       if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
103       ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr);
104 
105       ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
106       ierr = PetscMalloc(nasm->n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr);
107 
108       for (i=0; i<nasm->n; i++) {
109         ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr);
110         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr);
111         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr);
112         ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr);
113         ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr);
114         ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr);
115       }
116       ierr = PetscFree(subdms);CHKERRQ(ierr);
117     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!");
118   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!");
119   /* allocate the global vectors */
120   ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->x);CHKERRQ(ierr);
121   ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->xl);CHKERRQ(ierr);
122   ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->y);CHKERRQ(ierr);
123   ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->b);CHKERRQ(ierr);
124 
125   for (i=0; i<nasm->n; i++) {
126     DM subdm;
127     ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr);
128     ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);
129     ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);
130     ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);
131     ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr);
132     ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr);
133     ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr);
134   }
135   PetscFunctionReturn(0);
136 }
137 
138 #undef __FUNCT__
139 #define __FUNCT__ "SNESSetFromOptions_NASM"
140 PetscErrorCode SNESSetFromOptions_NASM(SNES snes)
141 {
142   PetscErrorCode    ierr;
143   PCASMType         asmtype;
144   PetscBool         flg,monflg;
145   SNES_NASM         *nasm = (SNES_NASM*)snes->data;
146 
147   PetscFunctionBegin;
148   ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr);
149   ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
150   if (flg) nasm->type = asmtype;
151   flg    = PETSC_FALSE;
152   monflg = PETSC_TRUE;
153   ierr   = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr);
154   if (flg) {
155     ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr);
156     ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr);
157   }
158   ierr = PetscOptionsTail();CHKERRQ(ierr);
159   PetscFunctionReturn(0);
160 }
161 
162 #undef __FUNCT__
163 #define __FUNCT__ "SNESView_NASM"
164 PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
165 {
166   PetscFunctionBegin;
167   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
168   PetscErrorCode ierr;
169   PetscMPIInt    rank;
170   PetscInt       i,N;
171   PetscBool      iascii,isstring;
172   PetscViewer    sviewer;
173   MPI_Comm       comm;
174 
175   PetscFunctionBegin;
176   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
177   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
178   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
179   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
180   ierr = MPI_Reduce(&nasm->n,&N,1,MPIU_INT,MPIU_SUM,0,comm);CHKERRQ(ierr);
181   if (iascii) {
182     ierr = PetscViewerASCIIPrintf(viewer, "  Nonlinear Additive Schwarz: total subdomain blocks = %D\n",N);CHKERRQ(ierr);
183     ierr = PetscViewerASCIIPrintf(viewer,"  Nonlinear Additive Schwarz: restriction/interpolation type - %s\n",SNESNASMTypes[nasm->type]);CHKERRQ(ierr);
184     ierr = PetscViewerASCIIPrintf(viewer,"  Nonlinear Additive Schwarz: subSNES iterations: %D subKSP iterations: %D\n",SNESNASMTypes[nasm->type]);CHKERRQ(ierr);
185     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
186     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr);
187     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
188     ierr = PetscViewerASCIIPrintf(viewer,"  Local SNES objects:\n");CHKERRQ(ierr);
189     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
190     if (!rank) {
191       for (i=0; i<nasm->n; i++) {
192         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
193         ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr);
194         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
195       }
196     }
197     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
198     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
199     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
200   } else if (isstring) {
201     ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr);
202     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
203     if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);}
204     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
205   }
206   PetscFunctionReturn(0);
207 }
208 
209 #undef __FUNCT__
210 #define __FUNCT__ "SNESNASMSetSubdomains"
211 PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
212 {
213   PetscErrorCode ierr;
214   PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*);
215 
216   PetscFunctionBegin;
217   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr);
218   ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);
219   PetscFunctionReturn(0);
220 }
221 
222 EXTERN_C_BEGIN
223 #undef __FUNCT__
224 #define __FUNCT__ "SNESNASMSetSubdomains_NASM"
225 PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
226 {
227   PetscInt       i;
228   PetscErrorCode ierr;
229   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
230 
231   PetscFunctionBegin;
232   if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
233 
234   /* tear down the previously set things */
235   ierr = SNESReset(snes);CHKERRQ(ierr);
236 
237   nasm->n = n;
238   if (oscatter) {
239     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);}
240   }
241   if (iscatter) {
242     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);}
243   }
244   if (gscatter) {
245     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);}
246   }
247   if (oscatter) {
248     ierr = PetscMalloc(n*sizeof(IS),&nasm->oscatter);CHKERRQ(ierr);
249     for (i=0; i<n; i++) {
250       nasm->oscatter[i] = oscatter[i];
251     }
252   }
253   if (iscatter) {
254     ierr = PetscMalloc(n*sizeof(IS),&nasm->iscatter);CHKERRQ(ierr);
255     for (i=0; i<n; i++) {
256       nasm->iscatter[i] = iscatter[i];
257     }
258   }
259   if (gscatter) {
260     ierr = PetscMalloc(n*sizeof(IS),&nasm->gscatter);CHKERRQ(ierr);
261     for (i=0; i<n; i++) {
262       nasm->gscatter[i] = gscatter[i];
263     }
264   }
265 
266   if (subsnes) {
267     ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr);
268     for (i=0; i<n; i++) {
269       nasm->subsnes[i] = subsnes[i];
270     }
271   }
272   PetscFunctionReturn(0);
273 }
274 EXTERN_C_END
275 
276 #undef __FUNCT__
277 #define __FUNCT__ "SNESNASMSolveLocal_Private"
278 PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X)
279 {
280   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
281   SNES           subsnes;
282   PetscInt       i;
283   PetscErrorCode ierr;
284   Vec            Xlloc,Xl,Bl,Yl;
285   VecScatter     iscat,oscat,gscat;
286   DM             dm,subdm;
287 
288   PetscFunctionBegin;
289   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
290   ierr = VecSet(Y,0);CHKERRQ(ierr);
291 
292   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
293   for (i=0; i<nasm->n; i++) {
294     /* scatter the solution to the local solution */
295     Xlloc = nasm->xl[i];
296     gscat   = nasm->gscatter[i];
297     oscat   = nasm->oscatter[i];
298     ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
299     if (B) {
300       /* scatter the RHS to the local RHS */
301       Bl   = nasm->b[i];
302       ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
303     }
304   }
305   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
306 
307   for (i=0; i<nasm->n; i++) {
308     Xlloc = nasm->xl[i];
309     gscat   = nasm->gscatter[i];
310     oscat   = nasm->oscatter[i];
311     ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
312     if (B) {
313       Bl   = nasm->b[i];
314       ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
315     }
316   }
317 
318   if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
319   for (i=0; i<nasm->n; i++) {
320     Xl    = nasm->x[i];
321     Xlloc = nasm->xl[i];
322     Yl    = nasm->y[i];
323     subsnes = nasm->subsnes[i];
324     ierr    = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
325     iscat   = nasm->iscatter[i];
326     oscat   = nasm->oscatter[i];
327     gscat   = nasm->gscatter[i];
328     ierr    = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
329     if (B) {
330       Bl = nasm->b[i];
331     } else {
332       Bl = NULL;
333     }
334     ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
335     ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
336     ierr = VecCopy(Xl,Yl);CHKERRQ(ierr);
337     ierr = SNESSolve(subsnes,Bl,Yl);CHKERRQ(ierr);
338     ierr = VecAXPY(Yl,-1.0,Xl);CHKERRQ(ierr);
339   }
340   if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
341 
342   ierr = MPI_Barrier(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
343 
344   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
345   for (i=0; i<nasm->n; i++) {
346     Yl    = nasm->y[i];
347     iscat   = nasm->iscatter[i];
348     oscat   = nasm->oscatter[i];
349     if (nasm->type == PC_ASM_BASIC) {
350       ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
351     } else if (nasm->type == PC_ASM_RESTRICT) {
352       ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
353     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
354   }
355 
356   for (i=0; i<nasm->n; i++) {
357     Yl    = nasm->y[i];
358     iscat   = nasm->iscatter[i];
359     oscat   = nasm->oscatter[i];
360     if (nasm->type == PC_ASM_BASIC) {
361       ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
362     } else if (nasm->type == PC_ASM_RESTRICT) {
363       ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
364     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
365   }
366   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
367 
368   ierr = MPI_Barrier(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
369 
370   ierr = VecAXPY(X,1.0,Y);CHKERRQ(ierr);
371   PetscFunctionReturn(0);
372 }
373 
374 #undef __FUNCT__
375 #define __FUNCT__ "SNESSolve_NASM"
376 PetscErrorCode SNESSolve_NASM(SNES snes)
377 {
378   Vec            F;
379   Vec            X;
380   Vec            B;
381   Vec            Y;
382   PetscInt       i;
383   PetscReal      fnorm;
384   PetscErrorCode ierr;
385   SNESNormType   normtype;
386 
387   PetscFunctionBegin;
388   X = snes->vec_sol;
389   Y = snes->vec_sol_update;
390   F = snes->vec_func;
391   B = snes->vec_rhs;
392 
393   ierr         = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
394   snes->iter   = 0;
395   snes->norm   = 0.;
396   ierr         = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
397   snes->reason = SNES_CONVERGED_ITERATING;
398   ierr         = SNESGetNormType(snes, &normtype);CHKERRQ(ierr);
399   if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
400     /* compute the initial function and preconditioned update delX */
401     if (!snes->vec_func_init_set) {
402       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
403       if (snes->domainerror) {
404         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
405         PetscFunctionReturn(0);
406       }
407     } else snes->vec_func_init_set = PETSC_FALSE;
408 
409     /* convergence test */
410     if (!snes->norm_init_set) {
411       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
412       if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
413     } else {
414       fnorm               = snes->norm_init;
415       snes->norm_init_set = PETSC_FALSE;
416     }
417     ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
418     snes->iter = 0;
419     snes->norm = fnorm;
420     ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
421     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
422     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
423 
424     /* set parameter for default relative tolerance convergence test */
425     snes->ttol = fnorm*snes->rtol;
426 
427     /* test convergence */
428     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
429     if (snes->reason) PetscFunctionReturn(0);
430   } else {
431     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
432     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
433     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
434   }
435 
436   /* Call general purpose update function */
437   if (snes->ops->update) {
438     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
439   }
440 
441   for (i = 0; i < snes->max_its; i++) {
442     ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr);
443     if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) {
444       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
445       if (snes->domainerror) {
446         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
447         PetscFunctionReturn(0);
448       }
449       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
450       if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
451     }
452     /* Monitor convergence */
453     ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
454     snes->iter = i+1;
455     snes->norm = fnorm;
456     ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
457     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
458     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
459     /* Test for convergence */
460     if (normtype == SNES_NORM_FUNCTION) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
461     if (snes->reason) PetscFunctionReturn(0);
462     /* Call general purpose update function */
463     if (snes->ops->update) {
464       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
465     }
466   }
467   if (normtype == SNES_NORM_FUNCTION) {
468     if (i == snes->max_its) {
469       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
470       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
471     }
472   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */
473   PetscFunctionReturn(0);
474 }
475 
476 /*MC
477   SNESNASM - Nonlinear Additive Schwartz
478 
479    Level: advanced
480 
481 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
482 M*/
483 
484 EXTERN_C_BEGIN
485 #undef __FUNCT__
486 #define __FUNCT__ "SNESCreate_NASM"
487 PetscErrorCode SNESCreate_NASM(SNES snes)
488 {
489   SNES_NASM      *nasm;
490   PetscErrorCode ierr;
491 
492   PetscFunctionBegin;
493   ierr       = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr);
494   snes->data = (void*)nasm;
495 
496   nasm->n        = PETSC_DECIDE;
497   nasm->subsnes  = 0;
498   nasm->x        = 0;
499   nasm->xl       = 0;
500   nasm->y        = 0;
501   nasm->b        = 0;
502   nasm->oscatter = 0;
503   nasm->iscatter = 0;
504   nasm->gscatter = 0;
505 
506   nasm->type = PC_ASM_BASIC;
507 
508   snes->ops->destroy        = SNESDestroy_NASM;
509   snes->ops->setup          = SNESSetUp_NASM;
510   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
511   snes->ops->view           = SNESView_NASM;
512   snes->ops->solve          = SNESSolve_NASM;
513   snes->ops->reset          = SNESReset_NASM;
514 
515   snes->usesksp = PETSC_FALSE;
516   snes->usespc  = PETSC_FALSE;
517 
518   nasm->eventrestrictinterp = 0;
519   nasm->eventsubsolve       = 0;
520 
521   if (!snes->tolerancesset) {
522     snes->max_its   = 10000;
523     snes->max_funcs = 10000;
524   }
525 
526   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNASMSetSubdomains_C","SNESNASMSetSubdomains_NASM",
527                                            SNESNASMSetSubdomains_NASM);CHKERRQ(ierr);
528   PetscFunctionReturn(0);
529 }
530 EXTERN_C_END
531