xref: /petsc/src/ksp/pc/impls/asm/asm.c (revision 265a2120fd39442a845a3b0193ce5532a9b543cd)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith   This file defines an additive Schwarz preconditioner for any Mat implementation.
44b9ad928SBarry Smith 
54b9ad928SBarry Smith   Note that each processor may have any number of subdomains. But in order to
64b9ad928SBarry Smith   deal easily with the VecScatter(), we treat each processor as if it has the
74b9ad928SBarry Smith   same number of subdomains.
84b9ad928SBarry Smith 
94b9ad928SBarry Smith        n - total number of true subdomains on all processors
104b9ad928SBarry Smith        n_local_true - actual number of subdomains on this processor
114b9ad928SBarry Smith        n_local = maximum over all processors of n_local_true
124b9ad928SBarry Smith */
13af0996ceSBarry Smith #include <petsc/private/pcimpl.h>     /*I "petscpc.h" I*/
141e25c274SJed Brown #include <petscdm.h>
154b9ad928SBarry Smith 
164b9ad928SBarry Smith typedef struct {
1713f74950SBarry Smith   PetscInt   n, n_local, n_local_true;
1813f74950SBarry Smith   PetscInt   overlap;             /* overlap requested by user */
194b9ad928SBarry Smith   KSP        *ksp;                /* linear solvers for each block */
20b0de9f23SBarry Smith   VecScatter restriction;         /* mapping from global to overlapping (process) subdomain*/
21b0de9f23SBarry Smith   VecScatter *lrestriction;       /* mapping from subregion to overlapping (process) subdomain */
22b0de9f23SBarry Smith   VecScatter *lprolongation;      /* mapping from non-overlapping subregion to overlapping (process) subdomain; used for restrict additive version of algorithms */
231dd8081eSeaulisa   Vec        lx, ly;              /* work vectors */
241dd8081eSeaulisa   Vec        *x,*y;  		  /* work vectors */
251dd8081eSeaulisa   IS         lis;                 /* index set that defines each overlapping multiplicative (process) subdomain */
262b691e39SMatthew Knepley   IS         *is;                 /* index set that defines each overlapping subdomain */
272b691e39SMatthew Knepley   IS         *is_local;           /* index set that defines each non-overlapping subdomain, may be NULL */
284b9ad928SBarry Smith   Mat        *mat,*pmat;          /* mat is not currently used */
294b9ad928SBarry Smith   PCASMType  type;                /* use reduced interpolation, restriction or both */
30ace3abfcSBarry Smith   PetscBool  type_set;            /* if user set this value (so won't change it for symmetric problems) */
31ace3abfcSBarry Smith   PetscBool  same_local_solves;   /* flag indicating whether all local solvers are same */
32ace3abfcSBarry Smith   PetscBool  sort_indices;        /* flag to sort subdomain indices */
33d709ea83SDmitry Karpeev   PetscBool  dm_subdomains;       /* whether DM is allowed to define subdomains */
3412cd4985SMatthew G. Knepley   PCCompositeType loctype;        /* the type of composition for local solves */
3580ec0b7dSPatrick Sanan   MatType    sub_mat_type;        /* the type of Mat used for subdomain solves (can be MATSAME or NULL) */
36fb745f2cSMatthew G. Knepley   /* For multiplicative solve */
37235cc4ceSMatthew G. Knepley   Mat       *lmats;               /* submatrices for overlapping multiplicative (process) subdomain */
384b9ad928SBarry Smith } PC_ASM;
394b9ad928SBarry Smith 
406849ba73SBarry Smith static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer)
414b9ad928SBarry Smith {
4292bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
436849ba73SBarry Smith   PetscErrorCode ierr;
4413f74950SBarry Smith   PetscMPIInt    rank;
454d219a6aSLisandro Dalcin   PetscInt       i,bsz;
46ace3abfcSBarry Smith   PetscBool      iascii,isstring;
474b9ad928SBarry Smith   PetscViewer    sviewer;
484b9ad928SBarry Smith 
494b9ad928SBarry Smith   PetscFunctionBegin;
50251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
51251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
5232077d6dSBarry Smith   if (iascii) {
533d03bb48SJed Brown     char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set";
548caf3d72SBarry Smith     if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof(overlaps),"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);}
558caf3d72SBarry Smith     if (osm->n > 0) {ierr = PetscSNPrintf(blocks,sizeof(blocks),"total subdomain blocks = %D",osm->n);CHKERRQ(ierr);}
56efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  %s, %s\n",blocks,overlaps);CHKERRQ(ierr);
57efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  restriction/interpolation type - %s\n",PCASMTypes[osm->type]);CHKERRQ(ierr);
58e64f0791SPatrick Sanan     if (osm->dm_subdomains) {ierr = PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: using DM to define subdomains\n");CHKERRQ(ierr);}
5912cd4985SMatthew G. Knepley     if (osm->loctype != PC_COMPOSITE_ADDITIVE) {ierr = PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: local solve composition type - %s\n",PCCompositeTypes[osm->loctype]);CHKERRQ(ierr);}
60ce94432eSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
6192bb6962SLisandro Dalcin     if (osm->same_local_solves) {
624d219a6aSLisandro Dalcin       if (osm->ksp) {
634b9ad928SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Local solve is same for all blocks, in the following KSP and PC objects:\n");CHKERRQ(ierr);
64c5e4d11fSDmitry Karpeev         ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
654d219a6aSLisandro Dalcin         if (!rank) {
664b9ad928SBarry Smith           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
6792bb6962SLisandro Dalcin           ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);
684b9ad928SBarry Smith           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
694b9ad928SBarry Smith         }
70c5e4d11fSDmitry Karpeev         ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
714d219a6aSLisandro Dalcin       }
724b9ad928SBarry Smith     } else {
73c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
744d219a6aSLisandro Dalcin       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,osm->n_local_true);CHKERRQ(ierr);
754d219a6aSLisandro Dalcin       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
764b9ad928SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr);
774b9ad928SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
784d219a6aSLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
79c5e4d11fSDmitry Karpeev       ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
80dd2fa690SBarry Smith       for (i=0; i<osm->n_local_true; i++) {
814d219a6aSLisandro Dalcin         ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr);
824d219a6aSLisandro Dalcin         ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr);
8392bb6962SLisandro Dalcin         ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr);
844d219a6aSLisandro Dalcin         ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
854b9ad928SBarry Smith       }
86c5e4d11fSDmitry Karpeev       ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
874b9ad928SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
884b9ad928SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
89c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
904b9ad928SBarry Smith     }
914b9ad928SBarry Smith   } else if (isstring) {
924d219a6aSLisandro Dalcin     ierr = PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCASMTypes[osm->type]);CHKERRQ(ierr);
93c5e4d11fSDmitry Karpeev     ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
9492bb6962SLisandro Dalcin     if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);}
95c5e4d11fSDmitry Karpeev     ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
964b9ad928SBarry Smith   }
974b9ad928SBarry Smith   PetscFunctionReturn(0);
984b9ad928SBarry Smith }
994b9ad928SBarry Smith 
10092bb6962SLisandro Dalcin static PetscErrorCode PCASMPrintSubdomains(PC pc)
10192bb6962SLisandro Dalcin {
10292bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
10392bb6962SLisandro Dalcin   const char     *prefix;
10492bb6962SLisandro Dalcin   char           fname[PETSC_MAX_PATH_LEN+1];
105643535aeSDmitry Karpeev   PetscViewer    viewer, sviewer;
106643535aeSDmitry Karpeev   char           *s;
10792bb6962SLisandro Dalcin   PetscInt       i,j,nidx;
10892bb6962SLisandro Dalcin   const PetscInt *idx;
109643535aeSDmitry Karpeev   PetscMPIInt    rank, size;
11092bb6962SLisandro Dalcin   PetscErrorCode ierr;
111846783a0SBarry Smith 
11292bb6962SLisandro Dalcin   PetscFunctionBegin;
113ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr);
114ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr);
11592bb6962SLisandro Dalcin   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
116c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
11792bb6962SLisandro Dalcin   if (fname[0] == 0) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
118ce94432eSBarry Smith   ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr);
119643535aeSDmitry Karpeev   for (i=0; i<osm->n_local; i++) {
120643535aeSDmitry Karpeev     if (i < osm->n_local_true) {
12192bb6962SLisandro Dalcin       ierr = ISGetLocalSize(osm->is[i],&nidx);CHKERRQ(ierr);
12292bb6962SLisandro Dalcin       ierr = ISGetIndices(osm->is[i],&idx);CHKERRQ(ierr);
123643535aeSDmitry Karpeev       /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
124854ce69bSBarry Smith       ierr = PetscMalloc1(16*(nidx+1)+512, &s);CHKERRQ(ierr);
125643535aeSDmitry Karpeev       ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, 16*(nidx+1)+512, &sviewer);CHKERRQ(ierr);
126643535aeSDmitry Karpeev       ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D with overlap:\n", rank, size, i);CHKERRQ(ierr);
12792bb6962SLisandro Dalcin       for (j=0; j<nidx; j++) {
128643535aeSDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr);
12992bb6962SLisandro Dalcin       }
13092bb6962SLisandro Dalcin       ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr);
131643535aeSDmitry Karpeev       ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr);
132643535aeSDmitry Karpeev       ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
133c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
134643535aeSDmitry Karpeev       ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr);
135643535aeSDmitry Karpeev       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
136c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
137643535aeSDmitry Karpeev       ierr = PetscFree(s);CHKERRQ(ierr);
1382b691e39SMatthew Knepley       if (osm->is_local) {
139643535aeSDmitry Karpeev         /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
140854ce69bSBarry Smith         ierr = PetscMalloc1(16*(nidx+1)+512, &s);CHKERRQ(ierr);
141643535aeSDmitry Karpeev         ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, 16*(nidx+1)+512, &sviewer);CHKERRQ(ierr);
14209d011a0SDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D without overlap:\n", rank, size, i);CHKERRQ(ierr);
1432b691e39SMatthew Knepley         ierr = ISGetLocalSize(osm->is_local[i],&nidx);CHKERRQ(ierr);
1442b691e39SMatthew Knepley         ierr = ISGetIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
1452b691e39SMatthew Knepley         for (j=0; j<nidx; j++) {
14609d011a0SDmitry Karpeev           ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr);
1472b691e39SMatthew Knepley         }
1482b691e39SMatthew Knepley         ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
14909d011a0SDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr);
150643535aeSDmitry Karpeev         ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
151c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
152643535aeSDmitry Karpeev         ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr);
153643535aeSDmitry Karpeev         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
154c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
155643535aeSDmitry Karpeev         ierr = PetscFree(s);CHKERRQ(ierr);
156643535aeSDmitry Karpeev       }
1572fa5cd67SKarl Rupp     } else {
158643535aeSDmitry Karpeev       /* Participate in collective viewer calls. */
159c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
160643535aeSDmitry Karpeev       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
161c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
162643535aeSDmitry Karpeev       /* Assume either all ranks have is_local or none do. */
163643535aeSDmitry Karpeev       if (osm->is_local) {
164c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
165643535aeSDmitry Karpeev         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
166c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
167643535aeSDmitry Karpeev       }
168fdc48646SDmitry Karpeev     }
16992bb6962SLisandro Dalcin   }
17092bb6962SLisandro Dalcin   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
171fcfd50ebSBarry Smith   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
17292bb6962SLisandro Dalcin   PetscFunctionReturn(0);
17392bb6962SLisandro Dalcin }
17492bb6962SLisandro Dalcin 
1756849ba73SBarry Smith static PetscErrorCode PCSetUp_ASM(PC pc)
1764b9ad928SBarry Smith {
1774b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
1786849ba73SBarry Smith   PetscErrorCode ierr;
179ace3abfcSBarry Smith   PetscBool      symset,flg;
18087e86712SBarry Smith   PetscInt       i,m,m_local;
1814b9ad928SBarry Smith   MatReuse       scall = MAT_REUSE_MATRIX;
1824b9ad928SBarry Smith   IS             isl;
1834b9ad928SBarry Smith   KSP            ksp;
1844b9ad928SBarry Smith   PC             subpc;
1852dcb1b2aSMatthew Knepley   const char     *prefix,*pprefix;
18623ce1328SBarry Smith   Vec            vec;
1870298fd71SBarry Smith   DM             *domain_dm = NULL;
1884b9ad928SBarry Smith 
1894b9ad928SBarry Smith   PetscFunctionBegin;
1904b9ad928SBarry Smith   if (!pc->setupcalled) {
191*265a2120SBarry Smith     PetscInt m;
19292bb6962SLisandro Dalcin 
19392bb6962SLisandro Dalcin     if (!osm->type_set) {
19492bb6962SLisandro Dalcin       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
1952fa5cd67SKarl Rupp       if (symset && flg) osm->type = PC_ASM_BASIC;
19692bb6962SLisandro Dalcin     }
19792bb6962SLisandro Dalcin 
1982b837212SDmitry Karpeev     /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */
1992b837212SDmitry Karpeev     if (osm->n_local_true == PETSC_DECIDE) {
20069ca1f37SDmitry Karpeev       /* no subdomains given */
20165db9045SDmitry Karpeev       /* try pc->dm first, if allowed */
202d709ea83SDmitry Karpeev       if (osm->dm_subdomains && pc->dm) {
203feb221f8SDmitry Karpeev         PetscInt  num_domains, d;
204feb221f8SDmitry Karpeev         char      **domain_names;
2058d4ac253SDmitry Karpeev         IS        *inner_domain_is, *outer_domain_is;
2068d4ac253SDmitry Karpeev         ierr = DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm);CHKERRQ(ierr);
207704f0395SPatrick Sanan         osm->overlap = -1; /* We do not want to increase the overlap of the IS.
208704f0395SPatrick Sanan                               A future improvement of this code might allow one to use
209704f0395SPatrick Sanan                               DM-defined subdomains and also increase the overlap,
210704f0395SPatrick Sanan                               but that is not currently supported */
21169ca1f37SDmitry Karpeev         if (num_domains) {
2128d4ac253SDmitry Karpeev           ierr = PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is);CHKERRQ(ierr);
21369ca1f37SDmitry Karpeev         }
214feb221f8SDmitry Karpeev         for (d = 0; d < num_domains; ++d) {
215a35b7b57SDmitry Karpeev           if (domain_names)    {ierr = PetscFree(domain_names[d]);CHKERRQ(ierr);}
216a35b7b57SDmitry Karpeev           if (inner_domain_is) {ierr = ISDestroy(&inner_domain_is[d]);CHKERRQ(ierr);}
217a35b7b57SDmitry Karpeev           if (outer_domain_is) {ierr = ISDestroy(&outer_domain_is[d]);CHKERRQ(ierr);}
218feb221f8SDmitry Karpeev         }
219feb221f8SDmitry Karpeev         ierr = PetscFree(domain_names);CHKERRQ(ierr);
2208d4ac253SDmitry Karpeev         ierr = PetscFree(inner_domain_is);CHKERRQ(ierr);
2218d4ac253SDmitry Karpeev         ierr = PetscFree(outer_domain_is);CHKERRQ(ierr);
222feb221f8SDmitry Karpeev       }
2232b837212SDmitry Karpeev       if (osm->n_local_true == PETSC_DECIDE) {
22469ca1f37SDmitry Karpeev         /* still no subdomains; use one subdomain per processor */
2252b837212SDmitry Karpeev         osm->n_local_true = 1;
226feb221f8SDmitry Karpeev       }
2272b837212SDmitry Karpeev     }
2282b837212SDmitry Karpeev     { /* determine the global and max number of subdomains */
2296ac3741eSJed Brown       struct {PetscInt max,sum;} inwork,outwork;
230c10200c1SHong Zhang       PetscMPIInt size;
231c10200c1SHong Zhang 
2326ac3741eSJed Brown       inwork.max   = osm->n_local_true;
2336ac3741eSJed Brown       inwork.sum   = osm->n_local_true;
234367daffbSBarry Smith       ierr         = MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
2356ac3741eSJed Brown       osm->n_local = outwork.max;
2366ac3741eSJed Brown       osm->n       = outwork.sum;
237c10200c1SHong Zhang 
238c10200c1SHong Zhang       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
239c10200c1SHong Zhang       if (outwork.max == 1 && outwork.sum == size) {
2407dae84e0SHong Zhang         /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */
241c10200c1SHong Zhang         ierr = MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE);CHKERRQ(ierr);
242c10200c1SHong Zhang       }
2434b9ad928SBarry Smith     }
24478904715SLisandro Dalcin     if (!osm->is) { /* create the index sets */
24578904715SLisandro Dalcin       ierr = PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);CHKERRQ(ierr);
2464b9ad928SBarry Smith     }
247f5234e35SJed Brown     if (osm->n_local_true > 1 && !osm->is_local) {
248785e854fSJed Brown       ierr = PetscMalloc1(osm->n_local_true,&osm->is_local);CHKERRQ(ierr);
249f5234e35SJed Brown       for (i=0; i<osm->n_local_true; i++) {
250f5234e35SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
251f5234e35SJed Brown           ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr);
252f5234e35SJed Brown           ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr);
253f5234e35SJed Brown         } else {
254f5234e35SJed Brown           ierr             = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr);
255f5234e35SJed Brown           osm->is_local[i] = osm->is[i];
256f5234e35SJed Brown         }
257f5234e35SJed Brown       }
258f5234e35SJed Brown     }
25992bb6962SLisandro Dalcin     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
26090d69ab7SBarry Smith     flg  = PETSC_FALSE;
261c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,prefix,"-pc_asm_print_subdomains",&flg,NULL);CHKERRQ(ierr);
26292bb6962SLisandro Dalcin     if (flg) { ierr = PCASMPrintSubdomains(pc);CHKERRQ(ierr); }
2634b9ad928SBarry Smith 
2643d03bb48SJed Brown     if (osm->overlap > 0) {
2654b9ad928SBarry Smith       /* Extend the "overlapping" regions by a number of steps */
26692bb6962SLisandro Dalcin       ierr = MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);CHKERRQ(ierr);
2673d03bb48SJed Brown     }
2686ed231c7SMatthew Knepley     if (osm->sort_indices) {
26992bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2704b9ad928SBarry Smith         ierr = ISSort(osm->is[i]);CHKERRQ(ierr);
2712b691e39SMatthew Knepley         if (osm->is_local) {
2722b691e39SMatthew Knepley           ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr);
2732b691e39SMatthew Knepley         }
2744b9ad928SBarry Smith       }
2756ed231c7SMatthew Knepley     }
2764b9ad928SBarry Smith 
277f6991133SBarry Smith     if (!osm->ksp) {
27878904715SLisandro Dalcin       /* Create the local solvers */
279785e854fSJed Brown       ierr = PetscMalloc1(osm->n_local_true,&osm->ksp);CHKERRQ(ierr);
280feb221f8SDmitry Karpeev       if (domain_dm) {
281feb221f8SDmitry Karpeev         ierr = PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n");CHKERRQ(ierr);
282feb221f8SDmitry Karpeev       }
28392bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2844b9ad928SBarry Smith         ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
285422a814eSBarry Smith         ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr);
2863bb1ff40SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr);
28792bb6962SLisandro Dalcin         ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
2884b9ad928SBarry Smith         ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
2894b9ad928SBarry Smith         ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
2904b9ad928SBarry Smith         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
2914b9ad928SBarry Smith         ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
2924b9ad928SBarry Smith         ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
293feb221f8SDmitry Karpeev         if (domain_dm) {
294feb221f8SDmitry Karpeev           ierr = KSPSetDM(ksp, domain_dm[i]);CHKERRQ(ierr);
295feb221f8SDmitry Karpeev           ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr);
296feb221f8SDmitry Karpeev           ierr = DMDestroy(&domain_dm[i]);CHKERRQ(ierr);
297feb221f8SDmitry Karpeev         }
2984b9ad928SBarry Smith         osm->ksp[i] = ksp;
2994b9ad928SBarry Smith       }
300feb221f8SDmitry Karpeev       if (domain_dm) {
301feb221f8SDmitry Karpeev         ierr = PetscFree(domain_dm);CHKERRQ(ierr);
302feb221f8SDmitry Karpeev       }
303f6991133SBarry Smith     }
3041dd8081eSeaulisa 
305fb745f2cSMatthew G. Knepley     ierr = ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis);CHKERRQ(ierr);
306fb745f2cSMatthew G. Knepley     ierr = ISSortRemoveDups(osm->lis);CHKERRQ(ierr);
307fb745f2cSMatthew G. Knepley     ierr = ISGetLocalSize(osm->lis, &m);CHKERRQ(ierr);
308fb745f2cSMatthew G. Knepley     ierr = VecCreateSeq(PETSC_COMM_SELF, m, &osm->lx);CHKERRQ(ierr);
309fb745f2cSMatthew G. Knepley     ierr = VecDuplicate(osm->lx, &osm->ly);CHKERRQ(ierr);
3101dd8081eSeaulisa 
3114b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
3124b9ad928SBarry Smith   } else {
3134b9ad928SBarry Smith     /*
3144b9ad928SBarry Smith        Destroy the blocks from the previous iteration
3154b9ad928SBarry Smith     */
316e09e08ccSBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
3174b9ad928SBarry Smith       ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
3184b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
3194b9ad928SBarry Smith     }
3204b9ad928SBarry Smith   }
3214b9ad928SBarry Smith 
32292bb6962SLisandro Dalcin   /*
32392bb6962SLisandro Dalcin      Extract out the submatrices
32492bb6962SLisandro Dalcin   */
3257dae84e0SHong Zhang   ierr = MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr);
32692bb6962SLisandro Dalcin   if (scall == MAT_INITIAL_MATRIX) {
32792bb6962SLisandro Dalcin     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
32892bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
3293bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr);
33078904715SLisandro Dalcin       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
33192bb6962SLisandro Dalcin     }
33292bb6962SLisandro Dalcin   }
33380ec0b7dSPatrick Sanan 
33480ec0b7dSPatrick Sanan   /* Convert the types of the submatrices (if needbe) */
33580ec0b7dSPatrick Sanan   if (osm->sub_mat_type) {
33680ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; i++) {
33780ec0b7dSPatrick Sanan       ierr = MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]));CHKERRQ(ierr);
33880ec0b7dSPatrick Sanan     }
33980ec0b7dSPatrick Sanan   }
34080ec0b7dSPatrick Sanan 
34180ec0b7dSPatrick Sanan   if(!pc->setupcalled){
34280ec0b7dSPatrick Sanan     /* Create the local work vectors (from the local matrices) and scatter contexts */
34380ec0b7dSPatrick Sanan     ierr = MatCreateVecs(pc->pmat,&vec,0);CHKERRQ(ierr);
3445b3c0d42Seaulisa 
3451dd8081eSeaulisa     if (osm->is_local && (osm->type == PC_ASM_INTERPOLATE || osm->type == PC_ASM_NONE )) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot use interpolate or none PCASMType if is_local was provided to PCASMSetLocalSubdomains()");
3461dd8081eSeaulisa     if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) {
3471dd8081eSeaulisa       ierr = PetscMalloc1(osm->n_local_true,&osm->lprolongation);CHKERRQ(ierr);
3481dd8081eSeaulisa     }
3491dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->lrestriction);CHKERRQ(ierr);
3501dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->x);CHKERRQ(ierr);
3511dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->y);CHKERRQ(ierr);
352347574c9Seaulisa 
353347574c9Seaulisa     ierr = ISGetLocalSize(osm->lis,&m);CHKERRQ(ierr);
354347574c9Seaulisa     ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
3551dd8081eSeaulisa     ierr = VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->restriction);CHKERRQ(ierr);
356347574c9Seaulisa     ierr = ISDestroy(&isl);CHKERRQ(ierr);
357347574c9Seaulisa 
358347574c9Seaulisa 
35980ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; ++i) {
36080ec0b7dSPatrick Sanan       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
36180ec0b7dSPatrick Sanan       ierr = MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);CHKERRQ(ierr);
36280ec0b7dSPatrick Sanan       ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr);
363d8b3b5e3Seaulisa 
3645b3c0d42Seaulisa       ISLocalToGlobalMapping ltog;
3655b3c0d42Seaulisa       IS                     isll;
3665b3c0d42Seaulisa       const PetscInt         *idx_is;
3675b3c0d42Seaulisa       PetscInt               *idx_lis,nout;
3685b3c0d42Seaulisa 
369b0de9f23SBarry Smith       /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */
3705b3c0d42Seaulisa       ierr = ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);CHKERRQ(ierr);
3715b3c0d42Seaulisa       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
3725b3c0d42Seaulisa       ierr = ISGetIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
3735b3c0d42Seaulisa       ierr = PetscMalloc1(m,&idx_lis);CHKERRQ(ierr);
3745b3c0d42Seaulisa       ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis);CHKERRQ(ierr);
3755b3c0d42Seaulisa       if (nout != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis");
3765b3c0d42Seaulisa       ierr = ISRestoreIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
3775b3c0d42Seaulisa       ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
378d8b3b5e3Seaulisa       ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
3795b3c0d42Seaulisa       ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
3801dd8081eSeaulisa       ierr = VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lrestriction[i]);CHKERRQ(ierr);
3815b3c0d42Seaulisa       ierr = ISDestroy(&isll);CHKERRQ(ierr);
3825b3c0d42Seaulisa       ierr = ISDestroy(&isl);CHKERRQ(ierr);
383b0de9f23SBarry Smith       if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overalapping is_local[i] entries */
384d8b3b5e3Seaulisa 	ISLocalToGlobalMapping ltog;
385d8b3b5e3Seaulisa         IS                     isll,isll_local;
386d8b3b5e3Seaulisa         const PetscInt         *idx_local;
387d8b3b5e3Seaulisa         PetscInt               *idx1, *idx2, nout;
388d8b3b5e3Seaulisa 
389d8b3b5e3Seaulisa         ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr);
390d8b3b5e3Seaulisa         ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
391d8b3b5e3Seaulisa 
392d8b3b5e3Seaulisa 	ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],&ltog);CHKERRQ(ierr);
393d8b3b5e3Seaulisa 	ierr = PetscMalloc1(m_local,&idx1);CHKERRQ(ierr);
394d8b3b5e3Seaulisa 	ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1);CHKERRQ(ierr);
395d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
396d8b3b5e3Seaulisa         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is");
397d8b3b5e3Seaulisa         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
398d8b3b5e3Seaulisa 
399d8b3b5e3Seaulisa 	ierr = ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);CHKERRQ(ierr);
400d8b3b5e3Seaulisa 	ierr = PetscMalloc1(m_local,&idx2);CHKERRQ(ierr);
401d8b3b5e3Seaulisa       	ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2);CHKERRQ(ierr);
402d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
403d8b3b5e3Seaulisa         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis");
404d8b3b5e3Seaulisa         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local);CHKERRQ(ierr);
405d8b3b5e3Seaulisa 
406d8b3b5e3Seaulisa 	ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
4071dd8081eSeaulisa 	ierr = VecScatterCreate(osm->y[i],isll,osm->ly,isll_local,&osm->lprolongation[i]);CHKERRQ(ierr);
408d8b3b5e3Seaulisa 
409d8b3b5e3Seaulisa         ierr = ISDestroy(&isll);CHKERRQ(ierr);
410d8b3b5e3Seaulisa 	ierr = ISDestroy(&isll_local);CHKERRQ(ierr);
411d8b3b5e3Seaulisa       }
41280ec0b7dSPatrick Sanan     }
41380ec0b7dSPatrick Sanan     ierr = VecDestroy(&vec);CHKERRQ(ierr);
41480ec0b7dSPatrick Sanan   }
41580ec0b7dSPatrick Sanan 
416fb745f2cSMatthew G. Knepley   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
417235cc4ceSMatthew G. Knepley     IS      *cis;
418235cc4ceSMatthew G. Knepley     PetscInt c;
419235cc4ceSMatthew G. Knepley 
420235cc4ceSMatthew G. Knepley     ierr = PetscMalloc1(osm->n_local_true, &cis);CHKERRQ(ierr);
421235cc4ceSMatthew G. Knepley     for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis;
4227dae84e0SHong Zhang     ierr = MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);CHKERRQ(ierr);
423235cc4ceSMatthew G. Knepley     ierr = PetscFree(cis);CHKERRQ(ierr);
424fb745f2cSMatthew G. Knepley   }
4254b9ad928SBarry Smith 
4264b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
4274b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
42892bb6962SLisandro Dalcin   ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
4294b9ad928SBarry Smith 
43092bb6962SLisandro Dalcin   /*
43192bb6962SLisandro Dalcin      Loop over subdomains putting them into local ksp
43292bb6962SLisandro Dalcin   */
43392bb6962SLisandro Dalcin   for (i=0; i<osm->n_local_true; i++) {
43423ee1639SBarry Smith     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr);
43592bb6962SLisandro Dalcin     if (!pc->setupcalled) {
436bf108f30SBarry Smith       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
4374b9ad928SBarry Smith     }
43892bb6962SLisandro Dalcin   }
4394b9ad928SBarry Smith   PetscFunctionReturn(0);
4404b9ad928SBarry Smith }
4414b9ad928SBarry Smith 
4426849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc)
4434b9ad928SBarry Smith {
4444b9ad928SBarry Smith   PC_ASM             *osm = (PC_ASM*)pc->data;
4456849ba73SBarry Smith   PetscErrorCode     ierr;
44613f74950SBarry Smith   PetscInt           i;
447539c167fSBarry Smith   KSPConvergedReason reason;
4484b9ad928SBarry Smith 
4494b9ad928SBarry Smith   PetscFunctionBegin;
4504b9ad928SBarry Smith   for (i=0; i<osm->n_local_true; i++) {
4514b9ad928SBarry Smith     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
452539c167fSBarry Smith     ierr = KSPGetConvergedReason(osm->ksp[i],&reason);CHKERRQ(ierr);
453539c167fSBarry Smith     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
454261222b2SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
455e0eafd54SHong Zhang     }
4564b9ad928SBarry Smith   }
4574b9ad928SBarry Smith   PetscFunctionReturn(0);
4584b9ad928SBarry Smith }
4594b9ad928SBarry Smith 
4606849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y)
4614b9ad928SBarry Smith {
4624b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
4636849ba73SBarry Smith   PetscErrorCode ierr;
4641dd8081eSeaulisa   PetscInt       i,n_local_true = osm->n_local_true;
4654b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
4664b9ad928SBarry Smith 
4674b9ad928SBarry Smith   PetscFunctionBegin;
4684b9ad928SBarry Smith   /*
4694b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
4704b9ad928SBarry Smith      subdomain values (leaving the other values 0).
4714b9ad928SBarry Smith   */
4724b9ad928SBarry Smith   if (!(osm->type & PC_ASM_RESTRICT)) {
4734b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
4744b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
4751dd8081eSeaulisa     ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr);
4764b9ad928SBarry Smith   }
477347574c9Seaulisa   if (!(osm->type & PC_ASM_INTERPOLATE)) {
478347574c9Seaulisa     reverse = SCATTER_REVERSE_LOCAL;
479347574c9Seaulisa   }
4804b9ad928SBarry Smith 
481347574c9Seaulisa   if(osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE){
482b0de9f23SBarry Smith     /* zero the global and the local solutions */
48312cd4985SMatthew G. Knepley     ierr = VecZeroEntries(y);CHKERRQ(ierr);
4845b3c0d42Seaulisa     ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
485347574c9Seaulisa 
486b0de9f23SBarry Smith     /* Copy the global RHS to local RHS including the ghost nodes */
4871dd8081eSeaulisa     ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
4881dd8081eSeaulisa     ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
489347574c9Seaulisa 
490b0de9f23SBarry Smith     /* Restrict local RHS to the overlapping 0-block RHS */
4911dd8081eSeaulisa     ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
4921dd8081eSeaulisa     ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0],  INSERT_VALUES, forward);CHKERRQ(ierr);
493d8b3b5e3Seaulisa 
49412cd4985SMatthew G. Knepley     /* do the local solves */
49512cd4985SMatthew G. Knepley     for (i = 0; i < n_local_true; ++i) {
496347574c9Seaulisa 
497b0de9f23SBarry Smith       /* solve the overlapping i-block */
49812cd4985SMatthew G. Knepley       ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
499d8b3b5e3Seaulisa 
500b0de9f23SBarry Smith       if (osm->lprolongation) { /* interpolate the non-overalapping i-block solution to the local solution (only for restrictive additive) */
5011dd8081eSeaulisa  	ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
5021dd8081eSeaulisa  	ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
503d8b3b5e3Seaulisa       }
504b0de9f23SBarry Smith       else{ /* interpolate the overalapping i-block solution to the local solution */
5051dd8081eSeaulisa         ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
5061dd8081eSeaulisa 	ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
507d8b3b5e3Seaulisa       }
508347574c9Seaulisa 
509347574c9Seaulisa       if (i < n_local_true-1) {
510b0de9f23SBarry Smith 	/* Restrict local RHS to the overlapping (i+1)-block RHS */
5111dd8081eSeaulisa 	ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
5121dd8081eSeaulisa 	ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
513347574c9Seaulisa 
514347574c9Seaulisa 	if ( osm->loctype == PC_COMPOSITE_MULTIPLICATIVE){
515b0de9f23SBarry Smith 	  /* udpdate the overlapping (i+1)-block RHS using the current local solution */
516347574c9Seaulisa 	  ierr = MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]);CHKERRQ(ierr);
517347574c9Seaulisa 	  ierr = VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1]); CHKERRQ(ierr);
5187c3d802fSMatthew G. Knepley 	}
51912cd4985SMatthew G. Knepley       }
52012cd4985SMatthew G. Knepley     }
521b0de9f23SBarry Smith     /* Add the local solution to the global solution including the ghost nodes */
5221dd8081eSeaulisa     ierr = VecScatterBegin(osm->restriction, osm->ly, y,  ADD_VALUES, reverse);CHKERRQ(ierr);
5231dd8081eSeaulisa     ierr = VecScatterEnd(osm->restriction,  osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
524347574c9Seaulisa   }else{
525347574c9Seaulisa     SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
52612cd4985SMatthew G. Knepley   }
5274b9ad928SBarry Smith   PetscFunctionReturn(0);
5284b9ad928SBarry Smith }
5294b9ad928SBarry Smith 
5306849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y)
5314b9ad928SBarry Smith {
5324b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
5336849ba73SBarry Smith   PetscErrorCode ierr;
5341dd8081eSeaulisa   PetscInt       i,n_local_true = osm->n_local_true;
5354b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
5364b9ad928SBarry Smith 
5374b9ad928SBarry Smith   PetscFunctionBegin;
5384b9ad928SBarry Smith   /*
5394b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
5404b9ad928SBarry Smith      subdomain values (leaving the other values 0).
5414b9ad928SBarry Smith 
5424b9ad928SBarry Smith      Note: these are reversed from the PCApply_ASM() because we are applying the
5434b9ad928SBarry Smith      transpose of the three terms
5444b9ad928SBarry Smith   */
545d8b3b5e3Seaulisa 
5464b9ad928SBarry Smith   if (!(osm->type & PC_ASM_INTERPOLATE)) {
5474b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
5484b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
5491dd8081eSeaulisa     ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr);
5504b9ad928SBarry Smith   }
5512fa5cd67SKarl Rupp   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;
5524b9ad928SBarry Smith 
553b0de9f23SBarry Smith   /* zero the global and the local solutions */
55410bd88b9SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
555d8b3b5e3Seaulisa   ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
556d8b3b5e3Seaulisa 
557b0de9f23SBarry Smith   /* Copy the global RHS to local RHS including the ghost nodes */
5581dd8081eSeaulisa   ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
5591dd8081eSeaulisa   ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
560d8b3b5e3Seaulisa 
561b0de9f23SBarry Smith   /* Restrict local RHS to the overlapping 0-block RHS */
5621dd8081eSeaulisa   ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
5631dd8081eSeaulisa   ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0],  INSERT_VALUES, forward);CHKERRQ(ierr);
564d8b3b5e3Seaulisa 
5654b9ad928SBarry Smith   /* do the local solves */
566d8b3b5e3Seaulisa   for (i = 0; i < n_local_true; ++i) {
567d8b3b5e3Seaulisa 
568b0de9f23SBarry Smith     /* solve the overlapping i-block */
569d8b3b5e3Seaulisa     ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
570d8b3b5e3Seaulisa 
571b0de9f23SBarry Smith     if (osm->lprolongation) { /* interpolate the non-overalapping i-block solution to the local solution */
5721dd8081eSeaulisa      ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
5731dd8081eSeaulisa      ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
574b9845e0eSMatthew Knepley     }
575b0de9f23SBarry Smith     else{ /* interpolate the overalapping i-block solution to the local solution */
5761dd8081eSeaulisa       ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
5771dd8081eSeaulisa       ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
5784b9ad928SBarry Smith     }
579d8b3b5e3Seaulisa 
580d8b3b5e3Seaulisa     if (i < n_local_true-1) {
581b0de9f23SBarry Smith       /* Restrict local RHS to the overlapping (i+1)-block RHS */
5821dd8081eSeaulisa       ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
5831dd8081eSeaulisa       ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
5844b9ad928SBarry Smith     }
5854b9ad928SBarry Smith   }
586b0de9f23SBarry Smith   /* Add the local solution to the global solution including the ghost nodes */
5871dd8081eSeaulisa   ierr = VecScatterBegin(osm->restriction, osm->ly, y,  ADD_VALUES, reverse);CHKERRQ(ierr);
5881dd8081eSeaulisa   ierr = VecScatterEnd(osm->restriction,  osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
589d8b3b5e3Seaulisa 
5904b9ad928SBarry Smith   PetscFunctionReturn(0);
591d8b3b5e3Seaulisa 
5924b9ad928SBarry Smith }
5934b9ad928SBarry Smith 
594e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc)
5954b9ad928SBarry Smith {
5964b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
5976849ba73SBarry Smith   PetscErrorCode ierr;
59813f74950SBarry Smith   PetscInt       i;
5994b9ad928SBarry Smith 
6004b9ad928SBarry Smith   PetscFunctionBegin;
60192bb6962SLisandro Dalcin   if (osm->ksp) {
60292bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
603e91c6855SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
60492bb6962SLisandro Dalcin     }
60592bb6962SLisandro Dalcin   }
606e09e08ccSBarry Smith   if (osm->pmat) {
60792bb6962SLisandro Dalcin     if (osm->n_local_true > 0) {
60830a70a9aSHong Zhang       ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
60992bb6962SLisandro Dalcin     }
61092bb6962SLisandro Dalcin   }
6111dd8081eSeaulisa   if (osm->lrestriction) {
6121dd8081eSeaulisa     ierr = VecScatterDestroy(&osm->restriction);CHKERRQ(ierr);
6131dd8081eSeaulisa     for (i=0; i<osm->n_local_true; i++) {
6141dd8081eSeaulisa       ierr = VecScatterDestroy(&osm->lrestriction[i]);CHKERRQ(ierr);
6151dd8081eSeaulisa       if (osm->lprolongation) {ierr = VecScatterDestroy(&osm->lprolongation[i]);CHKERRQ(ierr);}
616fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
617fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
6184b9ad928SBarry Smith     }
6191dd8081eSeaulisa     ierr = PetscFree(osm->lrestriction);CHKERRQ(ierr);
6201dd8081eSeaulisa     if (osm->lprolongation) {ierr = PetscFree(osm->lprolongation);CHKERRQ(ierr);}
62105b42c5fSBarry Smith     ierr = PetscFree(osm->x);CHKERRQ(ierr);
62278904715SLisandro Dalcin     ierr = PetscFree(osm->y);CHKERRQ(ierr);
6231dd8081eSeaulisa 
62492bb6962SLisandro Dalcin   }
6252b691e39SMatthew Knepley   ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
626fb745f2cSMatthew G. Knepley   ierr = ISDestroy(&osm->lis);CHKERRQ(ierr);
627fb745f2cSMatthew G. Knepley   ierr = VecDestroy(&osm->lx);CHKERRQ(ierr);
628fb745f2cSMatthew G. Knepley   ierr = VecDestroy(&osm->ly);CHKERRQ(ierr);
629347574c9Seaulisa   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
630347574c9Seaulisa     ierr = MatDestroyMatrices(osm->n_local_true, &osm->lmats);CHKERRQ(ierr);
631fb745f2cSMatthew G. Knepley   }
6322fa5cd67SKarl Rupp 
63380ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
63480ec0b7dSPatrick Sanan 
635e91c6855SBarry Smith   osm->is       = 0;
636e91c6855SBarry Smith   osm->is_local = 0;
637e91c6855SBarry Smith   PetscFunctionReturn(0);
638e91c6855SBarry Smith }
639e91c6855SBarry Smith 
640e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc)
641e91c6855SBarry Smith {
642e91c6855SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
643e91c6855SBarry Smith   PetscErrorCode ierr;
644e91c6855SBarry Smith   PetscInt       i;
645e91c6855SBarry Smith 
646e91c6855SBarry Smith   PetscFunctionBegin;
647e91c6855SBarry Smith   ierr = PCReset_ASM(pc);CHKERRQ(ierr);
648e91c6855SBarry Smith   if (osm->ksp) {
649e91c6855SBarry Smith     for (i=0; i<osm->n_local_true; i++) {
650fcfd50ebSBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
651e91c6855SBarry Smith     }
652e91c6855SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
653e91c6855SBarry Smith   }
654e91c6855SBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
6554b9ad928SBarry Smith   PetscFunctionReturn(0);
6564b9ad928SBarry Smith }
6574b9ad928SBarry Smith 
6584416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc)
6594b9ad928SBarry Smith {
6604b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
6616849ba73SBarry Smith   PetscErrorCode ierr;
6629dcbbd2bSBarry Smith   PetscInt       blocks,ovl;
663ace3abfcSBarry Smith   PetscBool      symset,flg;
66492bb6962SLisandro Dalcin   PCASMType      asmtype;
66512cd4985SMatthew G. Knepley   PCCompositeType loctype;
66680ec0b7dSPatrick Sanan   char           sub_mat_type[256];
6674b9ad928SBarry Smith 
6684b9ad928SBarry Smith   PetscFunctionBegin;
669bf108f30SBarry Smith   /* set the type to symmetric if matrix is symmetric */
67092bb6962SLisandro Dalcin   if (!osm->type_set && pc->pmat) {
67192bb6962SLisandro Dalcin     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
6722fa5cd67SKarl Rupp     if (symset && flg) osm->type = PC_ASM_BASIC;
673bf108f30SBarry Smith   }
674e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");CHKERRQ(ierr);
675d709ea83SDmitry Karpeev   ierr = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr);
67690d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
67765db9045SDmitry Karpeev   if (flg) {
678f77a5242SKarl Rupp     ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr);
679d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
68065db9045SDmitry Karpeev   }
68190d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
68265db9045SDmitry Karpeev   if (flg) {
68365db9045SDmitry Karpeev     ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr);
684d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
68565db9045SDmitry Karpeev   }
68690d69ab7SBarry Smith   flg  = PETSC_FALSE;
68790d69ab7SBarry Smith   ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
68892bb6962SLisandro Dalcin   if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); }
68912cd4985SMatthew G. Knepley   flg  = PETSC_FALSE;
69012cd4985SMatthew G. Knepley   ierr = PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr);
69112cd4985SMatthew G. Knepley   if (flg) {ierr = PCASMSetLocalType(pc,loctype);CHKERRQ(ierr); }
69280ec0b7dSPatrick Sanan   ierr = PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);CHKERRQ(ierr);
69380ec0b7dSPatrick Sanan   if(flg){
694459726d8SSatish Balay     ierr = PCASMSetSubMatType(pc,sub_mat_type);CHKERRQ(ierr);
69580ec0b7dSPatrick Sanan   }
6964b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
6974b9ad928SBarry Smith   PetscFunctionReturn(0);
6984b9ad928SBarry Smith }
6994b9ad928SBarry Smith 
7004b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/
7014b9ad928SBarry Smith 
7021e6b0712SBarry Smith static PetscErrorCode  PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[])
7034b9ad928SBarry Smith {
7044b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
70592bb6962SLisandro Dalcin   PetscErrorCode ierr;
70692bb6962SLisandro Dalcin   PetscInt       i;
7074b9ad928SBarry Smith 
7084b9ad928SBarry Smith   PetscFunctionBegin;
709e32f2f54SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
710ce94432eSBarry Smith   if (pc->setupcalled && (n != osm->n_local_true || is)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetLocalSubdomains() should be called before calling PCSetUp().");
711e7e72b3dSBarry Smith 
7124b9ad928SBarry Smith   if (!pc->setupcalled) {
71392bb6962SLisandro Dalcin     if (is) {
71492bb6962SLisandro Dalcin       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
71592bb6962SLisandro Dalcin     }
716832fc9a5SMatthew Knepley     if (is_local) {
717832fc9a5SMatthew Knepley       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
718832fc9a5SMatthew Knepley     }
7192b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
7202fa5cd67SKarl Rupp 
7214b9ad928SBarry Smith     osm->n_local_true = n;
72292bb6962SLisandro Dalcin     osm->is           = 0;
7232b691e39SMatthew Knepley     osm->is_local     = 0;
72492bb6962SLisandro Dalcin     if (is) {
725785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is);CHKERRQ(ierr);
7262fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is[i] = is[i];
7273d03bb48SJed Brown       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
7283d03bb48SJed Brown       osm->overlap = -1;
72992bb6962SLisandro Dalcin     }
7302b691e39SMatthew Knepley     if (is_local) {
731785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is_local);CHKERRQ(ierr);
7322fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is_local[i] = is_local[i];
733a35b7b57SDmitry Karpeev       if (!is) {
734785e854fSJed Brown         ierr = PetscMalloc1(osm->n_local_true,&osm->is);CHKERRQ(ierr);
735a35b7b57SDmitry Karpeev         for (i=0; i<osm->n_local_true; i++) {
736a35b7b57SDmitry Karpeev           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
737a35b7b57SDmitry Karpeev             ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr);
738a35b7b57SDmitry Karpeev             ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr);
739a35b7b57SDmitry Karpeev           } else {
740a35b7b57SDmitry Karpeev             ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr);
741a35b7b57SDmitry Karpeev             osm->is[i] = osm->is_local[i];
742a35b7b57SDmitry Karpeev           }
743a35b7b57SDmitry Karpeev         }
744a35b7b57SDmitry Karpeev       }
7452b691e39SMatthew Knepley     }
7464b9ad928SBarry Smith   }
7474b9ad928SBarry Smith   PetscFunctionReturn(0);
7484b9ad928SBarry Smith }
7494b9ad928SBarry Smith 
7501e6b0712SBarry Smith static PetscErrorCode  PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local)
7514b9ad928SBarry Smith {
7524b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
7536849ba73SBarry Smith   PetscErrorCode ierr;
75413f74950SBarry Smith   PetscMPIInt    rank,size;
75578904715SLisandro Dalcin   PetscInt       n;
7564b9ad928SBarry Smith 
7574b9ad928SBarry Smith   PetscFunctionBegin;
758ce94432eSBarry Smith   if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
759ce94432eSBarry Smith   if (is || is_local) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Use PCASMSetLocalSubdomains() to set specific index sets\n\they cannot be set globally yet.");
7604b9ad928SBarry Smith 
7614b9ad928SBarry Smith   /*
762880770e9SJed Brown      Split the subdomains equally among all processors
7634b9ad928SBarry Smith   */
764ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
765ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
7664b9ad928SBarry Smith   n    = N/size + ((N % size) > rank);
767e32f2f54SBarry Smith   if (!n) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Process %d must have at least one block: total processors %d total blocks %D",(int)rank,(int)size,N);
768e32f2f54SBarry Smith   if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp().");
7694b9ad928SBarry Smith   if (!pc->setupcalled) {
7702b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
7712fa5cd67SKarl Rupp 
7724b9ad928SBarry Smith     osm->n_local_true = n;
7734b9ad928SBarry Smith     osm->is           = 0;
7742b691e39SMatthew Knepley     osm->is_local     = 0;
7754b9ad928SBarry Smith   }
7764b9ad928SBarry Smith   PetscFunctionReturn(0);
7774b9ad928SBarry Smith }
7784b9ad928SBarry Smith 
7791e6b0712SBarry Smith static PetscErrorCode  PCASMSetOverlap_ASM(PC pc,PetscInt ovl)
7804b9ad928SBarry Smith {
78192bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
7824b9ad928SBarry Smith 
7834b9ad928SBarry Smith   PetscFunctionBegin;
784ce94432eSBarry Smith   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
785ce94432eSBarry Smith   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp().");
7862fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
7874b9ad928SBarry Smith   PetscFunctionReturn(0);
7884b9ad928SBarry Smith }
7894b9ad928SBarry Smith 
7901e6b0712SBarry Smith static PetscErrorCode  PCASMSetType_ASM(PC pc,PCASMType type)
7914b9ad928SBarry Smith {
79292bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
7934b9ad928SBarry Smith 
7944b9ad928SBarry Smith   PetscFunctionBegin;
7954b9ad928SBarry Smith   osm->type     = type;
796bf108f30SBarry Smith   osm->type_set = PETSC_TRUE;
7974b9ad928SBarry Smith   PetscFunctionReturn(0);
7984b9ad928SBarry Smith }
7994b9ad928SBarry Smith 
800c60c7ad4SBarry Smith static PetscErrorCode  PCASMGetType_ASM(PC pc,PCASMType *type)
801c60c7ad4SBarry Smith {
802c60c7ad4SBarry Smith   PC_ASM *osm = (PC_ASM*)pc->data;
803c60c7ad4SBarry Smith 
804c60c7ad4SBarry Smith   PetscFunctionBegin;
805c60c7ad4SBarry Smith   *type = osm->type;
806c60c7ad4SBarry Smith   PetscFunctionReturn(0);
807c60c7ad4SBarry Smith }
808c60c7ad4SBarry Smith 
80912cd4985SMatthew G. Knepley static PetscErrorCode  PCASMSetLocalType_ASM(PC pc, PCCompositeType type)
81012cd4985SMatthew G. Knepley {
81112cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
81212cd4985SMatthew G. Knepley 
81312cd4985SMatthew G. Knepley   PetscFunctionBegin;
814d0ecd4dfSBarry Smith   if (type != PC_COMPOSITE_ADDITIVE && type != PC_COMPOSITE_MULTIPLICATIVE) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Only supports additive or multiplicative as the local type");
81512cd4985SMatthew G. Knepley   osm->loctype = type;
81612cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
81712cd4985SMatthew G. Knepley }
81812cd4985SMatthew G. Knepley 
81912cd4985SMatthew G. Knepley static PetscErrorCode  PCASMGetLocalType_ASM(PC pc, PCCompositeType *type)
82012cd4985SMatthew G. Knepley {
82112cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
82212cd4985SMatthew G. Knepley 
82312cd4985SMatthew G. Knepley   PetscFunctionBegin;
82412cd4985SMatthew G. Knepley   *type = osm->loctype;
82512cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
82612cd4985SMatthew G. Knepley }
82712cd4985SMatthew G. Knepley 
8281e6b0712SBarry Smith static PetscErrorCode  PCASMSetSortIndices_ASM(PC pc,PetscBool  doSort)
8296ed231c7SMatthew Knepley {
8306ed231c7SMatthew Knepley   PC_ASM *osm = (PC_ASM*)pc->data;
8316ed231c7SMatthew Knepley 
8326ed231c7SMatthew Knepley   PetscFunctionBegin;
8336ed231c7SMatthew Knepley   osm->sort_indices = doSort;
8346ed231c7SMatthew Knepley   PetscFunctionReturn(0);
8356ed231c7SMatthew Knepley }
8366ed231c7SMatthew Knepley 
8371e6b0712SBarry Smith static PetscErrorCode  PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
8384b9ad928SBarry Smith {
83992bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
840dfbe8321SBarry Smith   PetscErrorCode ierr;
8414b9ad928SBarry Smith 
8424b9ad928SBarry Smith   PetscFunctionBegin;
843ce94432eSBarry Smith   if (osm->n_local_true < 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here");
8444b9ad928SBarry Smith 
8452fa5cd67SKarl Rupp   if (n_local) *n_local = osm->n_local_true;
84692bb6962SLisandro Dalcin   if (first_local) {
847ce94432eSBarry Smith     ierr          = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
84892bb6962SLisandro Dalcin     *first_local -= osm->n_local_true;
84992bb6962SLisandro Dalcin   }
85092bb6962SLisandro Dalcin   if (ksp) {
85192bb6962SLisandro Dalcin     /* Assume that local solves are now different; not necessarily
85292bb6962SLisandro Dalcin        true though!  This flag is used only for PCView_ASM() */
85392bb6962SLisandro Dalcin     *ksp                   = osm->ksp;
85492bb6962SLisandro Dalcin     osm->same_local_solves = PETSC_FALSE;
85592bb6962SLisandro Dalcin   }
8564b9ad928SBarry Smith   PetscFunctionReturn(0);
8574b9ad928SBarry Smith }
8584b9ad928SBarry Smith 
85980ec0b7dSPatrick Sanan static PetscErrorCode  PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type)
86080ec0b7dSPatrick Sanan {
86180ec0b7dSPatrick Sanan   PC_ASM         *osm = (PC_ASM*)pc->data;
86280ec0b7dSPatrick Sanan 
86380ec0b7dSPatrick Sanan   PetscFunctionBegin;
86480ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
86580ec0b7dSPatrick Sanan   PetscValidPointer(sub_mat_type,2);
86680ec0b7dSPatrick Sanan   *sub_mat_type = osm->sub_mat_type;
86780ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
86880ec0b7dSPatrick Sanan }
86980ec0b7dSPatrick Sanan 
87080ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type)
87180ec0b7dSPatrick Sanan {
87280ec0b7dSPatrick Sanan   PetscErrorCode    ierr;
87380ec0b7dSPatrick Sanan   PC_ASM            *osm = (PC_ASM*)pc->data;
87480ec0b7dSPatrick Sanan 
87580ec0b7dSPatrick Sanan   PetscFunctionBegin;
87680ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
87780ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
87880ec0b7dSPatrick Sanan   ierr = PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);CHKERRQ(ierr);
87980ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
88080ec0b7dSPatrick Sanan }
88180ec0b7dSPatrick Sanan 
8824b9ad928SBarry Smith /*@C
8831093a601SBarry Smith     PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner.
8844b9ad928SBarry Smith 
8854b9ad928SBarry Smith     Collective on PC
8864b9ad928SBarry Smith 
8874b9ad928SBarry Smith     Input Parameters:
8884b9ad928SBarry Smith +   pc - the preconditioner context
8894b9ad928SBarry Smith .   n - the number of subdomains for this processor (default value = 1)
8908c03b21aSDmitry Karpeev .   is - the index set that defines the subdomains for this processor
8910298fd71SBarry Smith          (or NULL for PETSc to determine subdomains)
892f1ee410cSBarry Smith -   is_local - the index sets that define the local part of the subdomains for this processor, not used unless PCASMType is PC_ASM_RESTRICT
893f1ee410cSBarry Smith          (or NULL to not provide these)
8944b9ad928SBarry Smith 
8954b9ad928SBarry Smith     Notes:
8961093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
8974b9ad928SBarry Smith 
8984b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
8994b9ad928SBarry Smith 
9004b9ad928SBarry Smith     Use PCASMSetTotalSubdomains() to set the subdomains for all processors.
9014b9ad928SBarry Smith 
902f1ee410cSBarry Smith     If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated
903f1ee410cSBarry Smith     back to form the global solution (this is the standard restricted additive Schwarz method)
904f1ee410cSBarry Smith 
905f1ee410cSBarry Smith     If the is_local is provided and PCASMType is PC_ASM_INTERPOLATE or PC_ASM_NONE then an error is generated since there is
906f1ee410cSBarry Smith     no code to handle that case.
907f1ee410cSBarry Smith 
9084b9ad928SBarry Smith     Level: advanced
9094b9ad928SBarry Smith 
9104b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz
9114b9ad928SBarry Smith 
9124b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
913f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType()
9144b9ad928SBarry Smith @*/
9157087cfbeSBarry Smith PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
9164b9ad928SBarry Smith {
9177bb14e67SBarry Smith   PetscErrorCode ierr;
9184b9ad928SBarry Smith 
9194b9ad928SBarry Smith   PetscFunctionBegin;
9200700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9217bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr);
9224b9ad928SBarry Smith   PetscFunctionReturn(0);
9234b9ad928SBarry Smith }
9244b9ad928SBarry Smith 
9254b9ad928SBarry Smith /*@C
926feb221f8SDmitry Karpeev     PCASMSetTotalSubdomains - Sets the subdomains for all processors for the
9274b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
9284b9ad928SBarry Smith     PC communicator must call this routine, with the same index sets.
9294b9ad928SBarry Smith 
9304b9ad928SBarry Smith     Collective on PC
9314b9ad928SBarry Smith 
9324b9ad928SBarry Smith     Input Parameters:
9334b9ad928SBarry Smith +   pc - the preconditioner context
934feb221f8SDmitry Karpeev .   N  - the number of subdomains for all processors
935feb221f8SDmitry Karpeev .   is - the index sets that define the subdomains for all processors
936dfaa0c5fSBarry Smith          (or NULL to ask PETSc to determine the subdomains)
9372b691e39SMatthew Knepley -   is_local - the index sets that define the local part of the subdomains for this processor
938f1ee410cSBarry Smith          (or NULL to not provide this information)
9394b9ad928SBarry Smith 
9404b9ad928SBarry Smith     Options Database Key:
9414b9ad928SBarry Smith     To set the total number of subdomain blocks rather than specify the
9424b9ad928SBarry Smith     index sets, use the option
9434b9ad928SBarry Smith .    -pc_asm_blocks <blks> - Sets total blocks
9444b9ad928SBarry Smith 
9454b9ad928SBarry Smith     Notes:
946f1ee410cSBarry Smith     Currently you cannot use this to set the actual subdomains with the argument is or is_local.
9474b9ad928SBarry Smith 
9484b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
9494b9ad928SBarry Smith 
9504b9ad928SBarry Smith     These index sets cannot be destroyed until after completion of the
9514b9ad928SBarry Smith     linear solves for which the ASM preconditioner is being used.
9524b9ad928SBarry Smith 
9534b9ad928SBarry Smith     Use PCASMSetLocalSubdomains() to set local subdomains.
9544b9ad928SBarry Smith 
9551093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
9561093a601SBarry Smith 
9574b9ad928SBarry Smith     Level: advanced
9584b9ad928SBarry Smith 
9594b9ad928SBarry Smith .keywords: PC, ASM, set, total, global, subdomains, additive Schwarz
9604b9ad928SBarry Smith 
9614b9ad928SBarry Smith .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
9624b9ad928SBarry Smith           PCASMCreateSubdomains2D()
9634b9ad928SBarry Smith @*/
9647087cfbeSBarry Smith PetscErrorCode  PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[])
9654b9ad928SBarry Smith {
9667bb14e67SBarry Smith   PetscErrorCode ierr;
9674b9ad928SBarry Smith 
9684b9ad928SBarry Smith   PetscFunctionBegin;
9690700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9707bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr);
9714b9ad928SBarry Smith   PetscFunctionReturn(0);
9724b9ad928SBarry Smith }
9734b9ad928SBarry Smith 
9744b9ad928SBarry Smith /*@
9754b9ad928SBarry Smith     PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
9764b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
977f1ee410cSBarry Smith     PC communicator must call this routine.
9784b9ad928SBarry Smith 
979ad4df100SBarry Smith     Logically Collective on PC
9804b9ad928SBarry Smith 
9814b9ad928SBarry Smith     Input Parameters:
9824b9ad928SBarry Smith +   pc  - the preconditioner context
9834b9ad928SBarry Smith -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
9844b9ad928SBarry Smith 
9854b9ad928SBarry Smith     Options Database Key:
9864b9ad928SBarry Smith .   -pc_asm_overlap <ovl> - Sets overlap
9874b9ad928SBarry Smith 
9884b9ad928SBarry Smith     Notes:
9894b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.  To use
9904b9ad928SBarry Smith     multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
9914b9ad928SBarry Smith     PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).
9924b9ad928SBarry Smith 
9934b9ad928SBarry Smith     The overlap defaults to 1, so if one desires that no additional
9944b9ad928SBarry Smith     overlap be computed beyond what may have been set with a call to
9954b9ad928SBarry Smith     PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
9964b9ad928SBarry Smith     must be set to be 0.  In particular, if one does not explicitly set
9974b9ad928SBarry Smith     the subdomains an application code, then all overlap would be computed
9984b9ad928SBarry Smith     internally by PETSc, and using an overlap of 0 would result in an ASM
9994b9ad928SBarry Smith     variant that is equivalent to the block Jacobi preconditioner.
10004b9ad928SBarry Smith 
1001f1ee410cSBarry Smith     The default algorithm used by PETSc to increase overlap is fast, but not scalable,
1002f1ee410cSBarry Smith     use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
1003f1ee410cSBarry Smith 
10044b9ad928SBarry Smith     Note that one can define initial index sets with any overlap via
1005f1ee410cSBarry Smith     PCASMSetLocalSubdomains(); the routine
10064b9ad928SBarry Smith     PCASMSetOverlap() merely allows PETSc to extend that overlap further
10074b9ad928SBarry Smith     if desired.
10084b9ad928SBarry Smith 
10094b9ad928SBarry Smith     Level: intermediate
10104b9ad928SBarry Smith 
10114b9ad928SBarry Smith .keywords: PC, ASM, set, overlap
10124b9ad928SBarry Smith 
10134b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
1014f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap()
10154b9ad928SBarry Smith @*/
10167087cfbeSBarry Smith PetscErrorCode  PCASMSetOverlap(PC pc,PetscInt ovl)
10174b9ad928SBarry Smith {
10187bb14e67SBarry Smith   PetscErrorCode ierr;
10194b9ad928SBarry Smith 
10204b9ad928SBarry Smith   PetscFunctionBegin;
10210700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1022c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,ovl,2);
10237bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
10244b9ad928SBarry Smith   PetscFunctionReturn(0);
10254b9ad928SBarry Smith }
10264b9ad928SBarry Smith 
10274b9ad928SBarry Smith /*@
10284b9ad928SBarry Smith     PCASMSetType - Sets the type of restriction and interpolation used
10294b9ad928SBarry Smith     for local problems in the additive Schwarz method.
10304b9ad928SBarry Smith 
1031ad4df100SBarry Smith     Logically Collective on PC
10324b9ad928SBarry Smith 
10334b9ad928SBarry Smith     Input Parameters:
10344b9ad928SBarry Smith +   pc  - the preconditioner context
10354b9ad928SBarry Smith -   type - variant of ASM, one of
10364b9ad928SBarry Smith .vb
10374b9ad928SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
10384b9ad928SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
10394b9ad928SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
10404b9ad928SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
10414b9ad928SBarry Smith .ve
10424b9ad928SBarry Smith 
10434b9ad928SBarry Smith     Options Database Key:
10444b9ad928SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
10454b9ad928SBarry Smith 
1046f1ee410cSBarry Smith     Notes: if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected
1047f1ee410cSBarry Smith     to limit the local processor interpolation
1048f1ee410cSBarry Smith 
10494b9ad928SBarry Smith     Level: intermediate
10504b9ad928SBarry Smith 
10514b9ad928SBarry Smith .keywords: PC, ASM, set, type
10524b9ad928SBarry Smith 
10534b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1054f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType()
10554b9ad928SBarry Smith @*/
10567087cfbeSBarry Smith PetscErrorCode  PCASMSetType(PC pc,PCASMType type)
10574b9ad928SBarry Smith {
10587bb14e67SBarry Smith   PetscErrorCode ierr;
10594b9ad928SBarry Smith 
10604b9ad928SBarry Smith   PetscFunctionBegin;
10610700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1062c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
10637bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr);
10644b9ad928SBarry Smith   PetscFunctionReturn(0);
10654b9ad928SBarry Smith }
10664b9ad928SBarry Smith 
1067c60c7ad4SBarry Smith /*@
1068c60c7ad4SBarry Smith     PCASMGetType - Gets the type of restriction and interpolation used
1069c60c7ad4SBarry Smith     for local problems in the additive Schwarz method.
1070c60c7ad4SBarry Smith 
1071c60c7ad4SBarry Smith     Logically Collective on PC
1072c60c7ad4SBarry Smith 
1073c60c7ad4SBarry Smith     Input Parameter:
1074c60c7ad4SBarry Smith .   pc  - the preconditioner context
1075c60c7ad4SBarry Smith 
1076c60c7ad4SBarry Smith     Output Parameter:
1077c60c7ad4SBarry Smith .   type - variant of ASM, one of
1078c60c7ad4SBarry Smith 
1079c60c7ad4SBarry Smith .vb
1080c60c7ad4SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
1081c60c7ad4SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1082c60c7ad4SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1083c60c7ad4SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
1084c60c7ad4SBarry Smith .ve
1085c60c7ad4SBarry Smith 
1086c60c7ad4SBarry Smith     Options Database Key:
1087c60c7ad4SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
1088c60c7ad4SBarry Smith 
1089c60c7ad4SBarry Smith     Level: intermediate
1090c60c7ad4SBarry Smith 
1091c60c7ad4SBarry Smith .keywords: PC, ASM, set, type
1092c60c7ad4SBarry Smith 
1093c60c7ad4SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1094f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType()
1095c60c7ad4SBarry Smith @*/
1096c60c7ad4SBarry Smith PetscErrorCode  PCASMGetType(PC pc,PCASMType *type)
1097c60c7ad4SBarry Smith {
1098c60c7ad4SBarry Smith   PetscErrorCode ierr;
1099c60c7ad4SBarry Smith 
1100c60c7ad4SBarry Smith   PetscFunctionBegin;
1101c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1102c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));CHKERRQ(ierr);
1103c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1104c60c7ad4SBarry Smith }
1105c60c7ad4SBarry Smith 
110612cd4985SMatthew G. Knepley /*@
110712cd4985SMatthew G. Knepley   PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method.
110812cd4985SMatthew G. Knepley 
110912cd4985SMatthew G. Knepley   Logically Collective on PC
111012cd4985SMatthew G. Knepley 
111112cd4985SMatthew G. Knepley   Input Parameters:
111212cd4985SMatthew G. Knepley + pc  - the preconditioner context
111312cd4985SMatthew G. Knepley - type - type of composition, one of
111412cd4985SMatthew G. Knepley .vb
111512cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
111612cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
111712cd4985SMatthew G. Knepley .ve
111812cd4985SMatthew G. Knepley 
111912cd4985SMatthew G. Knepley   Options Database Key:
112012cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
112112cd4985SMatthew G. Knepley 
112212cd4985SMatthew G. Knepley   Level: intermediate
112312cd4985SMatthew G. Knepley 
1124f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
112512cd4985SMatthew G. Knepley @*/
112612cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type)
112712cd4985SMatthew G. Knepley {
112812cd4985SMatthew G. Knepley   PetscErrorCode ierr;
112912cd4985SMatthew G. Knepley 
113012cd4985SMatthew G. Knepley   PetscFunctionBegin;
113112cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
113212cd4985SMatthew G. Knepley   PetscValidLogicalCollectiveEnum(pc, type, 2);
113312cd4985SMatthew G. Knepley   ierr = PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));CHKERRQ(ierr);
113412cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
113512cd4985SMatthew G. Knepley }
113612cd4985SMatthew G. Knepley 
113712cd4985SMatthew G. Knepley /*@
113812cd4985SMatthew G. Knepley   PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method.
113912cd4985SMatthew G. Knepley 
114012cd4985SMatthew G. Knepley   Logically Collective on PC
114112cd4985SMatthew G. Knepley 
114212cd4985SMatthew G. Knepley   Input Parameter:
114312cd4985SMatthew G. Knepley . pc  - the preconditioner context
114412cd4985SMatthew G. Knepley 
114512cd4985SMatthew G. Knepley   Output Parameter:
114612cd4985SMatthew G. Knepley . type - type of composition, one of
114712cd4985SMatthew G. Knepley .vb
114812cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
114912cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
115012cd4985SMatthew G. Knepley .ve
115112cd4985SMatthew G. Knepley 
115212cd4985SMatthew G. Knepley   Options Database Key:
115312cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
115412cd4985SMatthew G. Knepley 
115512cd4985SMatthew G. Knepley   Level: intermediate
115612cd4985SMatthew G. Knepley 
1157f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
115812cd4985SMatthew G. Knepley @*/
115912cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type)
116012cd4985SMatthew G. Knepley {
116112cd4985SMatthew G. Knepley   PetscErrorCode ierr;
116212cd4985SMatthew G. Knepley 
116312cd4985SMatthew G. Knepley   PetscFunctionBegin;
116412cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
116512cd4985SMatthew G. Knepley   PetscValidPointer(type, 2);
116612cd4985SMatthew G. Knepley   ierr = PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));CHKERRQ(ierr);
116712cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
116812cd4985SMatthew G. Knepley }
116912cd4985SMatthew G. Knepley 
11706ed231c7SMatthew Knepley /*@
11716ed231c7SMatthew Knepley     PCASMSetSortIndices - Determines whether subdomain indices are sorted.
11726ed231c7SMatthew Knepley 
1173ad4df100SBarry Smith     Logically Collective on PC
11746ed231c7SMatthew Knepley 
11756ed231c7SMatthew Knepley     Input Parameters:
11766ed231c7SMatthew Knepley +   pc  - the preconditioner context
11776ed231c7SMatthew Knepley -   doSort - sort the subdomain indices
11786ed231c7SMatthew Knepley 
11796ed231c7SMatthew Knepley     Level: intermediate
11806ed231c7SMatthew Knepley 
11816ed231c7SMatthew Knepley .keywords: PC, ASM, set, type
11826ed231c7SMatthew Knepley 
11836ed231c7SMatthew Knepley .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
11846ed231c7SMatthew Knepley           PCASMCreateSubdomains2D()
11856ed231c7SMatthew Knepley @*/
11867087cfbeSBarry Smith PetscErrorCode  PCASMSetSortIndices(PC pc,PetscBool doSort)
11876ed231c7SMatthew Knepley {
11887bb14e67SBarry Smith   PetscErrorCode ierr;
11896ed231c7SMatthew Knepley 
11906ed231c7SMatthew Knepley   PetscFunctionBegin;
11910700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1192acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
11937bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
11946ed231c7SMatthew Knepley   PetscFunctionReturn(0);
11956ed231c7SMatthew Knepley }
11966ed231c7SMatthew Knepley 
11974b9ad928SBarry Smith /*@C
11984b9ad928SBarry Smith    PCASMGetSubKSP - Gets the local KSP contexts for all blocks on
11994b9ad928SBarry Smith    this processor.
12004b9ad928SBarry Smith 
12014b9ad928SBarry Smith    Collective on PC iff first_local is requested
12024b9ad928SBarry Smith 
12034b9ad928SBarry Smith    Input Parameter:
12044b9ad928SBarry Smith .  pc - the preconditioner context
12054b9ad928SBarry Smith 
12064b9ad928SBarry Smith    Output Parameters:
12070298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
12080298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
12090298fd71SBarry Smith                  all processors must request or all must pass NULL
12104b9ad928SBarry Smith -  ksp - the array of KSP contexts
12114b9ad928SBarry Smith 
12124b9ad928SBarry Smith    Note:
1213d29017ddSJed Brown    After PCASMGetSubKSP() the array of KSPes is not to be freed.
12144b9ad928SBarry Smith 
12154b9ad928SBarry Smith    You must call KSPSetUp() before calling PCASMGetSubKSP().
12164b9ad928SBarry Smith 
1217d29017ddSJed Brown    Fortran note:
12182bf68e3eSBarry Smith    The output argument 'ksp' must be an array of sufficient length or PETSC_NULL_KSP. The latter can be used to learn the necessary length.
1219d29017ddSJed Brown 
12204b9ad928SBarry Smith    Level: advanced
12214b9ad928SBarry Smith 
12224b9ad928SBarry Smith .keywords: PC, ASM, additive Schwarz, get, sub, KSP, context
12234b9ad928SBarry Smith 
12244b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
12254b9ad928SBarry Smith           PCASMCreateSubdomains2D(),
12264b9ad928SBarry Smith @*/
12277087cfbeSBarry Smith PetscErrorCode  PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
12284b9ad928SBarry Smith {
12297bb14e67SBarry Smith   PetscErrorCode ierr;
12304b9ad928SBarry Smith 
12314b9ad928SBarry Smith   PetscFunctionBegin;
12320700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12337bb14e67SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
12344b9ad928SBarry Smith   PetscFunctionReturn(0);
12354b9ad928SBarry Smith }
12364b9ad928SBarry Smith 
12374b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
12384b9ad928SBarry Smith /*MC
12393b09bd56SBarry Smith    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
12404b9ad928SBarry Smith            its own KSP object.
12414b9ad928SBarry Smith 
12424b9ad928SBarry Smith    Options Database Keys:
124349517cdeSBarry Smith +  -pc_asm_blocks <blks> - Sets total blocks
12444b9ad928SBarry Smith .  -pc_asm_overlap <ovl> - Sets overlap
1245f1ee410cSBarry Smith .  -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict
1246f1ee410cSBarry Smith -  -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive
12474b9ad928SBarry Smith 
12483b09bd56SBarry Smith      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
12493b09bd56SBarry Smith       will get a different convergence rate due to the default option of -pc_asm_type restrict. Use
12503b09bd56SBarry Smith       -pc_asm_type basic to use the standard ASM.
12513b09bd56SBarry Smith 
12524b9ad928SBarry Smith    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
1253f1ee410cSBarry Smith      than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor.
12544b9ad928SBarry Smith 
12553b09bd56SBarry Smith      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1256d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
12574b9ad928SBarry Smith 
1258a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCASMGetSubKSP()
12594b9ad928SBarry Smith          and set the options directly on the resulting KSP object (you can access its PC
12604b9ad928SBarry Smith          with KSPGetPC())
12614b9ad928SBarry Smith 
12624b9ad928SBarry Smith    Level: beginner
12634b9ad928SBarry Smith 
12644b9ad928SBarry Smith    Concepts: additive Schwarz method
12654b9ad928SBarry Smith 
1266c582cd25SBarry Smith     References:
126796a0c994SBarry Smith +   1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
126896a0c994SBarry Smith      Courant Institute, New York University Technical report
126996a0c994SBarry Smith -   1. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
127096a0c994SBarry Smith     Cambridge University Press.
1271c582cd25SBarry Smith 
12724b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1273f1ee410cSBarry Smith            PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType()
1274f1ee410cSBarry Smith            PCASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType
1275e09e08ccSBarry Smith 
12764b9ad928SBarry Smith M*/
12774b9ad928SBarry Smith 
12788cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc)
12794b9ad928SBarry Smith {
1280dfbe8321SBarry Smith   PetscErrorCode ierr;
12814b9ad928SBarry Smith   PC_ASM         *osm;
12824b9ad928SBarry Smith 
12834b9ad928SBarry Smith   PetscFunctionBegin;
1284b00a9115SJed Brown   ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr);
12852fa5cd67SKarl Rupp 
12864b9ad928SBarry Smith   osm->n                 = PETSC_DECIDE;
12874b9ad928SBarry Smith   osm->n_local           = 0;
12882b837212SDmitry Karpeev   osm->n_local_true      = PETSC_DECIDE;
12894b9ad928SBarry Smith   osm->overlap           = 1;
12904b9ad928SBarry Smith   osm->ksp               = 0;
12912b691e39SMatthew Knepley   osm->restriction       = 0;
12925b3c0d42Seaulisa   osm->lprolongation     = 0;
12931dd8081eSeaulisa   osm->lrestriction      = 0;
129492bb6962SLisandro Dalcin   osm->x                 = 0;
129592bb6962SLisandro Dalcin   osm->y                 = 0;
12964b9ad928SBarry Smith   osm->is                = 0;
12972b691e39SMatthew Knepley   osm->is_local          = 0;
12984b9ad928SBarry Smith   osm->mat               = 0;
12994b9ad928SBarry Smith   osm->pmat              = 0;
13004b9ad928SBarry Smith   osm->type              = PC_ASM_RESTRICT;
130112cd4985SMatthew G. Knepley   osm->loctype           = PC_COMPOSITE_ADDITIVE;
13024b9ad928SBarry Smith   osm->same_local_solves = PETSC_TRUE;
13036ed231c7SMatthew Knepley   osm->sort_indices      = PETSC_TRUE;
1304d709ea83SDmitry Karpeev   osm->dm_subdomains     = PETSC_FALSE;
130580ec0b7dSPatrick Sanan   osm->sub_mat_type      = NULL;
13064b9ad928SBarry Smith 
130792bb6962SLisandro Dalcin   pc->data                 = (void*)osm;
13084b9ad928SBarry Smith   pc->ops->apply           = PCApply_ASM;
13094b9ad928SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_ASM;
13104b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_ASM;
1311e91c6855SBarry Smith   pc->ops->reset           = PCReset_ASM;
13124b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_ASM;
13134b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_ASM;
13144b9ad928SBarry Smith   pc->ops->setuponblocks   = PCSetUpOnBlocks_ASM;
13154b9ad928SBarry Smith   pc->ops->view            = PCView_ASM;
13164b9ad928SBarry Smith   pc->ops->applyrichardson = 0;
13174b9ad928SBarry Smith 
1318bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr);
1319bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr);
1320bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);CHKERRQ(ierr);
1321bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);CHKERRQ(ierr);
1322c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);CHKERRQ(ierr);
132312cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);CHKERRQ(ierr);
132412cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);CHKERRQ(ierr);
1325bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);CHKERRQ(ierr);
1326bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);CHKERRQ(ierr);
132780ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);CHKERRQ(ierr);
132880ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);CHKERRQ(ierr);
13294b9ad928SBarry Smith   PetscFunctionReturn(0);
13304b9ad928SBarry Smith }
13314b9ad928SBarry Smith 
133292bb6962SLisandro Dalcin /*@C
133392bb6962SLisandro Dalcin    PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
133492bb6962SLisandro Dalcin    preconditioner for a any problem on a general grid.
133592bb6962SLisandro Dalcin 
133692bb6962SLisandro Dalcin    Collective
133792bb6962SLisandro Dalcin 
133892bb6962SLisandro Dalcin    Input Parameters:
133992bb6962SLisandro Dalcin +  A - The global matrix operator
134092bb6962SLisandro Dalcin -  n - the number of local blocks
134192bb6962SLisandro Dalcin 
134292bb6962SLisandro Dalcin    Output Parameters:
134392bb6962SLisandro Dalcin .  outis - the array of index sets defining the subdomains
134492bb6962SLisandro Dalcin 
134592bb6962SLisandro Dalcin    Level: advanced
134692bb6962SLisandro Dalcin 
13477d6bfa3bSBarry Smith    Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap
13487d6bfa3bSBarry Smith     from these if you use PCASMSetLocalSubdomains()
13497d6bfa3bSBarry Smith 
13507d6bfa3bSBarry Smith     In the Fortran version you must provide the array outis[] already allocated of length n.
13517d6bfa3bSBarry Smith 
135292bb6962SLisandro Dalcin .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid
135392bb6962SLisandro Dalcin 
135492bb6962SLisandro Dalcin .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains()
135592bb6962SLisandro Dalcin @*/
13567087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
135792bb6962SLisandro Dalcin {
135892bb6962SLisandro Dalcin   MatPartitioning mpart;
135992bb6962SLisandro Dalcin   const char      *prefix;
1360c56a70eeSBarry Smith   void            (*f)(void);
136192bb6962SLisandro Dalcin   PetscInt        i,j,rstart,rend,bs;
136211bd1e4dSLisandro Dalcin   PetscBool       isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
13630298fd71SBarry Smith   Mat             Ad     = NULL, adj;
136492bb6962SLisandro Dalcin   IS              ispart,isnumb,*is;
136592bb6962SLisandro Dalcin   PetscErrorCode  ierr;
136692bb6962SLisandro Dalcin 
136792bb6962SLisandro Dalcin   PetscFunctionBegin;
13680700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
136992bb6962SLisandro Dalcin   PetscValidPointer(outis,3);
1370c1235816SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
137192bb6962SLisandro Dalcin 
137292bb6962SLisandro Dalcin   /* Get prefix, row distribution, and block size */
137392bb6962SLisandro Dalcin   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
137492bb6962SLisandro Dalcin   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
137592bb6962SLisandro Dalcin   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
137665e19b50SBarry Smith   if (rstart/bs*bs != rstart || rend/bs*bs != rend) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"bad row distribution [%D,%D) for matrix block size %D",rstart,rend,bs);
137765e19b50SBarry Smith 
137892bb6962SLisandro Dalcin   /* Get diagonal block from matrix if possible */
1379c56a70eeSBarry Smith   ierr = MatShellGetOperation(A,MATOP_GET_DIAGONAL_BLOCK,&f);CHKERRQ(ierr);
138092bb6962SLisandro Dalcin   if (f) {
138111bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
138292bb6962SLisandro Dalcin   }
138392bb6962SLisandro Dalcin   if (Ad) {
1384251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1385251f4c67SDmitry Karpeev     if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
138692bb6962SLisandro Dalcin   }
138792bb6962SLisandro Dalcin   if (Ad && n > 1) {
1388ace3abfcSBarry Smith     PetscBool match,done;
138992bb6962SLisandro Dalcin     /* Try to setup a good matrix partitioning if available */
139092bb6962SLisandro Dalcin     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
139192bb6962SLisandro Dalcin     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
139292bb6962SLisandro Dalcin     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1393251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
139492bb6962SLisandro Dalcin     if (!match) {
1395251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
139692bb6962SLisandro Dalcin     }
139792bb6962SLisandro Dalcin     if (!match) { /* assume a "good" partitioner is available */
13981a83f524SJed Brown       PetscInt       na;
13991a83f524SJed Brown       const PetscInt *ia,*ja;
140092bb6962SLisandro Dalcin       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
140192bb6962SLisandro Dalcin       if (done) {
140292bb6962SLisandro Dalcin         /* Build adjacency matrix by hand. Unfortunately a call to
140392bb6962SLisandro Dalcin            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
140492bb6962SLisandro Dalcin            remove the block-aij structure and we cannot expect
140592bb6962SLisandro Dalcin            MatPartitioning to split vertices as we need */
14061a83f524SJed Brown         PetscInt       i,j,len,nnz,cnt,*iia=0,*jja=0;
14071a83f524SJed Brown         const PetscInt *row;
140892bb6962SLisandro Dalcin         nnz = 0;
140992bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* count number of nonzeros */
141092bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
141192bb6962SLisandro Dalcin           row = ja + ia[i];
141292bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
141392bb6962SLisandro Dalcin             if (row[j] == i) { /* don't count diagonal */
141492bb6962SLisandro Dalcin               len--; break;
141592bb6962SLisandro Dalcin             }
141692bb6962SLisandro Dalcin           }
141792bb6962SLisandro Dalcin           nnz += len;
141892bb6962SLisandro Dalcin         }
1419854ce69bSBarry Smith         ierr   = PetscMalloc1(na+1,&iia);CHKERRQ(ierr);
1420854ce69bSBarry Smith         ierr   = PetscMalloc1(nnz,&jja);CHKERRQ(ierr);
142192bb6962SLisandro Dalcin         nnz    = 0;
142292bb6962SLisandro Dalcin         iia[0] = 0;
142392bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* fill adjacency */
142492bb6962SLisandro Dalcin           cnt = 0;
142592bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
142692bb6962SLisandro Dalcin           row = ja + ia[i];
142792bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
142892bb6962SLisandro Dalcin             if (row[j] != i) { /* if not diagonal */
142992bb6962SLisandro Dalcin               jja[nnz+cnt++] = row[j];
143092bb6962SLisandro Dalcin             }
143192bb6962SLisandro Dalcin           }
143292bb6962SLisandro Dalcin           nnz     += cnt;
143392bb6962SLisandro Dalcin           iia[i+1] = nnz;
143492bb6962SLisandro Dalcin         }
143592bb6962SLisandro Dalcin         /* Partitioning of the adjacency matrix */
14360298fd71SBarry Smith         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr);
143792bb6962SLisandro Dalcin         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
143892bb6962SLisandro Dalcin         ierr      = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
143992bb6962SLisandro Dalcin         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
144092bb6962SLisandro Dalcin         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
1441fcfd50ebSBarry Smith         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
144292bb6962SLisandro Dalcin         foundpart = PETSC_TRUE;
144392bb6962SLisandro Dalcin       }
144492bb6962SLisandro Dalcin       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
144592bb6962SLisandro Dalcin     }
1446fcfd50ebSBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
144792bb6962SLisandro Dalcin   }
144892bb6962SLisandro Dalcin 
1449785e854fSJed Brown   ierr   = PetscMalloc1(n,&is);CHKERRQ(ierr);
145092bb6962SLisandro Dalcin   *outis = is;
145192bb6962SLisandro Dalcin 
145292bb6962SLisandro Dalcin   if (!foundpart) {
145392bb6962SLisandro Dalcin 
145492bb6962SLisandro Dalcin     /* Partitioning by contiguous chunks of rows */
145592bb6962SLisandro Dalcin 
145692bb6962SLisandro Dalcin     PetscInt mbs   = (rend-rstart)/bs;
145792bb6962SLisandro Dalcin     PetscInt start = rstart;
145892bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
145992bb6962SLisandro Dalcin       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
146092bb6962SLisandro Dalcin       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
146192bb6962SLisandro Dalcin       start += count;
146292bb6962SLisandro Dalcin     }
146392bb6962SLisandro Dalcin 
146492bb6962SLisandro Dalcin   } else {
146592bb6962SLisandro Dalcin 
146692bb6962SLisandro Dalcin     /* Partitioning by adjacency of diagonal block  */
146792bb6962SLisandro Dalcin 
146892bb6962SLisandro Dalcin     const PetscInt *numbering;
146992bb6962SLisandro Dalcin     PetscInt       *count,nidx,*indices,*newidx,start=0;
147092bb6962SLisandro Dalcin     /* Get node count in each partition */
1471785e854fSJed Brown     ierr = PetscMalloc1(n,&count);CHKERRQ(ierr);
147292bb6962SLisandro Dalcin     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
147392bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
147492bb6962SLisandro Dalcin       for (i=0; i<n; i++) count[i] *= bs;
147592bb6962SLisandro Dalcin     }
147692bb6962SLisandro Dalcin     /* Build indices from node numbering */
147792bb6962SLisandro Dalcin     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1478785e854fSJed Brown     ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr);
147992bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
148092bb6962SLisandro Dalcin     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
148192bb6962SLisandro Dalcin     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
148292bb6962SLisandro Dalcin     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
148392bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1484785e854fSJed Brown       ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr);
14852fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
14862fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
14872fa5cd67SKarl Rupp       }
148892bb6962SLisandro Dalcin       ierr    = PetscFree(indices);CHKERRQ(ierr);
148992bb6962SLisandro Dalcin       nidx   *= bs;
149092bb6962SLisandro Dalcin       indices = newidx;
149192bb6962SLisandro Dalcin     }
149292bb6962SLisandro Dalcin     /* Shift to get global indices */
149392bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] += rstart;
149492bb6962SLisandro Dalcin 
149592bb6962SLisandro Dalcin     /* Build the index sets for each block */
149692bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
149770b3c8c7SBarry Smith       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
149892bb6962SLisandro Dalcin       ierr   = ISSort(is[i]);CHKERRQ(ierr);
149992bb6962SLisandro Dalcin       start += count[i];
150092bb6962SLisandro Dalcin     }
150192bb6962SLisandro Dalcin 
15023bf036e2SBarry Smith     ierr = PetscFree(count);CHKERRQ(ierr);
15033bf036e2SBarry Smith     ierr = PetscFree(indices);CHKERRQ(ierr);
1504fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1505fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
150692bb6962SLisandro Dalcin 
150792bb6962SLisandro Dalcin   }
150892bb6962SLisandro Dalcin   PetscFunctionReturn(0);
150992bb6962SLisandro Dalcin }
151092bb6962SLisandro Dalcin 
151192bb6962SLisandro Dalcin /*@C
151292bb6962SLisandro Dalcin    PCASMDestroySubdomains - Destroys the index sets created with
151392bb6962SLisandro Dalcin    PCASMCreateSubdomains(). Should be called after setting subdomains
151492bb6962SLisandro Dalcin    with PCASMSetLocalSubdomains().
151592bb6962SLisandro Dalcin 
151692bb6962SLisandro Dalcin    Collective
151792bb6962SLisandro Dalcin 
151892bb6962SLisandro Dalcin    Input Parameters:
151992bb6962SLisandro Dalcin +  n - the number of index sets
15202b691e39SMatthew Knepley .  is - the array of index sets
15210298fd71SBarry Smith -  is_local - the array of local index sets, can be NULL
152292bb6962SLisandro Dalcin 
152392bb6962SLisandro Dalcin    Level: advanced
152492bb6962SLisandro Dalcin 
152592bb6962SLisandro Dalcin .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid
152692bb6962SLisandro Dalcin 
152792bb6962SLisandro Dalcin .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains()
152892bb6962SLisandro Dalcin @*/
15297087cfbeSBarry Smith PetscErrorCode  PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
153092bb6962SLisandro Dalcin {
153192bb6962SLisandro Dalcin   PetscInt       i;
153292bb6962SLisandro Dalcin   PetscErrorCode ierr;
15335fd66863SKarl Rupp 
153492bb6962SLisandro Dalcin   PetscFunctionBegin;
1535a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1536a35b7b57SDmitry Karpeev   if (is) {
153792bb6962SLisandro Dalcin     PetscValidPointer(is,2);
1538fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); }
153992bb6962SLisandro Dalcin     ierr = PetscFree(is);CHKERRQ(ierr);
1540a35b7b57SDmitry Karpeev   }
15412b691e39SMatthew Knepley   if (is_local) {
15422b691e39SMatthew Knepley     PetscValidPointer(is_local,3);
1543fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); }
15442b691e39SMatthew Knepley     ierr = PetscFree(is_local);CHKERRQ(ierr);
15452b691e39SMatthew Knepley   }
154692bb6962SLisandro Dalcin   PetscFunctionReturn(0);
154792bb6962SLisandro Dalcin }
154892bb6962SLisandro Dalcin 
15494b9ad928SBarry Smith /*@
15504b9ad928SBarry Smith    PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
15514b9ad928SBarry Smith    preconditioner for a two-dimensional problem on a regular grid.
15524b9ad928SBarry Smith 
15534b9ad928SBarry Smith    Not Collective
15544b9ad928SBarry Smith 
15554b9ad928SBarry Smith    Input Parameters:
15564b9ad928SBarry Smith +  m, n - the number of mesh points in the x and y directions
15574b9ad928SBarry Smith .  M, N - the number of subdomains in the x and y directions
15584b9ad928SBarry Smith .  dof - degrees of freedom per node
15594b9ad928SBarry Smith -  overlap - overlap in mesh lines
15604b9ad928SBarry Smith 
15614b9ad928SBarry Smith    Output Parameters:
15624b9ad928SBarry Smith +  Nsub - the number of subdomains created
15633d03bb48SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
15643d03bb48SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
15654b9ad928SBarry Smith 
15664b9ad928SBarry Smith    Note:
15674b9ad928SBarry Smith    Presently PCAMSCreateSubdomains2d() is valid only for sequential
15684b9ad928SBarry Smith    preconditioners.  More general related routines are
15694b9ad928SBarry Smith    PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().
15704b9ad928SBarry Smith 
15714b9ad928SBarry Smith    Level: advanced
15724b9ad928SBarry Smith 
15734b9ad928SBarry Smith .keywords: PC, ASM, additive Schwarz, create, subdomains, 2D, regular grid
15744b9ad928SBarry Smith 
15754b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
15764b9ad928SBarry Smith           PCASMSetOverlap()
15774b9ad928SBarry Smith @*/
15787087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local)
15794b9ad928SBarry Smith {
15803d03bb48SJed Brown   PetscInt       i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer;
15816849ba73SBarry Smith   PetscErrorCode ierr;
158213f74950SBarry Smith   PetscInt       nidx,*idx,loc,ii,jj,count;
15834b9ad928SBarry Smith 
15844b9ad928SBarry Smith   PetscFunctionBegin;
1585e32f2f54SBarry Smith   if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," ");
15864b9ad928SBarry Smith 
15874b9ad928SBarry Smith   *Nsub     = N*M;
1588854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is);CHKERRQ(ierr);
1589854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is_local);CHKERRQ(ierr);
15904b9ad928SBarry Smith   ystart    = 0;
15913d03bb48SJed Brown   loc_outer = 0;
15924b9ad928SBarry Smith   for (i=0; i<N; i++) {
15934b9ad928SBarry Smith     height = n/N + ((n % N) > i); /* height of subdomain */
1594e32f2f54SBarry Smith     if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n");
15954b9ad928SBarry Smith     yleft  = ystart - overlap; if (yleft < 0) yleft = 0;
15964b9ad928SBarry Smith     yright = ystart + height + overlap; if (yright > n) yright = n;
15974b9ad928SBarry Smith     xstart = 0;
15984b9ad928SBarry Smith     for (j=0; j<M; j++) {
15994b9ad928SBarry Smith       width = m/M + ((m % M) > j); /* width of subdomain */
1600e32f2f54SBarry Smith       if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m");
16014b9ad928SBarry Smith       xleft  = xstart - overlap; if (xleft < 0) xleft = 0;
16024b9ad928SBarry Smith       xright = xstart + width + overlap; if (xright > m) xright = m;
16034b9ad928SBarry Smith       nidx   = (xright - xleft)*(yright - yleft);
1604785e854fSJed Brown       ierr   = PetscMalloc1(nidx,&idx);CHKERRQ(ierr);
16054b9ad928SBarry Smith       loc    = 0;
16064b9ad928SBarry Smith       for (ii=yleft; ii<yright; ii++) {
16074b9ad928SBarry Smith         count = m*ii + xleft;
16082fa5cd67SKarl Rupp         for (jj=xleft; jj<xright; jj++) idx[loc++] = count++;
16094b9ad928SBarry Smith       }
161070b3c8c7SBarry Smith       ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr);
16113d03bb48SJed Brown       if (overlap == 0) {
16123d03bb48SJed Brown         ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr);
16132fa5cd67SKarl Rupp 
16143d03bb48SJed Brown         (*is_local)[loc_outer] = (*is)[loc_outer];
16153d03bb48SJed Brown       } else {
16163d03bb48SJed Brown         for (loc=0,ii=ystart; ii<ystart+height; ii++) {
16173d03bb48SJed Brown           for (jj=xstart; jj<xstart+width; jj++) {
16183d03bb48SJed Brown             idx[loc++] = m*ii + jj;
16193d03bb48SJed Brown           }
16203d03bb48SJed Brown         }
162170b3c8c7SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr);
16223d03bb48SJed Brown       }
16234b9ad928SBarry Smith       ierr    = PetscFree(idx);CHKERRQ(ierr);
16244b9ad928SBarry Smith       xstart += width;
16253d03bb48SJed Brown       loc_outer++;
16264b9ad928SBarry Smith     }
16274b9ad928SBarry Smith     ystart += height;
16284b9ad928SBarry Smith   }
16294b9ad928SBarry Smith   for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); }
16304b9ad928SBarry Smith   PetscFunctionReturn(0);
16314b9ad928SBarry Smith }
16324b9ad928SBarry Smith 
16334b9ad928SBarry Smith /*@C
16344b9ad928SBarry Smith     PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
16354b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
16364b9ad928SBarry Smith 
1637ad4df100SBarry Smith     Not Collective
16384b9ad928SBarry Smith 
16394b9ad928SBarry Smith     Input Parameter:
16404b9ad928SBarry Smith .   pc - the preconditioner context
16414b9ad928SBarry Smith 
16424b9ad928SBarry Smith     Output Parameters:
16434b9ad928SBarry Smith +   n - the number of subdomains for this processor (default value = 1)
16442b691e39SMatthew Knepley .   is - the index sets that define the subdomains for this processor
16450298fd71SBarry Smith -   is_local - the index sets that define the local part of the subdomains for this processor (can be NULL)
16464b9ad928SBarry Smith 
16474b9ad928SBarry Smith 
16484b9ad928SBarry Smith     Notes:
16494b9ad928SBarry Smith     The IS numbering is in the parallel, global numbering of the vector.
16504b9ad928SBarry Smith 
16514b9ad928SBarry Smith     Level: advanced
16524b9ad928SBarry Smith 
16534b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz
16544b9ad928SBarry Smith 
16554b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
16564b9ad928SBarry Smith           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices()
16574b9ad928SBarry Smith @*/
16587087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
16594b9ad928SBarry Smith {
16602a808120SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
166192bb6962SLisandro Dalcin   PetscErrorCode ierr;
1662ace3abfcSBarry Smith   PetscBool      match;
16634b9ad928SBarry Smith 
16644b9ad928SBarry Smith   PetscFunctionBegin;
16650700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
16664482741eSBarry Smith   PetscValidIntPointer(n,2);
166792bb6962SLisandro Dalcin   if (is) PetscValidPointer(is,3);
1668251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
16692a808120SBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM");
16704b9ad928SBarry Smith   if (n) *n = osm->n_local_true;
16714b9ad928SBarry Smith   if (is) *is = osm->is;
16722b691e39SMatthew Knepley   if (is_local) *is_local = osm->is_local;
16734b9ad928SBarry Smith   PetscFunctionReturn(0);
16744b9ad928SBarry Smith }
16754b9ad928SBarry Smith 
16764b9ad928SBarry Smith /*@C
16774b9ad928SBarry Smith     PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
16784b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
16794b9ad928SBarry Smith 
1680ad4df100SBarry Smith     Not Collective
16814b9ad928SBarry Smith 
16824b9ad928SBarry Smith     Input Parameter:
16834b9ad928SBarry Smith .   pc - the preconditioner context
16844b9ad928SBarry Smith 
16854b9ad928SBarry Smith     Output Parameters:
16864b9ad928SBarry Smith +   n - the number of matrices for this processor (default value = 1)
16874b9ad928SBarry Smith -   mat - the matrices
16884b9ad928SBarry Smith 
16894b9ad928SBarry Smith     Level: advanced
16904b9ad928SBarry Smith 
1691cf739d55SBarry Smith     Notes: Call after PCSetUp() (or KSPSetUp()) but before PCApply() (or KSPApply()) and before PCSetUpOnBlocks())
1692cf739d55SBarry Smith 
1693cf739d55SBarry Smith            Usually one would use PCSetModifySubmatrices() to change the submatrices in building the preconditioner.
1694cf739d55SBarry Smith 
16954b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz, block Jacobi
16964b9ad928SBarry Smith 
16974b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1698cf739d55SBarry Smith           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubmatrices()
16994b9ad928SBarry Smith @*/
17007087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
17014b9ad928SBarry Smith {
17024b9ad928SBarry Smith   PC_ASM         *osm;
170392bb6962SLisandro Dalcin   PetscErrorCode ierr;
1704ace3abfcSBarry Smith   PetscBool      match;
17054b9ad928SBarry Smith 
17064b9ad928SBarry Smith   PetscFunctionBegin;
17070700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
170892bb6962SLisandro Dalcin   PetscValidIntPointer(n,2);
170992bb6962SLisandro Dalcin   if (mat) PetscValidPointer(mat,3);
1710ce94432eSBarry Smith   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1711251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
171292bb6962SLisandro Dalcin   if (!match) {
171392bb6962SLisandro Dalcin     if (n) *n = 0;
17140298fd71SBarry Smith     if (mat) *mat = NULL;
171592bb6962SLisandro Dalcin   } else {
17164b9ad928SBarry Smith     osm = (PC_ASM*)pc->data;
17174b9ad928SBarry Smith     if (n) *n = osm->n_local_true;
17184b9ad928SBarry Smith     if (mat) *mat = osm->pmat;
171992bb6962SLisandro Dalcin   }
17204b9ad928SBarry Smith   PetscFunctionReturn(0);
17214b9ad928SBarry Smith }
1722d709ea83SDmitry Karpeev 
1723d709ea83SDmitry Karpeev /*@
1724d709ea83SDmitry Karpeev     PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1725f1ee410cSBarry Smith 
1726d709ea83SDmitry Karpeev     Logically Collective
1727d709ea83SDmitry Karpeev 
1728d709ea83SDmitry Karpeev     Input Parameter:
1729d709ea83SDmitry Karpeev +   pc  - the preconditioner
1730d709ea83SDmitry Karpeev -   flg - boolean indicating whether to use subdomains defined by the DM
1731d709ea83SDmitry Karpeev 
1732d709ea83SDmitry Karpeev     Options Database Key:
1733d709ea83SDmitry Karpeev .   -pc_asm_dm_subdomains
1734d709ea83SDmitry Karpeev 
1735d709ea83SDmitry Karpeev     Level: intermediate
1736d709ea83SDmitry Karpeev 
1737d709ea83SDmitry Karpeev     Notes:
1738d709ea83SDmitry Karpeev     PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(),
1739d709ea83SDmitry Karpeev     so setting either of the first two effectively turns the latter off.
1740d709ea83SDmitry Karpeev 
1741d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1742d709ea83SDmitry Karpeev 
1743d709ea83SDmitry Karpeev .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1744d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1745d709ea83SDmitry Karpeev @*/
1746d709ea83SDmitry Karpeev PetscErrorCode  PCASMSetDMSubdomains(PC pc,PetscBool flg)
1747d709ea83SDmitry Karpeev {
1748d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1749d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1750d709ea83SDmitry Karpeev   PetscBool      match;
1751d709ea83SDmitry Karpeev 
1752d709ea83SDmitry Karpeev   PetscFunctionBegin;
1753d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1754d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
1755d709ea83SDmitry Karpeev   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1756d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1757d709ea83SDmitry Karpeev   if (match) {
1758d709ea83SDmitry Karpeev     osm->dm_subdomains = flg;
1759d709ea83SDmitry Karpeev   }
1760d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1761d709ea83SDmitry Karpeev }
1762d709ea83SDmitry Karpeev 
1763d709ea83SDmitry Karpeev /*@
1764d709ea83SDmitry Karpeev     PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1765d709ea83SDmitry Karpeev     Not Collective
1766d709ea83SDmitry Karpeev 
1767d709ea83SDmitry Karpeev     Input Parameter:
1768d709ea83SDmitry Karpeev .   pc  - the preconditioner
1769d709ea83SDmitry Karpeev 
1770d709ea83SDmitry Karpeev     Output Parameter:
1771d709ea83SDmitry Karpeev .   flg - boolean indicating whether to use subdomains defined by the DM
1772d709ea83SDmitry Karpeev 
1773d709ea83SDmitry Karpeev     Level: intermediate
1774d709ea83SDmitry Karpeev 
1775d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1776d709ea83SDmitry Karpeev 
1777d709ea83SDmitry Karpeev .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1778d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1779d709ea83SDmitry Karpeev @*/
1780d709ea83SDmitry Karpeev PetscErrorCode  PCASMGetDMSubdomains(PC pc,PetscBool* flg)
1781d709ea83SDmitry Karpeev {
1782d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1783d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1784d709ea83SDmitry Karpeev   PetscBool      match;
1785d709ea83SDmitry Karpeev 
1786d709ea83SDmitry Karpeev   PetscFunctionBegin;
1787d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1788d709ea83SDmitry Karpeev   PetscValidPointer(flg,2);
1789d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1790d709ea83SDmitry Karpeev   if (match) {
1791d709ea83SDmitry Karpeev     if (flg) *flg = osm->dm_subdomains;
1792d709ea83SDmitry Karpeev   }
1793d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1794d709ea83SDmitry Karpeev }
179580ec0b7dSPatrick Sanan 
179680ec0b7dSPatrick Sanan /*@
179780ec0b7dSPatrick Sanan      PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string.
179880ec0b7dSPatrick Sanan 
179980ec0b7dSPatrick Sanan    Not Collective
180080ec0b7dSPatrick Sanan 
180180ec0b7dSPatrick Sanan    Input Parameter:
180280ec0b7dSPatrick Sanan .  pc - the PC
180380ec0b7dSPatrick Sanan 
180480ec0b7dSPatrick Sanan    Output Parameter:
1805f1ee410cSBarry Smith .  -pc_asm_sub_mat_type - name of matrix type
180680ec0b7dSPatrick Sanan 
180780ec0b7dSPatrick Sanan    Level: advanced
180880ec0b7dSPatrick Sanan 
180980ec0b7dSPatrick Sanan .keywords: PC, PCASM, MatType, set
181080ec0b7dSPatrick Sanan 
181180ec0b7dSPatrick Sanan .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
181280ec0b7dSPatrick Sanan @*/
181380ec0b7dSPatrick Sanan PetscErrorCode  PCASMGetSubMatType(PC pc,MatType *sub_mat_type){
181480ec0b7dSPatrick Sanan   PetscErrorCode ierr;
181580ec0b7dSPatrick Sanan 
181680ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));CHKERRQ(ierr);
181780ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
181880ec0b7dSPatrick Sanan }
181980ec0b7dSPatrick Sanan 
182080ec0b7dSPatrick Sanan /*@
182180ec0b7dSPatrick Sanan      PCASMSetSubMatType - Set the type of matrix used for ASM subsolves
182280ec0b7dSPatrick Sanan 
182380ec0b7dSPatrick Sanan    Collective on Mat
182480ec0b7dSPatrick Sanan 
182580ec0b7dSPatrick Sanan    Input Parameters:
182680ec0b7dSPatrick Sanan +  pc             - the PC object
182780ec0b7dSPatrick Sanan -  sub_mat_type   - matrix type
182880ec0b7dSPatrick Sanan 
182980ec0b7dSPatrick Sanan    Options Database Key:
183080ec0b7dSPatrick Sanan .  -pc_asm_sub_mat_type  <sub_mat_type> - Sets the matrix type used for subsolves, for example, seqaijviennacl. If you specify a base name like aijviennacl, the corresponding sequential type is assumed.
183180ec0b7dSPatrick Sanan 
183280ec0b7dSPatrick Sanan    Notes:
183380ec0b7dSPatrick Sanan    See "${PETSC_DIR}/include/petscmat.h" for available types
183480ec0b7dSPatrick Sanan 
183580ec0b7dSPatrick Sanan   Level: advanced
183680ec0b7dSPatrick Sanan 
183780ec0b7dSPatrick Sanan .keywords: PC, PCASM, MatType, set
183880ec0b7dSPatrick Sanan 
183980ec0b7dSPatrick Sanan .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
184080ec0b7dSPatrick Sanan @*/
184180ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type)
184280ec0b7dSPatrick Sanan {
184380ec0b7dSPatrick Sanan   PetscErrorCode ierr;
184480ec0b7dSPatrick Sanan 
184580ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));CHKERRQ(ierr);
184680ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
184780ec0b7dSPatrick Sanan }
1848