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