xref: /petsc/src/snes/impls/nasm/nasm.c (revision 602bec5d85c6d2636ed66a5e71fca98e0385570f)
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 */
16d728fb7dSPeter Brune   PetscBool  finaljacobian;       /* compute the jacobian of the converged solution */
17b20c023fSPeter Brune 
18b20c023fSPeter Brune   /* logging events */
19b20c023fSPeter Brune   PetscLogEvent eventrestrictinterp;
20b20c023fSPeter Brune   PetscLogEvent eventsubsolve;
21*602bec5dSPeter Brune 
22*602bec5dSPeter Brune   PetscInt      fjtype;            /* type of computed jacobian */
23*602bec5dSPeter Brune   Vec           xinit;             /* initial solution in case the final jacobian type is computed as first */
24eaedb033SPeter Brune } SNES_NASM;
25eaedb033SPeter Brune 
26b20c023fSPeter Brune const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0};
27*602bec5dSPeter Brune const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"};
28b20c023fSPeter Brune 
29eaedb033SPeter Brune #undef __FUNCT__
30eaedb033SPeter Brune #define __FUNCT__ "SNESReset_NASM"
31eaedb033SPeter Brune PetscErrorCode SNESReset_NASM(SNES snes)
32eaedb033SPeter Brune {
33eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
34eaedb033SPeter Brune   PetscErrorCode ierr;
35eaedb033SPeter Brune   PetscInt       i;
366e111a19SKarl Rupp 
37eaedb033SPeter Brune   PetscFunctionBegin;
38eaedb033SPeter Brune   for (i=0; i<nasm->n; i++) {
39111ade9eSPeter Brune     if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); }
40f5f7c1b9SKarl Rupp     if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); }
41111ade9eSPeter Brune     if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); }
42bc8c1f72SJose Roman     if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); }
43eaedb033SPeter Brune 
44bc8c1f72SJose Roman     if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); }
45111ade9eSPeter Brune     if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); }
46111ade9eSPeter Brune     if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); }
47111ade9eSPeter Brune     if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); }
48eaedb033SPeter Brune   }
49111ade9eSPeter Brune 
50111ade9eSPeter Brune   if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);}
51111ade9eSPeter Brune   if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);}
52111ade9eSPeter Brune   if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);}
53111ade9eSPeter Brune   if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);}
54111ade9eSPeter Brune 
55*602bec5dSPeter Brune   if (nasm->xinit) {ierr = VecDestroy(&nasm->xinit);CHKERRQ(ierr);}
56*602bec5dSPeter Brune 
57111ade9eSPeter Brune   if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);}
58111ade9eSPeter Brune   if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);}
59111ade9eSPeter Brune   if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);}
60111ade9eSPeter Brune   if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);}
61b20c023fSPeter Brune 
62b20c023fSPeter Brune   nasm->eventrestrictinterp = 0;
63b20c023fSPeter Brune   nasm->eventsubsolve = 0;
64eaedb033SPeter Brune   PetscFunctionReturn(0);
65eaedb033SPeter Brune }
66eaedb033SPeter Brune 
67eaedb033SPeter Brune #undef __FUNCT__
68eaedb033SPeter Brune #define __FUNCT__ "SNESDestroy_NASM"
69eaedb033SPeter Brune PetscErrorCode SNESDestroy_NASM(SNES snes)
70eaedb033SPeter Brune {
71eaedb033SPeter Brune   PetscErrorCode ierr;
726e111a19SKarl Rupp 
73eaedb033SPeter Brune   PetscFunctionBegin;
74eaedb033SPeter Brune   ierr = SNESReset_NASM(snes);CHKERRQ(ierr);
7522d28d08SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
76eaedb033SPeter Brune   PetscFunctionReturn(0);
77eaedb033SPeter Brune }
78eaedb033SPeter Brune 
79eaedb033SPeter Brune #undef __FUNCT__
80111ade9eSPeter Brune #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private"
810adebc6cSBarry Smith PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx)
820adebc6cSBarry Smith {
83111ade9eSPeter Brune   PetscErrorCode ierr;
84111ade9eSPeter Brune   Vec            bcs = (Vec)ctx;
856e111a19SKarl Rupp 
86111ade9eSPeter Brune   PetscFunctionBegin;
87111ade9eSPeter Brune   ierr = VecCopy(bcs,l);CHKERRQ(ierr);
88111ade9eSPeter Brune   PetscFunctionReturn(0);
89111ade9eSPeter Brune }
90111ade9eSPeter Brune 
91111ade9eSPeter Brune #undef __FUNCT__
92eaedb033SPeter Brune #define __FUNCT__ "SNESSetUp_NASM"
93eaedb033SPeter Brune PetscErrorCode SNESSetUp_NASM(SNES snes)
94eaedb033SPeter Brune {
95eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
96eaedb033SPeter Brune   PetscErrorCode ierr;
9776857b2aSPeter Brune   DM             dm,subdm;
98111ade9eSPeter Brune   DM             *subdms;
99111ade9eSPeter Brune   PetscInt       i;
100eaedb033SPeter Brune   const char     *optionsprefix;
101111ade9eSPeter Brune   Vec            F;
102ed3c11a9SPeter Brune   PetscMPIInt    size;
103ed3c11a9SPeter Brune   KSP            ksp;
104ed3c11a9SPeter Brune   PC             pc;
105eaedb033SPeter Brune 
106eaedb033SPeter Brune   PetscFunctionBegin;
107eaedb033SPeter Brune   if (!nasm->subsnes) {
108eaedb033SPeter Brune     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1090a696f66SPeter Brune     if (dm) {
110eaedb033SPeter Brune       nasm->usesdm = PETSC_TRUE;
1110298fd71SBarry Smith       ierr         = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr);
112ce94432eSBarry Smith       if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
113111ade9eSPeter Brune       ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr);
114eaedb033SPeter Brune 
115eaedb033SPeter Brune       ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
116111ade9eSPeter Brune       ierr = PetscMalloc(nasm->n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr);
117111ade9eSPeter Brune       for (i=0; i<nasm->n; i++) {
118cdb298fcSPeter Brune         ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr);
119cdb298fcSPeter Brune         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr);
120cdb298fcSPeter Brune         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr);
121cdb298fcSPeter Brune         ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr);
122ed3c11a9SPeter Brune         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRQ(ierr);
123ed3c11a9SPeter Brune         if (size == 1) {
124ed3c11a9SPeter Brune           ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr);
125ed3c11a9SPeter Brune           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
126ed3c11a9SPeter Brune           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
127ed3c11a9SPeter Brune           ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
128ed3c11a9SPeter Brune         }
129cdb298fcSPeter Brune         ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr);
130111ade9eSPeter Brune         ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr);
131111ade9eSPeter Brune       }
132111ade9eSPeter Brune       ierr = PetscFree(subdms);CHKERRQ(ierr);
133ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!");
134ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!");
135111ade9eSPeter Brune   /* allocate the global vectors */
136e245e0aaSPeter Brune   if (!nasm->x) {
137111ade9eSPeter Brune     ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->x);CHKERRQ(ierr);
138e245e0aaSPeter Brune     ierr = PetscMemzero(nasm->x,nasm->n*sizeof(Vec));CHKERRQ(ierr);
139e245e0aaSPeter Brune   }
140e245e0aaSPeter Brune   if (!nasm->xl) {
141111ade9eSPeter Brune     ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->xl);CHKERRQ(ierr);
142e245e0aaSPeter Brune     ierr = PetscMemzero(nasm->xl,nasm->n*sizeof(Vec));CHKERRQ(ierr);
143e245e0aaSPeter Brune   }
144e245e0aaSPeter Brune   if (!nasm->y) {
145111ade9eSPeter Brune     ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->y);CHKERRQ(ierr);
146e245e0aaSPeter Brune     ierr = PetscMemzero(nasm->y,nasm->n*sizeof(Vec));CHKERRQ(ierr);
147e245e0aaSPeter Brune   }
148e245e0aaSPeter Brune   if (!nasm->b) {
149111ade9eSPeter Brune     ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->b);CHKERRQ(ierr);
150e245e0aaSPeter Brune     ierr = PetscMemzero(nasm->b,nasm->n*sizeof(Vec));CHKERRQ(ierr);
151e245e0aaSPeter Brune   }
152111ade9eSPeter Brune 
153111ade9eSPeter Brune   for (i=0; i<nasm->n; i++) {
1540298fd71SBarry Smith     ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr);
15576857b2aSPeter Brune     if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);}
15676857b2aSPeter Brune     if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);}
15776857b2aSPeter Brune     if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);}
15876857b2aSPeter Brune     if (!nasm->xl[i]) {
159111ade9eSPeter Brune       ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr);
160111ade9eSPeter Brune       ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr);
16176857b2aSPeter Brune     }
1620298fd71SBarry Smith     ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr);
163111ade9eSPeter Brune   }
164*602bec5dSPeter Brune   if (nasm->finaljacobian) {
165*602bec5dSPeter Brune     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
166*602bec5dSPeter Brune     if (nasm->fjtype == 2) {
167*602bec5dSPeter Brune       ierr = VecDuplicate(snes->vec_sol,&nasm->xinit);CHKERRQ(ierr);
168*602bec5dSPeter Brune     }
169*602bec5dSPeter Brune     for (i=0; i<nasm->n;i++) {
170*602bec5dSPeter Brune       ierr = SNESSetUpMatrices(nasm->subsnes[i]);CHKERRQ(ierr);
171*602bec5dSPeter Brune     }
172*602bec5dSPeter Brune   }
173eaedb033SPeter Brune   PetscFunctionReturn(0);
174eaedb033SPeter Brune }
175eaedb033SPeter Brune 
176eaedb033SPeter Brune #undef __FUNCT__
177eaedb033SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NASM"
178eaedb033SPeter Brune PetscErrorCode SNESSetFromOptions_NASM(SNES snes)
179eaedb033SPeter Brune {
180eaedb033SPeter Brune   PetscErrorCode    ierr;
181111ade9eSPeter Brune   PCASMType         asmtype;
182b20c023fSPeter Brune   PetscBool         flg,monflg;
183111ade9eSPeter Brune   SNES_NASM         *nasm = (SNES_NASM*)snes->data;
1846e111a19SKarl Rupp 
185eaedb033SPeter Brune   PetscFunctionBegin;
186111ade9eSPeter Brune   ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr);
187111ade9eSPeter Brune   ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
1881aa26658SKarl Rupp   if (flg) nasm->type = asmtype;
189b20c023fSPeter Brune   flg    = PETSC_FALSE;
190b20c023fSPeter Brune   monflg = PETSC_TRUE;
191b20c023fSPeter Brune   ierr   = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr);
192d728fb7dSPeter Brune   ierr   = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr);
193*602bec5dSPeter Brune   ierr   = PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL);CHKERRQ(ierr);
194b20c023fSPeter Brune   if (flg) {
195b20c023fSPeter Brune     ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr);
196b20c023fSPeter Brune     ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr);
197b20c023fSPeter Brune   }
198eaedb033SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
199eaedb033SPeter Brune   PetscFunctionReturn(0);
200eaedb033SPeter Brune }
201eaedb033SPeter Brune 
202eaedb033SPeter Brune #undef __FUNCT__
203eaedb033SPeter Brune #define __FUNCT__ "SNESView_NASM"
204eaedb033SPeter Brune PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
205eaedb033SPeter Brune {
206b20c023fSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
207b20c023fSPeter Brune   PetscErrorCode ierr;
208b20c023fSPeter Brune   PetscMPIInt    rank;
209b20c023fSPeter Brune   PetscInt       i,N;
210b20c023fSPeter Brune   PetscBool      iascii,isstring;
211b20c023fSPeter Brune   PetscViewer    sviewer;
212ce94432eSBarry Smith   MPI_Comm       comm;
213b20c023fSPeter Brune 
214b20c023fSPeter Brune   PetscFunctionBegin;
215ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
216b20c023fSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
217b20c023fSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
218b20c023fSPeter Brune   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
219b20c023fSPeter Brune   ierr = MPI_Reduce(&nasm->n,&N,1,MPIU_INT,MPIU_SUM,0,comm);CHKERRQ(ierr);
220b20c023fSPeter Brune   if (iascii) {
221b20c023fSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  Nonlinear Additive Schwarz: total subdomain blocks = %D\n",N);CHKERRQ(ierr);
222b20c023fSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  Nonlinear Additive Schwarz: restriction/interpolation type - %s\n",SNESNASMTypes[nasm->type]);CHKERRQ(ierr);
223b20c023fSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  Nonlinear Additive Schwarz: subSNES iterations: %D subKSP iterations: %D\n",SNESNASMTypes[nasm->type]);CHKERRQ(ierr);
224b20c023fSPeter Brune     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
225b20c023fSPeter Brune     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr);
226b20c023fSPeter Brune     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
227b20c023fSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  Local SNES objects:\n");CHKERRQ(ierr);
228b20c023fSPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
229b20c023fSPeter Brune     if (!rank) {
230b20c023fSPeter Brune       for (i=0; i<nasm->n; i++) {
231b20c023fSPeter Brune         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
232b20c023fSPeter Brune         ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr);
233b20c023fSPeter Brune         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
234b20c023fSPeter Brune       }
235b20c023fSPeter Brune     }
236b20c023fSPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
237b20c023fSPeter Brune     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
238b20c023fSPeter Brune     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
239b20c023fSPeter Brune   } else if (isstring) {
240b20c023fSPeter Brune     ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr);
241b20c023fSPeter Brune     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
242b20c023fSPeter Brune     if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);}
243b20c023fSPeter Brune     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
244b20c023fSPeter Brune   }
245eaedb033SPeter Brune   PetscFunctionReturn(0);
246eaedb033SPeter Brune }
247eaedb033SPeter Brune 
248eaedb033SPeter Brune #undef __FUNCT__
249eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains"
25076857b2aSPeter Brune /*@
25176857b2aSPeter Brune    SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems.
25276857b2aSPeter Brune 
25376857b2aSPeter Brune    Not Collective
25476857b2aSPeter Brune 
25576857b2aSPeter Brune    Input Parameters:
25676857b2aSPeter Brune +  SNES - the SNES context
25776857b2aSPeter Brune .  n - the number of local subdomains
25876857b2aSPeter Brune .  subsnes - solvers defined on the local subdomains
25976857b2aSPeter Brune .  iscatter - scatters into the nonoverlapping portions of the local subdomains
26076857b2aSPeter Brune .  oscatter - scatters into the overlapping portions of the local subdomains
26176857b2aSPeter Brune -  gscatter - scatters into the (ghosted) local vector of the local subdomain
26276857b2aSPeter Brune 
26376857b2aSPeter Brune    Level: intermediate
26476857b2aSPeter Brune 
26576857b2aSPeter Brune .keywords: SNES, NASM
26676857b2aSPeter Brune 
26776857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains()
26876857b2aSPeter Brune @*/
269a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
270a6dfd86eSKarl Rupp {
271eaedb033SPeter Brune   PetscErrorCode ierr;
272111ade9eSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*);
2736e111a19SKarl Rupp 
274eaedb033SPeter Brune   PetscFunctionBegin;
275eaedb033SPeter Brune   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr);
276111ade9eSPeter Brune   ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);
277eaedb033SPeter Brune   PetscFunctionReturn(0);
278eaedb033SPeter Brune }
279eaedb033SPeter Brune 
280eaedb033SPeter Brune #undef __FUNCT__
281eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains_NASM"
282a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
283a6dfd86eSKarl Rupp {
284eaedb033SPeter Brune   PetscInt       i;
285eaedb033SPeter Brune   PetscErrorCode ierr;
286eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
2876e111a19SKarl Rupp 
288eaedb033SPeter Brune   PetscFunctionBegin;
289ce94432eSBarry Smith   if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
290eaedb033SPeter Brune 
291111ade9eSPeter Brune   /* tear down the previously set things */
292111ade9eSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
293111ade9eSPeter Brune 
294eaedb033SPeter Brune   nasm->n = n;
295111ade9eSPeter Brune   if (oscatter) {
296111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);}
297eaedb033SPeter Brune   }
298111ade9eSPeter Brune   if (iscatter) {
299111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);}
300eaedb033SPeter Brune   }
301111ade9eSPeter Brune   if (gscatter) {
302111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);}
303111ade9eSPeter Brune   }
304111ade9eSPeter Brune   if (oscatter) {
305111ade9eSPeter Brune     ierr = PetscMalloc(n*sizeof(IS),&nasm->oscatter);CHKERRQ(ierr);
306eaedb033SPeter Brune     for (i=0; i<n; i++) {
307111ade9eSPeter Brune       nasm->oscatter[i] = oscatter[i];
308eaedb033SPeter Brune     }
309111ade9eSPeter Brune   }
310111ade9eSPeter Brune   if (iscatter) {
311111ade9eSPeter Brune     ierr = PetscMalloc(n*sizeof(IS),&nasm->iscatter);CHKERRQ(ierr);
312eaedb033SPeter Brune     for (i=0; i<n; i++) {
313111ade9eSPeter Brune       nasm->iscatter[i] = iscatter[i];
314eaedb033SPeter Brune     }
315eaedb033SPeter Brune   }
316111ade9eSPeter Brune   if (gscatter) {
317111ade9eSPeter Brune     ierr = PetscMalloc(n*sizeof(IS),&nasm->gscatter);CHKERRQ(ierr);
318eaedb033SPeter Brune     for (i=0; i<n; i++) {
319111ade9eSPeter Brune       nasm->gscatter[i] = gscatter[i];
320eaedb033SPeter Brune     }
321eaedb033SPeter Brune   }
322111ade9eSPeter Brune 
323eaedb033SPeter Brune   if (subsnes) {
324eaedb033SPeter Brune     ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr);
325eaedb033SPeter Brune     for (i=0; i<n; i++) {
326eaedb033SPeter Brune       nasm->subsnes[i] = subsnes[i];
327eaedb033SPeter Brune     }
328eaedb033SPeter Brune   }
329eaedb033SPeter Brune   PetscFunctionReturn(0);
330eaedb033SPeter Brune }
331eaedb033SPeter Brune 
332eaedb033SPeter Brune #undef __FUNCT__
33376857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains"
33476857b2aSPeter Brune /*@
33576857b2aSPeter Brune    SNESNASMGetSubdomains - Get the local subdomain context.
33676857b2aSPeter Brune 
33776857b2aSPeter Brune    Not Collective
33876857b2aSPeter Brune 
33976857b2aSPeter Brune    Input Parameters:
34076857b2aSPeter Brune .  SNES - the SNES context
34176857b2aSPeter Brune 
34276857b2aSPeter Brune    Output Parameters:
34376857b2aSPeter Brune +  n - the number of local subdomains
34476857b2aSPeter Brune .  subsnes - solvers defined on the local subdomains
34576857b2aSPeter Brune .  iscatter - scatters into the nonoverlapping portions of the local subdomains
34676857b2aSPeter Brune .  oscatter - scatters into the overlapping portions of the local subdomains
34776857b2aSPeter Brune -  gscatter - scatters into the (ghosted) local vector of the local subdomain
34876857b2aSPeter Brune 
34976857b2aSPeter Brune    Level: intermediate
35076857b2aSPeter Brune 
35176857b2aSPeter Brune .keywords: SNES, NASM
35276857b2aSPeter Brune 
35376857b2aSPeter Brune .seealso: SNESNASM, SNESNASMSetSubdomains()
35476857b2aSPeter Brune @*/
35576857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
35676857b2aSPeter Brune {
35776857b2aSPeter Brune   PetscErrorCode ierr;
35876857b2aSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**);
35976857b2aSPeter Brune 
36076857b2aSPeter Brune   PetscFunctionBegin;
36176857b2aSPeter Brune   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr);
36276857b2aSPeter Brune   ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);
36376857b2aSPeter Brune   PetscFunctionReturn(0);
36476857b2aSPeter Brune }
36576857b2aSPeter Brune 
36676857b2aSPeter Brune #undef __FUNCT__
36776857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains_NASM"
36876857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
36976857b2aSPeter Brune {
37076857b2aSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
37176857b2aSPeter Brune 
37276857b2aSPeter Brune   PetscFunctionBegin;
37376857b2aSPeter Brune   if (n) *n = nasm->n;
37476857b2aSPeter Brune   if (oscatter) *oscatter = nasm->oscatter;
37576857b2aSPeter Brune   if (iscatter) *iscatter = nasm->iscatter;
37676857b2aSPeter Brune   if (gscatter) *gscatter = nasm->gscatter;
37776857b2aSPeter Brune   if (subsnes)  *subsnes  = nasm->subsnes;
37876857b2aSPeter Brune   PetscFunctionReturn(0);
37976857b2aSPeter Brune }
38076857b2aSPeter Brune 
38176857b2aSPeter Brune #undef __FUNCT__
38276857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs"
38376857b2aSPeter Brune /*@
38476857b2aSPeter Brune    SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors
38576857b2aSPeter Brune 
38676857b2aSPeter Brune    Not Collective
38776857b2aSPeter Brune 
38876857b2aSPeter Brune    Input Parameters:
38976857b2aSPeter Brune .  SNES - the SNES context
39076857b2aSPeter Brune 
39176857b2aSPeter Brune    Output Parameters:
39276857b2aSPeter Brune +  n - the number of local subdomains
39376857b2aSPeter Brune .  x - The subdomain solution vector
39476857b2aSPeter Brune .  y - The subdomain step vector
39576857b2aSPeter Brune .  b - The subdomain RHS vector
39676857b2aSPeter Brune -  xl - The subdomain local vectors (ghosted)
39776857b2aSPeter Brune 
39876857b2aSPeter Brune    Level: developer
39976857b2aSPeter Brune 
40076857b2aSPeter Brune .keywords: SNES, NASM
40176857b2aSPeter Brune 
40276857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains()
40376857b2aSPeter Brune @*/
40476857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl)
40576857b2aSPeter Brune {
40676857b2aSPeter Brune   PetscErrorCode ierr;
40776857b2aSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**);
40876857b2aSPeter Brune 
40976857b2aSPeter Brune   PetscFunctionBegin;
41076857b2aSPeter Brune   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",(void (**)(void))&f);CHKERRQ(ierr);
41176857b2aSPeter Brune   ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);
41276857b2aSPeter Brune   PetscFunctionReturn(0);
41376857b2aSPeter Brune }
41476857b2aSPeter Brune 
41576857b2aSPeter Brune #undef __FUNCT__
41676857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs_NASM"
41776857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl)
41876857b2aSPeter Brune {
41976857b2aSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
42076857b2aSPeter Brune 
42176857b2aSPeter Brune   PetscFunctionBegin;
42276857b2aSPeter Brune   if (n)  *n  = nasm->n;
42376857b2aSPeter Brune   if (x)  *x  = nasm->x;
42476857b2aSPeter Brune   if (y)  *y  = nasm->y;
42576857b2aSPeter Brune   if (b)  *b  = nasm->b;
42676857b2aSPeter Brune   if (xl) *xl = nasm->xl;
42776857b2aSPeter Brune   PetscFunctionReturn(0);
42876857b2aSPeter Brune }
42976857b2aSPeter Brune 
430d728fb7dSPeter Brune #undef __FUNCT__
431d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian"
432d728fb7dSPeter Brune /*@
433d728fb7dSPeter Brune    SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence
434d728fb7dSPeter Brune 
435d728fb7dSPeter Brune    Collective on SNES
436d728fb7dSPeter Brune 
437d728fb7dSPeter Brune    Input Parameters:
438d728fb7dSPeter Brune +  SNES - the SNES context
439d728fb7dSPeter Brune -  flg - indication of whether to compute the jacobians or not
440d728fb7dSPeter Brune 
441d728fb7dSPeter Brune    Level: developer
442d728fb7dSPeter Brune 
443d728fb7dSPeter Brune    Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian
444d728fb7dSPeter Brune    is needed at each linear iteration.
445d728fb7dSPeter Brune 
446d728fb7dSPeter Brune .keywords: SNES, NASM, ASPIN
447d728fb7dSPeter Brune 
448d728fb7dSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains()
449d728fb7dSPeter Brune @*/
450d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg)
451d728fb7dSPeter Brune {
452d728fb7dSPeter Brune   PetscErrorCode (*f)(SNES,PetscBool);
453d728fb7dSPeter Brune   PetscErrorCode ierr;
454d728fb7dSPeter Brune 
455d728fb7dSPeter Brune   PetscFunctionBegin;
456d728fb7dSPeter Brune   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",(void (**)(void))&f);CHKERRQ(ierr);
457d728fb7dSPeter Brune   ierr = (f)(snes,flg);CHKERRQ(ierr);
458d728fb7dSPeter Brune   PetscFunctionReturn(0);
459d728fb7dSPeter Brune }
460d728fb7dSPeter Brune 
461d728fb7dSPeter Brune #undef __FUNCT__
462d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian_NASM"
463d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg)
464d728fb7dSPeter Brune {
465d728fb7dSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
466d728fb7dSPeter Brune 
467d728fb7dSPeter Brune   PetscFunctionBegin;
468d728fb7dSPeter Brune   nasm->finaljacobian = flg;
469*602bec5dSPeter Brune   if (flg) snes->usesksp = PETSC_TRUE;
470d728fb7dSPeter Brune   PetscFunctionReturn(0);
471d728fb7dSPeter Brune }
47276857b2aSPeter Brune 
47376857b2aSPeter Brune #undef __FUNCT__
474eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSolveLocal_Private"
4750adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X)
4760adebc6cSBarry Smith {
477eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
478258e1594SPeter Brune   SNES           subsnes;
479eaedb033SPeter Brune   PetscInt       i;
480eaedb033SPeter Brune   PetscErrorCode ierr;
481111ade9eSPeter Brune   Vec            Xlloc,Xl,Bl,Yl;
482111ade9eSPeter Brune   VecScatter     iscat,oscat,gscat;
483111ade9eSPeter Brune   DM             dm,subdm;
4840adebc6cSBarry Smith 
485eaedb033SPeter Brune   PetscFunctionBegin;
486eaedb033SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
487111ade9eSPeter Brune   ierr = VecSet(Y,0);CHKERRQ(ierr);
488b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
489eaedb033SPeter Brune   for (i=0; i<nasm->n; i++) {
49070c78f05SPeter Brune     /* scatter the solution to the local solution */
49170c78f05SPeter Brune     Xlloc = nasm->xl[i];
49270c78f05SPeter Brune     gscat   = nasm->gscatter[i];
49370c78f05SPeter Brune     oscat   = nasm->oscatter[i];
49470c78f05SPeter Brune     ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
49570c78f05SPeter Brune     if (B) {
49670c78f05SPeter Brune       /* scatter the RHS to the local RHS */
49770c78f05SPeter Brune       Bl   = nasm->b[i];
49870c78f05SPeter Brune       ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
49970c78f05SPeter Brune     }
50070c78f05SPeter Brune   }
501b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
502b20c023fSPeter Brune 
503b20c023fSPeter Brune 
504b20c023fSPeter Brune   if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
50570c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
50670c78f05SPeter Brune     Xl    = nasm->x[i];
50770c78f05SPeter Brune     Xlloc = nasm->xl[i];
50870c78f05SPeter Brune     Yl    = nasm->y[i];
509258e1594SPeter Brune     subsnes = nasm->subsnes[i];
510258e1594SPeter Brune     ierr    = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
511111ade9eSPeter Brune     iscat   = nasm->iscatter[i];
512111ade9eSPeter Brune     oscat   = nasm->oscatter[i];
513111ade9eSPeter Brune     gscat   = nasm->gscatter[i];
514ed3c11a9SPeter Brune     ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
51524b7f281SPeter Brune     if (B) {
51624b7f281SPeter Brune       Bl   = nasm->b[i];
517ed3c11a9SPeter Brune       ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
518ed3c11a9SPeter Brune     } else Bl = NULL;
519ed3c11a9SPeter Brune     ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
52070c78f05SPeter Brune     ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
52170c78f05SPeter Brune     ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
522111ade9eSPeter Brune     ierr = VecCopy(Xl,Yl);CHKERRQ(ierr);
523ed3c11a9SPeter Brune     ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr);
524ed3c11a9SPeter Brune     ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr);
525111ade9eSPeter Brune     if (nasm->type == PC_ASM_BASIC) {
526111ade9eSPeter Brune       ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
527111ade9eSPeter Brune     } else if (nasm->type == PC_ASM_RESTRICT) {
528111ade9eSPeter Brune       ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
529ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
530eaedb033SPeter Brune   }
531ed3c11a9SPeter Brune   if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
532ed3c11a9SPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
53370c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
53470c78f05SPeter Brune     Yl    = nasm->y[i];
53570c78f05SPeter Brune     iscat   = nasm->iscatter[i];
53670c78f05SPeter Brune     oscat   = nasm->oscatter[i];
53770c78f05SPeter Brune     if (nasm->type == PC_ASM_BASIC) {
53870c78f05SPeter Brune       ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
53970c78f05SPeter Brune     } else if (nasm->type == PC_ASM_RESTRICT) {
54070c78f05SPeter Brune       ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
541ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
54270c78f05SPeter Brune   }
543b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
544111ade9eSPeter Brune   ierr = VecAXPY(X,1.0,Y);CHKERRQ(ierr);
545eaedb033SPeter Brune   PetscFunctionReturn(0);
546eaedb033SPeter Brune }
547eaedb033SPeter Brune 
548eaedb033SPeter Brune #undef __FUNCT__
549d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMComputeFinalJacobian_Private"
550*602bec5dSPeter Brune PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal)
551d728fb7dSPeter Brune {
552*602bec5dSPeter Brune   Vec            X = Xfinal;
553d728fb7dSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
554d728fb7dSPeter Brune   SNES           subsnes;
555ca44f815SPeter Brune   PetscInt       i,lag = 1;
556d728fb7dSPeter Brune   PetscErrorCode ierr;
557e59f0a30SPeter Brune   Vec            Xlloc,Xl,Fl,F;
558d728fb7dSPeter Brune   VecScatter     oscat,gscat;
559d728fb7dSPeter Brune   DM             dm,subdm;
560d728fb7dSPeter Brune   MatStructure   flg = DIFFERENT_NONZERO_PATTERN;
561d728fb7dSPeter Brune   PetscFunctionBegin;
562*602bec5dSPeter Brune   if (nasm->fjtype == 2) X = nasm->xinit;
563e59f0a30SPeter Brune   F = snes->vec_func;
564e59f0a30SPeter Brune   if (snes->normtype == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);}
565d728fb7dSPeter Brune   ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
566d728fb7dSPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
567d728fb7dSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
568*602bec5dSPeter Brune   if (nasm->fjtype != 1) {
569d728fb7dSPeter Brune     for (i=0; i<nasm->n; i++) {
570d728fb7dSPeter Brune       Xlloc = nasm->xl[i];
571d728fb7dSPeter Brune       gscat = nasm->gscatter[i];
572d728fb7dSPeter Brune       oscat = nasm->oscatter[i];
573d728fb7dSPeter Brune       ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
574*602bec5dSPeter Brune     }
575d728fb7dSPeter Brune   }
576d728fb7dSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
577d728fb7dSPeter Brune   for (i=0; i<nasm->n; i++) {
578e59f0a30SPeter Brune     Fl      = nasm->subsnes[i]->vec_func;
579d728fb7dSPeter Brune     Xl      = nasm->x[i];
580d728fb7dSPeter Brune     Xlloc   = nasm->xl[i];
581d728fb7dSPeter Brune     subsnes = nasm->subsnes[i];
582d728fb7dSPeter Brune     oscat   = nasm->oscatter[i];
583d728fb7dSPeter Brune     gscat   = nasm->gscatter[i];
584*602bec5dSPeter Brune     if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);}
585ed3c11a9SPeter Brune     ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
586d728fb7dSPeter Brune     ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
587*602bec5dSPeter Brune     if (nasm->fjtype != 1) {
588d728fb7dSPeter Brune       ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
589d728fb7dSPeter Brune       ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
590*602bec5dSPeter Brune     }
591ca44f815SPeter Brune     if (subsnes->lagjacobian == -1)    subsnes->lagjacobian = -2;
592ca44f815SPeter Brune     else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian;
593*602bec5dSPeter Brune     ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr);
5942f543b25SPeter Brune     ierr = SNESComputeJacobian(subsnes,Xl,&subsnes->jacobian,&subsnes->jacobian_pre,&flg);CHKERRQ(ierr);
595ca44f815SPeter Brune     if (lag > 1) subsnes->lagjacobian = lag;
596d728fb7dSPeter Brune     ierr = KSPSetOperators(subsnes->ksp,subsnes->jacobian,subsnes->jacobian_pre,flg);CHKERRQ(ierr);
597d728fb7dSPeter Brune   }
598d728fb7dSPeter Brune   PetscFunctionReturn(0);
599d728fb7dSPeter Brune }
600d728fb7dSPeter Brune 
601d728fb7dSPeter Brune #undef __FUNCT__
602eaedb033SPeter Brune #define __FUNCT__ "SNESSolve_NASM"
603eaedb033SPeter Brune PetscErrorCode SNESSolve_NASM(SNES snes)
604eaedb033SPeter Brune {
605eaedb033SPeter Brune   Vec            F;
606eaedb033SPeter Brune   Vec            X;
607eaedb033SPeter Brune   Vec            B;
608111ade9eSPeter Brune   Vec            Y;
609eaedb033SPeter Brune   PetscInt       i;
610ed3c11a9SPeter Brune   PetscReal      fnorm = 0.0;
611eaedb033SPeter Brune   PetscErrorCode ierr;
612eaedb033SPeter Brune   SNESNormType   normtype;
613d728fb7dSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
614eaedb033SPeter Brune 
615eaedb033SPeter Brune   PetscFunctionBegin;
616eaedb033SPeter Brune   X = snes->vec_sol;
617111ade9eSPeter Brune   Y = snes->vec_sol_update;
618eaedb033SPeter Brune   F = snes->vec_func;
619eaedb033SPeter Brune   B = snes->vec_rhs;
620eaedb033SPeter Brune 
621eaedb033SPeter Brune   ierr         = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
622eaedb033SPeter Brune   snes->iter   = 0;
623eaedb033SPeter Brune   snes->norm   = 0.;
624eaedb033SPeter Brune   ierr         = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
625eaedb033SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
626eaedb033SPeter Brune   ierr         = SNESGetNormType(snes, &normtype);CHKERRQ(ierr);
627eaedb033SPeter Brune   if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
628eaedb033SPeter Brune     /* compute the initial function and preconditioned update delX */
629eaedb033SPeter Brune     if (!snes->vec_func_init_set) {
630eaedb033SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
631eaedb033SPeter Brune       if (snes->domainerror) {
632eaedb033SPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
633eaedb033SPeter Brune         PetscFunctionReturn(0);
634eaedb033SPeter Brune       }
6351aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
636eaedb033SPeter Brune 
637eaedb033SPeter Brune     /* convergence test */
638eaedb033SPeter Brune     if (!snes->norm_init_set) {
639eaedb033SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
640189a9710SBarry Smith       if (PetscIsInfOrNanReal(fnorm)) {
641189a9710SBarry Smith         snes->reason = SNES_DIVERGED_FNORM_NAN;
642189a9710SBarry Smith         PetscFunctionReturn(0);
643189a9710SBarry Smith       }
644eaedb033SPeter Brune     } else {
645eaedb033SPeter Brune       fnorm               = snes->norm_init;
646eaedb033SPeter Brune       snes->norm_init_set = PETSC_FALSE;
647eaedb033SPeter Brune     }
648eaedb033SPeter Brune     ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
649eaedb033SPeter Brune     snes->iter = 0;
650eaedb033SPeter Brune     snes->norm = fnorm;
651eaedb033SPeter Brune     ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
652a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
653eaedb033SPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
654eaedb033SPeter Brune 
655eaedb033SPeter Brune     /* set parameter for default relative tolerance convergence test */
656eaedb033SPeter Brune     snes->ttol = fnorm*snes->rtol;
657eaedb033SPeter Brune 
658eaedb033SPeter Brune     /* test convergence */
659eaedb033SPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
660eaedb033SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
661eaedb033SPeter Brune   } else {
662eaedb033SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
663a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
664eaedb033SPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
665eaedb033SPeter Brune   }
666eaedb033SPeter Brune 
667eaedb033SPeter Brune   /* Call general purpose update function */
668eaedb033SPeter Brune   if (snes->ops->update) {
669eaedb033SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
670eaedb033SPeter Brune   }
671*602bec5dSPeter Brune   /* copy the initial solution over for later */
672*602bec5dSPeter Brune   if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);}
673eaedb033SPeter Brune 
674eaedb033SPeter Brune   for (i = 0; i < snes->max_its; i++) {
675111ade9eSPeter Brune     ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr);
676eaedb033SPeter Brune     if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) {
677eaedb033SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
678eaedb033SPeter Brune       if (snes->domainerror) {
679eaedb033SPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
680d728fb7dSPeter Brune         break;
681eaedb033SPeter Brune       }
682eaedb033SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
683189a9710SBarry Smith       if (PetscIsInfOrNanReal(fnorm)) {
684189a9710SBarry Smith         snes->reason = SNES_DIVERGED_FNORM_NAN;
685189a9710SBarry Smith         break;
686189a9710SBarry Smith       }
687eaedb033SPeter Brune     }
688eaedb033SPeter Brune     /* Monitor convergence */
689eaedb033SPeter Brune     ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
690eaedb033SPeter Brune     snes->iter = i+1;
691eaedb033SPeter Brune     snes->norm = fnorm;
692eaedb033SPeter Brune     ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
693a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
694eaedb033SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
695eaedb033SPeter Brune     /* Test for convergence */
696e59f0a30SPeter Brune     if (normtype == SNES_NORM_FUNCTION) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
697d728fb7dSPeter Brune     if (snes->reason) break;
698eaedb033SPeter Brune     /* Call general purpose update function */
699e59f0a30SPeter Brune     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
700eaedb033SPeter Brune   }
701d728fb7dSPeter Brune   if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);}
702eaedb033SPeter Brune   if (normtype == SNES_NORM_FUNCTION) {
703eaedb033SPeter Brune     if (i == snes->max_its) {
704eaedb033SPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
705eaedb033SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
706eaedb033SPeter Brune     }
7071aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */
708eaedb033SPeter Brune   PetscFunctionReturn(0);
709eaedb033SPeter Brune }
710eaedb033SPeter Brune 
711eaedb033SPeter Brune /*MC
712eaedb033SPeter Brune   SNESNASM - Nonlinear Additive Schwartz
713eaedb033SPeter Brune 
71470c78603SPeter Brune    Options Database:
71570c78603SPeter Brune +  -snes_nasm_log - enable logging events for the communication and solve stages
71670c78603SPeter Brune .  -snes_nasm_type <basic,restrict> - type of subdomain update used
71770c78603SPeter Brune .  -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate
718*602bec5dSPeter Brune .  -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> pick state the jacobian is calculated at
71970c78603SPeter Brune .  -sub_snes_ - options prefix of the subdomain nonlinear solves
72070c78603SPeter Brune .  -sub_ksp_ - options prefix of the subdomain Krylov solver
72170c78603SPeter Brune -  -sub_pc_ - options prefix of the subdomain preconditioner
72270c78603SPeter Brune 
723eaedb033SPeter Brune    Level: advanced
724eaedb033SPeter Brune 
725eaedb033SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
726eaedb033SPeter Brune M*/
727eaedb033SPeter Brune 
728eaedb033SPeter Brune #undef __FUNCT__
729eaedb033SPeter Brune #define __FUNCT__ "SNESCreate_NASM"
7308cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes)
731eaedb033SPeter Brune {
732eaedb033SPeter Brune   SNES_NASM      *nasm;
733eaedb033SPeter Brune   PetscErrorCode ierr;
734eaedb033SPeter Brune 
735eaedb033SPeter Brune   PetscFunctionBegin;
736eaedb033SPeter Brune   ierr       = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr);
737eaedb033SPeter Brune   snes->data = (void*)nasm;
738eaedb033SPeter Brune 
739eaedb033SPeter Brune   nasm->n        = PETSC_DECIDE;
740eaedb033SPeter Brune   nasm->subsnes  = 0;
741eaedb033SPeter Brune   nasm->x        = 0;
742111ade9eSPeter Brune   nasm->xl       = 0;
743111ade9eSPeter Brune   nasm->y        = 0;
744eaedb033SPeter Brune   nasm->b        = 0;
745111ade9eSPeter Brune   nasm->oscatter = 0;
746111ade9eSPeter Brune   nasm->iscatter = 0;
747111ade9eSPeter Brune   nasm->gscatter = 0;
748111ade9eSPeter Brune 
749111ade9eSPeter Brune   nasm->type = PC_ASM_BASIC;
750d728fb7dSPeter Brune   nasm->finaljacobian = PETSC_FALSE;
751eaedb033SPeter Brune 
752eaedb033SPeter Brune   snes->ops->destroy        = SNESDestroy_NASM;
753eaedb033SPeter Brune   snes->ops->setup          = SNESSetUp_NASM;
754eaedb033SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
755eaedb033SPeter Brune   snes->ops->view           = SNESView_NASM;
756eaedb033SPeter Brune   snes->ops->solve          = SNESSolve_NASM;
757eaedb033SPeter Brune   snes->ops->reset          = SNESReset_NASM;
758eaedb033SPeter Brune 
759eaedb033SPeter Brune   snes->usesksp = PETSC_FALSE;
760eaedb033SPeter Brune   snes->usespc  = PETSC_FALSE;
761eaedb033SPeter Brune 
762*602bec5dSPeter Brune   nasm->fjtype              = 0;
763*602bec5dSPeter Brune   nasm->xinit               = NULL;
7640298fd71SBarry Smith   nasm->eventrestrictinterp = 0;
7650298fd71SBarry Smith   nasm->eventsubsolve       = 0;
766b20c023fSPeter Brune 
767eaedb033SPeter Brune   if (!snes->tolerancesset) {
768eaedb033SPeter Brune     snes->max_its   = 10000;
769eaedb033SPeter Brune     snes->max_funcs = 10000;
770eaedb033SPeter Brune   }
771eaedb033SPeter Brune 
77200de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C","SNESNASMSetSubdomains_NASM",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr);
77300de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C","SNESNASMGetSubdomains_NASM",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr);
77400de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C","SNESNASMGetSubdomainVecs_NASM",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr);
77500de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C","SNESNASMSetComputeFinalJacobian_NASM",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr);
776eaedb033SPeter Brune   PetscFunctionReturn(0);
777eaedb033SPeter Brune }
77899e0435eSBarry Smith 
779