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