xref: /petsc/src/snes/impls/nasm/nasm.c (revision 70c78f051b666a4db181aadf0c3f2dbd76338ef6)
1eaedb033SPeter Brune #include <petsc-private/snesimpl.h>             /*I   "petscsnes.h"   I*/
2111ade9eSPeter Brune #include <petscdm.h>
3eaedb033SPeter Brune 
4eaedb033SPeter Brune typedef struct {
5eaedb033SPeter Brune   PetscInt   n;                   /* local subdomains */
6eaedb033SPeter Brune   SNES       *subsnes;            /* nonlinear solvers for each subdomain */
7eaedb033SPeter Brune   Vec        *x;                  /* solution vectors */
8111ade9eSPeter Brune   Vec        *xl;                 /* solution local vectors */
9111ade9eSPeter Brune   Vec        *y;                  /* step vectors */
10eaedb033SPeter Brune   Vec        *b;                  /* rhs vectors */
11111ade9eSPeter Brune   VecScatter *oscatter;           /* scatter from global space to the subdomain global space */
12111ade9eSPeter Brune   VecScatter *iscatter;           /* scatter from global space to the nonoverlapping subdomain space */
13111ade9eSPeter Brune   VecScatter *gscatter;           /* scatter from global space to the subdomain local space */
14111ade9eSPeter Brune   PCASMType  type;                /* ASM type */
15111ade9eSPeter Brune   PetscBool  usesdm;              /* use the DM for setting up the subproblems */
16eaedb033SPeter Brune } SNES_NASM;
17eaedb033SPeter Brune 
18eaedb033SPeter Brune #undef __FUNCT__
19eaedb033SPeter Brune #define __FUNCT__ "SNESReset_NASM"
20eaedb033SPeter Brune PetscErrorCode SNESReset_NASM(SNES snes)
21eaedb033SPeter Brune {
22eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
23eaedb033SPeter Brune   PetscErrorCode ierr;
24eaedb033SPeter Brune   PetscInt       i;
256e111a19SKarl Rupp 
26eaedb033SPeter Brune   PetscFunctionBegin;
27eaedb033SPeter Brune   for (i=0; i<nasm->n; i++) {
28111ade9eSPeter Brune     if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); }
29f5f7c1b9SKarl Rupp     if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); }
30111ade9eSPeter Brune     if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); }
31bc8c1f72SJose Roman     if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); }
32eaedb033SPeter Brune 
33bc8c1f72SJose Roman     if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); }
34111ade9eSPeter Brune     if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); }
35111ade9eSPeter Brune     if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); }
36111ade9eSPeter Brune     if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); }
37eaedb033SPeter Brune   }
38111ade9eSPeter Brune 
39111ade9eSPeter Brune   if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);}
40111ade9eSPeter Brune   if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);}
41111ade9eSPeter Brune   if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);}
42111ade9eSPeter Brune   if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);}
43111ade9eSPeter Brune 
44111ade9eSPeter Brune   if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);}
45111ade9eSPeter Brune   if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);}
46111ade9eSPeter Brune   if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);}
47111ade9eSPeter Brune   if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);}
48eaedb033SPeter Brune   PetscFunctionReturn(0);
49eaedb033SPeter Brune }
50eaedb033SPeter Brune 
51eaedb033SPeter Brune #undef __FUNCT__
52eaedb033SPeter Brune #define __FUNCT__ "SNESDestroy_NASM"
53eaedb033SPeter Brune PetscErrorCode SNESDestroy_NASM(SNES snes)
54eaedb033SPeter Brune {
55eaedb033SPeter Brune   PetscErrorCode ierr;
566e111a19SKarl Rupp 
57eaedb033SPeter Brune   PetscFunctionBegin;
58eaedb033SPeter Brune   ierr = SNESReset_NASM(snes);CHKERRQ(ierr);
5922d28d08SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
60eaedb033SPeter Brune   PetscFunctionReturn(0);
61eaedb033SPeter Brune }
62eaedb033SPeter Brune 
63eaedb033SPeter Brune #undef __FUNCT__
64111ade9eSPeter Brune #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private"
650adebc6cSBarry Smith PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx)
660adebc6cSBarry Smith {
67111ade9eSPeter Brune   PetscErrorCode ierr;
68111ade9eSPeter Brune   Vec            bcs = (Vec)ctx;
696e111a19SKarl Rupp 
70111ade9eSPeter Brune   PetscFunctionBegin;
71111ade9eSPeter Brune   ierr = VecCopy(bcs,l);CHKERRQ(ierr);
72111ade9eSPeter Brune   PetscFunctionReturn(0);
73111ade9eSPeter Brune }
74111ade9eSPeter Brune 
75111ade9eSPeter Brune #undef __FUNCT__
76eaedb033SPeter Brune #define __FUNCT__ "SNESSetUp_NASM"
77eaedb033SPeter Brune PetscErrorCode SNESSetUp_NASM(SNES snes)
78eaedb033SPeter Brune {
79eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
80eaedb033SPeter Brune   PetscErrorCode ierr;
810a696f66SPeter Brune   DM             dm,ddm;
82111ade9eSPeter Brune   DM             *subdms;
83111ade9eSPeter Brune   PetscInt       i;
84eaedb033SPeter Brune   const char     *optionsprefix;
85111ade9eSPeter Brune   Vec            F;
86eaedb033SPeter Brune 
87eaedb033SPeter Brune   PetscFunctionBegin;
88eaedb033SPeter Brune   if (!nasm->subsnes) {
89eaedb033SPeter Brune     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
900a696f66SPeter Brune     if (dm) {
91eaedb033SPeter Brune       nasm->usesdm = PETSC_TRUE;
92111ade9eSPeter Brune       ierr         = DMCreateDomainDecomposition(dm,&nasm->n,PETSC_NULL,PETSC_NULL,PETSC_NULL,&subdms);CHKERRQ(ierr);
930a696f66SPeter Brune       if (!subdms) {
940a696f66SPeter Brune         ierr = DMCreateDomainDecompositionDM(dm,"default",&ddm);CHKERRQ(ierr);
95f23aa3ddSBarry Smith         if (!ddm) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
960a696f66SPeter Brune         ierr = SNESSetDM(snes,ddm);CHKERRQ(ierr);
970a696f66SPeter Brune         ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
980a696f66SPeter Brune         ierr = DMCreateDomainDecomposition(dm,&nasm->n,PETSC_NULL,PETSC_NULL,PETSC_NULL,&subdms);CHKERRQ(ierr);
990a696f66SPeter Brune       }
100f23aa3ddSBarry Smith       if (!subdms) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
101111ade9eSPeter Brune       ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr);
102eaedb033SPeter Brune 
103eaedb033SPeter Brune       ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
104111ade9eSPeter Brune       ierr = PetscMalloc(nasm->n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr);
105111ade9eSPeter Brune 
106111ade9eSPeter Brune       for (i=0; i<nasm->n; i++) {
107cdb298fcSPeter Brune         ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr);
108cdb298fcSPeter Brune         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr);
109cdb298fcSPeter Brune         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr);
110cdb298fcSPeter Brune         ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr);
111cdb298fcSPeter Brune         ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr);
112111ade9eSPeter Brune         ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr);
113111ade9eSPeter Brune       }
114111ade9eSPeter Brune       ierr = PetscFree(subdms);CHKERRQ(ierr);
115f23aa3ddSBarry Smith     } else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!");
116f23aa3ddSBarry Smith   } else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!");
117111ade9eSPeter Brune   /* allocate the global vectors */
118111ade9eSPeter Brune   ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->x);CHKERRQ(ierr);
119111ade9eSPeter Brune   ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->xl);CHKERRQ(ierr);
120111ade9eSPeter Brune   ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->y);CHKERRQ(ierr);
121111ade9eSPeter Brune   ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->b);CHKERRQ(ierr);
122111ade9eSPeter Brune 
123111ade9eSPeter Brune   for (i=0; i<nasm->n; i++) {
124111ade9eSPeter Brune     DM subdm;
125111ade9eSPeter Brune     ierr = SNESGetFunction(nasm->subsnes[i],&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
126111ade9eSPeter Brune     ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);
127111ade9eSPeter Brune     ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);
128111ade9eSPeter Brune     ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);
129111ade9eSPeter Brune     ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr);
130111ade9eSPeter Brune     ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr);
131111ade9eSPeter Brune     ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,PETSC_NULL,nasm->xl[i]);CHKERRQ(ierr);
132111ade9eSPeter Brune   }
133eaedb033SPeter Brune   PetscFunctionReturn(0);
134eaedb033SPeter Brune }
135eaedb033SPeter Brune 
136eaedb033SPeter Brune #undef __FUNCT__
137eaedb033SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NASM"
138eaedb033SPeter Brune PetscErrorCode SNESSetFromOptions_NASM(SNES snes)
139eaedb033SPeter Brune {
140eaedb033SPeter Brune   PetscErrorCode    ierr;
1410a696f66SPeter Brune   DM                dm,ddm;
1420a696f66SPeter Brune   char              ddm_name[1024];
143111ade9eSPeter Brune   PCASMType         asmtype;
144111ade9eSPeter Brune   const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0};
145111ade9eSPeter Brune   PetscBool         flg;
146111ade9eSPeter Brune   SNES_NASM         *nasm = (SNES_NASM*)snes->data;
1476e111a19SKarl Rupp 
148eaedb033SPeter Brune   PetscFunctionBegin;
149111ade9eSPeter Brune   ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr);
150111ade9eSPeter Brune   ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
1511aa26658SKarl Rupp   if (flg) nasm->type = asmtype;
1520a696f66SPeter Brune   ierr = PetscOptionsString("-snes_nasm_decomposition", "Name of the DM defining the composition", "SNESSetDM", ddm_name, ddm_name,1024,&flg);CHKERRQ(ierr);
1530a696f66SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1540a696f66SPeter Brune   if (flg) {
1550a696f66SPeter Brune     if (dm) {
1560a696f66SPeter Brune       ierr = DMCreateDomainDecompositionDM(dm, ddm_name, &ddm);CHKERRQ(ierr);
157f23aa3ddSBarry Smith       if (!ddm) SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "Unknown DM decomposition name %s", ddm_name);
1580a696f66SPeter Brune       ierr = PetscInfo(snes,"Using domain decomposition DM defined using options database\n");CHKERRQ(ierr);
1590a696f66SPeter Brune       ierr = SNESSetDM(snes,ddm);CHKERRQ(ierr);
160f23aa3ddSBarry Smith     } else SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "No DM to decompose");
1610a696f66SPeter Brune   }
162eaedb033SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
163eaedb033SPeter Brune   PetscFunctionReturn(0);
164eaedb033SPeter Brune }
165eaedb033SPeter Brune 
166eaedb033SPeter Brune #undef __FUNCT__
167eaedb033SPeter Brune #define __FUNCT__ "SNESView_NASM"
168eaedb033SPeter Brune PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
169eaedb033SPeter Brune {
170eaedb033SPeter Brune   PetscFunctionBegin;
171eaedb033SPeter Brune   PetscFunctionReturn(0);
172eaedb033SPeter Brune }
173eaedb033SPeter Brune 
174eaedb033SPeter Brune #undef __FUNCT__
175eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains"
176a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
177a6dfd86eSKarl Rupp {
178eaedb033SPeter Brune   PetscErrorCode ierr;
179111ade9eSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*);
1806e111a19SKarl Rupp 
181eaedb033SPeter Brune   PetscFunctionBegin;
182eaedb033SPeter Brune   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr);
183111ade9eSPeter Brune   ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);
184eaedb033SPeter Brune   PetscFunctionReturn(0);
185eaedb033SPeter Brune }
186eaedb033SPeter Brune 
187eaedb033SPeter Brune EXTERN_C_BEGIN
188eaedb033SPeter Brune #undef __FUNCT__
189eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains_NASM"
190a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
191a6dfd86eSKarl Rupp {
192eaedb033SPeter Brune   PetscInt       i;
193eaedb033SPeter Brune   PetscErrorCode ierr;
194eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
1956e111a19SKarl Rupp 
196eaedb033SPeter Brune   PetscFunctionBegin;
197eaedb033SPeter Brune   if (snes->setupcalled) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
198eaedb033SPeter Brune 
199111ade9eSPeter Brune   /* tear down the previously set things */
200111ade9eSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
201111ade9eSPeter Brune 
202eaedb033SPeter Brune   nasm->n = n;
203111ade9eSPeter Brune   if (oscatter) {
204111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);}
205eaedb033SPeter Brune   }
206111ade9eSPeter Brune   if (iscatter) {
207111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);}
208eaedb033SPeter Brune   }
209111ade9eSPeter Brune   if (gscatter) {
210111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);}
211111ade9eSPeter Brune   }
212111ade9eSPeter Brune   if (oscatter) {
213111ade9eSPeter Brune     ierr = PetscMalloc(n*sizeof(IS),&nasm->oscatter);CHKERRQ(ierr);
214eaedb033SPeter Brune     for (i=0; i<n; i++) {
215111ade9eSPeter Brune       nasm->oscatter[i] = oscatter[i];
216eaedb033SPeter Brune     }
217111ade9eSPeter Brune   }
218111ade9eSPeter Brune   if (iscatter) {
219111ade9eSPeter Brune     ierr = PetscMalloc(n*sizeof(IS),&nasm->iscatter);CHKERRQ(ierr);
220eaedb033SPeter Brune     for (i=0; i<n; i++) {
221111ade9eSPeter Brune       nasm->iscatter[i] = iscatter[i];
222eaedb033SPeter Brune     }
223eaedb033SPeter Brune   }
224111ade9eSPeter Brune   if (gscatter) {
225111ade9eSPeter Brune     ierr = PetscMalloc(n*sizeof(IS),&nasm->gscatter);CHKERRQ(ierr);
226eaedb033SPeter Brune     for (i=0; i<n; i++) {
227111ade9eSPeter Brune       nasm->gscatter[i] = gscatter[i];
228eaedb033SPeter Brune     }
229eaedb033SPeter Brune   }
230111ade9eSPeter Brune 
231eaedb033SPeter Brune   if (subsnes) {
232eaedb033SPeter Brune     ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr);
233eaedb033SPeter Brune     for (i=0; i<n; i++) {
234eaedb033SPeter Brune       nasm->subsnes[i] = subsnes[i];
235eaedb033SPeter Brune     }
236eaedb033SPeter Brune   }
237eaedb033SPeter Brune   PetscFunctionReturn(0);
238eaedb033SPeter Brune }
239eaedb033SPeter Brune EXTERN_C_END
240eaedb033SPeter Brune 
241eaedb033SPeter Brune #undef __FUNCT__
242eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSolveLocal_Private"
2430adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X)
2440adebc6cSBarry Smith {
245eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
246258e1594SPeter Brune   SNES           subsnes;
247eaedb033SPeter Brune   PetscInt       i;
248eaedb033SPeter Brune   PetscErrorCode ierr;
249111ade9eSPeter Brune   Vec            Xlloc,Xl,Bl,Yl;
250111ade9eSPeter Brune   VecScatter     iscat,oscat,gscat;
251111ade9eSPeter Brune   DM             dm,subdm;
2520adebc6cSBarry Smith 
253eaedb033SPeter Brune   PetscFunctionBegin;
254eaedb033SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
255111ade9eSPeter Brune   ierr = VecSet(Y,0);CHKERRQ(ierr);
256eaedb033SPeter Brune   for (i=0; i<nasm->n; i++) {
257*70c78f05SPeter Brune     /* scatter the solution to the local solution */
258*70c78f05SPeter Brune     Xlloc = nasm->xl[i];
259*70c78f05SPeter Brune     gscat   = nasm->gscatter[i];
260*70c78f05SPeter Brune     oscat   = nasm->oscatter[i];
261*70c78f05SPeter Brune     ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
262*70c78f05SPeter Brune     if (B) {
263*70c78f05SPeter Brune       /* scatter the RHS to the local RHS */
264*70c78f05SPeter Brune       Bl   = nasm->b[i];
265*70c78f05SPeter Brune       ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
266*70c78f05SPeter Brune     } else {
267*70c78f05SPeter Brune       Bl = PETSC_NULL;
268*70c78f05SPeter Brune     }
269*70c78f05SPeter Brune   }
270*70c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
271*70c78f05SPeter Brune     Xlloc = nasm->xl[i];
272*70c78f05SPeter Brune     ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
273*70c78f05SPeter Brune     if (B) {
274*70c78f05SPeter Brune       ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
275*70c78f05SPeter Brune     }
276*70c78f05SPeter Brune   }
277*70c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
278*70c78f05SPeter Brune     Xl    = nasm->x[i];
279*70c78f05SPeter Brune     Xlloc = nasm->xl[i];
280*70c78f05SPeter Brune     Yl    = nasm->y[i];
281258e1594SPeter Brune     subsnes = nasm->subsnes[i];
282258e1594SPeter Brune     ierr    = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
283111ade9eSPeter Brune     iscat   = nasm->iscatter[i];
284111ade9eSPeter Brune     oscat   = nasm->oscatter[i];
285111ade9eSPeter Brune     gscat   = nasm->gscatter[i];
286111ade9eSPeter Brune     ierr    = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
287111ade9eSPeter Brune 
288*70c78f05SPeter Brune     ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
289*70c78f05SPeter Brune     ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
290111ade9eSPeter Brune     ierr = VecCopy(Xl,Yl);CHKERRQ(ierr);
291258e1594SPeter Brune     ierr = SNESSolve(subsnes,Bl,Yl);CHKERRQ(ierr);
292111ade9eSPeter Brune     ierr = VecAXPY(Yl,-1.0,Xl);CHKERRQ(ierr);
293*70c78f05SPeter Brune   }
294111ade9eSPeter Brune 
295*70c78f05SPeter Brune   ierr = MPI_Barrier(((PetscObject)snes)->comm);CHKERRQ(ierr);
296*70c78f05SPeter Brune 
297*70c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
298*70c78f05SPeter Brune     Yl    = nasm->y[i];
299*70c78f05SPeter Brune     iscat   = nasm->iscatter[i];
300*70c78f05SPeter Brune     oscat   = nasm->oscatter[i];
301111ade9eSPeter Brune     if (nasm->type == PC_ASM_BASIC) {
302111ade9eSPeter Brune       ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
303111ade9eSPeter Brune     } else if (nasm->type == PC_ASM_RESTRICT) {
304111ade9eSPeter Brune       ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
305f23aa3ddSBarry Smith     } else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Only basic and interpolate types are supported for SNESNASM");
306eaedb033SPeter Brune   }
307eaedb033SPeter Brune 
308*70c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
309*70c78f05SPeter Brune     Yl    = nasm->y[i];
310*70c78f05SPeter Brune     iscat   = nasm->iscatter[i];
311*70c78f05SPeter Brune     oscat   = nasm->oscatter[i];
312*70c78f05SPeter Brune     if (nasm->type == PC_ASM_BASIC) {
313*70c78f05SPeter Brune       ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
314*70c78f05SPeter Brune     } else if (nasm->type == PC_ASM_RESTRICT) {
315*70c78f05SPeter Brune       ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
316*70c78f05SPeter Brune     } else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Only basic and interpolate types are supported for SNESNASM");
317*70c78f05SPeter Brune 
318*70c78f05SPeter Brune   }
319*70c78f05SPeter Brune 
320*70c78f05SPeter Brune   ierr = MPI_Barrier(((PetscObject)snes)->comm);CHKERRQ(ierr);
321cd939e56SPeter Brune 
322111ade9eSPeter Brune   ierr = VecAXPY(X,1.0,Y);CHKERRQ(ierr);
323eaedb033SPeter Brune   PetscFunctionReturn(0);
324eaedb033SPeter Brune }
325eaedb033SPeter Brune 
326eaedb033SPeter Brune #undef __FUNCT__
327eaedb033SPeter Brune #define __FUNCT__ "SNESSolve_NASM"
328eaedb033SPeter Brune PetscErrorCode SNESSolve_NASM(SNES snes)
329eaedb033SPeter Brune {
330eaedb033SPeter Brune   Vec            F;
331eaedb033SPeter Brune   Vec            X;
332eaedb033SPeter Brune   Vec            B;
333111ade9eSPeter Brune   Vec            Y;
334eaedb033SPeter Brune   PetscInt       i;
335eaedb033SPeter Brune   PetscReal      fnorm;
336eaedb033SPeter Brune   PetscErrorCode ierr;
337eaedb033SPeter Brune   SNESNormType   normtype;
338eaedb033SPeter Brune 
339eaedb033SPeter Brune   PetscFunctionBegin;
340eaedb033SPeter Brune   X = snes->vec_sol;
341111ade9eSPeter Brune   Y = snes->vec_sol_update;
342eaedb033SPeter Brune   F = snes->vec_func;
343eaedb033SPeter Brune   B = snes->vec_rhs;
344eaedb033SPeter Brune 
345eaedb033SPeter Brune   ierr         = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
346eaedb033SPeter Brune   snes->iter   = 0;
347eaedb033SPeter Brune   snes->norm   = 0.;
348eaedb033SPeter Brune   ierr         = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
349eaedb033SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
350eaedb033SPeter Brune   ierr         = SNESGetNormType(snes, &normtype);CHKERRQ(ierr);
351eaedb033SPeter Brune   if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
352eaedb033SPeter Brune     /* compute the initial function and preconditioned update delX */
353eaedb033SPeter Brune     if (!snes->vec_func_init_set) {
354eaedb033SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
355eaedb033SPeter Brune       if (snes->domainerror) {
356eaedb033SPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
357eaedb033SPeter Brune         PetscFunctionReturn(0);
358eaedb033SPeter Brune       }
3591aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
360eaedb033SPeter Brune 
361eaedb033SPeter Brune     /* convergence test */
362eaedb033SPeter Brune     if (!snes->norm_init_set) {
363eaedb033SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
364eaedb033SPeter Brune       if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
365eaedb033SPeter Brune     } else {
366eaedb033SPeter Brune       fnorm               = snes->norm_init;
367eaedb033SPeter Brune       snes->norm_init_set = PETSC_FALSE;
368eaedb033SPeter Brune     }
369eaedb033SPeter Brune     ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
370eaedb033SPeter Brune     snes->iter = 0;
371eaedb033SPeter Brune     snes->norm = fnorm;
372eaedb033SPeter Brune     ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
373eaedb033SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
374eaedb033SPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
375eaedb033SPeter Brune 
376eaedb033SPeter Brune     /* set parameter for default relative tolerance convergence test */
377eaedb033SPeter Brune     snes->ttol = fnorm*snes->rtol;
378eaedb033SPeter Brune 
379eaedb033SPeter Brune     /* test convergence */
380eaedb033SPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
381eaedb033SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
382eaedb033SPeter Brune   } else {
383eaedb033SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
384eaedb033SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
385eaedb033SPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
386eaedb033SPeter Brune   }
387eaedb033SPeter Brune 
388eaedb033SPeter Brune   /* Call general purpose update function */
389eaedb033SPeter Brune   if (snes->ops->update) {
390eaedb033SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
391eaedb033SPeter Brune   }
392eaedb033SPeter Brune 
393eaedb033SPeter Brune   for (i = 0; i < snes->max_its; i++) {
394111ade9eSPeter Brune     ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr);
395eaedb033SPeter Brune     if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) {
396eaedb033SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
397eaedb033SPeter Brune       if (snes->domainerror) {
398eaedb033SPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
399eaedb033SPeter Brune         PetscFunctionReturn(0);
400eaedb033SPeter Brune       }
401eaedb033SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
402eaedb033SPeter Brune       if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
403eaedb033SPeter Brune     }
404eaedb033SPeter Brune     /* Monitor convergence */
405eaedb033SPeter Brune     ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
406eaedb033SPeter Brune     snes->iter = i+1;
407eaedb033SPeter Brune     snes->norm = fnorm;
408eaedb033SPeter Brune     ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
409eaedb033SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
410eaedb033SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
411eaedb033SPeter Brune     /* Test for convergence */
412eaedb033SPeter Brune     if (normtype == SNES_NORM_FUNCTION) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
413eaedb033SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
414eaedb033SPeter Brune     /* Call general purpose update function */
415eaedb033SPeter Brune     if (snes->ops->update) {
416eaedb033SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
417eaedb033SPeter Brune     }
418eaedb033SPeter Brune   }
419eaedb033SPeter Brune   if (normtype == SNES_NORM_FUNCTION) {
420eaedb033SPeter Brune     if (i == snes->max_its) {
421eaedb033SPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
422eaedb033SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
423eaedb033SPeter Brune     }
4241aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */
425eaedb033SPeter Brune   PetscFunctionReturn(0);
426eaedb033SPeter Brune }
427eaedb033SPeter Brune 
428eaedb033SPeter Brune /*MC
429eaedb033SPeter Brune   SNESNASM - Nonlinear Additive Schwartz
430eaedb033SPeter Brune 
431eaedb033SPeter Brune    Level: advanced
432eaedb033SPeter Brune 
433eaedb033SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
434eaedb033SPeter Brune M*/
435eaedb033SPeter Brune 
436eaedb033SPeter Brune EXTERN_C_BEGIN
437eaedb033SPeter Brune #undef __FUNCT__
438eaedb033SPeter Brune #define __FUNCT__ "SNESCreate_NASM"
439eaedb033SPeter Brune PetscErrorCode SNESCreate_NASM(SNES snes)
440eaedb033SPeter Brune {
441eaedb033SPeter Brune   SNES_NASM      *nasm;
442eaedb033SPeter Brune   PetscErrorCode ierr;
443eaedb033SPeter Brune 
444eaedb033SPeter Brune   PetscFunctionBegin;
445eaedb033SPeter Brune   ierr       = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr);
446eaedb033SPeter Brune   snes->data = (void*)nasm;
447eaedb033SPeter Brune 
448eaedb033SPeter Brune   nasm->n        = PETSC_DECIDE;
449eaedb033SPeter Brune   nasm->subsnes  = 0;
450eaedb033SPeter Brune   nasm->x        = 0;
451111ade9eSPeter Brune   nasm->xl       = 0;
452111ade9eSPeter Brune   nasm->y        = 0;
453eaedb033SPeter Brune   nasm->b        = 0;
454111ade9eSPeter Brune   nasm->oscatter = 0;
455111ade9eSPeter Brune   nasm->iscatter = 0;
456111ade9eSPeter Brune   nasm->gscatter = 0;
457111ade9eSPeter Brune 
458111ade9eSPeter Brune   nasm->type = PC_ASM_BASIC;
459eaedb033SPeter Brune 
460eaedb033SPeter Brune   snes->ops->destroy        = SNESDestroy_NASM;
461eaedb033SPeter Brune   snes->ops->setup          = SNESSetUp_NASM;
462eaedb033SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
463eaedb033SPeter Brune   snes->ops->view           = SNESView_NASM;
464eaedb033SPeter Brune   snes->ops->solve          = SNESSolve_NASM;
465eaedb033SPeter Brune   snes->ops->reset          = SNESReset_NASM;
466eaedb033SPeter Brune 
467eaedb033SPeter Brune   snes->usesksp = PETSC_FALSE;
468eaedb033SPeter Brune   snes->usespc  = PETSC_FALSE;
469eaedb033SPeter Brune 
470eaedb033SPeter Brune   if (!snes->tolerancesset) {
471eaedb033SPeter Brune     snes->max_its   = 10000;
472eaedb033SPeter Brune     snes->max_funcs = 10000;
473eaedb033SPeter Brune   }
474eaedb033SPeter Brune 
475eaedb033SPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNASMSetSubdomains_C","SNESNASMSetSubdomains_NASM",
476eaedb033SPeter Brune                                            SNESNASMSetSubdomains_NASM);CHKERRQ(ierr);
477eaedb033SPeter Brune   PetscFunctionReturn(0);
478eaedb033SPeter Brune }
479eaedb033SPeter Brune EXTERN_C_END
480