xref: /petsc/src/snes/impls/nasm/nasm.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
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 */
20f10b3e88SPatrick Farrell   PetscBool  weight_set;          /* use a weight in the overlap updates */
21b20c023fSPeter Brune 
22b20c023fSPeter Brune   /* logging events */
23b20c023fSPeter Brune   PetscLogEvent eventrestrictinterp;
24b20c023fSPeter Brune   PetscLogEvent eventsubsolve;
25602bec5dSPeter Brune 
26602bec5dSPeter Brune   PetscInt      fjtype;            /* type of computed jacobian */
27602bec5dSPeter Brune   Vec           xinit;             /* initial solution in case the final jacobian type is computed as first */
28eaedb033SPeter Brune } SNES_NASM;
29eaedb033SPeter Brune 
309e5d0892SLisandro Dalcin const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",NULL};
31602bec5dSPeter Brune const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"};
32b20c023fSPeter Brune 
3340244768SBarry Smith static PetscErrorCode SNESReset_NASM(SNES snes)
34eaedb033SPeter Brune {
35eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
36eaedb033SPeter Brune   PetscInt       i;
376e111a19SKarl Rupp 
38eaedb033SPeter Brune   PetscFunctionBegin;
39eaedb033SPeter Brune   for (i=0; i<nasm->n; i++) {
409566063dSJacob Faibussowitsch     if (nasm->xl) PetscCall(VecDestroy(&nasm->xl[i]));
419566063dSJacob Faibussowitsch     if (nasm->x) PetscCall(VecDestroy(&nasm->x[i]));
429566063dSJacob Faibussowitsch     if (nasm->y) PetscCall(VecDestroy(&nasm->y[i]));
439566063dSJacob Faibussowitsch     if (nasm->b) PetscCall(VecDestroy(&nasm->b[i]));
44eaedb033SPeter Brune 
459566063dSJacob Faibussowitsch     if (nasm->subsnes) PetscCall(SNESDestroy(&nasm->subsnes[i]));
469566063dSJacob Faibussowitsch     if (nasm->oscatter) PetscCall(VecScatterDestroy(&nasm->oscatter[i]));
479566063dSJacob Faibussowitsch     if (nasm->oscatter_copy) PetscCall(VecScatterDestroy(&nasm->oscatter_copy[i]));
489566063dSJacob Faibussowitsch     if (nasm->iscatter) PetscCall(VecScatterDestroy(&nasm->iscatter[i]));
499566063dSJacob Faibussowitsch     if (nasm->gscatter) PetscCall(VecScatterDestroy(&nasm->gscatter[i]));
50eaedb033SPeter Brune   }
51111ade9eSPeter Brune 
529566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->x));
539566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->xl));
549566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->y));
559566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->b));
56111ade9eSPeter Brune 
579566063dSJacob Faibussowitsch   if (nasm->xinit) PetscCall(VecDestroy(&nasm->xinit));
58602bec5dSPeter Brune 
599566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->subsnes));
609566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->oscatter));
619566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->oscatter_copy));
629566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->iscatter));
639566063dSJacob Faibussowitsch   PetscCall(PetscFree(nasm->gscatter));
64b20c023fSPeter Brune 
65f10b3e88SPatrick Farrell   if (nasm->weight_set) {
669566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&nasm->weight));
67f10b3e88SPatrick Farrell   }
68f10b3e88SPatrick Farrell 
69b20c023fSPeter Brune   nasm->eventrestrictinterp = 0;
70b20c023fSPeter Brune   nasm->eventsubsolve = 0;
71eaedb033SPeter Brune   PetscFunctionReturn(0);
72eaedb033SPeter Brune }
73eaedb033SPeter Brune 
7440244768SBarry Smith static PetscErrorCode SNESDestroy_NASM(SNES snes)
75eaedb033SPeter Brune {
76eaedb033SPeter Brune   PetscFunctionBegin;
779566063dSJacob Faibussowitsch   PetscCall(SNESReset_NASM(snes));
789566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
79eaedb033SPeter Brune   PetscFunctionReturn(0);
80eaedb033SPeter Brune }
81eaedb033SPeter Brune 
8240244768SBarry Smith static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx)
830adebc6cSBarry Smith {
84111ade9eSPeter Brune   Vec            bcs = (Vec)ctx;
856e111a19SKarl Rupp 
86111ade9eSPeter Brune   PetscFunctionBegin;
879566063dSJacob Faibussowitsch   PetscCall(VecCopy(bcs,l));
88111ade9eSPeter Brune   PetscFunctionReturn(0);
89111ade9eSPeter Brune }
90111ade9eSPeter Brune 
9140244768SBarry Smith static PetscErrorCode SNESSetUp_NASM(SNES snes)
92eaedb033SPeter Brune {
93eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
9476857b2aSPeter Brune   DM             dm,subdm;
95111ade9eSPeter Brune   DM             *subdms;
96111ade9eSPeter Brune   PetscInt       i;
97eaedb033SPeter Brune   const char     *optionsprefix;
98111ade9eSPeter Brune   Vec            F;
99ed3c11a9SPeter Brune   PetscMPIInt    size;
100ed3c11a9SPeter Brune   KSP            ksp;
101ed3c11a9SPeter Brune   PC             pc;
102eaedb033SPeter Brune 
103eaedb033SPeter Brune   PetscFunctionBegin;
104eaedb033SPeter Brune   if (!nasm->subsnes) {
1059566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(snes,&dm));
1060a696f66SPeter Brune     if (dm) {
107eaedb033SPeter Brune       nasm->usesdm = PETSC_TRUE;
1089566063dSJacob Faibussowitsch       PetscCall(DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms));
10928b400f6SJacob Faibussowitsch       PetscCheck(subdms,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
1109566063dSJacob Faibussowitsch       PetscCall(DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter));
1119566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nasm->n, &nasm->oscatter_copy));
112f10b3e88SPatrick Farrell       for (i=0; i<nasm->n; i++) {
1139566063dSJacob Faibussowitsch         PetscCall(VecScatterCopy(nasm->oscatter[i], &nasm->oscatter_copy[i]));
114f10b3e88SPatrick Farrell       }
115eaedb033SPeter Brune 
1169566063dSJacob Faibussowitsch       PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
1179566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nasm->n,&nasm->subsnes));
118111ade9eSPeter Brune       for (i=0; i<nasm->n; i++) {
1199566063dSJacob Faibussowitsch         PetscCall(SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]));
1209566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)nasm->subsnes[i], (PetscObject)snes, 1));
1219566063dSJacob Faibussowitsch         PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix));
1229566063dSJacob Faibussowitsch         PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_"));
1239566063dSJacob Faibussowitsch         PetscCall(SNESSetDM(nasm->subsnes[i],subdms[i]));
1249566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size));
125ed3c11a9SPeter Brune         if (size == 1) {
1269566063dSJacob Faibussowitsch           PetscCall(SNESGetKSP(nasm->subsnes[i],&ksp));
1279566063dSJacob Faibussowitsch           PetscCall(KSPGetPC(ksp,&pc));
1289566063dSJacob Faibussowitsch           PetscCall(KSPSetType(ksp,KSPPREONLY));
1299566063dSJacob Faibussowitsch           PetscCall(PCSetType(pc,PCLU));
130ed3c11a9SPeter Brune         }
1319566063dSJacob Faibussowitsch         PetscCall(SNESSetFromOptions(nasm->subsnes[i]));
1329566063dSJacob Faibussowitsch         PetscCall(DMDestroy(&subdms[i]));
133111ade9eSPeter Brune       }
1349566063dSJacob Faibussowitsch       PetscCall(PetscFree(subdms));
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) {
1399566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(nasm->n,&nasm->x));
140e245e0aaSPeter Brune   }
141e245e0aaSPeter Brune   if (!nasm->xl) {
1429566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(nasm->n,&nasm->xl));
143e245e0aaSPeter Brune   }
144e245e0aaSPeter Brune   if (!nasm->y) {
1459566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(nasm->n,&nasm->y));
146e245e0aaSPeter Brune   }
147e245e0aaSPeter Brune   if (!nasm->b) {
1489566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(nasm->n,&nasm->b));
149e245e0aaSPeter Brune   }
150111ade9eSPeter Brune 
151111ade9eSPeter Brune   for (i=0; i<nasm->n; i++) {
1529566063dSJacob Faibussowitsch     PetscCall(SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL));
1539566063dSJacob Faibussowitsch     if (!nasm->x[i]) PetscCall(VecDuplicate(F,&nasm->x[i]));
1549566063dSJacob Faibussowitsch     if (!nasm->y[i]) PetscCall(VecDuplicate(F,&nasm->y[i]));
1559566063dSJacob Faibussowitsch     if (!nasm->b[i]) PetscCall(VecDuplicate(F,&nasm->b[i]));
15676857b2aSPeter Brune     if (!nasm->xl[i]) {
1579566063dSJacob Faibussowitsch       PetscCall(SNESGetDM(nasm->subsnes[i],&subdm));
1589566063dSJacob Faibussowitsch       PetscCall(DMCreateLocalVector(subdm,&nasm->xl[i]));
1599566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]));
160111ade9eSPeter Brune     }
16161ba4676SBarry Smith   }
162602bec5dSPeter Brune   if (nasm->finaljacobian) {
1639566063dSJacob Faibussowitsch     PetscCall(SNESSetUpMatrices(snes));
164602bec5dSPeter Brune     if (nasm->fjtype == 2) {
1659566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(snes->vec_sol,&nasm->xinit));
166602bec5dSPeter Brune     }
167602bec5dSPeter Brune     for (i=0; i<nasm->n;i++) {
1689566063dSJacob Faibussowitsch       PetscCall(SNESSetUpMatrices(nasm->subsnes[i]));
169602bec5dSPeter Brune     }
170602bec5dSPeter Brune   }
171eaedb033SPeter Brune   PetscFunctionReturn(0);
172eaedb033SPeter Brune }
173eaedb033SPeter Brune 
17440244768SBarry Smith static PetscErrorCode SNESSetFromOptions_NASM(PetscOptionItems *PetscOptionsObject,SNES snes)
175eaedb033SPeter Brune {
176111ade9eSPeter Brune   PCASMType         asmtype;
17783dc3634SPierre Jolivet   PetscBool         flg,monflg;
178111ade9eSPeter Brune   SNES_NASM         *nasm = (SNES_NASM*)snes->data;
1796e111a19SKarl Rupp 
180eaedb033SPeter Brune   PetscFunctionBegin;
1819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"Nonlinear Additive Schwarz options"));
1829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg));
1839566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESNASMSetType(snes,asmtype));
184b20c023fSPeter Brune   flg    = PETSC_FALSE;
185b20c023fSPeter Brune   monflg = PETSC_TRUE;
1869566063dSJacob Faibussowitsch   PetscCall(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));
1879566063dSJacob Faibussowitsch   if (flg) PetscCall(SNESNASMSetDamping(snes,nasm->damping));
1889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDeprecated("-snes_nasm_sub_view",NULL,"3.15","Use -snes_view ::ascii_info_detail"));
1899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL));
1909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL));
1919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg));
192b20c023fSPeter Brune   if (flg) {
1939566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve));
1949566063dSJacob Faibussowitsch     PetscCall(PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp));
195b20c023fSPeter Brune   }
1969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
197eaedb033SPeter Brune   PetscFunctionReturn(0);
198eaedb033SPeter Brune }
199eaedb033SPeter Brune 
20040244768SBarry Smith static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
201eaedb033SPeter Brune {
202b20c023fSPeter Brune   SNES_NASM         *nasm = (SNES_NASM*)snes->data;
203a4f17876SPeter Brune   PetscMPIInt       rank,size;
204dd2fa690SBarry Smith   PetscInt          i,N,bsz;
205b20c023fSPeter Brune   PetscBool         iascii,isstring;
206b20c023fSPeter Brune   PetscViewer       sviewer;
207ce94432eSBarry Smith   MPI_Comm          comm;
20883dc3634SPierre Jolivet   PetscViewerFormat format;
20983dc3634SPierre Jolivet   const char        *prefix;
210b20c023fSPeter Brune 
211b20c023fSPeter Brune   PetscFunctionBegin;
2129566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes,&comm));
2139566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
2149566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring));
2159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
2169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
2179566063dSJacob Faibussowitsch   PetscCallMPI(MPIU_Allreduce(&nasm->n,&N,1,MPIU_INT,MPI_SUM,comm));
218b20c023fSPeter Brune   if (iascii) {
2199566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  total subdomain blocks = %D\n",N));
2209566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer,&format));
22183dc3634SPierre Jolivet     if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
222a4f17876SPeter Brune       if (nasm->subsnes) {
2239566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  Local solver information for first block on rank 0:\n"));
2249566063dSJacob Faibussowitsch         PetscCall(SNESGetOptionsPrefix(snes,&prefix));
2259566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  Use -%ssnes_view ::ascii_info_detail to display information for all blocks\n",prefix?prefix:""));
2269566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
2279566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
228dd400576SPatrick Sanan         if (rank == 0) {
2299566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
2309566063dSJacob Faibussowitsch           PetscCall(SNESView(nasm->subsnes[0],sviewer));
2319566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
232a4f17876SPeter Brune         }
2339566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
2349566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
235a4f17876SPeter Brune       }
236a4f17876SPeter Brune     } else {
237a4f17876SPeter Brune       /* print the solver on each block */
2389566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
2399566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,nasm->n));
2409566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2419566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2429566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  Local solver information for each block is in the following SNES objects:\n"));
2439566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
2449566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n"));
2459566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
246a4f17876SPeter Brune       for (i=0; i<nasm->n; i++) {
2479566063dSJacob Faibussowitsch         PetscCall(VecGetLocalSize(nasm->x[i],&bsz));
2489566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz));
2499566063dSJacob Faibussowitsch         PetscCall(SNESView(nasm->subsnes[i],sviewer));
2509566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n"));
251b20c023fSPeter Brune       }
2529566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
2539566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2549566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
255a4f17876SPeter Brune     }
256b20c023fSPeter Brune   } else if (isstring) {
2579566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]));
2589566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
2599566063dSJacob Faibussowitsch     if (nasm->subsnes && rank == 0) PetscCall(SNESView(nasm->subsnes[0],sviewer));
2609566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
261b20c023fSPeter Brune   }
262eaedb033SPeter Brune   PetscFunctionReturn(0);
263eaedb033SPeter Brune }
264eaedb033SPeter Brune 
265e0331734SPeter Brune /*@
266e0331734SPeter Brune    SNESNASMSetType - Set the type of subdomain update used
267e0331734SPeter Brune 
268e0331734SPeter Brune    Logically Collective on SNES
269e0331734SPeter Brune 
270e0331734SPeter Brune    Input Parameters:
271e0331734SPeter Brune +  SNES - the SNES context
272e0331734SPeter Brune -  type - the type of update, PC_ASM_BASIC or PC_ASM_RESTRICT
273e0331734SPeter Brune 
274e0331734SPeter Brune    Level: intermediate
275e0331734SPeter Brune 
276e0331734SPeter Brune .seealso: SNESNASM, SNESNASMGetType(), PCASMSetType()
277e0331734SPeter Brune @*/
278e0331734SPeter Brune PetscErrorCode SNESNASMSetType(SNES snes,PCASMType type)
279e0331734SPeter Brune {
280e0331734SPeter Brune   PetscErrorCode (*f)(SNES,PCASMType);
281e0331734SPeter Brune 
282e0331734SPeter Brune   PetscFunctionBegin;
2839566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetType_C",&f));
2849566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes,type));
285e0331734SPeter Brune   PetscFunctionReturn(0);
286e0331734SPeter Brune }
287e0331734SPeter Brune 
28840244768SBarry Smith static PetscErrorCode SNESNASMSetType_NASM(SNES snes,PCASMType type)
289e0331734SPeter Brune {
290e0331734SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
291e0331734SPeter Brune 
292e0331734SPeter Brune   PetscFunctionBegin;
2932c71b3e2SJacob Faibussowitsch   PetscCheckFalse(type != PC_ASM_BASIC && type != PC_ASM_RESTRICT,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESNASM only supports basic and restrict types");
294e0331734SPeter Brune   nasm->type = type;
295e0331734SPeter Brune   PetscFunctionReturn(0);
296e0331734SPeter Brune }
297e0331734SPeter Brune 
298e0331734SPeter Brune /*@
299e0331734SPeter Brune    SNESNASMGetType - Get the type of subdomain update used
300e0331734SPeter Brune 
301e0331734SPeter Brune    Logically Collective on SNES
302e0331734SPeter Brune 
303e0331734SPeter Brune    Input Parameters:
304e0331734SPeter Brune .  SNES - the SNES context
305e0331734SPeter Brune 
306e0331734SPeter Brune    Output Parameters:
307e0331734SPeter Brune .  type - the type of update
308e0331734SPeter Brune 
309e0331734SPeter Brune    Level: intermediate
310e0331734SPeter Brune 
311e0331734SPeter Brune .seealso: SNESNASM, SNESNASMSetType(), PCASMGetType()
312e0331734SPeter Brune @*/
313e0331734SPeter Brune PetscErrorCode SNESNASMGetType(SNES snes,PCASMType *type)
314e0331734SPeter Brune {
315e0331734SPeter Brune   PetscFunctionBegin;
316*cac4c232SBarry Smith   PetscUseMethod(snes,"SNESNASMGetType_C",(SNES,PCASMType*),(snes,type));
317e0331734SPeter Brune   PetscFunctionReturn(0);
318e0331734SPeter Brune }
319e0331734SPeter Brune 
32040244768SBarry Smith static PetscErrorCode SNESNASMGetType_NASM(SNES snes,PCASMType *type)
321e0331734SPeter Brune {
322e0331734SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
323e0331734SPeter Brune 
324e0331734SPeter Brune   PetscFunctionBegin;
325e0331734SPeter Brune   *type = nasm->type;
326e0331734SPeter Brune   PetscFunctionReturn(0);
327e0331734SPeter Brune }
328e0331734SPeter Brune 
32976857b2aSPeter Brune /*@
33076857b2aSPeter Brune    SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems.
33176857b2aSPeter Brune 
33276857b2aSPeter Brune    Not Collective
33376857b2aSPeter Brune 
33476857b2aSPeter Brune    Input Parameters:
33576857b2aSPeter Brune +  SNES - the SNES context
33676857b2aSPeter Brune .  n - the number of local subdomains
33776857b2aSPeter Brune .  subsnes - solvers defined on the local subdomains
33876857b2aSPeter Brune .  iscatter - scatters into the nonoverlapping portions of the local subdomains
33976857b2aSPeter Brune .  oscatter - scatters into the overlapping portions of the local subdomains
34076857b2aSPeter Brune -  gscatter - scatters into the (ghosted) local vector of the local subdomain
34176857b2aSPeter Brune 
34276857b2aSPeter Brune    Level: intermediate
34376857b2aSPeter Brune 
34476857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains()
34576857b2aSPeter Brune @*/
346a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
347a6dfd86eSKarl Rupp {
348111ade9eSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*);
3496e111a19SKarl Rupp 
350eaedb033SPeter Brune   PetscFunctionBegin;
3519566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f));
3529566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes,n,subsnes,iscatter,oscatter,gscatter));
353eaedb033SPeter Brune   PetscFunctionReturn(0);
354eaedb033SPeter Brune }
355eaedb033SPeter Brune 
35640244768SBarry Smith static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
357a6dfd86eSKarl Rupp {
358eaedb033SPeter Brune   PetscInt       i;
359eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
3606e111a19SKarl Rupp 
361eaedb033SPeter Brune   PetscFunctionBegin;
36228b400f6SJacob Faibussowitsch   PetscCheck(!snes->setupcalled,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
363eaedb033SPeter Brune 
364111ade9eSPeter Brune   /* tear down the previously set things */
3659566063dSJacob Faibussowitsch   PetscCall(SNESReset(snes));
366111ade9eSPeter Brune 
367eaedb033SPeter Brune   nasm->n = n;
368111ade9eSPeter Brune   if (oscatter) {
3699566063dSJacob Faibussowitsch     for (i=0; i<n; i++) PetscCall(PetscObjectReference((PetscObject)oscatter[i]));
370eaedb033SPeter Brune   }
371111ade9eSPeter Brune   if (iscatter) {
3729566063dSJacob Faibussowitsch     for (i=0; i<n; i++) PetscCall(PetscObjectReference((PetscObject)iscatter[i]));
373eaedb033SPeter Brune   }
374111ade9eSPeter Brune   if (gscatter) {
3759566063dSJacob Faibussowitsch     for (i=0; i<n; i++) PetscCall(PetscObjectReference((PetscObject)gscatter[i]));
376111ade9eSPeter Brune   }
377111ade9eSPeter Brune   if (oscatter) {
3789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,&nasm->oscatter));
3799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,&nasm->oscatter_copy));
380eaedb033SPeter Brune     for (i=0; i<n; i++) {
381111ade9eSPeter Brune       nasm->oscatter[i] = oscatter[i];
3829566063dSJacob Faibussowitsch       PetscCall(VecScatterCopy(oscatter[i], &nasm->oscatter_copy[i]));
383eaedb033SPeter Brune     }
384111ade9eSPeter Brune   }
385111ade9eSPeter Brune   if (iscatter) {
3869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,&nasm->iscatter));
387eaedb033SPeter Brune     for (i=0; i<n; i++) {
388111ade9eSPeter Brune       nasm->iscatter[i] = iscatter[i];
389eaedb033SPeter Brune     }
390eaedb033SPeter Brune   }
391111ade9eSPeter Brune   if (gscatter) {
3929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,&nasm->gscatter));
393eaedb033SPeter Brune     for (i=0; i<n; i++) {
394111ade9eSPeter Brune       nasm->gscatter[i] = gscatter[i];
395eaedb033SPeter Brune     }
396eaedb033SPeter Brune   }
397111ade9eSPeter Brune 
398eaedb033SPeter Brune   if (subsnes) {
3999566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,&nasm->subsnes));
400eaedb033SPeter Brune     for (i=0; i<n; i++) {
401eaedb033SPeter Brune       nasm->subsnes[i] = subsnes[i];
402eaedb033SPeter Brune     }
403eaedb033SPeter Brune   }
404eaedb033SPeter Brune   PetscFunctionReturn(0);
405eaedb033SPeter Brune }
406eaedb033SPeter Brune 
40776857b2aSPeter Brune /*@
40876857b2aSPeter Brune    SNESNASMGetSubdomains - Get the local subdomain context.
40976857b2aSPeter Brune 
41076857b2aSPeter Brune    Not Collective
41176857b2aSPeter Brune 
412f899ff85SJose E. Roman    Input Parameter:
41376857b2aSPeter Brune .  SNES - the SNES context
41476857b2aSPeter Brune 
41576857b2aSPeter Brune    Output Parameters:
41676857b2aSPeter Brune +  n - the number of local subdomains
41776857b2aSPeter Brune .  subsnes - solvers defined on the local subdomains
41876857b2aSPeter Brune .  iscatter - scatters into the nonoverlapping portions of the local subdomains
41976857b2aSPeter Brune .  oscatter - scatters into the overlapping portions of the local subdomains
42076857b2aSPeter Brune -  gscatter - scatters into the (ghosted) local vector of the local subdomain
42176857b2aSPeter Brune 
42276857b2aSPeter Brune    Level: intermediate
42376857b2aSPeter Brune 
42476857b2aSPeter Brune .seealso: SNESNASM, SNESNASMSetSubdomains()
42576857b2aSPeter Brune @*/
42676857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
42776857b2aSPeter Brune {
42876857b2aSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**);
42976857b2aSPeter Brune 
43076857b2aSPeter Brune   PetscFunctionBegin;
4319566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f));
4329566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes,n,subsnes,iscatter,oscatter,gscatter));
43376857b2aSPeter Brune   PetscFunctionReturn(0);
43476857b2aSPeter Brune }
43576857b2aSPeter Brune 
43640244768SBarry Smith static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
43776857b2aSPeter Brune {
43876857b2aSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
43976857b2aSPeter Brune 
44076857b2aSPeter Brune   PetscFunctionBegin;
44176857b2aSPeter Brune   if (n) *n = nasm->n;
44276857b2aSPeter Brune   if (oscatter) *oscatter = nasm->oscatter;
44376857b2aSPeter Brune   if (iscatter) *iscatter = nasm->iscatter;
44476857b2aSPeter Brune   if (gscatter) *gscatter = nasm->gscatter;
44583dc3634SPierre Jolivet   if (subsnes)  *subsnes  = nasm->subsnes;
44676857b2aSPeter Brune   PetscFunctionReturn(0);
44776857b2aSPeter Brune }
44876857b2aSPeter Brune 
44976857b2aSPeter Brune /*@
45076857b2aSPeter Brune    SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors
45176857b2aSPeter Brune 
45276857b2aSPeter Brune    Not Collective
45376857b2aSPeter Brune 
454f899ff85SJose E. Roman    Input Parameter:
45576857b2aSPeter Brune .  SNES - the SNES context
45676857b2aSPeter Brune 
45776857b2aSPeter Brune    Output Parameters:
45876857b2aSPeter Brune +  n - the number of local subdomains
45976857b2aSPeter Brune .  x - The subdomain solution vector
46076857b2aSPeter Brune .  y - The subdomain step vector
46176857b2aSPeter Brune .  b - The subdomain RHS vector
46276857b2aSPeter Brune -  xl - The subdomain local vectors (ghosted)
46376857b2aSPeter Brune 
46476857b2aSPeter Brune    Level: developer
46576857b2aSPeter Brune 
46676857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains()
46776857b2aSPeter Brune @*/
46876857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl)
46976857b2aSPeter Brune {
47076857b2aSPeter Brune   PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**);
47176857b2aSPeter Brune 
47276857b2aSPeter Brune   PetscFunctionBegin;
4739566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f));
4749566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes,n,x,y,b,xl));
47576857b2aSPeter Brune   PetscFunctionReturn(0);
47676857b2aSPeter Brune }
47776857b2aSPeter Brune 
47840244768SBarry Smith static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl)
47976857b2aSPeter Brune {
48076857b2aSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
48176857b2aSPeter Brune 
48276857b2aSPeter Brune   PetscFunctionBegin;
48376857b2aSPeter Brune   if (n)  *n  = nasm->n;
48476857b2aSPeter Brune   if (x)  *x  = nasm->x;
48576857b2aSPeter Brune   if (y)  *y  = nasm->y;
48676857b2aSPeter Brune   if (b)  *b  = nasm->b;
48776857b2aSPeter Brune   if (xl) *xl = nasm->xl;
48876857b2aSPeter Brune   PetscFunctionReturn(0);
48976857b2aSPeter Brune }
49076857b2aSPeter Brune 
491d728fb7dSPeter Brune /*@
4928bf45196SRichard Tran Mills    SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain Jacobians upon convergence
493d728fb7dSPeter Brune 
494d728fb7dSPeter Brune    Collective on SNES
495d728fb7dSPeter Brune 
496d728fb7dSPeter Brune    Input Parameters:
497d728fb7dSPeter Brune +  SNES - the SNES context
4988bf45196SRichard Tran Mills -  flg - indication of whether to compute the Jacobians or not
499d728fb7dSPeter Brune 
500d728fb7dSPeter Brune    Level: developer
501d728fb7dSPeter Brune 
50295452b02SPatrick Sanan    Notes:
5038bf45196SRichard Tran Mills    This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global Jacobian
504d728fb7dSPeter Brune    is needed at each linear iteration.
505d728fb7dSPeter Brune 
506d728fb7dSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains()
507d728fb7dSPeter Brune @*/
508d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg)
509d728fb7dSPeter Brune {
510d728fb7dSPeter Brune   PetscErrorCode (*f)(SNES,PetscBool);
511d728fb7dSPeter Brune 
512d728fb7dSPeter Brune   PetscFunctionBegin;
5139566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f));
5149566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes,flg));
515d728fb7dSPeter Brune   PetscFunctionReturn(0);
516d728fb7dSPeter Brune }
517d728fb7dSPeter Brune 
51840244768SBarry Smith static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg)
519d728fb7dSPeter Brune {
520d728fb7dSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
521d728fb7dSPeter Brune 
522d728fb7dSPeter Brune   PetscFunctionBegin;
523d728fb7dSPeter Brune   nasm->finaljacobian = flg;
524d728fb7dSPeter Brune   PetscFunctionReturn(0);
525d728fb7dSPeter Brune }
52676857b2aSPeter Brune 
527610116beSPeter Brune /*@
528610116beSPeter Brune    SNESNASMSetDamping - Sets the update damping for NASM
529610116beSPeter Brune 
530610116beSPeter Brune    Logically collective on SNES
531610116beSPeter Brune 
532610116beSPeter Brune    Input Parameters:
533610116beSPeter Brune +  SNES - the SNES context
534610116beSPeter Brune -  dmp - damping
535610116beSPeter Brune 
536610116beSPeter Brune    Level: intermediate
537610116beSPeter Brune 
53895452b02SPatrick Sanan    Notes:
53995452b02SPatrick Sanan     The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
5405dfa0f3bSBarry Smith 
541610116beSPeter Brune .seealso: SNESNASM, SNESNASMGetDamping()
542610116beSPeter Brune @*/
543610116beSPeter Brune PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp)
544610116beSPeter Brune {
545610116beSPeter Brune   PetscErrorCode (*f)(SNES,PetscReal);
546610116beSPeter Brune 
547610116beSPeter Brune   PetscFunctionBegin;
5489566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f));
5499566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes,dmp));
550610116beSPeter Brune   PetscFunctionReturn(0);
551610116beSPeter Brune }
552610116beSPeter Brune 
55340244768SBarry Smith static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp)
554610116beSPeter Brune {
555610116beSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
556610116beSPeter Brune 
557610116beSPeter Brune   PetscFunctionBegin;
558610116beSPeter Brune   nasm->damping = dmp;
559610116beSPeter Brune   PetscFunctionReturn(0);
560610116beSPeter Brune }
561610116beSPeter Brune 
562610116beSPeter Brune /*@
563610116beSPeter Brune    SNESNASMGetDamping - Gets the update damping for NASM
564610116beSPeter Brune 
565610116beSPeter Brune    Not Collective
566610116beSPeter Brune 
567610116beSPeter Brune    Input Parameters:
568610116beSPeter Brune +  SNES - the SNES context
569610116beSPeter Brune -  dmp - damping
570610116beSPeter Brune 
571610116beSPeter Brune    Level: intermediate
572610116beSPeter Brune 
573610116beSPeter Brune .seealso: SNESNASM, SNESNASMSetDamping()
574610116beSPeter Brune @*/
575610116beSPeter Brune PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp)
576610116beSPeter Brune {
577610116beSPeter Brune   PetscFunctionBegin;
578*cac4c232SBarry Smith   PetscUseMethod(snes,"SNESNASMGetDamping_C",(SNES,PetscReal*),(snes,dmp));
579610116beSPeter Brune   PetscFunctionReturn(0);
580610116beSPeter Brune }
581610116beSPeter Brune 
58240244768SBarry Smith static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp)
583610116beSPeter Brune {
584610116beSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
585610116beSPeter Brune 
586610116beSPeter Brune   PetscFunctionBegin;
587610116beSPeter Brune   *dmp = nasm->damping;
588610116beSPeter Brune   PetscFunctionReturn(0);
589610116beSPeter Brune }
590610116beSPeter Brune 
59114eb1c5cSMatthew G. Knepley /*
59214eb1c5cSMatthew G. Knepley   Input Parameters:
59314eb1c5cSMatthew G. Knepley + snes - The solver
59414eb1c5cSMatthew G. Knepley . B - The RHS vector
59514eb1c5cSMatthew G. Knepley - X - The initial guess
59614eb1c5cSMatthew G. Knepley 
59714eb1c5cSMatthew G. Knepley   Output Parameters:
59814eb1c5cSMatthew G. Knepley . Y - The solution update
59914eb1c5cSMatthew G. Knepley 
60014eb1c5cSMatthew G. Knepley   TODO: All scatters should be packed into one
60114eb1c5cSMatthew G. Knepley */
6020adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X)
6030adebc6cSBarry Smith {
604eaedb033SPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
605258e1594SPeter Brune   SNES           subsnes;
606eaedb033SPeter Brune   PetscInt       i;
607610116beSPeter Brune   PetscReal      dmp;
608f10b3e88SPatrick Farrell   Vec            Xl,Bl,Yl,Xlloc;
609f10b3e88SPatrick Farrell   VecScatter     iscat,oscat,gscat,oscat_copy;
610111ade9eSPeter Brune   DM             dm,subdm;
611e0331734SPeter Brune   PCASMType      type;
6120adebc6cSBarry Smith 
613eaedb033SPeter Brune   PetscFunctionBegin;
6149566063dSJacob Faibussowitsch   PetscCall(SNESNASMGetType(snes,&type));
6159566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
6169566063dSJacob Faibussowitsch   PetscCall(VecSet(Y,0));
6179566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0));
618eaedb033SPeter Brune   for (i=0; i<nasm->n; i++) {
619f10b3e88SPatrick Farrell     /* scatter the solution to the global solution and the local solution */
620f10b3e88SPatrick Farrell     Xl      = nasm->x[i];
62170c78f05SPeter Brune     Xlloc   = nasm->xl[i];
62270c78f05SPeter Brune     oscat   = nasm->oscatter[i];
623f10b3e88SPatrick Farrell     oscat_copy = nasm->oscatter_copy[i];
624f10b3e88SPatrick Farrell     gscat   = nasm->gscatter[i];
6259566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD));
6269566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD));
62770c78f05SPeter Brune     if (B) {
62870c78f05SPeter Brune       /* scatter the RHS to the local RHS */
62970c78f05SPeter Brune       Bl   = nasm->b[i];
6309566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(oscat_copy,B,Bl,INSERT_VALUES,SCATTER_FORWARD));
63170c78f05SPeter Brune     }
63270c78f05SPeter Brune   }
6339566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0));
634b20c023fSPeter Brune 
6359566063dSJacob Faibussowitsch   if (nasm->eventsubsolve) PetscCall(PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0));
63670c78f05SPeter Brune   for (i=0; i<nasm->n; i++) {
63770c78f05SPeter Brune     Xl    = nasm->x[i];
63870c78f05SPeter Brune     Xlloc = nasm->xl[i];
63970c78f05SPeter Brune     Yl    = nasm->y[i];
640258e1594SPeter Brune     subsnes = nasm->subsnes[i];
6419566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(subsnes,&subdm));
642111ade9eSPeter Brune     iscat   = nasm->iscatter[i];
643111ade9eSPeter Brune     oscat   = nasm->oscatter[i];
644f10b3e88SPatrick Farrell     oscat_copy = nasm->oscatter_copy[i];
645111ade9eSPeter Brune     gscat   = nasm->gscatter[i];
6469566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD));
6479566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD));
64824b7f281SPeter Brune     if (B) {
64924b7f281SPeter Brune       Bl   = nasm->b[i];
6509566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(oscat_copy,B,Bl,INSERT_VALUES,SCATTER_FORWARD));
651ed3c11a9SPeter Brune     } else Bl = NULL;
652f10b3e88SPatrick Farrell 
6539566063dSJacob Faibussowitsch     PetscCall(DMSubDomainRestrict(dm,oscat,gscat,subdm));
6549566063dSJacob Faibussowitsch     PetscCall(VecCopy(Xl,Yl));
6559566063dSJacob Faibussowitsch     PetscCall(SNESSolve(subsnes,Bl,Xl));
6569566063dSJacob Faibussowitsch     PetscCall(VecAYPX(Yl,-1.0,Xl));
6579566063dSJacob Faibussowitsch     PetscCall(VecScale(Yl, nasm->damping));
658e0331734SPeter Brune     if (type == PC_ASM_BASIC) {
6599566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE));
6609566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE));
661e0331734SPeter Brune     } else if (type == PC_ASM_RESTRICT) {
6629566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE));
6639566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE));
664ce94432eSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
665eaedb033SPeter Brune   }
6669566063dSJacob Faibussowitsch   if (nasm->eventsubsolve) PetscCall(PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0));
6679566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0));
668f10b3e88SPatrick Farrell   if (nasm->weight_set) {
6699566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Y,Y,nasm->weight));
670f10b3e88SPatrick Farrell   }
6719566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0));
6729566063dSJacob Faibussowitsch   PetscCall(SNESNASMGetDamping(snes,&dmp));
6739566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X,dmp,Y));
674eaedb033SPeter Brune   PetscFunctionReturn(0);
675eaedb033SPeter Brune }
676eaedb033SPeter Brune 
67740244768SBarry Smith static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal)
678d728fb7dSPeter Brune {
679602bec5dSPeter Brune   Vec            X = Xfinal;
680d728fb7dSPeter Brune   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
681d728fb7dSPeter Brune   SNES           subsnes;
682ca44f815SPeter Brune   PetscInt       i,lag = 1;
683e59f0a30SPeter Brune   Vec            Xlloc,Xl,Fl,F;
684d728fb7dSPeter Brune   VecScatter     oscat,gscat;
685d728fb7dSPeter Brune   DM             dm,subdm;
686d1e9a80fSBarry Smith 
687d728fb7dSPeter Brune   PetscFunctionBegin;
688602bec5dSPeter Brune   if (nasm->fjtype == 2) X = nasm->xinit;
689e59f0a30SPeter Brune   F = snes->vec_func;
6909566063dSJacob Faibussowitsch   if (snes->normschedule == SNES_NORM_NONE) PetscCall(SNESComputeFunction(snes,X,F));
6919566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre));
6929566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
6939566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0));
694602bec5dSPeter Brune   if (nasm->fjtype != 1) {
695d728fb7dSPeter Brune     for (i=0; i<nasm->n; i++) {
696d728fb7dSPeter Brune       Xlloc = nasm->xl[i];
697d728fb7dSPeter Brune       gscat = nasm->gscatter[i];
6989566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD));
699602bec5dSPeter Brune     }
700d728fb7dSPeter Brune   }
7019566063dSJacob Faibussowitsch   if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0));
702d728fb7dSPeter Brune   for (i=0; i<nasm->n; i++) {
703e59f0a30SPeter Brune     Fl      = nasm->subsnes[i]->vec_func;
704d728fb7dSPeter Brune     Xl      = nasm->x[i];
705d728fb7dSPeter Brune     Xlloc   = nasm->xl[i];
706d728fb7dSPeter Brune     subsnes = nasm->subsnes[i];
707d728fb7dSPeter Brune     oscat   = nasm->oscatter[i];
708d728fb7dSPeter Brune     gscat   = nasm->gscatter[i];
7099566063dSJacob Faibussowitsch     if (nasm->fjtype != 1) PetscCall(VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD));
7109566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(subsnes,&subdm));
7119566063dSJacob Faibussowitsch     PetscCall(DMSubDomainRestrict(dm,oscat,gscat,subdm));
712602bec5dSPeter Brune     if (nasm->fjtype != 1) {
7139566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl));
7149566063dSJacob Faibussowitsch       PetscCall(DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl));
715602bec5dSPeter Brune     }
716ca44f815SPeter Brune     if (subsnes->lagjacobian == -1)    subsnes->lagjacobian = -2;
717ca44f815SPeter Brune     else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian;
7189566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(subsnes,Xl,Fl));
7199566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(subsnes,Xl,subsnes->jacobian,subsnes->jacobian_pre));
720ca44f815SPeter Brune     if (lag > 1) subsnes->lagjacobian = lag;
721d728fb7dSPeter Brune   }
722d728fb7dSPeter Brune   PetscFunctionReturn(0);
723d728fb7dSPeter Brune }
724d728fb7dSPeter Brune 
72540244768SBarry Smith static PetscErrorCode SNESSolve_NASM(SNES snes)
726eaedb033SPeter Brune {
727eaedb033SPeter Brune   Vec              F;
728eaedb033SPeter Brune   Vec              X;
729eaedb033SPeter Brune   Vec              B;
730111ade9eSPeter Brune   Vec              Y;
731eaedb033SPeter Brune   PetscInt         i;
732ed3c11a9SPeter Brune   PetscReal        fnorm = 0.0;
733365a6726SPeter Brune   SNESNormSchedule normschedule;
734d728fb7dSPeter Brune   SNES_NASM        *nasm = (SNES_NASM*)snes->data;
735eaedb033SPeter Brune 
736eaedb033SPeter Brune   PetscFunctionBegin;
737c579b300SPatrick Farrell 
7382c71b3e2SJacob Faibussowitsch   PetscCheckFalse(snes->xl || snes->xu || snes->ops->computevariablebounds,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
739c579b300SPatrick Farrell 
7409566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite));
741eaedb033SPeter Brune   X = snes->vec_sol;
742111ade9eSPeter Brune   Y = snes->vec_sol_update;
743eaedb033SPeter Brune   F = snes->vec_func;
744eaedb033SPeter Brune   B = snes->vec_rhs;
745eaedb033SPeter Brune 
7469566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
747eaedb033SPeter Brune   snes->iter   = 0;
748eaedb033SPeter Brune   snes->norm   = 0.;
7499566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
750eaedb033SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
7519566063dSJacob Faibussowitsch   PetscCall(SNESGetNormSchedule(snes, &normschedule));
752365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
753eaedb033SPeter Brune     /* compute the initial function and preconditioned update delX */
754eaedb033SPeter Brune     if (!snes->vec_func_init_set) {
7559566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes,X,F));
7561aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
757eaedb033SPeter Brune 
7589566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
759422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
7609566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
761eaedb033SPeter Brune     snes->iter = 0;
762eaedb033SPeter Brune     snes->norm = fnorm;
7639566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7649566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0));
7659566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes,0,snes->norm));
766eaedb033SPeter Brune 
767eaedb033SPeter Brune     /* test convergence */
7689566063dSJacob Faibussowitsch     PetscCall((*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP));
769eaedb033SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
770eaedb033SPeter Brune   } else {
7719566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7729566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0));
7739566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes,0,snes->norm));
774eaedb033SPeter Brune   }
775eaedb033SPeter Brune 
776eaedb033SPeter Brune   /* Call general purpose update function */
777eaedb033SPeter Brune   if (snes->ops->update) {
7789566063dSJacob Faibussowitsch     PetscCall((*snes->ops->update)(snes, snes->iter));
779eaedb033SPeter Brune   }
780602bec5dSPeter Brune   /* copy the initial solution over for later */
7819566063dSJacob Faibussowitsch   if (nasm->fjtype == 2) PetscCall(VecCopy(X,nasm->xinit));
782eaedb033SPeter Brune 
783eaedb033SPeter Brune   for (i=0; i < snes->max_its; i++) {
7849566063dSJacob Faibussowitsch     PetscCall(SNESNASMSolveLocal_Private(snes,B,Y,X));
785365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
7869566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes,X,F));
7879566063dSJacob Faibussowitsch       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
788422a814eSBarry Smith       SNESCheckFunctionNorm(snes,fnorm);
789eaedb033SPeter Brune     }
790eaedb033SPeter Brune     /* Monitor convergence */
7919566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
792eaedb033SPeter Brune     snes->iter = i+1;
793eaedb033SPeter Brune     snes->norm = fnorm;
7949566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
7959566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0));
7969566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes,snes->iter,snes->norm));
797eaedb033SPeter Brune     /* Test for convergence */
7989566063dSJacob Faibussowitsch     if (normschedule == SNES_NORM_ALWAYS) PetscCall((*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP));
799d728fb7dSPeter Brune     if (snes->reason) break;
800eaedb033SPeter Brune     /* Call general purpose update function */
8019566063dSJacob Faibussowitsch     if (snes->ops->update) PetscCall((*snes->ops->update)(snes, snes->iter));
802eaedb033SPeter Brune   }
80307b62357SFande Kong   if (nasm->finaljacobian) {
8049566063dSJacob Faibussowitsch     PetscCall(SNESNASMComputeFinalJacobian_Private(snes,X));
80507b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
80607b62357SFande Kong   }
807365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS) {
808eaedb033SPeter Brune     if (i == snes->max_its) {
8099566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its));
810eaedb033SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
811eaedb033SPeter Brune     }
8121aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */
813eaedb033SPeter Brune   PetscFunctionReturn(0);
814eaedb033SPeter Brune }
815eaedb033SPeter Brune 
816eaedb033SPeter Brune /*MC
817f3f228e0SPierre Jolivet   SNESNASM - Nonlinear Additive Schwarz
818eaedb033SPeter Brune 
81970c78603SPeter Brune    Options Database:
82070c78603SPeter Brune +  -snes_nasm_log - enable logging events for the communication and solve stages
82170c78603SPeter Brune .  -snes_nasm_type <basic,restrict> - type of subdomain update used
8225dfa0f3bSBarry Smith .  -snes_asm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)
82370c78603SPeter Brune .  -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate
824150d49b7SHong Zhang .  -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the jacobian is calculated at
82570c78603SPeter Brune .  -sub_snes_ - options prefix of the subdomain nonlinear solves
82670c78603SPeter Brune .  -sub_ksp_ - options prefix of the subdomain Krylov solver
82770c78603SPeter Brune -  -sub_pc_ - options prefix of the subdomain preconditioner
82870c78603SPeter Brune 
829eaedb033SPeter Brune    Level: advanced
830eaedb033SPeter Brune 
831951fe5abSBarry Smith    Developer Note: This is a non-Newton based nonlinear solver that does not directly require a Jacobian; hence the flag snes->usesksp is set to
8327f851d37SRichard Tran Mills        false and SNESView() and -snes_view do not display a KSP object. However, if the flag nasm->finaljacobian is set (for example, if
833951fe5abSBarry Smith        NASM is used as a nonlinear preconditioner for  KSPASPIN) then SNESSetUpMatrices() is called to generate the Jacobian (needed by KSPASPIN)
834951fe5abSBarry Smith        and this utilizes the KSP for storing the matrices, but the KSP is never used for solving a linear system. Note that when SNESNASM is
835951fe5abSBarry Smith        used by SNESASPIN they share the same Jacobian matrices because SNESSetUp() (called on the outer SNES KSPASPIN) causes the inner SNES
836951fe5abSBarry Smith        object (in this case SNESNASM) to inherit the outer Jacobian matrices.
837951fe5abSBarry Smith 
8384f02bc6aSBarry Smith    References:
839606c0280SSatish Balay .  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
8404f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
8414f02bc6aSBarry Smith 
8425dfa0f3bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types), SNESNASMSetType(), SNESNASMGetType(), SNESNASMSetSubdomains(), SNESNASMGetSubdomains(), SNESNASMGetSubdomainVecs(), SNESNASMSetComputeFinalJacobian(), SNESNASMSetDamping(), SNESNASMGetDamping()
843eaedb033SPeter Brune M*/
844eaedb033SPeter Brune 
8458cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes)
846eaedb033SPeter Brune {
847eaedb033SPeter Brune   SNES_NASM      *nasm;
848eaedb033SPeter Brune 
849eaedb033SPeter Brune   PetscFunctionBegin;
8509566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes,&nasm));
851eaedb033SPeter Brune   snes->data = (void*)nasm;
852eaedb033SPeter Brune 
853eaedb033SPeter Brune   nasm->n        = PETSC_DECIDE;
8549e5d0892SLisandro Dalcin   nasm->subsnes  = NULL;
8559e5d0892SLisandro Dalcin   nasm->x        = NULL;
8569e5d0892SLisandro Dalcin   nasm->xl       = NULL;
8579e5d0892SLisandro Dalcin   nasm->y        = NULL;
8589e5d0892SLisandro Dalcin   nasm->b        = NULL;
8599e5d0892SLisandro Dalcin   nasm->oscatter = NULL;
8609e5d0892SLisandro Dalcin   nasm->oscatter_copy = NULL;
8619e5d0892SLisandro Dalcin   nasm->iscatter = NULL;
8629e5d0892SLisandro Dalcin   nasm->gscatter = NULL;
863610116beSPeter Brune   nasm->damping  = 1.;
864111ade9eSPeter Brune 
865111ade9eSPeter Brune   nasm->type              = PC_ASM_BASIC;
866d728fb7dSPeter Brune   nasm->finaljacobian     = PETSC_FALSE;
867f10b3e88SPatrick Farrell   nasm->weight_set        = PETSC_FALSE;
868eaedb033SPeter Brune 
869eaedb033SPeter Brune   snes->ops->destroy        = SNESDestroy_NASM;
870eaedb033SPeter Brune   snes->ops->setup          = SNESSetUp_NASM;
871eaedb033SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
872eaedb033SPeter Brune   snes->ops->view           = SNESView_NASM;
873eaedb033SPeter Brune   snes->ops->solve          = SNESSolve_NASM;
874eaedb033SPeter Brune   snes->ops->reset          = SNESReset_NASM;
875eaedb033SPeter Brune 
876eaedb033SPeter Brune   snes->usesksp = PETSC_FALSE;
877efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
878eaedb033SPeter Brune 
8794fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
8804fc747eaSLawrence Mitchell 
881602bec5dSPeter Brune   nasm->fjtype              = 0;
882602bec5dSPeter Brune   nasm->xinit               = NULL;
8830298fd71SBarry Smith   nasm->eventrestrictinterp = 0;
8840298fd71SBarry Smith   nasm->eventsubsolve       = 0;
885b20c023fSPeter Brune 
886eaedb033SPeter Brune   if (!snes->tolerancesset) {
887eaedb033SPeter Brune     snes->max_its   = 10000;
888eaedb033SPeter Brune     snes->max_funcs = 10000;
889eaedb033SPeter Brune   }
890eaedb033SPeter Brune 
8919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetType_C",SNESNASMSetType_NASM));
8929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetType_C",SNESNASMGetType_NASM));
8939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM));
8949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM));
8959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",SNESNASMSetDamping_NASM));
8969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",SNESNASMGetDamping_NASM));
8979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM));
8989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM));
899eaedb033SPeter Brune   PetscFunctionReturn(0);
900eaedb033SPeter Brune }
90199e0435eSBarry Smith 
902448b6425SPatrick Farrell /*@
903448b6425SPatrick Farrell    SNESNASMGetSNES - Gets a subsolver
904448b6425SPatrick Farrell 
905448b6425SPatrick Farrell    Not collective
906448b6425SPatrick Farrell 
907448b6425SPatrick Farrell    Input Parameters:
908448b6425SPatrick Farrell +  snes - the SNES context
909448b6425SPatrick Farrell -  i - the number of the subsnes to get
910448b6425SPatrick Farrell 
911448b6425SPatrick Farrell    Output Parameters:
912448b6425SPatrick Farrell .  subsnes - the subsolver context
913448b6425SPatrick Farrell 
914448b6425SPatrick Farrell    Level: intermediate
915448b6425SPatrick Farrell 
916448b6425SPatrick Farrell .seealso: SNESNASM, SNESNASMGetNumber()
917448b6425SPatrick Farrell @*/
918448b6425SPatrick Farrell PetscErrorCode SNESNASMGetSNES(SNES snes,PetscInt i,SNES *subsnes)
919448b6425SPatrick Farrell {
920448b6425SPatrick Farrell   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
921448b6425SPatrick Farrell 
922448b6425SPatrick Farrell   PetscFunctionBegin;
9232c71b3e2SJacob Faibussowitsch   PetscCheckFalse(i < 0 || i >= nasm->n,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"No such subsolver");
924448b6425SPatrick Farrell   *subsnes = nasm->subsnes[i];
925448b6425SPatrick Farrell   PetscFunctionReturn(0);
926448b6425SPatrick Farrell }
927448b6425SPatrick Farrell 
928448b6425SPatrick Farrell /*@
929448b6425SPatrick Farrell    SNESNASMGetNumber - Gets number of subsolvers
930448b6425SPatrick Farrell 
931448b6425SPatrick Farrell    Not collective
932448b6425SPatrick Farrell 
933448b6425SPatrick Farrell    Input Parameters:
934448b6425SPatrick Farrell .  snes - the SNES context
935448b6425SPatrick Farrell 
936448b6425SPatrick Farrell    Output Parameters:
937448b6425SPatrick Farrell .  n - the number of subsolvers
938448b6425SPatrick Farrell 
939448b6425SPatrick Farrell    Level: intermediate
940448b6425SPatrick Farrell 
941448b6425SPatrick Farrell .seealso: SNESNASM, SNESNASMGetSNES()
942448b6425SPatrick Farrell @*/
943448b6425SPatrick Farrell PetscErrorCode SNESNASMGetNumber(SNES snes,PetscInt *n)
944448b6425SPatrick Farrell {
945448b6425SPatrick Farrell   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
946448b6425SPatrick Farrell 
947448b6425SPatrick Farrell   PetscFunctionBegin;
948448b6425SPatrick Farrell   *n = nasm->n;
949448b6425SPatrick Farrell   PetscFunctionReturn(0);
950448b6425SPatrick Farrell }
951448b6425SPatrick Farrell 
952f10b3e88SPatrick Farrell /*@
953f10b3e88SPatrick Farrell    SNESNASMSetWeight - Sets weight to use when adding overlapping updates
954f10b3e88SPatrick Farrell 
955f10b3e88SPatrick Farrell    Collective
956f10b3e88SPatrick Farrell 
957f10b3e88SPatrick Farrell    Input Parameters:
958f10b3e88SPatrick Farrell +  snes - the SNES context
959f10b3e88SPatrick Farrell -  weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in)
960f10b3e88SPatrick Farrell 
961f10b3e88SPatrick Farrell    Level: intermediate
962f10b3e88SPatrick Farrell 
963f10b3e88SPatrick Farrell .seealso: SNESNASM
964f10b3e88SPatrick Farrell @*/
965f10b3e88SPatrick Farrell PetscErrorCode SNESNASMSetWeight(SNES snes,Vec weight)
966f10b3e88SPatrick Farrell {
967f10b3e88SPatrick Farrell   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
968f10b3e88SPatrick Farrell 
969f10b3e88SPatrick Farrell   PetscFunctionBegin;
970f10b3e88SPatrick Farrell 
9719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&nasm->weight));
972f10b3e88SPatrick Farrell   nasm->weight_set = PETSC_TRUE;
973f10b3e88SPatrick Farrell   nasm->weight     = weight;
9749566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)nasm->weight));
975f10b3e88SPatrick Farrell 
976f10b3e88SPatrick Farrell   PetscFunctionReturn(0);
977f10b3e88SPatrick Farrell }
978