xref: /petsc/src/ksp/pc/impls/asm/asm.c (revision 5b3c0d426ef84c76d57947eaebec531d6acb0f6f)
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 */
38*5b3c0d42Seaulisa   VecScatter *lprolongation;      /* mapping from subregion to overlapping process subdomain */
394b9ad928SBarry Smith } PC_ASM;
404b9ad928SBarry Smith 
416849ba73SBarry Smith static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer)
424b9ad928SBarry Smith {
4392bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
446849ba73SBarry Smith   PetscErrorCode ierr;
4513f74950SBarry Smith   PetscMPIInt    rank;
464d219a6aSLisandro Dalcin   PetscInt       i,bsz;
47ace3abfcSBarry Smith   PetscBool      iascii,isstring;
484b9ad928SBarry Smith   PetscViewer    sviewer;
494b9ad928SBarry Smith 
504b9ad928SBarry Smith   PetscFunctionBegin;
51251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
52251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
5332077d6dSBarry Smith   if (iascii) {
543d03bb48SJed Brown     char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set";
558caf3d72SBarry Smith     if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof(overlaps),"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);}
568caf3d72SBarry Smith     if (osm->n > 0) {ierr = PetscSNPrintf(blocks,sizeof(blocks),"total subdomain blocks = %D",osm->n);CHKERRQ(ierr);}
57efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  %s, %s\n",blocks,overlaps);CHKERRQ(ierr);
58efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  restriction/interpolation type - %s\n",PCASMTypes[osm->type]);CHKERRQ(ierr);
59e64f0791SPatrick Sanan     if (osm->dm_subdomains) {ierr = PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: using DM to define subdomains\n");CHKERRQ(ierr);}
6012cd4985SMatthew G. Knepley     if (osm->loctype != PC_COMPOSITE_ADDITIVE) {ierr = PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: local solve composition type - %s\n",PCCompositeTypes[osm->loctype]);CHKERRQ(ierr);}
61ce94432eSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
6292bb6962SLisandro Dalcin     if (osm->same_local_solves) {
634d219a6aSLisandro Dalcin       if (osm->ksp) {
644b9ad928SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Local solve is same for all blocks, in the following KSP and PC objects:\n");CHKERRQ(ierr);
65c5e4d11fSDmitry Karpeev         ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
664d219a6aSLisandro Dalcin         if (!rank) {
674b9ad928SBarry Smith           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
6892bb6962SLisandro Dalcin           ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);
694b9ad928SBarry Smith           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
704b9ad928SBarry Smith         }
71c5e4d11fSDmitry Karpeev         ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
724d219a6aSLisandro Dalcin       }
734b9ad928SBarry Smith     } else {
74c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
754d219a6aSLisandro Dalcin       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,osm->n_local_true);CHKERRQ(ierr);
764d219a6aSLisandro Dalcin       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
774b9ad928SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr);
784b9ad928SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
794d219a6aSLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
80c5e4d11fSDmitry Karpeev       ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
81dd2fa690SBarry Smith       for (i=0; i<osm->n_local_true; i++) {
824d219a6aSLisandro Dalcin         ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr);
834d219a6aSLisandro Dalcin         ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr);
8492bb6962SLisandro Dalcin         ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr);
854d219a6aSLisandro Dalcin         ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
864b9ad928SBarry Smith       }
87c5e4d11fSDmitry Karpeev       ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
884b9ad928SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
894b9ad928SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
90c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
914b9ad928SBarry Smith     }
924b9ad928SBarry Smith   } else if (isstring) {
934d219a6aSLisandro Dalcin     ierr = PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCASMTypes[osm->type]);CHKERRQ(ierr);
94c5e4d11fSDmitry Karpeev     ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
9592bb6962SLisandro Dalcin     if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);}
96c5e4d11fSDmitry Karpeev     ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
974b9ad928SBarry Smith   }
984b9ad928SBarry Smith   PetscFunctionReturn(0);
994b9ad928SBarry Smith }
1004b9ad928SBarry Smith 
10192bb6962SLisandro Dalcin static PetscErrorCode PCASMPrintSubdomains(PC pc)
10292bb6962SLisandro Dalcin {
10392bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
10492bb6962SLisandro Dalcin   const char     *prefix;
10592bb6962SLisandro Dalcin   char           fname[PETSC_MAX_PATH_LEN+1];
106643535aeSDmitry Karpeev   PetscViewer    viewer, sviewer;
107643535aeSDmitry Karpeev   char           *s;
10892bb6962SLisandro Dalcin   PetscInt       i,j,nidx;
10992bb6962SLisandro Dalcin   const PetscInt *idx;
110643535aeSDmitry Karpeev   PetscMPIInt    rank, size;
11192bb6962SLisandro Dalcin   PetscErrorCode ierr;
112846783a0SBarry Smith 
11392bb6962SLisandro Dalcin   PetscFunctionBegin;
114ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr);
115ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr);
11692bb6962SLisandro Dalcin   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
117c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
11892bb6962SLisandro Dalcin   if (fname[0] == 0) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
119ce94432eSBarry Smith   ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr);
120643535aeSDmitry Karpeev   for (i=0; i<osm->n_local; i++) {
121643535aeSDmitry Karpeev     if (i < osm->n_local_true) {
12292bb6962SLisandro Dalcin       ierr = ISGetLocalSize(osm->is[i],&nidx);CHKERRQ(ierr);
12392bb6962SLisandro Dalcin       ierr = ISGetIndices(osm->is[i],&idx);CHKERRQ(ierr);
124643535aeSDmitry Karpeev       /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
125854ce69bSBarry Smith       ierr = PetscMalloc1(16*(nidx+1)+512, &s);CHKERRQ(ierr);
126643535aeSDmitry Karpeev       ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, 16*(nidx+1)+512, &sviewer);CHKERRQ(ierr);
127643535aeSDmitry Karpeev       ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D with overlap:\n", rank, size, i);CHKERRQ(ierr);
12892bb6962SLisandro Dalcin       for (j=0; j<nidx; j++) {
129643535aeSDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr);
13092bb6962SLisandro Dalcin       }
13192bb6962SLisandro Dalcin       ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr);
132643535aeSDmitry Karpeev       ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr);
133643535aeSDmitry Karpeev       ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
134c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
135643535aeSDmitry Karpeev       ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr);
136643535aeSDmitry Karpeev       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
137c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
138643535aeSDmitry Karpeev       ierr = PetscFree(s);CHKERRQ(ierr);
1392b691e39SMatthew Knepley       if (osm->is_local) {
140643535aeSDmitry Karpeev         /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
141854ce69bSBarry Smith         ierr = PetscMalloc1(16*(nidx+1)+512, &s);CHKERRQ(ierr);
142643535aeSDmitry Karpeev         ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, 16*(nidx+1)+512, &sviewer);CHKERRQ(ierr);
14309d011a0SDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D without overlap:\n", rank, size, i);CHKERRQ(ierr);
1442b691e39SMatthew Knepley         ierr = ISGetLocalSize(osm->is_local[i],&nidx);CHKERRQ(ierr);
1452b691e39SMatthew Knepley         ierr = ISGetIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
1462b691e39SMatthew Knepley         for (j=0; j<nidx; j++) {
14709d011a0SDmitry Karpeev           ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr);
1482b691e39SMatthew Knepley         }
1492b691e39SMatthew Knepley         ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
15009d011a0SDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr);
151643535aeSDmitry Karpeev         ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
152c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
153643535aeSDmitry Karpeev         ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr);
154643535aeSDmitry Karpeev         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
155c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
156643535aeSDmitry Karpeev         ierr = PetscFree(s);CHKERRQ(ierr);
157643535aeSDmitry Karpeev       }
1582fa5cd67SKarl Rupp     } else {
159643535aeSDmitry Karpeev       /* Participate in collective viewer calls. */
160c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
161643535aeSDmitry Karpeev       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
162c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
163643535aeSDmitry Karpeev       /* Assume either all ranks have is_local or none do. */
164643535aeSDmitry Karpeev       if (osm->is_local) {
165c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
166643535aeSDmitry Karpeev         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
167c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
168643535aeSDmitry Karpeev       }
169fdc48646SDmitry Karpeev     }
17092bb6962SLisandro Dalcin   }
17192bb6962SLisandro Dalcin   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
172fcfd50ebSBarry Smith   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
17392bb6962SLisandro Dalcin   PetscFunctionReturn(0);
17492bb6962SLisandro Dalcin }
17592bb6962SLisandro Dalcin 
1766849ba73SBarry Smith static PetscErrorCode PCSetUp_ASM(PC pc)
1774b9ad928SBarry Smith {
1784b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
1796849ba73SBarry Smith   PetscErrorCode ierr;
180ace3abfcSBarry Smith   PetscBool      symset,flg;
18187e86712SBarry Smith   PetscInt       i,m,m_local;
1824b9ad928SBarry Smith   MatReuse       scall = MAT_REUSE_MATRIX;
1834b9ad928SBarry Smith   IS             isl;
1844b9ad928SBarry Smith   KSP            ksp;
1854b9ad928SBarry Smith   PC             subpc;
1862dcb1b2aSMatthew Knepley   const char     *prefix,*pprefix;
18723ce1328SBarry Smith   Vec            vec;
1880298fd71SBarry Smith   DM             *domain_dm = NULL;
1894b9ad928SBarry Smith 
1904b9ad928SBarry Smith   PetscFunctionBegin;
1914b9ad928SBarry Smith   if (!pc->setupcalled) {
19292bb6962SLisandro Dalcin 
19392bb6962SLisandro Dalcin     if (!osm->type_set) {
19492bb6962SLisandro Dalcin       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
1952fa5cd67SKarl Rupp       if (symset && flg) osm->type = PC_ASM_BASIC;
19692bb6962SLisandro Dalcin     }
19792bb6962SLisandro Dalcin 
1982b837212SDmitry Karpeev     /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */
1992b837212SDmitry Karpeev     if (osm->n_local_true == PETSC_DECIDE) {
20069ca1f37SDmitry Karpeev       /* no subdomains given */
20165db9045SDmitry Karpeev       /* try pc->dm first, if allowed */
202d709ea83SDmitry Karpeev       if (osm->dm_subdomains && pc->dm) {
203feb221f8SDmitry Karpeev         PetscInt  num_domains, d;
204feb221f8SDmitry Karpeev         char      **domain_names;
2058d4ac253SDmitry Karpeev         IS        *inner_domain_is, *outer_domain_is;
2068d4ac253SDmitry Karpeev         ierr = DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm);CHKERRQ(ierr);
207704f0395SPatrick Sanan         osm->overlap = -1; /* We do not want to increase the overlap of the IS.
208704f0395SPatrick Sanan                               A future improvement of this code might allow one to use
209704f0395SPatrick Sanan                               DM-defined subdomains and also increase the overlap,
210704f0395SPatrick Sanan                               but that is not currently supported */
21169ca1f37SDmitry Karpeev         if (num_domains) {
2128d4ac253SDmitry Karpeev           ierr = PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is);CHKERRQ(ierr);
21369ca1f37SDmitry Karpeev         }
214feb221f8SDmitry Karpeev         for (d = 0; d < num_domains; ++d) {
215a35b7b57SDmitry Karpeev           if (domain_names)    {ierr = PetscFree(domain_names[d]);CHKERRQ(ierr);}
216a35b7b57SDmitry Karpeev           if (inner_domain_is) {ierr = ISDestroy(&inner_domain_is[d]);CHKERRQ(ierr);}
217a35b7b57SDmitry Karpeev           if (outer_domain_is) {ierr = ISDestroy(&outer_domain_is[d]);CHKERRQ(ierr);}
218feb221f8SDmitry Karpeev         }
219feb221f8SDmitry Karpeev         ierr = PetscFree(domain_names);CHKERRQ(ierr);
2208d4ac253SDmitry Karpeev         ierr = PetscFree(inner_domain_is);CHKERRQ(ierr);
2218d4ac253SDmitry Karpeev         ierr = PetscFree(outer_domain_is);CHKERRQ(ierr);
222feb221f8SDmitry Karpeev       }
2232b837212SDmitry Karpeev       if (osm->n_local_true == PETSC_DECIDE) {
22469ca1f37SDmitry Karpeev         /* still no subdomains; use one subdomain per processor */
2252b837212SDmitry Karpeev         osm->n_local_true = 1;
226feb221f8SDmitry Karpeev       }
2272b837212SDmitry Karpeev     }
2282b837212SDmitry Karpeev     { /* determine the global and max number of subdomains */
2296ac3741eSJed Brown       struct {PetscInt max,sum;} inwork,outwork;
230c10200c1SHong Zhang       PetscMPIInt size;
231c10200c1SHong Zhang 
2326ac3741eSJed Brown       inwork.max   = osm->n_local_true;
2336ac3741eSJed Brown       inwork.sum   = osm->n_local_true;
234367daffbSBarry Smith       ierr         = MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
2356ac3741eSJed Brown       osm->n_local = outwork.max;
2366ac3741eSJed Brown       osm->n       = outwork.sum;
237c10200c1SHong Zhang 
238c10200c1SHong Zhang       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
239c10200c1SHong Zhang       if (outwork.max == 1 && outwork.sum == size) {
2407dae84e0SHong Zhang         /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */
241c10200c1SHong Zhang         ierr = MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE);CHKERRQ(ierr);
242c10200c1SHong Zhang       }
2434b9ad928SBarry Smith     }
24478904715SLisandro Dalcin     if (!osm->is) { /* create the index sets */
24578904715SLisandro Dalcin       ierr = PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);CHKERRQ(ierr);
2464b9ad928SBarry Smith     }
247f5234e35SJed Brown     if (osm->n_local_true > 1 && !osm->is_local) {
248785e854fSJed Brown       ierr = PetscMalloc1(osm->n_local_true,&osm->is_local);CHKERRQ(ierr);
249f5234e35SJed Brown       for (i=0; i<osm->n_local_true; i++) {
250f5234e35SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
251f5234e35SJed Brown           ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr);
252f5234e35SJed Brown           ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr);
253f5234e35SJed Brown         } else {
254f5234e35SJed Brown           ierr             = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr);
255f5234e35SJed Brown           osm->is_local[i] = osm->is[i];
256f5234e35SJed Brown         }
257f5234e35SJed Brown       }
258f5234e35SJed Brown     }
25992bb6962SLisandro Dalcin     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
26090d69ab7SBarry Smith     flg  = PETSC_FALSE;
261c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,prefix,"-pc_asm_print_subdomains",&flg,NULL);CHKERRQ(ierr);
26292bb6962SLisandro Dalcin     if (flg) { ierr = PCASMPrintSubdomains(pc);CHKERRQ(ierr); }
2634b9ad928SBarry Smith 
2643d03bb48SJed Brown     if (osm->overlap > 0) {
2654b9ad928SBarry Smith       /* Extend the "overlapping" regions by a number of steps */
26692bb6962SLisandro Dalcin       ierr = MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);CHKERRQ(ierr);
2673d03bb48SJed Brown     }
2686ed231c7SMatthew Knepley     if (osm->sort_indices) {
26992bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2704b9ad928SBarry Smith         ierr = ISSort(osm->is[i]);CHKERRQ(ierr);
2712b691e39SMatthew Knepley         if (osm->is_local) {
2722b691e39SMatthew Knepley           ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr);
2732b691e39SMatthew Knepley         }
2744b9ad928SBarry Smith       }
2756ed231c7SMatthew Knepley     }
2764b9ad928SBarry Smith 
277f6991133SBarry Smith     if (!osm->ksp) {
27878904715SLisandro Dalcin       /* Create the local solvers */
279785e854fSJed Brown       ierr = PetscMalloc1(osm->n_local_true,&osm->ksp);CHKERRQ(ierr);
280feb221f8SDmitry Karpeev       if (domain_dm) {
281feb221f8SDmitry Karpeev         ierr = PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n");CHKERRQ(ierr);
282feb221f8SDmitry Karpeev       }
28392bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2844b9ad928SBarry Smith         ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
285422a814eSBarry Smith         ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr);
2863bb1ff40SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr);
28792bb6962SLisandro Dalcin         ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
2884b9ad928SBarry Smith         ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
2894b9ad928SBarry Smith         ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
2904b9ad928SBarry Smith         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
2914b9ad928SBarry Smith         ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
2924b9ad928SBarry Smith         ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
293feb221f8SDmitry Karpeev         if (domain_dm) {
294feb221f8SDmitry Karpeev           ierr = KSPSetDM(ksp, domain_dm[i]);CHKERRQ(ierr);
295feb221f8SDmitry Karpeev           ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr);
296feb221f8SDmitry Karpeev           ierr = DMDestroy(&domain_dm[i]);CHKERRQ(ierr);
297feb221f8SDmitry Karpeev         }
2984b9ad928SBarry Smith         osm->ksp[i] = ksp;
2994b9ad928SBarry Smith       }
300feb221f8SDmitry Karpeev       if (domain_dm) {
301feb221f8SDmitry Karpeev         ierr = PetscFree(domain_dm);CHKERRQ(ierr);
302feb221f8SDmitry Karpeev       }
303f6991133SBarry Smith     }
304fb745f2cSMatthew G. Knepley     if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
305fb745f2cSMatthew G. Knepley       PetscInt m;
306fb745f2cSMatthew G. Knepley       ierr = ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &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);
349*5b3c0d42Seaulisa     if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE)
350*5b3c0d42Seaulisa       ierr = PetscMalloc1(osm->n_local,&osm->lprolongation);CHKERRQ(ierr);
351*5b3c0d42Seaulisa 
35280ec0b7dSPatrick Sanan     ierr = PetscMalloc1(osm->n_local,&osm->x);CHKERRQ(ierr);
35380ec0b7dSPatrick Sanan     ierr = PetscMalloc1(osm->n_local,&osm->y);CHKERRQ(ierr);
35480ec0b7dSPatrick Sanan     ierr = PetscMalloc1(osm->n_local,&osm->y_local);CHKERRQ(ierr);
35580ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; ++i) {
35680ec0b7dSPatrick Sanan       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
35780ec0b7dSPatrick Sanan       ierr = MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);CHKERRQ(ierr);
35880ec0b7dSPatrick Sanan       ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
35980ec0b7dSPatrick Sanan       ierr = VecScatterCreate(vec,osm->is[i],osm->x[i],isl,&osm->restriction[i]);CHKERRQ(ierr);
36080ec0b7dSPatrick Sanan       ierr = ISDestroy(&isl);CHKERRQ(ierr);
36180ec0b7dSPatrick Sanan       ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr);
362f1ee410cSBarry Smith       if (osm->localization) {
36380ec0b7dSPatrick Sanan         ISLocalToGlobalMapping ltog;
36480ec0b7dSPatrick Sanan         IS                     isll;
36580ec0b7dSPatrick Sanan         const PetscInt         *idx_local;
36680ec0b7dSPatrick Sanan         PetscInt               *idx,nout;
36780ec0b7dSPatrick Sanan 
36880ec0b7dSPatrick Sanan         ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],&ltog);CHKERRQ(ierr);
36980ec0b7dSPatrick Sanan         ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr);
37080ec0b7dSPatrick Sanan         ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
37180ec0b7dSPatrick Sanan         ierr = PetscMalloc1(m_local,&idx);CHKERRQ(ierr);
37280ec0b7dSPatrick Sanan         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx);CHKERRQ(ierr);
37380ec0b7dSPatrick Sanan         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
37480ec0b7dSPatrick Sanan         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is");
37580ec0b7dSPatrick Sanan         ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
37680ec0b7dSPatrick Sanan         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
37780ec0b7dSPatrick Sanan         ierr = ISCreateStride(PETSC_COMM_SELF,m_local,0,1,&isl);CHKERRQ(ierr);
37880ec0b7dSPatrick Sanan         ierr = VecCreateSeq(PETSC_COMM_SELF,m_local,&osm->y_local[i]);CHKERRQ(ierr);
37980ec0b7dSPatrick Sanan         ierr = VecScatterCreate(osm->y[i],isll,osm->y_local[i],isl,&osm->localization[i]);CHKERRQ(ierr);
38080ec0b7dSPatrick Sanan         ierr = ISDestroy(&isll);CHKERRQ(ierr);
38180ec0b7dSPatrick Sanan 
38280ec0b7dSPatrick Sanan         ierr = VecScatterCreate(vec,osm->is_local[i],osm->y_local[i],isl,&osm->prolongation[i]);CHKERRQ(ierr);
38380ec0b7dSPatrick Sanan         ierr = ISDestroy(&isl);CHKERRQ(ierr);
38480ec0b7dSPatrick Sanan       } else {
38580ec0b7dSPatrick Sanan         osm->y_local[i] = osm->y[i];
38680ec0b7dSPatrick Sanan         ierr = PetscObjectReference((PetscObject) osm->y[i]);CHKERRQ(ierr);
38780ec0b7dSPatrick Sanan         osm->prolongation[i] = osm->restriction[i];
38880ec0b7dSPatrick Sanan         ierr = PetscObjectReference((PetscObject) osm->restriction[i]);CHKERRQ(ierr);
38980ec0b7dSPatrick Sanan       }
390*5b3c0d42Seaulisa 
391*5b3c0d42Seaulisa       if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE && i < osm->n_local_true - 1) {
392*5b3c0d42Seaulisa         ISLocalToGlobalMapping ltog;
393*5b3c0d42Seaulisa         IS                     isll;
394*5b3c0d42Seaulisa         const PetscInt         *idx_is;
395*5b3c0d42Seaulisa         PetscInt               *idx_lis,nout;
396*5b3c0d42Seaulisa 
397*5b3c0d42Seaulisa         ierr = ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);CHKERRQ(ierr);
398*5b3c0d42Seaulisa         ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
399*5b3c0d42Seaulisa         ierr = ISGetIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
400*5b3c0d42Seaulisa         ierr = PetscMalloc1(m,&idx_lis);CHKERRQ(ierr);
401*5b3c0d42Seaulisa         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis);CHKERRQ(ierr);
402*5b3c0d42Seaulisa         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
403*5b3c0d42Seaulisa         if (nout != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis");
404*5b3c0d42Seaulisa         ierr = ISRestoreIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
405*5b3c0d42Seaulisa         ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
406*5b3c0d42Seaulisa         ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
407*5b3c0d42Seaulisa         ierr = VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lprolongation[i]);CHKERRQ(ierr);
408*5b3c0d42Seaulisa         ierr = ISDestroy(&isll);CHKERRQ(ierr);
409*5b3c0d42Seaulisa         ierr = ISDestroy(&isl);CHKERRQ(ierr);
410*5b3c0d42Seaulisa       }
411*5b3c0d42Seaulisa 
41280ec0b7dSPatrick Sanan     }
41380ec0b7dSPatrick Sanan     for (i=osm->n_local_true; i<osm->n_local; i++) {
41480ec0b7dSPatrick Sanan       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&osm->x[i]);CHKERRQ(ierr);
41580ec0b7dSPatrick Sanan       ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr);
41680ec0b7dSPatrick Sanan       ierr = VecDuplicate(osm->x[i],&osm->y_local[i]);CHKERRQ(ierr);
41780ec0b7dSPatrick Sanan       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&isl);CHKERRQ(ierr);
41880ec0b7dSPatrick Sanan       ierr = VecScatterCreate(vec,isl,osm->x[i],isl,&osm->restriction[i]);CHKERRQ(ierr);
419f1ee410cSBarry Smith       if (osm->localization) {
42080ec0b7dSPatrick Sanan         ierr = VecScatterCreate(osm->y[i],isl,osm->y_local[i],isl,&osm->localization[i]);CHKERRQ(ierr);
42180ec0b7dSPatrick Sanan         ierr = VecScatterCreate(vec,isl,osm->x[i],isl,&osm->prolongation[i]);CHKERRQ(ierr);
42280ec0b7dSPatrick Sanan       } else {
42380ec0b7dSPatrick Sanan         osm->prolongation[i] = osm->restriction[i];
42480ec0b7dSPatrick Sanan         ierr                 = PetscObjectReference((PetscObject) osm->restriction[i]);CHKERRQ(ierr);
42580ec0b7dSPatrick Sanan       }
42680ec0b7dSPatrick Sanan       ierr = ISDestroy(&isl);CHKERRQ(ierr);
42780ec0b7dSPatrick Sanan     }
42880ec0b7dSPatrick Sanan     ierr = VecDestroy(&vec);CHKERRQ(ierr);
42980ec0b7dSPatrick Sanan   }
43080ec0b7dSPatrick Sanan 
431fb745f2cSMatthew G. Knepley   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
432235cc4ceSMatthew G. Knepley     IS      *cis;
433235cc4ceSMatthew G. Knepley     PetscInt c;
434235cc4ceSMatthew G. Knepley 
435235cc4ceSMatthew G. Knepley     ierr = PetscMalloc1(osm->n_local_true, &cis);CHKERRQ(ierr);
436235cc4ceSMatthew G. Knepley     for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis;
4377dae84e0SHong Zhang     ierr = MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);CHKERRQ(ierr);
438235cc4ceSMatthew G. Knepley     ierr = PetscFree(cis);CHKERRQ(ierr);
439fb745f2cSMatthew G. Knepley   }
4404b9ad928SBarry Smith 
4414b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
4424b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
44392bb6962SLisandro Dalcin   ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
4444b9ad928SBarry Smith 
44592bb6962SLisandro Dalcin   /*
44692bb6962SLisandro Dalcin      Loop over subdomains putting them into local ksp
44792bb6962SLisandro Dalcin   */
44892bb6962SLisandro Dalcin   for (i=0; i<osm->n_local_true; i++) {
44923ee1639SBarry Smith     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr);
45092bb6962SLisandro Dalcin     if (!pc->setupcalled) {
451bf108f30SBarry Smith       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
4524b9ad928SBarry Smith     }
45392bb6962SLisandro Dalcin   }
4544b9ad928SBarry Smith   PetscFunctionReturn(0);
4554b9ad928SBarry Smith }
4564b9ad928SBarry Smith 
4576849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc)
4584b9ad928SBarry Smith {
4594b9ad928SBarry Smith   PC_ASM             *osm = (PC_ASM*)pc->data;
4606849ba73SBarry Smith   PetscErrorCode     ierr;
46113f74950SBarry Smith   PetscInt           i;
462539c167fSBarry Smith   KSPConvergedReason reason;
4634b9ad928SBarry Smith 
4644b9ad928SBarry Smith   PetscFunctionBegin;
4654b9ad928SBarry Smith   for (i=0; i<osm->n_local_true; i++) {
4664b9ad928SBarry Smith     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
467539c167fSBarry Smith     ierr = KSPGetConvergedReason(osm->ksp[i],&reason);CHKERRQ(ierr);
468539c167fSBarry Smith     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
469261222b2SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
470e0eafd54SHong Zhang     }
4714b9ad928SBarry Smith   }
4724b9ad928SBarry Smith   PetscFunctionReturn(0);
4734b9ad928SBarry Smith }
4744b9ad928SBarry Smith 
4756849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y)
4764b9ad928SBarry Smith {
4774b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
4786849ba73SBarry Smith   PetscErrorCode ierr;
47913f74950SBarry Smith   PetscInt       i,n_local = osm->n_local,n_local_true = osm->n_local_true;
4804b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
4814b9ad928SBarry Smith 
4824b9ad928SBarry Smith   PetscFunctionBegin;
4834b9ad928SBarry Smith   /*
4844b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
4854b9ad928SBarry Smith      subdomain values (leaving the other values 0).
4864b9ad928SBarry Smith   */
4874b9ad928SBarry Smith   if (!(osm->type & PC_ASM_RESTRICT)) {
4884b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
4894b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
490f21eeb63SBarry Smith     for (i=0; i<n_local_true; i++) {
4913d03bb48SJed Brown       ierr = VecZeroEntries(osm->x[i]);CHKERRQ(ierr);
4924b9ad928SBarry Smith     }
4934b9ad928SBarry Smith   }
4942fa5cd67SKarl Rupp   if (!(osm->type & PC_ASM_INTERPOLATE)) reverse = SCATTER_REVERSE_LOCAL;
4954b9ad928SBarry Smith 
49612cd4985SMatthew G. Knepley   switch (osm->loctype)
49712cd4985SMatthew G. Knepley   {
49812cd4985SMatthew G. Knepley   case PC_COMPOSITE_ADDITIVE:
4994b9ad928SBarry Smith     for (i=0; i<n_local; i++) {
5002b691e39SMatthew Knepley       ierr = VecScatterBegin(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
5014b9ad928SBarry Smith     }
50210bd88b9SJed Brown     ierr = VecZeroEntries(y);CHKERRQ(ierr);
5034b9ad928SBarry Smith     /* do the local solves */
5044b9ad928SBarry Smith     for (i=0; i<n_local_true; i++) {
5052b691e39SMatthew Knepley       ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
50623ce1328SBarry Smith       ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
507b9845e0eSMatthew Knepley       if (osm->localization) {
50853888de8SMatthew Knepley         ierr = VecScatterBegin(osm->localization[i],osm->y[i],osm->y_local[i],INSERT_VALUES,forward);CHKERRQ(ierr);
50953888de8SMatthew Knepley         ierr = VecScatterEnd(osm->localization[i],osm->y[i],osm->y_local[i],INSERT_VALUES,forward);CHKERRQ(ierr);
510b9845e0eSMatthew Knepley       }
51153888de8SMatthew Knepley       ierr = VecScatterBegin(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
5124b9ad928SBarry Smith     }
5134b9ad928SBarry Smith     /* handle the rest of the scatters that do not have local solves */
5144b9ad928SBarry Smith     for (i=n_local_true; i<n_local; i++) {
5152b691e39SMatthew Knepley       ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
51653888de8SMatthew Knepley       ierr = VecScatterBegin(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
5174b9ad928SBarry Smith     }
5184b9ad928SBarry Smith     for (i=0; i<n_local; i++) {
51953888de8SMatthew Knepley       ierr = VecScatterEnd(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
5204b9ad928SBarry Smith     }
52112cd4985SMatthew G. Knepley     break;
52212cd4985SMatthew G. Knepley   case PC_COMPOSITE_MULTIPLICATIVE:
52312cd4985SMatthew G. Knepley     ierr = VecZeroEntries(y);CHKERRQ(ierr);
524*5b3c0d42Seaulisa     ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
52512cd4985SMatthew G. Knepley     /* do the local solves */
52612cd4985SMatthew G. Knepley     for (i = 0; i < n_local_true; ++i) {
527*5b3c0d42Seaulisa       if (i == 0) {
52812cd4985SMatthew G. Knepley         ierr = VecZeroEntries(osm->x[i]);CHKERRQ(ierr);
52912cd4985SMatthew G. Knepley       }
53012cd4985SMatthew G. Knepley       ierr = VecScatterBegin(osm->restriction[i], x, osm->x[i], ADD_VALUES, forward);CHKERRQ(ierr);
53112cd4985SMatthew G. Knepley       ierr = VecScatterEnd(osm->restriction[i], x, osm->x[i], ADD_VALUES, forward);CHKERRQ(ierr);
53212cd4985SMatthew G. Knepley       ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
53312cd4985SMatthew G. Knepley       if (osm->localization) {
53412cd4985SMatthew G. Knepley         ierr = VecScatterBegin(osm->localization[i], osm->y[i], osm->y_local[i], INSERT_VALUES, forward);CHKERRQ(ierr);
53512cd4985SMatthew G. Knepley         ierr = VecScatterEnd(osm->localization[i], osm->y[i], osm->y_local[i], INSERT_VALUES, forward);CHKERRQ(ierr);
53612cd4985SMatthew G. Knepley       }
53712cd4985SMatthew G. Knepley       ierr = VecScatterBegin(osm->prolongation[i], osm->y_local[i], y, ADD_VALUES, reverse);CHKERRQ(ierr);
53812cd4985SMatthew G. Knepley       ierr = VecScatterEnd(osm->prolongation[i], osm->y_local[i], y, ADD_VALUES, reverse);CHKERRQ(ierr);
5397c3d802fSMatthew G. Knepley       if (i < n_local_true-1) {
540*5b3c0d42Seaulisa 	ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
541*5b3c0d42Seaulisa         ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
542*5b3c0d42Seaulisa 	ierr = VecCopy(osm->ly, osm->lx);CHKERRQ(ierr);
543*5b3c0d42Seaulisa         ierr = VecScale(osm->lx, -1.0);CHKERRQ(ierr);
544*5b3c0d42Seaulisa         ierr = MatMult(osm->lmats[i+1], osm->lx, osm->x[i+1]);CHKERRQ(ierr);
5457c3d802fSMatthew G. Knepley       }
54612cd4985SMatthew G. Knepley     }
54712cd4985SMatthew G. Knepley     /* handle the rest of the scatters that do not have local solves */
54812cd4985SMatthew G. Knepley     for (i = n_local_true; i < n_local; ++i) {
54912cd4985SMatthew G. Knepley       ierr = VecScatterBegin(osm->restriction[i], x, osm->x[i], INSERT_VALUES, forward);CHKERRQ(ierr);
55012cd4985SMatthew G. Knepley       ierr = VecScatterEnd(osm->restriction[i], x, osm->x[i], INSERT_VALUES, forward);CHKERRQ(ierr);
55112cd4985SMatthew G. Knepley       ierr = VecScatterBegin(osm->prolongation[i], osm->y_local[i], y, ADD_VALUES, reverse);CHKERRQ(ierr);
55212cd4985SMatthew G. Knepley       ierr = VecScatterEnd(osm->prolongation[i], osm->y_local[i], y, ADD_VALUES, reverse);CHKERRQ(ierr);
55312cd4985SMatthew G. Knepley     }
55412cd4985SMatthew G. Knepley     break;
55512cd4985SMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
55612cd4985SMatthew G. Knepley   }
5574b9ad928SBarry Smith   PetscFunctionReturn(0);
5584b9ad928SBarry Smith }
5594b9ad928SBarry Smith 
5606849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y)
5614b9ad928SBarry Smith {
5624b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
5636849ba73SBarry Smith   PetscErrorCode ierr;
56413f74950SBarry Smith   PetscInt       i,n_local = osm->n_local,n_local_true = osm->n_local_true;
5654b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
5664b9ad928SBarry Smith 
5674b9ad928SBarry Smith   PetscFunctionBegin;
5684b9ad928SBarry Smith   /*
5694b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
5704b9ad928SBarry Smith      subdomain values (leaving the other values 0).
5714b9ad928SBarry Smith 
5724b9ad928SBarry Smith      Note: these are reversed from the PCApply_ASM() because we are applying the
5734b9ad928SBarry Smith      transpose of the three terms
5744b9ad928SBarry Smith   */
5754b9ad928SBarry Smith   if (!(osm->type & PC_ASM_INTERPOLATE)) {
5764b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
5774b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
578f21eeb63SBarry Smith     for (i=0; i<n_local_true; i++) {
57910bd88b9SJed Brown       ierr = VecZeroEntries(osm->x[i]);CHKERRQ(ierr);
5804b9ad928SBarry Smith     }
5814b9ad928SBarry Smith   }
5822fa5cd67SKarl Rupp   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;
5834b9ad928SBarry Smith 
5844b9ad928SBarry Smith   for (i=0; i<n_local; i++) {
5852b691e39SMatthew Knepley     ierr = VecScatterBegin(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
5864b9ad928SBarry Smith   }
58710bd88b9SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
5884b9ad928SBarry Smith   /* do the local solves */
5894b9ad928SBarry Smith   for (i=0; i<n_local_true; i++) {
5902b691e39SMatthew Knepley     ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
59123ce1328SBarry Smith     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
592b9845e0eSMatthew Knepley     if (osm->localization) {
59353888de8SMatthew Knepley       ierr = VecScatterBegin(osm->localization[i],osm->y[i],osm->y_local[i],INSERT_VALUES,forward);CHKERRQ(ierr);
59453888de8SMatthew Knepley       ierr = VecScatterEnd(osm->localization[i],osm->y[i],osm->y_local[i],INSERT_VALUES,forward);CHKERRQ(ierr);
595b9845e0eSMatthew Knepley     }
59653888de8SMatthew Knepley     ierr = VecScatterBegin(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
5974b9ad928SBarry Smith   }
5984b9ad928SBarry Smith   /* handle the rest of the scatters that do not have local solves */
5994b9ad928SBarry Smith   for (i=n_local_true; i<n_local; i++) {
6002b691e39SMatthew Knepley     ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
60153888de8SMatthew Knepley     ierr = VecScatterBegin(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
6024b9ad928SBarry Smith   }
6034b9ad928SBarry Smith   for (i=0; i<n_local; i++) {
60453888de8SMatthew Knepley     ierr = VecScatterEnd(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
6054b9ad928SBarry Smith   }
6064b9ad928SBarry Smith   PetscFunctionReturn(0);
6074b9ad928SBarry Smith }
6084b9ad928SBarry Smith 
609e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc)
6104b9ad928SBarry Smith {
6114b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
6126849ba73SBarry Smith   PetscErrorCode ierr;
61313f74950SBarry Smith   PetscInt       i;
6144b9ad928SBarry Smith 
6154b9ad928SBarry Smith   PetscFunctionBegin;
61692bb6962SLisandro Dalcin   if (osm->ksp) {
61792bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
618e91c6855SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
61992bb6962SLisandro Dalcin     }
62092bb6962SLisandro Dalcin   }
621e09e08ccSBarry Smith   if (osm->pmat) {
62292bb6962SLisandro Dalcin     if (osm->n_local_true > 0) {
62330a70a9aSHong Zhang       ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
62492bb6962SLisandro Dalcin     }
62592bb6962SLisandro Dalcin   }
6262b691e39SMatthew Knepley   if (osm->restriction) {
6274b9ad928SBarry Smith     for (i=0; i<osm->n_local; i++) {
628fcfd50ebSBarry Smith       ierr = VecScatterDestroy(&osm->restriction[i]);CHKERRQ(ierr);
629fcfd50ebSBarry Smith       if (osm->localization) {ierr = VecScatterDestroy(&osm->localization[i]);CHKERRQ(ierr);}
630fcfd50ebSBarry Smith       ierr = VecScatterDestroy(&osm->prolongation[i]);CHKERRQ(ierr);
631*5b3c0d42Seaulisa       if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE && i < osm->n_local_true - 1){
632*5b3c0d42Seaulisa 	ierr = VecScatterDestroy(&osm->lprolongation[i]);CHKERRQ(ierr);
633*5b3c0d42Seaulisa       }
634fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
635fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
636fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y_local[i]);CHKERRQ(ierr);
6374b9ad928SBarry Smith     }
6382b691e39SMatthew Knepley     ierr = PetscFree(osm->restriction);CHKERRQ(ierr);
639b9845e0eSMatthew Knepley     if (osm->localization) {ierr = PetscFree(osm->localization);CHKERRQ(ierr);}
6402b691e39SMatthew Knepley     ierr = PetscFree(osm->prolongation);CHKERRQ(ierr);
641*5b3c0d42Seaulisa     if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE){
642*5b3c0d42Seaulisa       ierr = PetscFree(osm->lprolongation);CHKERRQ(ierr);
643*5b3c0d42Seaulisa     }
64405b42c5fSBarry Smith     ierr = PetscFree(osm->x);CHKERRQ(ierr);
64578904715SLisandro Dalcin     ierr = PetscFree(osm->y);CHKERRQ(ierr);
64653888de8SMatthew Knepley     ierr = PetscFree(osm->y_local);CHKERRQ(ierr);
64792bb6962SLisandro Dalcin   }
6482b691e39SMatthew Knepley   ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
649fb745f2cSMatthew G. Knepley   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
650fb745f2cSMatthew G. Knepley     ierr = ISDestroy(&osm->lis);CHKERRQ(ierr);
651235cc4ceSMatthew G. Knepley     ierr = MatDestroyMatrices(osm->n_local_true, &osm->lmats);CHKERRQ(ierr);
652fb745f2cSMatthew G. Knepley     ierr = VecDestroy(&osm->lx);CHKERRQ(ierr);
653fb745f2cSMatthew G. Knepley     ierr = VecDestroy(&osm->ly);CHKERRQ(ierr);
654fb745f2cSMatthew G. Knepley   }
6552fa5cd67SKarl Rupp 
65680ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
65780ec0b7dSPatrick Sanan 
658e91c6855SBarry Smith   osm->is       = 0;
659e91c6855SBarry Smith   osm->is_local = 0;
660e91c6855SBarry Smith   PetscFunctionReturn(0);
661e91c6855SBarry Smith }
662e91c6855SBarry Smith 
663e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc)
664e91c6855SBarry Smith {
665e91c6855SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
666e91c6855SBarry Smith   PetscErrorCode ierr;
667e91c6855SBarry Smith   PetscInt       i;
668e91c6855SBarry Smith 
669e91c6855SBarry Smith   PetscFunctionBegin;
670e91c6855SBarry Smith   ierr = PCReset_ASM(pc);CHKERRQ(ierr);
671e91c6855SBarry Smith   if (osm->ksp) {
672e91c6855SBarry Smith     for (i=0; i<osm->n_local_true; i++) {
673fcfd50ebSBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
674e91c6855SBarry Smith     }
675e91c6855SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
676e91c6855SBarry Smith   }
677e91c6855SBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
6784b9ad928SBarry Smith   PetscFunctionReturn(0);
6794b9ad928SBarry Smith }
6804b9ad928SBarry Smith 
6814416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc)
6824b9ad928SBarry Smith {
6834b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
6846849ba73SBarry Smith   PetscErrorCode ierr;
6859dcbbd2bSBarry Smith   PetscInt       blocks,ovl;
686ace3abfcSBarry Smith   PetscBool      symset,flg;
68792bb6962SLisandro Dalcin   PCASMType      asmtype;
68812cd4985SMatthew G. Knepley   PCCompositeType loctype;
68980ec0b7dSPatrick Sanan   char           sub_mat_type[256];
6904b9ad928SBarry Smith 
6914b9ad928SBarry Smith   PetscFunctionBegin;
692bf108f30SBarry Smith   /* set the type to symmetric if matrix is symmetric */
69392bb6962SLisandro Dalcin   if (!osm->type_set && pc->pmat) {
69492bb6962SLisandro Dalcin     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
6952fa5cd67SKarl Rupp     if (symset && flg) osm->type = PC_ASM_BASIC;
696bf108f30SBarry Smith   }
697e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");CHKERRQ(ierr);
698d709ea83SDmitry Karpeev   ierr = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr);
69990d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
70065db9045SDmitry Karpeev   if (flg) {
701f77a5242SKarl Rupp     ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr);
702d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
70365db9045SDmitry Karpeev   }
70490d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
70565db9045SDmitry Karpeev   if (flg) {
70665db9045SDmitry Karpeev     ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr);
707d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
70865db9045SDmitry Karpeev   }
70990d69ab7SBarry Smith   flg  = PETSC_FALSE;
71090d69ab7SBarry Smith   ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
71192bb6962SLisandro Dalcin   if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); }
71212cd4985SMatthew G. Knepley   flg  = PETSC_FALSE;
71312cd4985SMatthew G. Knepley   ierr = PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr);
71412cd4985SMatthew G. Knepley   if (flg) {ierr = PCASMSetLocalType(pc,loctype);CHKERRQ(ierr); }
71580ec0b7dSPatrick Sanan   ierr = PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);CHKERRQ(ierr);
71680ec0b7dSPatrick Sanan   if(flg){
717459726d8SSatish Balay     ierr = PCASMSetSubMatType(pc,sub_mat_type);CHKERRQ(ierr);
71880ec0b7dSPatrick Sanan   }
7194b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7204b9ad928SBarry Smith   PetscFunctionReturn(0);
7214b9ad928SBarry Smith }
7224b9ad928SBarry Smith 
7234b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/
7244b9ad928SBarry Smith 
7251e6b0712SBarry Smith static PetscErrorCode  PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[])
7264b9ad928SBarry Smith {
7274b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
72892bb6962SLisandro Dalcin   PetscErrorCode ierr;
72992bb6962SLisandro Dalcin   PetscInt       i;
7304b9ad928SBarry Smith 
7314b9ad928SBarry Smith   PetscFunctionBegin;
732e32f2f54SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
733ce94432eSBarry 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().");
734e7e72b3dSBarry Smith 
7354b9ad928SBarry Smith   if (!pc->setupcalled) {
73692bb6962SLisandro Dalcin     if (is) {
73792bb6962SLisandro Dalcin       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
73892bb6962SLisandro Dalcin     }
739832fc9a5SMatthew Knepley     if (is_local) {
740832fc9a5SMatthew Knepley       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
741832fc9a5SMatthew Knepley     }
7422b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
7432fa5cd67SKarl Rupp 
7444b9ad928SBarry Smith     osm->n_local_true = n;
74592bb6962SLisandro Dalcin     osm->is           = 0;
7462b691e39SMatthew Knepley     osm->is_local     = 0;
74792bb6962SLisandro Dalcin     if (is) {
748785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is);CHKERRQ(ierr);
7492fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is[i] = is[i];
7503d03bb48SJed Brown       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
7513d03bb48SJed Brown       osm->overlap = -1;
75292bb6962SLisandro Dalcin     }
7532b691e39SMatthew Knepley     if (is_local) {
754785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is_local);CHKERRQ(ierr);
7552fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is_local[i] = is_local[i];
756a35b7b57SDmitry Karpeev       if (!is) {
757785e854fSJed Brown         ierr = PetscMalloc1(osm->n_local_true,&osm->is);CHKERRQ(ierr);
758a35b7b57SDmitry Karpeev         for (i=0; i<osm->n_local_true; i++) {
759a35b7b57SDmitry Karpeev           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
760a35b7b57SDmitry Karpeev             ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr);
761a35b7b57SDmitry Karpeev             ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr);
762a35b7b57SDmitry Karpeev           } else {
763a35b7b57SDmitry Karpeev             ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr);
764a35b7b57SDmitry Karpeev             osm->is[i] = osm->is_local[i];
765a35b7b57SDmitry Karpeev           }
766a35b7b57SDmitry Karpeev         }
767a35b7b57SDmitry Karpeev       }
7682b691e39SMatthew Knepley     }
7694b9ad928SBarry Smith   }
7704b9ad928SBarry Smith   PetscFunctionReturn(0);
7714b9ad928SBarry Smith }
7724b9ad928SBarry Smith 
7731e6b0712SBarry Smith static PetscErrorCode  PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local)
7744b9ad928SBarry Smith {
7754b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
7766849ba73SBarry Smith   PetscErrorCode ierr;
77713f74950SBarry Smith   PetscMPIInt    rank,size;
77878904715SLisandro Dalcin   PetscInt       n;
7794b9ad928SBarry Smith 
7804b9ad928SBarry Smith   PetscFunctionBegin;
781ce94432eSBarry Smith   if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
782ce94432eSBarry 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.");
7834b9ad928SBarry Smith 
7844b9ad928SBarry Smith   /*
785880770e9SJed Brown      Split the subdomains equally among all processors
7864b9ad928SBarry Smith   */
787ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
788ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
7894b9ad928SBarry Smith   n    = N/size + ((N % size) > rank);
790e32f2f54SBarry 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);
791e32f2f54SBarry Smith   if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp().");
7924b9ad928SBarry Smith   if (!pc->setupcalled) {
7932b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
7942fa5cd67SKarl Rupp 
7954b9ad928SBarry Smith     osm->n_local_true = n;
7964b9ad928SBarry Smith     osm->is           = 0;
7972b691e39SMatthew Knepley     osm->is_local     = 0;
7984b9ad928SBarry Smith   }
7994b9ad928SBarry Smith   PetscFunctionReturn(0);
8004b9ad928SBarry Smith }
8014b9ad928SBarry Smith 
8021e6b0712SBarry Smith static PetscErrorCode  PCASMSetOverlap_ASM(PC pc,PetscInt ovl)
8034b9ad928SBarry Smith {
80492bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
8054b9ad928SBarry Smith 
8064b9ad928SBarry Smith   PetscFunctionBegin;
807ce94432eSBarry Smith   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
808ce94432eSBarry Smith   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp().");
8092fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
8104b9ad928SBarry Smith   PetscFunctionReturn(0);
8114b9ad928SBarry Smith }
8124b9ad928SBarry Smith 
8131e6b0712SBarry Smith static PetscErrorCode  PCASMSetType_ASM(PC pc,PCASMType type)
8144b9ad928SBarry Smith {
81592bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
8164b9ad928SBarry Smith 
8174b9ad928SBarry Smith   PetscFunctionBegin;
8184b9ad928SBarry Smith   osm->type     = type;
819bf108f30SBarry Smith   osm->type_set = PETSC_TRUE;
8204b9ad928SBarry Smith   PetscFunctionReturn(0);
8214b9ad928SBarry Smith }
8224b9ad928SBarry Smith 
823c60c7ad4SBarry Smith static PetscErrorCode  PCASMGetType_ASM(PC pc,PCASMType *type)
824c60c7ad4SBarry Smith {
825c60c7ad4SBarry Smith   PC_ASM *osm = (PC_ASM*)pc->data;
826c60c7ad4SBarry Smith 
827c60c7ad4SBarry Smith   PetscFunctionBegin;
828c60c7ad4SBarry Smith   *type = osm->type;
829c60c7ad4SBarry Smith   PetscFunctionReturn(0);
830c60c7ad4SBarry Smith }
831c60c7ad4SBarry Smith 
83212cd4985SMatthew G. Knepley static PetscErrorCode  PCASMSetLocalType_ASM(PC pc, PCCompositeType type)
83312cd4985SMatthew G. Knepley {
83412cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
83512cd4985SMatthew G. Knepley 
83612cd4985SMatthew G. Knepley   PetscFunctionBegin;
83712cd4985SMatthew G. Knepley   osm->loctype = type;
83812cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
83912cd4985SMatthew G. Knepley }
84012cd4985SMatthew G. Knepley 
84112cd4985SMatthew G. Knepley static PetscErrorCode  PCASMGetLocalType_ASM(PC pc, PCCompositeType *type)
84212cd4985SMatthew G. Knepley {
84312cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
84412cd4985SMatthew G. Knepley 
84512cd4985SMatthew G. Knepley   PetscFunctionBegin;
84612cd4985SMatthew G. Knepley   *type = osm->loctype;
84712cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
84812cd4985SMatthew G. Knepley }
84912cd4985SMatthew G. Knepley 
8501e6b0712SBarry Smith static PetscErrorCode  PCASMSetSortIndices_ASM(PC pc,PetscBool  doSort)
8516ed231c7SMatthew Knepley {
8526ed231c7SMatthew Knepley   PC_ASM *osm = (PC_ASM*)pc->data;
8536ed231c7SMatthew Knepley 
8546ed231c7SMatthew Knepley   PetscFunctionBegin;
8556ed231c7SMatthew Knepley   osm->sort_indices = doSort;
8566ed231c7SMatthew Knepley   PetscFunctionReturn(0);
8576ed231c7SMatthew Knepley }
8586ed231c7SMatthew Knepley 
8591e6b0712SBarry Smith static PetscErrorCode  PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
8604b9ad928SBarry Smith {
86192bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
862dfbe8321SBarry Smith   PetscErrorCode ierr;
8634b9ad928SBarry Smith 
8644b9ad928SBarry Smith   PetscFunctionBegin;
865ce94432eSBarry 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");
8664b9ad928SBarry Smith 
8672fa5cd67SKarl Rupp   if (n_local) *n_local = osm->n_local_true;
86892bb6962SLisandro Dalcin   if (first_local) {
869ce94432eSBarry Smith     ierr          = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
87092bb6962SLisandro Dalcin     *first_local -= osm->n_local_true;
87192bb6962SLisandro Dalcin   }
87292bb6962SLisandro Dalcin   if (ksp) {
87392bb6962SLisandro Dalcin     /* Assume that local solves are now different; not necessarily
87492bb6962SLisandro Dalcin        true though!  This flag is used only for PCView_ASM() */
87592bb6962SLisandro Dalcin     *ksp                   = osm->ksp;
87692bb6962SLisandro Dalcin     osm->same_local_solves = PETSC_FALSE;
87792bb6962SLisandro Dalcin   }
8784b9ad928SBarry Smith   PetscFunctionReturn(0);
8794b9ad928SBarry Smith }
8804b9ad928SBarry Smith 
88180ec0b7dSPatrick Sanan static PetscErrorCode  PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type)
88280ec0b7dSPatrick Sanan {
88380ec0b7dSPatrick Sanan   PC_ASM         *osm = (PC_ASM*)pc->data;
88480ec0b7dSPatrick Sanan 
88580ec0b7dSPatrick Sanan   PetscFunctionBegin;
88680ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
88780ec0b7dSPatrick Sanan   PetscValidPointer(sub_mat_type,2);
88880ec0b7dSPatrick Sanan   *sub_mat_type = osm->sub_mat_type;
88980ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
89080ec0b7dSPatrick Sanan }
89180ec0b7dSPatrick Sanan 
89280ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type)
89380ec0b7dSPatrick Sanan {
89480ec0b7dSPatrick Sanan   PetscErrorCode    ierr;
89580ec0b7dSPatrick Sanan   PC_ASM            *osm = (PC_ASM*)pc->data;
89680ec0b7dSPatrick Sanan 
89780ec0b7dSPatrick Sanan   PetscFunctionBegin;
89880ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
89980ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
90080ec0b7dSPatrick Sanan   ierr = PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);CHKERRQ(ierr);
90180ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
90280ec0b7dSPatrick Sanan }
90380ec0b7dSPatrick Sanan 
9044b9ad928SBarry Smith /*@C
9051093a601SBarry Smith     PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner.
9064b9ad928SBarry Smith 
9074b9ad928SBarry Smith     Collective on PC
9084b9ad928SBarry Smith 
9094b9ad928SBarry Smith     Input Parameters:
9104b9ad928SBarry Smith +   pc - the preconditioner context
9114b9ad928SBarry Smith .   n - the number of subdomains for this processor (default value = 1)
9128c03b21aSDmitry Karpeev .   is - the index set that defines the subdomains for this processor
9130298fd71SBarry Smith          (or NULL for PETSc to determine subdomains)
914f1ee410cSBarry 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
915f1ee410cSBarry Smith          (or NULL to not provide these)
9164b9ad928SBarry Smith 
9174b9ad928SBarry Smith     Notes:
9181093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
9194b9ad928SBarry Smith 
9204b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
9214b9ad928SBarry Smith 
9224b9ad928SBarry Smith     Use PCASMSetTotalSubdomains() to set the subdomains for all processors.
9234b9ad928SBarry Smith 
924f1ee410cSBarry Smith     If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated
925f1ee410cSBarry Smith     back to form the global solution (this is the standard restricted additive Schwarz method)
926f1ee410cSBarry Smith 
927f1ee410cSBarry 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
928f1ee410cSBarry Smith     no code to handle that case.
929f1ee410cSBarry Smith 
9304b9ad928SBarry Smith     Level: advanced
9314b9ad928SBarry Smith 
9324b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz
9334b9ad928SBarry Smith 
9344b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
935f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType()
9364b9ad928SBarry Smith @*/
9377087cfbeSBarry Smith PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
9384b9ad928SBarry Smith {
9397bb14e67SBarry Smith   PetscErrorCode ierr;
9404b9ad928SBarry Smith 
9414b9ad928SBarry Smith   PetscFunctionBegin;
9420700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9437bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr);
9444b9ad928SBarry Smith   PetscFunctionReturn(0);
9454b9ad928SBarry Smith }
9464b9ad928SBarry Smith 
9474b9ad928SBarry Smith /*@C
948feb221f8SDmitry Karpeev     PCASMSetTotalSubdomains - Sets the subdomains for all processors for the
9494b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
9504b9ad928SBarry Smith     PC communicator must call this routine, with the same index sets.
9514b9ad928SBarry Smith 
9524b9ad928SBarry Smith     Collective on PC
9534b9ad928SBarry Smith 
9544b9ad928SBarry Smith     Input Parameters:
9554b9ad928SBarry Smith +   pc - the preconditioner context
956feb221f8SDmitry Karpeev .   N  - the number of subdomains for all processors
957feb221f8SDmitry Karpeev .   is - the index sets that define the subdomains for all processors
958dfaa0c5fSBarry Smith          (or NULL to ask PETSc to determine the subdomains)
9592b691e39SMatthew Knepley -   is_local - the index sets that define the local part of the subdomains for this processor
960f1ee410cSBarry Smith          (or NULL to not provide this information)
9614b9ad928SBarry Smith 
9624b9ad928SBarry Smith     Options Database Key:
9634b9ad928SBarry Smith     To set the total number of subdomain blocks rather than specify the
9644b9ad928SBarry Smith     index sets, use the option
9654b9ad928SBarry Smith .    -pc_asm_blocks <blks> - Sets total blocks
9664b9ad928SBarry Smith 
9674b9ad928SBarry Smith     Notes:
968f1ee410cSBarry Smith     Currently you cannot use this to set the actual subdomains with the argument is or is_local.
9694b9ad928SBarry Smith 
9704b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
9714b9ad928SBarry Smith 
9724b9ad928SBarry Smith     These index sets cannot be destroyed until after completion of the
9734b9ad928SBarry Smith     linear solves for which the ASM preconditioner is being used.
9744b9ad928SBarry Smith 
9754b9ad928SBarry Smith     Use PCASMSetLocalSubdomains() to set local subdomains.
9764b9ad928SBarry Smith 
9771093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
9781093a601SBarry Smith 
9794b9ad928SBarry Smith     Level: advanced
9804b9ad928SBarry Smith 
9814b9ad928SBarry Smith .keywords: PC, ASM, set, total, global, subdomains, additive Schwarz
9824b9ad928SBarry Smith 
9834b9ad928SBarry Smith .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
9844b9ad928SBarry Smith           PCASMCreateSubdomains2D()
9854b9ad928SBarry Smith @*/
9867087cfbeSBarry Smith PetscErrorCode  PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[])
9874b9ad928SBarry Smith {
9887bb14e67SBarry Smith   PetscErrorCode ierr;
9894b9ad928SBarry Smith 
9904b9ad928SBarry Smith   PetscFunctionBegin;
9910700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9927bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr);
9934b9ad928SBarry Smith   PetscFunctionReturn(0);
9944b9ad928SBarry Smith }
9954b9ad928SBarry Smith 
9964b9ad928SBarry Smith /*@
9974b9ad928SBarry Smith     PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
9984b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
999f1ee410cSBarry Smith     PC communicator must call this routine.
10004b9ad928SBarry Smith 
1001ad4df100SBarry Smith     Logically Collective on PC
10024b9ad928SBarry Smith 
10034b9ad928SBarry Smith     Input Parameters:
10044b9ad928SBarry Smith +   pc  - the preconditioner context
10054b9ad928SBarry Smith -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
10064b9ad928SBarry Smith 
10074b9ad928SBarry Smith     Options Database Key:
10084b9ad928SBarry Smith .   -pc_asm_overlap <ovl> - Sets overlap
10094b9ad928SBarry Smith 
10104b9ad928SBarry Smith     Notes:
10114b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.  To use
10124b9ad928SBarry Smith     multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
10134b9ad928SBarry Smith     PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).
10144b9ad928SBarry Smith 
10154b9ad928SBarry Smith     The overlap defaults to 1, so if one desires that no additional
10164b9ad928SBarry Smith     overlap be computed beyond what may have been set with a call to
10174b9ad928SBarry Smith     PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
10184b9ad928SBarry Smith     must be set to be 0.  In particular, if one does not explicitly set
10194b9ad928SBarry Smith     the subdomains an application code, then all overlap would be computed
10204b9ad928SBarry Smith     internally by PETSc, and using an overlap of 0 would result in an ASM
10214b9ad928SBarry Smith     variant that is equivalent to the block Jacobi preconditioner.
10224b9ad928SBarry Smith 
1023f1ee410cSBarry Smith     The default algorithm used by PETSc to increase overlap is fast, but not scalable,
1024f1ee410cSBarry Smith     use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
1025f1ee410cSBarry Smith 
10264b9ad928SBarry Smith     Note that one can define initial index sets with any overlap via
1027f1ee410cSBarry Smith     PCASMSetLocalSubdomains(); the routine
10284b9ad928SBarry Smith     PCASMSetOverlap() merely allows PETSc to extend that overlap further
10294b9ad928SBarry Smith     if desired.
10304b9ad928SBarry Smith 
10314b9ad928SBarry Smith     Level: intermediate
10324b9ad928SBarry Smith 
10334b9ad928SBarry Smith .keywords: PC, ASM, set, overlap
10344b9ad928SBarry Smith 
10354b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
1036f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap()
10374b9ad928SBarry Smith @*/
10387087cfbeSBarry Smith PetscErrorCode  PCASMSetOverlap(PC pc,PetscInt ovl)
10394b9ad928SBarry Smith {
10407bb14e67SBarry Smith   PetscErrorCode ierr;
10414b9ad928SBarry Smith 
10424b9ad928SBarry Smith   PetscFunctionBegin;
10430700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1044c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,ovl,2);
10457bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
10464b9ad928SBarry Smith   PetscFunctionReturn(0);
10474b9ad928SBarry Smith }
10484b9ad928SBarry Smith 
10494b9ad928SBarry Smith /*@
10504b9ad928SBarry Smith     PCASMSetType - Sets the type of restriction and interpolation used
10514b9ad928SBarry Smith     for local problems in the additive Schwarz method.
10524b9ad928SBarry Smith 
1053ad4df100SBarry Smith     Logically Collective on PC
10544b9ad928SBarry Smith 
10554b9ad928SBarry Smith     Input Parameters:
10564b9ad928SBarry Smith +   pc  - the preconditioner context
10574b9ad928SBarry Smith -   type - variant of ASM, one of
10584b9ad928SBarry Smith .vb
10594b9ad928SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
10604b9ad928SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
10614b9ad928SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
10624b9ad928SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
10634b9ad928SBarry Smith .ve
10644b9ad928SBarry Smith 
10654b9ad928SBarry Smith     Options Database Key:
10664b9ad928SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
10674b9ad928SBarry Smith 
1068f1ee410cSBarry Smith     Notes: if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected
1069f1ee410cSBarry Smith     to limit the local processor interpolation
1070f1ee410cSBarry Smith 
10714b9ad928SBarry Smith     Level: intermediate
10724b9ad928SBarry Smith 
10734b9ad928SBarry Smith .keywords: PC, ASM, set, type
10744b9ad928SBarry Smith 
10754b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1076f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType()
10774b9ad928SBarry Smith @*/
10787087cfbeSBarry Smith PetscErrorCode  PCASMSetType(PC pc,PCASMType type)
10794b9ad928SBarry Smith {
10807bb14e67SBarry Smith   PetscErrorCode ierr;
10814b9ad928SBarry Smith 
10824b9ad928SBarry Smith   PetscFunctionBegin;
10830700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1084c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
10857bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr);
10864b9ad928SBarry Smith   PetscFunctionReturn(0);
10874b9ad928SBarry Smith }
10884b9ad928SBarry Smith 
1089c60c7ad4SBarry Smith /*@
1090c60c7ad4SBarry Smith     PCASMGetType - Gets the type of restriction and interpolation used
1091c60c7ad4SBarry Smith     for local problems in the additive Schwarz method.
1092c60c7ad4SBarry Smith 
1093c60c7ad4SBarry Smith     Logically Collective on PC
1094c60c7ad4SBarry Smith 
1095c60c7ad4SBarry Smith     Input Parameter:
1096c60c7ad4SBarry Smith .   pc  - the preconditioner context
1097c60c7ad4SBarry Smith 
1098c60c7ad4SBarry Smith     Output Parameter:
1099c60c7ad4SBarry Smith .   type - variant of ASM, one of
1100c60c7ad4SBarry Smith 
1101c60c7ad4SBarry Smith .vb
1102c60c7ad4SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
1103c60c7ad4SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1104c60c7ad4SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1105c60c7ad4SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
1106c60c7ad4SBarry Smith .ve
1107c60c7ad4SBarry Smith 
1108c60c7ad4SBarry Smith     Options Database Key:
1109c60c7ad4SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
1110c60c7ad4SBarry Smith 
1111c60c7ad4SBarry Smith     Level: intermediate
1112c60c7ad4SBarry Smith 
1113c60c7ad4SBarry Smith .keywords: PC, ASM, set, type
1114c60c7ad4SBarry Smith 
1115c60c7ad4SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1116f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType()
1117c60c7ad4SBarry Smith @*/
1118c60c7ad4SBarry Smith PetscErrorCode  PCASMGetType(PC pc,PCASMType *type)
1119c60c7ad4SBarry Smith {
1120c60c7ad4SBarry Smith   PetscErrorCode ierr;
1121c60c7ad4SBarry Smith 
1122c60c7ad4SBarry Smith   PetscFunctionBegin;
1123c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1124c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));CHKERRQ(ierr);
1125c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1126c60c7ad4SBarry Smith }
1127c60c7ad4SBarry Smith 
112812cd4985SMatthew G. Knepley /*@
112912cd4985SMatthew G. Knepley   PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method.
113012cd4985SMatthew G. Knepley 
113112cd4985SMatthew G. Knepley   Logically Collective on PC
113212cd4985SMatthew G. Knepley 
113312cd4985SMatthew G. Knepley   Input Parameters:
113412cd4985SMatthew G. Knepley + pc  - the preconditioner context
113512cd4985SMatthew G. Knepley - type - type of composition, one of
113612cd4985SMatthew G. Knepley .vb
113712cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
113812cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
113912cd4985SMatthew G. Knepley .ve
114012cd4985SMatthew G. Knepley 
114112cd4985SMatthew G. Knepley   Options Database Key:
114212cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
114312cd4985SMatthew G. Knepley 
114412cd4985SMatthew G. Knepley   Level: intermediate
114512cd4985SMatthew G. Knepley 
1146f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
114712cd4985SMatthew G. Knepley @*/
114812cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type)
114912cd4985SMatthew G. Knepley {
115012cd4985SMatthew G. Knepley   PetscErrorCode ierr;
115112cd4985SMatthew G. Knepley 
115212cd4985SMatthew G. Knepley   PetscFunctionBegin;
115312cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
115412cd4985SMatthew G. Knepley   PetscValidLogicalCollectiveEnum(pc, type, 2);
115512cd4985SMatthew G. Knepley   ierr = PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));CHKERRQ(ierr);
115612cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
115712cd4985SMatthew G. Knepley }
115812cd4985SMatthew G. Knepley 
115912cd4985SMatthew G. Knepley /*@
116012cd4985SMatthew G. Knepley   PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method.
116112cd4985SMatthew G. Knepley 
116212cd4985SMatthew G. Knepley   Logically Collective on PC
116312cd4985SMatthew G. Knepley 
116412cd4985SMatthew G. Knepley   Input Parameter:
116512cd4985SMatthew G. Knepley . pc  - the preconditioner context
116612cd4985SMatthew G. Knepley 
116712cd4985SMatthew G. Knepley   Output Parameter:
116812cd4985SMatthew G. Knepley . type - type of composition, one of
116912cd4985SMatthew G. Knepley .vb
117012cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
117112cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
117212cd4985SMatthew G. Knepley .ve
117312cd4985SMatthew G. Knepley 
117412cd4985SMatthew G. Knepley   Options Database Key:
117512cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
117612cd4985SMatthew G. Knepley 
117712cd4985SMatthew G. Knepley   Level: intermediate
117812cd4985SMatthew G. Knepley 
1179f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
118012cd4985SMatthew G. Knepley @*/
118112cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type)
118212cd4985SMatthew G. Knepley {
118312cd4985SMatthew G. Knepley   PetscErrorCode ierr;
118412cd4985SMatthew G. Knepley 
118512cd4985SMatthew G. Knepley   PetscFunctionBegin;
118612cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
118712cd4985SMatthew G. Knepley   PetscValidPointer(type, 2);
118812cd4985SMatthew G. Knepley   ierr = PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));CHKERRQ(ierr);
118912cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
119012cd4985SMatthew G. Knepley }
119112cd4985SMatthew G. Knepley 
11926ed231c7SMatthew Knepley /*@
11936ed231c7SMatthew Knepley     PCASMSetSortIndices - Determines whether subdomain indices are sorted.
11946ed231c7SMatthew Knepley 
1195ad4df100SBarry Smith     Logically Collective on PC
11966ed231c7SMatthew Knepley 
11976ed231c7SMatthew Knepley     Input Parameters:
11986ed231c7SMatthew Knepley +   pc  - the preconditioner context
11996ed231c7SMatthew Knepley -   doSort - sort the subdomain indices
12006ed231c7SMatthew Knepley 
12016ed231c7SMatthew Knepley     Level: intermediate
12026ed231c7SMatthew Knepley 
12036ed231c7SMatthew Knepley .keywords: PC, ASM, set, type
12046ed231c7SMatthew Knepley 
12056ed231c7SMatthew Knepley .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
12066ed231c7SMatthew Knepley           PCASMCreateSubdomains2D()
12076ed231c7SMatthew Knepley @*/
12087087cfbeSBarry Smith PetscErrorCode  PCASMSetSortIndices(PC pc,PetscBool doSort)
12096ed231c7SMatthew Knepley {
12107bb14e67SBarry Smith   PetscErrorCode ierr;
12116ed231c7SMatthew Knepley 
12126ed231c7SMatthew Knepley   PetscFunctionBegin;
12130700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1214acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
12157bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
12166ed231c7SMatthew Knepley   PetscFunctionReturn(0);
12176ed231c7SMatthew Knepley }
12186ed231c7SMatthew Knepley 
12194b9ad928SBarry Smith /*@C
12204b9ad928SBarry Smith    PCASMGetSubKSP - Gets the local KSP contexts for all blocks on
12214b9ad928SBarry Smith    this processor.
12224b9ad928SBarry Smith 
12234b9ad928SBarry Smith    Collective on PC iff first_local is requested
12244b9ad928SBarry Smith 
12254b9ad928SBarry Smith    Input Parameter:
12264b9ad928SBarry Smith .  pc - the preconditioner context
12274b9ad928SBarry Smith 
12284b9ad928SBarry Smith    Output Parameters:
12290298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
12300298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
12310298fd71SBarry Smith                  all processors must request or all must pass NULL
12324b9ad928SBarry Smith -  ksp - the array of KSP contexts
12334b9ad928SBarry Smith 
12344b9ad928SBarry Smith    Note:
1235d29017ddSJed Brown    After PCASMGetSubKSP() the array of KSPes is not to be freed.
12364b9ad928SBarry Smith 
12374b9ad928SBarry Smith    You must call KSPSetUp() before calling PCASMGetSubKSP().
12384b9ad928SBarry Smith 
1239d29017ddSJed Brown    Fortran note:
12402bf68e3eSBarry 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.
1241d29017ddSJed Brown 
12424b9ad928SBarry Smith    Level: advanced
12434b9ad928SBarry Smith 
12444b9ad928SBarry Smith .keywords: PC, ASM, additive Schwarz, get, sub, KSP, context
12454b9ad928SBarry Smith 
12464b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
12474b9ad928SBarry Smith           PCASMCreateSubdomains2D(),
12484b9ad928SBarry Smith @*/
12497087cfbeSBarry Smith PetscErrorCode  PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
12504b9ad928SBarry Smith {
12517bb14e67SBarry Smith   PetscErrorCode ierr;
12524b9ad928SBarry Smith 
12534b9ad928SBarry Smith   PetscFunctionBegin;
12540700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12557bb14e67SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
12564b9ad928SBarry Smith   PetscFunctionReturn(0);
12574b9ad928SBarry Smith }
12584b9ad928SBarry Smith 
12594b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
12604b9ad928SBarry Smith /*MC
12613b09bd56SBarry Smith    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
12624b9ad928SBarry Smith            its own KSP object.
12634b9ad928SBarry Smith 
12644b9ad928SBarry Smith    Options Database Keys:
126549517cdeSBarry Smith +  -pc_asm_blocks <blks> - Sets total blocks
12664b9ad928SBarry Smith .  -pc_asm_overlap <ovl> - Sets overlap
1267f1ee410cSBarry Smith .  -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict
1268f1ee410cSBarry Smith -  -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive
12694b9ad928SBarry Smith 
12703b09bd56SBarry Smith      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
12713b09bd56SBarry Smith       will get a different convergence rate due to the default option of -pc_asm_type restrict. Use
12723b09bd56SBarry Smith       -pc_asm_type basic to use the standard ASM.
12733b09bd56SBarry Smith 
12744b9ad928SBarry Smith    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
1275f1ee410cSBarry Smith      than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor.
12764b9ad928SBarry Smith 
12773b09bd56SBarry Smith      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1278d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
12794b9ad928SBarry Smith 
1280a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCASMGetSubKSP()
12814b9ad928SBarry Smith          and set the options directly on the resulting KSP object (you can access its PC
12824b9ad928SBarry Smith          with KSPGetPC())
12834b9ad928SBarry Smith 
12844b9ad928SBarry Smith    Level: beginner
12854b9ad928SBarry Smith 
12864b9ad928SBarry Smith    Concepts: additive Schwarz method
12874b9ad928SBarry Smith 
1288c582cd25SBarry Smith     References:
128996a0c994SBarry Smith +   1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
129096a0c994SBarry Smith      Courant Institute, New York University Technical report
129196a0c994SBarry Smith -   1. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
129296a0c994SBarry Smith     Cambridge University Press.
1293c582cd25SBarry Smith 
12944b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1295f1ee410cSBarry Smith            PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType()
1296f1ee410cSBarry Smith            PCASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType
1297e09e08ccSBarry Smith 
12984b9ad928SBarry Smith M*/
12994b9ad928SBarry Smith 
13008cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc)
13014b9ad928SBarry Smith {
1302dfbe8321SBarry Smith   PetscErrorCode ierr;
13034b9ad928SBarry Smith   PC_ASM         *osm;
13044b9ad928SBarry Smith 
13054b9ad928SBarry Smith   PetscFunctionBegin;
1306b00a9115SJed Brown   ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr);
13072fa5cd67SKarl Rupp 
13084b9ad928SBarry Smith   osm->n                 = PETSC_DECIDE;
13094b9ad928SBarry Smith   osm->n_local           = 0;
13102b837212SDmitry Karpeev   osm->n_local_true      = PETSC_DECIDE;
13114b9ad928SBarry Smith   osm->overlap           = 1;
13124b9ad928SBarry Smith   osm->ksp               = 0;
13132b691e39SMatthew Knepley   osm->restriction       = 0;
1314b9845e0eSMatthew Knepley   osm->localization      = 0;
13152b691e39SMatthew Knepley   osm->prolongation      = 0;
1316*5b3c0d42Seaulisa   osm->lprolongation     = 0;
131792bb6962SLisandro Dalcin   osm->x                 = 0;
131892bb6962SLisandro Dalcin   osm->y                 = 0;
131953888de8SMatthew Knepley   osm->y_local           = 0;
13204b9ad928SBarry Smith   osm->is                = 0;
13212b691e39SMatthew Knepley   osm->is_local          = 0;
13224b9ad928SBarry Smith   osm->mat               = 0;
13234b9ad928SBarry Smith   osm->pmat              = 0;
13244b9ad928SBarry Smith   osm->type              = PC_ASM_RESTRICT;
132512cd4985SMatthew G. Knepley   osm->loctype           = PC_COMPOSITE_ADDITIVE;
13264b9ad928SBarry Smith   osm->same_local_solves = PETSC_TRUE;
13276ed231c7SMatthew Knepley   osm->sort_indices      = PETSC_TRUE;
1328d709ea83SDmitry Karpeev   osm->dm_subdomains     = PETSC_FALSE;
132980ec0b7dSPatrick Sanan   osm->sub_mat_type      = NULL;
13304b9ad928SBarry Smith 
133192bb6962SLisandro Dalcin   pc->data                 = (void*)osm;
13324b9ad928SBarry Smith   pc->ops->apply           = PCApply_ASM;
13334b9ad928SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_ASM;
13344b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_ASM;
1335e91c6855SBarry Smith   pc->ops->reset           = PCReset_ASM;
13364b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_ASM;
13374b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_ASM;
13384b9ad928SBarry Smith   pc->ops->setuponblocks   = PCSetUpOnBlocks_ASM;
13394b9ad928SBarry Smith   pc->ops->view            = PCView_ASM;
13404b9ad928SBarry Smith   pc->ops->applyrichardson = 0;
13414b9ad928SBarry Smith 
1342bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr);
1343bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr);
1344bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);CHKERRQ(ierr);
1345bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);CHKERRQ(ierr);
1346c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);CHKERRQ(ierr);
134712cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);CHKERRQ(ierr);
134812cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);CHKERRQ(ierr);
1349bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);CHKERRQ(ierr);
1350bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);CHKERRQ(ierr);
135180ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);CHKERRQ(ierr);
135280ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);CHKERRQ(ierr);
13534b9ad928SBarry Smith   PetscFunctionReturn(0);
13544b9ad928SBarry Smith }
13554b9ad928SBarry Smith 
135692bb6962SLisandro Dalcin /*@C
135792bb6962SLisandro Dalcin    PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
135892bb6962SLisandro Dalcin    preconditioner for a any problem on a general grid.
135992bb6962SLisandro Dalcin 
136092bb6962SLisandro Dalcin    Collective
136192bb6962SLisandro Dalcin 
136292bb6962SLisandro Dalcin    Input Parameters:
136392bb6962SLisandro Dalcin +  A - The global matrix operator
136492bb6962SLisandro Dalcin -  n - the number of local blocks
136592bb6962SLisandro Dalcin 
136692bb6962SLisandro Dalcin    Output Parameters:
136792bb6962SLisandro Dalcin .  outis - the array of index sets defining the subdomains
136892bb6962SLisandro Dalcin 
136992bb6962SLisandro Dalcin    Level: advanced
137092bb6962SLisandro Dalcin 
13717d6bfa3bSBarry Smith    Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap
13727d6bfa3bSBarry Smith     from these if you use PCASMSetLocalSubdomains()
13737d6bfa3bSBarry Smith 
13747d6bfa3bSBarry Smith     In the Fortran version you must provide the array outis[] already allocated of length n.
13757d6bfa3bSBarry Smith 
137692bb6962SLisandro Dalcin .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid
137792bb6962SLisandro Dalcin 
137892bb6962SLisandro Dalcin .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains()
137992bb6962SLisandro Dalcin @*/
13807087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
138192bb6962SLisandro Dalcin {
138292bb6962SLisandro Dalcin   MatPartitioning mpart;
138392bb6962SLisandro Dalcin   const char      *prefix;
1384c56a70eeSBarry Smith   void            (*f)(void);
138592bb6962SLisandro Dalcin   PetscInt        i,j,rstart,rend,bs;
138611bd1e4dSLisandro Dalcin   PetscBool       isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
13870298fd71SBarry Smith   Mat             Ad     = NULL, adj;
138892bb6962SLisandro Dalcin   IS              ispart,isnumb,*is;
138992bb6962SLisandro Dalcin   PetscErrorCode  ierr;
139092bb6962SLisandro Dalcin 
139192bb6962SLisandro Dalcin   PetscFunctionBegin;
13920700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
139392bb6962SLisandro Dalcin   PetscValidPointer(outis,3);
1394c1235816SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
139592bb6962SLisandro Dalcin 
139692bb6962SLisandro Dalcin   /* Get prefix, row distribution, and block size */
139792bb6962SLisandro Dalcin   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
139892bb6962SLisandro Dalcin   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
139992bb6962SLisandro Dalcin   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
140065e19b50SBarry 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);
140165e19b50SBarry Smith 
140292bb6962SLisandro Dalcin   /* Get diagonal block from matrix if possible */
1403c56a70eeSBarry Smith   ierr = MatShellGetOperation(A,MATOP_GET_DIAGONAL_BLOCK,&f);CHKERRQ(ierr);
140492bb6962SLisandro Dalcin   if (f) {
140511bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
140692bb6962SLisandro Dalcin   }
140792bb6962SLisandro Dalcin   if (Ad) {
1408251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1409251f4c67SDmitry Karpeev     if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
141092bb6962SLisandro Dalcin   }
141192bb6962SLisandro Dalcin   if (Ad && n > 1) {
1412ace3abfcSBarry Smith     PetscBool match,done;
141392bb6962SLisandro Dalcin     /* Try to setup a good matrix partitioning if available */
141492bb6962SLisandro Dalcin     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
141592bb6962SLisandro Dalcin     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
141692bb6962SLisandro Dalcin     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1417251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
141892bb6962SLisandro Dalcin     if (!match) {
1419251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
142092bb6962SLisandro Dalcin     }
142192bb6962SLisandro Dalcin     if (!match) { /* assume a "good" partitioner is available */
14221a83f524SJed Brown       PetscInt       na;
14231a83f524SJed Brown       const PetscInt *ia,*ja;
142492bb6962SLisandro Dalcin       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
142592bb6962SLisandro Dalcin       if (done) {
142692bb6962SLisandro Dalcin         /* Build adjacency matrix by hand. Unfortunately a call to
142792bb6962SLisandro Dalcin            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
142892bb6962SLisandro Dalcin            remove the block-aij structure and we cannot expect
142992bb6962SLisandro Dalcin            MatPartitioning to split vertices as we need */
14301a83f524SJed Brown         PetscInt       i,j,len,nnz,cnt,*iia=0,*jja=0;
14311a83f524SJed Brown         const PetscInt *row;
143292bb6962SLisandro Dalcin         nnz = 0;
143392bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* count number of nonzeros */
143492bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
143592bb6962SLisandro Dalcin           row = ja + ia[i];
143692bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
143792bb6962SLisandro Dalcin             if (row[j] == i) { /* don't count diagonal */
143892bb6962SLisandro Dalcin               len--; break;
143992bb6962SLisandro Dalcin             }
144092bb6962SLisandro Dalcin           }
144192bb6962SLisandro Dalcin           nnz += len;
144292bb6962SLisandro Dalcin         }
1443854ce69bSBarry Smith         ierr   = PetscMalloc1(na+1,&iia);CHKERRQ(ierr);
1444854ce69bSBarry Smith         ierr   = PetscMalloc1(nnz,&jja);CHKERRQ(ierr);
144592bb6962SLisandro Dalcin         nnz    = 0;
144692bb6962SLisandro Dalcin         iia[0] = 0;
144792bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* fill adjacency */
144892bb6962SLisandro Dalcin           cnt = 0;
144992bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
145092bb6962SLisandro Dalcin           row = ja + ia[i];
145192bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
145292bb6962SLisandro Dalcin             if (row[j] != i) { /* if not diagonal */
145392bb6962SLisandro Dalcin               jja[nnz+cnt++] = row[j];
145492bb6962SLisandro Dalcin             }
145592bb6962SLisandro Dalcin           }
145692bb6962SLisandro Dalcin           nnz     += cnt;
145792bb6962SLisandro Dalcin           iia[i+1] = nnz;
145892bb6962SLisandro Dalcin         }
145992bb6962SLisandro Dalcin         /* Partitioning of the adjacency matrix */
14600298fd71SBarry Smith         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr);
146192bb6962SLisandro Dalcin         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
146292bb6962SLisandro Dalcin         ierr      = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
146392bb6962SLisandro Dalcin         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
146492bb6962SLisandro Dalcin         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
1465fcfd50ebSBarry Smith         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
146692bb6962SLisandro Dalcin         foundpart = PETSC_TRUE;
146792bb6962SLisandro Dalcin       }
146892bb6962SLisandro Dalcin       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
146992bb6962SLisandro Dalcin     }
1470fcfd50ebSBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
147192bb6962SLisandro Dalcin   }
147292bb6962SLisandro Dalcin 
1473785e854fSJed Brown   ierr   = PetscMalloc1(n,&is);CHKERRQ(ierr);
147492bb6962SLisandro Dalcin   *outis = is;
147592bb6962SLisandro Dalcin 
147692bb6962SLisandro Dalcin   if (!foundpart) {
147792bb6962SLisandro Dalcin 
147892bb6962SLisandro Dalcin     /* Partitioning by contiguous chunks of rows */
147992bb6962SLisandro Dalcin 
148092bb6962SLisandro Dalcin     PetscInt mbs   = (rend-rstart)/bs;
148192bb6962SLisandro Dalcin     PetscInt start = rstart;
148292bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
148392bb6962SLisandro Dalcin       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
148492bb6962SLisandro Dalcin       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
148592bb6962SLisandro Dalcin       start += count;
148692bb6962SLisandro Dalcin     }
148792bb6962SLisandro Dalcin 
148892bb6962SLisandro Dalcin   } else {
148992bb6962SLisandro Dalcin 
149092bb6962SLisandro Dalcin     /* Partitioning by adjacency of diagonal block  */
149192bb6962SLisandro Dalcin 
149292bb6962SLisandro Dalcin     const PetscInt *numbering;
149392bb6962SLisandro Dalcin     PetscInt       *count,nidx,*indices,*newidx,start=0;
149492bb6962SLisandro Dalcin     /* Get node count in each partition */
1495785e854fSJed Brown     ierr = PetscMalloc1(n,&count);CHKERRQ(ierr);
149692bb6962SLisandro Dalcin     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
149792bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
149892bb6962SLisandro Dalcin       for (i=0; i<n; i++) count[i] *= bs;
149992bb6962SLisandro Dalcin     }
150092bb6962SLisandro Dalcin     /* Build indices from node numbering */
150192bb6962SLisandro Dalcin     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1502785e854fSJed Brown     ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr);
150392bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
150492bb6962SLisandro Dalcin     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
150592bb6962SLisandro Dalcin     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
150692bb6962SLisandro Dalcin     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
150792bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1508785e854fSJed Brown       ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr);
15092fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
15102fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
15112fa5cd67SKarl Rupp       }
151292bb6962SLisandro Dalcin       ierr    = PetscFree(indices);CHKERRQ(ierr);
151392bb6962SLisandro Dalcin       nidx   *= bs;
151492bb6962SLisandro Dalcin       indices = newidx;
151592bb6962SLisandro Dalcin     }
151692bb6962SLisandro Dalcin     /* Shift to get global indices */
151792bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] += rstart;
151892bb6962SLisandro Dalcin 
151992bb6962SLisandro Dalcin     /* Build the index sets for each block */
152092bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
152170b3c8c7SBarry Smith       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
152292bb6962SLisandro Dalcin       ierr   = ISSort(is[i]);CHKERRQ(ierr);
152392bb6962SLisandro Dalcin       start += count[i];
152492bb6962SLisandro Dalcin     }
152592bb6962SLisandro Dalcin 
15263bf036e2SBarry Smith     ierr = PetscFree(count);CHKERRQ(ierr);
15273bf036e2SBarry Smith     ierr = PetscFree(indices);CHKERRQ(ierr);
1528fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1529fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
153092bb6962SLisandro Dalcin 
153192bb6962SLisandro Dalcin   }
153292bb6962SLisandro Dalcin   PetscFunctionReturn(0);
153392bb6962SLisandro Dalcin }
153492bb6962SLisandro Dalcin 
153592bb6962SLisandro Dalcin /*@C
153692bb6962SLisandro Dalcin    PCASMDestroySubdomains - Destroys the index sets created with
153792bb6962SLisandro Dalcin    PCASMCreateSubdomains(). Should be called after setting subdomains
153892bb6962SLisandro Dalcin    with PCASMSetLocalSubdomains().
153992bb6962SLisandro Dalcin 
154092bb6962SLisandro Dalcin    Collective
154192bb6962SLisandro Dalcin 
154292bb6962SLisandro Dalcin    Input Parameters:
154392bb6962SLisandro Dalcin +  n - the number of index sets
15442b691e39SMatthew Knepley .  is - the array of index sets
15450298fd71SBarry Smith -  is_local - the array of local index sets, can be NULL
154692bb6962SLisandro Dalcin 
154792bb6962SLisandro Dalcin    Level: advanced
154892bb6962SLisandro Dalcin 
154992bb6962SLisandro Dalcin .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid
155092bb6962SLisandro Dalcin 
155192bb6962SLisandro Dalcin .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains()
155292bb6962SLisandro Dalcin @*/
15537087cfbeSBarry Smith PetscErrorCode  PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
155492bb6962SLisandro Dalcin {
155592bb6962SLisandro Dalcin   PetscInt       i;
155692bb6962SLisandro Dalcin   PetscErrorCode ierr;
15575fd66863SKarl Rupp 
155892bb6962SLisandro Dalcin   PetscFunctionBegin;
1559a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1560a35b7b57SDmitry Karpeev   if (is) {
156192bb6962SLisandro Dalcin     PetscValidPointer(is,2);
1562fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); }
156392bb6962SLisandro Dalcin     ierr = PetscFree(is);CHKERRQ(ierr);
1564a35b7b57SDmitry Karpeev   }
15652b691e39SMatthew Knepley   if (is_local) {
15662b691e39SMatthew Knepley     PetscValidPointer(is_local,3);
1567fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); }
15682b691e39SMatthew Knepley     ierr = PetscFree(is_local);CHKERRQ(ierr);
15692b691e39SMatthew Knepley   }
157092bb6962SLisandro Dalcin   PetscFunctionReturn(0);
157192bb6962SLisandro Dalcin }
157292bb6962SLisandro Dalcin 
15734b9ad928SBarry Smith /*@
15744b9ad928SBarry Smith    PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
15754b9ad928SBarry Smith    preconditioner for a two-dimensional problem on a regular grid.
15764b9ad928SBarry Smith 
15774b9ad928SBarry Smith    Not Collective
15784b9ad928SBarry Smith 
15794b9ad928SBarry Smith    Input Parameters:
15804b9ad928SBarry Smith +  m, n - the number of mesh points in the x and y directions
15814b9ad928SBarry Smith .  M, N - the number of subdomains in the x and y directions
15824b9ad928SBarry Smith .  dof - degrees of freedom per node
15834b9ad928SBarry Smith -  overlap - overlap in mesh lines
15844b9ad928SBarry Smith 
15854b9ad928SBarry Smith    Output Parameters:
15864b9ad928SBarry Smith +  Nsub - the number of subdomains created
15873d03bb48SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
15883d03bb48SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
15894b9ad928SBarry Smith 
15904b9ad928SBarry Smith    Note:
15914b9ad928SBarry Smith    Presently PCAMSCreateSubdomains2d() is valid only for sequential
15924b9ad928SBarry Smith    preconditioners.  More general related routines are
15934b9ad928SBarry Smith    PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().
15944b9ad928SBarry Smith 
15954b9ad928SBarry Smith    Level: advanced
15964b9ad928SBarry Smith 
15974b9ad928SBarry Smith .keywords: PC, ASM, additive Schwarz, create, subdomains, 2D, regular grid
15984b9ad928SBarry Smith 
15994b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
16004b9ad928SBarry Smith           PCASMSetOverlap()
16014b9ad928SBarry Smith @*/
16027087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local)
16034b9ad928SBarry Smith {
16043d03bb48SJed Brown   PetscInt       i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer;
16056849ba73SBarry Smith   PetscErrorCode ierr;
160613f74950SBarry Smith   PetscInt       nidx,*idx,loc,ii,jj,count;
16074b9ad928SBarry Smith 
16084b9ad928SBarry Smith   PetscFunctionBegin;
1609e32f2f54SBarry Smith   if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," ");
16104b9ad928SBarry Smith 
16114b9ad928SBarry Smith   *Nsub     = N*M;
1612854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is);CHKERRQ(ierr);
1613854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is_local);CHKERRQ(ierr);
16144b9ad928SBarry Smith   ystart    = 0;
16153d03bb48SJed Brown   loc_outer = 0;
16164b9ad928SBarry Smith   for (i=0; i<N; i++) {
16174b9ad928SBarry Smith     height = n/N + ((n % N) > i); /* height of subdomain */
1618e32f2f54SBarry Smith     if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n");
16194b9ad928SBarry Smith     yleft  = ystart - overlap; if (yleft < 0) yleft = 0;
16204b9ad928SBarry Smith     yright = ystart + height + overlap; if (yright > n) yright = n;
16214b9ad928SBarry Smith     xstart = 0;
16224b9ad928SBarry Smith     for (j=0; j<M; j++) {
16234b9ad928SBarry Smith       width = m/M + ((m % M) > j); /* width of subdomain */
1624e32f2f54SBarry Smith       if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m");
16254b9ad928SBarry Smith       xleft  = xstart - overlap; if (xleft < 0) xleft = 0;
16264b9ad928SBarry Smith       xright = xstart + width + overlap; if (xright > m) xright = m;
16274b9ad928SBarry Smith       nidx   = (xright - xleft)*(yright - yleft);
1628785e854fSJed Brown       ierr   = PetscMalloc1(nidx,&idx);CHKERRQ(ierr);
16294b9ad928SBarry Smith       loc    = 0;
16304b9ad928SBarry Smith       for (ii=yleft; ii<yright; ii++) {
16314b9ad928SBarry Smith         count = m*ii + xleft;
16322fa5cd67SKarl Rupp         for (jj=xleft; jj<xright; jj++) idx[loc++] = count++;
16334b9ad928SBarry Smith       }
163470b3c8c7SBarry Smith       ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr);
16353d03bb48SJed Brown       if (overlap == 0) {
16363d03bb48SJed Brown         ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr);
16372fa5cd67SKarl Rupp 
16383d03bb48SJed Brown         (*is_local)[loc_outer] = (*is)[loc_outer];
16393d03bb48SJed Brown       } else {
16403d03bb48SJed Brown         for (loc=0,ii=ystart; ii<ystart+height; ii++) {
16413d03bb48SJed Brown           for (jj=xstart; jj<xstart+width; jj++) {
16423d03bb48SJed Brown             idx[loc++] = m*ii + jj;
16433d03bb48SJed Brown           }
16443d03bb48SJed Brown         }
164570b3c8c7SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr);
16463d03bb48SJed Brown       }
16474b9ad928SBarry Smith       ierr    = PetscFree(idx);CHKERRQ(ierr);
16484b9ad928SBarry Smith       xstart += width;
16493d03bb48SJed Brown       loc_outer++;
16504b9ad928SBarry Smith     }
16514b9ad928SBarry Smith     ystart += height;
16524b9ad928SBarry Smith   }
16534b9ad928SBarry Smith   for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); }
16544b9ad928SBarry Smith   PetscFunctionReturn(0);
16554b9ad928SBarry Smith }
16564b9ad928SBarry Smith 
16574b9ad928SBarry Smith /*@C
16584b9ad928SBarry Smith     PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
16594b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
16604b9ad928SBarry Smith 
1661ad4df100SBarry Smith     Not Collective
16624b9ad928SBarry Smith 
16634b9ad928SBarry Smith     Input Parameter:
16644b9ad928SBarry Smith .   pc - the preconditioner context
16654b9ad928SBarry Smith 
16664b9ad928SBarry Smith     Output Parameters:
16674b9ad928SBarry Smith +   n - the number of subdomains for this processor (default value = 1)
16682b691e39SMatthew Knepley .   is - the index sets that define the subdomains for this processor
16690298fd71SBarry Smith -   is_local - the index sets that define the local part of the subdomains for this processor (can be NULL)
16704b9ad928SBarry Smith 
16714b9ad928SBarry Smith 
16724b9ad928SBarry Smith     Notes:
16734b9ad928SBarry Smith     The IS numbering is in the parallel, global numbering of the vector.
16744b9ad928SBarry Smith 
16754b9ad928SBarry Smith     Level: advanced
16764b9ad928SBarry Smith 
16774b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz
16784b9ad928SBarry Smith 
16794b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
16804b9ad928SBarry Smith           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices()
16814b9ad928SBarry Smith @*/
16827087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
16834b9ad928SBarry Smith {
16842a808120SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
168592bb6962SLisandro Dalcin   PetscErrorCode ierr;
1686ace3abfcSBarry Smith   PetscBool      match;
16874b9ad928SBarry Smith 
16884b9ad928SBarry Smith   PetscFunctionBegin;
16890700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
16904482741eSBarry Smith   PetscValidIntPointer(n,2);
169192bb6962SLisandro Dalcin   if (is) PetscValidPointer(is,3);
1692251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
16932a808120SBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM");
16944b9ad928SBarry Smith   if (n) *n = osm->n_local_true;
16954b9ad928SBarry Smith   if (is) *is = osm->is;
16962b691e39SMatthew Knepley   if (is_local) *is_local = osm->is_local;
16974b9ad928SBarry Smith   PetscFunctionReturn(0);
16984b9ad928SBarry Smith }
16994b9ad928SBarry Smith 
17004b9ad928SBarry Smith /*@C
17014b9ad928SBarry Smith     PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
17024b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
17034b9ad928SBarry Smith 
1704ad4df100SBarry Smith     Not Collective
17054b9ad928SBarry Smith 
17064b9ad928SBarry Smith     Input Parameter:
17074b9ad928SBarry Smith .   pc - the preconditioner context
17084b9ad928SBarry Smith 
17094b9ad928SBarry Smith     Output Parameters:
17104b9ad928SBarry Smith +   n - the number of matrices for this processor (default value = 1)
17114b9ad928SBarry Smith -   mat - the matrices
17124b9ad928SBarry Smith 
17134b9ad928SBarry Smith     Level: advanced
17144b9ad928SBarry Smith 
1715cf739d55SBarry Smith     Notes: Call after PCSetUp() (or KSPSetUp()) but before PCApply() (or KSPApply()) and before PCSetUpOnBlocks())
1716cf739d55SBarry Smith 
1717cf739d55SBarry Smith            Usually one would use PCSetModifySubmatrices() to change the submatrices in building the preconditioner.
1718cf739d55SBarry Smith 
17194b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz, block Jacobi
17204b9ad928SBarry Smith 
17214b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1722cf739d55SBarry Smith           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubmatrices()
17234b9ad928SBarry Smith @*/
17247087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
17254b9ad928SBarry Smith {
17264b9ad928SBarry Smith   PC_ASM         *osm;
172792bb6962SLisandro Dalcin   PetscErrorCode ierr;
1728ace3abfcSBarry Smith   PetscBool      match;
17294b9ad928SBarry Smith 
17304b9ad928SBarry Smith   PetscFunctionBegin;
17310700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
173292bb6962SLisandro Dalcin   PetscValidIntPointer(n,2);
173392bb6962SLisandro Dalcin   if (mat) PetscValidPointer(mat,3);
1734ce94432eSBarry Smith   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1735251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
173692bb6962SLisandro Dalcin   if (!match) {
173792bb6962SLisandro Dalcin     if (n) *n = 0;
17380298fd71SBarry Smith     if (mat) *mat = NULL;
173992bb6962SLisandro Dalcin   } else {
17404b9ad928SBarry Smith     osm = (PC_ASM*)pc->data;
17414b9ad928SBarry Smith     if (n) *n = osm->n_local_true;
17424b9ad928SBarry Smith     if (mat) *mat = osm->pmat;
174392bb6962SLisandro Dalcin   }
17444b9ad928SBarry Smith   PetscFunctionReturn(0);
17454b9ad928SBarry Smith }
1746d709ea83SDmitry Karpeev 
1747d709ea83SDmitry Karpeev /*@
1748d709ea83SDmitry Karpeev     PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1749f1ee410cSBarry Smith 
1750d709ea83SDmitry Karpeev     Logically Collective
1751d709ea83SDmitry Karpeev 
1752d709ea83SDmitry Karpeev     Input Parameter:
1753d709ea83SDmitry Karpeev +   pc  - the preconditioner
1754d709ea83SDmitry Karpeev -   flg - boolean indicating whether to use subdomains defined by the DM
1755d709ea83SDmitry Karpeev 
1756d709ea83SDmitry Karpeev     Options Database Key:
1757d709ea83SDmitry Karpeev .   -pc_asm_dm_subdomains
1758d709ea83SDmitry Karpeev 
1759d709ea83SDmitry Karpeev     Level: intermediate
1760d709ea83SDmitry Karpeev 
1761d709ea83SDmitry Karpeev     Notes:
1762d709ea83SDmitry Karpeev     PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(),
1763d709ea83SDmitry Karpeev     so setting either of the first two effectively turns the latter off.
1764d709ea83SDmitry Karpeev 
1765d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1766d709ea83SDmitry Karpeev 
1767d709ea83SDmitry Karpeev .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1768d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1769d709ea83SDmitry Karpeev @*/
1770d709ea83SDmitry Karpeev PetscErrorCode  PCASMSetDMSubdomains(PC pc,PetscBool flg)
1771d709ea83SDmitry Karpeev {
1772d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1773d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1774d709ea83SDmitry Karpeev   PetscBool      match;
1775d709ea83SDmitry Karpeev 
1776d709ea83SDmitry Karpeev   PetscFunctionBegin;
1777d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1778d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
1779d709ea83SDmitry Karpeev   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1780d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1781d709ea83SDmitry Karpeev   if (match) {
1782d709ea83SDmitry Karpeev     osm->dm_subdomains = flg;
1783d709ea83SDmitry Karpeev   }
1784d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1785d709ea83SDmitry Karpeev }
1786d709ea83SDmitry Karpeev 
1787d709ea83SDmitry Karpeev /*@
1788d709ea83SDmitry Karpeev     PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1789d709ea83SDmitry Karpeev     Not Collective
1790d709ea83SDmitry Karpeev 
1791d709ea83SDmitry Karpeev     Input Parameter:
1792d709ea83SDmitry Karpeev .   pc  - the preconditioner
1793d709ea83SDmitry Karpeev 
1794d709ea83SDmitry Karpeev     Output Parameter:
1795d709ea83SDmitry Karpeev .   flg - boolean indicating whether to use subdomains defined by the DM
1796d709ea83SDmitry Karpeev 
1797d709ea83SDmitry Karpeev     Level: intermediate
1798d709ea83SDmitry Karpeev 
1799d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1800d709ea83SDmitry Karpeev 
1801d709ea83SDmitry Karpeev .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1802d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1803d709ea83SDmitry Karpeev @*/
1804d709ea83SDmitry Karpeev PetscErrorCode  PCASMGetDMSubdomains(PC pc,PetscBool* flg)
1805d709ea83SDmitry Karpeev {
1806d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1807d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1808d709ea83SDmitry Karpeev   PetscBool      match;
1809d709ea83SDmitry Karpeev 
1810d709ea83SDmitry Karpeev   PetscFunctionBegin;
1811d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1812d709ea83SDmitry Karpeev   PetscValidPointer(flg,2);
1813d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1814d709ea83SDmitry Karpeev   if (match) {
1815d709ea83SDmitry Karpeev     if (flg) *flg = osm->dm_subdomains;
1816d709ea83SDmitry Karpeev   }
1817d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1818d709ea83SDmitry Karpeev }
181980ec0b7dSPatrick Sanan 
182080ec0b7dSPatrick Sanan /*@
182180ec0b7dSPatrick Sanan      PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string.
182280ec0b7dSPatrick Sanan 
182380ec0b7dSPatrick Sanan    Not Collective
182480ec0b7dSPatrick Sanan 
182580ec0b7dSPatrick Sanan    Input Parameter:
182680ec0b7dSPatrick Sanan .  pc - the PC
182780ec0b7dSPatrick Sanan 
182880ec0b7dSPatrick Sanan    Output Parameter:
1829f1ee410cSBarry Smith .  -pc_asm_sub_mat_type - name of matrix type
183080ec0b7dSPatrick Sanan 
183180ec0b7dSPatrick Sanan    Level: advanced
183280ec0b7dSPatrick Sanan 
183380ec0b7dSPatrick Sanan .keywords: PC, PCASM, MatType, set
183480ec0b7dSPatrick Sanan 
183580ec0b7dSPatrick Sanan .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
183680ec0b7dSPatrick Sanan @*/
183780ec0b7dSPatrick Sanan PetscErrorCode  PCASMGetSubMatType(PC pc,MatType *sub_mat_type){
183880ec0b7dSPatrick Sanan   PetscErrorCode ierr;
183980ec0b7dSPatrick Sanan 
184080ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));CHKERRQ(ierr);
184180ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
184280ec0b7dSPatrick Sanan }
184380ec0b7dSPatrick Sanan 
184480ec0b7dSPatrick Sanan /*@
184580ec0b7dSPatrick Sanan      PCASMSetSubMatType - Set the type of matrix used for ASM subsolves
184680ec0b7dSPatrick Sanan 
184780ec0b7dSPatrick Sanan    Collective on Mat
184880ec0b7dSPatrick Sanan 
184980ec0b7dSPatrick Sanan    Input Parameters:
185080ec0b7dSPatrick Sanan +  pc             - the PC object
185180ec0b7dSPatrick Sanan -  sub_mat_type   - matrix type
185280ec0b7dSPatrick Sanan 
185380ec0b7dSPatrick Sanan    Options Database Key:
185480ec0b7dSPatrick 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.
185580ec0b7dSPatrick Sanan 
185680ec0b7dSPatrick Sanan    Notes:
185780ec0b7dSPatrick Sanan    See "${PETSC_DIR}/include/petscmat.h" for available types
185880ec0b7dSPatrick Sanan 
185980ec0b7dSPatrick Sanan   Level: advanced
186080ec0b7dSPatrick Sanan 
186180ec0b7dSPatrick Sanan .keywords: PC, PCASM, MatType, set
186280ec0b7dSPatrick Sanan 
186380ec0b7dSPatrick Sanan .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
186480ec0b7dSPatrick Sanan @*/
186580ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type)
186680ec0b7dSPatrick Sanan {
186780ec0b7dSPatrick Sanan   PetscErrorCode ierr;
186880ec0b7dSPatrick Sanan 
186980ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));CHKERRQ(ierr);
187080ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
187180ec0b7dSPatrick Sanan }
1872