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