xref: /petsc/src/ksp/pc/impls/asm/asm.c (revision d0ecd4df9cacc91dc5dd638d339d9c6d5b321098)
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 */
203c4dc4c4SMatthew Knepley   VecScatter *restriction;        /* mapping from global to subregion */
21f1ee410cSBarry Smith   VecScatter *localization;       /* mapping from overlapping to non-overlapping subregion; used for restrict version of algorithms */
223c4dc4c4SMatthew Knepley   VecScatter *prolongation;       /* mapping from subregion to global */
2353888de8SMatthew Knepley   Vec        *x,*y,*y_local;      /* work vectors */
242b691e39SMatthew Knepley   IS         *is;                 /* index set that defines each overlapping subdomain */
252b691e39SMatthew Knepley   IS         *is_local;           /* index set that defines each non-overlapping subdomain, may be NULL */
264b9ad928SBarry Smith   Mat        *mat,*pmat;          /* mat is not currently used */
274b9ad928SBarry Smith   PCASMType  type;                /* use reduced interpolation, restriction or both */
28ace3abfcSBarry Smith   PetscBool  type_set;            /* if user set this value (so won't change it for symmetric problems) */
29ace3abfcSBarry Smith   PetscBool  same_local_solves;   /* flag indicating whether all local solvers are same */
30ace3abfcSBarry Smith   PetscBool  sort_indices;        /* flag to sort subdomain indices */
31d709ea83SDmitry Karpeev   PetscBool  dm_subdomains;       /* whether DM is allowed to define subdomains */
3212cd4985SMatthew G. Knepley   PCCompositeType loctype;        /* the type of composition for local solves */
3380ec0b7dSPatrick Sanan   MatType    sub_mat_type;        /* the type of Mat used for subdomain solves (can be MATSAME or NULL) */
34fb745f2cSMatthew G. Knepley   /* For multiplicative solve */
35235cc4ceSMatthew G. Knepley   Mat       *lmats;               /* submatrices for overlapping multiplicative (process) subdomain */
36fb745f2cSMatthew G. Knepley   Vec        lx, ly;              /* work vectors */
37fb745f2cSMatthew G. Knepley   IS         lis;                 /* index set that defines each 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) {
19192bb6962SLisandro Dalcin 
19292bb6962SLisandro Dalcin     if (!osm->type_set) {
19392bb6962SLisandro Dalcin       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
1942fa5cd67SKarl Rupp       if (symset && flg) osm->type = PC_ASM_BASIC;
19592bb6962SLisandro Dalcin     }
19692bb6962SLisandro Dalcin 
1972b837212SDmitry Karpeev     /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */
1982b837212SDmitry Karpeev     if (osm->n_local_true == PETSC_DECIDE) {
19969ca1f37SDmitry Karpeev       /* no subdomains given */
20065db9045SDmitry Karpeev       /* try pc->dm first, if allowed */
201d709ea83SDmitry Karpeev       if (osm->dm_subdomains && pc->dm) {
202feb221f8SDmitry Karpeev         PetscInt  num_domains, d;
203feb221f8SDmitry Karpeev         char      **domain_names;
2048d4ac253SDmitry Karpeev         IS        *inner_domain_is, *outer_domain_is;
2058d4ac253SDmitry Karpeev         ierr = DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm);CHKERRQ(ierr);
206704f0395SPatrick Sanan         osm->overlap = -1; /* We do not want to increase the overlap of the IS.
207704f0395SPatrick Sanan                               A future improvement of this code might allow one to use
208704f0395SPatrick Sanan                               DM-defined subdomains and also increase the overlap,
209704f0395SPatrick Sanan                               but that is not currently supported */
21069ca1f37SDmitry Karpeev         if (num_domains) {
2118d4ac253SDmitry Karpeev           ierr = PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is);CHKERRQ(ierr);
21269ca1f37SDmitry Karpeev         }
213feb221f8SDmitry Karpeev         for (d = 0; d < num_domains; ++d) {
214a35b7b57SDmitry Karpeev           if (domain_names)    {ierr = PetscFree(domain_names[d]);CHKERRQ(ierr);}
215a35b7b57SDmitry Karpeev           if (inner_domain_is) {ierr = ISDestroy(&inner_domain_is[d]);CHKERRQ(ierr);}
216a35b7b57SDmitry Karpeev           if (outer_domain_is) {ierr = ISDestroy(&outer_domain_is[d]);CHKERRQ(ierr);}
217feb221f8SDmitry Karpeev         }
218feb221f8SDmitry Karpeev         ierr = PetscFree(domain_names);CHKERRQ(ierr);
2198d4ac253SDmitry Karpeev         ierr = PetscFree(inner_domain_is);CHKERRQ(ierr);
2208d4ac253SDmitry Karpeev         ierr = PetscFree(outer_domain_is);CHKERRQ(ierr);
221feb221f8SDmitry Karpeev       }
2222b837212SDmitry Karpeev       if (osm->n_local_true == PETSC_DECIDE) {
22369ca1f37SDmitry Karpeev         /* still no subdomains; use one subdomain per processor */
2242b837212SDmitry Karpeev         osm->n_local_true = 1;
225feb221f8SDmitry Karpeev       }
2262b837212SDmitry Karpeev     }
2272b837212SDmitry Karpeev     { /* determine the global and max number of subdomains */
2286ac3741eSJed Brown       struct {PetscInt max,sum;} inwork,outwork;
229c10200c1SHong Zhang       PetscMPIInt size;
230c10200c1SHong Zhang 
2316ac3741eSJed Brown       inwork.max   = osm->n_local_true;
2326ac3741eSJed Brown       inwork.sum   = osm->n_local_true;
233367daffbSBarry Smith       ierr         = MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
2346ac3741eSJed Brown       osm->n_local = outwork.max;
2356ac3741eSJed Brown       osm->n       = outwork.sum;
236c10200c1SHong Zhang 
237c10200c1SHong Zhang       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
238c10200c1SHong Zhang       if (outwork.max == 1 && outwork.sum == size) {
2397dae84e0SHong Zhang         /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */
240c10200c1SHong Zhang         ierr = MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE);CHKERRQ(ierr);
241c10200c1SHong Zhang       }
2424b9ad928SBarry Smith     }
24378904715SLisandro Dalcin     if (!osm->is) { /* create the index sets */
24478904715SLisandro Dalcin       ierr = PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);CHKERRQ(ierr);
2454b9ad928SBarry Smith     }
246f5234e35SJed Brown     if (osm->n_local_true > 1 && !osm->is_local) {
247785e854fSJed Brown       ierr = PetscMalloc1(osm->n_local_true,&osm->is_local);CHKERRQ(ierr);
248f5234e35SJed Brown       for (i=0; i<osm->n_local_true; i++) {
249f5234e35SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
250f5234e35SJed Brown           ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr);
251f5234e35SJed Brown           ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr);
252f5234e35SJed Brown         } else {
253f5234e35SJed Brown           ierr             = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr);
254f5234e35SJed Brown           osm->is_local[i] = osm->is[i];
255f5234e35SJed Brown         }
256f5234e35SJed Brown       }
257f5234e35SJed Brown     }
25892bb6962SLisandro Dalcin     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
25990d69ab7SBarry Smith     flg  = PETSC_FALSE;
260c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,prefix,"-pc_asm_print_subdomains",&flg,NULL);CHKERRQ(ierr);
26192bb6962SLisandro Dalcin     if (flg) { ierr = PCASMPrintSubdomains(pc);CHKERRQ(ierr); }
2624b9ad928SBarry Smith 
2633d03bb48SJed Brown     if (osm->overlap > 0) {
2644b9ad928SBarry Smith       /* Extend the "overlapping" regions by a number of steps */
26592bb6962SLisandro Dalcin       ierr = MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);CHKERRQ(ierr);
2663d03bb48SJed Brown     }
2676ed231c7SMatthew Knepley     if (osm->sort_indices) {
26892bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2694b9ad928SBarry Smith         ierr = ISSort(osm->is[i]);CHKERRQ(ierr);
2702b691e39SMatthew Knepley         if (osm->is_local) {
2712b691e39SMatthew Knepley           ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr);
2722b691e39SMatthew Knepley         }
2734b9ad928SBarry Smith       }
2746ed231c7SMatthew Knepley     }
2754b9ad928SBarry Smith 
276f6991133SBarry Smith     if (!osm->ksp) {
27778904715SLisandro Dalcin       /* Create the local solvers */
278785e854fSJed Brown       ierr = PetscMalloc1(osm->n_local_true,&osm->ksp);CHKERRQ(ierr);
279feb221f8SDmitry Karpeev       if (domain_dm) {
280feb221f8SDmitry Karpeev         ierr = PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n");CHKERRQ(ierr);
281feb221f8SDmitry Karpeev       }
28292bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2834b9ad928SBarry Smith         ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
284422a814eSBarry Smith         ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr);
2853bb1ff40SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr);
28692bb6962SLisandro Dalcin         ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
2874b9ad928SBarry Smith         ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
2884b9ad928SBarry Smith         ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
2894b9ad928SBarry Smith         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
2904b9ad928SBarry Smith         ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
2914b9ad928SBarry Smith         ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
292feb221f8SDmitry Karpeev         if (domain_dm) {
293feb221f8SDmitry Karpeev           ierr = KSPSetDM(ksp, domain_dm[i]);CHKERRQ(ierr);
294feb221f8SDmitry Karpeev           ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr);
295feb221f8SDmitry Karpeev           ierr = DMDestroy(&domain_dm[i]);CHKERRQ(ierr);
296feb221f8SDmitry Karpeev         }
2974b9ad928SBarry Smith         osm->ksp[i] = ksp;
2984b9ad928SBarry Smith       }
299feb221f8SDmitry Karpeev       if (domain_dm) {
300feb221f8SDmitry Karpeev         ierr = PetscFree(domain_dm);CHKERRQ(ierr);
301feb221f8SDmitry Karpeev       }
302f6991133SBarry Smith     }
303fb745f2cSMatthew G. Knepley     if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
304fb745f2cSMatthew G. Knepley       PetscInt m;
305fb745f2cSMatthew G. Knepley 
306*d0ecd4dfSBarry Smith       ierr = ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is_local, &osm->lis);CHKERRQ(ierr);
307fb745f2cSMatthew G. Knepley       ierr = ISSortRemoveDups(osm->lis);CHKERRQ(ierr);
308fb745f2cSMatthew G. Knepley       ierr = ISGetLocalSize(osm->lis, &m);CHKERRQ(ierr);
309fb745f2cSMatthew G. Knepley       ierr = VecCreateSeq(PETSC_COMM_SELF, m, &osm->lx);CHKERRQ(ierr);
310fb745f2cSMatthew G. Knepley       ierr = VecDuplicate(osm->lx, &osm->ly);CHKERRQ(ierr);
311fb745f2cSMatthew G. Knepley     }
3124b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
3134b9ad928SBarry Smith   } else {
3144b9ad928SBarry Smith     /*
3154b9ad928SBarry Smith        Destroy the blocks from the previous iteration
3164b9ad928SBarry Smith     */
317e09e08ccSBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
3184b9ad928SBarry Smith       ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
3194b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
3204b9ad928SBarry Smith     }
3214b9ad928SBarry Smith   }
3224b9ad928SBarry Smith 
32392bb6962SLisandro Dalcin   /*
32492bb6962SLisandro Dalcin      Extract out the submatrices
32592bb6962SLisandro Dalcin   */
3267dae84e0SHong Zhang   ierr = MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr);
32792bb6962SLisandro Dalcin   if (scall == MAT_INITIAL_MATRIX) {
32892bb6962SLisandro Dalcin     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
32992bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
3303bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr);
33178904715SLisandro Dalcin       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
33292bb6962SLisandro Dalcin     }
33392bb6962SLisandro Dalcin   }
33480ec0b7dSPatrick Sanan 
33580ec0b7dSPatrick Sanan   /* Convert the types of the submatrices (if needbe) */
33680ec0b7dSPatrick Sanan   if (osm->sub_mat_type) {
33780ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; i++) {
33880ec0b7dSPatrick Sanan       ierr = MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]));CHKERRQ(ierr);
33980ec0b7dSPatrick Sanan     }
34080ec0b7dSPatrick Sanan   }
34180ec0b7dSPatrick Sanan 
34280ec0b7dSPatrick Sanan   if(!pc->setupcalled){
34380ec0b7dSPatrick Sanan     /* Create the local work vectors (from the local matrices) and scatter contexts */
34480ec0b7dSPatrick Sanan     ierr = MatCreateVecs(pc->pmat,&vec,0);CHKERRQ(ierr);
34580ec0b7dSPatrick Sanan     ierr = PetscMalloc1(osm->n_local,&osm->restriction);CHKERRQ(ierr);
346f1ee410cSBarry Smith     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()");
347f1ee410cSBarry Smith     if (osm->is_local && osm->type == PC_ASM_RESTRICT) {ierr = PetscMalloc1(osm->n_local,&osm->localization);CHKERRQ(ierr);}
34880ec0b7dSPatrick Sanan     ierr = PetscMalloc1(osm->n_local,&osm->prolongation);CHKERRQ(ierr);
34980ec0b7dSPatrick Sanan     ierr = PetscMalloc1(osm->n_local,&osm->x);CHKERRQ(ierr);
35080ec0b7dSPatrick Sanan     ierr = PetscMalloc1(osm->n_local,&osm->y);CHKERRQ(ierr);
35180ec0b7dSPatrick Sanan     ierr = PetscMalloc1(osm->n_local,&osm->y_local);CHKERRQ(ierr);
35280ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; ++i) {
35380ec0b7dSPatrick Sanan       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
35480ec0b7dSPatrick Sanan       ierr = MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);CHKERRQ(ierr);
35580ec0b7dSPatrick Sanan       ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
35680ec0b7dSPatrick Sanan       ierr = VecScatterCreate(vec,osm->is[i],osm->x[i],isl,&osm->restriction[i]);CHKERRQ(ierr);
35780ec0b7dSPatrick Sanan       ierr = ISDestroy(&isl);CHKERRQ(ierr);
35880ec0b7dSPatrick Sanan       ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr);
359f1ee410cSBarry Smith       if (osm->localization) {
36080ec0b7dSPatrick Sanan         ISLocalToGlobalMapping ltog;
36180ec0b7dSPatrick Sanan         IS                     isll;
36280ec0b7dSPatrick Sanan         const PetscInt         *idx_local;
36380ec0b7dSPatrick Sanan         PetscInt               *idx,nout;
36480ec0b7dSPatrick Sanan 
36580ec0b7dSPatrick Sanan         ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],&ltog);CHKERRQ(ierr);
36680ec0b7dSPatrick Sanan         ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr);
36780ec0b7dSPatrick Sanan         ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
36880ec0b7dSPatrick Sanan         ierr = PetscMalloc1(m_local,&idx);CHKERRQ(ierr);
36980ec0b7dSPatrick Sanan         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx);CHKERRQ(ierr);
37080ec0b7dSPatrick Sanan         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
37180ec0b7dSPatrick Sanan         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is");
37280ec0b7dSPatrick Sanan         ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
37380ec0b7dSPatrick Sanan         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
37480ec0b7dSPatrick Sanan         ierr = ISCreateStride(PETSC_COMM_SELF,m_local,0,1,&isl);CHKERRQ(ierr);
37580ec0b7dSPatrick Sanan         ierr = VecCreateSeq(PETSC_COMM_SELF,m_local,&osm->y_local[i]);CHKERRQ(ierr);
37680ec0b7dSPatrick Sanan         ierr = VecScatterCreate(osm->y[i],isll,osm->y_local[i],isl,&osm->localization[i]);CHKERRQ(ierr);
37780ec0b7dSPatrick Sanan         ierr = ISDestroy(&isll);CHKERRQ(ierr);
37880ec0b7dSPatrick Sanan 
37980ec0b7dSPatrick Sanan         ierr = VecScatterCreate(vec,osm->is_local[i],osm->y_local[i],isl,&osm->prolongation[i]);CHKERRQ(ierr);
38080ec0b7dSPatrick Sanan         ierr = ISDestroy(&isl);CHKERRQ(ierr);
38180ec0b7dSPatrick Sanan       } else {
38280ec0b7dSPatrick Sanan         osm->y_local[i] = osm->y[i];
38380ec0b7dSPatrick Sanan         ierr = PetscObjectReference((PetscObject) osm->y[i]);CHKERRQ(ierr);
38480ec0b7dSPatrick Sanan         osm->prolongation[i] = osm->restriction[i];
38580ec0b7dSPatrick Sanan         ierr = PetscObjectReference((PetscObject) osm->restriction[i]);CHKERRQ(ierr);
38680ec0b7dSPatrick Sanan       }
38780ec0b7dSPatrick Sanan     }
38880ec0b7dSPatrick Sanan     for (i=osm->n_local_true; i<osm->n_local; i++) {
38980ec0b7dSPatrick Sanan       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&osm->x[i]);CHKERRQ(ierr);
39080ec0b7dSPatrick Sanan       ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr);
39180ec0b7dSPatrick Sanan       ierr = VecDuplicate(osm->x[i],&osm->y_local[i]);CHKERRQ(ierr);
39280ec0b7dSPatrick Sanan       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&isl);CHKERRQ(ierr);
39380ec0b7dSPatrick Sanan       ierr = VecScatterCreate(vec,isl,osm->x[i],isl,&osm->restriction[i]);CHKERRQ(ierr);
394f1ee410cSBarry Smith       if (osm->localization) {
39580ec0b7dSPatrick Sanan         ierr = VecScatterCreate(osm->y[i],isl,osm->y_local[i],isl,&osm->localization[i]);CHKERRQ(ierr);
39680ec0b7dSPatrick Sanan         ierr = VecScatterCreate(vec,isl,osm->x[i],isl,&osm->prolongation[i]);CHKERRQ(ierr);
39780ec0b7dSPatrick Sanan       } else {
39880ec0b7dSPatrick Sanan         osm->prolongation[i] = osm->restriction[i];
39980ec0b7dSPatrick Sanan         ierr                 = PetscObjectReference((PetscObject) osm->restriction[i]);CHKERRQ(ierr);
40080ec0b7dSPatrick Sanan       }
40180ec0b7dSPatrick Sanan       ierr = ISDestroy(&isl);CHKERRQ(ierr);
40280ec0b7dSPatrick Sanan     }
40380ec0b7dSPatrick Sanan     ierr = VecDestroy(&vec);CHKERRQ(ierr);
40480ec0b7dSPatrick Sanan   }
40580ec0b7dSPatrick Sanan 
406fb745f2cSMatthew G. Knepley   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
407235cc4ceSMatthew G. Knepley     IS      *cis;
408235cc4ceSMatthew G. Knepley     PetscInt c;
409235cc4ceSMatthew G. Knepley 
410235cc4ceSMatthew G. Knepley     ierr = PetscMalloc1(osm->n_local_true, &cis);CHKERRQ(ierr);
411235cc4ceSMatthew G. Knepley     for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis;
4127dae84e0SHong Zhang     ierr = MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);CHKERRQ(ierr);
413235cc4ceSMatthew G. Knepley     ierr = PetscFree(cis);CHKERRQ(ierr);
414fb745f2cSMatthew G. Knepley   }
4154b9ad928SBarry Smith 
4164b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
4174b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
41892bb6962SLisandro Dalcin   ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
4194b9ad928SBarry Smith 
42092bb6962SLisandro Dalcin   /*
42192bb6962SLisandro Dalcin      Loop over subdomains putting them into local ksp
42292bb6962SLisandro Dalcin   */
42392bb6962SLisandro Dalcin   for (i=0; i<osm->n_local_true; i++) {
42423ee1639SBarry Smith     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr);
42592bb6962SLisandro Dalcin     if (!pc->setupcalled) {
426bf108f30SBarry Smith       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
4274b9ad928SBarry Smith     }
42892bb6962SLisandro Dalcin   }
4294b9ad928SBarry Smith   PetscFunctionReturn(0);
4304b9ad928SBarry Smith }
4314b9ad928SBarry Smith 
4326849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc)
4334b9ad928SBarry Smith {
4344b9ad928SBarry Smith   PC_ASM             *osm = (PC_ASM*)pc->data;
4356849ba73SBarry Smith   PetscErrorCode     ierr;
43613f74950SBarry Smith   PetscInt           i;
437539c167fSBarry Smith   KSPConvergedReason reason;
4384b9ad928SBarry Smith 
4394b9ad928SBarry Smith   PetscFunctionBegin;
4404b9ad928SBarry Smith   for (i=0; i<osm->n_local_true; i++) {
4414b9ad928SBarry Smith     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
442539c167fSBarry Smith     ierr = KSPGetConvergedReason(osm->ksp[i],&reason);CHKERRQ(ierr);
443539c167fSBarry Smith     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
444261222b2SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
445e0eafd54SHong Zhang     }
4464b9ad928SBarry Smith   }
4474b9ad928SBarry Smith   PetscFunctionReturn(0);
4484b9ad928SBarry Smith }
4494b9ad928SBarry Smith 
4506849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y)
4514b9ad928SBarry Smith {
4524b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
4536849ba73SBarry Smith   PetscErrorCode ierr;
45413f74950SBarry Smith   PetscInt       i,n_local = osm->n_local,n_local_true = osm->n_local_true;
4554b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
4564b9ad928SBarry Smith 
4574b9ad928SBarry Smith   PetscFunctionBegin;
4584b9ad928SBarry Smith   /*
4594b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
4604b9ad928SBarry Smith      subdomain values (leaving the other values 0).
4614b9ad928SBarry Smith   */
4624b9ad928SBarry Smith   if (!(osm->type & PC_ASM_RESTRICT)) {
4634b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
4644b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
465f21eeb63SBarry Smith     for (i=0; i<n_local_true; i++) {
4663d03bb48SJed Brown       ierr = VecZeroEntries(osm->x[i]);CHKERRQ(ierr);
4674b9ad928SBarry Smith     }
4684b9ad928SBarry Smith   }
4692fa5cd67SKarl Rupp   if (!(osm->type & PC_ASM_INTERPOLATE)) reverse = SCATTER_REVERSE_LOCAL;
4704b9ad928SBarry Smith 
47112cd4985SMatthew G. Knepley   switch (osm->loctype)
47212cd4985SMatthew G. Knepley   {
47312cd4985SMatthew G. Knepley   case PC_COMPOSITE_ADDITIVE:
4744b9ad928SBarry Smith     for (i=0; i<n_local; i++) {
4752b691e39SMatthew Knepley       ierr = VecScatterBegin(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
4764b9ad928SBarry Smith     }
47710bd88b9SJed Brown     ierr = VecZeroEntries(y);CHKERRQ(ierr);
4784b9ad928SBarry Smith     /* do the local solves */
4794b9ad928SBarry Smith     for (i=0; i<n_local_true; i++) {
4802b691e39SMatthew Knepley       ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
48123ce1328SBarry Smith       ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
482b9845e0eSMatthew Knepley       if (osm->localization) {
48353888de8SMatthew Knepley         ierr = VecScatterBegin(osm->localization[i],osm->y[i],osm->y_local[i],INSERT_VALUES,forward);CHKERRQ(ierr);
48453888de8SMatthew Knepley         ierr = VecScatterEnd(osm->localization[i],osm->y[i],osm->y_local[i],INSERT_VALUES,forward);CHKERRQ(ierr);
485b9845e0eSMatthew Knepley       }
48653888de8SMatthew Knepley       ierr = VecScatterBegin(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
4874b9ad928SBarry Smith     }
4884b9ad928SBarry Smith     /* handle the rest of the scatters that do not have local solves */
4894b9ad928SBarry Smith     for (i=n_local_true; i<n_local; i++) {
4902b691e39SMatthew Knepley       ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
49153888de8SMatthew Knepley       ierr = VecScatterBegin(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
4924b9ad928SBarry Smith     }
4934b9ad928SBarry Smith     for (i=0; i<n_local; i++) {
49453888de8SMatthew Knepley       ierr = VecScatterEnd(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
4954b9ad928SBarry Smith     }
49612cd4985SMatthew G. Knepley     break;
49712cd4985SMatthew G. Knepley   case PC_COMPOSITE_MULTIPLICATIVE:
49812cd4985SMatthew G. Knepley     ierr = VecZeroEntries(y);CHKERRQ(ierr);
49912cd4985SMatthew G. Knepley     /* do the local solves */
50012cd4985SMatthew G. Knepley     for (i = 0; i < n_local_true; ++i) {
50112cd4985SMatthew G. Knepley       if (i > 0) {
5027c3d802fSMatthew G. Knepley         /* Update rhs */
5037c3d802fSMatthew G. Knepley         ierr = VecScatterBegin(osm->restriction[i], osm->lx, osm->x[i], INSERT_VALUES, forward);CHKERRQ(ierr);
5047c3d802fSMatthew G. Knepley         ierr = VecScatterEnd(osm->restriction[i], osm->lx, osm->x[i], INSERT_VALUES, forward);CHKERRQ(ierr);
50512cd4985SMatthew G. Knepley       } else {
50612cd4985SMatthew G. Knepley         ierr = VecZeroEntries(osm->x[i]);CHKERRQ(ierr);
50712cd4985SMatthew G. Knepley       }
50812cd4985SMatthew G. Knepley       ierr = VecScatterBegin(osm->restriction[i], x, osm->x[i], ADD_VALUES, forward);CHKERRQ(ierr);
50912cd4985SMatthew G. Knepley       ierr = VecScatterEnd(osm->restriction[i], x, osm->x[i], ADD_VALUES, forward);CHKERRQ(ierr);
51012cd4985SMatthew G. Knepley       ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
51112cd4985SMatthew G. Knepley       if (osm->localization) {
51212cd4985SMatthew G. Knepley         ierr = VecScatterBegin(osm->localization[i], osm->y[i], osm->y_local[i], INSERT_VALUES, forward);CHKERRQ(ierr);
51312cd4985SMatthew G. Knepley         ierr = VecScatterEnd(osm->localization[i], osm->y[i], osm->y_local[i], INSERT_VALUES, forward);CHKERRQ(ierr);
51412cd4985SMatthew G. Knepley       }
51512cd4985SMatthew G. Knepley       ierr = VecScatterBegin(osm->prolongation[i], osm->y_local[i], y, ADD_VALUES, reverse);CHKERRQ(ierr);
51612cd4985SMatthew G. Knepley       ierr = VecScatterEnd(osm->prolongation[i], osm->y_local[i], y, ADD_VALUES, reverse);CHKERRQ(ierr);
5177c3d802fSMatthew G. Knepley       if (i < n_local_true-1) {
5187c3d802fSMatthew G. Knepley         ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
5197c3d802fSMatthew G. Knepley         ierr = VecScatterBegin(osm->prolongation[i], osm->y_local[i], osm->ly, INSERT_VALUES, reverse);CHKERRQ(ierr);
5207c3d802fSMatthew G. Knepley         ierr = VecScatterEnd(osm->prolongation[i], osm->y_local[i], osm->ly, INSERT_VALUES, reverse);CHKERRQ(ierr);
5217c3d802fSMatthew G. Knepley         ierr = VecScale(osm->ly, -1.0);CHKERRQ(ierr);
522235cc4ceSMatthew G. Knepley         ierr = MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]);CHKERRQ(ierr);
523235cc4ceSMatthew G. Knepley         ierr = VecScatterBegin(osm->restriction[i+1], osm->y[i+1], osm->lx, INSERT_VALUES, reverse);CHKERRQ(ierr);
524235cc4ceSMatthew G. Knepley         ierr = VecScatterEnd(osm->restriction[i+1], osm->y[i+1], osm->lx, INSERT_VALUES, reverse);CHKERRQ(ierr);
5257c3d802fSMatthew G. Knepley       }
52612cd4985SMatthew G. Knepley     }
52712cd4985SMatthew G. Knepley     /* handle the rest of the scatters that do not have local solves */
52812cd4985SMatthew G. Knepley     for (i = n_local_true; i < n_local; ++i) {
52912cd4985SMatthew G. Knepley       ierr = VecScatterBegin(osm->restriction[i], x, osm->x[i], INSERT_VALUES, forward);CHKERRQ(ierr);
53012cd4985SMatthew G. Knepley       ierr = VecScatterEnd(osm->restriction[i], x, osm->x[i], INSERT_VALUES, forward);CHKERRQ(ierr);
53112cd4985SMatthew G. Knepley       ierr = VecScatterBegin(osm->prolongation[i], osm->y_local[i], y, ADD_VALUES, reverse);CHKERRQ(ierr);
53212cd4985SMatthew G. Knepley       ierr = VecScatterEnd(osm->prolongation[i], osm->y_local[i], y, ADD_VALUES, reverse);CHKERRQ(ierr);
53312cd4985SMatthew G. Knepley     }
53412cd4985SMatthew G. Knepley     break;
53512cd4985SMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
53612cd4985SMatthew G. Knepley   }
5374b9ad928SBarry Smith   PetscFunctionReturn(0);
5384b9ad928SBarry Smith }
5394b9ad928SBarry Smith 
5406849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y)
5414b9ad928SBarry Smith {
5424b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
5436849ba73SBarry Smith   PetscErrorCode ierr;
54413f74950SBarry Smith   PetscInt       i,n_local = osm->n_local,n_local_true = osm->n_local_true;
5454b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
5464b9ad928SBarry Smith 
5474b9ad928SBarry Smith   PetscFunctionBegin;
5484b9ad928SBarry Smith   /*
5494b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
5504b9ad928SBarry Smith      subdomain values (leaving the other values 0).
5514b9ad928SBarry Smith 
5524b9ad928SBarry Smith      Note: these are reversed from the PCApply_ASM() because we are applying the
5534b9ad928SBarry Smith      transpose of the three terms
5544b9ad928SBarry Smith   */
5554b9ad928SBarry Smith   if (!(osm->type & PC_ASM_INTERPOLATE)) {
5564b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
5574b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
558f21eeb63SBarry Smith     for (i=0; i<n_local_true; i++) {
55910bd88b9SJed Brown       ierr = VecZeroEntries(osm->x[i]);CHKERRQ(ierr);
5604b9ad928SBarry Smith     }
5614b9ad928SBarry Smith   }
5622fa5cd67SKarl Rupp   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;
5634b9ad928SBarry Smith 
5644b9ad928SBarry Smith   for (i=0; i<n_local; i++) {
5652b691e39SMatthew Knepley     ierr = VecScatterBegin(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
5664b9ad928SBarry Smith   }
56710bd88b9SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
5684b9ad928SBarry Smith   /* do the local solves */
5694b9ad928SBarry Smith   for (i=0; i<n_local_true; i++) {
5702b691e39SMatthew Knepley     ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
57123ce1328SBarry Smith     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
572b9845e0eSMatthew Knepley     if (osm->localization) {
57353888de8SMatthew Knepley       ierr = VecScatterBegin(osm->localization[i],osm->y[i],osm->y_local[i],INSERT_VALUES,forward);CHKERRQ(ierr);
57453888de8SMatthew Knepley       ierr = VecScatterEnd(osm->localization[i],osm->y[i],osm->y_local[i],INSERT_VALUES,forward);CHKERRQ(ierr);
575b9845e0eSMatthew Knepley     }
57653888de8SMatthew Knepley     ierr = VecScatterBegin(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
5774b9ad928SBarry Smith   }
5784b9ad928SBarry Smith   /* handle the rest of the scatters that do not have local solves */
5794b9ad928SBarry Smith   for (i=n_local_true; i<n_local; i++) {
5802b691e39SMatthew Knepley     ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
58153888de8SMatthew Knepley     ierr = VecScatterBegin(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
5824b9ad928SBarry Smith   }
5834b9ad928SBarry Smith   for (i=0; i<n_local; i++) {
58453888de8SMatthew Knepley     ierr = VecScatterEnd(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
5854b9ad928SBarry Smith   }
5864b9ad928SBarry Smith   PetscFunctionReturn(0);
5874b9ad928SBarry Smith }
5884b9ad928SBarry Smith 
589e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc)
5904b9ad928SBarry Smith {
5914b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
5926849ba73SBarry Smith   PetscErrorCode ierr;
59313f74950SBarry Smith   PetscInt       i;
5944b9ad928SBarry Smith 
5954b9ad928SBarry Smith   PetscFunctionBegin;
59692bb6962SLisandro Dalcin   if (osm->ksp) {
59792bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
598e91c6855SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
59992bb6962SLisandro Dalcin     }
60092bb6962SLisandro Dalcin   }
601e09e08ccSBarry Smith   if (osm->pmat) {
60292bb6962SLisandro Dalcin     if (osm->n_local_true > 0) {
60330a70a9aSHong Zhang       ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
60492bb6962SLisandro Dalcin     }
60592bb6962SLisandro Dalcin   }
6062b691e39SMatthew Knepley   if (osm->restriction) {
6074b9ad928SBarry Smith     for (i=0; i<osm->n_local; i++) {
608fcfd50ebSBarry Smith       ierr = VecScatterDestroy(&osm->restriction[i]);CHKERRQ(ierr);
609fcfd50ebSBarry Smith       if (osm->localization) {ierr = VecScatterDestroy(&osm->localization[i]);CHKERRQ(ierr);}
610fcfd50ebSBarry Smith       ierr = VecScatterDestroy(&osm->prolongation[i]);CHKERRQ(ierr);
611fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
612fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
613fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y_local[i]);CHKERRQ(ierr);
6144b9ad928SBarry Smith     }
6152b691e39SMatthew Knepley     ierr = PetscFree(osm->restriction);CHKERRQ(ierr);
616b9845e0eSMatthew Knepley     if (osm->localization) {ierr = PetscFree(osm->localization);CHKERRQ(ierr);}
6172b691e39SMatthew Knepley     ierr = PetscFree(osm->prolongation);CHKERRQ(ierr);
61805b42c5fSBarry Smith     ierr = PetscFree(osm->x);CHKERRQ(ierr);
61978904715SLisandro Dalcin     ierr = PetscFree(osm->y);CHKERRQ(ierr);
62053888de8SMatthew Knepley     ierr = PetscFree(osm->y_local);CHKERRQ(ierr);
62192bb6962SLisandro Dalcin   }
6222b691e39SMatthew Knepley   ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
623fb745f2cSMatthew G. Knepley   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
624fb745f2cSMatthew G. Knepley     ierr = ISDestroy(&osm->lis);CHKERRQ(ierr);
625235cc4ceSMatthew G. Knepley     ierr = MatDestroyMatrices(osm->n_local_true, &osm->lmats);CHKERRQ(ierr);
626fb745f2cSMatthew G. Knepley     ierr = VecDestroy(&osm->lx);CHKERRQ(ierr);
627fb745f2cSMatthew G. Knepley     ierr = VecDestroy(&osm->ly);CHKERRQ(ierr);
628fb745f2cSMatthew G. Knepley   }
6292fa5cd67SKarl Rupp 
63080ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
63180ec0b7dSPatrick Sanan 
632e91c6855SBarry Smith   osm->is       = 0;
633e91c6855SBarry Smith   osm->is_local = 0;
634e91c6855SBarry Smith   PetscFunctionReturn(0);
635e91c6855SBarry Smith }
636e91c6855SBarry Smith 
637e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc)
638e91c6855SBarry Smith {
639e91c6855SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
640e91c6855SBarry Smith   PetscErrorCode ierr;
641e91c6855SBarry Smith   PetscInt       i;
642e91c6855SBarry Smith 
643e91c6855SBarry Smith   PetscFunctionBegin;
644e91c6855SBarry Smith   ierr = PCReset_ASM(pc);CHKERRQ(ierr);
645e91c6855SBarry Smith   if (osm->ksp) {
646e91c6855SBarry Smith     for (i=0; i<osm->n_local_true; i++) {
647fcfd50ebSBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
648e91c6855SBarry Smith     }
649e91c6855SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
650e91c6855SBarry Smith   }
651e91c6855SBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
6524b9ad928SBarry Smith   PetscFunctionReturn(0);
6534b9ad928SBarry Smith }
6544b9ad928SBarry Smith 
6554416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc)
6564b9ad928SBarry Smith {
6574b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
6586849ba73SBarry Smith   PetscErrorCode ierr;
6599dcbbd2bSBarry Smith   PetscInt       blocks,ovl;
660ace3abfcSBarry Smith   PetscBool      symset,flg;
66192bb6962SLisandro Dalcin   PCASMType      asmtype;
66212cd4985SMatthew G. Knepley   PCCompositeType loctype;
66380ec0b7dSPatrick Sanan   char           sub_mat_type[256];
6644b9ad928SBarry Smith 
6654b9ad928SBarry Smith   PetscFunctionBegin;
666bf108f30SBarry Smith   /* set the type to symmetric if matrix is symmetric */
66792bb6962SLisandro Dalcin   if (!osm->type_set && pc->pmat) {
66892bb6962SLisandro Dalcin     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
6692fa5cd67SKarl Rupp     if (symset && flg) osm->type = PC_ASM_BASIC;
670bf108f30SBarry Smith   }
671e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");CHKERRQ(ierr);
672d709ea83SDmitry Karpeev   ierr = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr);
67390d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
67465db9045SDmitry Karpeev   if (flg) {
675f77a5242SKarl Rupp     ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr);
676d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
67765db9045SDmitry Karpeev   }
67890d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
67965db9045SDmitry Karpeev   if (flg) {
68065db9045SDmitry Karpeev     ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr);
681d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
68265db9045SDmitry Karpeev   }
68390d69ab7SBarry Smith   flg  = PETSC_FALSE;
68490d69ab7SBarry Smith   ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
68592bb6962SLisandro Dalcin   if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); }
68612cd4985SMatthew G. Knepley   flg  = PETSC_FALSE;
68712cd4985SMatthew G. Knepley   ierr = PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr);
68812cd4985SMatthew G. Knepley   if (flg) {ierr = PCASMSetLocalType(pc,loctype);CHKERRQ(ierr); }
68980ec0b7dSPatrick Sanan   ierr = PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);CHKERRQ(ierr);
69080ec0b7dSPatrick Sanan   if(flg){
691459726d8SSatish Balay     ierr = PCASMSetSubMatType(pc,sub_mat_type);CHKERRQ(ierr);
69280ec0b7dSPatrick Sanan   }
6934b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
6944b9ad928SBarry Smith   PetscFunctionReturn(0);
6954b9ad928SBarry Smith }
6964b9ad928SBarry Smith 
6974b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/
6984b9ad928SBarry Smith 
6991e6b0712SBarry Smith static PetscErrorCode  PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[])
7004b9ad928SBarry Smith {
7014b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
70292bb6962SLisandro Dalcin   PetscErrorCode ierr;
70392bb6962SLisandro Dalcin   PetscInt       i;
7044b9ad928SBarry Smith 
7054b9ad928SBarry Smith   PetscFunctionBegin;
706e32f2f54SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
707ce94432eSBarry 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().");
708e7e72b3dSBarry Smith 
7094b9ad928SBarry Smith   if (!pc->setupcalled) {
71092bb6962SLisandro Dalcin     if (is) {
71192bb6962SLisandro Dalcin       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
71292bb6962SLisandro Dalcin     }
713832fc9a5SMatthew Knepley     if (is_local) {
714832fc9a5SMatthew Knepley       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
715832fc9a5SMatthew Knepley     }
7162b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
7172fa5cd67SKarl Rupp 
7184b9ad928SBarry Smith     osm->n_local_true = n;
71992bb6962SLisandro Dalcin     osm->is           = 0;
7202b691e39SMatthew Knepley     osm->is_local     = 0;
72192bb6962SLisandro Dalcin     if (is) {
722785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is);CHKERRQ(ierr);
7232fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is[i] = is[i];
7243d03bb48SJed Brown       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
7253d03bb48SJed Brown       osm->overlap = -1;
72692bb6962SLisandro Dalcin     }
7272b691e39SMatthew Knepley     if (is_local) {
728785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is_local);CHKERRQ(ierr);
7292fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is_local[i] = is_local[i];
730a35b7b57SDmitry Karpeev       if (!is) {
731785e854fSJed Brown         ierr = PetscMalloc1(osm->n_local_true,&osm->is);CHKERRQ(ierr);
732a35b7b57SDmitry Karpeev         for (i=0; i<osm->n_local_true; i++) {
733a35b7b57SDmitry Karpeev           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
734a35b7b57SDmitry Karpeev             ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr);
735a35b7b57SDmitry Karpeev             ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr);
736a35b7b57SDmitry Karpeev           } else {
737a35b7b57SDmitry Karpeev             ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr);
738a35b7b57SDmitry Karpeev             osm->is[i] = osm->is_local[i];
739a35b7b57SDmitry Karpeev           }
740a35b7b57SDmitry Karpeev         }
741a35b7b57SDmitry Karpeev       }
7422b691e39SMatthew Knepley     }
7434b9ad928SBarry Smith   }
7444b9ad928SBarry Smith   PetscFunctionReturn(0);
7454b9ad928SBarry Smith }
7464b9ad928SBarry Smith 
7471e6b0712SBarry Smith static PetscErrorCode  PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local)
7484b9ad928SBarry Smith {
7494b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
7506849ba73SBarry Smith   PetscErrorCode ierr;
75113f74950SBarry Smith   PetscMPIInt    rank,size;
75278904715SLisandro Dalcin   PetscInt       n;
7534b9ad928SBarry Smith 
7544b9ad928SBarry Smith   PetscFunctionBegin;
755ce94432eSBarry Smith   if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
756ce94432eSBarry 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.");
7574b9ad928SBarry Smith 
7584b9ad928SBarry Smith   /*
759880770e9SJed Brown      Split the subdomains equally among all processors
7604b9ad928SBarry Smith   */
761ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
762ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
7634b9ad928SBarry Smith   n    = N/size + ((N % size) > rank);
764e32f2f54SBarry 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);
765e32f2f54SBarry Smith   if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp().");
7664b9ad928SBarry Smith   if (!pc->setupcalled) {
7672b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
7682fa5cd67SKarl Rupp 
7694b9ad928SBarry Smith     osm->n_local_true = n;
7704b9ad928SBarry Smith     osm->is           = 0;
7712b691e39SMatthew Knepley     osm->is_local     = 0;
7724b9ad928SBarry Smith   }
7734b9ad928SBarry Smith   PetscFunctionReturn(0);
7744b9ad928SBarry Smith }
7754b9ad928SBarry Smith 
7761e6b0712SBarry Smith static PetscErrorCode  PCASMSetOverlap_ASM(PC pc,PetscInt ovl)
7774b9ad928SBarry Smith {
77892bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
7794b9ad928SBarry Smith 
7804b9ad928SBarry Smith   PetscFunctionBegin;
781ce94432eSBarry Smith   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
782ce94432eSBarry Smith   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp().");
7832fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
7844b9ad928SBarry Smith   PetscFunctionReturn(0);
7854b9ad928SBarry Smith }
7864b9ad928SBarry Smith 
7871e6b0712SBarry Smith static PetscErrorCode  PCASMSetType_ASM(PC pc,PCASMType type)
7884b9ad928SBarry Smith {
78992bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
7904b9ad928SBarry Smith 
7914b9ad928SBarry Smith   PetscFunctionBegin;
7924b9ad928SBarry Smith   osm->type     = type;
793bf108f30SBarry Smith   osm->type_set = PETSC_TRUE;
7944b9ad928SBarry Smith   PetscFunctionReturn(0);
7954b9ad928SBarry Smith }
7964b9ad928SBarry Smith 
797c60c7ad4SBarry Smith static PetscErrorCode  PCASMGetType_ASM(PC pc,PCASMType *type)
798c60c7ad4SBarry Smith {
799c60c7ad4SBarry Smith   PC_ASM *osm = (PC_ASM*)pc->data;
800c60c7ad4SBarry Smith 
801c60c7ad4SBarry Smith   PetscFunctionBegin;
802c60c7ad4SBarry Smith   *type = osm->type;
803c60c7ad4SBarry Smith   PetscFunctionReturn(0);
804c60c7ad4SBarry Smith }
805c60c7ad4SBarry Smith 
80612cd4985SMatthew G. Knepley static PetscErrorCode  PCASMSetLocalType_ASM(PC pc, PCCompositeType type)
80712cd4985SMatthew G. Knepley {
80812cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
80912cd4985SMatthew G. Knepley 
81012cd4985SMatthew G. Knepley   PetscFunctionBegin;
811*d0ecd4dfSBarry 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");
81212cd4985SMatthew G. Knepley   osm->loctype = type;
81312cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
81412cd4985SMatthew G. Knepley }
81512cd4985SMatthew G. Knepley 
81612cd4985SMatthew G. Knepley static PetscErrorCode  PCASMGetLocalType_ASM(PC pc, PCCompositeType *type)
81712cd4985SMatthew G. Knepley {
81812cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
81912cd4985SMatthew G. Knepley 
82012cd4985SMatthew G. Knepley   PetscFunctionBegin;
82112cd4985SMatthew G. Knepley   *type = osm->loctype;
82212cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
82312cd4985SMatthew G. Knepley }
82412cd4985SMatthew G. Knepley 
8251e6b0712SBarry Smith static PetscErrorCode  PCASMSetSortIndices_ASM(PC pc,PetscBool  doSort)
8266ed231c7SMatthew Knepley {
8276ed231c7SMatthew Knepley   PC_ASM *osm = (PC_ASM*)pc->data;
8286ed231c7SMatthew Knepley 
8296ed231c7SMatthew Knepley   PetscFunctionBegin;
8306ed231c7SMatthew Knepley   osm->sort_indices = doSort;
8316ed231c7SMatthew Knepley   PetscFunctionReturn(0);
8326ed231c7SMatthew Knepley }
8336ed231c7SMatthew Knepley 
8341e6b0712SBarry Smith static PetscErrorCode  PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
8354b9ad928SBarry Smith {
83692bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
837dfbe8321SBarry Smith   PetscErrorCode ierr;
8384b9ad928SBarry Smith 
8394b9ad928SBarry Smith   PetscFunctionBegin;
840ce94432eSBarry 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");
8414b9ad928SBarry Smith 
8422fa5cd67SKarl Rupp   if (n_local) *n_local = osm->n_local_true;
84392bb6962SLisandro Dalcin   if (first_local) {
844ce94432eSBarry Smith     ierr          = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
84592bb6962SLisandro Dalcin     *first_local -= osm->n_local_true;
84692bb6962SLisandro Dalcin   }
84792bb6962SLisandro Dalcin   if (ksp) {
84892bb6962SLisandro Dalcin     /* Assume that local solves are now different; not necessarily
84992bb6962SLisandro Dalcin        true though!  This flag is used only for PCView_ASM() */
85092bb6962SLisandro Dalcin     *ksp                   = osm->ksp;
85192bb6962SLisandro Dalcin     osm->same_local_solves = PETSC_FALSE;
85292bb6962SLisandro Dalcin   }
8534b9ad928SBarry Smith   PetscFunctionReturn(0);
8544b9ad928SBarry Smith }
8554b9ad928SBarry Smith 
85680ec0b7dSPatrick Sanan static PetscErrorCode  PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type)
85780ec0b7dSPatrick Sanan {
85880ec0b7dSPatrick Sanan   PC_ASM         *osm = (PC_ASM*)pc->data;
85980ec0b7dSPatrick Sanan 
86080ec0b7dSPatrick Sanan   PetscFunctionBegin;
86180ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
86280ec0b7dSPatrick Sanan   PetscValidPointer(sub_mat_type,2);
86380ec0b7dSPatrick Sanan   *sub_mat_type = osm->sub_mat_type;
86480ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
86580ec0b7dSPatrick Sanan }
86680ec0b7dSPatrick Sanan 
86780ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type)
86880ec0b7dSPatrick Sanan {
86980ec0b7dSPatrick Sanan   PetscErrorCode    ierr;
87080ec0b7dSPatrick Sanan   PC_ASM            *osm = (PC_ASM*)pc->data;
87180ec0b7dSPatrick Sanan 
87280ec0b7dSPatrick Sanan   PetscFunctionBegin;
87380ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
87480ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
87580ec0b7dSPatrick Sanan   ierr = PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);CHKERRQ(ierr);
87680ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
87780ec0b7dSPatrick Sanan }
87880ec0b7dSPatrick Sanan 
8794b9ad928SBarry Smith /*@C
8801093a601SBarry Smith     PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner.
8814b9ad928SBarry Smith 
8824b9ad928SBarry Smith     Collective on PC
8834b9ad928SBarry Smith 
8844b9ad928SBarry Smith     Input Parameters:
8854b9ad928SBarry Smith +   pc - the preconditioner context
8864b9ad928SBarry Smith .   n - the number of subdomains for this processor (default value = 1)
8878c03b21aSDmitry Karpeev .   is - the index set that defines the subdomains for this processor
8880298fd71SBarry Smith          (or NULL for PETSc to determine subdomains)
889f1ee410cSBarry 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
890f1ee410cSBarry Smith          (or NULL to not provide these)
8914b9ad928SBarry Smith 
8924b9ad928SBarry Smith     Notes:
8931093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
8944b9ad928SBarry Smith 
8954b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
8964b9ad928SBarry Smith 
8974b9ad928SBarry Smith     Use PCASMSetTotalSubdomains() to set the subdomains for all processors.
8984b9ad928SBarry Smith 
899f1ee410cSBarry Smith     If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated
900f1ee410cSBarry Smith     back to form the global solution (this is the standard restricted additive Schwarz method)
901f1ee410cSBarry Smith 
902f1ee410cSBarry 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
903f1ee410cSBarry Smith     no code to handle that case.
904f1ee410cSBarry Smith 
9054b9ad928SBarry Smith     Level: advanced
9064b9ad928SBarry Smith 
9074b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz
9084b9ad928SBarry Smith 
9094b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
910f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType()
9114b9ad928SBarry Smith @*/
9127087cfbeSBarry Smith PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
9134b9ad928SBarry Smith {
9147bb14e67SBarry Smith   PetscErrorCode ierr;
9154b9ad928SBarry Smith 
9164b9ad928SBarry Smith   PetscFunctionBegin;
9170700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9187bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr);
9194b9ad928SBarry Smith   PetscFunctionReturn(0);
9204b9ad928SBarry Smith }
9214b9ad928SBarry Smith 
9224b9ad928SBarry Smith /*@C
923feb221f8SDmitry Karpeev     PCASMSetTotalSubdomains - Sets the subdomains for all processors for the
9244b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
9254b9ad928SBarry Smith     PC communicator must call this routine, with the same index sets.
9264b9ad928SBarry Smith 
9274b9ad928SBarry Smith     Collective on PC
9284b9ad928SBarry Smith 
9294b9ad928SBarry Smith     Input Parameters:
9304b9ad928SBarry Smith +   pc - the preconditioner context
931feb221f8SDmitry Karpeev .   N  - the number of subdomains for all processors
932feb221f8SDmitry Karpeev .   is - the index sets that define the subdomains for all processors
933dfaa0c5fSBarry Smith          (or NULL to ask PETSc to determine the subdomains)
9342b691e39SMatthew Knepley -   is_local - the index sets that define the local part of the subdomains for this processor
935f1ee410cSBarry Smith          (or NULL to not provide this information)
9364b9ad928SBarry Smith 
9374b9ad928SBarry Smith     Options Database Key:
9384b9ad928SBarry Smith     To set the total number of subdomain blocks rather than specify the
9394b9ad928SBarry Smith     index sets, use the option
9404b9ad928SBarry Smith .    -pc_asm_blocks <blks> - Sets total blocks
9414b9ad928SBarry Smith 
9424b9ad928SBarry Smith     Notes:
943f1ee410cSBarry Smith     Currently you cannot use this to set the actual subdomains with the argument is or is_local.
9444b9ad928SBarry Smith 
9454b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
9464b9ad928SBarry Smith 
9474b9ad928SBarry Smith     These index sets cannot be destroyed until after completion of the
9484b9ad928SBarry Smith     linear solves for which the ASM preconditioner is being used.
9494b9ad928SBarry Smith 
9504b9ad928SBarry Smith     Use PCASMSetLocalSubdomains() to set local subdomains.
9514b9ad928SBarry Smith 
9521093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
9531093a601SBarry Smith 
9544b9ad928SBarry Smith     Level: advanced
9554b9ad928SBarry Smith 
9564b9ad928SBarry Smith .keywords: PC, ASM, set, total, global, subdomains, additive Schwarz
9574b9ad928SBarry Smith 
9584b9ad928SBarry Smith .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
9594b9ad928SBarry Smith           PCASMCreateSubdomains2D()
9604b9ad928SBarry Smith @*/
9617087cfbeSBarry Smith PetscErrorCode  PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[])
9624b9ad928SBarry Smith {
9637bb14e67SBarry Smith   PetscErrorCode ierr;
9644b9ad928SBarry Smith 
9654b9ad928SBarry Smith   PetscFunctionBegin;
9660700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9677bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr);
9684b9ad928SBarry Smith   PetscFunctionReturn(0);
9694b9ad928SBarry Smith }
9704b9ad928SBarry Smith 
9714b9ad928SBarry Smith /*@
9724b9ad928SBarry Smith     PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
9734b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
974f1ee410cSBarry Smith     PC communicator must call this routine.
9754b9ad928SBarry Smith 
976ad4df100SBarry Smith     Logically Collective on PC
9774b9ad928SBarry Smith 
9784b9ad928SBarry Smith     Input Parameters:
9794b9ad928SBarry Smith +   pc  - the preconditioner context
9804b9ad928SBarry Smith -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
9814b9ad928SBarry Smith 
9824b9ad928SBarry Smith     Options Database Key:
9834b9ad928SBarry Smith .   -pc_asm_overlap <ovl> - Sets overlap
9844b9ad928SBarry Smith 
9854b9ad928SBarry Smith     Notes:
9864b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.  To use
9874b9ad928SBarry Smith     multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
9884b9ad928SBarry Smith     PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).
9894b9ad928SBarry Smith 
9904b9ad928SBarry Smith     The overlap defaults to 1, so if one desires that no additional
9914b9ad928SBarry Smith     overlap be computed beyond what may have been set with a call to
9924b9ad928SBarry Smith     PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
9934b9ad928SBarry Smith     must be set to be 0.  In particular, if one does not explicitly set
9944b9ad928SBarry Smith     the subdomains an application code, then all overlap would be computed
9954b9ad928SBarry Smith     internally by PETSc, and using an overlap of 0 would result in an ASM
9964b9ad928SBarry Smith     variant that is equivalent to the block Jacobi preconditioner.
9974b9ad928SBarry Smith 
998f1ee410cSBarry Smith     The default algorithm used by PETSc to increase overlap is fast, but not scalable,
999f1ee410cSBarry Smith     use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
1000f1ee410cSBarry Smith 
10014b9ad928SBarry Smith     Note that one can define initial index sets with any overlap via
1002f1ee410cSBarry Smith     PCASMSetLocalSubdomains(); the routine
10034b9ad928SBarry Smith     PCASMSetOverlap() merely allows PETSc to extend that overlap further
10044b9ad928SBarry Smith     if desired.
10054b9ad928SBarry Smith 
10064b9ad928SBarry Smith     Level: intermediate
10074b9ad928SBarry Smith 
10084b9ad928SBarry Smith .keywords: PC, ASM, set, overlap
10094b9ad928SBarry Smith 
10104b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
1011f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap()
10124b9ad928SBarry Smith @*/
10137087cfbeSBarry Smith PetscErrorCode  PCASMSetOverlap(PC pc,PetscInt ovl)
10144b9ad928SBarry Smith {
10157bb14e67SBarry Smith   PetscErrorCode ierr;
10164b9ad928SBarry Smith 
10174b9ad928SBarry Smith   PetscFunctionBegin;
10180700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1019c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,ovl,2);
10207bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
10214b9ad928SBarry Smith   PetscFunctionReturn(0);
10224b9ad928SBarry Smith }
10234b9ad928SBarry Smith 
10244b9ad928SBarry Smith /*@
10254b9ad928SBarry Smith     PCASMSetType - Sets the type of restriction and interpolation used
10264b9ad928SBarry Smith     for local problems in the additive Schwarz method.
10274b9ad928SBarry Smith 
1028ad4df100SBarry Smith     Logically Collective on PC
10294b9ad928SBarry Smith 
10304b9ad928SBarry Smith     Input Parameters:
10314b9ad928SBarry Smith +   pc  - the preconditioner context
10324b9ad928SBarry Smith -   type - variant of ASM, one of
10334b9ad928SBarry Smith .vb
10344b9ad928SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
10354b9ad928SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
10364b9ad928SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
10374b9ad928SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
10384b9ad928SBarry Smith .ve
10394b9ad928SBarry Smith 
10404b9ad928SBarry Smith     Options Database Key:
10414b9ad928SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
10424b9ad928SBarry Smith 
1043f1ee410cSBarry Smith     Notes: if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected
1044f1ee410cSBarry Smith     to limit the local processor interpolation
1045f1ee410cSBarry Smith 
10464b9ad928SBarry Smith     Level: intermediate
10474b9ad928SBarry Smith 
10484b9ad928SBarry Smith .keywords: PC, ASM, set, type
10494b9ad928SBarry Smith 
10504b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1051f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType()
10524b9ad928SBarry Smith @*/
10537087cfbeSBarry Smith PetscErrorCode  PCASMSetType(PC pc,PCASMType type)
10544b9ad928SBarry Smith {
10557bb14e67SBarry Smith   PetscErrorCode ierr;
10564b9ad928SBarry Smith 
10574b9ad928SBarry Smith   PetscFunctionBegin;
10580700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1059c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
10607bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr);
10614b9ad928SBarry Smith   PetscFunctionReturn(0);
10624b9ad928SBarry Smith }
10634b9ad928SBarry Smith 
1064c60c7ad4SBarry Smith /*@
1065c60c7ad4SBarry Smith     PCASMGetType - Gets the type of restriction and interpolation used
1066c60c7ad4SBarry Smith     for local problems in the additive Schwarz method.
1067c60c7ad4SBarry Smith 
1068c60c7ad4SBarry Smith     Logically Collective on PC
1069c60c7ad4SBarry Smith 
1070c60c7ad4SBarry Smith     Input Parameter:
1071c60c7ad4SBarry Smith .   pc  - the preconditioner context
1072c60c7ad4SBarry Smith 
1073c60c7ad4SBarry Smith     Output Parameter:
1074c60c7ad4SBarry Smith .   type - variant of ASM, one of
1075c60c7ad4SBarry Smith 
1076c60c7ad4SBarry Smith .vb
1077c60c7ad4SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
1078c60c7ad4SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1079c60c7ad4SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1080c60c7ad4SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
1081c60c7ad4SBarry Smith .ve
1082c60c7ad4SBarry Smith 
1083c60c7ad4SBarry Smith     Options Database Key:
1084c60c7ad4SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
1085c60c7ad4SBarry Smith 
1086c60c7ad4SBarry Smith     Level: intermediate
1087c60c7ad4SBarry Smith 
1088c60c7ad4SBarry Smith .keywords: PC, ASM, set, type
1089c60c7ad4SBarry Smith 
1090c60c7ad4SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1091f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType()
1092c60c7ad4SBarry Smith @*/
1093c60c7ad4SBarry Smith PetscErrorCode  PCASMGetType(PC pc,PCASMType *type)
1094c60c7ad4SBarry Smith {
1095c60c7ad4SBarry Smith   PetscErrorCode ierr;
1096c60c7ad4SBarry Smith 
1097c60c7ad4SBarry Smith   PetscFunctionBegin;
1098c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1099c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));CHKERRQ(ierr);
1100c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1101c60c7ad4SBarry Smith }
1102c60c7ad4SBarry Smith 
110312cd4985SMatthew G. Knepley /*@
110412cd4985SMatthew G. Knepley   PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method.
110512cd4985SMatthew G. Knepley 
110612cd4985SMatthew G. Knepley   Logically Collective on PC
110712cd4985SMatthew G. Knepley 
110812cd4985SMatthew G. Knepley   Input Parameters:
110912cd4985SMatthew G. Knepley + pc  - the preconditioner context
111012cd4985SMatthew G. Knepley - type - type of composition, one of
111112cd4985SMatthew G. Knepley .vb
111212cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
111312cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
111412cd4985SMatthew G. Knepley .ve
111512cd4985SMatthew G. Knepley 
111612cd4985SMatthew G. Knepley   Options Database Key:
111712cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
111812cd4985SMatthew G. Knepley 
111912cd4985SMatthew G. Knepley   Level: intermediate
112012cd4985SMatthew G. Knepley 
1121f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
112212cd4985SMatthew G. Knepley @*/
112312cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type)
112412cd4985SMatthew G. Knepley {
112512cd4985SMatthew G. Knepley   PetscErrorCode ierr;
112612cd4985SMatthew G. Knepley 
112712cd4985SMatthew G. Knepley   PetscFunctionBegin;
112812cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
112912cd4985SMatthew G. Knepley   PetscValidLogicalCollectiveEnum(pc, type, 2);
113012cd4985SMatthew G. Knepley   ierr = PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));CHKERRQ(ierr);
113112cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
113212cd4985SMatthew G. Knepley }
113312cd4985SMatthew G. Knepley 
113412cd4985SMatthew G. Knepley /*@
113512cd4985SMatthew G. Knepley   PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method.
113612cd4985SMatthew G. Knepley 
113712cd4985SMatthew G. Knepley   Logically Collective on PC
113812cd4985SMatthew G. Knepley 
113912cd4985SMatthew G. Knepley   Input Parameter:
114012cd4985SMatthew G. Knepley . pc  - the preconditioner context
114112cd4985SMatthew G. Knepley 
114212cd4985SMatthew G. Knepley   Output Parameter:
114312cd4985SMatthew G. Knepley . type - type of composition, one of
114412cd4985SMatthew G. Knepley .vb
114512cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
114612cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
114712cd4985SMatthew G. Knepley .ve
114812cd4985SMatthew G. Knepley 
114912cd4985SMatthew G. Knepley   Options Database Key:
115012cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
115112cd4985SMatthew G. Knepley 
115212cd4985SMatthew G. Knepley   Level: intermediate
115312cd4985SMatthew G. Knepley 
1154f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
115512cd4985SMatthew G. Knepley @*/
115612cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type)
115712cd4985SMatthew G. Knepley {
115812cd4985SMatthew G. Knepley   PetscErrorCode ierr;
115912cd4985SMatthew G. Knepley 
116012cd4985SMatthew G. Knepley   PetscFunctionBegin;
116112cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
116212cd4985SMatthew G. Knepley   PetscValidPointer(type, 2);
116312cd4985SMatthew G. Knepley   ierr = PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));CHKERRQ(ierr);
116412cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
116512cd4985SMatthew G. Knepley }
116612cd4985SMatthew G. Knepley 
11676ed231c7SMatthew Knepley /*@
11686ed231c7SMatthew Knepley     PCASMSetSortIndices - Determines whether subdomain indices are sorted.
11696ed231c7SMatthew Knepley 
1170ad4df100SBarry Smith     Logically Collective on PC
11716ed231c7SMatthew Knepley 
11726ed231c7SMatthew Knepley     Input Parameters:
11736ed231c7SMatthew Knepley +   pc  - the preconditioner context
11746ed231c7SMatthew Knepley -   doSort - sort the subdomain indices
11756ed231c7SMatthew Knepley 
11766ed231c7SMatthew Knepley     Level: intermediate
11776ed231c7SMatthew Knepley 
11786ed231c7SMatthew Knepley .keywords: PC, ASM, set, type
11796ed231c7SMatthew Knepley 
11806ed231c7SMatthew Knepley .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
11816ed231c7SMatthew Knepley           PCASMCreateSubdomains2D()
11826ed231c7SMatthew Knepley @*/
11837087cfbeSBarry Smith PetscErrorCode  PCASMSetSortIndices(PC pc,PetscBool doSort)
11846ed231c7SMatthew Knepley {
11857bb14e67SBarry Smith   PetscErrorCode ierr;
11866ed231c7SMatthew Knepley 
11876ed231c7SMatthew Knepley   PetscFunctionBegin;
11880700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1189acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
11907bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
11916ed231c7SMatthew Knepley   PetscFunctionReturn(0);
11926ed231c7SMatthew Knepley }
11936ed231c7SMatthew Knepley 
11944b9ad928SBarry Smith /*@C
11954b9ad928SBarry Smith    PCASMGetSubKSP - Gets the local KSP contexts for all blocks on
11964b9ad928SBarry Smith    this processor.
11974b9ad928SBarry Smith 
11984b9ad928SBarry Smith    Collective on PC iff first_local is requested
11994b9ad928SBarry Smith 
12004b9ad928SBarry Smith    Input Parameter:
12014b9ad928SBarry Smith .  pc - the preconditioner context
12024b9ad928SBarry Smith 
12034b9ad928SBarry Smith    Output Parameters:
12040298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
12050298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
12060298fd71SBarry Smith                  all processors must request or all must pass NULL
12074b9ad928SBarry Smith -  ksp - the array of KSP contexts
12084b9ad928SBarry Smith 
12094b9ad928SBarry Smith    Note:
1210d29017ddSJed Brown    After PCASMGetSubKSP() the array of KSPes is not to be freed.
12114b9ad928SBarry Smith 
12124b9ad928SBarry Smith    You must call KSPSetUp() before calling PCASMGetSubKSP().
12134b9ad928SBarry Smith 
1214d29017ddSJed Brown    Fortran note:
12152bf68e3eSBarry 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.
1216d29017ddSJed Brown 
12174b9ad928SBarry Smith    Level: advanced
12184b9ad928SBarry Smith 
12194b9ad928SBarry Smith .keywords: PC, ASM, additive Schwarz, get, sub, KSP, context
12204b9ad928SBarry Smith 
12214b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
12224b9ad928SBarry Smith           PCASMCreateSubdomains2D(),
12234b9ad928SBarry Smith @*/
12247087cfbeSBarry Smith PetscErrorCode  PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
12254b9ad928SBarry Smith {
12267bb14e67SBarry Smith   PetscErrorCode ierr;
12274b9ad928SBarry Smith 
12284b9ad928SBarry Smith   PetscFunctionBegin;
12290700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12307bb14e67SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
12314b9ad928SBarry Smith   PetscFunctionReturn(0);
12324b9ad928SBarry Smith }
12334b9ad928SBarry Smith 
12344b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
12354b9ad928SBarry Smith /*MC
12363b09bd56SBarry Smith    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
12374b9ad928SBarry Smith            its own KSP object.
12384b9ad928SBarry Smith 
12394b9ad928SBarry Smith    Options Database Keys:
124049517cdeSBarry Smith +  -pc_asm_blocks <blks> - Sets total blocks
12414b9ad928SBarry Smith .  -pc_asm_overlap <ovl> - Sets overlap
1242f1ee410cSBarry Smith .  -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict
1243f1ee410cSBarry Smith -  -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive
12444b9ad928SBarry Smith 
12453b09bd56SBarry Smith      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
12463b09bd56SBarry Smith       will get a different convergence rate due to the default option of -pc_asm_type restrict. Use
12473b09bd56SBarry Smith       -pc_asm_type basic to use the standard ASM.
12483b09bd56SBarry Smith 
12494b9ad928SBarry Smith    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
1250f1ee410cSBarry Smith      than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor.
12514b9ad928SBarry Smith 
12523b09bd56SBarry Smith      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1253d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
12544b9ad928SBarry Smith 
1255a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCASMGetSubKSP()
12564b9ad928SBarry Smith          and set the options directly on the resulting KSP object (you can access its PC
12574b9ad928SBarry Smith          with KSPGetPC())
12584b9ad928SBarry Smith 
12594b9ad928SBarry Smith    Level: beginner
12604b9ad928SBarry Smith 
12614b9ad928SBarry Smith    Concepts: additive Schwarz method
12624b9ad928SBarry Smith 
1263c582cd25SBarry Smith     References:
126496a0c994SBarry Smith +   1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
126596a0c994SBarry Smith      Courant Institute, New York University Technical report
126696a0c994SBarry Smith -   1. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
126796a0c994SBarry Smith     Cambridge University Press.
1268c582cd25SBarry Smith 
12694b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1270f1ee410cSBarry Smith            PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType()
1271f1ee410cSBarry Smith            PCASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType
1272e09e08ccSBarry Smith 
12734b9ad928SBarry Smith M*/
12744b9ad928SBarry Smith 
12758cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc)
12764b9ad928SBarry Smith {
1277dfbe8321SBarry Smith   PetscErrorCode ierr;
12784b9ad928SBarry Smith   PC_ASM         *osm;
12794b9ad928SBarry Smith 
12804b9ad928SBarry Smith   PetscFunctionBegin;
1281b00a9115SJed Brown   ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr);
12822fa5cd67SKarl Rupp 
12834b9ad928SBarry Smith   osm->n                 = PETSC_DECIDE;
12844b9ad928SBarry Smith   osm->n_local           = 0;
12852b837212SDmitry Karpeev   osm->n_local_true      = PETSC_DECIDE;
12864b9ad928SBarry Smith   osm->overlap           = 1;
12874b9ad928SBarry Smith   osm->ksp               = 0;
12882b691e39SMatthew Knepley   osm->restriction       = 0;
1289b9845e0eSMatthew Knepley   osm->localization      = 0;
12902b691e39SMatthew Knepley   osm->prolongation      = 0;
129192bb6962SLisandro Dalcin   osm->x                 = 0;
129292bb6962SLisandro Dalcin   osm->y                 = 0;
129353888de8SMatthew Knepley   osm->y_local           = 0;
12944b9ad928SBarry Smith   osm->is                = 0;
12952b691e39SMatthew Knepley   osm->is_local          = 0;
12964b9ad928SBarry Smith   osm->mat               = 0;
12974b9ad928SBarry Smith   osm->pmat              = 0;
12984b9ad928SBarry Smith   osm->type              = PC_ASM_RESTRICT;
129912cd4985SMatthew G. Knepley   osm->loctype           = PC_COMPOSITE_ADDITIVE;
13004b9ad928SBarry Smith   osm->same_local_solves = PETSC_TRUE;
13016ed231c7SMatthew Knepley   osm->sort_indices      = PETSC_TRUE;
1302d709ea83SDmitry Karpeev   osm->dm_subdomains     = PETSC_FALSE;
130380ec0b7dSPatrick Sanan   osm->sub_mat_type      = NULL;
13044b9ad928SBarry Smith 
130592bb6962SLisandro Dalcin   pc->data                 = (void*)osm;
13064b9ad928SBarry Smith   pc->ops->apply           = PCApply_ASM;
13074b9ad928SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_ASM;
13084b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_ASM;
1309e91c6855SBarry Smith   pc->ops->reset           = PCReset_ASM;
13104b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_ASM;
13114b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_ASM;
13124b9ad928SBarry Smith   pc->ops->setuponblocks   = PCSetUpOnBlocks_ASM;
13134b9ad928SBarry Smith   pc->ops->view            = PCView_ASM;
13144b9ad928SBarry Smith   pc->ops->applyrichardson = 0;
13154b9ad928SBarry Smith 
1316bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr);
1317bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr);
1318bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);CHKERRQ(ierr);
1319bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);CHKERRQ(ierr);
1320c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);CHKERRQ(ierr);
132112cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);CHKERRQ(ierr);
132212cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);CHKERRQ(ierr);
1323bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);CHKERRQ(ierr);
1324bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);CHKERRQ(ierr);
132580ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);CHKERRQ(ierr);
132680ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);CHKERRQ(ierr);
13274b9ad928SBarry Smith   PetscFunctionReturn(0);
13284b9ad928SBarry Smith }
13294b9ad928SBarry Smith 
133092bb6962SLisandro Dalcin /*@C
133192bb6962SLisandro Dalcin    PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
133292bb6962SLisandro Dalcin    preconditioner for a any problem on a general grid.
133392bb6962SLisandro Dalcin 
133492bb6962SLisandro Dalcin    Collective
133592bb6962SLisandro Dalcin 
133692bb6962SLisandro Dalcin    Input Parameters:
133792bb6962SLisandro Dalcin +  A - The global matrix operator
133892bb6962SLisandro Dalcin -  n - the number of local blocks
133992bb6962SLisandro Dalcin 
134092bb6962SLisandro Dalcin    Output Parameters:
134192bb6962SLisandro Dalcin .  outis - the array of index sets defining the subdomains
134292bb6962SLisandro Dalcin 
134392bb6962SLisandro Dalcin    Level: advanced
134492bb6962SLisandro Dalcin 
13457d6bfa3bSBarry Smith    Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap
13467d6bfa3bSBarry Smith     from these if you use PCASMSetLocalSubdomains()
13477d6bfa3bSBarry Smith 
13487d6bfa3bSBarry Smith     In the Fortran version you must provide the array outis[] already allocated of length n.
13497d6bfa3bSBarry Smith 
135092bb6962SLisandro Dalcin .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid
135192bb6962SLisandro Dalcin 
135292bb6962SLisandro Dalcin .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains()
135392bb6962SLisandro Dalcin @*/
13547087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
135592bb6962SLisandro Dalcin {
135692bb6962SLisandro Dalcin   MatPartitioning mpart;
135792bb6962SLisandro Dalcin   const char      *prefix;
1358c56a70eeSBarry Smith   void            (*f)(void);
135992bb6962SLisandro Dalcin   PetscInt        i,j,rstart,rend,bs;
136011bd1e4dSLisandro Dalcin   PetscBool       isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
13610298fd71SBarry Smith   Mat             Ad     = NULL, adj;
136292bb6962SLisandro Dalcin   IS              ispart,isnumb,*is;
136392bb6962SLisandro Dalcin   PetscErrorCode  ierr;
136492bb6962SLisandro Dalcin 
136592bb6962SLisandro Dalcin   PetscFunctionBegin;
13660700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
136792bb6962SLisandro Dalcin   PetscValidPointer(outis,3);
1368c1235816SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
136992bb6962SLisandro Dalcin 
137092bb6962SLisandro Dalcin   /* Get prefix, row distribution, and block size */
137192bb6962SLisandro Dalcin   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
137292bb6962SLisandro Dalcin   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
137392bb6962SLisandro Dalcin   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
137465e19b50SBarry 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);
137565e19b50SBarry Smith 
137692bb6962SLisandro Dalcin   /* Get diagonal block from matrix if possible */
1377c56a70eeSBarry Smith   ierr = MatShellGetOperation(A,MATOP_GET_DIAGONAL_BLOCK,&f);CHKERRQ(ierr);
137892bb6962SLisandro Dalcin   if (f) {
137911bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
138092bb6962SLisandro Dalcin   }
138192bb6962SLisandro Dalcin   if (Ad) {
1382251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1383251f4c67SDmitry Karpeev     if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
138492bb6962SLisandro Dalcin   }
138592bb6962SLisandro Dalcin   if (Ad && n > 1) {
1386ace3abfcSBarry Smith     PetscBool match,done;
138792bb6962SLisandro Dalcin     /* Try to setup a good matrix partitioning if available */
138892bb6962SLisandro Dalcin     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
138992bb6962SLisandro Dalcin     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
139092bb6962SLisandro Dalcin     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1391251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
139292bb6962SLisandro Dalcin     if (!match) {
1393251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
139492bb6962SLisandro Dalcin     }
139592bb6962SLisandro Dalcin     if (!match) { /* assume a "good" partitioner is available */
13961a83f524SJed Brown       PetscInt       na;
13971a83f524SJed Brown       const PetscInt *ia,*ja;
139892bb6962SLisandro Dalcin       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
139992bb6962SLisandro Dalcin       if (done) {
140092bb6962SLisandro Dalcin         /* Build adjacency matrix by hand. Unfortunately a call to
140192bb6962SLisandro Dalcin            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
140292bb6962SLisandro Dalcin            remove the block-aij structure and we cannot expect
140392bb6962SLisandro Dalcin            MatPartitioning to split vertices as we need */
14041a83f524SJed Brown         PetscInt       i,j,len,nnz,cnt,*iia=0,*jja=0;
14051a83f524SJed Brown         const PetscInt *row;
140692bb6962SLisandro Dalcin         nnz = 0;
140792bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* count number of nonzeros */
140892bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
140992bb6962SLisandro Dalcin           row = ja + ia[i];
141092bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
141192bb6962SLisandro Dalcin             if (row[j] == i) { /* don't count diagonal */
141292bb6962SLisandro Dalcin               len--; break;
141392bb6962SLisandro Dalcin             }
141492bb6962SLisandro Dalcin           }
141592bb6962SLisandro Dalcin           nnz += len;
141692bb6962SLisandro Dalcin         }
1417854ce69bSBarry Smith         ierr   = PetscMalloc1(na+1,&iia);CHKERRQ(ierr);
1418854ce69bSBarry Smith         ierr   = PetscMalloc1(nnz,&jja);CHKERRQ(ierr);
141992bb6962SLisandro Dalcin         nnz    = 0;
142092bb6962SLisandro Dalcin         iia[0] = 0;
142192bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* fill adjacency */
142292bb6962SLisandro Dalcin           cnt = 0;
142392bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
142492bb6962SLisandro Dalcin           row = ja + ia[i];
142592bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
142692bb6962SLisandro Dalcin             if (row[j] != i) { /* if not diagonal */
142792bb6962SLisandro Dalcin               jja[nnz+cnt++] = row[j];
142892bb6962SLisandro Dalcin             }
142992bb6962SLisandro Dalcin           }
143092bb6962SLisandro Dalcin           nnz     += cnt;
143192bb6962SLisandro Dalcin           iia[i+1] = nnz;
143292bb6962SLisandro Dalcin         }
143392bb6962SLisandro Dalcin         /* Partitioning of the adjacency matrix */
14340298fd71SBarry Smith         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr);
143592bb6962SLisandro Dalcin         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
143692bb6962SLisandro Dalcin         ierr      = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
143792bb6962SLisandro Dalcin         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
143892bb6962SLisandro Dalcin         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
1439fcfd50ebSBarry Smith         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
144092bb6962SLisandro Dalcin         foundpart = PETSC_TRUE;
144192bb6962SLisandro Dalcin       }
144292bb6962SLisandro Dalcin       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
144392bb6962SLisandro Dalcin     }
1444fcfd50ebSBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
144592bb6962SLisandro Dalcin   }
144692bb6962SLisandro Dalcin 
1447785e854fSJed Brown   ierr   = PetscMalloc1(n,&is);CHKERRQ(ierr);
144892bb6962SLisandro Dalcin   *outis = is;
144992bb6962SLisandro Dalcin 
145092bb6962SLisandro Dalcin   if (!foundpart) {
145192bb6962SLisandro Dalcin 
145292bb6962SLisandro Dalcin     /* Partitioning by contiguous chunks of rows */
145392bb6962SLisandro Dalcin 
145492bb6962SLisandro Dalcin     PetscInt mbs   = (rend-rstart)/bs;
145592bb6962SLisandro Dalcin     PetscInt start = rstart;
145692bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
145792bb6962SLisandro Dalcin       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
145892bb6962SLisandro Dalcin       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
145992bb6962SLisandro Dalcin       start += count;
146092bb6962SLisandro Dalcin     }
146192bb6962SLisandro Dalcin 
146292bb6962SLisandro Dalcin   } else {
146392bb6962SLisandro Dalcin 
146492bb6962SLisandro Dalcin     /* Partitioning by adjacency of diagonal block  */
146592bb6962SLisandro Dalcin 
146692bb6962SLisandro Dalcin     const PetscInt *numbering;
146792bb6962SLisandro Dalcin     PetscInt       *count,nidx,*indices,*newidx,start=0;
146892bb6962SLisandro Dalcin     /* Get node count in each partition */
1469785e854fSJed Brown     ierr = PetscMalloc1(n,&count);CHKERRQ(ierr);
147092bb6962SLisandro Dalcin     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
147192bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
147292bb6962SLisandro Dalcin       for (i=0; i<n; i++) count[i] *= bs;
147392bb6962SLisandro Dalcin     }
147492bb6962SLisandro Dalcin     /* Build indices from node numbering */
147592bb6962SLisandro Dalcin     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1476785e854fSJed Brown     ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr);
147792bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
147892bb6962SLisandro Dalcin     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
147992bb6962SLisandro Dalcin     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
148092bb6962SLisandro Dalcin     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
148192bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1482785e854fSJed Brown       ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr);
14832fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
14842fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
14852fa5cd67SKarl Rupp       }
148692bb6962SLisandro Dalcin       ierr    = PetscFree(indices);CHKERRQ(ierr);
148792bb6962SLisandro Dalcin       nidx   *= bs;
148892bb6962SLisandro Dalcin       indices = newidx;
148992bb6962SLisandro Dalcin     }
149092bb6962SLisandro Dalcin     /* Shift to get global indices */
149192bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] += rstart;
149292bb6962SLisandro Dalcin 
149392bb6962SLisandro Dalcin     /* Build the index sets for each block */
149492bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
149570b3c8c7SBarry Smith       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
149692bb6962SLisandro Dalcin       ierr   = ISSort(is[i]);CHKERRQ(ierr);
149792bb6962SLisandro Dalcin       start += count[i];
149892bb6962SLisandro Dalcin     }
149992bb6962SLisandro Dalcin 
15003bf036e2SBarry Smith     ierr = PetscFree(count);CHKERRQ(ierr);
15013bf036e2SBarry Smith     ierr = PetscFree(indices);CHKERRQ(ierr);
1502fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1503fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
150492bb6962SLisandro Dalcin 
150592bb6962SLisandro Dalcin   }
150692bb6962SLisandro Dalcin   PetscFunctionReturn(0);
150792bb6962SLisandro Dalcin }
150892bb6962SLisandro Dalcin 
150992bb6962SLisandro Dalcin /*@C
151092bb6962SLisandro Dalcin    PCASMDestroySubdomains - Destroys the index sets created with
151192bb6962SLisandro Dalcin    PCASMCreateSubdomains(). Should be called after setting subdomains
151292bb6962SLisandro Dalcin    with PCASMSetLocalSubdomains().
151392bb6962SLisandro Dalcin 
151492bb6962SLisandro Dalcin    Collective
151592bb6962SLisandro Dalcin 
151692bb6962SLisandro Dalcin    Input Parameters:
151792bb6962SLisandro Dalcin +  n - the number of index sets
15182b691e39SMatthew Knepley .  is - the array of index sets
15190298fd71SBarry Smith -  is_local - the array of local index sets, can be NULL
152092bb6962SLisandro Dalcin 
152192bb6962SLisandro Dalcin    Level: advanced
152292bb6962SLisandro Dalcin 
152392bb6962SLisandro Dalcin .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid
152492bb6962SLisandro Dalcin 
152592bb6962SLisandro Dalcin .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains()
152692bb6962SLisandro Dalcin @*/
15277087cfbeSBarry Smith PetscErrorCode  PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
152892bb6962SLisandro Dalcin {
152992bb6962SLisandro Dalcin   PetscInt       i;
153092bb6962SLisandro Dalcin   PetscErrorCode ierr;
15315fd66863SKarl Rupp 
153292bb6962SLisandro Dalcin   PetscFunctionBegin;
1533a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1534a35b7b57SDmitry Karpeev   if (is) {
153592bb6962SLisandro Dalcin     PetscValidPointer(is,2);
1536fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); }
153792bb6962SLisandro Dalcin     ierr = PetscFree(is);CHKERRQ(ierr);
1538a35b7b57SDmitry Karpeev   }
15392b691e39SMatthew Knepley   if (is_local) {
15402b691e39SMatthew Knepley     PetscValidPointer(is_local,3);
1541fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); }
15422b691e39SMatthew Knepley     ierr = PetscFree(is_local);CHKERRQ(ierr);
15432b691e39SMatthew Knepley   }
154492bb6962SLisandro Dalcin   PetscFunctionReturn(0);
154592bb6962SLisandro Dalcin }
154692bb6962SLisandro Dalcin 
15474b9ad928SBarry Smith /*@
15484b9ad928SBarry Smith    PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
15494b9ad928SBarry Smith    preconditioner for a two-dimensional problem on a regular grid.
15504b9ad928SBarry Smith 
15514b9ad928SBarry Smith    Not Collective
15524b9ad928SBarry Smith 
15534b9ad928SBarry Smith    Input Parameters:
15544b9ad928SBarry Smith +  m, n - the number of mesh points in the x and y directions
15554b9ad928SBarry Smith .  M, N - the number of subdomains in the x and y directions
15564b9ad928SBarry Smith .  dof - degrees of freedom per node
15574b9ad928SBarry Smith -  overlap - overlap in mesh lines
15584b9ad928SBarry Smith 
15594b9ad928SBarry Smith    Output Parameters:
15604b9ad928SBarry Smith +  Nsub - the number of subdomains created
15613d03bb48SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
15623d03bb48SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
15634b9ad928SBarry Smith 
15644b9ad928SBarry Smith    Note:
15654b9ad928SBarry Smith    Presently PCAMSCreateSubdomains2d() is valid only for sequential
15664b9ad928SBarry Smith    preconditioners.  More general related routines are
15674b9ad928SBarry Smith    PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().
15684b9ad928SBarry Smith 
15694b9ad928SBarry Smith    Level: advanced
15704b9ad928SBarry Smith 
15714b9ad928SBarry Smith .keywords: PC, ASM, additive Schwarz, create, subdomains, 2D, regular grid
15724b9ad928SBarry Smith 
15734b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
15744b9ad928SBarry Smith           PCASMSetOverlap()
15754b9ad928SBarry Smith @*/
15767087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local)
15774b9ad928SBarry Smith {
15783d03bb48SJed Brown   PetscInt       i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer;
15796849ba73SBarry Smith   PetscErrorCode ierr;
158013f74950SBarry Smith   PetscInt       nidx,*idx,loc,ii,jj,count;
15814b9ad928SBarry Smith 
15824b9ad928SBarry Smith   PetscFunctionBegin;
1583e32f2f54SBarry Smith   if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," ");
15844b9ad928SBarry Smith 
15854b9ad928SBarry Smith   *Nsub     = N*M;
1586854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is);CHKERRQ(ierr);
1587854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is_local);CHKERRQ(ierr);
15884b9ad928SBarry Smith   ystart    = 0;
15893d03bb48SJed Brown   loc_outer = 0;
15904b9ad928SBarry Smith   for (i=0; i<N; i++) {
15914b9ad928SBarry Smith     height = n/N + ((n % N) > i); /* height of subdomain */
1592e32f2f54SBarry Smith     if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n");
15934b9ad928SBarry Smith     yleft  = ystart - overlap; if (yleft < 0) yleft = 0;
15944b9ad928SBarry Smith     yright = ystart + height + overlap; if (yright > n) yright = n;
15954b9ad928SBarry Smith     xstart = 0;
15964b9ad928SBarry Smith     for (j=0; j<M; j++) {
15974b9ad928SBarry Smith       width = m/M + ((m % M) > j); /* width of subdomain */
1598e32f2f54SBarry Smith       if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m");
15994b9ad928SBarry Smith       xleft  = xstart - overlap; if (xleft < 0) xleft = 0;
16004b9ad928SBarry Smith       xright = xstart + width + overlap; if (xright > m) xright = m;
16014b9ad928SBarry Smith       nidx   = (xright - xleft)*(yright - yleft);
1602785e854fSJed Brown       ierr   = PetscMalloc1(nidx,&idx);CHKERRQ(ierr);
16034b9ad928SBarry Smith       loc    = 0;
16044b9ad928SBarry Smith       for (ii=yleft; ii<yright; ii++) {
16054b9ad928SBarry Smith         count = m*ii + xleft;
16062fa5cd67SKarl Rupp         for (jj=xleft; jj<xright; jj++) idx[loc++] = count++;
16074b9ad928SBarry Smith       }
160870b3c8c7SBarry Smith       ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr);
16093d03bb48SJed Brown       if (overlap == 0) {
16103d03bb48SJed Brown         ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr);
16112fa5cd67SKarl Rupp 
16123d03bb48SJed Brown         (*is_local)[loc_outer] = (*is)[loc_outer];
16133d03bb48SJed Brown       } else {
16143d03bb48SJed Brown         for (loc=0,ii=ystart; ii<ystart+height; ii++) {
16153d03bb48SJed Brown           for (jj=xstart; jj<xstart+width; jj++) {
16163d03bb48SJed Brown             idx[loc++] = m*ii + jj;
16173d03bb48SJed Brown           }
16183d03bb48SJed Brown         }
161970b3c8c7SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr);
16203d03bb48SJed Brown       }
16214b9ad928SBarry Smith       ierr    = PetscFree(idx);CHKERRQ(ierr);
16224b9ad928SBarry Smith       xstart += width;
16233d03bb48SJed Brown       loc_outer++;
16244b9ad928SBarry Smith     }
16254b9ad928SBarry Smith     ystart += height;
16264b9ad928SBarry Smith   }
16274b9ad928SBarry Smith   for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); }
16284b9ad928SBarry Smith   PetscFunctionReturn(0);
16294b9ad928SBarry Smith }
16304b9ad928SBarry Smith 
16314b9ad928SBarry Smith /*@C
16324b9ad928SBarry Smith     PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
16334b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
16344b9ad928SBarry Smith 
1635ad4df100SBarry Smith     Not Collective
16364b9ad928SBarry Smith 
16374b9ad928SBarry Smith     Input Parameter:
16384b9ad928SBarry Smith .   pc - the preconditioner context
16394b9ad928SBarry Smith 
16404b9ad928SBarry Smith     Output Parameters:
16414b9ad928SBarry Smith +   n - the number of subdomains for this processor (default value = 1)
16422b691e39SMatthew Knepley .   is - the index sets that define the subdomains for this processor
16430298fd71SBarry Smith -   is_local - the index sets that define the local part of the subdomains for this processor (can be NULL)
16444b9ad928SBarry Smith 
16454b9ad928SBarry Smith 
16464b9ad928SBarry Smith     Notes:
16474b9ad928SBarry Smith     The IS numbering is in the parallel, global numbering of the vector.
16484b9ad928SBarry Smith 
16494b9ad928SBarry Smith     Level: advanced
16504b9ad928SBarry Smith 
16514b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz
16524b9ad928SBarry Smith 
16534b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
16544b9ad928SBarry Smith           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices()
16554b9ad928SBarry Smith @*/
16567087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
16574b9ad928SBarry Smith {
16582a808120SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
165992bb6962SLisandro Dalcin   PetscErrorCode ierr;
1660ace3abfcSBarry Smith   PetscBool      match;
16614b9ad928SBarry Smith 
16624b9ad928SBarry Smith   PetscFunctionBegin;
16630700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
16644482741eSBarry Smith   PetscValidIntPointer(n,2);
166592bb6962SLisandro Dalcin   if (is) PetscValidPointer(is,3);
1666251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
16672a808120SBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM");
16684b9ad928SBarry Smith   if (n) *n = osm->n_local_true;
16694b9ad928SBarry Smith   if (is) *is = osm->is;
16702b691e39SMatthew Knepley   if (is_local) *is_local = osm->is_local;
16714b9ad928SBarry Smith   PetscFunctionReturn(0);
16724b9ad928SBarry Smith }
16734b9ad928SBarry Smith 
16744b9ad928SBarry Smith /*@C
16754b9ad928SBarry Smith     PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
16764b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
16774b9ad928SBarry Smith 
1678ad4df100SBarry Smith     Not Collective
16794b9ad928SBarry Smith 
16804b9ad928SBarry Smith     Input Parameter:
16814b9ad928SBarry Smith .   pc - the preconditioner context
16824b9ad928SBarry Smith 
16834b9ad928SBarry Smith     Output Parameters:
16844b9ad928SBarry Smith +   n - the number of matrices for this processor (default value = 1)
16854b9ad928SBarry Smith -   mat - the matrices
16864b9ad928SBarry Smith 
16874b9ad928SBarry Smith     Level: advanced
16884b9ad928SBarry Smith 
1689cf739d55SBarry Smith     Notes: Call after PCSetUp() (or KSPSetUp()) but before PCApply() (or KSPApply()) and before PCSetUpOnBlocks())
1690cf739d55SBarry Smith 
1691cf739d55SBarry Smith            Usually one would use PCSetModifySubmatrices() to change the submatrices in building the preconditioner.
1692cf739d55SBarry Smith 
16934b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz, block Jacobi
16944b9ad928SBarry Smith 
16954b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1696cf739d55SBarry Smith           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubmatrices()
16974b9ad928SBarry Smith @*/
16987087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
16994b9ad928SBarry Smith {
17004b9ad928SBarry Smith   PC_ASM         *osm;
170192bb6962SLisandro Dalcin   PetscErrorCode ierr;
1702ace3abfcSBarry Smith   PetscBool      match;
17034b9ad928SBarry Smith 
17044b9ad928SBarry Smith   PetscFunctionBegin;
17050700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
170692bb6962SLisandro Dalcin   PetscValidIntPointer(n,2);
170792bb6962SLisandro Dalcin   if (mat) PetscValidPointer(mat,3);
1708ce94432eSBarry Smith   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1709251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
171092bb6962SLisandro Dalcin   if (!match) {
171192bb6962SLisandro Dalcin     if (n) *n = 0;
17120298fd71SBarry Smith     if (mat) *mat = NULL;
171392bb6962SLisandro Dalcin   } else {
17144b9ad928SBarry Smith     osm = (PC_ASM*)pc->data;
17154b9ad928SBarry Smith     if (n) *n = osm->n_local_true;
17164b9ad928SBarry Smith     if (mat) *mat = osm->pmat;
171792bb6962SLisandro Dalcin   }
17184b9ad928SBarry Smith   PetscFunctionReturn(0);
17194b9ad928SBarry Smith }
1720d709ea83SDmitry Karpeev 
1721d709ea83SDmitry Karpeev /*@
1722d709ea83SDmitry Karpeev     PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1723f1ee410cSBarry Smith 
1724d709ea83SDmitry Karpeev     Logically Collective
1725d709ea83SDmitry Karpeev 
1726d709ea83SDmitry Karpeev     Input Parameter:
1727d709ea83SDmitry Karpeev +   pc  - the preconditioner
1728d709ea83SDmitry Karpeev -   flg - boolean indicating whether to use subdomains defined by the DM
1729d709ea83SDmitry Karpeev 
1730d709ea83SDmitry Karpeev     Options Database Key:
1731d709ea83SDmitry Karpeev .   -pc_asm_dm_subdomains
1732d709ea83SDmitry Karpeev 
1733d709ea83SDmitry Karpeev     Level: intermediate
1734d709ea83SDmitry Karpeev 
1735d709ea83SDmitry Karpeev     Notes:
1736d709ea83SDmitry Karpeev     PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(),
1737d709ea83SDmitry Karpeev     so setting either of the first two effectively turns the latter off.
1738d709ea83SDmitry Karpeev 
1739d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1740d709ea83SDmitry Karpeev 
1741d709ea83SDmitry Karpeev .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1742d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1743d709ea83SDmitry Karpeev @*/
1744d709ea83SDmitry Karpeev PetscErrorCode  PCASMSetDMSubdomains(PC pc,PetscBool flg)
1745d709ea83SDmitry Karpeev {
1746d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1747d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1748d709ea83SDmitry Karpeev   PetscBool      match;
1749d709ea83SDmitry Karpeev 
1750d709ea83SDmitry Karpeev   PetscFunctionBegin;
1751d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1752d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
1753d709ea83SDmitry Karpeev   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1754d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1755d709ea83SDmitry Karpeev   if (match) {
1756d709ea83SDmitry Karpeev     osm->dm_subdomains = flg;
1757d709ea83SDmitry Karpeev   }
1758d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1759d709ea83SDmitry Karpeev }
1760d709ea83SDmitry Karpeev 
1761d709ea83SDmitry Karpeev /*@
1762d709ea83SDmitry Karpeev     PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1763d709ea83SDmitry Karpeev     Not Collective
1764d709ea83SDmitry Karpeev 
1765d709ea83SDmitry Karpeev     Input Parameter:
1766d709ea83SDmitry Karpeev .   pc  - the preconditioner
1767d709ea83SDmitry Karpeev 
1768d709ea83SDmitry Karpeev     Output Parameter:
1769d709ea83SDmitry Karpeev .   flg - boolean indicating whether to use subdomains defined by the DM
1770d709ea83SDmitry Karpeev 
1771d709ea83SDmitry Karpeev     Level: intermediate
1772d709ea83SDmitry Karpeev 
1773d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1774d709ea83SDmitry Karpeev 
1775d709ea83SDmitry Karpeev .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1776d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1777d709ea83SDmitry Karpeev @*/
1778d709ea83SDmitry Karpeev PetscErrorCode  PCASMGetDMSubdomains(PC pc,PetscBool* flg)
1779d709ea83SDmitry Karpeev {
1780d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1781d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1782d709ea83SDmitry Karpeev   PetscBool      match;
1783d709ea83SDmitry Karpeev 
1784d709ea83SDmitry Karpeev   PetscFunctionBegin;
1785d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1786d709ea83SDmitry Karpeev   PetscValidPointer(flg,2);
1787d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1788d709ea83SDmitry Karpeev   if (match) {
1789d709ea83SDmitry Karpeev     if (flg) *flg = osm->dm_subdomains;
1790d709ea83SDmitry Karpeev   }
1791d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1792d709ea83SDmitry Karpeev }
179380ec0b7dSPatrick Sanan 
179480ec0b7dSPatrick Sanan /*@
179580ec0b7dSPatrick Sanan      PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string.
179680ec0b7dSPatrick Sanan 
179780ec0b7dSPatrick Sanan    Not Collective
179880ec0b7dSPatrick Sanan 
179980ec0b7dSPatrick Sanan    Input Parameter:
180080ec0b7dSPatrick Sanan .  pc - the PC
180180ec0b7dSPatrick Sanan 
180280ec0b7dSPatrick Sanan    Output Parameter:
1803f1ee410cSBarry Smith .  -pc_asm_sub_mat_type - name of matrix type
180480ec0b7dSPatrick Sanan 
180580ec0b7dSPatrick Sanan    Level: advanced
180680ec0b7dSPatrick Sanan 
180780ec0b7dSPatrick Sanan .keywords: PC, PCASM, MatType, set
180880ec0b7dSPatrick Sanan 
180980ec0b7dSPatrick Sanan .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
181080ec0b7dSPatrick Sanan @*/
181180ec0b7dSPatrick Sanan PetscErrorCode  PCASMGetSubMatType(PC pc,MatType *sub_mat_type){
181280ec0b7dSPatrick Sanan   PetscErrorCode ierr;
181380ec0b7dSPatrick Sanan 
181480ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));CHKERRQ(ierr);
181580ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
181680ec0b7dSPatrick Sanan }
181780ec0b7dSPatrick Sanan 
181880ec0b7dSPatrick Sanan /*@
181980ec0b7dSPatrick Sanan      PCASMSetSubMatType - Set the type of matrix used for ASM subsolves
182080ec0b7dSPatrick Sanan 
182180ec0b7dSPatrick Sanan    Collective on Mat
182280ec0b7dSPatrick Sanan 
182380ec0b7dSPatrick Sanan    Input Parameters:
182480ec0b7dSPatrick Sanan +  pc             - the PC object
182580ec0b7dSPatrick Sanan -  sub_mat_type   - matrix type
182680ec0b7dSPatrick Sanan 
182780ec0b7dSPatrick Sanan    Options Database Key:
182880ec0b7dSPatrick 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.
182980ec0b7dSPatrick Sanan 
183080ec0b7dSPatrick Sanan    Notes:
183180ec0b7dSPatrick Sanan    See "${PETSC_DIR}/include/petscmat.h" for available types
183280ec0b7dSPatrick Sanan 
183380ec0b7dSPatrick Sanan   Level: advanced
183480ec0b7dSPatrick Sanan 
183580ec0b7dSPatrick Sanan .keywords: PC, PCASM, MatType, set
183680ec0b7dSPatrick Sanan 
183780ec0b7dSPatrick Sanan .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
183880ec0b7dSPatrick Sanan @*/
183980ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type)
184080ec0b7dSPatrick Sanan {
184180ec0b7dSPatrick Sanan   PetscErrorCode ierr;
184280ec0b7dSPatrick Sanan 
184380ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));CHKERRQ(ierr);
184480ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
184580ec0b7dSPatrick Sanan }
1846