xref: /petsc/src/snes/impls/nasm/nasm.c (revision eaedb03373539b57b5e56e72e3c75b64208da092)
1*eaedb033SPeter Brune #include <petsc-private/snesimpl.h>             /*I   "petscsnes.h"   I*/
2*eaedb033SPeter Brune #include <petscdmshell.h>
3*eaedb033SPeter Brune 
4*eaedb033SPeter Brune typedef struct {
5*eaedb033SPeter Brune   PetscInt   n;                   /* local subdomains */
6*eaedb033SPeter Brune   SNES       *subsnes;            /* nonlinear solvers for each subdomain */
7*eaedb033SPeter Brune 
8*eaedb033SPeter Brune   Vec        *r;                  /* function vectors */
9*eaedb033SPeter Brune   Vec        *x;                  /* solution vectors */
10*eaedb033SPeter Brune   Vec        *b;                  /* rhs vectors */
11*eaedb033SPeter Brune 
12*eaedb033SPeter Brune   IS         *ois;
13*eaedb033SPeter Brune   IS         *iis;
14*eaedb033SPeter Brune   PetscBool  usesdm;               /* use the outer DM's communication facilities rather than ISes */
15*eaedb033SPeter Brune } SNES_NASM;
16*eaedb033SPeter Brune 
17*eaedb033SPeter Brune #undef __FUNCT__
18*eaedb033SPeter Brune #define __FUNCT__ "SNESReset_NASM"
19*eaedb033SPeter Brune PetscErrorCode SNESReset_NASM(SNES snes)
20*eaedb033SPeter Brune {
21*eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM *)snes->data;
22*eaedb033SPeter Brune   PetscErrorCode ierr;
23*eaedb033SPeter Brune   PetscInt       i;
24*eaedb033SPeter Brune   PetscFunctionBegin;
25*eaedb033SPeter Brune   for (i=0;i<nasm->n;i++) {
26*eaedb033SPeter Brune     if (nasm->x) ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr);
27*eaedb033SPeter Brune     if (nasm->r) ierr = VecDestroy(&nasm->r[i]);CHKERRQ(ierr);
28*eaedb033SPeter Brune     if (nasm->b) ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr);
29*eaedb033SPeter Brune 
30*eaedb033SPeter Brune     if (nasm->subsnes) ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr);
31*eaedb033SPeter Brune     if (nasm->ois) ierr = ISDestroy(&nasm->ois[i]);CHKERRQ(ierr);
32*eaedb033SPeter Brune     if (nasm->iis) ierr = ISDestroy(&nasm->iis[i]);CHKERRQ(ierr);
33*eaedb033SPeter Brune   }
34*eaedb033SPeter Brune   PetscFunctionReturn(0);
35*eaedb033SPeter Brune }
36*eaedb033SPeter Brune 
37*eaedb033SPeter Brune #undef __FUNCT__
38*eaedb033SPeter Brune #define __FUNCT__ "SNESDestroy_NASM"
39*eaedb033SPeter Brune PetscErrorCode SNESDestroy_NASM(SNES snes)
40*eaedb033SPeter Brune {
41*eaedb033SPeter Brune   PetscErrorCode ierr;
42*eaedb033SPeter Brune   PetscFunctionBegin;
43*eaedb033SPeter Brune   ierr = SNESReset_NASM(snes);CHKERRQ(ierr);
44*eaedb033SPeter Brune   ierr = PetscFree(snes->data);
45*eaedb033SPeter Brune   PetscFunctionReturn(0);
46*eaedb033SPeter Brune }
47*eaedb033SPeter Brune 
48*eaedb033SPeter Brune #undef __FUNCT__
49*eaedb033SPeter Brune #define __FUNCT__ "SNESSetUp_NASM"
50*eaedb033SPeter Brune PetscErrorCode SNESSetUp_NASM(SNES snes)
51*eaedb033SPeter Brune {
52*eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM *)snes->data;
53*eaedb033SPeter Brune   PetscErrorCode ierr;
54*eaedb033SPeter Brune   DM             dm;
55*eaedb033SPeter Brune   DM             subdm;
56*eaedb033SPeter Brune   PetscErrorCode (*f)(SNES,Vec,Vec,void*);
57*eaedb033SPeter Brune 
58*eaedb033SPeter Brune   Mat            A;
59*eaedb033SPeter Brune   PetscErrorCode (*fj)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
60*eaedb033SPeter Brune 
61*eaedb033SPeter Brune   PetscInt       n;
62*eaedb033SPeter Brune   void           *ctx;
63*eaedb033SPeter Brune   const char     *optionsprefix;
64*eaedb033SPeter Brune 
65*eaedb033SPeter Brune   PetscFunctionBegin;
66*eaedb033SPeter Brune 
67*eaedb033SPeter Brune   if (!nasm->subsnes) {
68*eaedb033SPeter Brune     if (snes->dm) {
69*eaedb033SPeter Brune       ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
70*eaedb033SPeter Brune       nasm->n      = 1;
71*eaedb033SPeter Brune       nasm->usesdm = PETSC_TRUE;
72*eaedb033SPeter Brune       /* create a single solver */
73*eaedb033SPeter Brune       ierr = PetscMalloc(sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr);
74*eaedb033SPeter Brune       ierr = PetscMalloc(sizeof(Vec),&nasm->r);CHKERRQ(ierr);
75*eaedb033SPeter Brune       ierr = PetscMalloc(sizeof(Vec),&nasm->x);CHKERRQ(ierr);
76*eaedb033SPeter Brune       ierr = PetscMalloc(sizeof(Vec),&nasm->b);CHKERRQ(ierr);
77*eaedb033SPeter Brune       ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[0]);CHKERRQ(ierr);
78*eaedb033SPeter Brune 
79*eaedb033SPeter Brune       ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
80*eaedb033SPeter Brune       ierr = SNESAppendOptionsPrefix(nasm->subsnes[0],optionsprefix);CHKERRQ(ierr);
81*eaedb033SPeter Brune       ierr = SNESAppendOptionsPrefix(nasm->subsnes[0],"sub_");CHKERRQ(ierr);
82*eaedb033SPeter Brune 
83*eaedb033SPeter Brune       ierr = SNESGetDM(nasm->subsnes[0],&subdm);CHKERRQ(ierr);
84*eaedb033SPeter Brune 
85*eaedb033SPeter Brune       /* set up the local function */
86*eaedb033SPeter Brune       ierr = DMSNESGetBlockFunction(dm,&f,&ctx);CHKERRQ(ierr);
87*eaedb033SPeter Brune       ierr = DMCreateLocalVector(dm,&nasm->r[0]);CHKERRQ(ierr);
88*eaedb033SPeter Brune       ierr = DMShellSetGlobalVector(dm,nasm->r[0]);CHKERRQ(ierr);
89*eaedb033SPeter Brune       ierr = DMCreateLocalVector(dm,&nasm->b[0]);CHKERRQ(ierr);
90*eaedb033SPeter Brune       ierr = DMCreateLocalVector(dm,&nasm->x[0]);CHKERRQ(ierr);
91*eaedb033SPeter Brune 
92*eaedb033SPeter Brune       if (!f) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide a block solve!");CHKERRQ(ierr);
93*eaedb033SPeter Brune 
94*eaedb033SPeter Brune       ierr = SNESSetFunction(nasm->subsnes[0],nasm->r[0],f,ctx);CHKERRQ(ierr);
95*eaedb033SPeter Brune 
96*eaedb033SPeter Brune       /* set up the local jacobian -- TODO: do this correctly */
97*eaedb033SPeter Brune       ierr = DMSNESGetBlockJacobian(dm,&fj,&ctx);CHKERRQ(ierr);
98*eaedb033SPeter Brune       ierr = VecGetSize(nasm->r[0],&n);CHKERRQ(ierr);
99*eaedb033SPeter Brune       ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,PETSC_DEFAULT,PETSC_NULL,&A);CHKERRQ(ierr);
100*eaedb033SPeter Brune       ierr = DMShellSetMatrix(subdm,A);CHKERRQ(ierr);
101*eaedb033SPeter Brune       ierr = SNESSetJacobian(nasm->subsnes[0],A,A,fj,ctx);CHKERRQ(ierr);
102*eaedb033SPeter Brune 
103*eaedb033SPeter Brune       ierr = SNESSetFromOptions(nasm->subsnes[0]);CHKERRQ(ierr);
104*eaedb033SPeter Brune     } else {
105*eaedb033SPeter Brune       SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!");
106*eaedb033SPeter Brune     }
107*eaedb033SPeter Brune   } else {
108*eaedb033SPeter Brune     SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!");
109*eaedb033SPeter Brune   }
110*eaedb033SPeter Brune   PetscFunctionReturn(0);
111*eaedb033SPeter Brune }
112*eaedb033SPeter Brune 
113*eaedb033SPeter Brune #undef __FUNCT__
114*eaedb033SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NASM"
115*eaedb033SPeter Brune PetscErrorCode SNESSetFromOptions_NASM(SNES snes)
116*eaedb033SPeter Brune {
117*eaedb033SPeter Brune   PetscErrorCode ierr;
118*eaedb033SPeter Brune 
119*eaedb033SPeter Brune   PetscFunctionBegin;
120*eaedb033SPeter Brune   ierr = PetscOptionsHead("SNES NASM options");CHKERRQ(ierr);
121*eaedb033SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
122*eaedb033SPeter Brune   PetscFunctionReturn(0);
123*eaedb033SPeter Brune }
124*eaedb033SPeter Brune 
125*eaedb033SPeter Brune #undef __FUNCT__
126*eaedb033SPeter Brune #define __FUNCT__ "SNESView_NASM"
127*eaedb033SPeter Brune PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
128*eaedb033SPeter Brune {
129*eaedb033SPeter Brune   PetscFunctionBegin;
130*eaedb033SPeter Brune   PetscFunctionReturn(0);
131*eaedb033SPeter Brune }
132*eaedb033SPeter Brune 
133*eaedb033SPeter Brune #undef __FUNCT__
134*eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains"
135*eaedb033SPeter Brune PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],IS iis[],IS ois[]) {
136*eaedb033SPeter Brune   PetscErrorCode ierr;
137*eaedb033SPeter Brune   PetscErrorCode (*f)(SNES,PetscInt,SNES*,IS*,IS*);
138*eaedb033SPeter Brune   PetscFunctionBegin;
139*eaedb033SPeter Brune   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr);
140*eaedb033SPeter Brune   ierr = (f)(snes,n,subsnes,iis,ois);CHKERRQ(ierr);
141*eaedb033SPeter Brune   PetscFunctionReturn(0);
142*eaedb033SPeter Brune }
143*eaedb033SPeter Brune 
144*eaedb033SPeter Brune EXTERN_C_BEGIN
145*eaedb033SPeter Brune #undef __FUNCT__
146*eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains_NASM"
147*eaedb033SPeter Brune PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],IS iis[],IS ois[]) {
148*eaedb033SPeter Brune   PetscInt       i;
149*eaedb033SPeter Brune   PetscErrorCode ierr;
150*eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM *)snes->data;
151*eaedb033SPeter Brune   PetscFunctionBegin;
152*eaedb033SPeter Brune   if (snes->setupcalled)SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
153*eaedb033SPeter Brune 
154*eaedb033SPeter Brune   if (!snes->setupcalled) {
155*eaedb033SPeter Brune     nasm->n       = n;
156*eaedb033SPeter Brune     nasm->ois     = 0;
157*eaedb033SPeter Brune     nasm->iis     = 0;
158*eaedb033SPeter Brune     if (ois) {
159*eaedb033SPeter Brune       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);}
160*eaedb033SPeter Brune     }
161*eaedb033SPeter Brune     if (iis) {
162*eaedb033SPeter Brune       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);}
163*eaedb033SPeter Brune     }
164*eaedb033SPeter Brune     if (ois) {
165*eaedb033SPeter Brune       ierr = PetscMalloc(n*sizeof(IS),&nasm->ois);CHKERRQ(ierr);
166*eaedb033SPeter Brune       for (i=0; i<n; i++) {
167*eaedb033SPeter Brune         nasm->ois[i] = ois[i];
168*eaedb033SPeter Brune       }
169*eaedb033SPeter Brune       if (!iis) {
170*eaedb033SPeter Brune         ierr = PetscMalloc(n*sizeof(IS),&nasm->iis);CHKERRQ(ierr);
171*eaedb033SPeter Brune         for (i=0; i<n; i++) {
172*eaedb033SPeter Brune           for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);}
173*eaedb033SPeter Brune           nasm->iis[i] = ois[i];
174*eaedb033SPeter Brune         }
175*eaedb033SPeter Brune       }
176*eaedb033SPeter Brune     }
177*eaedb033SPeter Brune     if (iis) {
178*eaedb033SPeter Brune       ierr = PetscMalloc(n*sizeof(IS),&nasm->iis);CHKERRQ(ierr);
179*eaedb033SPeter Brune       for (i=0; i<n; i++) {
180*eaedb033SPeter Brune         nasm->iis[i] = iis[i];
181*eaedb033SPeter Brune       }
182*eaedb033SPeter Brune       if (!ois) {
183*eaedb033SPeter Brune         ierr = PetscMalloc(n*sizeof(IS),&nasm->ois);CHKERRQ(ierr);
184*eaedb033SPeter Brune         for (i=0; i<n; i++) {
185*eaedb033SPeter Brune           for (i=0; i<n; i++) {
186*eaedb033SPeter Brune             ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);
187*eaedb033SPeter Brune             nasm->ois[i] = iis[i];
188*eaedb033SPeter Brune           }
189*eaedb033SPeter Brune         }
190*eaedb033SPeter Brune       }
191*eaedb033SPeter Brune     }
192*eaedb033SPeter Brune   }
193*eaedb033SPeter Brune   if (subsnes) {
194*eaedb033SPeter Brune     ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr);
195*eaedb033SPeter Brune     for (i=0; i<n; i++) {
196*eaedb033SPeter Brune       nasm->subsnes[i] = subsnes[i];
197*eaedb033SPeter Brune     }
198*eaedb033SPeter Brune   }
199*eaedb033SPeter Brune   PetscFunctionReturn(0);
200*eaedb033SPeter Brune }
201*eaedb033SPeter Brune EXTERN_C_END
202*eaedb033SPeter Brune 
203*eaedb033SPeter Brune #undef __FUNCT__
204*eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSolveLocal_Private"
205*eaedb033SPeter Brune PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec X) {
206*eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM *)snes->data;
207*eaedb033SPeter Brune   PetscInt       i;
208*eaedb033SPeter Brune   PetscErrorCode ierr;
209*eaedb033SPeter Brune   Vec            Xl,Bl;
210*eaedb033SPeter Brune   DM             dm;
211*eaedb033SPeter Brune   PetscFunctionBegin;
212*eaedb033SPeter Brune     /* restrict to the local system */
213*eaedb033SPeter Brune   if (nasm->usesdm) {
214*eaedb033SPeter Brune     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
215*eaedb033SPeter Brune     ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,nasm->x[0]);CHKERRQ(ierr);
216*eaedb033SPeter Brune     ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,nasm->x[0]);CHKERRQ(ierr);
217*eaedb033SPeter Brune     if (B) {
218*eaedb033SPeter Brune       ierr = DMGlobalToLocalBegin(dm,B,INSERT_VALUES,nasm->b[0]);CHKERRQ(ierr);
219*eaedb033SPeter Brune       ierr = DMGlobalToLocalEnd(dm,B,INSERT_VALUES,nasm->b[0]);CHKERRQ(ierr);
220*eaedb033SPeter Brune     }
221*eaedb033SPeter Brune   } else {
222*eaedb033SPeter Brune     /*
223*eaedb033SPeter Brune     for (i = 0;i < nasm->n;i++) {
224*eaedb033SPeter Brune       ierr = VecScatterBegin(nasm->gorestriction,X,nasm->gx,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
225*eaedb033SPeter Brune       ierr = VecScatterEnd(nasm->gorestriction,X,nasm->gx,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
226*eaedb033SPeter Brune       if (B) {
227*eaedb033SPeter Brune         ierr = VecScatterBegin(nasm->gorestriction,B,nasm->gb,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
228*eaedb033SPeter Brune         ierr = VecScatterEnd(nasm->gorestriction,B,nasm->gb,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
229*eaedb033SPeter Brune      } else {
230*eaedb033SPeter Brune      }
231*eaedb033SPeter Brune      }
232*eaedb033SPeter Brune      */
233*eaedb033SPeter Brune   }
234*eaedb033SPeter Brune   for(i=0;i<nasm->n;i++) {
235*eaedb033SPeter Brune     Xl = nasm->x[i];
236*eaedb033SPeter Brune     if (B) {
237*eaedb033SPeter Brune       Bl = nasm->b[i];
238*eaedb033SPeter Brune     } else {
239*eaedb033SPeter Brune       Bl = PETSC_NULL;
240*eaedb033SPeter Brune     }
241*eaedb033SPeter Brune     ierr = SNESSolve(nasm->subsnes[i],Bl,Xl);CHKERRQ(ierr);
242*eaedb033SPeter Brune   }
243*eaedb033SPeter Brune 
244*eaedb033SPeter Brune   if (nasm->usesdm) {
245*eaedb033SPeter Brune     ierr = DMLocalToGlobalBegin(dm,Xl,INSERT_VALUES,X);CHKERRQ(ierr);
246*eaedb033SPeter Brune     ierr = DMLocalToGlobalEnd(dm,Xl,INSERT_VALUES,X);CHKERRQ(ierr);
247*eaedb033SPeter Brune   } else {
248*eaedb033SPeter Brune     /*
249*eaedb033SPeter Brune     ierr = VecScatterBegin(nasm->girestriction,nasm->gx,X,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
250*eaedb033SPeter Brune     ierr = VecScatterEnd(nasm->girestriction,nasm->gx,X,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
251*eaedb033SPeter Brune      */
252*eaedb033SPeter Brune   }
253*eaedb033SPeter Brune   PetscFunctionReturn(0);
254*eaedb033SPeter Brune }
255*eaedb033SPeter Brune 
256*eaedb033SPeter Brune #undef __FUNCT__
257*eaedb033SPeter Brune #define __FUNCT__ "SNESSolve_NASM"
258*eaedb033SPeter Brune PetscErrorCode SNESSolve_NASM(SNES snes)
259*eaedb033SPeter Brune {
260*eaedb033SPeter Brune   Vec            F;
261*eaedb033SPeter Brune   Vec            X;
262*eaedb033SPeter Brune   Vec            B;
263*eaedb033SPeter Brune   PetscInt       i;
264*eaedb033SPeter Brune   PetscReal      fnorm;
265*eaedb033SPeter Brune   PetscErrorCode ierr;
266*eaedb033SPeter Brune   SNESNormType   normtype;
267*eaedb033SPeter Brune 
268*eaedb033SPeter Brune   PetscFunctionBegin;
269*eaedb033SPeter Brune 
270*eaedb033SPeter Brune   X = snes->vec_sol;
271*eaedb033SPeter Brune   F = snes->vec_func;
272*eaedb033SPeter Brune   B = snes->vec_rhs;
273*eaedb033SPeter Brune 
274*eaedb033SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
275*eaedb033SPeter Brune   snes->iter = 0;
276*eaedb033SPeter Brune   snes->norm = 0.;
277*eaedb033SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
278*eaedb033SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
279*eaedb033SPeter Brune   ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr);
280*eaedb033SPeter Brune   if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
281*eaedb033SPeter Brune     /* compute the initial function and preconditioned update delX */
282*eaedb033SPeter Brune     if (!snes->vec_func_init_set) {
283*eaedb033SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
284*eaedb033SPeter Brune       if (snes->domainerror) {
285*eaedb033SPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
286*eaedb033SPeter Brune         PetscFunctionReturn(0);
287*eaedb033SPeter Brune       }
288*eaedb033SPeter Brune     } else {
289*eaedb033SPeter Brune       snes->vec_func_init_set = PETSC_FALSE;
290*eaedb033SPeter Brune     }
291*eaedb033SPeter Brune 
292*eaedb033SPeter Brune     /* convergence test */
293*eaedb033SPeter Brune     if (!snes->norm_init_set) {
294*eaedb033SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
295*eaedb033SPeter Brune       if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
296*eaedb033SPeter Brune     } else {
297*eaedb033SPeter Brune       fnorm = snes->norm_init;
298*eaedb033SPeter Brune       snes->norm_init_set = PETSC_FALSE;
299*eaedb033SPeter Brune     }
300*eaedb033SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
301*eaedb033SPeter Brune     snes->iter = 0;
302*eaedb033SPeter Brune     snes->norm = fnorm;
303*eaedb033SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
304*eaedb033SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
305*eaedb033SPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
306*eaedb033SPeter Brune 
307*eaedb033SPeter Brune     /* set parameter for default relative tolerance convergence test */
308*eaedb033SPeter Brune     snes->ttol = fnorm*snes->rtol;
309*eaedb033SPeter Brune 
310*eaedb033SPeter Brune     /* test convergence */
311*eaedb033SPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
312*eaedb033SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
313*eaedb033SPeter Brune   } else {
314*eaedb033SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
315*eaedb033SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
316*eaedb033SPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
317*eaedb033SPeter Brune   }
318*eaedb033SPeter Brune 
319*eaedb033SPeter Brune   /* Call general purpose update function */
320*eaedb033SPeter Brune   if (snes->ops->update) {
321*eaedb033SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
322*eaedb033SPeter Brune   }
323*eaedb033SPeter Brune 
324*eaedb033SPeter Brune   for (i = 0; i < snes->max_its; i++) {
325*eaedb033SPeter Brune     ierr = SNESNASMSolveLocal_Private(snes,B,X);CHKERRQ(ierr);
326*eaedb033SPeter Brune     if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) {
327*eaedb033SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
328*eaedb033SPeter Brune       if (snes->domainerror) {
329*eaedb033SPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
330*eaedb033SPeter Brune         PetscFunctionReturn(0);
331*eaedb033SPeter Brune       }
332*eaedb033SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
333*eaedb033SPeter Brune       if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
334*eaedb033SPeter Brune     }
335*eaedb033SPeter Brune     /* Monitor convergence */
336*eaedb033SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
337*eaedb033SPeter Brune     snes->iter = i+1;
338*eaedb033SPeter Brune     snes->norm = fnorm;
339*eaedb033SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
340*eaedb033SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
341*eaedb033SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
342*eaedb033SPeter Brune     /* Test for convergence */
343*eaedb033SPeter Brune     if (normtype == SNES_NORM_FUNCTION) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
344*eaedb033SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
345*eaedb033SPeter Brune     /* Call general purpose update function */
346*eaedb033SPeter Brune     if (snes->ops->update) {
347*eaedb033SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
348*eaedb033SPeter Brune     }
349*eaedb033SPeter Brune   }
350*eaedb033SPeter Brune   if (normtype == SNES_NORM_FUNCTION) {
351*eaedb033SPeter Brune     if (i == snes->max_its) {
352*eaedb033SPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
353*eaedb033SPeter Brune       if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
354*eaedb033SPeter Brune     }
355*eaedb033SPeter Brune   } else {
356*eaedb033SPeter Brune     if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */
357*eaedb033SPeter Brune   }
358*eaedb033SPeter Brune   PetscFunctionReturn(0);
359*eaedb033SPeter Brune }
360*eaedb033SPeter Brune 
361*eaedb033SPeter Brune /*MC
362*eaedb033SPeter Brune   SNESNASM - Nonlinear Additive Schwartz
363*eaedb033SPeter Brune 
364*eaedb033SPeter Brune    Level: advanced
365*eaedb033SPeter Brune 
366*eaedb033SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
367*eaedb033SPeter Brune M*/
368*eaedb033SPeter Brune 
369*eaedb033SPeter Brune EXTERN_C_BEGIN
370*eaedb033SPeter Brune #undef __FUNCT__
371*eaedb033SPeter Brune #define __FUNCT__ "SNESCreate_NASM"
372*eaedb033SPeter Brune PetscErrorCode SNESCreate_NASM(SNES snes)
373*eaedb033SPeter Brune {
374*eaedb033SPeter Brune   SNES_NASM        *nasm;
375*eaedb033SPeter Brune   PetscErrorCode ierr;
376*eaedb033SPeter Brune 
377*eaedb033SPeter Brune   PetscFunctionBegin;
378*eaedb033SPeter Brune 
379*eaedb033SPeter Brune 
380*eaedb033SPeter Brune   ierr = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr);
381*eaedb033SPeter Brune   snes->data = (void*)nasm;
382*eaedb033SPeter Brune 
383*eaedb033SPeter Brune   nasm->n                 = PETSC_DECIDE;
384*eaedb033SPeter Brune   nasm->subsnes           = 0;
385*eaedb033SPeter Brune   nasm->x                 = 0;
386*eaedb033SPeter Brune   nasm->b                 = 0;
387*eaedb033SPeter Brune   nasm->ois               = 0;
388*eaedb033SPeter Brune   nasm->iis               = 0;
389*eaedb033SPeter Brune 
390*eaedb033SPeter Brune   snes->ops->destroy        = SNESDestroy_NASM;
391*eaedb033SPeter Brune   snes->ops->setup          = SNESSetUp_NASM;
392*eaedb033SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
393*eaedb033SPeter Brune   snes->ops->view           = SNESView_NASM;
394*eaedb033SPeter Brune   snes->ops->solve          = SNESSolve_NASM;
395*eaedb033SPeter Brune   snes->ops->reset          = SNESReset_NASM;
396*eaedb033SPeter Brune 
397*eaedb033SPeter Brune   snes->usesksp             = PETSC_FALSE;
398*eaedb033SPeter Brune   snes->usespc              = PETSC_FALSE;
399*eaedb033SPeter Brune 
400*eaedb033SPeter Brune   if (!snes->tolerancesset) {
401*eaedb033SPeter Brune     snes->max_its             = 10000;
402*eaedb033SPeter Brune     snes->max_funcs           = 10000;
403*eaedb033SPeter Brune   }
404*eaedb033SPeter Brune 
405*eaedb033SPeter Brune     ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNASMSetSubdomains_C","SNESNASMSetSubdomains_NASM",
406*eaedb033SPeter Brune                     SNESNASMSetSubdomains_NASM);CHKERRQ(ierr);
407*eaedb033SPeter Brune 
408*eaedb033SPeter Brune   PetscFunctionReturn(0);
409*eaedb033SPeter Brune }
410*eaedb033SPeter Brune EXTERN_C_END
411