xref: /petsc/src/snes/impls/nasm/nasm.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
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 */
16b20c023fSPeter Brune 
17b20c023fSPeter Brune   /* logging events */
18b20c023fSPeter Brune   PetscLogEvent eventrestrictinterp;
19b20c023fSPeter Brune   PetscLogEvent eventsubsolve;
20eaedb033SPeter Brune } SNES_NASM;
21eaedb033SPeter Brune 
22b20c023fSPeter Brune const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0};
23b20c023fSPeter Brune 
24eaedb033SPeter Brune #undef __FUNCT__
25eaedb033SPeter Brune #define __FUNCT__ "SNESReset_NASM"
26eaedb033SPeter Brune PetscErrorCode SNESReset_NASM(SNES snes)
27eaedb033SPeter Brune {
28eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
29eaedb033SPeter Brune   PetscErrorCode ierr;
30eaedb033SPeter Brune   PetscInt       i;
316e111a19SKarl Rupp 
32eaedb033SPeter Brune   PetscFunctionBegin;
33eaedb033SPeter Brune   for (i=0; i<nasm->n; i++) {
34111ade9eSPeter Brune     if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); }
35f5f7c1b9SKarl Rupp     if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); }
36111ade9eSPeter Brune     if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); }
37bc8c1f72SJose Roman     if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); }
38eaedb033SPeter Brune 
39bc8c1f72SJose Roman     if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); }
40111ade9eSPeter Brune     if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); }
41111ade9eSPeter Brune     if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); }
42111ade9eSPeter Brune     if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); }
43eaedb033SPeter Brune   }
44111ade9eSPeter Brune 
45111ade9eSPeter Brune   if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);}
46111ade9eSPeter Brune   if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);}
47111ade9eSPeter Brune   if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);}
48111ade9eSPeter Brune   if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);}
49111ade9eSPeter Brune 
50111ade9eSPeter Brune   if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);}
51111ade9eSPeter Brune   if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);}
52111ade9eSPeter Brune   if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);}
53111ade9eSPeter Brune   if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);}
54b20c023fSPeter Brune 
55b20c023fSPeter Brune   nasm->eventrestrictinterp = 0;
56b20c023fSPeter Brune   nasm->eventsubsolve = 0;
57eaedb033SPeter Brune   PetscFunctionReturn(0);
58eaedb033SPeter Brune }
59eaedb033SPeter Brune 
60eaedb033SPeter Brune #undef __FUNCT__
61eaedb033SPeter Brune #define __FUNCT__ "SNESDestroy_NASM"
62eaedb033SPeter Brune PetscErrorCode SNESDestroy_NASM(SNES snes)
63eaedb033SPeter Brune {
64eaedb033SPeter Brune   PetscErrorCode ierr;
656e111a19SKarl Rupp 
66eaedb033SPeter Brune   PetscFunctionBegin;
67eaedb033SPeter Brune   ierr = SNESReset_NASM(snes);CHKERRQ(ierr);
6822d28d08SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
69eaedb033SPeter Brune   PetscFunctionReturn(0);
70eaedb033SPeter Brune }
71eaedb033SPeter Brune 
72eaedb033SPeter Brune #undef __FUNCT__
73111ade9eSPeter Brune #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private"
740adebc6cSBarry Smith PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx)
750adebc6cSBarry Smith {
76111ade9eSPeter Brune   PetscErrorCode ierr;
77111ade9eSPeter Brune   Vec            bcs = (Vec)ctx;
786e111a19SKarl Rupp 
79111ade9eSPeter Brune   PetscFunctionBegin;
80111ade9eSPeter Brune   ierr = VecCopy(bcs,l);CHKERRQ(ierr);
81111ade9eSPeter Brune   PetscFunctionReturn(0);
82111ade9eSPeter Brune }
83111ade9eSPeter Brune 
84111ade9eSPeter Brune #undef __FUNCT__
85eaedb033SPeter Brune #define __FUNCT__ "SNESSetUp_NASM"
86eaedb033SPeter Brune PetscErrorCode SNESSetUp_NASM(SNES snes)
87eaedb033SPeter Brune {
88eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
89eaedb033SPeter Brune   PetscErrorCode ierr;
90d387c056SBarry Smith   DM             dm;
91111ade9eSPeter Brune   DM             *subdms;
92111ade9eSPeter Brune   PetscInt       i;
93eaedb033SPeter Brune   const char     *optionsprefix;
94111ade9eSPeter Brune   Vec            F;
95eaedb033SPeter Brune 
96eaedb033SPeter Brune   PetscFunctionBegin;
97eaedb033SPeter Brune   if (!nasm->subsnes) {
98eaedb033SPeter Brune     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
990a696f66SPeter Brune     if (dm) {
100eaedb033SPeter Brune       nasm->usesdm = PETSC_TRUE;
1010298fd71SBarry Smith       ierr         = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr);
102*ce94432eSBarry Smith       if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
103111ade9eSPeter Brune       ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr);
104eaedb033SPeter Brune 
105eaedb033SPeter Brune       ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
106111ade9eSPeter Brune       ierr = PetscMalloc(nasm->n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr);
107111ade9eSPeter Brune 
108111ade9eSPeter Brune       for (i=0; i<nasm->n; i++) {
109cdb298fcSPeter Brune         ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr);
110cdb298fcSPeter Brune         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr);
111cdb298fcSPeter Brune         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr);
112cdb298fcSPeter Brune         ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr);
113cdb298fcSPeter Brune         ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr);
114111ade9eSPeter Brune         ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr);
115111ade9eSPeter Brune       }
116111ade9eSPeter Brune       ierr = PetscFree(subdms);CHKERRQ(ierr);
117*ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!");
118*ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!");
119111ade9eSPeter Brune   /* allocate the global vectors */
120111ade9eSPeter Brune   ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->x);CHKERRQ(ierr);
121111ade9eSPeter Brune   ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->xl);CHKERRQ(ierr);
122111ade9eSPeter Brune   ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->y);CHKERRQ(ierr);
123111ade9eSPeter Brune   ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->b);CHKERRQ(ierr);
124111ade9eSPeter Brune 
125111ade9eSPeter Brune   for (i=0; i<nasm->n; i++) {
126111ade9eSPeter Brune     DM subdm;
1270298fd71SBarry Smith     ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr);
128111ade9eSPeter Brune     ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);
129111ade9eSPeter Brune     ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);
130111ade9eSPeter Brune     ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);
131111ade9eSPeter Brune     ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr);
132111ade9eSPeter Brune     ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr);
1330298fd71SBarry Smith     ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr);
134111ade9eSPeter Brune   }
135eaedb033SPeter Brune   PetscFunctionReturn(0);
136eaedb033SPeter Brune }
137eaedb033SPeter Brune 
138eaedb033SPeter Brune #undef __FUNCT__
139eaedb033SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NASM"
140eaedb033SPeter Brune PetscErrorCode SNESSetFromOptions_NASM(SNES snes)
141eaedb033SPeter Brune {
142eaedb033SPeter Brune   PetscErrorCode    ierr;
143111ade9eSPeter Brune   PCASMType         asmtype;
144b20c023fSPeter Brune   PetscBool         flg,monflg;
145111ade9eSPeter Brune   SNES_NASM         *nasm = (SNES_NASM*)snes->data;
1466e111a19SKarl Rupp 
147eaedb033SPeter Brune   PetscFunctionBegin;
148111ade9eSPeter Brune   ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr);
149111ade9eSPeter Brune   ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
1501aa26658SKarl Rupp   if (flg) nasm->type = asmtype;
151b20c023fSPeter Brune   flg    = PETSC_FALSE;
152b20c023fSPeter Brune   monflg = PETSC_TRUE;
153b20c023fSPeter Brune   ierr   = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr);
154b20c023fSPeter Brune   if (flg) {
155b20c023fSPeter Brune     ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr);
156b20c023fSPeter Brune     ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr);
157b20c023fSPeter Brune   }
158eaedb033SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
159eaedb033SPeter Brune   PetscFunctionReturn(0);
160eaedb033SPeter Brune }
161eaedb033SPeter Brune 
162eaedb033SPeter Brune #undef __FUNCT__
163eaedb033SPeter Brune #define __FUNCT__ "SNESView_NASM"
164eaedb033SPeter Brune PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
165eaedb033SPeter Brune {
166eaedb033SPeter Brune   PetscFunctionBegin;
167b20c023fSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
168b20c023fSPeter Brune   PetscErrorCode ierr;
169b20c023fSPeter Brune   PetscMPIInt    rank;
170b20c023fSPeter Brune   PetscInt       i,N;
171b20c023fSPeter Brune   PetscBool      iascii,isstring;
172b20c023fSPeter Brune   PetscViewer    sviewer;
173*ce94432eSBarry Smith   MPI_Comm       comm;
174b20c023fSPeter Brune 
175b20c023fSPeter Brune   PetscFunctionBegin;
176*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
177b20c023fSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
178b20c023fSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
179b20c023fSPeter Brune   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
180b20c023fSPeter Brune   ierr = MPI_Reduce(&nasm->n,&N,1,MPIU_INT,MPIU_SUM,0,comm);CHKERRQ(ierr);
181b20c023fSPeter Brune   if (iascii) {
182b20c023fSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  Nonlinear Additive Schwarz: total subdomain blocks = %D\n",N);CHKERRQ(ierr);
183b20c023fSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  Nonlinear Additive Schwarz: restriction/interpolation type - %s\n",SNESNASMTypes[nasm->type]);CHKERRQ(ierr);
184b20c023fSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  Nonlinear Additive Schwarz: subSNES iterations: %D subKSP iterations: %D\n",SNESNASMTypes[nasm->type]);CHKERRQ(ierr);
185b20c023fSPeter Brune     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
186b20c023fSPeter Brune     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr);
187b20c023fSPeter Brune     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
188b20c023fSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  Local SNES objects:\n");CHKERRQ(ierr);
189b20c023fSPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
190b20c023fSPeter Brune     if (!rank) {
191b20c023fSPeter Brune       for (i=0; i<nasm->n; i++) {
192b20c023fSPeter Brune         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
193b20c023fSPeter Brune         ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr);
194b20c023fSPeter Brune         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
195b20c023fSPeter Brune       }
196b20c023fSPeter Brune     }
197b20c023fSPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
198b20c023fSPeter Brune     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
199b20c023fSPeter Brune     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
200b20c023fSPeter Brune   } else if (isstring) {
201b20c023fSPeter Brune     ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr);
202b20c023fSPeter Brune     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
203b20c023fSPeter Brune     if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);}
204b20c023fSPeter Brune     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
205b20c023fSPeter Brune   }
206eaedb033SPeter Brune   PetscFunctionReturn(0);
207eaedb033SPeter Brune }
208eaedb033SPeter Brune 
209eaedb033SPeter Brune #undef __FUNCT__
210eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains"
211a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
212a6dfd86eSKarl Rupp {
213eaedb033SPeter Brune   PetscErrorCode ierr;
214111ade9eSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*);
2156e111a19SKarl Rupp 
216eaedb033SPeter Brune   PetscFunctionBegin;
217eaedb033SPeter Brune   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr);
218111ade9eSPeter Brune   ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);
219eaedb033SPeter Brune   PetscFunctionReturn(0);
220eaedb033SPeter Brune }
221eaedb033SPeter Brune 
222eaedb033SPeter Brune EXTERN_C_BEGIN
223eaedb033SPeter Brune #undef __FUNCT__
224eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains_NASM"
225a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
226a6dfd86eSKarl Rupp {
227eaedb033SPeter Brune   PetscInt       i;
228eaedb033SPeter Brune   PetscErrorCode ierr;
229eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
2306e111a19SKarl Rupp 
231eaedb033SPeter Brune   PetscFunctionBegin;
232*ce94432eSBarry Smith   if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
233eaedb033SPeter Brune 
234111ade9eSPeter Brune   /* tear down the previously set things */
235111ade9eSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
236111ade9eSPeter Brune 
237eaedb033SPeter Brune   nasm->n = n;
238111ade9eSPeter Brune   if (oscatter) {
239111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);}
240eaedb033SPeter Brune   }
241111ade9eSPeter Brune   if (iscatter) {
242111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);}
243eaedb033SPeter Brune   }
244111ade9eSPeter Brune   if (gscatter) {
245111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);}
246111ade9eSPeter Brune   }
247111ade9eSPeter Brune   if (oscatter) {
248111ade9eSPeter Brune     ierr = PetscMalloc(n*sizeof(IS),&nasm->oscatter);CHKERRQ(ierr);
249eaedb033SPeter Brune     for (i=0; i<n; i++) {
250111ade9eSPeter Brune       nasm->oscatter[i] = oscatter[i];
251eaedb033SPeter Brune     }
252111ade9eSPeter Brune   }
253111ade9eSPeter Brune   if (iscatter) {
254111ade9eSPeter Brune     ierr = PetscMalloc(n*sizeof(IS),&nasm->iscatter);CHKERRQ(ierr);
255eaedb033SPeter Brune     for (i=0; i<n; i++) {
256111ade9eSPeter Brune       nasm->iscatter[i] = iscatter[i];
257eaedb033SPeter Brune     }
258eaedb033SPeter Brune   }
259111ade9eSPeter Brune   if (gscatter) {
260111ade9eSPeter Brune     ierr = PetscMalloc(n*sizeof(IS),&nasm->gscatter);CHKERRQ(ierr);
261eaedb033SPeter Brune     for (i=0; i<n; i++) {
262111ade9eSPeter Brune       nasm->gscatter[i] = gscatter[i];
263eaedb033SPeter Brune     }
264eaedb033SPeter Brune   }
265111ade9eSPeter Brune 
266eaedb033SPeter Brune   if (subsnes) {
267eaedb033SPeter Brune     ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr);
268eaedb033SPeter Brune     for (i=0; i<n; i++) {
269eaedb033SPeter Brune       nasm->subsnes[i] = subsnes[i];
270eaedb033SPeter Brune     }
271eaedb033SPeter Brune   }
272eaedb033SPeter Brune   PetscFunctionReturn(0);
273eaedb033SPeter Brune }
274eaedb033SPeter Brune EXTERN_C_END
275eaedb033SPeter Brune 
276eaedb033SPeter Brune #undef __FUNCT__
277eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSolveLocal_Private"
2780adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X)
2790adebc6cSBarry Smith {
280eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
281258e1594SPeter Brune   SNES           subsnes;
282eaedb033SPeter Brune   PetscInt       i;
283eaedb033SPeter Brune   PetscErrorCode ierr;
284111ade9eSPeter Brune   Vec            Xlloc,Xl,Bl,Yl;
285111ade9eSPeter Brune   VecScatter     iscat,oscat,gscat;
286111ade9eSPeter Brune   DM             dm,subdm;
2870adebc6cSBarry Smith 
288eaedb033SPeter Brune   PetscFunctionBegin;
289eaedb033SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
290111ade9eSPeter Brune   ierr = VecSet(Y,0);CHKERRQ(ierr);
291b20c023fSPeter Brune 
292b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
293eaedb033SPeter Brune   for (i=0; i<nasm->n; i++) {
29470c78f05SPeter Brune     /* scatter the solution to the local solution */
29570c78f05SPeter Brune     Xlloc = nasm->xl[i];
29670c78f05SPeter Brune     gscat   = nasm->gscatter[i];
29770c78f05SPeter Brune     oscat   = nasm->oscatter[i];
29870c78f05SPeter Brune     ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
29970c78f05SPeter Brune     if (B) {
30070c78f05SPeter Brune       /* scatter the RHS to the local RHS */
30170c78f05SPeter Brune       Bl   = nasm->b[i];
30270c78f05SPeter Brune       ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
30370c78f05SPeter Brune     }
30470c78f05SPeter Brune   }
305b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
306b20c023fSPeter Brune 
30770c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
30870c78f05SPeter Brune     Xlloc = nasm->xl[i];
309d590fa63SPeter Brune     gscat   = nasm->gscatter[i];
310d590fa63SPeter Brune     oscat   = nasm->oscatter[i];
31170c78f05SPeter Brune     ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31270c78f05SPeter Brune     if (B) {
31324b7f281SPeter Brune       Bl   = nasm->b[i];
31470c78f05SPeter Brune       ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
31570c78f05SPeter Brune     }
31670c78f05SPeter Brune   }
317b20c023fSPeter Brune 
318b20c023fSPeter Brune   if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
31970c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
32070c78f05SPeter Brune     Xl    = nasm->x[i];
32170c78f05SPeter Brune     Xlloc = nasm->xl[i];
32270c78f05SPeter Brune     Yl    = nasm->y[i];
323258e1594SPeter Brune     subsnes = nasm->subsnes[i];
324258e1594SPeter Brune     ierr    = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
325111ade9eSPeter Brune     iscat   = nasm->iscatter[i];
326111ade9eSPeter Brune     oscat   = nasm->oscatter[i];
327111ade9eSPeter Brune     gscat   = nasm->gscatter[i];
328111ade9eSPeter Brune     ierr    = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
32924b7f281SPeter Brune     if (B) {
33024b7f281SPeter Brune       Bl = nasm->b[i];
33124b7f281SPeter Brune     } else {
33224b7f281SPeter Brune       Bl = NULL;
33324b7f281SPeter Brune     }
33470c78f05SPeter Brune     ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
33570c78f05SPeter Brune     ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
336111ade9eSPeter Brune     ierr = VecCopy(Xl,Yl);CHKERRQ(ierr);
337258e1594SPeter Brune     ierr = SNESSolve(subsnes,Bl,Yl);CHKERRQ(ierr);
338111ade9eSPeter Brune     ierr = VecAXPY(Yl,-1.0,Xl);CHKERRQ(ierr);
33970c78f05SPeter Brune   }
340b20c023fSPeter Brune   if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
341111ade9eSPeter Brune 
342*ce94432eSBarry Smith   ierr = MPI_Barrier(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
34370c78f05SPeter Brune 
344b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
34570c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
34670c78f05SPeter Brune     Yl    = nasm->y[i];
34770c78f05SPeter Brune     iscat   = nasm->iscatter[i];
34870c78f05SPeter Brune     oscat   = nasm->oscatter[i];
349111ade9eSPeter Brune     if (nasm->type == PC_ASM_BASIC) {
350111ade9eSPeter Brune       ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
351111ade9eSPeter Brune     } else if (nasm->type == PC_ASM_RESTRICT) {
352111ade9eSPeter Brune       ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
353*ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
354eaedb033SPeter Brune   }
355eaedb033SPeter Brune 
35670c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
35770c78f05SPeter Brune     Yl    = nasm->y[i];
35870c78f05SPeter Brune     iscat   = nasm->iscatter[i];
35970c78f05SPeter Brune     oscat   = nasm->oscatter[i];
36070c78f05SPeter Brune     if (nasm->type == PC_ASM_BASIC) {
36170c78f05SPeter Brune       ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
36270c78f05SPeter Brune     } else if (nasm->type == PC_ASM_RESTRICT) {
36370c78f05SPeter Brune       ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
364*ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
36570c78f05SPeter Brune   }
366b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
36770c78f05SPeter Brune 
368*ce94432eSBarry Smith   ierr = MPI_Barrier(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
369cd939e56SPeter Brune 
370111ade9eSPeter Brune   ierr = VecAXPY(X,1.0,Y);CHKERRQ(ierr);
371eaedb033SPeter Brune   PetscFunctionReturn(0);
372eaedb033SPeter Brune }
373eaedb033SPeter Brune 
374eaedb033SPeter Brune #undef __FUNCT__
375eaedb033SPeter Brune #define __FUNCT__ "SNESSolve_NASM"
376eaedb033SPeter Brune PetscErrorCode SNESSolve_NASM(SNES snes)
377eaedb033SPeter Brune {
378eaedb033SPeter Brune   Vec            F;
379eaedb033SPeter Brune   Vec            X;
380eaedb033SPeter Brune   Vec            B;
381111ade9eSPeter Brune   Vec            Y;
382eaedb033SPeter Brune   PetscInt       i;
383eaedb033SPeter Brune   PetscReal      fnorm;
384eaedb033SPeter Brune   PetscErrorCode ierr;
385eaedb033SPeter Brune   SNESNormType   normtype;
386eaedb033SPeter Brune 
387eaedb033SPeter Brune   PetscFunctionBegin;
388eaedb033SPeter Brune   X = snes->vec_sol;
389111ade9eSPeter Brune   Y = snes->vec_sol_update;
390eaedb033SPeter Brune   F = snes->vec_func;
391eaedb033SPeter Brune   B = snes->vec_rhs;
392eaedb033SPeter Brune 
393eaedb033SPeter Brune   ierr         = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
394eaedb033SPeter Brune   snes->iter   = 0;
395eaedb033SPeter Brune   snes->norm   = 0.;
396eaedb033SPeter Brune   ierr         = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
397eaedb033SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
398eaedb033SPeter Brune   ierr         = SNESGetNormType(snes, &normtype);CHKERRQ(ierr);
399eaedb033SPeter Brune   if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
400eaedb033SPeter Brune     /* compute the initial function and preconditioned update delX */
401eaedb033SPeter Brune     if (!snes->vec_func_init_set) {
402eaedb033SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
403eaedb033SPeter Brune       if (snes->domainerror) {
404eaedb033SPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
405eaedb033SPeter Brune         PetscFunctionReturn(0);
406eaedb033SPeter Brune       }
4071aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
408eaedb033SPeter Brune 
409eaedb033SPeter Brune     /* convergence test */
410eaedb033SPeter Brune     if (!snes->norm_init_set) {
411eaedb033SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
412eaedb033SPeter Brune       if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
413eaedb033SPeter Brune     } else {
414eaedb033SPeter Brune       fnorm               = snes->norm_init;
415eaedb033SPeter Brune       snes->norm_init_set = PETSC_FALSE;
416eaedb033SPeter Brune     }
417eaedb033SPeter Brune     ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
418eaedb033SPeter Brune     snes->iter = 0;
419eaedb033SPeter Brune     snes->norm = fnorm;
420eaedb033SPeter Brune     ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
421a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
422eaedb033SPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
423eaedb033SPeter Brune 
424eaedb033SPeter Brune     /* set parameter for default relative tolerance convergence test */
425eaedb033SPeter Brune     snes->ttol = fnorm*snes->rtol;
426eaedb033SPeter Brune 
427eaedb033SPeter Brune     /* test convergence */
428eaedb033SPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
429eaedb033SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
430eaedb033SPeter Brune   } else {
431eaedb033SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
432a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
433eaedb033SPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
434eaedb033SPeter Brune   }
435eaedb033SPeter Brune 
436eaedb033SPeter Brune   /* Call general purpose update function */
437eaedb033SPeter Brune   if (snes->ops->update) {
438eaedb033SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
439eaedb033SPeter Brune   }
440eaedb033SPeter Brune 
441eaedb033SPeter Brune   for (i = 0; i < snes->max_its; i++) {
442111ade9eSPeter Brune     ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr);
443eaedb033SPeter Brune     if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) {
444eaedb033SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
445eaedb033SPeter Brune       if (snes->domainerror) {
446eaedb033SPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
447eaedb033SPeter Brune         PetscFunctionReturn(0);
448eaedb033SPeter Brune       }
449eaedb033SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
450eaedb033SPeter Brune       if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
451eaedb033SPeter Brune     }
452eaedb033SPeter Brune     /* Monitor convergence */
453eaedb033SPeter Brune     ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
454eaedb033SPeter Brune     snes->iter = i+1;
455eaedb033SPeter Brune     snes->norm = fnorm;
456eaedb033SPeter Brune     ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
457a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
458eaedb033SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
459eaedb033SPeter Brune     /* Test for convergence */
460eaedb033SPeter Brune     if (normtype == SNES_NORM_FUNCTION) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
461eaedb033SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
462eaedb033SPeter Brune     /* Call general purpose update function */
463eaedb033SPeter Brune     if (snes->ops->update) {
464eaedb033SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
465eaedb033SPeter Brune     }
466eaedb033SPeter Brune   }
467eaedb033SPeter Brune   if (normtype == SNES_NORM_FUNCTION) {
468eaedb033SPeter Brune     if (i == snes->max_its) {
469eaedb033SPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
470eaedb033SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
471eaedb033SPeter Brune     }
4721aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */
473eaedb033SPeter Brune   PetscFunctionReturn(0);
474eaedb033SPeter Brune }
475eaedb033SPeter Brune 
476eaedb033SPeter Brune /*MC
477eaedb033SPeter Brune   SNESNASM - Nonlinear Additive Schwartz
478eaedb033SPeter Brune 
479eaedb033SPeter Brune    Level: advanced
480eaedb033SPeter Brune 
481eaedb033SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
482eaedb033SPeter Brune M*/
483eaedb033SPeter Brune 
484eaedb033SPeter Brune EXTERN_C_BEGIN
485eaedb033SPeter Brune #undef __FUNCT__
486eaedb033SPeter Brune #define __FUNCT__ "SNESCreate_NASM"
487eaedb033SPeter Brune PetscErrorCode SNESCreate_NASM(SNES snes)
488eaedb033SPeter Brune {
489eaedb033SPeter Brune   SNES_NASM      *nasm;
490eaedb033SPeter Brune   PetscErrorCode ierr;
491eaedb033SPeter Brune 
492eaedb033SPeter Brune   PetscFunctionBegin;
493eaedb033SPeter Brune   ierr       = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr);
494eaedb033SPeter Brune   snes->data = (void*)nasm;
495eaedb033SPeter Brune 
496eaedb033SPeter Brune   nasm->n        = PETSC_DECIDE;
497eaedb033SPeter Brune   nasm->subsnes  = 0;
498eaedb033SPeter Brune   nasm->x        = 0;
499111ade9eSPeter Brune   nasm->xl       = 0;
500111ade9eSPeter Brune   nasm->y        = 0;
501eaedb033SPeter Brune   nasm->b        = 0;
502111ade9eSPeter Brune   nasm->oscatter = 0;
503111ade9eSPeter Brune   nasm->iscatter = 0;
504111ade9eSPeter Brune   nasm->gscatter = 0;
505111ade9eSPeter Brune 
506111ade9eSPeter Brune   nasm->type = PC_ASM_BASIC;
507eaedb033SPeter Brune 
508eaedb033SPeter Brune   snes->ops->destroy        = SNESDestroy_NASM;
509eaedb033SPeter Brune   snes->ops->setup          = SNESSetUp_NASM;
510eaedb033SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
511eaedb033SPeter Brune   snes->ops->view           = SNESView_NASM;
512eaedb033SPeter Brune   snes->ops->solve          = SNESSolve_NASM;
513eaedb033SPeter Brune   snes->ops->reset          = SNESReset_NASM;
514eaedb033SPeter Brune 
515eaedb033SPeter Brune   snes->usesksp = PETSC_FALSE;
516eaedb033SPeter Brune   snes->usespc  = PETSC_FALSE;
517eaedb033SPeter Brune 
5180298fd71SBarry Smith   nasm->eventrestrictinterp = 0;
5190298fd71SBarry Smith   nasm->eventsubsolve       = 0;
520b20c023fSPeter Brune 
521eaedb033SPeter Brune   if (!snes->tolerancesset) {
522eaedb033SPeter Brune     snes->max_its   = 10000;
523eaedb033SPeter Brune     snes->max_funcs = 10000;
524eaedb033SPeter Brune   }
525eaedb033SPeter Brune 
526eaedb033SPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNASMSetSubdomains_C","SNESNASMSetSubdomains_NASM",
527eaedb033SPeter Brune                                            SNESNASMSetSubdomains_NASM);CHKERRQ(ierr);
528eaedb033SPeter Brune   PetscFunctionReturn(0);
529eaedb033SPeter Brune }
530eaedb033SPeter Brune EXTERN_C_END
531