xref: /petsc/src/snes/impls/nasm/nasm.c (revision 4fc747eaadbeca11629f314a99edccbc2ed7b3d3)
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 */
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 */
17610116beSPeter Brune   PetscReal  damping;             /* damping parameter for updates from the blocks */
18a4f17876SPeter Brune   PetscBool  same_local_solves;   /* flag to determine if the solvers have been individually modified */
19b20c023fSPeter Brune 
20b20c023fSPeter Brune   /* logging events */
21b20c023fSPeter Brune   PetscLogEvent eventrestrictinterp;
22b20c023fSPeter Brune   PetscLogEvent eventsubsolve;
23602bec5dSPeter Brune 
24602bec5dSPeter Brune   PetscInt      fjtype;            /* type of computed jacobian */
25602bec5dSPeter Brune   Vec           xinit;             /* initial solution in case the final jacobian type is computed as first */
26eaedb033SPeter Brune } SNES_NASM;
27eaedb033SPeter Brune 
28b20c023fSPeter Brune const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0};
29602bec5dSPeter Brune const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"};
30b20c023fSPeter Brune 
31eaedb033SPeter Brune #undef __FUNCT__
32eaedb033SPeter Brune #define __FUNCT__ "SNESReset_NASM"
3340244768SBarry Smith static PetscErrorCode SNESReset_NASM(SNES snes)
34eaedb033SPeter Brune {
35eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
36eaedb033SPeter Brune   PetscErrorCode ierr;
37eaedb033SPeter Brune   PetscInt       i;
386e111a19SKarl Rupp 
39eaedb033SPeter Brune   PetscFunctionBegin;
40eaedb033SPeter Brune   for (i=0; i<nasm->n; i++) {
41111ade9eSPeter Brune     if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); }
42f5f7c1b9SKarl Rupp     if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); }
43111ade9eSPeter Brune     if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); }
44bc8c1f72SJose Roman     if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); }
45eaedb033SPeter Brune 
46bc8c1f72SJose Roman     if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); }
47111ade9eSPeter Brune     if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); }
48111ade9eSPeter Brune     if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); }
49111ade9eSPeter Brune     if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); }
50eaedb033SPeter Brune   }
51111ade9eSPeter Brune 
52111ade9eSPeter Brune   if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);}
53111ade9eSPeter Brune   if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);}
54111ade9eSPeter Brune   if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);}
55111ade9eSPeter Brune   if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);}
56111ade9eSPeter Brune 
57602bec5dSPeter Brune   if (nasm->xinit) {ierr = VecDestroy(&nasm->xinit);CHKERRQ(ierr);}
58602bec5dSPeter Brune 
59111ade9eSPeter Brune   if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);}
60111ade9eSPeter Brune   if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);}
61111ade9eSPeter Brune   if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);}
62111ade9eSPeter Brune   if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);}
63b20c023fSPeter Brune 
64b20c023fSPeter Brune   nasm->eventrestrictinterp = 0;
65b20c023fSPeter Brune   nasm->eventsubsolve = 0;
66eaedb033SPeter Brune   PetscFunctionReturn(0);
67eaedb033SPeter Brune }
68eaedb033SPeter Brune 
69eaedb033SPeter Brune #undef __FUNCT__
70eaedb033SPeter Brune #define __FUNCT__ "SNESDestroy_NASM"
7140244768SBarry Smith static PetscErrorCode SNESDestroy_NASM(SNES snes)
72eaedb033SPeter Brune {
73eaedb033SPeter Brune   PetscErrorCode ierr;
746e111a19SKarl Rupp 
75eaedb033SPeter Brune   PetscFunctionBegin;
76eaedb033SPeter Brune   ierr = SNESReset_NASM(snes);CHKERRQ(ierr);
7722d28d08SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
78eaedb033SPeter Brune   PetscFunctionReturn(0);
79eaedb033SPeter Brune }
80eaedb033SPeter Brune 
81eaedb033SPeter Brune #undef __FUNCT__
82111ade9eSPeter Brune #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private"
8340244768SBarry Smith static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx)
840adebc6cSBarry Smith {
85111ade9eSPeter Brune   PetscErrorCode ierr;
86111ade9eSPeter Brune   Vec            bcs = (Vec)ctx;
876e111a19SKarl Rupp 
88111ade9eSPeter Brune   PetscFunctionBegin;
89111ade9eSPeter Brune   ierr = VecCopy(bcs,l);CHKERRQ(ierr);
90111ade9eSPeter Brune   PetscFunctionReturn(0);
91111ade9eSPeter Brune }
92111ade9eSPeter Brune 
93111ade9eSPeter Brune #undef __FUNCT__
94eaedb033SPeter Brune #define __FUNCT__ "SNESSetUp_NASM"
9540244768SBarry Smith static PetscErrorCode SNESSetUp_NASM(SNES snes)
96eaedb033SPeter Brune {
97eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
98eaedb033SPeter Brune   PetscErrorCode ierr;
9976857b2aSPeter Brune   DM             dm,subdm;
100111ade9eSPeter Brune   DM             *subdms;
101111ade9eSPeter Brune   PetscInt       i;
102eaedb033SPeter Brune   const char     *optionsprefix;
103111ade9eSPeter Brune   Vec            F;
104ed3c11a9SPeter Brune   PetscMPIInt    size;
105ed3c11a9SPeter Brune   KSP            ksp;
106ed3c11a9SPeter Brune   PC             pc;
107eaedb033SPeter Brune 
108eaedb033SPeter Brune   PetscFunctionBegin;
109eaedb033SPeter Brune   if (!nasm->subsnes) {
110eaedb033SPeter Brune     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1110a696f66SPeter Brune     if (dm) {
112eaedb033SPeter Brune       nasm->usesdm = PETSC_TRUE;
1130298fd71SBarry Smith       ierr         = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr);
114ce94432eSBarry Smith       if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
115111ade9eSPeter Brune       ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr);
116eaedb033SPeter Brune 
117eaedb033SPeter Brune       ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
118785e854fSJed Brown       ierr = PetscMalloc1(nasm->n,&nasm->subsnes);CHKERRQ(ierr);
119111ade9eSPeter Brune       for (i=0; i<nasm->n; i++) {
120cdb298fcSPeter Brune         ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr);
121cdb298fcSPeter Brune         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr);
122cdb298fcSPeter Brune         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr);
123cdb298fcSPeter Brune         ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr);
124ed3c11a9SPeter Brune         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRQ(ierr);
125ed3c11a9SPeter Brune         if (size == 1) {
126ed3c11a9SPeter Brune           ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr);
127ed3c11a9SPeter Brune           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
128ed3c11a9SPeter Brune           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
129ed3c11a9SPeter Brune           ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
130ed3c11a9SPeter Brune         }
131cdb298fcSPeter Brune         ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr);
132111ade9eSPeter Brune         ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr);
133111ade9eSPeter Brune       }
134111ade9eSPeter Brune       ierr = PetscFree(subdms);CHKERRQ(ierr);
135ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!");
136ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!");
137111ade9eSPeter Brune   /* allocate the global vectors */
138e245e0aaSPeter Brune   if (!nasm->x) {
1391795a4d1SJed Brown     ierr = PetscCalloc1(nasm->n,&nasm->x);CHKERRQ(ierr);
140e245e0aaSPeter Brune   }
141e245e0aaSPeter Brune   if (!nasm->xl) {
1421795a4d1SJed Brown     ierr = PetscCalloc1(nasm->n,&nasm->xl);CHKERRQ(ierr);
143e245e0aaSPeter Brune   }
144e245e0aaSPeter Brune   if (!nasm->y) {
1451795a4d1SJed Brown     ierr = PetscCalloc1(nasm->n,&nasm->y);CHKERRQ(ierr);
146e245e0aaSPeter Brune   }
147e245e0aaSPeter Brune   if (!nasm->b) {
1481795a4d1SJed Brown     ierr = PetscCalloc1(nasm->n,&nasm->b);CHKERRQ(ierr);
149e245e0aaSPeter Brune   }
150111ade9eSPeter Brune 
151111ade9eSPeter Brune   for (i=0; i<nasm->n; i++) {
1520298fd71SBarry Smith     ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr);
15376857b2aSPeter Brune     if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);}
15476857b2aSPeter Brune     if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);}
15576857b2aSPeter Brune     if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);}
15676857b2aSPeter Brune     if (!nasm->xl[i]) {
157111ade9eSPeter Brune       ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr);
158111ade9eSPeter Brune       ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr);
1590298fd71SBarry Smith       ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr);
160111ade9eSPeter Brune     }
16161ba4676SBarry Smith   }
162602bec5dSPeter Brune   if (nasm->finaljacobian) {
163602bec5dSPeter Brune     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
164602bec5dSPeter Brune     if (nasm->fjtype == 2) {
165602bec5dSPeter Brune       ierr = VecDuplicate(snes->vec_sol,&nasm->xinit);CHKERRQ(ierr);
166602bec5dSPeter Brune     }
167602bec5dSPeter Brune     for (i=0; i<nasm->n;i++) {
168602bec5dSPeter Brune       ierr = SNESSetUpMatrices(nasm->subsnes[i]);CHKERRQ(ierr);
169602bec5dSPeter Brune     }
170602bec5dSPeter Brune   }
171eaedb033SPeter Brune   PetscFunctionReturn(0);
172eaedb033SPeter Brune }
173eaedb033SPeter Brune 
174eaedb033SPeter Brune #undef __FUNCT__
175eaedb033SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NASM"
17640244768SBarry Smith static PetscErrorCode SNESSetFromOptions_NASM(PetscOptionItems *PetscOptionsObject,SNES snes)
177eaedb033SPeter Brune {
178eaedb033SPeter Brune   PetscErrorCode    ierr;
179111ade9eSPeter Brune   PCASMType         asmtype;
180a4f17876SPeter Brune   PetscBool         flg,monflg,subviewflg;
181111ade9eSPeter Brune   SNES_NASM         *nasm = (SNES_NASM*)snes->data;
1826e111a19SKarl Rupp 
183eaedb033SPeter Brune   PetscFunctionBegin;
184e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Nonlinear Additive Schwartz options");CHKERRQ(ierr);
185111ade9eSPeter Brune   ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
186e0331734SPeter Brune   if (flg) {ierr = SNESNASMSetType(snes,asmtype);CHKERRQ(ierr);}
187b20c023fSPeter Brune   flg    = PETSC_FALSE;
188b20c023fSPeter Brune   monflg = PETSC_TRUE;
1895dfa0f3bSBarry 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);
190610116beSPeter Brune   if (flg) {ierr = SNESNASMSetDamping(snes,nasm->damping);CHKERRQ(ierr);}
191a4f17876SPeter Brune   subviewflg = PETSC_FALSE;
192a4f17876SPeter Brune   ierr   = PetscOptionsBool("-snes_nasm_sub_view","Print detailed information for every processor when using -snes_view","",subviewflg,&subviewflg,&flg);CHKERRQ(ierr);
193a4f17876SPeter Brune   if (flg) {
194a4f17876SPeter Brune     nasm->same_local_solves = PETSC_FALSE;
195a4f17876SPeter Brune     if (!subviewflg) {
196a4f17876SPeter Brune       nasm->same_local_solves = PETSC_TRUE;
197a4f17876SPeter Brune     }
198a4f17876SPeter Brune   }
199d728fb7dSPeter Brune   ierr   = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr);
200602bec5dSPeter Brune   ierr   = PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL);CHKERRQ(ierr);
201a4f17876SPeter Brune   ierr   = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr);
202b20c023fSPeter Brune   if (flg) {
203b20c023fSPeter Brune     ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr);
204b20c023fSPeter Brune     ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr);
205b20c023fSPeter Brune   }
206eaedb033SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
207eaedb033SPeter Brune   PetscFunctionReturn(0);
208eaedb033SPeter Brune }
209eaedb033SPeter Brune 
210eaedb033SPeter Brune #undef __FUNCT__
211eaedb033SPeter Brune #define __FUNCT__ "SNESView_NASM"
21240244768SBarry Smith static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
213eaedb033SPeter Brune {
214b20c023fSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
215b20c023fSPeter Brune   PetscErrorCode ierr;
216a4f17876SPeter Brune   PetscMPIInt    rank,size;
217dd2fa690SBarry Smith   PetscInt       i,N,bsz;
218b20c023fSPeter Brune   PetscBool      iascii,isstring;
219b20c023fSPeter Brune   PetscViewer    sviewer;
220ce94432eSBarry Smith   MPI_Comm       comm;
221b20c023fSPeter Brune 
222b20c023fSPeter Brune   PetscFunctionBegin;
223ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
224b20c023fSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
225b20c023fSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
226b20c023fSPeter Brune   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
227a4f17876SPeter Brune   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
228b2566f29SBarry Smith   ierr = MPIU_Allreduce(&nasm->n,&N,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
229b20c023fSPeter Brune   if (iascii) {
230b20c023fSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "  Nonlinear Additive Schwarz: total subdomain blocks = %D\n",N);CHKERRQ(ierr);
231a4f17876SPeter Brune     if (nasm->same_local_solves) {
232a4f17876SPeter Brune       if (nasm->subsnes) {
233a4f17876SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  Local solve is the same for all blocks:\n");CHKERRQ(ierr);
234a4f17876SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
2353f08860eSBarry Smith         ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
236a4f17876SPeter Brune         if (!rank) {
237a4f17876SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
238a4f17876SPeter Brune           ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);
239a4f17876SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
240a4f17876SPeter Brune         }
2413f08860eSBarry Smith         ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
242a4f17876SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
243a4f17876SPeter Brune       }
244a4f17876SPeter Brune     } else {
245a4f17876SPeter Brune       /* print the solver on each block */
2461575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
247b20c023fSPeter Brune       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr);
248b20c023fSPeter Brune       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
2491575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
250a4f17876SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following SNES objects:\n");CHKERRQ(ierr);
251b20c023fSPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
252a4f17876SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
2533f08860eSBarry Smith       ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
254a4f17876SPeter Brune       for (i=0; i<nasm->n; i++) {
255a4f17876SPeter Brune         ierr = VecGetLocalSize(nasm->x[i],&bsz);CHKERRQ(ierr);
256a4f17876SPeter Brune         ierr = PetscViewerASCIIPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr);
257b20c023fSPeter Brune         ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr);
258a4f17876SPeter Brune         ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
259b20c023fSPeter Brune       }
2603f08860eSBarry Smith       ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
261a4f17876SPeter Brune       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
262b20c023fSPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
263a4f17876SPeter Brune     }
264b20c023fSPeter Brune   } else if (isstring) {
265b20c023fSPeter Brune     ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr);
2663f08860eSBarry Smith     ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
267b20c023fSPeter Brune     if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);}
2683f08860eSBarry Smith     ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
269b20c023fSPeter Brune   }
270eaedb033SPeter Brune   PetscFunctionReturn(0);
271eaedb033SPeter Brune }
272eaedb033SPeter Brune 
273eaedb033SPeter Brune #undef __FUNCT__
274e0331734SPeter Brune #define __FUNCT__ "SNESNASMSetType"
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 
301e0331734SPeter Brune #undef __FUNCT__
302e0331734SPeter Brune #define __FUNCT__ "SNESNASMSetType_NASM"
30340244768SBarry Smith static PetscErrorCode SNESNASMSetType_NASM(SNES snes,PCASMType type)
304e0331734SPeter Brune {
305e0331734SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
306e0331734SPeter Brune 
307e0331734SPeter Brune   PetscFunctionBegin;
308e0331734SPeter 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");
309e0331734SPeter Brune   nasm->type = type;
310e0331734SPeter Brune   PetscFunctionReturn(0);
311e0331734SPeter Brune }
312e0331734SPeter Brune 
313e0331734SPeter Brune #undef __FUNCT__
314e0331734SPeter Brune #define __FUNCT__ "SNESNASMGetType"
315e0331734SPeter Brune /*@
316e0331734SPeter Brune    SNESNASMGetType - Get the type of subdomain update used
317e0331734SPeter Brune 
318e0331734SPeter Brune    Logically Collective on SNES
319e0331734SPeter Brune 
320e0331734SPeter Brune    Input Parameters:
321e0331734SPeter Brune .  SNES - the SNES context
322e0331734SPeter Brune 
323e0331734SPeter Brune    Output Parameters:
324e0331734SPeter Brune .  type - the type of update
325e0331734SPeter Brune 
326e0331734SPeter Brune    Level: intermediate
327e0331734SPeter Brune 
328e0331734SPeter Brune .keywords: SNES, NASM
329e0331734SPeter Brune 
330e0331734SPeter Brune .seealso: SNESNASM, SNESNASMSetType(), PCASMGetType()
331e0331734SPeter Brune @*/
332e0331734SPeter Brune PetscErrorCode SNESNASMGetType(SNES snes,PCASMType *type)
333e0331734SPeter Brune {
334e0331734SPeter Brune   PetscErrorCode ierr;
335e0331734SPeter Brune 
336e0331734SPeter Brune   PetscFunctionBegin;
3372a808120SBarry Smith   ierr = PetscUseMethod(snes,"SNESNASMGetType_C",(SNES,PCASMType*),(snes,type));CHKERRQ(ierr);
338e0331734SPeter Brune   PetscFunctionReturn(0);
339e0331734SPeter Brune }
340e0331734SPeter Brune 
341e0331734SPeter Brune #undef __FUNCT__
342e0331734SPeter Brune #define __FUNCT__ "SNESNASMGetType_NASM"
34340244768SBarry Smith static PetscErrorCode SNESNASMGetType_NASM(SNES snes,PCASMType *type)
344e0331734SPeter Brune {
345e0331734SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
346e0331734SPeter Brune 
347e0331734SPeter Brune   PetscFunctionBegin;
348e0331734SPeter Brune   *type = nasm->type;
349e0331734SPeter Brune   PetscFunctionReturn(0);
350e0331734SPeter Brune }
351e0331734SPeter Brune 
352e0331734SPeter Brune #undef __FUNCT__
353eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains"
35476857b2aSPeter Brune /*@
35576857b2aSPeter Brune    SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems.
35676857b2aSPeter Brune 
35776857b2aSPeter Brune    Not Collective
35876857b2aSPeter Brune 
35976857b2aSPeter Brune    Input Parameters:
36076857b2aSPeter Brune +  SNES - the SNES context
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, SNESNASMGetSubdomains()
37276857b2aSPeter Brune @*/
373a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
374a6dfd86eSKarl Rupp {
375eaedb033SPeter Brune   PetscErrorCode ierr;
376111ade9eSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*);
3776e111a19SKarl Rupp 
378eaedb033SPeter Brune   PetscFunctionBegin;
3790005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f);CHKERRQ(ierr);
3804b838e8fSPeter Brune   if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);}
381eaedb033SPeter Brune   PetscFunctionReturn(0);
382eaedb033SPeter Brune }
383eaedb033SPeter Brune 
384eaedb033SPeter Brune #undef __FUNCT__
385eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains_NASM"
38640244768SBarry Smith static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
387a6dfd86eSKarl Rupp {
388eaedb033SPeter Brune   PetscInt       i;
389eaedb033SPeter Brune   PetscErrorCode ierr;
390eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
3916e111a19SKarl Rupp 
392eaedb033SPeter Brune   PetscFunctionBegin;
393ce94432eSBarry Smith   if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
394eaedb033SPeter Brune 
395111ade9eSPeter Brune   /* tear down the previously set things */
396111ade9eSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
397111ade9eSPeter Brune 
398eaedb033SPeter Brune   nasm->n = n;
399111ade9eSPeter Brune   if (oscatter) {
400111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);}
401eaedb033SPeter Brune   }
402111ade9eSPeter Brune   if (iscatter) {
403111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);}
404eaedb033SPeter Brune   }
405111ade9eSPeter Brune   if (gscatter) {
406111ade9eSPeter Brune     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);}
407111ade9eSPeter Brune   }
408111ade9eSPeter Brune   if (oscatter) {
409785e854fSJed Brown     ierr = PetscMalloc1(n,&nasm->oscatter);CHKERRQ(ierr);
410eaedb033SPeter Brune     for (i=0; i<n; i++) {
411111ade9eSPeter Brune       nasm->oscatter[i] = oscatter[i];
412eaedb033SPeter Brune     }
413111ade9eSPeter Brune   }
414111ade9eSPeter Brune   if (iscatter) {
415785e854fSJed Brown     ierr = PetscMalloc1(n,&nasm->iscatter);CHKERRQ(ierr);
416eaedb033SPeter Brune     for (i=0; i<n; i++) {
417111ade9eSPeter Brune       nasm->iscatter[i] = iscatter[i];
418eaedb033SPeter Brune     }
419eaedb033SPeter Brune   }
420111ade9eSPeter Brune   if (gscatter) {
421785e854fSJed Brown     ierr = PetscMalloc1(n,&nasm->gscatter);CHKERRQ(ierr);
422eaedb033SPeter Brune     for (i=0; i<n; i++) {
423111ade9eSPeter Brune       nasm->gscatter[i] = gscatter[i];
424eaedb033SPeter Brune     }
425eaedb033SPeter Brune   }
426111ade9eSPeter Brune 
427eaedb033SPeter Brune   if (subsnes) {
428785e854fSJed Brown     ierr = PetscMalloc1(n,&nasm->subsnes);CHKERRQ(ierr);
429eaedb033SPeter Brune     for (i=0; i<n; i++) {
430eaedb033SPeter Brune       nasm->subsnes[i] = subsnes[i];
431eaedb033SPeter Brune     }
432a4f17876SPeter Brune     nasm->same_local_solves = PETSC_FALSE;
433eaedb033SPeter Brune   }
434eaedb033SPeter Brune   PetscFunctionReturn(0);
435eaedb033SPeter Brune }
436eaedb033SPeter Brune 
437eaedb033SPeter Brune #undef __FUNCT__
43876857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains"
43976857b2aSPeter Brune /*@
44076857b2aSPeter Brune    SNESNASMGetSubdomains - Get the local subdomain context.
44176857b2aSPeter Brune 
44276857b2aSPeter Brune    Not Collective
44376857b2aSPeter Brune 
44476857b2aSPeter Brune    Input Parameters:
44576857b2aSPeter Brune .  SNES - the SNES context
44676857b2aSPeter Brune 
44776857b2aSPeter Brune    Output Parameters:
44876857b2aSPeter Brune +  n - the number of local subdomains
44976857b2aSPeter Brune .  subsnes - solvers defined on the local subdomains
45076857b2aSPeter Brune .  iscatter - scatters into the nonoverlapping portions of the local subdomains
45176857b2aSPeter Brune .  oscatter - scatters into the overlapping portions of the local subdomains
45276857b2aSPeter Brune -  gscatter - scatters into the (ghosted) local vector of the local subdomain
45376857b2aSPeter Brune 
45476857b2aSPeter Brune    Level: intermediate
45576857b2aSPeter Brune 
45676857b2aSPeter Brune .keywords: SNES, NASM
45776857b2aSPeter Brune 
45876857b2aSPeter Brune .seealso: SNESNASM, SNESNASMSetSubdomains()
45976857b2aSPeter Brune @*/
46076857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
46176857b2aSPeter Brune {
46276857b2aSPeter Brune   PetscErrorCode ierr;
46376857b2aSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**);
46476857b2aSPeter Brune 
46576857b2aSPeter Brune   PetscFunctionBegin;
4660005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f);CHKERRQ(ierr);
4674b838e8fSPeter Brune   if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);}
46876857b2aSPeter Brune   PetscFunctionReturn(0);
46976857b2aSPeter Brune }
47076857b2aSPeter Brune 
47176857b2aSPeter Brune #undef __FUNCT__
47276857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains_NASM"
47340244768SBarry Smith static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
47476857b2aSPeter Brune {
47576857b2aSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
47676857b2aSPeter Brune 
47776857b2aSPeter Brune   PetscFunctionBegin;
47876857b2aSPeter Brune   if (n) *n = nasm->n;
47976857b2aSPeter Brune   if (oscatter) *oscatter = nasm->oscatter;
48076857b2aSPeter Brune   if (iscatter) *iscatter = nasm->iscatter;
48176857b2aSPeter Brune   if (gscatter) *gscatter = nasm->gscatter;
482a4f17876SPeter Brune   if (subsnes)  {
483a4f17876SPeter Brune     *subsnes  = nasm->subsnes;
484a4f17876SPeter Brune     nasm->same_local_solves = PETSC_FALSE;
485a4f17876SPeter Brune   }
48676857b2aSPeter Brune   PetscFunctionReturn(0);
48776857b2aSPeter Brune }
48876857b2aSPeter Brune 
48976857b2aSPeter Brune #undef __FUNCT__
49076857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs"
49176857b2aSPeter Brune /*@
49276857b2aSPeter Brune    SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors
49376857b2aSPeter Brune 
49476857b2aSPeter Brune    Not Collective
49576857b2aSPeter Brune 
49676857b2aSPeter Brune    Input Parameters:
49776857b2aSPeter Brune .  SNES - the SNES context
49876857b2aSPeter Brune 
49976857b2aSPeter Brune    Output Parameters:
50076857b2aSPeter Brune +  n - the number of local subdomains
50176857b2aSPeter Brune .  x - The subdomain solution vector
50276857b2aSPeter Brune .  y - The subdomain step vector
50376857b2aSPeter Brune .  b - The subdomain RHS vector
50476857b2aSPeter Brune -  xl - The subdomain local vectors (ghosted)
50576857b2aSPeter Brune 
50676857b2aSPeter Brune    Level: developer
50776857b2aSPeter Brune 
50876857b2aSPeter Brune .keywords: SNES, NASM
50976857b2aSPeter Brune 
51076857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains()
51176857b2aSPeter Brune @*/
51276857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl)
51376857b2aSPeter Brune {
51476857b2aSPeter Brune   PetscErrorCode ierr;
51576857b2aSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**);
51676857b2aSPeter Brune 
51776857b2aSPeter Brune   PetscFunctionBegin;
5180005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f);CHKERRQ(ierr);
5194b838e8fSPeter Brune   if (f) {ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);}
52076857b2aSPeter Brune   PetscFunctionReturn(0);
52176857b2aSPeter Brune }
52276857b2aSPeter Brune 
52376857b2aSPeter Brune #undef __FUNCT__
52476857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs_NASM"
52540244768SBarry Smith static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl)
52676857b2aSPeter Brune {
52776857b2aSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
52876857b2aSPeter Brune 
52976857b2aSPeter Brune   PetscFunctionBegin;
53076857b2aSPeter Brune   if (n)  *n  = nasm->n;
53176857b2aSPeter Brune   if (x)  *x  = nasm->x;
53276857b2aSPeter Brune   if (y)  *y  = nasm->y;
53376857b2aSPeter Brune   if (b)  *b  = nasm->b;
53476857b2aSPeter Brune   if (xl) *xl = nasm->xl;
53576857b2aSPeter Brune   PetscFunctionReturn(0);
53676857b2aSPeter Brune }
53776857b2aSPeter Brune 
538d728fb7dSPeter Brune #undef __FUNCT__
539d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian"
540d728fb7dSPeter Brune /*@
541d728fb7dSPeter Brune    SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence
542d728fb7dSPeter Brune 
543d728fb7dSPeter Brune    Collective on SNES
544d728fb7dSPeter Brune 
545d728fb7dSPeter Brune    Input Parameters:
546d728fb7dSPeter Brune +  SNES - the SNES context
547d728fb7dSPeter Brune -  flg - indication of whether to compute the jacobians or not
548d728fb7dSPeter Brune 
549d728fb7dSPeter Brune    Level: developer
550d728fb7dSPeter Brune 
551d728fb7dSPeter Brune    Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian
552d728fb7dSPeter Brune    is needed at each linear iteration.
553d728fb7dSPeter Brune 
554d728fb7dSPeter Brune .keywords: SNES, NASM, ASPIN
555d728fb7dSPeter Brune 
556d728fb7dSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains()
557d728fb7dSPeter Brune @*/
558d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg)
559d728fb7dSPeter Brune {
560d728fb7dSPeter Brune   PetscErrorCode (*f)(SNES,PetscBool);
561d728fb7dSPeter Brune   PetscErrorCode ierr;
562d728fb7dSPeter Brune 
563d728fb7dSPeter Brune   PetscFunctionBegin;
5640005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f);CHKERRQ(ierr);
5654b838e8fSPeter Brune   if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);}
566d728fb7dSPeter Brune   PetscFunctionReturn(0);
567d728fb7dSPeter Brune }
568d728fb7dSPeter Brune 
569d728fb7dSPeter Brune #undef __FUNCT__
570d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian_NASM"
57140244768SBarry Smith static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg)
572d728fb7dSPeter Brune {
573d728fb7dSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
574d728fb7dSPeter Brune 
575d728fb7dSPeter Brune   PetscFunctionBegin;
576d728fb7dSPeter Brune   nasm->finaljacobian = flg;
577602bec5dSPeter Brune   if (flg) snes->usesksp = PETSC_TRUE;
578d728fb7dSPeter Brune   PetscFunctionReturn(0);
579d728fb7dSPeter Brune }
58076857b2aSPeter Brune 
58176857b2aSPeter Brune #undef __FUNCT__
582610116beSPeter Brune #define __FUNCT__ "SNESNASMSetDamping"
583610116beSPeter Brune /*@
584610116beSPeter Brune    SNESNASMSetDamping - Sets the update damping for NASM
585610116beSPeter Brune 
586610116beSPeter Brune    Logically collective on SNES
587610116beSPeter Brune 
588610116beSPeter Brune    Input Parameters:
589610116beSPeter Brune +  SNES - the SNES context
590610116beSPeter Brune -  dmp - damping
591610116beSPeter Brune 
592610116beSPeter Brune    Level: intermediate
593610116beSPeter Brune 
5945dfa0f3bSBarry Smith    Notes: The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
5955dfa0f3bSBarry Smith 
596610116beSPeter Brune .keywords: SNES, NASM, damping
597610116beSPeter Brune 
598610116beSPeter Brune .seealso: SNESNASM, SNESNASMGetDamping()
599610116beSPeter Brune @*/
600610116beSPeter Brune PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp)
601610116beSPeter Brune {
602610116beSPeter Brune   PetscErrorCode (*f)(SNES,PetscReal);
603610116beSPeter Brune   PetscErrorCode ierr;
604610116beSPeter Brune 
605610116beSPeter Brune   PetscFunctionBegin;
606610116beSPeter Brune   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr);
607e24b0d94SPeter Brune   if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);}
608610116beSPeter Brune   PetscFunctionReturn(0);
609610116beSPeter Brune }
610610116beSPeter Brune 
611610116beSPeter Brune #undef __FUNCT__
612610116beSPeter Brune #define __FUNCT__ "SNESNASMSetDamping_NASM"
61340244768SBarry Smith static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp)
614610116beSPeter Brune {
615610116beSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
616610116beSPeter Brune 
617610116beSPeter Brune   PetscFunctionBegin;
618610116beSPeter Brune   nasm->damping = dmp;
619610116beSPeter Brune   PetscFunctionReturn(0);
620610116beSPeter Brune }
621610116beSPeter Brune 
622610116beSPeter Brune #undef __FUNCT__
623610116beSPeter Brune #define __FUNCT__ "SNESNASMGetDamping"
624610116beSPeter Brune /*@
625610116beSPeter Brune    SNESNASMGetDamping - Gets the update damping for NASM
626610116beSPeter Brune 
627610116beSPeter Brune    Not Collective
628610116beSPeter Brune 
629610116beSPeter Brune    Input Parameters:
630610116beSPeter Brune +  SNES - the SNES context
631610116beSPeter Brune -  dmp - damping
632610116beSPeter Brune 
633610116beSPeter Brune    Level: intermediate
634610116beSPeter Brune 
635610116beSPeter Brune .keywords: SNES, NASM, damping
636610116beSPeter Brune 
637610116beSPeter Brune .seealso: SNESNASM, SNESNASMSetDamping()
638610116beSPeter Brune @*/
639610116beSPeter Brune PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp)
640610116beSPeter Brune {
641610116beSPeter Brune   PetscErrorCode ierr;
642610116beSPeter Brune 
643610116beSPeter Brune   PetscFunctionBegin;
6444a2f8832SBarry Smith   ierr = PetscUseMethod(snes,"SNESNASMGetDamping_C",(SNES,PetscReal*),(snes,dmp));CHKERRQ(ierr);
645610116beSPeter Brune   PetscFunctionReturn(0);
646610116beSPeter Brune }
647610116beSPeter Brune 
648610116beSPeter Brune #undef __FUNCT__
649610116beSPeter Brune #define __FUNCT__ "SNESNASMGetDamping_NASM"
65040244768SBarry Smith static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp)
651610116beSPeter Brune {
652610116beSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
653610116beSPeter Brune 
654610116beSPeter Brune   PetscFunctionBegin;
655610116beSPeter Brune   *dmp = nasm->damping;
656610116beSPeter Brune   PetscFunctionReturn(0);
657610116beSPeter Brune }
658610116beSPeter Brune 
659610116beSPeter Brune 
660610116beSPeter Brune #undef __FUNCT__
661eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSolveLocal_Private"
66214eb1c5cSMatthew G. Knepley /*
66314eb1c5cSMatthew G. Knepley   Input Parameters:
66414eb1c5cSMatthew G. Knepley + snes - The solver
66514eb1c5cSMatthew G. Knepley . B - The RHS vector
66614eb1c5cSMatthew G. Knepley - X - The initial guess
66714eb1c5cSMatthew G. Knepley 
66814eb1c5cSMatthew G. Knepley   Output Parameters:
66914eb1c5cSMatthew G. Knepley . Y - The solution update
67014eb1c5cSMatthew G. Knepley 
67114eb1c5cSMatthew G. Knepley   TODO: All scatters should be packed into one
67214eb1c5cSMatthew G. Knepley */
6730adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X)
6740adebc6cSBarry Smith {
675eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
676258e1594SPeter Brune   SNES           subsnes;
677eaedb033SPeter Brune   PetscInt       i;
678610116beSPeter Brune   PetscReal      dmp;
679eaedb033SPeter Brune   PetscErrorCode ierr;
680111ade9eSPeter Brune   Vec            Xlloc,Xl,Bl,Yl;
681111ade9eSPeter Brune   VecScatter     iscat,oscat,gscat;
682111ade9eSPeter Brune   DM             dm,subdm;
683e0331734SPeter Brune   PCASMType      type;
6840adebc6cSBarry Smith 
685eaedb033SPeter Brune   PetscFunctionBegin;
686e0331734SPeter Brune   ierr = SNESNASMGetType(snes,&type);CHKERRQ(ierr);
687eaedb033SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
688111ade9eSPeter Brune   ierr = VecSet(Y,0);CHKERRQ(ierr);
689b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
690eaedb033SPeter Brune   for (i=0; i<nasm->n; i++) {
69170c78f05SPeter Brune     /* scatter the solution to the local solution */
69270c78f05SPeter Brune     Xlloc = nasm->xl[i];
69370c78f05SPeter Brune     gscat   = nasm->gscatter[i];
69470c78f05SPeter Brune     oscat   = nasm->oscatter[i];
69570c78f05SPeter Brune     ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
69670c78f05SPeter Brune     if (B) {
69770c78f05SPeter Brune       /* scatter the RHS to the local RHS */
69870c78f05SPeter Brune       Bl   = nasm->b[i];
69970c78f05SPeter Brune       ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
70070c78f05SPeter Brune     }
70170c78f05SPeter Brune   }
702b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
703b20c023fSPeter Brune 
704b20c023fSPeter Brune 
705b20c023fSPeter Brune   if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
70670c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
70770c78f05SPeter Brune     Xl    = nasm->x[i];
70870c78f05SPeter Brune     Xlloc = nasm->xl[i];
70970c78f05SPeter Brune     Yl    = nasm->y[i];
710258e1594SPeter Brune     subsnes = nasm->subsnes[i];
711258e1594SPeter Brune     ierr    = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
712111ade9eSPeter Brune     iscat   = nasm->iscatter[i];
713111ade9eSPeter Brune     oscat   = nasm->oscatter[i];
714111ade9eSPeter Brune     gscat   = nasm->gscatter[i];
715ed3c11a9SPeter Brune     ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
71624b7f281SPeter Brune     if (B) {
71724b7f281SPeter Brune       Bl   = nasm->b[i];
718ed3c11a9SPeter Brune       ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
719ed3c11a9SPeter Brune     } else Bl = NULL;
720ed3c11a9SPeter Brune     ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
72114eb1c5cSMatthew G. Knepley     /* Could scatter directly from X */
72270c78f05SPeter Brune     ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
72370c78f05SPeter Brune     ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
724111ade9eSPeter Brune     ierr = VecCopy(Xl,Yl);CHKERRQ(ierr);
725ed3c11a9SPeter Brune     ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr);
726ed3c11a9SPeter Brune     ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr);
727e0331734SPeter Brune     if (type == PC_ASM_BASIC) {
728111ade9eSPeter Brune       ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
729e0331734SPeter Brune     } else if (type == PC_ASM_RESTRICT) {
730111ade9eSPeter Brune       ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
731ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
732eaedb033SPeter Brune   }
733ed3c11a9SPeter Brune   if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
734ed3c11a9SPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
73570c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
73670c78f05SPeter Brune     Yl    = nasm->y[i];
73770c78f05SPeter Brune     iscat   = nasm->iscatter[i];
73870c78f05SPeter Brune     oscat   = nasm->oscatter[i];
739e0331734SPeter Brune     if (type == PC_ASM_BASIC) {
74070c78f05SPeter Brune       ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
741e0331734SPeter Brune     } else if (type == PC_ASM_RESTRICT) {
74270c78f05SPeter Brune       ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
743ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
74470c78f05SPeter Brune   }
745b20c023fSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
7465dfa0f3bSBarry Smith   ierr = SNESNASMGetDamping(snes,&dmp);CHKERRQ(ierr);
747610116beSPeter Brune   ierr = VecAXPY(X,dmp,Y);CHKERRQ(ierr);
748eaedb033SPeter Brune   PetscFunctionReturn(0);
749eaedb033SPeter Brune }
750eaedb033SPeter Brune 
751eaedb033SPeter Brune #undef __FUNCT__
752d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMComputeFinalJacobian_Private"
75340244768SBarry Smith static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal)
754d728fb7dSPeter Brune {
755602bec5dSPeter Brune   Vec            X = Xfinal;
756d728fb7dSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
757d728fb7dSPeter Brune   SNES           subsnes;
758ca44f815SPeter Brune   PetscInt       i,lag = 1;
759d728fb7dSPeter Brune   PetscErrorCode ierr;
760e59f0a30SPeter Brune   Vec            Xlloc,Xl,Fl,F;
761d728fb7dSPeter Brune   VecScatter     oscat,gscat;
762d728fb7dSPeter Brune   DM             dm,subdm;
763d1e9a80fSBarry Smith 
764d728fb7dSPeter Brune   PetscFunctionBegin;
765602bec5dSPeter Brune   if (nasm->fjtype == 2) X = nasm->xinit;
766e59f0a30SPeter Brune   F = snes->vec_func;
767365a6726SPeter Brune   if (snes->normschedule == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);}
768d1e9a80fSBarry Smith   ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
769d728fb7dSPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
770d728fb7dSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
771602bec5dSPeter Brune   if (nasm->fjtype != 1) {
772d728fb7dSPeter Brune     for (i=0; i<nasm->n; i++) {
773d728fb7dSPeter Brune       Xlloc = nasm->xl[i];
774d728fb7dSPeter Brune       gscat = nasm->gscatter[i];
775d728fb7dSPeter Brune       ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
776602bec5dSPeter Brune     }
777d728fb7dSPeter Brune   }
778d728fb7dSPeter Brune   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
779d728fb7dSPeter Brune   for (i=0; i<nasm->n; i++) {
780e59f0a30SPeter Brune     Fl      = nasm->subsnes[i]->vec_func;
781d728fb7dSPeter Brune     Xl      = nasm->x[i];
782d728fb7dSPeter Brune     Xlloc   = nasm->xl[i];
783d728fb7dSPeter Brune     subsnes = nasm->subsnes[i];
784d728fb7dSPeter Brune     oscat   = nasm->oscatter[i];
785d728fb7dSPeter Brune     gscat   = nasm->gscatter[i];
786602bec5dSPeter Brune     if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);}
787ed3c11a9SPeter Brune     ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
788d728fb7dSPeter Brune     ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
789602bec5dSPeter Brune     if (nasm->fjtype != 1) {
790d728fb7dSPeter Brune       ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
791d728fb7dSPeter Brune       ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
792602bec5dSPeter Brune     }
793ca44f815SPeter Brune     if (subsnes->lagjacobian == -1)    subsnes->lagjacobian = -2;
794ca44f815SPeter Brune     else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian;
795602bec5dSPeter Brune     ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr);
796d1e9a80fSBarry Smith     ierr = SNESComputeJacobian(subsnes,Xl,subsnes->jacobian,subsnes->jacobian_pre);CHKERRQ(ierr);
797ca44f815SPeter Brune     if (lag > 1) subsnes->lagjacobian = lag;
798d728fb7dSPeter Brune   }
799d728fb7dSPeter Brune   PetscFunctionReturn(0);
800d728fb7dSPeter Brune }
801d728fb7dSPeter Brune 
802d728fb7dSPeter Brune #undef __FUNCT__
803eaedb033SPeter Brune #define __FUNCT__ "SNESSolve_NASM"
80440244768SBarry Smith static PetscErrorCode SNESSolve_NASM(SNES snes)
805eaedb033SPeter Brune {
806eaedb033SPeter Brune   Vec              F;
807eaedb033SPeter Brune   Vec              X;
808eaedb033SPeter Brune   Vec              B;
809111ade9eSPeter Brune   Vec              Y;
810eaedb033SPeter Brune   PetscInt         i;
811ed3c11a9SPeter Brune   PetscReal        fnorm = 0.0;
812eaedb033SPeter Brune   PetscErrorCode   ierr;
813365a6726SPeter Brune   SNESNormSchedule normschedule;
814d728fb7dSPeter Brune   SNES_NASM        *nasm = (SNES_NASM*)snes->data;
815eaedb033SPeter Brune 
816eaedb033SPeter Brune   PetscFunctionBegin;
817c579b300SPatrick Farrell 
8186c4ed002SBarry 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);
819c579b300SPatrick Farrell 
820fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
821eaedb033SPeter Brune   X = snes->vec_sol;
822111ade9eSPeter Brune   Y = snes->vec_sol_update;
823eaedb033SPeter Brune   F = snes->vec_func;
824eaedb033SPeter Brune   B = snes->vec_rhs;
825eaedb033SPeter Brune 
826e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
827eaedb033SPeter Brune   snes->iter   = 0;
828eaedb033SPeter Brune   snes->norm   = 0.;
829e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
830eaedb033SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
831365a6726SPeter Brune   ierr         = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr);
832365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
833eaedb033SPeter Brune     /* compute the initial function and preconditioned update delX */
834eaedb033SPeter Brune     if (!snes->vec_func_init_set) {
835eaedb033SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
8361aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
837eaedb033SPeter Brune 
838eaedb033SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
839422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
840e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
841eaedb033SPeter Brune     snes->iter = 0;
842eaedb033SPeter Brune     snes->norm = fnorm;
843e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
844a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
845eaedb033SPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
846eaedb033SPeter Brune 
847eaedb033SPeter Brune     /* test convergence */
848eaedb033SPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
849eaedb033SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
850eaedb033SPeter Brune   } else {
851e04113cfSBarry Smith     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
852a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
853eaedb033SPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
854eaedb033SPeter Brune   }
855eaedb033SPeter Brune 
856eaedb033SPeter Brune   /* Call general purpose update function */
857eaedb033SPeter Brune   if (snes->ops->update) {
858eaedb033SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
859eaedb033SPeter Brune   }
860602bec5dSPeter Brune   /* copy the initial solution over for later */
861602bec5dSPeter Brune   if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);}
862eaedb033SPeter Brune 
863eaedb033SPeter Brune   for (i = 0; i < snes->max_its; i++) {
864111ade9eSPeter Brune     ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr);
865365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
866eaedb033SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
867eaedb033SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
868422a814eSBarry Smith       SNESCheckFunctionNorm(snes,fnorm);
869eaedb033SPeter Brune     }
870eaedb033SPeter Brune     /* Monitor convergence */
871e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
872eaedb033SPeter Brune     snes->iter = i+1;
873eaedb033SPeter Brune     snes->norm = fnorm;
874e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
875a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
876eaedb033SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
877eaedb033SPeter Brune     /* Test for convergence */
878365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
879d728fb7dSPeter Brune     if (snes->reason) break;
880eaedb033SPeter Brune     /* Call general purpose update function */
881e59f0a30SPeter Brune     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
882eaedb033SPeter Brune   }
883d728fb7dSPeter Brune   if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);}
884365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS) {
885eaedb033SPeter Brune     if (i == snes->max_its) {
886eaedb033SPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
887eaedb033SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
888eaedb033SPeter Brune     }
8891aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */
890eaedb033SPeter Brune   PetscFunctionReturn(0);
891eaedb033SPeter Brune }
892eaedb033SPeter Brune 
893eaedb033SPeter Brune /*MC
894eaedb033SPeter Brune   SNESNASM - Nonlinear Additive Schwartz
895eaedb033SPeter Brune 
89670c78603SPeter Brune    Options Database:
89770c78603SPeter Brune +  -snes_nasm_log - enable logging events for the communication and solve stages
89870c78603SPeter Brune .  -snes_nasm_type <basic,restrict> - type of subdomain update used
8995dfa0f3bSBarry Smith .  -snes_asm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
90070c78603SPeter Brune .  -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate
901150d49b7SHong Zhang .  -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the jacobian is calculated at
90270c78603SPeter Brune .  -sub_snes_ - options prefix of the subdomain nonlinear solves
90370c78603SPeter Brune .  -sub_ksp_ - options prefix of the subdomain Krylov solver
90470c78603SPeter Brune -  -sub_pc_ - options prefix of the subdomain preconditioner
90570c78603SPeter Brune 
906eaedb033SPeter Brune    Level: advanced
907eaedb033SPeter Brune 
9084f02bc6aSBarry Smith    References:
90996a0c994SBarry Smith .  1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
9104f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
9114f02bc6aSBarry Smith 
9125dfa0f3bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types), SNESNASMSetType(), SNESNASMGetType(), SNESNASMSetSubdomains(), SNESNASMGetSubdomains(), SNESNASMGetSubdomainVecs(), SNESNASMSetComputeFinalJacobian(), SNESNASMSetDamping(), SNESNASMGetDamping()
913eaedb033SPeter Brune M*/
914eaedb033SPeter Brune 
915eaedb033SPeter Brune #undef __FUNCT__
916eaedb033SPeter Brune #define __FUNCT__ "SNESCreate_NASM"
9178cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes)
918eaedb033SPeter Brune {
919eaedb033SPeter Brune   SNES_NASM      *nasm;
920eaedb033SPeter Brune   PetscErrorCode ierr;
921eaedb033SPeter Brune 
922eaedb033SPeter Brune   PetscFunctionBegin;
923b00a9115SJed Brown   ierr       = PetscNewLog(snes,&nasm);CHKERRQ(ierr);
924eaedb033SPeter Brune   snes->data = (void*)nasm;
925eaedb033SPeter Brune 
926eaedb033SPeter Brune   nasm->n        = PETSC_DECIDE;
927eaedb033SPeter Brune   nasm->subsnes  = 0;
928eaedb033SPeter Brune   nasm->x        = 0;
929111ade9eSPeter Brune   nasm->xl       = 0;
930111ade9eSPeter Brune   nasm->y        = 0;
931eaedb033SPeter Brune   nasm->b        = 0;
932111ade9eSPeter Brune   nasm->oscatter = 0;
933111ade9eSPeter Brune   nasm->iscatter = 0;
934111ade9eSPeter Brune   nasm->gscatter = 0;
935610116beSPeter Brune   nasm->damping  = 1.;
936111ade9eSPeter Brune 
937111ade9eSPeter Brune   nasm->type = PC_ASM_BASIC;
938d728fb7dSPeter Brune   nasm->finaljacobian = PETSC_FALSE;
939a4f17876SPeter Brune   nasm->same_local_solves = PETSC_TRUE;
940eaedb033SPeter Brune 
941eaedb033SPeter Brune   snes->ops->destroy        = SNESDestroy_NASM;
942eaedb033SPeter Brune   snes->ops->setup          = SNESSetUp_NASM;
943eaedb033SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
944eaedb033SPeter Brune   snes->ops->view           = SNESView_NASM;
945eaedb033SPeter Brune   snes->ops->solve          = SNESSolve_NASM;
946eaedb033SPeter Brune   snes->ops->reset          = SNESReset_NASM;
947eaedb033SPeter Brune 
948eaedb033SPeter Brune   snes->usesksp = PETSC_FALSE;
949eaedb033SPeter Brune   snes->usespc  = PETSC_FALSE;
950eaedb033SPeter Brune 
951*4fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
952*4fc747eaSLawrence Mitchell 
953602bec5dSPeter Brune   nasm->fjtype              = 0;
954602bec5dSPeter Brune   nasm->xinit               = NULL;
9550298fd71SBarry Smith   nasm->eventrestrictinterp = 0;
9560298fd71SBarry Smith   nasm->eventsubsolve       = 0;
957b20c023fSPeter Brune 
958eaedb033SPeter Brune   if (!snes->tolerancesset) {
959eaedb033SPeter Brune     snes->max_its   = 10000;
960eaedb033SPeter Brune     snes->max_funcs = 10000;
961eaedb033SPeter Brune   }
962eaedb033SPeter Brune 
963e0331734SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetType_C",SNESNASMSetType_NASM);CHKERRQ(ierr);
964e0331734SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetType_C",SNESNASMGetType_NASM);CHKERRQ(ierr);
965bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr);
966bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr);
96790bcee39SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",SNESNASMSetDamping_NASM);CHKERRQ(ierr);
96890bcee39SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",SNESNASMGetDamping_NASM);CHKERRQ(ierr);
969bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr);
970bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr);
971eaedb033SPeter Brune   PetscFunctionReturn(0);
972eaedb033SPeter Brune }
97399e0435eSBarry Smith 
974