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