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