xref: /petsc/src/snes/impls/nasm/nasm.c (revision 95452b02e12c0ee11232c7ff2b24b568a8e07e43)
1af0996ceSBarry Smith #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 */
11f10b3e88SPatrick Farrell   Vec        weight;              /* weighting for adding updates on overlaps, in global space */
12111ade9eSPeter Brune   VecScatter *oscatter;           /* scatter from global space to the subdomain global space */
13f10b3e88SPatrick Farrell   VecScatter *oscatter_copy;      /* copy of the above */
14111ade9eSPeter Brune   VecScatter *iscatter;           /* scatter from global space to the nonoverlapping subdomain space */
15111ade9eSPeter Brune   VecScatter *gscatter;           /* scatter from global space to the subdomain local space */
16111ade9eSPeter Brune   PCASMType  type;                /* ASM type */
17111ade9eSPeter Brune   PetscBool  usesdm;              /* use the DM for setting up the subproblems */
18d728fb7dSPeter Brune   PetscBool  finaljacobian;       /* compute the jacobian of the converged solution */
19610116beSPeter Brune   PetscReal  damping;             /* damping parameter for updates from the blocks */
20a4f17876SPeter Brune   PetscBool  same_local_solves;   /* flag to determine if the solvers have been individually modified */
21f10b3e88SPatrick Farrell   PetscBool  weight_set;          /* use a weight in the overlap updates */
22b20c023fSPeter Brune 
23b20c023fSPeter Brune   /* logging events */
24b20c023fSPeter Brune   PetscLogEvent eventrestrictinterp;
25b20c023fSPeter Brune   PetscLogEvent eventsubsolve;
26602bec5dSPeter Brune 
27602bec5dSPeter Brune   PetscInt      fjtype;            /* type of computed jacobian */
28602bec5dSPeter Brune   Vec           xinit;             /* initial solution in case the final jacobian type is computed as first */
29eaedb033SPeter Brune } SNES_NASM;
30eaedb033SPeter Brune 
31b20c023fSPeter Brune const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0};
32602bec5dSPeter Brune const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"};
33b20c023fSPeter Brune 
3440244768SBarry Smith static PetscErrorCode SNESReset_NASM(SNES snes)
35eaedb033SPeter Brune {
36eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
37eaedb033SPeter Brune   PetscErrorCode ierr;
38eaedb033SPeter Brune   PetscInt       i;
396e111a19SKarl Rupp 
40eaedb033SPeter Brune   PetscFunctionBegin;
41eaedb033SPeter Brune   for (i=0; i<nasm->n; i++) {
42111ade9eSPeter Brune     if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); }
43f5f7c1b9SKarl Rupp     if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); }
44111ade9eSPeter Brune     if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); }
45bc8c1f72SJose Roman     if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); }
46eaedb033SPeter Brune 
47bc8c1f72SJose Roman     if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); }
48111ade9eSPeter Brune     if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); }
49f10b3e88SPatrick Farrell     if (nasm->oscatter_copy) { ierr = VecScatterDestroy(&nasm->oscatter_copy[i]);CHKERRQ(ierr); }
50111ade9eSPeter Brune     if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); }
51111ade9eSPeter Brune     if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); }
52eaedb033SPeter Brune   }
53111ade9eSPeter Brune 
540677ede6SBarry Smith   ierr = PetscFree(nasm->x);CHKERRQ(ierr);
550677ede6SBarry Smith   ierr = PetscFree(nasm->xl);CHKERRQ(ierr);
560677ede6SBarry Smith   ierr = PetscFree(nasm->y);CHKERRQ(ierr);
570677ede6SBarry Smith   ierr = PetscFree(nasm->b);CHKERRQ(ierr);
58111ade9eSPeter Brune 
59602bec5dSPeter Brune   if (nasm->xinit) {ierr = VecDestroy(&nasm->xinit);CHKERRQ(ierr);}
60602bec5dSPeter Brune 
610677ede6SBarry Smith   ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);
620677ede6SBarry Smith   ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);
63f10b3e88SPatrick Farrell   ierr = PetscFree(nasm->oscatter_copy);CHKERRQ(ierr);
640677ede6SBarry Smith   ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);
650677ede6SBarry Smith   ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);
66b20c023fSPeter Brune 
67f10b3e88SPatrick Farrell   if (nasm->weight_set) {
68f10b3e88SPatrick Farrell     ierr = VecDestroy(&nasm->weight);CHKERRQ(ierr);
69f10b3e88SPatrick Farrell   }
70f10b3e88SPatrick Farrell 
71b20c023fSPeter Brune   nasm->eventrestrictinterp = 0;
72b20c023fSPeter Brune   nasm->eventsubsolve = 0;
73eaedb033SPeter Brune   PetscFunctionReturn(0);
74eaedb033SPeter Brune }
75eaedb033SPeter Brune 
7640244768SBarry Smith static PetscErrorCode SNESDestroy_NASM(SNES snes)
77eaedb033SPeter Brune {
78eaedb033SPeter Brune   PetscErrorCode ierr;
796e111a19SKarl Rupp 
80eaedb033SPeter Brune   PetscFunctionBegin;
81eaedb033SPeter Brune   ierr = SNESReset_NASM(snes);CHKERRQ(ierr);
8222d28d08SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
83eaedb033SPeter Brune   PetscFunctionReturn(0);
84eaedb033SPeter Brune }
85eaedb033SPeter Brune 
8640244768SBarry Smith static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx)
870adebc6cSBarry Smith {
88111ade9eSPeter Brune   PetscErrorCode ierr;
89111ade9eSPeter Brune   Vec            bcs = (Vec)ctx;
906e111a19SKarl Rupp 
91111ade9eSPeter Brune   PetscFunctionBegin;
92111ade9eSPeter Brune   ierr = VecCopy(bcs,l);CHKERRQ(ierr);
93111ade9eSPeter Brune   PetscFunctionReturn(0);
94111ade9eSPeter Brune }
95111ade9eSPeter Brune 
9640244768SBarry Smith static PetscErrorCode SNESSetUp_NASM(SNES snes)
97eaedb033SPeter Brune {
98eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
99eaedb033SPeter Brune   PetscErrorCode ierr;
10076857b2aSPeter Brune   DM             dm,subdm;
101111ade9eSPeter Brune   DM             *subdms;
102111ade9eSPeter Brune   PetscInt       i;
103eaedb033SPeter Brune   const char     *optionsprefix;
104111ade9eSPeter Brune   Vec            F;
105ed3c11a9SPeter Brune   PetscMPIInt    size;
106ed3c11a9SPeter Brune   KSP            ksp;
107ed3c11a9SPeter Brune   PC             pc;
108eaedb033SPeter Brune 
109eaedb033SPeter Brune   PetscFunctionBegin;
110eaedb033SPeter Brune   if (!nasm->subsnes) {
111eaedb033SPeter Brune     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1120a696f66SPeter Brune     if (dm) {
113eaedb033SPeter Brune       nasm->usesdm = PETSC_TRUE;
1140298fd71SBarry Smith       ierr         = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr);
115ce94432eSBarry Smith       if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
116111ade9eSPeter Brune       ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr);
117f10b3e88SPatrick Farrell       ierr = PetscMalloc1(nasm->n, &nasm->oscatter_copy);CHKERRQ(ierr);
118f10b3e88SPatrick Farrell       for (i=0; i<nasm->n; i++) {
119f10b3e88SPatrick Farrell         ierr = VecScatterCopy(nasm->oscatter[i], &nasm->oscatter_copy[i]);CHKERRQ(ierr);
120f10b3e88SPatrick Farrell       }
121eaedb033SPeter Brune 
122eaedb033SPeter Brune       ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
123785e854fSJed Brown       ierr = PetscMalloc1(nasm->n,&nasm->subsnes);CHKERRQ(ierr);
124111ade9eSPeter Brune       for (i=0; i<nasm->n; i++) {
125cdb298fcSPeter Brune         ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr);
1267a2f0ea0SPatrick Farrell         ierr = PetscObjectIncrementTabLevel((PetscObject)nasm->subsnes[i], (PetscObject)snes, 1);CHKERRQ(ierr);
127cdb298fcSPeter Brune         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr);
128cdb298fcSPeter Brune         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr);
129cdb298fcSPeter Brune         ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr);
130ed3c11a9SPeter Brune         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRQ(ierr);
131ed3c11a9SPeter Brune         if (size == 1) {
132ed3c11a9SPeter Brune           ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr);
133ed3c11a9SPeter Brune           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
134ed3c11a9SPeter Brune           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
135ed3c11a9SPeter Brune           ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
136ed3c11a9SPeter Brune         }
137cdb298fcSPeter Brune         ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr);
138111ade9eSPeter Brune         ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr);
139111ade9eSPeter Brune       }
140111ade9eSPeter Brune       ierr = PetscFree(subdms);CHKERRQ(ierr);
141ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!");
142ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!");
143111ade9eSPeter Brune   /* allocate the global vectors */
144e245e0aaSPeter Brune   if (!nasm->x) {
1451795a4d1SJed Brown     ierr = PetscCalloc1(nasm->n,&nasm->x);CHKERRQ(ierr);
146e245e0aaSPeter Brune   }
147e245e0aaSPeter Brune   if (!nasm->xl) {
1481795a4d1SJed Brown     ierr = PetscCalloc1(nasm->n,&nasm->xl);CHKERRQ(ierr);
149e245e0aaSPeter Brune   }
150e245e0aaSPeter Brune   if (!nasm->y) {
1511795a4d1SJed Brown     ierr = PetscCalloc1(nasm->n,&nasm->y);CHKERRQ(ierr);
152e245e0aaSPeter Brune   }
153e245e0aaSPeter Brune   if (!nasm->b) {
1541795a4d1SJed Brown     ierr = PetscCalloc1(nasm->n,&nasm->b);CHKERRQ(ierr);
155e245e0aaSPeter Brune   }
156111ade9eSPeter Brune 
157111ade9eSPeter Brune   for (i=0; i<nasm->n; i++) {
1580298fd71SBarry Smith     ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr);
15976857b2aSPeter Brune     if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);}
16076857b2aSPeter Brune     if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);}
16176857b2aSPeter Brune     if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);}
16276857b2aSPeter Brune     if (!nasm->xl[i]) {
163111ade9eSPeter Brune       ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr);
164111ade9eSPeter Brune       ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr);
1650298fd71SBarry Smith       ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr);
166111ade9eSPeter Brune     }
16761ba4676SBarry Smith   }
168602bec5dSPeter Brune   if (nasm->finaljacobian) {
169602bec5dSPeter Brune     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
170602bec5dSPeter Brune     if (nasm->fjtype == 2) {
171602bec5dSPeter Brune       ierr = VecDuplicate(snes->vec_sol,&nasm->xinit);CHKERRQ(ierr);
172602bec5dSPeter Brune     }
173602bec5dSPeter Brune     for (i=0; i<nasm->n;i++) {
174602bec5dSPeter Brune       ierr = SNESSetUpMatrices(nasm->subsnes[i]);CHKERRQ(ierr);
175602bec5dSPeter Brune     }
176602bec5dSPeter Brune   }
177eaedb033SPeter Brune   PetscFunctionReturn(0);
178eaedb033SPeter Brune }
179eaedb033SPeter Brune 
18040244768SBarry Smith static PetscErrorCode SNESSetFromOptions_NASM(PetscOptionItems *PetscOptionsObject,SNES snes)
181eaedb033SPeter Brune {
182eaedb033SPeter Brune   PetscErrorCode    ierr;
183111ade9eSPeter Brune   PCASMType         asmtype;
184a4f17876SPeter Brune   PetscBool         flg,monflg,subviewflg;
185111ade9eSPeter Brune   SNES_NASM         *nasm = (SNES_NASM*)snes->data;
1866e111a19SKarl Rupp 
187eaedb033SPeter Brune   PetscFunctionBegin;
188e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Nonlinear Additive Schwartz options");CHKERRQ(ierr);
189111ade9eSPeter Brune   ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
190e0331734SPeter Brune   if (flg) {ierr = SNESNASMSetType(snes,asmtype);CHKERRQ(ierr);}
191b20c023fSPeter Brune   flg    = PETSC_FALSE;
192b20c023fSPeter Brune   monflg = PETSC_TRUE;
1935dfa0f3bSBarry Smith   ierr   = PetscOptionsReal("-snes_nasm_damping","The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)","SNESNASMSetDamping",nasm->damping,&nasm->damping,&flg);CHKERRQ(ierr);
194610116beSPeter Brune   if (flg) {ierr = SNESNASMSetDamping(snes,nasm->damping);CHKERRQ(ierr);}
195a4f17876SPeter Brune   subviewflg = PETSC_FALSE;
196a4f17876SPeter Brune   ierr   = PetscOptionsBool("-snes_nasm_sub_view","Print detailed information for every processor when using -snes_view","",subviewflg,&subviewflg,&flg);CHKERRQ(ierr);
197a4f17876SPeter Brune   if (flg) {
198a4f17876SPeter Brune     nasm->same_local_solves = PETSC_FALSE;
199a4f17876SPeter Brune     if (!subviewflg) {
200a4f17876SPeter Brune       nasm->same_local_solves = PETSC_TRUE;
201a4f17876SPeter Brune     }
202a4f17876SPeter Brune   }
203d728fb7dSPeter Brune   ierr   = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr);
204602bec5dSPeter Brune   ierr   = PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL);CHKERRQ(ierr);
205a4f17876SPeter Brune   ierr   = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr);
206b20c023fSPeter Brune   if (flg) {
207b20c023fSPeter Brune     ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr);
208b20c023fSPeter Brune     ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr);
209b20c023fSPeter Brune   }
210eaedb033SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
211eaedb033SPeter Brune   PetscFunctionReturn(0);
212eaedb033SPeter Brune }
213eaedb033SPeter Brune 
21440244768SBarry Smith static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
215eaedb033SPeter Brune {
216b20c023fSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
217b20c023fSPeter Brune   PetscErrorCode ierr;
218a4f17876SPeter Brune   PetscMPIInt    rank,size;
219dd2fa690SBarry Smith   PetscInt       i,N,bsz;
220b20c023fSPeter Brune   PetscBool      iascii,isstring;
221b20c023fSPeter Brune   PetscViewer    sviewer;
222ce94432eSBarry Smith   MPI_Comm       comm;
223b20c023fSPeter Brune 
224b20c023fSPeter Brune   PetscFunctionBegin;
225ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
226b20c023fSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
227b20c023fSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
228b20c023fSPeter Brune   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
229a4f17876SPeter Brune   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
230b2566f29SBarry Smith   ierr = MPIU_Allreduce(&nasm->n,&N,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
231b20c023fSPeter Brune   if (iascii) {
232efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer, "  total subdomain blocks = %D\n",N);CHKERRQ(ierr);
233a4f17876SPeter Brune     if (nasm->same_local_solves) {
234a4f17876SPeter Brune       if (nasm->subsnes) {
235a4f17876SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  Local solve is the same for all blocks:\n");CHKERRQ(ierr);
236a4f17876SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
2373f08860eSBarry Smith         ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
238a4f17876SPeter Brune         if (!rank) {
239a4f17876SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
240a4f17876SPeter Brune           ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);
241a4f17876SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
242a4f17876SPeter Brune         }
2433f08860eSBarry Smith         ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
244a4f17876SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
245a4f17876SPeter Brune       }
246a4f17876SPeter Brune     } else {
247a4f17876SPeter Brune       /* print the solver on each block */
2481575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
249b20c023fSPeter Brune       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr);
250b20c023fSPeter Brune       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
2511575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
252a4f17876SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following SNES objects:\n");CHKERRQ(ierr);
253b20c023fSPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
254a4f17876SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
2553f08860eSBarry Smith       ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
256a4f17876SPeter Brune       for (i=0; i<nasm->n; i++) {
257a4f17876SPeter Brune         ierr = VecGetLocalSize(nasm->x[i],&bsz);CHKERRQ(ierr);
258a4f17876SPeter Brune         ierr = PetscViewerASCIIPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr);
259b20c023fSPeter Brune         ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr);
260a4f17876SPeter Brune         ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
261b20c023fSPeter Brune       }
2623f08860eSBarry Smith       ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
263a4f17876SPeter Brune       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
264b20c023fSPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
265a4f17876SPeter Brune     }
266b20c023fSPeter Brune   } else if (isstring) {
267b20c023fSPeter Brune     ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr);
2683f08860eSBarry Smith     ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
269b20c023fSPeter Brune     if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);}
2703f08860eSBarry Smith     ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
271b20c023fSPeter Brune   }
272eaedb033SPeter Brune   PetscFunctionReturn(0);
273eaedb033SPeter Brune }
274eaedb033SPeter Brune 
275e0331734SPeter Brune /*@
276e0331734SPeter Brune    SNESNASMSetType - Set the type of subdomain update used
277e0331734SPeter Brune 
278e0331734SPeter Brune    Logically Collective on SNES
279e0331734SPeter Brune 
280e0331734SPeter Brune    Input Parameters:
281e0331734SPeter Brune +  SNES - the SNES context
282e0331734SPeter Brune -  type - the type of update, PC_ASM_BASIC or PC_ASM_RESTRICT
283e0331734SPeter Brune 
284e0331734SPeter Brune    Level: intermediate
285e0331734SPeter Brune 
286e0331734SPeter Brune .keywords: SNES, NASM
287e0331734SPeter Brune 
288e0331734SPeter Brune .seealso: SNESNASM, SNESNASMGetType(), PCASMSetType()
289e0331734SPeter Brune @*/
290e0331734SPeter Brune PetscErrorCode SNESNASMSetType(SNES snes,PCASMType type)
291e0331734SPeter Brune {
292e0331734SPeter Brune   PetscErrorCode ierr;
293e0331734SPeter Brune   PetscErrorCode (*f)(SNES,PCASMType);
294e0331734SPeter Brune 
295e0331734SPeter Brune   PetscFunctionBegin;
296e0331734SPeter Brune   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetType_C",&f);CHKERRQ(ierr);
297e0331734SPeter Brune   if (f) {ierr = (f)(snes,type);CHKERRQ(ierr);}
298e0331734SPeter Brune   PetscFunctionReturn(0);
299e0331734SPeter Brune }
300e0331734SPeter Brune 
30140244768SBarry Smith static PetscErrorCode SNESNASMSetType_NASM(SNES snes,PCASMType type)
302e0331734SPeter Brune {
303e0331734SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
304e0331734SPeter Brune 
305e0331734SPeter Brune   PetscFunctionBegin;
306e0331734SPeter Brune   if (type != PC_ASM_BASIC && type != PC_ASM_RESTRICT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESNASM only supports basic and restrict types");
307e0331734SPeter Brune   nasm->type = type;
308e0331734SPeter Brune   PetscFunctionReturn(0);
309e0331734SPeter Brune }
310e0331734SPeter Brune 
311e0331734SPeter Brune /*@
312e0331734SPeter Brune    SNESNASMGetType - Get the type of subdomain update used
313e0331734SPeter Brune 
314e0331734SPeter Brune    Logically Collective on SNES
315e0331734SPeter Brune 
316e0331734SPeter Brune    Input Parameters:
317e0331734SPeter Brune .  SNES - the SNES context
318e0331734SPeter Brune 
319e0331734SPeter Brune    Output Parameters:
320e0331734SPeter Brune .  type - the type of update
321e0331734SPeter Brune 
322e0331734SPeter Brune    Level: intermediate
323e0331734SPeter Brune 
324e0331734SPeter Brune .keywords: SNES, NASM
325e0331734SPeter Brune 
326e0331734SPeter Brune .seealso: SNESNASM, SNESNASMSetType(), PCASMGetType()
327e0331734SPeter Brune @*/
328e0331734SPeter Brune PetscErrorCode SNESNASMGetType(SNES snes,PCASMType *type)
329e0331734SPeter Brune {
330e0331734SPeter Brune   PetscErrorCode ierr;
331e0331734SPeter Brune 
332e0331734SPeter Brune   PetscFunctionBegin;
3332a808120SBarry Smith   ierr = PetscUseMethod(snes,"SNESNASMGetType_C",(SNES,PCASMType*),(snes,type));CHKERRQ(ierr);
334e0331734SPeter Brune   PetscFunctionReturn(0);
335e0331734SPeter Brune }
336e0331734SPeter Brune 
33740244768SBarry Smith static PetscErrorCode SNESNASMGetType_NASM(SNES snes,PCASMType *type)
338e0331734SPeter Brune {
339e0331734SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
340e0331734SPeter Brune 
341e0331734SPeter Brune   PetscFunctionBegin;
342e0331734SPeter Brune   *type = nasm->type;
343e0331734SPeter Brune   PetscFunctionReturn(0);
344e0331734SPeter Brune }
345e0331734SPeter Brune 
34676857b2aSPeter Brune /*@
34776857b2aSPeter Brune    SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems.
34876857b2aSPeter Brune 
34976857b2aSPeter Brune    Not Collective
35076857b2aSPeter Brune 
35176857b2aSPeter Brune    Input Parameters:
35276857b2aSPeter Brune +  SNES - the SNES context
35376857b2aSPeter Brune .  n - the number of local subdomains
35476857b2aSPeter Brune .  subsnes - solvers defined on the local subdomains
35576857b2aSPeter Brune .  iscatter - scatters into the nonoverlapping portions of the local subdomains
35676857b2aSPeter Brune .  oscatter - scatters into the overlapping portions of the local subdomains
35776857b2aSPeter Brune -  gscatter - scatters into the (ghosted) local vector of the local subdomain
35876857b2aSPeter Brune 
35976857b2aSPeter Brune    Level: intermediate
36076857b2aSPeter Brune 
36176857b2aSPeter Brune .keywords: SNES, NASM
36276857b2aSPeter Brune 
36376857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains()
36476857b2aSPeter Brune @*/
365a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
366a6dfd86eSKarl Rupp {
367eaedb033SPeter Brune   PetscErrorCode ierr;
368111ade9eSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*);
3696e111a19SKarl Rupp 
370eaedb033SPeter Brune   PetscFunctionBegin;
3710005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f);CHKERRQ(ierr);
3724b838e8fSPeter Brune   if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);}
373eaedb033SPeter Brune   PetscFunctionReturn(0);
374eaedb033SPeter Brune }
375eaedb033SPeter Brune 
37640244768SBarry Smith static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
377a6dfd86eSKarl Rupp {
378eaedb033SPeter Brune   PetscInt       i;
379eaedb033SPeter Brune   PetscErrorCode ierr;
380eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
3816e111a19SKarl Rupp 
382eaedb033SPeter Brune   PetscFunctionBegin;
383ce94432eSBarry Smith   if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
384eaedb033SPeter Brune 
385111ade9eSPeter Brune   /* tear down the previously set things */
386111ade9eSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
387111ade9eSPeter Brune 
388eaedb033SPeter Brune   nasm->n = n;
389111ade9eSPeter Brune   if (oscatter) {
390111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);}
391eaedb033SPeter Brune   }
392111ade9eSPeter Brune   if (iscatter) {
393111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);}
394eaedb033SPeter Brune   }
395111ade9eSPeter Brune   if (gscatter) {
396111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);}
397111ade9eSPeter Brune   }
398111ade9eSPeter Brune   if (oscatter) {
399785e854fSJed Brown     ierr = PetscMalloc1(n,&nasm->oscatter);CHKERRQ(ierr);
400f10b3e88SPatrick Farrell     ierr = PetscMalloc1(n,&nasm->oscatter_copy);CHKERRQ(ierr);
401eaedb033SPeter Brune     for (i=0; i<n; i++) {
402111ade9eSPeter Brune       nasm->oscatter[i] = oscatter[i];
403f10b3e88SPatrick Farrell       ierr = VecScatterCopy(oscatter[i], &nasm->oscatter_copy[i]);CHKERRQ(ierr);
404eaedb033SPeter Brune     }
405111ade9eSPeter Brune   }
406111ade9eSPeter Brune   if (iscatter) {
407785e854fSJed Brown     ierr = PetscMalloc1(n,&nasm->iscatter);CHKERRQ(ierr);
408eaedb033SPeter Brune     for (i=0; i<n; i++) {
409111ade9eSPeter Brune       nasm->iscatter[i] = iscatter[i];
410eaedb033SPeter Brune     }
411eaedb033SPeter Brune   }
412111ade9eSPeter Brune   if (gscatter) {
413785e854fSJed Brown     ierr = PetscMalloc1(n,&nasm->gscatter);CHKERRQ(ierr);
414eaedb033SPeter Brune     for (i=0; i<n; i++) {
415111ade9eSPeter Brune       nasm->gscatter[i] = gscatter[i];
416eaedb033SPeter Brune     }
417eaedb033SPeter Brune   }
418111ade9eSPeter Brune 
419eaedb033SPeter Brune   if (subsnes) {
420785e854fSJed Brown     ierr = PetscMalloc1(n,&nasm->subsnes);CHKERRQ(ierr);
421eaedb033SPeter Brune     for (i=0; i<n; i++) {
422eaedb033SPeter Brune       nasm->subsnes[i] = subsnes[i];
423eaedb033SPeter Brune     }
424a4f17876SPeter Brune     nasm->same_local_solves = PETSC_FALSE;
425eaedb033SPeter Brune   }
426eaedb033SPeter Brune   PetscFunctionReturn(0);
427eaedb033SPeter Brune }
428eaedb033SPeter Brune 
42976857b2aSPeter Brune /*@
43076857b2aSPeter Brune    SNESNASMGetSubdomains - Get the local subdomain context.
43176857b2aSPeter Brune 
43276857b2aSPeter Brune    Not Collective
43376857b2aSPeter Brune 
43476857b2aSPeter Brune    Input Parameters:
43576857b2aSPeter Brune .  SNES - the SNES context
43676857b2aSPeter Brune 
43776857b2aSPeter Brune    Output Parameters:
43876857b2aSPeter Brune +  n - the number of local subdomains
43976857b2aSPeter Brune .  subsnes - solvers defined on the local subdomains
44076857b2aSPeter Brune .  iscatter - scatters into the nonoverlapping portions of the local subdomains
44176857b2aSPeter Brune .  oscatter - scatters into the overlapping portions of the local subdomains
44276857b2aSPeter Brune -  gscatter - scatters into the (ghosted) local vector of the local subdomain
44376857b2aSPeter Brune 
44476857b2aSPeter Brune    Level: intermediate
44576857b2aSPeter Brune 
44676857b2aSPeter Brune .keywords: SNES, NASM
44776857b2aSPeter Brune 
44876857b2aSPeter Brune .seealso: SNESNASM, SNESNASMSetSubdomains()
44976857b2aSPeter Brune @*/
45076857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
45176857b2aSPeter Brune {
45276857b2aSPeter Brune   PetscErrorCode ierr;
45376857b2aSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**);
45476857b2aSPeter Brune 
45576857b2aSPeter Brune   PetscFunctionBegin;
4560005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f);CHKERRQ(ierr);
4574b838e8fSPeter Brune   if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);}
45876857b2aSPeter Brune   PetscFunctionReturn(0);
45976857b2aSPeter Brune }
46076857b2aSPeter Brune 
46140244768SBarry Smith static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
46276857b2aSPeter Brune {
46376857b2aSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
46476857b2aSPeter Brune 
46576857b2aSPeter Brune   PetscFunctionBegin;
46676857b2aSPeter Brune   if (n) *n = nasm->n;
46776857b2aSPeter Brune   if (oscatter) *oscatter = nasm->oscatter;
46876857b2aSPeter Brune   if (iscatter) *iscatter = nasm->iscatter;
46976857b2aSPeter Brune   if (gscatter) *gscatter = nasm->gscatter;
470a4f17876SPeter Brune   if (subsnes)  {
471a4f17876SPeter Brune     *subsnes  = nasm->subsnes;
472a4f17876SPeter Brune     nasm->same_local_solves = PETSC_FALSE;
473a4f17876SPeter Brune   }
47476857b2aSPeter Brune   PetscFunctionReturn(0);
47576857b2aSPeter Brune }
47676857b2aSPeter Brune 
47776857b2aSPeter Brune /*@
47876857b2aSPeter Brune    SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors
47976857b2aSPeter Brune 
48076857b2aSPeter Brune    Not Collective
48176857b2aSPeter Brune 
48276857b2aSPeter Brune    Input Parameters:
48376857b2aSPeter Brune .  SNES - the SNES context
48476857b2aSPeter Brune 
48576857b2aSPeter Brune    Output Parameters:
48676857b2aSPeter Brune +  n - the number of local subdomains
48776857b2aSPeter Brune .  x - The subdomain solution vector
48876857b2aSPeter Brune .  y - The subdomain step vector
48976857b2aSPeter Brune .  b - The subdomain RHS vector
49076857b2aSPeter Brune -  xl - The subdomain local vectors (ghosted)
49176857b2aSPeter Brune 
49276857b2aSPeter Brune    Level: developer
49376857b2aSPeter Brune 
49476857b2aSPeter Brune .keywords: SNES, NASM
49576857b2aSPeter Brune 
49676857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains()
49776857b2aSPeter Brune @*/
49876857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl)
49976857b2aSPeter Brune {
50076857b2aSPeter Brune   PetscErrorCode ierr;
50176857b2aSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**);
50276857b2aSPeter Brune 
50376857b2aSPeter Brune   PetscFunctionBegin;
5040005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f);CHKERRQ(ierr);
5054b838e8fSPeter Brune   if (f) {ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);}
50676857b2aSPeter Brune   PetscFunctionReturn(0);
50776857b2aSPeter Brune }
50876857b2aSPeter Brune 
50940244768SBarry Smith static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl)
51076857b2aSPeter Brune {
51176857b2aSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
51276857b2aSPeter Brune 
51376857b2aSPeter Brune   PetscFunctionBegin;
51476857b2aSPeter Brune   if (n)  *n  = nasm->n;
51576857b2aSPeter Brune   if (x)  *x  = nasm->x;
51676857b2aSPeter Brune   if (y)  *y  = nasm->y;
51776857b2aSPeter Brune   if (b)  *b  = nasm->b;
51876857b2aSPeter Brune   if (xl) *xl = nasm->xl;
51976857b2aSPeter Brune   PetscFunctionReturn(0);
52076857b2aSPeter Brune }
52176857b2aSPeter Brune 
522d728fb7dSPeter Brune /*@
523d728fb7dSPeter Brune    SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence
524d728fb7dSPeter Brune 
525d728fb7dSPeter Brune    Collective on SNES
526d728fb7dSPeter Brune 
527d728fb7dSPeter Brune    Input Parameters:
528d728fb7dSPeter Brune +  SNES - the SNES context
529d728fb7dSPeter Brune -  flg - indication of whether to compute the jacobians or not
530d728fb7dSPeter Brune 
531d728fb7dSPeter Brune    Level: developer
532d728fb7dSPeter Brune 
533*95452b02SPatrick Sanan    Notes:
534*95452b02SPatrick Sanan     This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian
535d728fb7dSPeter Brune    is needed at each linear iteration.
536d728fb7dSPeter Brune 
537d728fb7dSPeter Brune .keywords: SNES, NASM, ASPIN
538d728fb7dSPeter Brune 
539d728fb7dSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains()
540d728fb7dSPeter Brune @*/
541d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg)
542d728fb7dSPeter Brune {
543d728fb7dSPeter Brune   PetscErrorCode (*f)(SNES,PetscBool);
544d728fb7dSPeter Brune   PetscErrorCode ierr;
545d728fb7dSPeter Brune 
546d728fb7dSPeter Brune   PetscFunctionBegin;
5470005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f);CHKERRQ(ierr);
5484b838e8fSPeter Brune   if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);}
549d728fb7dSPeter Brune   PetscFunctionReturn(0);
550d728fb7dSPeter Brune }
551d728fb7dSPeter Brune 
55240244768SBarry Smith static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg)
553d728fb7dSPeter Brune {
554d728fb7dSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
555d728fb7dSPeter Brune 
556d728fb7dSPeter Brune   PetscFunctionBegin;
557d728fb7dSPeter Brune   nasm->finaljacobian = flg;
558602bec5dSPeter Brune   if (flg) snes->usesksp = PETSC_TRUE;
559d728fb7dSPeter Brune   PetscFunctionReturn(0);
560d728fb7dSPeter Brune }
56176857b2aSPeter Brune 
562610116beSPeter Brune /*@
563610116beSPeter Brune    SNESNASMSetDamping - Sets the update damping for NASM
564610116beSPeter Brune 
565610116beSPeter Brune    Logically collective on SNES
566610116beSPeter Brune 
567610116beSPeter Brune    Input Parameters:
568610116beSPeter Brune +  SNES - the SNES context
569610116beSPeter Brune -  dmp - damping
570610116beSPeter Brune 
571610116beSPeter Brune    Level: intermediate
572610116beSPeter Brune 
573*95452b02SPatrick Sanan    Notes:
574*95452b02SPatrick Sanan     The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
5755dfa0f3bSBarry Smith 
576610116beSPeter Brune .keywords: SNES, NASM, damping
577610116beSPeter Brune 
578610116beSPeter Brune .seealso: SNESNASM, SNESNASMGetDamping()
579610116beSPeter Brune @*/
580610116beSPeter Brune PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp)
581610116beSPeter Brune {
582610116beSPeter Brune   PetscErrorCode (*f)(SNES,PetscReal);
583610116beSPeter Brune   PetscErrorCode ierr;
584610116beSPeter Brune 
585610116beSPeter Brune   PetscFunctionBegin;
586610116beSPeter Brune   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr);
587e24b0d94SPeter Brune   if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);}
588610116beSPeter Brune   PetscFunctionReturn(0);
589610116beSPeter Brune }
590610116beSPeter Brune 
59140244768SBarry Smith static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp)
592610116beSPeter Brune {
593610116beSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
594610116beSPeter Brune 
595610116beSPeter Brune   PetscFunctionBegin;
596610116beSPeter Brune   nasm->damping = dmp;
597610116beSPeter Brune   PetscFunctionReturn(0);
598610116beSPeter Brune }
599610116beSPeter Brune 
600610116beSPeter Brune /*@
601610116beSPeter Brune    SNESNASMGetDamping - Gets the update damping for NASM
602610116beSPeter Brune 
603610116beSPeter Brune    Not Collective
604610116beSPeter Brune 
605610116beSPeter Brune    Input Parameters:
606610116beSPeter Brune +  SNES - the SNES context
607610116beSPeter Brune -  dmp - damping
608610116beSPeter Brune 
609610116beSPeter Brune    Level: intermediate
610610116beSPeter Brune 
611610116beSPeter Brune .keywords: SNES, NASM, damping
612610116beSPeter Brune 
613610116beSPeter Brune .seealso: SNESNASM, SNESNASMSetDamping()
614610116beSPeter Brune @*/
615610116beSPeter Brune PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp)
616610116beSPeter Brune {
617610116beSPeter Brune   PetscErrorCode ierr;
618610116beSPeter Brune 
619610116beSPeter Brune   PetscFunctionBegin;
6204a2f8832SBarry Smith   ierr = PetscUseMethod(snes,"SNESNASMGetDamping_C",(SNES,PetscReal*),(snes,dmp));CHKERRQ(ierr);
621610116beSPeter Brune   PetscFunctionReturn(0);
622610116beSPeter Brune }
623610116beSPeter Brune 
62440244768SBarry Smith static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp)
625610116beSPeter Brune {
626610116beSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
627610116beSPeter Brune 
628610116beSPeter Brune   PetscFunctionBegin;
629610116beSPeter Brune   *dmp = nasm->damping;
630610116beSPeter Brune   PetscFunctionReturn(0);
631610116beSPeter Brune }
632610116beSPeter Brune 
633610116beSPeter Brune 
63414eb1c5cSMatthew G. Knepley /*
63514eb1c5cSMatthew G. Knepley   Input Parameters:
63614eb1c5cSMatthew G. Knepley + snes - The solver
63714eb1c5cSMatthew G. Knepley . B - The RHS vector
63814eb1c5cSMatthew G. Knepley - X - The initial guess
63914eb1c5cSMatthew G. Knepley 
64014eb1c5cSMatthew G. Knepley   Output Parameters:
64114eb1c5cSMatthew G. Knepley . Y - The solution update
64214eb1c5cSMatthew G. Knepley 
64314eb1c5cSMatthew G. Knepley   TODO: All scatters should be packed into one
64414eb1c5cSMatthew G. Knepley */
6450adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X)
6460adebc6cSBarry Smith {
647eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
648258e1594SPeter Brune   SNES           subsnes;
649eaedb033SPeter Brune   PetscInt       i;
650610116beSPeter Brune   PetscReal      dmp;
651eaedb033SPeter Brune   PetscErrorCode ierr;
652f10b3e88SPatrick Farrell   Vec            Xl,Bl,Yl,Xlloc;
653f10b3e88SPatrick Farrell   VecScatter     iscat,oscat,gscat,oscat_copy;
654111ade9eSPeter Brune   DM             dm,subdm;
655e0331734SPeter Brune   PCASMType      type;
6560adebc6cSBarry Smith 
657eaedb033SPeter Brune   PetscFunctionBegin;
658e0331734SPeter Brune   ierr = SNESNASMGetType(snes,&type);CHKERRQ(ierr);
659eaedb033SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
660111ade9eSPeter Brune   ierr = VecSet(Y,0);CHKERRQ(ierr);
661b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
662eaedb033SPeter Brune   for (i=0; i<nasm->n; i++) {
663f10b3e88SPatrick Farrell     /* scatter the solution to the global solution and the local solution */
664f10b3e88SPatrick Farrell     Xl      = nasm->x[i];
66570c78f05SPeter Brune     Xlloc   = nasm->xl[i];
66670c78f05SPeter Brune     oscat   = nasm->oscatter[i];
667f10b3e88SPatrick Farrell     oscat_copy = nasm->oscatter_copy[i];
668f10b3e88SPatrick Farrell     gscat   = nasm->gscatter[i];
669f10b3e88SPatrick Farrell     ierr = VecScatterBegin(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
67070c78f05SPeter Brune     ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
67170c78f05SPeter Brune     if (B) {
67270c78f05SPeter Brune       /* scatter the RHS to the local RHS */
67370c78f05SPeter Brune       Bl   = nasm->b[i];
674f10b3e88SPatrick Farrell       ierr = VecScatterBegin(oscat_copy,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
67570c78f05SPeter Brune     }
67670c78f05SPeter Brune   }
677b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
678b20c023fSPeter Brune 
679b20c023fSPeter Brune 
680b20c023fSPeter Brune   if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
68170c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
68270c78f05SPeter Brune     Xl    = nasm->x[i];
68370c78f05SPeter Brune     Xlloc = nasm->xl[i];
68470c78f05SPeter Brune     Yl    = nasm->y[i];
685258e1594SPeter Brune     subsnes = nasm->subsnes[i];
686258e1594SPeter Brune     ierr    = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
687111ade9eSPeter Brune     iscat   = nasm->iscatter[i];
688111ade9eSPeter Brune     oscat   = nasm->oscatter[i];
689f10b3e88SPatrick Farrell     oscat_copy = nasm->oscatter_copy[i];
690111ade9eSPeter Brune     gscat   = nasm->gscatter[i];
691f10b3e88SPatrick Farrell     ierr = VecScatterEnd(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
692ed3c11a9SPeter Brune     ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
69324b7f281SPeter Brune     if (B) {
69424b7f281SPeter Brune       Bl   = nasm->b[i];
695f10b3e88SPatrick Farrell       ierr = VecScatterEnd(oscat_copy,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
696ed3c11a9SPeter Brune     } else Bl = NULL;
697f10b3e88SPatrick Farrell 
698ed3c11a9SPeter Brune     ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
699111ade9eSPeter Brune     ierr = VecCopy(Xl,Yl);CHKERRQ(ierr);
700ed3c11a9SPeter Brune     ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr);
701ed3c11a9SPeter Brune     ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr);
702f10b3e88SPatrick Farrell     ierr = VecScale(Yl, nasm->damping);CHKERRQ(ierr);
703e0331734SPeter Brune     if (type == PC_ASM_BASIC) {
704111ade9eSPeter Brune       ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
705e0331734SPeter Brune     } else if (type == PC_ASM_RESTRICT) {
706111ade9eSPeter Brune       ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
707ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
708eaedb033SPeter Brune   }
709ed3c11a9SPeter Brune   if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
710ed3c11a9SPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
71170c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
71270c78f05SPeter Brune     Yl    = nasm->y[i];
71370c78f05SPeter Brune     iscat   = nasm->iscatter[i];
71470c78f05SPeter Brune     oscat   = nasm->oscatter[i];
715e0331734SPeter Brune     if (type == PC_ASM_BASIC) {
71670c78f05SPeter Brune       ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
717e0331734SPeter Brune     } else if (type == PC_ASM_RESTRICT) {
71870c78f05SPeter Brune       ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
719ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
72070c78f05SPeter Brune   }
721f10b3e88SPatrick Farrell   if (nasm->weight_set) {
722f10b3e88SPatrick Farrell     ierr = VecPointwiseMult(Y,Y,nasm->weight);CHKERRQ(ierr);
723f10b3e88SPatrick Farrell   }
724b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
7255dfa0f3bSBarry Smith   ierr = SNESNASMGetDamping(snes,&dmp);CHKERRQ(ierr);
726610116beSPeter Brune   ierr = VecAXPY(X,dmp,Y);CHKERRQ(ierr);
727eaedb033SPeter Brune   PetscFunctionReturn(0);
728eaedb033SPeter Brune }
729eaedb033SPeter Brune 
73040244768SBarry Smith static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal)
731d728fb7dSPeter Brune {
732602bec5dSPeter Brune   Vec            X = Xfinal;
733d728fb7dSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
734d728fb7dSPeter Brune   SNES           subsnes;
735ca44f815SPeter Brune   PetscInt       i,lag = 1;
736d728fb7dSPeter Brune   PetscErrorCode ierr;
737e59f0a30SPeter Brune   Vec            Xlloc,Xl,Fl,F;
738d728fb7dSPeter Brune   VecScatter     oscat,gscat;
739d728fb7dSPeter Brune   DM             dm,subdm;
740d1e9a80fSBarry Smith 
741d728fb7dSPeter Brune   PetscFunctionBegin;
742602bec5dSPeter Brune   if (nasm->fjtype == 2) X = nasm->xinit;
743e59f0a30SPeter Brune   F = snes->vec_func;
744365a6726SPeter Brune   if (snes->normschedule == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);}
745d1e9a80fSBarry Smith   ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
746d728fb7dSPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
747d728fb7dSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
748602bec5dSPeter Brune   if (nasm->fjtype != 1) {
749d728fb7dSPeter Brune     for (i=0; i<nasm->n; i++) {
750d728fb7dSPeter Brune       Xlloc = nasm->xl[i];
751d728fb7dSPeter Brune       gscat = nasm->gscatter[i];
752d728fb7dSPeter Brune       ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
753602bec5dSPeter Brune     }
754d728fb7dSPeter Brune   }
755d728fb7dSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
756d728fb7dSPeter Brune   for (i=0; i<nasm->n; i++) {
757e59f0a30SPeter Brune     Fl      = nasm->subsnes[i]->vec_func;
758d728fb7dSPeter Brune     Xl      = nasm->x[i];
759d728fb7dSPeter Brune     Xlloc   = nasm->xl[i];
760d728fb7dSPeter Brune     subsnes = nasm->subsnes[i];
761d728fb7dSPeter Brune     oscat   = nasm->oscatter[i];
762d728fb7dSPeter Brune     gscat   = nasm->gscatter[i];
763602bec5dSPeter Brune     if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);}
764ed3c11a9SPeter Brune     ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
765d728fb7dSPeter Brune     ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
766602bec5dSPeter Brune     if (nasm->fjtype != 1) {
767d728fb7dSPeter Brune       ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
768d728fb7dSPeter Brune       ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
769602bec5dSPeter Brune     }
770ca44f815SPeter Brune     if (subsnes->lagjacobian == -1)    subsnes->lagjacobian = -2;
771ca44f815SPeter Brune     else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian;
772602bec5dSPeter Brune     ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr);
773d1e9a80fSBarry Smith     ierr = SNESComputeJacobian(subsnes,Xl,subsnes->jacobian,subsnes->jacobian_pre);CHKERRQ(ierr);
774ca44f815SPeter Brune     if (lag > 1) subsnes->lagjacobian = lag;
775d728fb7dSPeter Brune   }
776d728fb7dSPeter Brune   PetscFunctionReturn(0);
777d728fb7dSPeter Brune }
778d728fb7dSPeter Brune 
77940244768SBarry Smith static PetscErrorCode SNESSolve_NASM(SNES snes)
780eaedb033SPeter Brune {
781eaedb033SPeter Brune   Vec              F;
782eaedb033SPeter Brune   Vec              X;
783eaedb033SPeter Brune   Vec              B;
784111ade9eSPeter Brune   Vec              Y;
785eaedb033SPeter Brune   PetscInt         i;
786ed3c11a9SPeter Brune   PetscReal        fnorm = 0.0;
787eaedb033SPeter Brune   PetscErrorCode   ierr;
788365a6726SPeter Brune   SNESNormSchedule normschedule;
789d728fb7dSPeter Brune   SNES_NASM        *nasm = (SNES_NASM*)snes->data;
790eaedb033SPeter Brune 
791eaedb033SPeter Brune   PetscFunctionBegin;
792c579b300SPatrick Farrell 
7936c4ed002SBarry Smith   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
794c579b300SPatrick Farrell 
795fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
796eaedb033SPeter Brune   X = snes->vec_sol;
797111ade9eSPeter Brune   Y = snes->vec_sol_update;
798eaedb033SPeter Brune   F = snes->vec_func;
799eaedb033SPeter Brune   B = snes->vec_rhs;
800eaedb033SPeter Brune 
801e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
802eaedb033SPeter Brune   snes->iter   = 0;
803eaedb033SPeter Brune   snes->norm   = 0.;
804e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
805eaedb033SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
806365a6726SPeter Brune   ierr         = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr);
807365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
808eaedb033SPeter Brune     /* compute the initial function and preconditioned update delX */
809eaedb033SPeter Brune     if (!snes->vec_func_init_set) {
810eaedb033SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
8111aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
812eaedb033SPeter Brune 
813eaedb033SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
814422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
815e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
816eaedb033SPeter Brune     snes->iter = 0;
817eaedb033SPeter Brune     snes->norm = fnorm;
818e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
819a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
820eaedb033SPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
821eaedb033SPeter Brune 
822eaedb033SPeter Brune     /* test convergence */
823eaedb033SPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
824eaedb033SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
825eaedb033SPeter Brune   } else {
826e04113cfSBarry Smith     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
827a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
828eaedb033SPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
829eaedb033SPeter Brune   }
830eaedb033SPeter Brune 
831eaedb033SPeter Brune   /* Call general purpose update function */
832eaedb033SPeter Brune   if (snes->ops->update) {
833eaedb033SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
834eaedb033SPeter Brune   }
835602bec5dSPeter Brune   /* copy the initial solution over for later */
836602bec5dSPeter Brune   if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);}
837eaedb033SPeter Brune 
838eaedb033SPeter Brune   for (i=0; i < snes->max_its; i++) {
839111ade9eSPeter Brune     ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr);
840365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
841eaedb033SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
842eaedb033SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
843422a814eSBarry Smith       SNESCheckFunctionNorm(snes,fnorm);
844eaedb033SPeter Brune     }
845eaedb033SPeter Brune     /* Monitor convergence */
846e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
847eaedb033SPeter Brune     snes->iter = i+1;
848eaedb033SPeter Brune     snes->norm = fnorm;
849e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
850a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
851eaedb033SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
852eaedb033SPeter Brune     /* Test for convergence */
853365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
854d728fb7dSPeter Brune     if (snes->reason) break;
855eaedb033SPeter Brune     /* Call general purpose update function */
856e59f0a30SPeter Brune     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
857eaedb033SPeter Brune   }
858d728fb7dSPeter Brune   if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);}
859365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS) {
860eaedb033SPeter Brune     if (i == snes->max_its) {
861eaedb033SPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
862eaedb033SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
863eaedb033SPeter Brune     }
8641aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */
865eaedb033SPeter Brune   PetscFunctionReturn(0);
866eaedb033SPeter Brune }
867eaedb033SPeter Brune 
868eaedb033SPeter Brune /*MC
869eaedb033SPeter Brune   SNESNASM - Nonlinear Additive Schwartz
870eaedb033SPeter Brune 
87170c78603SPeter Brune    Options Database:
87270c78603SPeter Brune +  -snes_nasm_log - enable logging events for the communication and solve stages
87370c78603SPeter Brune .  -snes_nasm_type <basic,restrict> - type of subdomain update used
8745dfa0f3bSBarry Smith .  -snes_asm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
87570c78603SPeter Brune .  -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate
876150d49b7SHong Zhang .  -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the jacobian is calculated at
87770c78603SPeter Brune .  -sub_snes_ - options prefix of the subdomain nonlinear solves
87870c78603SPeter Brune .  -sub_ksp_ - options prefix of the subdomain Krylov solver
87970c78603SPeter Brune -  -sub_pc_ - options prefix of the subdomain preconditioner
88070c78603SPeter Brune 
881eaedb033SPeter Brune    Level: advanced
882eaedb033SPeter Brune 
8834f02bc6aSBarry Smith    References:
88496a0c994SBarry Smith .  1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
8854f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
8864f02bc6aSBarry Smith 
8875dfa0f3bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types), SNESNASMSetType(), SNESNASMGetType(), SNESNASMSetSubdomains(), SNESNASMGetSubdomains(), SNESNASMGetSubdomainVecs(), SNESNASMSetComputeFinalJacobian(), SNESNASMSetDamping(), SNESNASMGetDamping()
888eaedb033SPeter Brune M*/
889eaedb033SPeter Brune 
8908cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes)
891eaedb033SPeter Brune {
892eaedb033SPeter Brune   SNES_NASM      *nasm;
893eaedb033SPeter Brune   PetscErrorCode ierr;
894eaedb033SPeter Brune 
895eaedb033SPeter Brune   PetscFunctionBegin;
896b00a9115SJed Brown   ierr       = PetscNewLog(snes,&nasm);CHKERRQ(ierr);
897eaedb033SPeter Brune   snes->data = (void*)nasm;
898eaedb033SPeter Brune 
899eaedb033SPeter Brune   nasm->n        = PETSC_DECIDE;
900eaedb033SPeter Brune   nasm->subsnes  = 0;
901eaedb033SPeter Brune   nasm->x        = 0;
902111ade9eSPeter Brune   nasm->xl       = 0;
903111ade9eSPeter Brune   nasm->y        = 0;
904eaedb033SPeter Brune   nasm->b        = 0;
905111ade9eSPeter Brune   nasm->oscatter = 0;
906f10b3e88SPatrick Farrell   nasm->oscatter_copy = 0;
907111ade9eSPeter Brune   nasm->iscatter = 0;
908111ade9eSPeter Brune   nasm->gscatter = 0;
909610116beSPeter Brune   nasm->damping  = 1.;
910111ade9eSPeter Brune 
911111ade9eSPeter Brune   nasm->type              = PC_ASM_BASIC;
912d728fb7dSPeter Brune   nasm->finaljacobian     = PETSC_FALSE;
913a4f17876SPeter Brune   nasm->same_local_solves = PETSC_TRUE;
914f10b3e88SPatrick Farrell   nasm->weight_set        = PETSC_FALSE;
915eaedb033SPeter Brune 
916eaedb033SPeter Brune   snes->ops->destroy        = SNESDestroy_NASM;
917eaedb033SPeter Brune   snes->ops->setup          = SNESSetUp_NASM;
918eaedb033SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
919eaedb033SPeter Brune   snes->ops->view           = SNESView_NASM;
920eaedb033SPeter Brune   snes->ops->solve          = SNESSolve_NASM;
921eaedb033SPeter Brune   snes->ops->reset          = SNESReset_NASM;
922eaedb033SPeter Brune 
923eaedb033SPeter Brune   snes->usesksp = PETSC_FALSE;
924efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
925eaedb033SPeter Brune 
9264fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
9274fc747eaSLawrence Mitchell 
928602bec5dSPeter Brune   nasm->fjtype              = 0;
929602bec5dSPeter Brune   nasm->xinit               = NULL;
9300298fd71SBarry Smith   nasm->eventrestrictinterp = 0;
9310298fd71SBarry Smith   nasm->eventsubsolve       = 0;
932b20c023fSPeter Brune 
933eaedb033SPeter Brune   if (!snes->tolerancesset) {
934eaedb033SPeter Brune     snes->max_its   = 10000;
935eaedb033SPeter Brune     snes->max_funcs = 10000;
936eaedb033SPeter Brune   }
937eaedb033SPeter Brune 
938e0331734SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetType_C",SNESNASMSetType_NASM);CHKERRQ(ierr);
939e0331734SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetType_C",SNESNASMGetType_NASM);CHKERRQ(ierr);
940bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr);
941bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr);
94290bcee39SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",SNESNASMSetDamping_NASM);CHKERRQ(ierr);
94390bcee39SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",SNESNASMGetDamping_NASM);CHKERRQ(ierr);
944bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr);
945bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr);
946eaedb033SPeter Brune   PetscFunctionReturn(0);
947eaedb033SPeter Brune }
94899e0435eSBarry Smith 
949448b6425SPatrick Farrell /*@
950448b6425SPatrick Farrell    SNESNASMGetSNES - Gets a subsolver
951448b6425SPatrick Farrell 
952448b6425SPatrick Farrell    Not collective
953448b6425SPatrick Farrell 
954448b6425SPatrick Farrell    Input Parameters:
955448b6425SPatrick Farrell +  snes - the SNES context
956448b6425SPatrick Farrell -  i - the number of the subsnes to get
957448b6425SPatrick Farrell 
958448b6425SPatrick Farrell    Output Parameters:
959448b6425SPatrick Farrell .  subsnes - the subsolver context
960448b6425SPatrick Farrell 
961448b6425SPatrick Farrell    Level: intermediate
962448b6425SPatrick Farrell 
963448b6425SPatrick Farrell .keywords: SNES, NASM
964448b6425SPatrick Farrell 
965448b6425SPatrick Farrell .seealso: SNESNASM, SNESNASMGetNumber()
966448b6425SPatrick Farrell @*/
967448b6425SPatrick Farrell PetscErrorCode SNESNASMGetSNES(SNES snes,PetscInt i,SNES *subsnes)
968448b6425SPatrick Farrell {
969448b6425SPatrick Farrell   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
970448b6425SPatrick Farrell 
971448b6425SPatrick Farrell   PetscFunctionBegin;
972448b6425SPatrick Farrell   if (i < 0 || i >= nasm->n) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"No such subsolver");
973448b6425SPatrick Farrell   *subsnes = nasm->subsnes[i];
974448b6425SPatrick Farrell   PetscFunctionReturn(0);
975448b6425SPatrick Farrell }
976448b6425SPatrick Farrell 
977448b6425SPatrick Farrell /*@
978448b6425SPatrick Farrell    SNESNASMGetNumber - Gets number of subsolvers
979448b6425SPatrick Farrell 
980448b6425SPatrick Farrell    Not collective
981448b6425SPatrick Farrell 
982448b6425SPatrick Farrell    Input Parameters:
983448b6425SPatrick Farrell .  snes - the SNES context
984448b6425SPatrick Farrell 
985448b6425SPatrick Farrell    Output Parameters:
986448b6425SPatrick Farrell .  n - the number of subsolvers
987448b6425SPatrick Farrell 
988448b6425SPatrick Farrell    Level: intermediate
989448b6425SPatrick Farrell 
990448b6425SPatrick Farrell .keywords: SNES, NASM
991448b6425SPatrick Farrell 
992448b6425SPatrick Farrell .seealso: SNESNASM, SNESNASMGetSNES()
993448b6425SPatrick Farrell @*/
994448b6425SPatrick Farrell PetscErrorCode SNESNASMGetNumber(SNES snes,PetscInt *n)
995448b6425SPatrick Farrell {
996448b6425SPatrick Farrell   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
997448b6425SPatrick Farrell 
998448b6425SPatrick Farrell   PetscFunctionBegin;
999448b6425SPatrick Farrell   *n = nasm->n;
1000448b6425SPatrick Farrell   PetscFunctionReturn(0);
1001448b6425SPatrick Farrell }
1002448b6425SPatrick Farrell 
1003f10b3e88SPatrick Farrell /*@
1004f10b3e88SPatrick Farrell    SNESNASMSetWeight - Sets weight to use when adding overlapping updates
1005f10b3e88SPatrick Farrell 
1006f10b3e88SPatrick Farrell    Collective
1007f10b3e88SPatrick Farrell 
1008f10b3e88SPatrick Farrell    Input Parameters:
1009f10b3e88SPatrick Farrell +  snes - the SNES context
1010f10b3e88SPatrick Farrell -  weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in)
1011f10b3e88SPatrick Farrell 
1012f10b3e88SPatrick Farrell    Level: intermediate
1013f10b3e88SPatrick Farrell 
1014f10b3e88SPatrick Farrell .keywords: SNES, NASM
1015f10b3e88SPatrick Farrell 
1016f10b3e88SPatrick Farrell .seealso: SNESNASM
1017f10b3e88SPatrick Farrell @*/
1018f10b3e88SPatrick Farrell PetscErrorCode SNESNASMSetWeight(SNES snes,Vec weight)
1019f10b3e88SPatrick Farrell {
1020f10b3e88SPatrick Farrell   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
1021f10b3e88SPatrick Farrell   PetscErrorCode ierr;
1022f10b3e88SPatrick Farrell 
1023f10b3e88SPatrick Farrell   PetscFunctionBegin;
1024f10b3e88SPatrick Farrell 
1025f10b3e88SPatrick Farrell   ierr = VecDestroy(&nasm->weight);CHKERRQ(ierr);
1026f10b3e88SPatrick Farrell   nasm->weight_set = PETSC_TRUE;
1027f10b3e88SPatrick Farrell   nasm->weight     = weight;
1028f10b3e88SPatrick Farrell   ierr = PetscObjectReference((PetscObject)nasm->weight);CHKERRQ(ierr);
1029f10b3e88SPatrick Farrell 
1030f10b3e88SPatrick Farrell   PetscFunctionReturn(0);
1031f10b3e88SPatrick Farrell }
1032