xref: /petsc/src/snes/impls/nasm/nasm.c (revision a4f1787623b00ac75b859ba8d8c035f42daafb83)
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 */
17*a4f17876SPeter Brune   PetscBool  same_local_solves;   /* flag to determine if the solvers have been individually modified */
18b20c023fSPeter Brune 
19b20c023fSPeter Brune   /* logging events */
20b20c023fSPeter Brune   PetscLogEvent eventrestrictinterp;
21b20c023fSPeter Brune   PetscLogEvent eventsubsolve;
22eaedb033SPeter Brune } SNES_NASM;
23eaedb033SPeter Brune 
24b20c023fSPeter Brune const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0};
25b20c023fSPeter Brune 
26eaedb033SPeter Brune #undef __FUNCT__
27eaedb033SPeter Brune #define __FUNCT__ "SNESReset_NASM"
28eaedb033SPeter Brune PetscErrorCode SNESReset_NASM(SNES snes)
29eaedb033SPeter Brune {
30eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
31eaedb033SPeter Brune   PetscErrorCode ierr;
32eaedb033SPeter Brune   PetscInt       i;
336e111a19SKarl Rupp 
34eaedb033SPeter Brune   PetscFunctionBegin;
35eaedb033SPeter Brune   for (i=0; i<nasm->n; i++) {
36111ade9eSPeter Brune     if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); }
37f5f7c1b9SKarl Rupp     if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); }
38111ade9eSPeter Brune     if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); }
39bc8c1f72SJose Roman     if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); }
40eaedb033SPeter Brune 
41bc8c1f72SJose Roman     if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); }
42111ade9eSPeter Brune     if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); }
43111ade9eSPeter Brune     if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); }
44111ade9eSPeter Brune     if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); }
45eaedb033SPeter Brune   }
46111ade9eSPeter Brune 
47111ade9eSPeter Brune   if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);}
48111ade9eSPeter Brune   if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);}
49111ade9eSPeter Brune   if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);}
50111ade9eSPeter Brune   if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);}
51111ade9eSPeter Brune 
52111ade9eSPeter Brune   if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);}
53111ade9eSPeter Brune   if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);}
54111ade9eSPeter Brune   if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);}
55111ade9eSPeter Brune   if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);}
56b20c023fSPeter Brune 
57b20c023fSPeter Brune   nasm->eventrestrictinterp = 0;
58b20c023fSPeter Brune   nasm->eventsubsolve = 0;
59eaedb033SPeter Brune   PetscFunctionReturn(0);
60eaedb033SPeter Brune }
61eaedb033SPeter Brune 
62eaedb033SPeter Brune #undef __FUNCT__
63eaedb033SPeter Brune #define __FUNCT__ "SNESDestroy_NASM"
64eaedb033SPeter Brune PetscErrorCode SNESDestroy_NASM(SNES snes)
65eaedb033SPeter Brune {
66eaedb033SPeter Brune   PetscErrorCode ierr;
676e111a19SKarl Rupp 
68eaedb033SPeter Brune   PetscFunctionBegin;
69eaedb033SPeter Brune   ierr = SNESReset_NASM(snes);CHKERRQ(ierr);
7022d28d08SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
71eaedb033SPeter Brune   PetscFunctionReturn(0);
72eaedb033SPeter Brune }
73eaedb033SPeter Brune 
74eaedb033SPeter Brune #undef __FUNCT__
75111ade9eSPeter Brune #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private"
760adebc6cSBarry Smith PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx)
770adebc6cSBarry Smith {
78111ade9eSPeter Brune   PetscErrorCode ierr;
79111ade9eSPeter Brune   Vec            bcs = (Vec)ctx;
806e111a19SKarl Rupp 
81111ade9eSPeter Brune   PetscFunctionBegin;
82111ade9eSPeter Brune   ierr = VecCopy(bcs,l);CHKERRQ(ierr);
83111ade9eSPeter Brune   PetscFunctionReturn(0);
84111ade9eSPeter Brune }
85111ade9eSPeter Brune 
86111ade9eSPeter Brune #undef __FUNCT__
87eaedb033SPeter Brune #define __FUNCT__ "SNESSetUp_NASM"
88eaedb033SPeter Brune PetscErrorCode SNESSetUp_NASM(SNES snes)
89eaedb033SPeter Brune {
90eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
91eaedb033SPeter Brune   PetscErrorCode ierr;
9276857b2aSPeter Brune   DM             dm,subdm;
93111ade9eSPeter Brune   DM             *subdms;
94111ade9eSPeter Brune   PetscInt       i;
95eaedb033SPeter Brune   const char     *optionsprefix;
96111ade9eSPeter Brune   Vec            F;
97ed3c11a9SPeter Brune   PetscMPIInt    size;
98ed3c11a9SPeter Brune   KSP            ksp;
99ed3c11a9SPeter Brune   PC             pc;
100eaedb033SPeter Brune 
101eaedb033SPeter Brune   PetscFunctionBegin;
102eaedb033SPeter Brune   if (!nasm->subsnes) {
103eaedb033SPeter Brune     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1040a696f66SPeter Brune     if (dm) {
105eaedb033SPeter Brune       nasm->usesdm = PETSC_TRUE;
1060298fd71SBarry Smith       ierr         = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr);
107ce94432eSBarry Smith       if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
108111ade9eSPeter Brune       ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr);
109eaedb033SPeter Brune 
110eaedb033SPeter Brune       ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
111111ade9eSPeter Brune       ierr = PetscMalloc(nasm->n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr);
112111ade9eSPeter Brune       for (i=0; i<nasm->n; i++) {
113cdb298fcSPeter Brune         ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr);
114cdb298fcSPeter Brune         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr);
115cdb298fcSPeter Brune         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr);
116cdb298fcSPeter Brune         ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr);
117ed3c11a9SPeter Brune         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRQ(ierr);
118ed3c11a9SPeter Brune         if (size == 1) {
119ed3c11a9SPeter Brune           ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr);
120ed3c11a9SPeter Brune           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
121ed3c11a9SPeter Brune           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
122ed3c11a9SPeter Brune           ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
123ed3c11a9SPeter Brune         }
124cdb298fcSPeter Brune         ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr);
125111ade9eSPeter Brune         ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr);
126111ade9eSPeter Brune       }
127111ade9eSPeter Brune       ierr = PetscFree(subdms);CHKERRQ(ierr);
128ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!");
129ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!");
130111ade9eSPeter Brune   /* allocate the global vectors */
131e245e0aaSPeter Brune   if (!nasm->x) {
132111ade9eSPeter Brune     ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->x);CHKERRQ(ierr);
133e245e0aaSPeter Brune     ierr = PetscMemzero(nasm->x,nasm->n*sizeof(Vec));CHKERRQ(ierr);
134e245e0aaSPeter Brune   }
135e245e0aaSPeter Brune   if (!nasm->xl) {
136111ade9eSPeter Brune     ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->xl);CHKERRQ(ierr);
137e245e0aaSPeter Brune     ierr = PetscMemzero(nasm->xl,nasm->n*sizeof(Vec));CHKERRQ(ierr);
138e245e0aaSPeter Brune   }
139e245e0aaSPeter Brune   if (!nasm->y) {
140111ade9eSPeter Brune     ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->y);CHKERRQ(ierr);
141e245e0aaSPeter Brune     ierr = PetscMemzero(nasm->y,nasm->n*sizeof(Vec));CHKERRQ(ierr);
142e245e0aaSPeter Brune   }
143e245e0aaSPeter Brune   if (!nasm->b) {
144111ade9eSPeter Brune     ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->b);CHKERRQ(ierr);
145e245e0aaSPeter Brune     ierr = PetscMemzero(nasm->b,nasm->n*sizeof(Vec));CHKERRQ(ierr);
146e245e0aaSPeter Brune   }
147111ade9eSPeter Brune 
148111ade9eSPeter Brune   for (i=0; i<nasm->n; i++) {
1490298fd71SBarry Smith     ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr);
15076857b2aSPeter Brune     if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);}
15176857b2aSPeter Brune     if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);}
15276857b2aSPeter Brune     if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);}
15376857b2aSPeter Brune     if (!nasm->xl[i]) {
154111ade9eSPeter Brune       ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr);
155111ade9eSPeter Brune       ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr);
15676857b2aSPeter Brune     }
1570298fd71SBarry Smith     ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr);
158111ade9eSPeter Brune   }
159d728fb7dSPeter Brune   if (nasm->finaljacobian) {ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);}
160eaedb033SPeter Brune   PetscFunctionReturn(0);
161eaedb033SPeter Brune }
162eaedb033SPeter Brune 
163eaedb033SPeter Brune #undef __FUNCT__
164eaedb033SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NASM"
165eaedb033SPeter Brune PetscErrorCode SNESSetFromOptions_NASM(SNES snes)
166eaedb033SPeter Brune {
167eaedb033SPeter Brune   PetscErrorCode    ierr;
168111ade9eSPeter Brune   PCASMType         asmtype;
169*a4f17876SPeter Brune   PetscBool         flg,monflg,subviewflg;
170111ade9eSPeter Brune   SNES_NASM         *nasm = (SNES_NASM*)snes->data;
1716e111a19SKarl Rupp 
172eaedb033SPeter Brune   PetscFunctionBegin;
173111ade9eSPeter Brune   ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr);
174111ade9eSPeter Brune   ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
1751aa26658SKarl Rupp   if (flg) nasm->type = asmtype;
176b20c023fSPeter Brune   flg    = PETSC_FALSE;
177b20c023fSPeter Brune   monflg = PETSC_TRUE;
178*a4f17876SPeter Brune   subviewflg = PETSC_FALSE;
179*a4f17876SPeter Brune   ierr   = PetscOptionsBool("-snes_nasm_sub_view","Print detailed information for every processor when using -snes_view","",subviewflg,&subviewflg,&flg);CHKERRQ(ierr);
180*a4f17876SPeter Brune   if (flg) {
181*a4f17876SPeter Brune     nasm->same_local_solves = PETSC_FALSE;
182*a4f17876SPeter Brune     if (!subviewflg) {
183*a4f17876SPeter Brune       nasm->same_local_solves = PETSC_TRUE;
184*a4f17876SPeter Brune     }
185*a4f17876SPeter Brune   }
186d728fb7dSPeter Brune   ierr   = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr);
187*a4f17876SPeter Brune   ierr   = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr);
188b20c023fSPeter Brune   if (flg) {
189b20c023fSPeter Brune     ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr);
190b20c023fSPeter Brune     ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr);
191b20c023fSPeter Brune   }
192eaedb033SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
193eaedb033SPeter Brune   PetscFunctionReturn(0);
194eaedb033SPeter Brune }
195eaedb033SPeter Brune 
196eaedb033SPeter Brune #undef __FUNCT__
197eaedb033SPeter Brune #define __FUNCT__ "SNESView_NASM"
198eaedb033SPeter Brune PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
199eaedb033SPeter Brune {
200b20c023fSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
201b20c023fSPeter Brune   PetscErrorCode ierr;
202*a4f17876SPeter Brune   PetscMPIInt    rank,size;
203*a4f17876SPeter Brune   PetscInt       i,j,N,bsz,nmax,nmin;
204b20c023fSPeter Brune   PetscBool      iascii,isstring;
205b20c023fSPeter Brune   PetscViewer    sviewer;
206ce94432eSBarry Smith   MPI_Comm       comm;
207b20c023fSPeter Brune 
208b20c023fSPeter Brune   PetscFunctionBegin;
209ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
210b20c023fSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
211b20c023fSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
212b20c023fSPeter Brune   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
213*a4f17876SPeter Brune   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
214b20c023fSPeter Brune   ierr = MPI_Reduce(&nasm->n,&N,1,MPIU_INT,MPIU_SUM,0,comm);CHKERRQ(ierr);
215b20c023fSPeter Brune   if (iascii) {
216b20c023fSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  Nonlinear Additive Schwarz: total subdomain blocks = %D\n",N);CHKERRQ(ierr);
217*a4f17876SPeter Brune     ierr = MPI_Reduce(&nasm->n,&nmax,1,MPIU_INT,MPIU_MAX,0,comm);CHKERRQ(ierr);
218*a4f17876SPeter Brune     ierr = MPI_Reduce(&nasm->n,&nmin,1,MPIU_INT,MPIU_MIN,0,comm);CHKERRQ(ierr);
219*a4f17876SPeter Brune     if (nasm->same_local_solves) {
220*a4f17876SPeter Brune       if (nasm->subsnes) {
221*a4f17876SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  Local solve is the same for all blocks:\n");CHKERRQ(ierr);
222*a4f17876SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
223*a4f17876SPeter Brune         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
224*a4f17876SPeter Brune         if (!rank) {
225*a4f17876SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
226*a4f17876SPeter Brune           ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);
227*a4f17876SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
228*a4f17876SPeter Brune         }
229*a4f17876SPeter Brune         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
230*a4f17876SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
231*a4f17876SPeter Brune       }
232*a4f17876SPeter Brune     } else {
233*a4f17876SPeter Brune       /* print the solver on each block */
234b20c023fSPeter Brune       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
235b20c023fSPeter Brune       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr);
236b20c023fSPeter Brune       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
237*a4f17876SPeter Brune       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
238*a4f17876SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following SNES objects:\n");CHKERRQ(ierr);
239b20c023fSPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
240*a4f17876SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
241*a4f17876SPeter Brune       for (j=0; j<size; j++) {
242b20c023fSPeter Brune         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
243*a4f17876SPeter Brune         if (rank == j) {
244*a4f17876SPeter Brune           for (i=0; i<nasm->n; i++) {
245*a4f17876SPeter Brune             ierr = VecGetLocalSize(nasm->x[i],&bsz);CHKERRQ(ierr);
246*a4f17876SPeter Brune             ierr = PetscViewerASCIIPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr);
247b20c023fSPeter Brune             ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr);
248*a4f17876SPeter Brune             ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
249b20c023fSPeter Brune           }
250b20c023fSPeter Brune         }
251*a4f17876SPeter Brune         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
252*a4f17876SPeter Brune         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
253*a4f17876SPeter Brune       }
254b20c023fSPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
255*a4f17876SPeter Brune     }
256b20c023fSPeter Brune   } else if (isstring) {
257b20c023fSPeter Brune     ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr);
258b20c023fSPeter Brune     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
259b20c023fSPeter Brune     if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);}
260b20c023fSPeter Brune     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
261b20c023fSPeter Brune   }
262eaedb033SPeter Brune   PetscFunctionReturn(0);
263eaedb033SPeter Brune }
264eaedb033SPeter Brune 
265eaedb033SPeter Brune #undef __FUNCT__
266eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains"
26776857b2aSPeter Brune /*@
26876857b2aSPeter Brune    SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems.
26976857b2aSPeter Brune 
27076857b2aSPeter Brune    Not Collective
27176857b2aSPeter Brune 
27276857b2aSPeter Brune    Input Parameters:
27376857b2aSPeter Brune +  SNES - the SNES context
27476857b2aSPeter Brune .  n - the number of local subdomains
27576857b2aSPeter Brune .  subsnes - solvers defined on the local subdomains
27676857b2aSPeter Brune .  iscatter - scatters into the nonoverlapping portions of the local subdomains
27776857b2aSPeter Brune .  oscatter - scatters into the overlapping portions of the local subdomains
27876857b2aSPeter Brune -  gscatter - scatters into the (ghosted) local vector of the local subdomain
27976857b2aSPeter Brune 
28076857b2aSPeter Brune    Level: intermediate
28176857b2aSPeter Brune 
28276857b2aSPeter Brune .keywords: SNES, NASM
28376857b2aSPeter Brune 
28476857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains()
28576857b2aSPeter Brune @*/
286a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
287a6dfd86eSKarl Rupp {
288eaedb033SPeter Brune   PetscErrorCode ierr;
289111ade9eSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*);
2906e111a19SKarl Rupp 
291eaedb033SPeter Brune   PetscFunctionBegin;
2920005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f);CHKERRQ(ierr);
293111ade9eSPeter Brune   ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);
294eaedb033SPeter Brune   PetscFunctionReturn(0);
295eaedb033SPeter Brune }
296eaedb033SPeter Brune 
297eaedb033SPeter Brune #undef __FUNCT__
298eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains_NASM"
299a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
300a6dfd86eSKarl Rupp {
301eaedb033SPeter Brune   PetscInt       i;
302eaedb033SPeter Brune   PetscErrorCode ierr;
303eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
3046e111a19SKarl Rupp 
305eaedb033SPeter Brune   PetscFunctionBegin;
306ce94432eSBarry Smith   if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
307eaedb033SPeter Brune 
308111ade9eSPeter Brune   /* tear down the previously set things */
309111ade9eSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
310111ade9eSPeter Brune 
311eaedb033SPeter Brune   nasm->n = n;
312111ade9eSPeter Brune   if (oscatter) {
313111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);}
314eaedb033SPeter Brune   }
315111ade9eSPeter Brune   if (iscatter) {
316111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);}
317eaedb033SPeter Brune   }
318111ade9eSPeter Brune   if (gscatter) {
319111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);}
320111ade9eSPeter Brune   }
321111ade9eSPeter Brune   if (oscatter) {
322111ade9eSPeter Brune     ierr = PetscMalloc(n*sizeof(IS),&nasm->oscatter);CHKERRQ(ierr);
323eaedb033SPeter Brune     for (i=0; i<n; i++) {
324111ade9eSPeter Brune       nasm->oscatter[i] = oscatter[i];
325eaedb033SPeter Brune     }
326111ade9eSPeter Brune   }
327111ade9eSPeter Brune   if (iscatter) {
328111ade9eSPeter Brune     ierr = PetscMalloc(n*sizeof(IS),&nasm->iscatter);CHKERRQ(ierr);
329eaedb033SPeter Brune     for (i=0; i<n; i++) {
330111ade9eSPeter Brune       nasm->iscatter[i] = iscatter[i];
331eaedb033SPeter Brune     }
332eaedb033SPeter Brune   }
333111ade9eSPeter Brune   if (gscatter) {
334111ade9eSPeter Brune     ierr = PetscMalloc(n*sizeof(IS),&nasm->gscatter);CHKERRQ(ierr);
335eaedb033SPeter Brune     for (i=0; i<n; i++) {
336111ade9eSPeter Brune       nasm->gscatter[i] = gscatter[i];
337eaedb033SPeter Brune     }
338eaedb033SPeter Brune   }
339111ade9eSPeter Brune 
340eaedb033SPeter Brune   if (subsnes) {
341eaedb033SPeter Brune     ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr);
342eaedb033SPeter Brune     for (i=0; i<n; i++) {
343eaedb033SPeter Brune       nasm->subsnes[i] = subsnes[i];
344eaedb033SPeter Brune     }
345*a4f17876SPeter Brune     nasm->same_local_solves = PETSC_FALSE;
346eaedb033SPeter Brune   }
347eaedb033SPeter Brune   PetscFunctionReturn(0);
348eaedb033SPeter Brune }
349eaedb033SPeter Brune 
350eaedb033SPeter Brune #undef __FUNCT__
35176857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains"
35276857b2aSPeter Brune /*@
35376857b2aSPeter Brune    SNESNASMGetSubdomains - Get the local subdomain context.
35476857b2aSPeter Brune 
35576857b2aSPeter Brune    Not Collective
35676857b2aSPeter Brune 
35776857b2aSPeter Brune    Input Parameters:
35876857b2aSPeter Brune .  SNES - the SNES context
35976857b2aSPeter Brune 
36076857b2aSPeter Brune    Output Parameters:
36176857b2aSPeter Brune +  n - the number of local subdomains
36276857b2aSPeter Brune .  subsnes - solvers defined on the local subdomains
36376857b2aSPeter Brune .  iscatter - scatters into the nonoverlapping portions of the local subdomains
36476857b2aSPeter Brune .  oscatter - scatters into the overlapping portions of the local subdomains
36576857b2aSPeter Brune -  gscatter - scatters into the (ghosted) local vector of the local subdomain
36676857b2aSPeter Brune 
36776857b2aSPeter Brune    Level: intermediate
36876857b2aSPeter Brune 
36976857b2aSPeter Brune .keywords: SNES, NASM
37076857b2aSPeter Brune 
37176857b2aSPeter Brune .seealso: SNESNASM, SNESNASMSetSubdomains()
37276857b2aSPeter Brune @*/
37376857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
37476857b2aSPeter Brune {
37576857b2aSPeter Brune   PetscErrorCode ierr;
37676857b2aSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**);
37776857b2aSPeter Brune 
37876857b2aSPeter Brune   PetscFunctionBegin;
3790005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f);CHKERRQ(ierr);
38076857b2aSPeter Brune   ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);
38176857b2aSPeter Brune   PetscFunctionReturn(0);
38276857b2aSPeter Brune }
38376857b2aSPeter Brune 
38476857b2aSPeter Brune #undef __FUNCT__
38576857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains_NASM"
38676857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
38776857b2aSPeter Brune {
38876857b2aSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
38976857b2aSPeter Brune 
39076857b2aSPeter Brune   PetscFunctionBegin;
39176857b2aSPeter Brune   if (n) *n = nasm->n;
39276857b2aSPeter Brune   if (oscatter) *oscatter = nasm->oscatter;
39376857b2aSPeter Brune   if (iscatter) *iscatter = nasm->iscatter;
39476857b2aSPeter Brune   if (gscatter) *gscatter = nasm->gscatter;
395*a4f17876SPeter Brune   if (subsnes)  {
396*a4f17876SPeter Brune     *subsnes  = nasm->subsnes;
397*a4f17876SPeter Brune     nasm->same_local_solves = PETSC_FALSE;
398*a4f17876SPeter Brune   }
39976857b2aSPeter Brune   PetscFunctionReturn(0);
40076857b2aSPeter Brune }
40176857b2aSPeter Brune 
40276857b2aSPeter Brune #undef __FUNCT__
40376857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs"
40476857b2aSPeter Brune /*@
40576857b2aSPeter Brune    SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors
40676857b2aSPeter Brune 
40776857b2aSPeter Brune    Not Collective
40876857b2aSPeter Brune 
40976857b2aSPeter Brune    Input Parameters:
41076857b2aSPeter Brune .  SNES - the SNES context
41176857b2aSPeter Brune 
41276857b2aSPeter Brune    Output Parameters:
41376857b2aSPeter Brune +  n - the number of local subdomains
41476857b2aSPeter Brune .  x - The subdomain solution vector
41576857b2aSPeter Brune .  y - The subdomain step vector
41676857b2aSPeter Brune .  b - The subdomain RHS vector
41776857b2aSPeter Brune -  xl - The subdomain local vectors (ghosted)
41876857b2aSPeter Brune 
41976857b2aSPeter Brune    Level: developer
42076857b2aSPeter Brune 
42176857b2aSPeter Brune .keywords: SNES, NASM
42276857b2aSPeter Brune 
42376857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains()
42476857b2aSPeter Brune @*/
42576857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl)
42676857b2aSPeter Brune {
42776857b2aSPeter Brune   PetscErrorCode ierr;
42876857b2aSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**);
42976857b2aSPeter Brune 
43076857b2aSPeter Brune   PetscFunctionBegin;
4310005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f);CHKERRQ(ierr);
43276857b2aSPeter Brune   ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);
43376857b2aSPeter Brune   PetscFunctionReturn(0);
43476857b2aSPeter Brune }
43576857b2aSPeter Brune 
43676857b2aSPeter Brune #undef __FUNCT__
43776857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs_NASM"
43876857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl)
43976857b2aSPeter Brune {
44076857b2aSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
44176857b2aSPeter Brune 
44276857b2aSPeter Brune   PetscFunctionBegin;
44376857b2aSPeter Brune   if (n)  *n  = nasm->n;
44476857b2aSPeter Brune   if (x)  *x  = nasm->x;
44576857b2aSPeter Brune   if (y)  *y  = nasm->y;
44676857b2aSPeter Brune   if (b)  *b  = nasm->b;
44776857b2aSPeter Brune   if (xl) *xl = nasm->xl;
44876857b2aSPeter Brune   PetscFunctionReturn(0);
44976857b2aSPeter Brune }
45076857b2aSPeter Brune 
451d728fb7dSPeter Brune #undef __FUNCT__
452d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian"
453d728fb7dSPeter Brune /*@
454d728fb7dSPeter Brune    SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence
455d728fb7dSPeter Brune 
456d728fb7dSPeter Brune    Collective on SNES
457d728fb7dSPeter Brune 
458d728fb7dSPeter Brune    Input Parameters:
459d728fb7dSPeter Brune +  SNES - the SNES context
460d728fb7dSPeter Brune -  flg - indication of whether to compute the jacobians or not
461d728fb7dSPeter Brune 
462d728fb7dSPeter Brune    Level: developer
463d728fb7dSPeter Brune 
464d728fb7dSPeter Brune    Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian
465d728fb7dSPeter Brune    is needed at each linear iteration.
466d728fb7dSPeter Brune 
467d728fb7dSPeter Brune .keywords: SNES, NASM, ASPIN
468d728fb7dSPeter Brune 
469d728fb7dSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains()
470d728fb7dSPeter Brune @*/
471d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg)
472d728fb7dSPeter Brune {
473d728fb7dSPeter Brune   PetscErrorCode (*f)(SNES,PetscBool);
474d728fb7dSPeter Brune   PetscErrorCode ierr;
475d728fb7dSPeter Brune 
476d728fb7dSPeter Brune   PetscFunctionBegin;
4770005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f);CHKERRQ(ierr);
478d728fb7dSPeter Brune   ierr = (f)(snes,flg);CHKERRQ(ierr);
479d728fb7dSPeter Brune   PetscFunctionReturn(0);
480d728fb7dSPeter Brune }
481d728fb7dSPeter Brune 
482d728fb7dSPeter Brune #undef __FUNCT__
483d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian_NASM"
484d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg)
485d728fb7dSPeter Brune {
486d728fb7dSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
487d728fb7dSPeter Brune 
488d728fb7dSPeter Brune   PetscFunctionBegin;
489d728fb7dSPeter Brune   nasm->finaljacobian = flg;
490d728fb7dSPeter Brune   PetscFunctionReturn(0);
491d728fb7dSPeter Brune }
49276857b2aSPeter Brune 
49376857b2aSPeter Brune #undef __FUNCT__
494eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSolveLocal_Private"
4950adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X)
4960adebc6cSBarry Smith {
497eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
498258e1594SPeter Brune   SNES           subsnes;
499eaedb033SPeter Brune   PetscInt       i;
500eaedb033SPeter Brune   PetscErrorCode ierr;
501111ade9eSPeter Brune   Vec            Xlloc,Xl,Bl,Yl;
502111ade9eSPeter Brune   VecScatter     iscat,oscat,gscat;
503111ade9eSPeter Brune   DM             dm,subdm;
5040adebc6cSBarry Smith 
505eaedb033SPeter Brune   PetscFunctionBegin;
506eaedb033SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
507111ade9eSPeter Brune   ierr = VecSet(Y,0);CHKERRQ(ierr);
508b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
509eaedb033SPeter Brune   for (i=0; i<nasm->n; i++) {
51070c78f05SPeter Brune     /* scatter the solution to the local solution */
51170c78f05SPeter Brune     Xlloc = nasm->xl[i];
51270c78f05SPeter Brune     gscat   = nasm->gscatter[i];
51370c78f05SPeter Brune     oscat   = nasm->oscatter[i];
51470c78f05SPeter Brune     ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
51570c78f05SPeter Brune     if (B) {
51670c78f05SPeter Brune       /* scatter the RHS to the local RHS */
51770c78f05SPeter Brune       Bl   = nasm->b[i];
51870c78f05SPeter Brune       ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
51970c78f05SPeter Brune     }
52070c78f05SPeter Brune   }
521b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
522b20c023fSPeter Brune 
523b20c023fSPeter Brune 
524b20c023fSPeter Brune   if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
52570c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
52670c78f05SPeter Brune     Xl    = nasm->x[i];
52770c78f05SPeter Brune     Xlloc = nasm->xl[i];
52870c78f05SPeter Brune     Yl    = nasm->y[i];
529258e1594SPeter Brune     subsnes = nasm->subsnes[i];
530258e1594SPeter Brune     ierr    = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
531111ade9eSPeter Brune     iscat   = nasm->iscatter[i];
532111ade9eSPeter Brune     oscat   = nasm->oscatter[i];
533111ade9eSPeter Brune     gscat   = nasm->gscatter[i];
534ed3c11a9SPeter Brune     ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
53524b7f281SPeter Brune     if (B) {
53624b7f281SPeter Brune       Bl   = nasm->b[i];
537ed3c11a9SPeter Brune       ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
538ed3c11a9SPeter Brune     } else Bl = NULL;
539ed3c11a9SPeter Brune     ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
54070c78f05SPeter Brune     ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
54170c78f05SPeter Brune     ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
542111ade9eSPeter Brune     ierr = VecCopy(Xl,Yl);CHKERRQ(ierr);
543ed3c11a9SPeter Brune     ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr);
544ed3c11a9SPeter Brune     ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr);
545111ade9eSPeter Brune     if (nasm->type == PC_ASM_BASIC) {
546111ade9eSPeter Brune       ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
547111ade9eSPeter Brune     } else if (nasm->type == PC_ASM_RESTRICT) {
548111ade9eSPeter Brune       ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
549ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
550eaedb033SPeter Brune   }
551ed3c11a9SPeter Brune   if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
552ed3c11a9SPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
55370c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
55470c78f05SPeter Brune     Yl    = nasm->y[i];
55570c78f05SPeter Brune     iscat   = nasm->iscatter[i];
55670c78f05SPeter Brune     oscat   = nasm->oscatter[i];
55770c78f05SPeter Brune     if (nasm->type == PC_ASM_BASIC) {
55870c78f05SPeter Brune       ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
55970c78f05SPeter Brune     } else if (nasm->type == PC_ASM_RESTRICT) {
56070c78f05SPeter Brune       ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
561ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
56270c78f05SPeter Brune   }
563b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
564111ade9eSPeter Brune   ierr = VecAXPY(X,1.0,Y);CHKERRQ(ierr);
565eaedb033SPeter Brune   PetscFunctionReturn(0);
566eaedb033SPeter Brune }
567eaedb033SPeter Brune 
568eaedb033SPeter Brune #undef __FUNCT__
569d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMComputeFinalJacobian_Private"
570d728fb7dSPeter Brune PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec X)
571d728fb7dSPeter Brune {
572d728fb7dSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
573d728fb7dSPeter Brune   SNES           subsnes;
574ca44f815SPeter Brune   PetscInt       i,lag = 1;
575d728fb7dSPeter Brune   PetscErrorCode ierr;
576e59f0a30SPeter Brune   Vec            Xlloc,Xl,Fl,F;
577d728fb7dSPeter Brune   VecScatter     oscat,gscat;
578d728fb7dSPeter Brune   DM             dm,subdm;
579d728fb7dSPeter Brune   MatStructure   flg = DIFFERENT_NONZERO_PATTERN;
580d728fb7dSPeter Brune   PetscFunctionBegin;
581e59f0a30SPeter Brune   F = snes->vec_func;
582e59f0a30SPeter Brune   if (snes->normtype == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);}
583d728fb7dSPeter Brune   ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
584d728fb7dSPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
585d728fb7dSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
586d728fb7dSPeter Brune   for (i=0; i<nasm->n; i++) {
587d728fb7dSPeter Brune     Xlloc = nasm->xl[i];
588e59f0a30SPeter Brune     Fl    = nasm->subsnes[i]->vec_func;
589d728fb7dSPeter Brune     gscat = nasm->gscatter[i];
590d728fb7dSPeter Brune     oscat = nasm->oscatter[i];
591d728fb7dSPeter Brune     ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
592e59f0a30SPeter Brune     ierr = VecScatterBegin(oscat,F,Fl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
593d728fb7dSPeter Brune   }
594d728fb7dSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
595d728fb7dSPeter Brune   for (i=0; i<nasm->n; i++) {
596e59f0a30SPeter Brune     Fl      = nasm->subsnes[i]->vec_func;
597d728fb7dSPeter Brune     Xl      = nasm->x[i];
598d728fb7dSPeter Brune     Xlloc   = nasm->xl[i];
599d728fb7dSPeter Brune     subsnes = nasm->subsnes[i];
600d728fb7dSPeter Brune     oscat   = nasm->oscatter[i];
601d728fb7dSPeter Brune     gscat   = nasm->gscatter[i];
602ed3c11a9SPeter Brune     ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
603e59f0a30SPeter Brune     ierr = VecScatterEnd(oscat,F,Fl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
604ed3c11a9SPeter Brune     ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
605d728fb7dSPeter Brune     ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
606d728fb7dSPeter Brune     ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
607d728fb7dSPeter Brune     ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
608ca44f815SPeter Brune 
609ca44f815SPeter Brune     if (subsnes->lagjacobian == -1)    subsnes->lagjacobian = -2;
610ca44f815SPeter Brune     else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian;
6112f543b25SPeter Brune     ierr = SNESComputeJacobian(subsnes,Xl,&subsnes->jacobian,&subsnes->jacobian_pre,&flg);CHKERRQ(ierr);
612ca44f815SPeter Brune     if (lag > 1) subsnes->lagjacobian = lag;
613d728fb7dSPeter Brune     ierr = KSPSetOperators(subsnes->ksp,subsnes->jacobian,subsnes->jacobian_pre,flg);CHKERRQ(ierr);
614d728fb7dSPeter Brune   }
615d728fb7dSPeter Brune   PetscFunctionReturn(0);
616d728fb7dSPeter Brune }
617d728fb7dSPeter Brune 
618d728fb7dSPeter Brune #undef __FUNCT__
619eaedb033SPeter Brune #define __FUNCT__ "SNESSolve_NASM"
620eaedb033SPeter Brune PetscErrorCode SNESSolve_NASM(SNES snes)
621eaedb033SPeter Brune {
622eaedb033SPeter Brune   Vec            F;
623eaedb033SPeter Brune   Vec            X;
624eaedb033SPeter Brune   Vec            B;
625111ade9eSPeter Brune   Vec            Y;
626eaedb033SPeter Brune   PetscInt       i;
627ed3c11a9SPeter Brune   PetscReal      fnorm = 0.0;
628eaedb033SPeter Brune   PetscErrorCode ierr;
629eaedb033SPeter Brune   SNESNormType   normtype;
630d728fb7dSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
631eaedb033SPeter Brune 
632eaedb033SPeter Brune   PetscFunctionBegin;
633eaedb033SPeter Brune   X = snes->vec_sol;
634111ade9eSPeter Brune   Y = snes->vec_sol_update;
635eaedb033SPeter Brune   F = snes->vec_func;
636eaedb033SPeter Brune   B = snes->vec_rhs;
637eaedb033SPeter Brune 
638ce8c27fbSBarry Smith   ierr         = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
639eaedb033SPeter Brune   snes->iter   = 0;
640eaedb033SPeter Brune   snes->norm   = 0.;
641ce8c27fbSBarry Smith   ierr         = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
642eaedb033SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
643eaedb033SPeter Brune   ierr         = SNESGetNormType(snes, &normtype);CHKERRQ(ierr);
644eaedb033SPeter Brune   if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
645eaedb033SPeter Brune     /* compute the initial function and preconditioned update delX */
646eaedb033SPeter Brune     if (!snes->vec_func_init_set) {
647eaedb033SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
648eaedb033SPeter Brune       if (snes->domainerror) {
649eaedb033SPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
650eaedb033SPeter Brune         PetscFunctionReturn(0);
651eaedb033SPeter Brune       }
6521aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
653eaedb033SPeter Brune 
654eaedb033SPeter Brune     /* convergence test */
655eaedb033SPeter Brune     if (!snes->norm_init_set) {
656eaedb033SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
657189a9710SBarry Smith       if (PetscIsInfOrNanReal(fnorm)) {
658189a9710SBarry Smith         snes->reason = SNES_DIVERGED_FNORM_NAN;
659189a9710SBarry Smith         PetscFunctionReturn(0);
660189a9710SBarry Smith       }
661eaedb033SPeter Brune     } else {
662eaedb033SPeter Brune       fnorm               = snes->norm_init;
663eaedb033SPeter Brune       snes->norm_init_set = PETSC_FALSE;
664eaedb033SPeter Brune     }
665ce8c27fbSBarry Smith     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
666eaedb033SPeter Brune     snes->iter = 0;
667eaedb033SPeter Brune     snes->norm = fnorm;
668ce8c27fbSBarry Smith     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
669a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
670eaedb033SPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
671eaedb033SPeter Brune 
672eaedb033SPeter Brune     /* set parameter for default relative tolerance convergence test */
673eaedb033SPeter Brune     snes->ttol = fnorm*snes->rtol;
674eaedb033SPeter Brune 
675eaedb033SPeter Brune     /* test convergence */
676eaedb033SPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
677eaedb033SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
678eaedb033SPeter Brune   } else {
679ce8c27fbSBarry Smith     ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
680a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
681eaedb033SPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
682eaedb033SPeter Brune   }
683eaedb033SPeter Brune 
684eaedb033SPeter Brune   /* Call general purpose update function */
685eaedb033SPeter Brune   if (snes->ops->update) {
686eaedb033SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
687eaedb033SPeter Brune   }
688eaedb033SPeter Brune 
689eaedb033SPeter Brune   for (i = 0; i < snes->max_its; i++) {
690111ade9eSPeter Brune     ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr);
691eaedb033SPeter Brune     if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) {
692eaedb033SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
693eaedb033SPeter Brune       if (snes->domainerror) {
694eaedb033SPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
695d728fb7dSPeter Brune         break;
696eaedb033SPeter Brune       }
697eaedb033SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
698189a9710SBarry Smith       if (PetscIsInfOrNanReal(fnorm)) {
699189a9710SBarry Smith         snes->reason = SNES_DIVERGED_FNORM_NAN;
700189a9710SBarry Smith         break;
701189a9710SBarry Smith       }
702eaedb033SPeter Brune     }
703eaedb033SPeter Brune     /* Monitor convergence */
704ce8c27fbSBarry Smith     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
705eaedb033SPeter Brune     snes->iter = i+1;
706eaedb033SPeter Brune     snes->norm = fnorm;
707ce8c27fbSBarry Smith     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
708a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
709eaedb033SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
710eaedb033SPeter Brune     /* Test for convergence */
711e59f0a30SPeter Brune     if (normtype == SNES_NORM_FUNCTION) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
712d728fb7dSPeter Brune     if (snes->reason) break;
713eaedb033SPeter Brune     /* Call general purpose update function */
714e59f0a30SPeter Brune     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
715eaedb033SPeter Brune   }
716d728fb7dSPeter Brune   if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);}
717eaedb033SPeter Brune   if (normtype == SNES_NORM_FUNCTION) {
718eaedb033SPeter Brune     if (i == snes->max_its) {
719eaedb033SPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
720eaedb033SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
721eaedb033SPeter Brune     }
7221aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */
723eaedb033SPeter Brune   PetscFunctionReturn(0);
724eaedb033SPeter Brune }
725eaedb033SPeter Brune 
726eaedb033SPeter Brune /*MC
727eaedb033SPeter Brune   SNESNASM - Nonlinear Additive Schwartz
728eaedb033SPeter Brune 
72970c78603SPeter Brune    Options Database:
73070c78603SPeter Brune +  -snes_nasm_log - enable logging events for the communication and solve stages
73170c78603SPeter Brune .  -snes_nasm_type <basic,restrict> - type of subdomain update used
73270c78603SPeter Brune .  -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate
73370c78603SPeter Brune .  -sub_snes_ - options prefix of the subdomain nonlinear solves
73470c78603SPeter Brune .  -sub_ksp_ - options prefix of the subdomain Krylov solver
73570c78603SPeter Brune -  -sub_pc_ - options prefix of the subdomain preconditioner
73670c78603SPeter Brune 
737eaedb033SPeter Brune    Level: advanced
738eaedb033SPeter Brune 
739eaedb033SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
740eaedb033SPeter Brune M*/
741eaedb033SPeter Brune 
742eaedb033SPeter Brune #undef __FUNCT__
743eaedb033SPeter Brune #define __FUNCT__ "SNESCreate_NASM"
7448cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes)
745eaedb033SPeter Brune {
746eaedb033SPeter Brune   SNES_NASM      *nasm;
747eaedb033SPeter Brune   PetscErrorCode ierr;
748eaedb033SPeter Brune 
749eaedb033SPeter Brune   PetscFunctionBegin;
750eaedb033SPeter Brune   ierr       = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr);
751eaedb033SPeter Brune   snes->data = (void*)nasm;
752eaedb033SPeter Brune 
753eaedb033SPeter Brune   nasm->n        = PETSC_DECIDE;
754eaedb033SPeter Brune   nasm->subsnes  = 0;
755eaedb033SPeter Brune   nasm->x        = 0;
756111ade9eSPeter Brune   nasm->xl       = 0;
757111ade9eSPeter Brune   nasm->y        = 0;
758eaedb033SPeter Brune   nasm->b        = 0;
759111ade9eSPeter Brune   nasm->oscatter = 0;
760111ade9eSPeter Brune   nasm->iscatter = 0;
761111ade9eSPeter Brune   nasm->gscatter = 0;
762111ade9eSPeter Brune 
763111ade9eSPeter Brune   nasm->type = PC_ASM_BASIC;
764d728fb7dSPeter Brune   nasm->finaljacobian = PETSC_FALSE;
765*a4f17876SPeter Brune   nasm->same_local_solves = PETSC_TRUE;
766eaedb033SPeter Brune 
767eaedb033SPeter Brune   snes->ops->destroy        = SNESDestroy_NASM;
768eaedb033SPeter Brune   snes->ops->setup          = SNESSetUp_NASM;
769eaedb033SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
770eaedb033SPeter Brune   snes->ops->view           = SNESView_NASM;
771eaedb033SPeter Brune   snes->ops->solve          = SNESSolve_NASM;
772eaedb033SPeter Brune   snes->ops->reset          = SNESReset_NASM;
773eaedb033SPeter Brune 
774eaedb033SPeter Brune   snes->usesksp = PETSC_FALSE;
775eaedb033SPeter Brune   snes->usespc  = PETSC_FALSE;
776eaedb033SPeter Brune 
7770298fd71SBarry Smith   nasm->eventrestrictinterp = 0;
7780298fd71SBarry Smith   nasm->eventsubsolve       = 0;
779b20c023fSPeter Brune 
780eaedb033SPeter Brune   if (!snes->tolerancesset) {
781eaedb033SPeter Brune     snes->max_its   = 10000;
782eaedb033SPeter Brune     snes->max_funcs = 10000;
783eaedb033SPeter Brune   }
784eaedb033SPeter Brune 
785bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr);
786bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr);
787bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr);
788bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr);
789eaedb033SPeter Brune   PetscFunctionReturn(0);
790eaedb033SPeter Brune }
79199e0435eSBarry Smith 
792