xref: /petsc/src/ksp/pc/impls/asm/asm.c (revision 347574c96490b43456ea59747b83d5ed1c0f64e5)
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 */
385b3c0d42Seaulisa   VecScatter *lprolongation;      /* mapping from subregion to overlapping process subdomain */
39*347574c9Seaulisa   VecScatter lrestriction;        /* mapping from global to overlapping process subdomain*/
404b9ad928SBarry Smith } PC_ASM;
414b9ad928SBarry Smith 
426849ba73SBarry Smith static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer)
434b9ad928SBarry Smith {
4492bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
456849ba73SBarry Smith   PetscErrorCode ierr;
4613f74950SBarry Smith   PetscMPIInt    rank;
474d219a6aSLisandro Dalcin   PetscInt       i,bsz;
48ace3abfcSBarry Smith   PetscBool      iascii,isstring;
494b9ad928SBarry Smith   PetscViewer    sviewer;
504b9ad928SBarry Smith 
514b9ad928SBarry Smith   PetscFunctionBegin;
52251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
53251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
5432077d6dSBarry Smith   if (iascii) {
553d03bb48SJed Brown     char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set";
568caf3d72SBarry Smith     if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof(overlaps),"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);}
578caf3d72SBarry Smith     if (osm->n > 0) {ierr = PetscSNPrintf(blocks,sizeof(blocks),"total subdomain blocks = %D",osm->n);CHKERRQ(ierr);}
58efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  %s, %s\n",blocks,overlaps);CHKERRQ(ierr);
59efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  restriction/interpolation type - %s\n",PCASMTypes[osm->type]);CHKERRQ(ierr);
60e64f0791SPatrick Sanan     if (osm->dm_subdomains) {ierr = PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: using DM to define subdomains\n");CHKERRQ(ierr);}
6112cd4985SMatthew G. Knepley     if (osm->loctype != PC_COMPOSITE_ADDITIVE) {ierr = PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: local solve composition type - %s\n",PCCompositeTypes[osm->loctype]);CHKERRQ(ierr);}
62ce94432eSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
6392bb6962SLisandro Dalcin     if (osm->same_local_solves) {
644d219a6aSLisandro Dalcin       if (osm->ksp) {
654b9ad928SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Local solve is same for all blocks, in the following KSP and PC objects:\n");CHKERRQ(ierr);
66c5e4d11fSDmitry Karpeev         ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
674d219a6aSLisandro Dalcin         if (!rank) {
684b9ad928SBarry Smith           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
6992bb6962SLisandro Dalcin           ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);
704b9ad928SBarry Smith           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
714b9ad928SBarry Smith         }
72c5e4d11fSDmitry Karpeev         ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
734d219a6aSLisandro Dalcin       }
744b9ad928SBarry Smith     } else {
75c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
764d219a6aSLisandro Dalcin       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,osm->n_local_true);CHKERRQ(ierr);
774d219a6aSLisandro Dalcin       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
784b9ad928SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr);
794b9ad928SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
804d219a6aSLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
81c5e4d11fSDmitry Karpeev       ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
82dd2fa690SBarry Smith       for (i=0; i<osm->n_local_true; i++) {
834d219a6aSLisandro Dalcin         ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr);
844d219a6aSLisandro Dalcin         ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr);
8592bb6962SLisandro Dalcin         ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr);
864d219a6aSLisandro Dalcin         ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
874b9ad928SBarry Smith       }
88c5e4d11fSDmitry Karpeev       ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
894b9ad928SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
904b9ad928SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
91c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
924b9ad928SBarry Smith     }
934b9ad928SBarry Smith   } else if (isstring) {
944d219a6aSLisandro Dalcin     ierr = PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCASMTypes[osm->type]);CHKERRQ(ierr);
95c5e4d11fSDmitry Karpeev     ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
9692bb6962SLisandro Dalcin     if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);}
97c5e4d11fSDmitry Karpeev     ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
984b9ad928SBarry Smith   }
994b9ad928SBarry Smith   PetscFunctionReturn(0);
1004b9ad928SBarry Smith }
1014b9ad928SBarry Smith 
10292bb6962SLisandro Dalcin static PetscErrorCode PCASMPrintSubdomains(PC pc)
10392bb6962SLisandro Dalcin {
10492bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
10592bb6962SLisandro Dalcin   const char     *prefix;
10692bb6962SLisandro Dalcin   char           fname[PETSC_MAX_PATH_LEN+1];
107643535aeSDmitry Karpeev   PetscViewer    viewer, sviewer;
108643535aeSDmitry Karpeev   char           *s;
10992bb6962SLisandro Dalcin   PetscInt       i,j,nidx;
11092bb6962SLisandro Dalcin   const PetscInt *idx;
111643535aeSDmitry Karpeev   PetscMPIInt    rank, size;
11292bb6962SLisandro Dalcin   PetscErrorCode ierr;
113846783a0SBarry Smith 
11492bb6962SLisandro Dalcin   PetscFunctionBegin;
115ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr);
116ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr);
11792bb6962SLisandro Dalcin   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
118c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
11992bb6962SLisandro Dalcin   if (fname[0] == 0) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
120ce94432eSBarry Smith   ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr);
121643535aeSDmitry Karpeev   for (i=0; i<osm->n_local; i++) {
122643535aeSDmitry Karpeev     if (i < osm->n_local_true) {
12392bb6962SLisandro Dalcin       ierr = ISGetLocalSize(osm->is[i],&nidx);CHKERRQ(ierr);
12492bb6962SLisandro Dalcin       ierr = ISGetIndices(osm->is[i],&idx);CHKERRQ(ierr);
125643535aeSDmitry Karpeev       /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
126854ce69bSBarry Smith       ierr = PetscMalloc1(16*(nidx+1)+512, &s);CHKERRQ(ierr);
127643535aeSDmitry Karpeev       ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, 16*(nidx+1)+512, &sviewer);CHKERRQ(ierr);
128643535aeSDmitry Karpeev       ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D with overlap:\n", rank, size, i);CHKERRQ(ierr);
12992bb6962SLisandro Dalcin       for (j=0; j<nidx; j++) {
130643535aeSDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr);
13192bb6962SLisandro Dalcin       }
13292bb6962SLisandro Dalcin       ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr);
133643535aeSDmitry Karpeev       ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr);
134643535aeSDmitry Karpeev       ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
135c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
136643535aeSDmitry Karpeev       ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr);
137643535aeSDmitry Karpeev       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
138c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
139643535aeSDmitry Karpeev       ierr = PetscFree(s);CHKERRQ(ierr);
1402b691e39SMatthew Knepley       if (osm->is_local) {
141643535aeSDmitry Karpeev         /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
142854ce69bSBarry Smith         ierr = PetscMalloc1(16*(nidx+1)+512, &s);CHKERRQ(ierr);
143643535aeSDmitry Karpeev         ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, 16*(nidx+1)+512, &sviewer);CHKERRQ(ierr);
14409d011a0SDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D without overlap:\n", rank, size, i);CHKERRQ(ierr);
1452b691e39SMatthew Knepley         ierr = ISGetLocalSize(osm->is_local[i],&nidx);CHKERRQ(ierr);
1462b691e39SMatthew Knepley         ierr = ISGetIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
1472b691e39SMatthew Knepley         for (j=0; j<nidx; j++) {
14809d011a0SDmitry Karpeev           ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr);
1492b691e39SMatthew Knepley         }
1502b691e39SMatthew Knepley         ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
15109d011a0SDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr);
152643535aeSDmitry Karpeev         ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
153c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
154643535aeSDmitry Karpeev         ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr);
155643535aeSDmitry Karpeev         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
156c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
157643535aeSDmitry Karpeev         ierr = PetscFree(s);CHKERRQ(ierr);
158643535aeSDmitry Karpeev       }
1592fa5cd67SKarl Rupp     } else {
160643535aeSDmitry Karpeev       /* Participate in collective viewer calls. */
161c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
162643535aeSDmitry Karpeev       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
163c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
164643535aeSDmitry Karpeev       /* Assume either all ranks have is_local or none do. */
165643535aeSDmitry Karpeev       if (osm->is_local) {
166c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
167643535aeSDmitry Karpeev         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
168c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
169643535aeSDmitry Karpeev       }
170fdc48646SDmitry Karpeev     }
17192bb6962SLisandro Dalcin   }
17292bb6962SLisandro Dalcin   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
173fcfd50ebSBarry Smith   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
17492bb6962SLisandro Dalcin   PetscFunctionReturn(0);
17592bb6962SLisandro Dalcin }
17692bb6962SLisandro Dalcin 
1776849ba73SBarry Smith static PetscErrorCode PCSetUp_ASM(PC pc)
1784b9ad928SBarry Smith {
1794b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
1806849ba73SBarry Smith   PetscErrorCode ierr;
181ace3abfcSBarry Smith   PetscBool      symset,flg;
18287e86712SBarry Smith   PetscInt       i,m,m_local;
1834b9ad928SBarry Smith   MatReuse       scall = MAT_REUSE_MATRIX;
1844b9ad928SBarry Smith   IS             isl;
1854b9ad928SBarry Smith   KSP            ksp;
1864b9ad928SBarry Smith   PC             subpc;
1872dcb1b2aSMatthew Knepley   const char     *prefix,*pprefix;
18823ce1328SBarry Smith   Vec            vec;
1890298fd71SBarry Smith   DM             *domain_dm = NULL;
1904b9ad928SBarry Smith 
1914b9ad928SBarry Smith   PetscFunctionBegin;
1924b9ad928SBarry Smith   if (!pc->setupcalled) {
19392bb6962SLisandro Dalcin 
19492bb6962SLisandro Dalcin     if (!osm->type_set) {
19592bb6962SLisandro Dalcin       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
1962fa5cd67SKarl Rupp       if (symset && flg) osm->type = PC_ASM_BASIC;
19792bb6962SLisandro Dalcin     }
19892bb6962SLisandro Dalcin 
1992b837212SDmitry Karpeev     /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */
2002b837212SDmitry Karpeev     if (osm->n_local_true == PETSC_DECIDE) {
20169ca1f37SDmitry Karpeev       /* no subdomains given */
20265db9045SDmitry Karpeev       /* try pc->dm first, if allowed */
203d709ea83SDmitry Karpeev       if (osm->dm_subdomains && pc->dm) {
204feb221f8SDmitry Karpeev         PetscInt  num_domains, d;
205feb221f8SDmitry Karpeev         char      **domain_names;
2068d4ac253SDmitry Karpeev         IS        *inner_domain_is, *outer_domain_is;
2078d4ac253SDmitry Karpeev         ierr = DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm);CHKERRQ(ierr);
208704f0395SPatrick Sanan         osm->overlap = -1; /* We do not want to increase the overlap of the IS.
209704f0395SPatrick Sanan                               A future improvement of this code might allow one to use
210704f0395SPatrick Sanan                               DM-defined subdomains and also increase the overlap,
211704f0395SPatrick Sanan                               but that is not currently supported */
21269ca1f37SDmitry Karpeev         if (num_domains) {
2138d4ac253SDmitry Karpeev           ierr = PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is);CHKERRQ(ierr);
21469ca1f37SDmitry Karpeev         }
215feb221f8SDmitry Karpeev         for (d = 0; d < num_domains; ++d) {
216a35b7b57SDmitry Karpeev           if (domain_names)    {ierr = PetscFree(domain_names[d]);CHKERRQ(ierr);}
217a35b7b57SDmitry Karpeev           if (inner_domain_is) {ierr = ISDestroy(&inner_domain_is[d]);CHKERRQ(ierr);}
218a35b7b57SDmitry Karpeev           if (outer_domain_is) {ierr = ISDestroy(&outer_domain_is[d]);CHKERRQ(ierr);}
219feb221f8SDmitry Karpeev         }
220feb221f8SDmitry Karpeev         ierr = PetscFree(domain_names);CHKERRQ(ierr);
2218d4ac253SDmitry Karpeev         ierr = PetscFree(inner_domain_is);CHKERRQ(ierr);
2228d4ac253SDmitry Karpeev         ierr = PetscFree(outer_domain_is);CHKERRQ(ierr);
223feb221f8SDmitry Karpeev       }
2242b837212SDmitry Karpeev       if (osm->n_local_true == PETSC_DECIDE) {
22569ca1f37SDmitry Karpeev         /* still no subdomains; use one subdomain per processor */
2262b837212SDmitry Karpeev         osm->n_local_true = 1;
227feb221f8SDmitry Karpeev       }
2282b837212SDmitry Karpeev     }
2292b837212SDmitry Karpeev     { /* determine the global and max number of subdomains */
2306ac3741eSJed Brown       struct {PetscInt max,sum;} inwork,outwork;
231c10200c1SHong Zhang       PetscMPIInt size;
232c10200c1SHong Zhang 
2336ac3741eSJed Brown       inwork.max   = osm->n_local_true;
2346ac3741eSJed Brown       inwork.sum   = osm->n_local_true;
235367daffbSBarry Smith       ierr         = MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
2366ac3741eSJed Brown       osm->n_local = outwork.max;
2376ac3741eSJed Brown       osm->n       = outwork.sum;
238c10200c1SHong Zhang 
239c10200c1SHong Zhang       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
240c10200c1SHong Zhang       if (outwork.max == 1 && outwork.sum == size) {
2417dae84e0SHong Zhang         /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */
242c10200c1SHong Zhang         ierr = MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE);CHKERRQ(ierr);
243c10200c1SHong Zhang       }
2444b9ad928SBarry Smith     }
24578904715SLisandro Dalcin     if (!osm->is) { /* create the index sets */
24678904715SLisandro Dalcin       ierr = PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);CHKERRQ(ierr);
2474b9ad928SBarry Smith     }
248f5234e35SJed Brown     if (osm->n_local_true > 1 && !osm->is_local) {
249785e854fSJed Brown       ierr = PetscMalloc1(osm->n_local_true,&osm->is_local);CHKERRQ(ierr);
250f5234e35SJed Brown       for (i=0; i<osm->n_local_true; i++) {
251f5234e35SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
252f5234e35SJed Brown           ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr);
253f5234e35SJed Brown           ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr);
254f5234e35SJed Brown         } else {
255f5234e35SJed Brown           ierr             = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr);
256f5234e35SJed Brown           osm->is_local[i] = osm->is[i];
257f5234e35SJed Brown         }
258f5234e35SJed Brown       }
259f5234e35SJed Brown     }
26092bb6962SLisandro Dalcin     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
26190d69ab7SBarry Smith     flg  = PETSC_FALSE;
262c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,prefix,"-pc_asm_print_subdomains",&flg,NULL);CHKERRQ(ierr);
26392bb6962SLisandro Dalcin     if (flg) { ierr = PCASMPrintSubdomains(pc);CHKERRQ(ierr); }
2644b9ad928SBarry Smith 
2653d03bb48SJed Brown     if (osm->overlap > 0) {
2664b9ad928SBarry Smith       /* Extend the "overlapping" regions by a number of steps */
26792bb6962SLisandro Dalcin       ierr = MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);CHKERRQ(ierr);
2683d03bb48SJed Brown     }
2696ed231c7SMatthew Knepley     if (osm->sort_indices) {
27092bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2714b9ad928SBarry Smith         ierr = ISSort(osm->is[i]);CHKERRQ(ierr);
2722b691e39SMatthew Knepley         if (osm->is_local) {
2732b691e39SMatthew Knepley           ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr);
2742b691e39SMatthew Knepley         }
2754b9ad928SBarry Smith       }
2766ed231c7SMatthew Knepley     }
2774b9ad928SBarry Smith 
278f6991133SBarry Smith     if (!osm->ksp) {
27978904715SLisandro Dalcin       /* Create the local solvers */
280785e854fSJed Brown       ierr = PetscMalloc1(osm->n_local_true,&osm->ksp);CHKERRQ(ierr);
281feb221f8SDmitry Karpeev       if (domain_dm) {
282feb221f8SDmitry Karpeev         ierr = PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n");CHKERRQ(ierr);
283feb221f8SDmitry Karpeev       }
28492bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2854b9ad928SBarry Smith         ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
286422a814eSBarry Smith         ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr);
2873bb1ff40SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr);
28892bb6962SLisandro Dalcin         ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
2894b9ad928SBarry Smith         ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
2904b9ad928SBarry Smith         ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
2914b9ad928SBarry Smith         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
2924b9ad928SBarry Smith         ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
2934b9ad928SBarry Smith         ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
294feb221f8SDmitry Karpeev         if (domain_dm) {
295feb221f8SDmitry Karpeev           ierr = KSPSetDM(ksp, domain_dm[i]);CHKERRQ(ierr);
296feb221f8SDmitry Karpeev           ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr);
297feb221f8SDmitry Karpeev           ierr = DMDestroy(&domain_dm[i]);CHKERRQ(ierr);
298feb221f8SDmitry Karpeev         }
2994b9ad928SBarry Smith         osm->ksp[i] = ksp;
3004b9ad928SBarry Smith       }
301feb221f8SDmitry Karpeev       if (domain_dm) {
302feb221f8SDmitry Karpeev         ierr = PetscFree(domain_dm);CHKERRQ(ierr);
303feb221f8SDmitry Karpeev       }
304f6991133SBarry Smith     }
305*347574c9Seaulisa     ierr = ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis);CHKERRQ(ierr);
306*347574c9Seaulisa     ierr = ISSortRemoveDups(osm->lis);CHKERRQ(ierr);
307*347574c9Seaulisa     //if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
308fb745f2cSMatthew G. Knepley       PetscInt m;
309fb745f2cSMatthew G. Knepley       ierr = ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis);CHKERRQ(ierr);
310fb745f2cSMatthew G. Knepley       ierr = ISSortRemoveDups(osm->lis);CHKERRQ(ierr);
311fb745f2cSMatthew G. Knepley       ierr = ISGetLocalSize(osm->lis, &m);CHKERRQ(ierr);
312fb745f2cSMatthew G. Knepley       ierr = VecCreateSeq(PETSC_COMM_SELF, m, &osm->lx);CHKERRQ(ierr);
313fb745f2cSMatthew G. Knepley       ierr = VecDuplicate(osm->lx, &osm->ly);CHKERRQ(ierr);
314*347574c9Seaulisa     //}
3154b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
3164b9ad928SBarry Smith   } else {
3174b9ad928SBarry Smith     /*
3184b9ad928SBarry Smith        Destroy the blocks from the previous iteration
3194b9ad928SBarry Smith     */
320e09e08ccSBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
3214b9ad928SBarry Smith       ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
3224b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
3234b9ad928SBarry Smith     }
3244b9ad928SBarry Smith   }
3254b9ad928SBarry Smith 
32692bb6962SLisandro Dalcin   /*
32792bb6962SLisandro Dalcin      Extract out the submatrices
32892bb6962SLisandro Dalcin   */
3297dae84e0SHong Zhang   ierr = MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr);
33092bb6962SLisandro Dalcin   if (scall == MAT_INITIAL_MATRIX) {
33192bb6962SLisandro Dalcin     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
33292bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
3333bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr);
33478904715SLisandro Dalcin       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
33592bb6962SLisandro Dalcin     }
33692bb6962SLisandro Dalcin   }
33780ec0b7dSPatrick Sanan 
33880ec0b7dSPatrick Sanan   /* Convert the types of the submatrices (if needbe) */
33980ec0b7dSPatrick Sanan   if (osm->sub_mat_type) {
34080ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; i++) {
34180ec0b7dSPatrick Sanan       ierr = MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]));CHKERRQ(ierr);
34280ec0b7dSPatrick Sanan     }
34380ec0b7dSPatrick Sanan   }
34480ec0b7dSPatrick Sanan 
34580ec0b7dSPatrick Sanan   if(!pc->setupcalled){
34680ec0b7dSPatrick Sanan     /* Create the local work vectors (from the local matrices) and scatter contexts */
34780ec0b7dSPatrick Sanan     ierr = MatCreateVecs(pc->pmat,&vec,0);CHKERRQ(ierr);
34880ec0b7dSPatrick Sanan     ierr = PetscMalloc1(osm->n_local,&osm->restriction);CHKERRQ(ierr);
349f1ee410cSBarry 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()");
350f1ee410cSBarry Smith     if (osm->is_local && osm->type == PC_ASM_RESTRICT) {ierr = PetscMalloc1(osm->n_local,&osm->localization);CHKERRQ(ierr);}
35180ec0b7dSPatrick Sanan     ierr = PetscMalloc1(osm->n_local,&osm->prolongation);CHKERRQ(ierr);
352*347574c9Seaulisa     //if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE)
3535b3c0d42Seaulisa       ierr = PetscMalloc1(osm->n_local,&osm->lprolongation);CHKERRQ(ierr);
3545b3c0d42Seaulisa 
35580ec0b7dSPatrick Sanan     ierr = PetscMalloc1(osm->n_local,&osm->x);CHKERRQ(ierr);
35680ec0b7dSPatrick Sanan     ierr = PetscMalloc1(osm->n_local,&osm->y);CHKERRQ(ierr);
35780ec0b7dSPatrick Sanan     ierr = PetscMalloc1(osm->n_local,&osm->y_local);CHKERRQ(ierr);
358*347574c9Seaulisa 
359*347574c9Seaulisa     ierr = ISGetLocalSize(osm->lis,&m);CHKERRQ(ierr);
360*347574c9Seaulisa     ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
361*347574c9Seaulisa     ierr = VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->lrestriction);CHKERRQ(ierr);
362*347574c9Seaulisa     ierr = ISDestroy(&isl);CHKERRQ(ierr);
363*347574c9Seaulisa 
364*347574c9Seaulisa 
36580ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; ++i) {
36680ec0b7dSPatrick Sanan       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
36780ec0b7dSPatrick Sanan       ierr = MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);CHKERRQ(ierr);
36880ec0b7dSPatrick Sanan       ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
36980ec0b7dSPatrick Sanan       ierr = VecScatterCreate(vec,osm->is[i],osm->x[i],isl,&osm->restriction[i]);CHKERRQ(ierr);
37080ec0b7dSPatrick Sanan       ierr = ISDestroy(&isl);CHKERRQ(ierr);
37180ec0b7dSPatrick Sanan       ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr);
372f1ee410cSBarry Smith       if (osm->localization) {
37380ec0b7dSPatrick Sanan         ISLocalToGlobalMapping ltog;
37480ec0b7dSPatrick Sanan         IS                     isll;
37580ec0b7dSPatrick Sanan         const PetscInt         *idx_local;
37680ec0b7dSPatrick Sanan         PetscInt               *idx,nout;
37780ec0b7dSPatrick Sanan 
37880ec0b7dSPatrick Sanan         ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],&ltog);CHKERRQ(ierr);
37980ec0b7dSPatrick Sanan         ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr);
38080ec0b7dSPatrick Sanan         ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
38180ec0b7dSPatrick Sanan         ierr = PetscMalloc1(m_local,&idx);CHKERRQ(ierr);
38280ec0b7dSPatrick Sanan         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx);CHKERRQ(ierr);
38380ec0b7dSPatrick Sanan         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
38480ec0b7dSPatrick Sanan         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is");
38580ec0b7dSPatrick Sanan         ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
38680ec0b7dSPatrick Sanan         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
38780ec0b7dSPatrick Sanan         ierr = ISCreateStride(PETSC_COMM_SELF,m_local,0,1,&isl);CHKERRQ(ierr);
38880ec0b7dSPatrick Sanan         ierr = VecCreateSeq(PETSC_COMM_SELF,m_local,&osm->y_local[i]);CHKERRQ(ierr);
38980ec0b7dSPatrick Sanan         ierr = VecScatterCreate(osm->y[i],isll,osm->y_local[i],isl,&osm->localization[i]);CHKERRQ(ierr);
39080ec0b7dSPatrick Sanan         ierr = ISDestroy(&isll);CHKERRQ(ierr);
39180ec0b7dSPatrick Sanan 
39280ec0b7dSPatrick Sanan         ierr = VecScatterCreate(vec,osm->is_local[i],osm->y_local[i],isl,&osm->prolongation[i]);CHKERRQ(ierr);
39380ec0b7dSPatrick Sanan         ierr = ISDestroy(&isl);CHKERRQ(ierr);
39480ec0b7dSPatrick Sanan       } else {
39580ec0b7dSPatrick Sanan         osm->y_local[i] = osm->y[i];
39680ec0b7dSPatrick Sanan         ierr = PetscObjectReference((PetscObject) osm->y[i]);CHKERRQ(ierr);
39780ec0b7dSPatrick Sanan         osm->prolongation[i] = osm->restriction[i];
39880ec0b7dSPatrick Sanan         ierr = PetscObjectReference((PetscObject) osm->restriction[i]);CHKERRQ(ierr);
39980ec0b7dSPatrick Sanan       }
4005b3c0d42Seaulisa 
401*347574c9Seaulisa       //if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE){// && i < osm->n_local_true - 1) {
4025b3c0d42Seaulisa         ISLocalToGlobalMapping ltog;
4035b3c0d42Seaulisa         IS                     isll;
4045b3c0d42Seaulisa         const PetscInt         *idx_is;
4055b3c0d42Seaulisa         PetscInt               *idx_lis,nout;
4065b3c0d42Seaulisa 
4075b3c0d42Seaulisa         ierr = ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);CHKERRQ(ierr);
4085b3c0d42Seaulisa         ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
4095b3c0d42Seaulisa         ierr = ISGetIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
4105b3c0d42Seaulisa         ierr = PetscMalloc1(m,&idx_lis);CHKERRQ(ierr);
4115b3c0d42Seaulisa         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis);CHKERRQ(ierr);
4125b3c0d42Seaulisa         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
4135b3c0d42Seaulisa         if (nout != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis");
4145b3c0d42Seaulisa         ierr = ISRestoreIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
4155b3c0d42Seaulisa         ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
4165b3c0d42Seaulisa         ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
4175b3c0d42Seaulisa         ierr = VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lprolongation[i]);CHKERRQ(ierr);
4185b3c0d42Seaulisa         ierr = ISDestroy(&isll);CHKERRQ(ierr);
4195b3c0d42Seaulisa         ierr = ISDestroy(&isl);CHKERRQ(ierr);
420*347574c9Seaulisa       //}
4215b3c0d42Seaulisa 
42280ec0b7dSPatrick Sanan     }
42380ec0b7dSPatrick Sanan     for (i=osm->n_local_true; i<osm->n_local; i++) {
42480ec0b7dSPatrick Sanan       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&osm->x[i]);CHKERRQ(ierr);
42580ec0b7dSPatrick Sanan       ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr);
42680ec0b7dSPatrick Sanan       ierr = VecDuplicate(osm->x[i],&osm->y_local[i]);CHKERRQ(ierr);
42780ec0b7dSPatrick Sanan       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&isl);CHKERRQ(ierr);
42880ec0b7dSPatrick Sanan       ierr = VecScatterCreate(vec,isl,osm->x[i],isl,&osm->restriction[i]);CHKERRQ(ierr);
429f1ee410cSBarry Smith       if (osm->localization) {
43080ec0b7dSPatrick Sanan         ierr = VecScatterCreate(osm->y[i],isl,osm->y_local[i],isl,&osm->localization[i]);CHKERRQ(ierr);
43180ec0b7dSPatrick Sanan         ierr = VecScatterCreate(vec,isl,osm->x[i],isl,&osm->prolongation[i]);CHKERRQ(ierr);
43280ec0b7dSPatrick Sanan       } else {
43380ec0b7dSPatrick Sanan         osm->prolongation[i] = osm->restriction[i];
43480ec0b7dSPatrick Sanan         ierr                 = PetscObjectReference((PetscObject) osm->restriction[i]);CHKERRQ(ierr);
43580ec0b7dSPatrick Sanan       }
43680ec0b7dSPatrick Sanan       ierr = ISDestroy(&isl);CHKERRQ(ierr);
43780ec0b7dSPatrick Sanan     }
43880ec0b7dSPatrick Sanan     ierr = VecDestroy(&vec);CHKERRQ(ierr);
43980ec0b7dSPatrick Sanan   }
44080ec0b7dSPatrick Sanan 
441fb745f2cSMatthew G. Knepley   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
442235cc4ceSMatthew G. Knepley     IS      *cis;
443235cc4ceSMatthew G. Knepley     PetscInt c;
444235cc4ceSMatthew G. Knepley 
445235cc4ceSMatthew G. Knepley     ierr = PetscMalloc1(osm->n_local_true, &cis);CHKERRQ(ierr);
446235cc4ceSMatthew G. Knepley     for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis;
4477dae84e0SHong Zhang     ierr = MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);CHKERRQ(ierr);
448235cc4ceSMatthew G. Knepley     ierr = PetscFree(cis);CHKERRQ(ierr);
449fb745f2cSMatthew G. Knepley   }
4504b9ad928SBarry Smith 
4514b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
4524b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
45392bb6962SLisandro Dalcin   ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
4544b9ad928SBarry Smith 
45592bb6962SLisandro Dalcin   /*
45692bb6962SLisandro Dalcin      Loop over subdomains putting them into local ksp
45792bb6962SLisandro Dalcin   */
45892bb6962SLisandro Dalcin   for (i=0; i<osm->n_local_true; i++) {
45923ee1639SBarry Smith     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr);
46092bb6962SLisandro Dalcin     if (!pc->setupcalled) {
461bf108f30SBarry Smith       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
4624b9ad928SBarry Smith     }
46392bb6962SLisandro Dalcin   }
4644b9ad928SBarry Smith   PetscFunctionReturn(0);
4654b9ad928SBarry Smith }
4664b9ad928SBarry Smith 
4676849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc)
4684b9ad928SBarry Smith {
4694b9ad928SBarry Smith   PC_ASM             *osm = (PC_ASM*)pc->data;
4706849ba73SBarry Smith   PetscErrorCode     ierr;
47113f74950SBarry Smith   PetscInt           i;
472539c167fSBarry Smith   KSPConvergedReason reason;
4734b9ad928SBarry Smith 
4744b9ad928SBarry Smith   PetscFunctionBegin;
4754b9ad928SBarry Smith   for (i=0; i<osm->n_local_true; i++) {
4764b9ad928SBarry Smith     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
477539c167fSBarry Smith     ierr = KSPGetConvergedReason(osm->ksp[i],&reason);CHKERRQ(ierr);
478539c167fSBarry Smith     if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
479261222b2SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
480e0eafd54SHong Zhang     }
4814b9ad928SBarry Smith   }
4824b9ad928SBarry Smith   PetscFunctionReturn(0);
4834b9ad928SBarry Smith }
4844b9ad928SBarry Smith 
4856849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y)
4864b9ad928SBarry Smith {
4874b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
4886849ba73SBarry Smith   PetscErrorCode ierr;
48913f74950SBarry Smith   PetscInt       i,n_local = osm->n_local,n_local_true = osm->n_local_true;
4904b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
4914b9ad928SBarry Smith 
4924b9ad928SBarry Smith   PetscFunctionBegin;
4934b9ad928SBarry Smith   /*
4944b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
4954b9ad928SBarry Smith      subdomain values (leaving the other values 0).
4964b9ad928SBarry Smith   */
4974b9ad928SBarry Smith   if (!(osm->type & PC_ASM_RESTRICT)) {
4984b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
4994b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
500f21eeb63SBarry Smith     for (i=0; i<n_local_true; i++) {
5013d03bb48SJed Brown       ierr = VecZeroEntries(osm->x[i]);CHKERRQ(ierr);
5024b9ad928SBarry Smith     }
5034b9ad928SBarry Smith   }
504*347574c9Seaulisa   if (!(osm->type & PC_ASM_INTERPOLATE)) {
505*347574c9Seaulisa     reverse = SCATTER_REVERSE_LOCAL;
506*347574c9Seaulisa   }
5074b9ad928SBarry Smith 
508*347574c9Seaulisa //   switch (osm->loctype)
509*347574c9Seaulisa //   {
510*347574c9Seaulisa //   case PC_COMPOSITE_ADDITIVE:
511*347574c9Seaulisa //     for (i=0; i<n_local; i++) {
512*347574c9Seaulisa //       ierr = VecScatterBegin(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
513*347574c9Seaulisa //     }
514*347574c9Seaulisa //     ierr = VecZeroEntries(y);CHKERRQ(ierr);
515*347574c9Seaulisa //     /* do the local solves */
516*347574c9Seaulisa //     for (i=0; i<n_local_true; i++) {
517*347574c9Seaulisa //       ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
518*347574c9Seaulisa //       ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
519*347574c9Seaulisa //       if (osm->localization) {
520*347574c9Seaulisa //         ierr = VecScatterBegin(osm->localization[i],osm->y[i],osm->y_local[i],INSERT_VALUES,forward);CHKERRQ(ierr);
521*347574c9Seaulisa //         ierr = VecScatterEnd(osm->localization[i],osm->y[i],osm->y_local[i],INSERT_VALUES,forward);CHKERRQ(ierr);
522*347574c9Seaulisa //       }
523*347574c9Seaulisa //       ierr = VecScatterBegin(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
524*347574c9Seaulisa //     }
525*347574c9Seaulisa //     /* handle the rest of the scatters that do not have local solves */
526*347574c9Seaulisa //     for (i=n_local_true; i<n_local; i++) {
527*347574c9Seaulisa //       ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
528*347574c9Seaulisa //       ierr = VecScatterBegin(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
529*347574c9Seaulisa //     }
530*347574c9Seaulisa //     for (i=0; i<n_local; i++) {
531*347574c9Seaulisa //       ierr = VecScatterEnd(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
532*347574c9Seaulisa //     }
533*347574c9Seaulisa //     break;
534*347574c9Seaulisa //   case PC_COMPOSITE_MULTIPLICATIVE:
535*347574c9Seaulisa //     ierr = VecZeroEntries(y);CHKERRQ(ierr);
536*347574c9Seaulisa //     ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
537*347574c9Seaulisa //     /* do the local solves */
538*347574c9Seaulisa //     for (i = 0; i < n_local_true; ++i) {
539*347574c9Seaulisa //       if (i == 0) {
540*347574c9Seaulisa //         ierr = VecZeroEntries(osm->x[i]);CHKERRQ(ierr);
541*347574c9Seaulisa //       }
542*347574c9Seaulisa //       ierr = VecScatterBegin(osm->restriction[i], x, osm->x[i], ADD_VALUES, forward);CHKERRQ(ierr);
543*347574c9Seaulisa //       ierr = VecScatterEnd(osm->restriction[i], x, osm->x[i], ADD_VALUES, forward);CHKERRQ(ierr);
544*347574c9Seaulisa //       ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
545*347574c9Seaulisa //       if (osm->localization) {
546*347574c9Seaulisa //         ierr = VecScatterBegin(osm->localization[i], osm->y[i], osm->y_local[i], INSERT_VALUES, forward);CHKERRQ(ierr);
547*347574c9Seaulisa //         ierr = VecScatterEnd(osm->localization[i], osm->y[i], osm->y_local[i], INSERT_VALUES, forward);CHKERRQ(ierr);
548*347574c9Seaulisa //       }
549*347574c9Seaulisa //       ierr = VecScatterBegin(osm->prolongation[i], osm->y_local[i], y, ADD_VALUES, reverse);CHKERRQ(ierr);
550*347574c9Seaulisa //       ierr = VecScatterEnd(osm->prolongation[i], osm->y_local[i], y, ADD_VALUES, reverse);CHKERRQ(ierr);
551*347574c9Seaulisa //       if (i < n_local_true-1) {
552*347574c9Seaulisa // 	ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
553*347574c9Seaulisa //         ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
554*347574c9Seaulisa // 	ierr = VecCopy(osm->ly, osm->lx);CHKERRQ(ierr);
555*347574c9Seaulisa //         ierr = VecScale(osm->lx, -1.0);CHKERRQ(ierr);
556*347574c9Seaulisa //         ierr = MatMult(osm->lmats[i+1], osm->lx, osm->x[i+1]);CHKERRQ(ierr);
557*347574c9Seaulisa //       }
558*347574c9Seaulisa //     }
559*347574c9Seaulisa //     /* handle the rest of the scatters that do not have local solves */
560*347574c9Seaulisa //     for (i = n_local_true; i < n_local; ++i) {
561*347574c9Seaulisa //       ierr = VecScatterBegin(osm->restriction[i], x, osm->x[i], INSERT_VALUES, forward);CHKERRQ(ierr);
562*347574c9Seaulisa //       ierr = VecScatterEnd(osm->restriction[i], x, osm->x[i], INSERT_VALUES, forward);CHKERRQ(ierr);
563*347574c9Seaulisa //       ierr = VecScatterBegin(osm->prolongation[i], osm->y_local[i], y, ADD_VALUES, reverse);CHKERRQ(ierr);
564*347574c9Seaulisa //       ierr = VecScatterEnd(osm->prolongation[i], osm->y_local[i], y, ADD_VALUES, reverse);CHKERRQ(ierr);
565*347574c9Seaulisa //     }
566*347574c9Seaulisa   if(osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE){
56712cd4985SMatthew G. Knepley     ierr = VecZeroEntries(y);CHKERRQ(ierr);
5685b3c0d42Seaulisa     ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
569*347574c9Seaulisa 
570*347574c9Seaulisa     ierr = VecScatterBegin(osm->lrestriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
571*347574c9Seaulisa     ierr = VecScatterEnd(osm->lrestriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
572*347574c9Seaulisa 
573*347574c9Seaulisa     ierr = VecScatterBegin(osm->lprolongation[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
574*347574c9Seaulisa     ierr = VecScatterEnd(osm->lprolongation[0], osm->lx, osm->x[0],  INSERT_VALUES, forward);CHKERRQ(ierr);
57512cd4985SMatthew G. Knepley     /* do the local solves */
57612cd4985SMatthew G. Knepley     for (i = 0; i < n_local_true; ++i) {
577*347574c9Seaulisa 
57812cd4985SMatthew G. Knepley       ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
579*347574c9Seaulisa //       if (osm->localization) {
580*347574c9Seaulisa //         ierr = VecScatterBegin(osm->localization[i], osm->y[i], osm->y_local[i], INSERT_VALUES, forward);CHKERRQ(ierr);
581*347574c9Seaulisa //         ierr = VecScatterEnd(osm->localization[i], osm->y[i], osm->y_local[i], INSERT_VALUES, forward);CHKERRQ(ierr);
582*347574c9Seaulisa //       }
5835b3c0d42Seaulisa       ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
5845b3c0d42Seaulisa       ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
585*347574c9Seaulisa 
586*347574c9Seaulisa       if (i < n_local_true-1) {
587*347574c9Seaulisa 	ierr = VecScatterBegin(osm->lprolongation[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
588*347574c9Seaulisa 	ierr = VecScatterEnd(osm->lprolongation[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
589*347574c9Seaulisa 
590*347574c9Seaulisa 	if ( osm->loctype == PC_COMPOSITE_MULTIPLICATIVE){
591*347574c9Seaulisa 	  ierr = MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]);CHKERRQ(ierr);
592*347574c9Seaulisa 	  ierr = VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1]); CHKERRQ(ierr);
5937c3d802fSMatthew G. Knepley 	}
59412cd4985SMatthew G. Knepley       }
59512cd4985SMatthew G. Knepley     }
596*347574c9Seaulisa     ierr = VecScatterBegin(osm->lrestriction, osm->ly, y,  ADD_VALUES, reverse);CHKERRQ(ierr);
597*347574c9Seaulisa     ierr = VecScatterEnd(osm->lrestriction,  osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
598*347574c9Seaulisa 
599*347574c9Seaulisa     //break;
600*347574c9Seaulisa   //default: SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
601*347574c9Seaulisa   }else{
602*347574c9Seaulisa     SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
60312cd4985SMatthew G. Knepley   }
6044b9ad928SBarry Smith   PetscFunctionReturn(0);
6054b9ad928SBarry Smith }
6064b9ad928SBarry Smith 
6076849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y)
6084b9ad928SBarry Smith {
6094b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
6106849ba73SBarry Smith   PetscErrorCode ierr;
61113f74950SBarry Smith   PetscInt       i,n_local = osm->n_local,n_local_true = osm->n_local_true;
6124b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
6134b9ad928SBarry Smith 
6144b9ad928SBarry Smith   PetscFunctionBegin;
6154b9ad928SBarry Smith   /*
6164b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
6174b9ad928SBarry Smith      subdomain values (leaving the other values 0).
6184b9ad928SBarry Smith 
6194b9ad928SBarry Smith      Note: these are reversed from the PCApply_ASM() because we are applying the
6204b9ad928SBarry Smith      transpose of the three terms
6214b9ad928SBarry Smith   */
6224b9ad928SBarry Smith   if (!(osm->type & PC_ASM_INTERPOLATE)) {
6234b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
6244b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
625f21eeb63SBarry Smith     for (i=0; i<n_local_true; i++) {
62610bd88b9SJed Brown       ierr = VecZeroEntries(osm->x[i]);CHKERRQ(ierr);
6274b9ad928SBarry Smith     }
6284b9ad928SBarry Smith   }
6292fa5cd67SKarl Rupp   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;
6304b9ad928SBarry Smith 
6314b9ad928SBarry Smith   for (i=0; i<n_local; i++) {
6322b691e39SMatthew Knepley     ierr = VecScatterBegin(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
6334b9ad928SBarry Smith   }
63410bd88b9SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
6354b9ad928SBarry Smith   /* do the local solves */
6364b9ad928SBarry Smith   for (i=0; i<n_local_true; i++) {
6372b691e39SMatthew Knepley     ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
63823ce1328SBarry Smith     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
639b9845e0eSMatthew Knepley     if (osm->localization) {
64053888de8SMatthew Knepley       ierr = VecScatterBegin(osm->localization[i],osm->y[i],osm->y_local[i],INSERT_VALUES,forward);CHKERRQ(ierr);
64153888de8SMatthew Knepley       ierr = VecScatterEnd(osm->localization[i],osm->y[i],osm->y_local[i],INSERT_VALUES,forward);CHKERRQ(ierr);
642b9845e0eSMatthew Knepley     }
64353888de8SMatthew Knepley     ierr = VecScatterBegin(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
6444b9ad928SBarry Smith   }
6454b9ad928SBarry Smith   /* handle the rest of the scatters that do not have local solves */
6464b9ad928SBarry Smith   for (i=n_local_true; i<n_local; i++) {
6472b691e39SMatthew Knepley     ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr);
64853888de8SMatthew Knepley     ierr = VecScatterBegin(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
6494b9ad928SBarry Smith   }
6504b9ad928SBarry Smith   for (i=0; i<n_local; i++) {
65153888de8SMatthew Knepley     ierr = VecScatterEnd(osm->prolongation[i],osm->y_local[i],y,ADD_VALUES,reverse);CHKERRQ(ierr);
6524b9ad928SBarry Smith   }
6534b9ad928SBarry Smith   PetscFunctionReturn(0);
6544b9ad928SBarry Smith }
6554b9ad928SBarry Smith 
656e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc)
6574b9ad928SBarry Smith {
6584b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
6596849ba73SBarry Smith   PetscErrorCode ierr;
66013f74950SBarry Smith   PetscInt       i;
6614b9ad928SBarry Smith 
6624b9ad928SBarry Smith   PetscFunctionBegin;
66392bb6962SLisandro Dalcin   if (osm->ksp) {
66492bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
665e91c6855SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
66692bb6962SLisandro Dalcin     }
66792bb6962SLisandro Dalcin   }
668e09e08ccSBarry Smith   if (osm->pmat) {
66992bb6962SLisandro Dalcin     if (osm->n_local_true > 0) {
67030a70a9aSHong Zhang       ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
67192bb6962SLisandro Dalcin     }
67292bb6962SLisandro Dalcin   }
6732b691e39SMatthew Knepley   if (osm->restriction) {
674*347574c9Seaulisa     ierr = VecScatterDestroy(&osm->lrestriction);CHKERRQ(ierr);
6754b9ad928SBarry Smith     for (i=0; i<osm->n_local; i++) {
676fcfd50ebSBarry Smith       ierr = VecScatterDestroy(&osm->restriction[i]);CHKERRQ(ierr);
677fcfd50ebSBarry Smith       if (osm->localization) {ierr = VecScatterDestroy(&osm->localization[i]);CHKERRQ(ierr);}
678fcfd50ebSBarry Smith       ierr = VecScatterDestroy(&osm->prolongation[i]);CHKERRQ(ierr);
679*347574c9Seaulisa //       if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE && i < osm->n_local_true){// - 1){
680*347574c9Seaulisa       if (i < osm->n_local_true){// - 1){
6815b3c0d42Seaulisa 	ierr = VecScatterDestroy(&osm->lprolongation[i]);CHKERRQ(ierr);
6825b3c0d42Seaulisa       }
683*347574c9Seaulisa 
684fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
685fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
686fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y_local[i]);CHKERRQ(ierr);
6874b9ad928SBarry Smith     }
6882b691e39SMatthew Knepley     ierr = PetscFree(osm->restriction);CHKERRQ(ierr);
689b9845e0eSMatthew Knepley     if (osm->localization) {ierr = PetscFree(osm->localization);CHKERRQ(ierr);}
6902b691e39SMatthew Knepley     ierr = PetscFree(osm->prolongation);CHKERRQ(ierr);
691*347574c9Seaulisa     //if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE){
6925b3c0d42Seaulisa       ierr = PetscFree(osm->lprolongation);CHKERRQ(ierr);
693*347574c9Seaulisa     //}
69405b42c5fSBarry Smith     ierr = PetscFree(osm->x);CHKERRQ(ierr);
69578904715SLisandro Dalcin     ierr = PetscFree(osm->y);CHKERRQ(ierr);
69653888de8SMatthew Knepley     ierr = PetscFree(osm->y_local);CHKERRQ(ierr);
69792bb6962SLisandro Dalcin   }
6982b691e39SMatthew Knepley   ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
699fb745f2cSMatthew G. Knepley   ierr = ISDestroy(&osm->lis);CHKERRQ(ierr);
700fb745f2cSMatthew G. Knepley   ierr = VecDestroy(&osm->lx);CHKERRQ(ierr);
701fb745f2cSMatthew G. Knepley   ierr = VecDestroy(&osm->ly);CHKERRQ(ierr);
702*347574c9Seaulisa   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
703*347574c9Seaulisa     //ierr = ISDestroy(&osm->lis);CHKERRQ(ierr);
704*347574c9Seaulisa     ierr = MatDestroyMatrices(osm->n_local_true, &osm->lmats);CHKERRQ(ierr);
705*347574c9Seaulisa //     ierr = VecDestroy(&osm->lx);CHKERRQ(ierr);
706*347574c9Seaulisa //     ierr = VecDestroy(&osm->ly);CHKERRQ(ierr);
707fb745f2cSMatthew G. Knepley   }
7082fa5cd67SKarl Rupp 
70980ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
71080ec0b7dSPatrick Sanan 
711e91c6855SBarry Smith   osm->is       = 0;
712e91c6855SBarry Smith   osm->is_local = 0;
713e91c6855SBarry Smith   PetscFunctionReturn(0);
714e91c6855SBarry Smith }
715e91c6855SBarry Smith 
716e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc)
717e91c6855SBarry Smith {
718e91c6855SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
719e91c6855SBarry Smith   PetscErrorCode ierr;
720e91c6855SBarry Smith   PetscInt       i;
721e91c6855SBarry Smith 
722e91c6855SBarry Smith   PetscFunctionBegin;
723e91c6855SBarry Smith   ierr = PCReset_ASM(pc);CHKERRQ(ierr);
724e91c6855SBarry Smith   if (osm->ksp) {
725e91c6855SBarry Smith     for (i=0; i<osm->n_local_true; i++) {
726fcfd50ebSBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
727e91c6855SBarry Smith     }
728e91c6855SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
729e91c6855SBarry Smith   }
730e91c6855SBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
7314b9ad928SBarry Smith   PetscFunctionReturn(0);
7324b9ad928SBarry Smith }
7334b9ad928SBarry Smith 
7344416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc)
7354b9ad928SBarry Smith {
7364b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
7376849ba73SBarry Smith   PetscErrorCode ierr;
7389dcbbd2bSBarry Smith   PetscInt       blocks,ovl;
739ace3abfcSBarry Smith   PetscBool      symset,flg;
74092bb6962SLisandro Dalcin   PCASMType      asmtype;
74112cd4985SMatthew G. Knepley   PCCompositeType loctype;
74280ec0b7dSPatrick Sanan   char           sub_mat_type[256];
7434b9ad928SBarry Smith 
7444b9ad928SBarry Smith   PetscFunctionBegin;
745bf108f30SBarry Smith   /* set the type to symmetric if matrix is symmetric */
74692bb6962SLisandro Dalcin   if (!osm->type_set && pc->pmat) {
74792bb6962SLisandro Dalcin     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
7482fa5cd67SKarl Rupp     if (symset && flg) osm->type = PC_ASM_BASIC;
749bf108f30SBarry Smith   }
750e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");CHKERRQ(ierr);
751d709ea83SDmitry Karpeev   ierr = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr);
75290d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
75365db9045SDmitry Karpeev   if (flg) {
754f77a5242SKarl Rupp     ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr);
755d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
75665db9045SDmitry Karpeev   }
75790d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
75865db9045SDmitry Karpeev   if (flg) {
75965db9045SDmitry Karpeev     ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr);
760d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
76165db9045SDmitry Karpeev   }
76290d69ab7SBarry Smith   flg  = PETSC_FALSE;
76390d69ab7SBarry Smith   ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
76492bb6962SLisandro Dalcin   if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); }
76512cd4985SMatthew G. Knepley   flg  = PETSC_FALSE;
76612cd4985SMatthew G. Knepley   ierr = PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr);
76712cd4985SMatthew G. Knepley   if (flg) {ierr = PCASMSetLocalType(pc,loctype);CHKERRQ(ierr); }
76880ec0b7dSPatrick Sanan   ierr = PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);CHKERRQ(ierr);
76980ec0b7dSPatrick Sanan   if(flg){
770459726d8SSatish Balay     ierr = PCASMSetSubMatType(pc,sub_mat_type);CHKERRQ(ierr);
77180ec0b7dSPatrick Sanan   }
7724b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7734b9ad928SBarry Smith   PetscFunctionReturn(0);
7744b9ad928SBarry Smith }
7754b9ad928SBarry Smith 
7764b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/
7774b9ad928SBarry Smith 
7781e6b0712SBarry Smith static PetscErrorCode  PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[])
7794b9ad928SBarry Smith {
7804b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
78192bb6962SLisandro Dalcin   PetscErrorCode ierr;
78292bb6962SLisandro Dalcin   PetscInt       i;
7834b9ad928SBarry Smith 
7844b9ad928SBarry Smith   PetscFunctionBegin;
785e32f2f54SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
786ce94432eSBarry 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().");
787e7e72b3dSBarry Smith 
7884b9ad928SBarry Smith   if (!pc->setupcalled) {
78992bb6962SLisandro Dalcin     if (is) {
79092bb6962SLisandro Dalcin       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
79192bb6962SLisandro Dalcin     }
792832fc9a5SMatthew Knepley     if (is_local) {
793832fc9a5SMatthew Knepley       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
794832fc9a5SMatthew Knepley     }
7952b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
7962fa5cd67SKarl Rupp 
7974b9ad928SBarry Smith     osm->n_local_true = n;
79892bb6962SLisandro Dalcin     osm->is           = 0;
7992b691e39SMatthew Knepley     osm->is_local     = 0;
80092bb6962SLisandro Dalcin     if (is) {
801785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is);CHKERRQ(ierr);
8022fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is[i] = is[i];
8033d03bb48SJed Brown       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
8043d03bb48SJed Brown       osm->overlap = -1;
80592bb6962SLisandro Dalcin     }
8062b691e39SMatthew Knepley     if (is_local) {
807785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is_local);CHKERRQ(ierr);
8082fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is_local[i] = is_local[i];
809a35b7b57SDmitry Karpeev       if (!is) {
810785e854fSJed Brown         ierr = PetscMalloc1(osm->n_local_true,&osm->is);CHKERRQ(ierr);
811a35b7b57SDmitry Karpeev         for (i=0; i<osm->n_local_true; i++) {
812a35b7b57SDmitry Karpeev           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
813a35b7b57SDmitry Karpeev             ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr);
814a35b7b57SDmitry Karpeev             ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr);
815a35b7b57SDmitry Karpeev           } else {
816a35b7b57SDmitry Karpeev             ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr);
817a35b7b57SDmitry Karpeev             osm->is[i] = osm->is_local[i];
818a35b7b57SDmitry Karpeev           }
819a35b7b57SDmitry Karpeev         }
820a35b7b57SDmitry Karpeev       }
8212b691e39SMatthew Knepley     }
8224b9ad928SBarry Smith   }
8234b9ad928SBarry Smith   PetscFunctionReturn(0);
8244b9ad928SBarry Smith }
8254b9ad928SBarry Smith 
8261e6b0712SBarry Smith static PetscErrorCode  PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local)
8274b9ad928SBarry Smith {
8284b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
8296849ba73SBarry Smith   PetscErrorCode ierr;
83013f74950SBarry Smith   PetscMPIInt    rank,size;
83178904715SLisandro Dalcin   PetscInt       n;
8324b9ad928SBarry Smith 
8334b9ad928SBarry Smith   PetscFunctionBegin;
834ce94432eSBarry Smith   if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
835ce94432eSBarry 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.");
8364b9ad928SBarry Smith 
8374b9ad928SBarry Smith   /*
838880770e9SJed Brown      Split the subdomains equally among all processors
8394b9ad928SBarry Smith   */
840ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
841ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
8424b9ad928SBarry Smith   n    = N/size + ((N % size) > rank);
843e32f2f54SBarry 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);
844e32f2f54SBarry Smith   if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp().");
8454b9ad928SBarry Smith   if (!pc->setupcalled) {
8462b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
8472fa5cd67SKarl Rupp 
8484b9ad928SBarry Smith     osm->n_local_true = n;
8494b9ad928SBarry Smith     osm->is           = 0;
8502b691e39SMatthew Knepley     osm->is_local     = 0;
8514b9ad928SBarry Smith   }
8524b9ad928SBarry Smith   PetscFunctionReturn(0);
8534b9ad928SBarry Smith }
8544b9ad928SBarry Smith 
8551e6b0712SBarry Smith static PetscErrorCode  PCASMSetOverlap_ASM(PC pc,PetscInt ovl)
8564b9ad928SBarry Smith {
85792bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
8584b9ad928SBarry Smith 
8594b9ad928SBarry Smith   PetscFunctionBegin;
860ce94432eSBarry Smith   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
861ce94432eSBarry Smith   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp().");
8622fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
8634b9ad928SBarry Smith   PetscFunctionReturn(0);
8644b9ad928SBarry Smith }
8654b9ad928SBarry Smith 
8661e6b0712SBarry Smith static PetscErrorCode  PCASMSetType_ASM(PC pc,PCASMType type)
8674b9ad928SBarry Smith {
86892bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
8694b9ad928SBarry Smith 
8704b9ad928SBarry Smith   PetscFunctionBegin;
8714b9ad928SBarry Smith   osm->type     = type;
872bf108f30SBarry Smith   osm->type_set = PETSC_TRUE;
8734b9ad928SBarry Smith   PetscFunctionReturn(0);
8744b9ad928SBarry Smith }
8754b9ad928SBarry Smith 
876c60c7ad4SBarry Smith static PetscErrorCode  PCASMGetType_ASM(PC pc,PCASMType *type)
877c60c7ad4SBarry Smith {
878c60c7ad4SBarry Smith   PC_ASM *osm = (PC_ASM*)pc->data;
879c60c7ad4SBarry Smith 
880c60c7ad4SBarry Smith   PetscFunctionBegin;
881c60c7ad4SBarry Smith   *type = osm->type;
882c60c7ad4SBarry Smith   PetscFunctionReturn(0);
883c60c7ad4SBarry Smith }
884c60c7ad4SBarry Smith 
88512cd4985SMatthew G. Knepley static PetscErrorCode  PCASMSetLocalType_ASM(PC pc, PCCompositeType type)
88612cd4985SMatthew G. Knepley {
88712cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
88812cd4985SMatthew G. Knepley 
88912cd4985SMatthew G. Knepley   PetscFunctionBegin;
890d0ecd4dfSBarry Smith   if (type != PC_COMPOSITE_ADDITIVE && type != PC_COMPOSITE_MULTIPLICATIVE) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Only supports additive or multiplicative as the local type");
89112cd4985SMatthew G. Knepley   osm->loctype = type;
89212cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
89312cd4985SMatthew G. Knepley }
89412cd4985SMatthew G. Knepley 
89512cd4985SMatthew G. Knepley static PetscErrorCode  PCASMGetLocalType_ASM(PC pc, PCCompositeType *type)
89612cd4985SMatthew G. Knepley {
89712cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
89812cd4985SMatthew G. Knepley 
89912cd4985SMatthew G. Knepley   PetscFunctionBegin;
90012cd4985SMatthew G. Knepley   *type = osm->loctype;
90112cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
90212cd4985SMatthew G. Knepley }
90312cd4985SMatthew G. Knepley 
9041e6b0712SBarry Smith static PetscErrorCode  PCASMSetSortIndices_ASM(PC pc,PetscBool  doSort)
9056ed231c7SMatthew Knepley {
9066ed231c7SMatthew Knepley   PC_ASM *osm = (PC_ASM*)pc->data;
9076ed231c7SMatthew Knepley 
9086ed231c7SMatthew Knepley   PetscFunctionBegin;
9096ed231c7SMatthew Knepley   osm->sort_indices = doSort;
9106ed231c7SMatthew Knepley   PetscFunctionReturn(0);
9116ed231c7SMatthew Knepley }
9126ed231c7SMatthew Knepley 
9131e6b0712SBarry Smith static PetscErrorCode  PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
9144b9ad928SBarry Smith {
91592bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
916dfbe8321SBarry Smith   PetscErrorCode ierr;
9174b9ad928SBarry Smith 
9184b9ad928SBarry Smith   PetscFunctionBegin;
919ce94432eSBarry 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");
9204b9ad928SBarry Smith 
9212fa5cd67SKarl Rupp   if (n_local) *n_local = osm->n_local_true;
92292bb6962SLisandro Dalcin   if (first_local) {
923ce94432eSBarry Smith     ierr          = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
92492bb6962SLisandro Dalcin     *first_local -= osm->n_local_true;
92592bb6962SLisandro Dalcin   }
92692bb6962SLisandro Dalcin   if (ksp) {
92792bb6962SLisandro Dalcin     /* Assume that local solves are now different; not necessarily
92892bb6962SLisandro Dalcin        true though!  This flag is used only for PCView_ASM() */
92992bb6962SLisandro Dalcin     *ksp                   = osm->ksp;
93092bb6962SLisandro Dalcin     osm->same_local_solves = PETSC_FALSE;
93192bb6962SLisandro Dalcin   }
9324b9ad928SBarry Smith   PetscFunctionReturn(0);
9334b9ad928SBarry Smith }
9344b9ad928SBarry Smith 
93580ec0b7dSPatrick Sanan static PetscErrorCode  PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type)
93680ec0b7dSPatrick Sanan {
93780ec0b7dSPatrick Sanan   PC_ASM         *osm = (PC_ASM*)pc->data;
93880ec0b7dSPatrick Sanan 
93980ec0b7dSPatrick Sanan   PetscFunctionBegin;
94080ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
94180ec0b7dSPatrick Sanan   PetscValidPointer(sub_mat_type,2);
94280ec0b7dSPatrick Sanan   *sub_mat_type = osm->sub_mat_type;
94380ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
94480ec0b7dSPatrick Sanan }
94580ec0b7dSPatrick Sanan 
94680ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type)
94780ec0b7dSPatrick Sanan {
94880ec0b7dSPatrick Sanan   PetscErrorCode    ierr;
94980ec0b7dSPatrick Sanan   PC_ASM            *osm = (PC_ASM*)pc->data;
95080ec0b7dSPatrick Sanan 
95180ec0b7dSPatrick Sanan   PetscFunctionBegin;
95280ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
95380ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
95480ec0b7dSPatrick Sanan   ierr = PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);CHKERRQ(ierr);
95580ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
95680ec0b7dSPatrick Sanan }
95780ec0b7dSPatrick Sanan 
9584b9ad928SBarry Smith /*@C
9591093a601SBarry Smith     PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner.
9604b9ad928SBarry Smith 
9614b9ad928SBarry Smith     Collective on PC
9624b9ad928SBarry Smith 
9634b9ad928SBarry Smith     Input Parameters:
9644b9ad928SBarry Smith +   pc - the preconditioner context
9654b9ad928SBarry Smith .   n - the number of subdomains for this processor (default value = 1)
9668c03b21aSDmitry Karpeev .   is - the index set that defines the subdomains for this processor
9670298fd71SBarry Smith          (or NULL for PETSc to determine subdomains)
968f1ee410cSBarry 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
969f1ee410cSBarry Smith          (or NULL to not provide these)
9704b9ad928SBarry Smith 
9714b9ad928SBarry Smith     Notes:
9721093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
9734b9ad928SBarry Smith 
9744b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
9754b9ad928SBarry Smith 
9764b9ad928SBarry Smith     Use PCASMSetTotalSubdomains() to set the subdomains for all processors.
9774b9ad928SBarry Smith 
978f1ee410cSBarry Smith     If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated
979f1ee410cSBarry Smith     back to form the global solution (this is the standard restricted additive Schwarz method)
980f1ee410cSBarry Smith 
981f1ee410cSBarry 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
982f1ee410cSBarry Smith     no code to handle that case.
983f1ee410cSBarry Smith 
9844b9ad928SBarry Smith     Level: advanced
9854b9ad928SBarry Smith 
9864b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz
9874b9ad928SBarry Smith 
9884b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
989f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType()
9904b9ad928SBarry Smith @*/
9917087cfbeSBarry Smith PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
9924b9ad928SBarry Smith {
9937bb14e67SBarry Smith   PetscErrorCode ierr;
9944b9ad928SBarry Smith 
9954b9ad928SBarry Smith   PetscFunctionBegin;
9960700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9977bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr);
9984b9ad928SBarry Smith   PetscFunctionReturn(0);
9994b9ad928SBarry Smith }
10004b9ad928SBarry Smith 
10014b9ad928SBarry Smith /*@C
1002feb221f8SDmitry Karpeev     PCASMSetTotalSubdomains - Sets the subdomains for all processors for the
10034b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
10044b9ad928SBarry Smith     PC communicator must call this routine, with the same index sets.
10054b9ad928SBarry Smith 
10064b9ad928SBarry Smith     Collective on PC
10074b9ad928SBarry Smith 
10084b9ad928SBarry Smith     Input Parameters:
10094b9ad928SBarry Smith +   pc - the preconditioner context
1010feb221f8SDmitry Karpeev .   N  - the number of subdomains for all processors
1011feb221f8SDmitry Karpeev .   is - the index sets that define the subdomains for all processors
1012dfaa0c5fSBarry Smith          (or NULL to ask PETSc to determine the subdomains)
10132b691e39SMatthew Knepley -   is_local - the index sets that define the local part of the subdomains for this processor
1014f1ee410cSBarry Smith          (or NULL to not provide this information)
10154b9ad928SBarry Smith 
10164b9ad928SBarry Smith     Options Database Key:
10174b9ad928SBarry Smith     To set the total number of subdomain blocks rather than specify the
10184b9ad928SBarry Smith     index sets, use the option
10194b9ad928SBarry Smith .    -pc_asm_blocks <blks> - Sets total blocks
10204b9ad928SBarry Smith 
10214b9ad928SBarry Smith     Notes:
1022f1ee410cSBarry Smith     Currently you cannot use this to set the actual subdomains with the argument is or is_local.
10234b9ad928SBarry Smith 
10244b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
10254b9ad928SBarry Smith 
10264b9ad928SBarry Smith     These index sets cannot be destroyed until after completion of the
10274b9ad928SBarry Smith     linear solves for which the ASM preconditioner is being used.
10284b9ad928SBarry Smith 
10294b9ad928SBarry Smith     Use PCASMSetLocalSubdomains() to set local subdomains.
10304b9ad928SBarry Smith 
10311093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
10321093a601SBarry Smith 
10334b9ad928SBarry Smith     Level: advanced
10344b9ad928SBarry Smith 
10354b9ad928SBarry Smith .keywords: PC, ASM, set, total, global, subdomains, additive Schwarz
10364b9ad928SBarry Smith 
10374b9ad928SBarry Smith .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
10384b9ad928SBarry Smith           PCASMCreateSubdomains2D()
10394b9ad928SBarry Smith @*/
10407087cfbeSBarry Smith PetscErrorCode  PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[])
10414b9ad928SBarry Smith {
10427bb14e67SBarry Smith   PetscErrorCode ierr;
10434b9ad928SBarry Smith 
10444b9ad928SBarry Smith   PetscFunctionBegin;
10450700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10467bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr);
10474b9ad928SBarry Smith   PetscFunctionReturn(0);
10484b9ad928SBarry Smith }
10494b9ad928SBarry Smith 
10504b9ad928SBarry Smith /*@
10514b9ad928SBarry Smith     PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
10524b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
1053f1ee410cSBarry Smith     PC communicator must call this routine.
10544b9ad928SBarry Smith 
1055ad4df100SBarry Smith     Logically Collective on PC
10564b9ad928SBarry Smith 
10574b9ad928SBarry Smith     Input Parameters:
10584b9ad928SBarry Smith +   pc  - the preconditioner context
10594b9ad928SBarry Smith -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
10604b9ad928SBarry Smith 
10614b9ad928SBarry Smith     Options Database Key:
10624b9ad928SBarry Smith .   -pc_asm_overlap <ovl> - Sets overlap
10634b9ad928SBarry Smith 
10644b9ad928SBarry Smith     Notes:
10654b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.  To use
10664b9ad928SBarry Smith     multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
10674b9ad928SBarry Smith     PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).
10684b9ad928SBarry Smith 
10694b9ad928SBarry Smith     The overlap defaults to 1, so if one desires that no additional
10704b9ad928SBarry Smith     overlap be computed beyond what may have been set with a call to
10714b9ad928SBarry Smith     PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
10724b9ad928SBarry Smith     must be set to be 0.  In particular, if one does not explicitly set
10734b9ad928SBarry Smith     the subdomains an application code, then all overlap would be computed
10744b9ad928SBarry Smith     internally by PETSc, and using an overlap of 0 would result in an ASM
10754b9ad928SBarry Smith     variant that is equivalent to the block Jacobi preconditioner.
10764b9ad928SBarry Smith 
1077f1ee410cSBarry Smith     The default algorithm used by PETSc to increase overlap is fast, but not scalable,
1078f1ee410cSBarry Smith     use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
1079f1ee410cSBarry Smith 
10804b9ad928SBarry Smith     Note that one can define initial index sets with any overlap via
1081f1ee410cSBarry Smith     PCASMSetLocalSubdomains(); the routine
10824b9ad928SBarry Smith     PCASMSetOverlap() merely allows PETSc to extend that overlap further
10834b9ad928SBarry Smith     if desired.
10844b9ad928SBarry Smith 
10854b9ad928SBarry Smith     Level: intermediate
10864b9ad928SBarry Smith 
10874b9ad928SBarry Smith .keywords: PC, ASM, set, overlap
10884b9ad928SBarry Smith 
10894b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
1090f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap()
10914b9ad928SBarry Smith @*/
10927087cfbeSBarry Smith PetscErrorCode  PCASMSetOverlap(PC pc,PetscInt ovl)
10934b9ad928SBarry Smith {
10947bb14e67SBarry Smith   PetscErrorCode ierr;
10954b9ad928SBarry Smith 
10964b9ad928SBarry Smith   PetscFunctionBegin;
10970700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1098c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,ovl,2);
10997bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
11004b9ad928SBarry Smith   PetscFunctionReturn(0);
11014b9ad928SBarry Smith }
11024b9ad928SBarry Smith 
11034b9ad928SBarry Smith /*@
11044b9ad928SBarry Smith     PCASMSetType - Sets the type of restriction and interpolation used
11054b9ad928SBarry Smith     for local problems in the additive Schwarz method.
11064b9ad928SBarry Smith 
1107ad4df100SBarry Smith     Logically Collective on PC
11084b9ad928SBarry Smith 
11094b9ad928SBarry Smith     Input Parameters:
11104b9ad928SBarry Smith +   pc  - the preconditioner context
11114b9ad928SBarry Smith -   type - variant of ASM, one of
11124b9ad928SBarry Smith .vb
11134b9ad928SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
11144b9ad928SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
11154b9ad928SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
11164b9ad928SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
11174b9ad928SBarry Smith .ve
11184b9ad928SBarry Smith 
11194b9ad928SBarry Smith     Options Database Key:
11204b9ad928SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
11214b9ad928SBarry Smith 
1122f1ee410cSBarry Smith     Notes: if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected
1123f1ee410cSBarry Smith     to limit the local processor interpolation
1124f1ee410cSBarry Smith 
11254b9ad928SBarry Smith     Level: intermediate
11264b9ad928SBarry Smith 
11274b9ad928SBarry Smith .keywords: PC, ASM, set, type
11284b9ad928SBarry Smith 
11294b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1130f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType()
11314b9ad928SBarry Smith @*/
11327087cfbeSBarry Smith PetscErrorCode  PCASMSetType(PC pc,PCASMType type)
11334b9ad928SBarry Smith {
11347bb14e67SBarry Smith   PetscErrorCode ierr;
11354b9ad928SBarry Smith 
11364b9ad928SBarry Smith   PetscFunctionBegin;
11370700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1138c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
11397bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr);
11404b9ad928SBarry Smith   PetscFunctionReturn(0);
11414b9ad928SBarry Smith }
11424b9ad928SBarry Smith 
1143c60c7ad4SBarry Smith /*@
1144c60c7ad4SBarry Smith     PCASMGetType - Gets the type of restriction and interpolation used
1145c60c7ad4SBarry Smith     for local problems in the additive Schwarz method.
1146c60c7ad4SBarry Smith 
1147c60c7ad4SBarry Smith     Logically Collective on PC
1148c60c7ad4SBarry Smith 
1149c60c7ad4SBarry Smith     Input Parameter:
1150c60c7ad4SBarry Smith .   pc  - the preconditioner context
1151c60c7ad4SBarry Smith 
1152c60c7ad4SBarry Smith     Output Parameter:
1153c60c7ad4SBarry Smith .   type - variant of ASM, one of
1154c60c7ad4SBarry Smith 
1155c60c7ad4SBarry Smith .vb
1156c60c7ad4SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
1157c60c7ad4SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1158c60c7ad4SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1159c60c7ad4SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
1160c60c7ad4SBarry Smith .ve
1161c60c7ad4SBarry Smith 
1162c60c7ad4SBarry Smith     Options Database Key:
1163c60c7ad4SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
1164c60c7ad4SBarry Smith 
1165c60c7ad4SBarry Smith     Level: intermediate
1166c60c7ad4SBarry Smith 
1167c60c7ad4SBarry Smith .keywords: PC, ASM, set, type
1168c60c7ad4SBarry Smith 
1169c60c7ad4SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1170f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType()
1171c60c7ad4SBarry Smith @*/
1172c60c7ad4SBarry Smith PetscErrorCode  PCASMGetType(PC pc,PCASMType *type)
1173c60c7ad4SBarry Smith {
1174c60c7ad4SBarry Smith   PetscErrorCode ierr;
1175c60c7ad4SBarry Smith 
1176c60c7ad4SBarry Smith   PetscFunctionBegin;
1177c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1178c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));CHKERRQ(ierr);
1179c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1180c60c7ad4SBarry Smith }
1181c60c7ad4SBarry Smith 
118212cd4985SMatthew G. Knepley /*@
118312cd4985SMatthew G. Knepley   PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method.
118412cd4985SMatthew G. Knepley 
118512cd4985SMatthew G. Knepley   Logically Collective on PC
118612cd4985SMatthew G. Knepley 
118712cd4985SMatthew G. Knepley   Input Parameters:
118812cd4985SMatthew G. Knepley + pc  - the preconditioner context
118912cd4985SMatthew G. Knepley - type - type of composition, one of
119012cd4985SMatthew G. Knepley .vb
119112cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
119212cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
119312cd4985SMatthew G. Knepley .ve
119412cd4985SMatthew G. Knepley 
119512cd4985SMatthew G. Knepley   Options Database Key:
119612cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
119712cd4985SMatthew G. Knepley 
119812cd4985SMatthew G. Knepley   Level: intermediate
119912cd4985SMatthew G. Knepley 
1200f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
120112cd4985SMatthew G. Knepley @*/
120212cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type)
120312cd4985SMatthew G. Knepley {
120412cd4985SMatthew G. Knepley   PetscErrorCode ierr;
120512cd4985SMatthew G. Knepley 
120612cd4985SMatthew G. Knepley   PetscFunctionBegin;
120712cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
120812cd4985SMatthew G. Knepley   PetscValidLogicalCollectiveEnum(pc, type, 2);
120912cd4985SMatthew G. Knepley   ierr = PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));CHKERRQ(ierr);
121012cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
121112cd4985SMatthew G. Knepley }
121212cd4985SMatthew G. Knepley 
121312cd4985SMatthew G. Knepley /*@
121412cd4985SMatthew G. Knepley   PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method.
121512cd4985SMatthew G. Knepley 
121612cd4985SMatthew G. Knepley   Logically Collective on PC
121712cd4985SMatthew G. Knepley 
121812cd4985SMatthew G. Knepley   Input Parameter:
121912cd4985SMatthew G. Knepley . pc  - the preconditioner context
122012cd4985SMatthew G. Knepley 
122112cd4985SMatthew G. Knepley   Output Parameter:
122212cd4985SMatthew G. Knepley . type - type of composition, one of
122312cd4985SMatthew G. Knepley .vb
122412cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
122512cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
122612cd4985SMatthew G. Knepley .ve
122712cd4985SMatthew G. Knepley 
122812cd4985SMatthew G. Knepley   Options Database Key:
122912cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
123012cd4985SMatthew G. Knepley 
123112cd4985SMatthew G. Knepley   Level: intermediate
123212cd4985SMatthew G. Knepley 
1233f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
123412cd4985SMatthew G. Knepley @*/
123512cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type)
123612cd4985SMatthew G. Knepley {
123712cd4985SMatthew G. Knepley   PetscErrorCode ierr;
123812cd4985SMatthew G. Knepley 
123912cd4985SMatthew G. Knepley   PetscFunctionBegin;
124012cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
124112cd4985SMatthew G. Knepley   PetscValidPointer(type, 2);
124212cd4985SMatthew G. Knepley   ierr = PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));CHKERRQ(ierr);
124312cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
124412cd4985SMatthew G. Knepley }
124512cd4985SMatthew G. Knepley 
12466ed231c7SMatthew Knepley /*@
12476ed231c7SMatthew Knepley     PCASMSetSortIndices - Determines whether subdomain indices are sorted.
12486ed231c7SMatthew Knepley 
1249ad4df100SBarry Smith     Logically Collective on PC
12506ed231c7SMatthew Knepley 
12516ed231c7SMatthew Knepley     Input Parameters:
12526ed231c7SMatthew Knepley +   pc  - the preconditioner context
12536ed231c7SMatthew Knepley -   doSort - sort the subdomain indices
12546ed231c7SMatthew Knepley 
12556ed231c7SMatthew Knepley     Level: intermediate
12566ed231c7SMatthew Knepley 
12576ed231c7SMatthew Knepley .keywords: PC, ASM, set, type
12586ed231c7SMatthew Knepley 
12596ed231c7SMatthew Knepley .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
12606ed231c7SMatthew Knepley           PCASMCreateSubdomains2D()
12616ed231c7SMatthew Knepley @*/
12627087cfbeSBarry Smith PetscErrorCode  PCASMSetSortIndices(PC pc,PetscBool doSort)
12636ed231c7SMatthew Knepley {
12647bb14e67SBarry Smith   PetscErrorCode ierr;
12656ed231c7SMatthew Knepley 
12666ed231c7SMatthew Knepley   PetscFunctionBegin;
12670700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1268acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
12697bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
12706ed231c7SMatthew Knepley   PetscFunctionReturn(0);
12716ed231c7SMatthew Knepley }
12726ed231c7SMatthew Knepley 
12734b9ad928SBarry Smith /*@C
12744b9ad928SBarry Smith    PCASMGetSubKSP - Gets the local KSP contexts for all blocks on
12754b9ad928SBarry Smith    this processor.
12764b9ad928SBarry Smith 
12774b9ad928SBarry Smith    Collective on PC iff first_local is requested
12784b9ad928SBarry Smith 
12794b9ad928SBarry Smith    Input Parameter:
12804b9ad928SBarry Smith .  pc - the preconditioner context
12814b9ad928SBarry Smith 
12824b9ad928SBarry Smith    Output Parameters:
12830298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
12840298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
12850298fd71SBarry Smith                  all processors must request or all must pass NULL
12864b9ad928SBarry Smith -  ksp - the array of KSP contexts
12874b9ad928SBarry Smith 
12884b9ad928SBarry Smith    Note:
1289d29017ddSJed Brown    After PCASMGetSubKSP() the array of KSPes is not to be freed.
12904b9ad928SBarry Smith 
12914b9ad928SBarry Smith    You must call KSPSetUp() before calling PCASMGetSubKSP().
12924b9ad928SBarry Smith 
1293d29017ddSJed Brown    Fortran note:
12942bf68e3eSBarry 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.
1295d29017ddSJed Brown 
12964b9ad928SBarry Smith    Level: advanced
12974b9ad928SBarry Smith 
12984b9ad928SBarry Smith .keywords: PC, ASM, additive Schwarz, get, sub, KSP, context
12994b9ad928SBarry Smith 
13004b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
13014b9ad928SBarry Smith           PCASMCreateSubdomains2D(),
13024b9ad928SBarry Smith @*/
13037087cfbeSBarry Smith PetscErrorCode  PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
13044b9ad928SBarry Smith {
13057bb14e67SBarry Smith   PetscErrorCode ierr;
13064b9ad928SBarry Smith 
13074b9ad928SBarry Smith   PetscFunctionBegin;
13080700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13097bb14e67SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
13104b9ad928SBarry Smith   PetscFunctionReturn(0);
13114b9ad928SBarry Smith }
13124b9ad928SBarry Smith 
13134b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
13144b9ad928SBarry Smith /*MC
13153b09bd56SBarry Smith    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
13164b9ad928SBarry Smith            its own KSP object.
13174b9ad928SBarry Smith 
13184b9ad928SBarry Smith    Options Database Keys:
131949517cdeSBarry Smith +  -pc_asm_blocks <blks> - Sets total blocks
13204b9ad928SBarry Smith .  -pc_asm_overlap <ovl> - Sets overlap
1321f1ee410cSBarry Smith .  -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict
1322f1ee410cSBarry Smith -  -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive
13234b9ad928SBarry Smith 
13243b09bd56SBarry Smith      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
13253b09bd56SBarry Smith       will get a different convergence rate due to the default option of -pc_asm_type restrict. Use
13263b09bd56SBarry Smith       -pc_asm_type basic to use the standard ASM.
13273b09bd56SBarry Smith 
13284b9ad928SBarry Smith    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
1329f1ee410cSBarry Smith      than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor.
13304b9ad928SBarry Smith 
13313b09bd56SBarry Smith      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1332d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
13334b9ad928SBarry Smith 
1334a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCASMGetSubKSP()
13354b9ad928SBarry Smith          and set the options directly on the resulting KSP object (you can access its PC
13364b9ad928SBarry Smith          with KSPGetPC())
13374b9ad928SBarry Smith 
13384b9ad928SBarry Smith    Level: beginner
13394b9ad928SBarry Smith 
13404b9ad928SBarry Smith    Concepts: additive Schwarz method
13414b9ad928SBarry Smith 
1342c582cd25SBarry Smith     References:
134396a0c994SBarry Smith +   1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
134496a0c994SBarry Smith      Courant Institute, New York University Technical report
134596a0c994SBarry Smith -   1. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
134696a0c994SBarry Smith     Cambridge University Press.
1347c582cd25SBarry Smith 
13484b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1349f1ee410cSBarry Smith            PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType()
1350f1ee410cSBarry Smith            PCASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType
1351e09e08ccSBarry Smith 
13524b9ad928SBarry Smith M*/
13534b9ad928SBarry Smith 
13548cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc)
13554b9ad928SBarry Smith {
1356dfbe8321SBarry Smith   PetscErrorCode ierr;
13574b9ad928SBarry Smith   PC_ASM         *osm;
13584b9ad928SBarry Smith 
13594b9ad928SBarry Smith   PetscFunctionBegin;
1360b00a9115SJed Brown   ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr);
13612fa5cd67SKarl Rupp 
13624b9ad928SBarry Smith   osm->n                 = PETSC_DECIDE;
13634b9ad928SBarry Smith   osm->n_local           = 0;
13642b837212SDmitry Karpeev   osm->n_local_true      = PETSC_DECIDE;
13654b9ad928SBarry Smith   osm->overlap           = 1;
13664b9ad928SBarry Smith   osm->ksp               = 0;
13672b691e39SMatthew Knepley   osm->restriction       = 0;
1368*347574c9Seaulisa   osm->lrestriction      = 0;
1369b9845e0eSMatthew Knepley   osm->localization      = 0;
13702b691e39SMatthew Knepley   osm->prolongation      = 0;
13715b3c0d42Seaulisa   osm->lprolongation     = 0;
137292bb6962SLisandro Dalcin   osm->x                 = 0;
137392bb6962SLisandro Dalcin   osm->y                 = 0;
137453888de8SMatthew Knepley   osm->y_local           = 0;
13754b9ad928SBarry Smith   osm->is                = 0;
13762b691e39SMatthew Knepley   osm->is_local          = 0;
13774b9ad928SBarry Smith   osm->mat               = 0;
13784b9ad928SBarry Smith   osm->pmat              = 0;
13794b9ad928SBarry Smith   osm->type              = PC_ASM_RESTRICT;
138012cd4985SMatthew G. Knepley   osm->loctype           = PC_COMPOSITE_ADDITIVE;
13814b9ad928SBarry Smith   osm->same_local_solves = PETSC_TRUE;
13826ed231c7SMatthew Knepley   osm->sort_indices      = PETSC_TRUE;
1383d709ea83SDmitry Karpeev   osm->dm_subdomains     = PETSC_FALSE;
138480ec0b7dSPatrick Sanan   osm->sub_mat_type      = NULL;
13854b9ad928SBarry Smith 
138692bb6962SLisandro Dalcin   pc->data                 = (void*)osm;
13874b9ad928SBarry Smith   pc->ops->apply           = PCApply_ASM;
13884b9ad928SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_ASM;
13894b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_ASM;
1390e91c6855SBarry Smith   pc->ops->reset           = PCReset_ASM;
13914b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_ASM;
13924b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_ASM;
13934b9ad928SBarry Smith   pc->ops->setuponblocks   = PCSetUpOnBlocks_ASM;
13944b9ad928SBarry Smith   pc->ops->view            = PCView_ASM;
13954b9ad928SBarry Smith   pc->ops->applyrichardson = 0;
13964b9ad928SBarry Smith 
1397bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr);
1398bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr);
1399bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);CHKERRQ(ierr);
1400bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);CHKERRQ(ierr);
1401c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);CHKERRQ(ierr);
140212cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);CHKERRQ(ierr);
140312cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);CHKERRQ(ierr);
1404bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);CHKERRQ(ierr);
1405bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);CHKERRQ(ierr);
140680ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);CHKERRQ(ierr);
140780ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);CHKERRQ(ierr);
14084b9ad928SBarry Smith   PetscFunctionReturn(0);
14094b9ad928SBarry Smith }
14104b9ad928SBarry Smith 
141192bb6962SLisandro Dalcin /*@C
141292bb6962SLisandro Dalcin    PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
141392bb6962SLisandro Dalcin    preconditioner for a any problem on a general grid.
141492bb6962SLisandro Dalcin 
141592bb6962SLisandro Dalcin    Collective
141692bb6962SLisandro Dalcin 
141792bb6962SLisandro Dalcin    Input Parameters:
141892bb6962SLisandro Dalcin +  A - The global matrix operator
141992bb6962SLisandro Dalcin -  n - the number of local blocks
142092bb6962SLisandro Dalcin 
142192bb6962SLisandro Dalcin    Output Parameters:
142292bb6962SLisandro Dalcin .  outis - the array of index sets defining the subdomains
142392bb6962SLisandro Dalcin 
142492bb6962SLisandro Dalcin    Level: advanced
142592bb6962SLisandro Dalcin 
14267d6bfa3bSBarry Smith    Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap
14277d6bfa3bSBarry Smith     from these if you use PCASMSetLocalSubdomains()
14287d6bfa3bSBarry Smith 
14297d6bfa3bSBarry Smith     In the Fortran version you must provide the array outis[] already allocated of length n.
14307d6bfa3bSBarry Smith 
143192bb6962SLisandro Dalcin .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid
143292bb6962SLisandro Dalcin 
143392bb6962SLisandro Dalcin .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains()
143492bb6962SLisandro Dalcin @*/
14357087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
143692bb6962SLisandro Dalcin {
143792bb6962SLisandro Dalcin   MatPartitioning mpart;
143892bb6962SLisandro Dalcin   const char      *prefix;
1439c56a70eeSBarry Smith   void            (*f)(void);
144092bb6962SLisandro Dalcin   PetscInt        i,j,rstart,rend,bs;
144111bd1e4dSLisandro Dalcin   PetscBool       isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
14420298fd71SBarry Smith   Mat             Ad     = NULL, adj;
144392bb6962SLisandro Dalcin   IS              ispart,isnumb,*is;
144492bb6962SLisandro Dalcin   PetscErrorCode  ierr;
144592bb6962SLisandro Dalcin 
144692bb6962SLisandro Dalcin   PetscFunctionBegin;
14470700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
144892bb6962SLisandro Dalcin   PetscValidPointer(outis,3);
1449c1235816SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
145092bb6962SLisandro Dalcin 
145192bb6962SLisandro Dalcin   /* Get prefix, row distribution, and block size */
145292bb6962SLisandro Dalcin   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
145392bb6962SLisandro Dalcin   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
145492bb6962SLisandro Dalcin   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
145565e19b50SBarry 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);
145665e19b50SBarry Smith 
145792bb6962SLisandro Dalcin   /* Get diagonal block from matrix if possible */
1458c56a70eeSBarry Smith   ierr = MatShellGetOperation(A,MATOP_GET_DIAGONAL_BLOCK,&f);CHKERRQ(ierr);
145992bb6962SLisandro Dalcin   if (f) {
146011bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
146192bb6962SLisandro Dalcin   }
146292bb6962SLisandro Dalcin   if (Ad) {
1463251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1464251f4c67SDmitry Karpeev     if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
146592bb6962SLisandro Dalcin   }
146692bb6962SLisandro Dalcin   if (Ad && n > 1) {
1467ace3abfcSBarry Smith     PetscBool match,done;
146892bb6962SLisandro Dalcin     /* Try to setup a good matrix partitioning if available */
146992bb6962SLisandro Dalcin     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
147092bb6962SLisandro Dalcin     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
147192bb6962SLisandro Dalcin     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1472251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
147392bb6962SLisandro Dalcin     if (!match) {
1474251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
147592bb6962SLisandro Dalcin     }
147692bb6962SLisandro Dalcin     if (!match) { /* assume a "good" partitioner is available */
14771a83f524SJed Brown       PetscInt       na;
14781a83f524SJed Brown       const PetscInt *ia,*ja;
147992bb6962SLisandro Dalcin       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
148092bb6962SLisandro Dalcin       if (done) {
148192bb6962SLisandro Dalcin         /* Build adjacency matrix by hand. Unfortunately a call to
148292bb6962SLisandro Dalcin            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
148392bb6962SLisandro Dalcin            remove the block-aij structure and we cannot expect
148492bb6962SLisandro Dalcin            MatPartitioning to split vertices as we need */
14851a83f524SJed Brown         PetscInt       i,j,len,nnz,cnt,*iia=0,*jja=0;
14861a83f524SJed Brown         const PetscInt *row;
148792bb6962SLisandro Dalcin         nnz = 0;
148892bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* count number of nonzeros */
148992bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
149092bb6962SLisandro Dalcin           row = ja + ia[i];
149192bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
149292bb6962SLisandro Dalcin             if (row[j] == i) { /* don't count diagonal */
149392bb6962SLisandro Dalcin               len--; break;
149492bb6962SLisandro Dalcin             }
149592bb6962SLisandro Dalcin           }
149692bb6962SLisandro Dalcin           nnz += len;
149792bb6962SLisandro Dalcin         }
1498854ce69bSBarry Smith         ierr   = PetscMalloc1(na+1,&iia);CHKERRQ(ierr);
1499854ce69bSBarry Smith         ierr   = PetscMalloc1(nnz,&jja);CHKERRQ(ierr);
150092bb6962SLisandro Dalcin         nnz    = 0;
150192bb6962SLisandro Dalcin         iia[0] = 0;
150292bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* fill adjacency */
150392bb6962SLisandro Dalcin           cnt = 0;
150492bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
150592bb6962SLisandro Dalcin           row = ja + ia[i];
150692bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
150792bb6962SLisandro Dalcin             if (row[j] != i) { /* if not diagonal */
150892bb6962SLisandro Dalcin               jja[nnz+cnt++] = row[j];
150992bb6962SLisandro Dalcin             }
151092bb6962SLisandro Dalcin           }
151192bb6962SLisandro Dalcin           nnz     += cnt;
151292bb6962SLisandro Dalcin           iia[i+1] = nnz;
151392bb6962SLisandro Dalcin         }
151492bb6962SLisandro Dalcin         /* Partitioning of the adjacency matrix */
15150298fd71SBarry Smith         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr);
151692bb6962SLisandro Dalcin         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
151792bb6962SLisandro Dalcin         ierr      = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
151892bb6962SLisandro Dalcin         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
151992bb6962SLisandro Dalcin         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
1520fcfd50ebSBarry Smith         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
152192bb6962SLisandro Dalcin         foundpart = PETSC_TRUE;
152292bb6962SLisandro Dalcin       }
152392bb6962SLisandro Dalcin       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
152492bb6962SLisandro Dalcin     }
1525fcfd50ebSBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
152692bb6962SLisandro Dalcin   }
152792bb6962SLisandro Dalcin 
1528785e854fSJed Brown   ierr   = PetscMalloc1(n,&is);CHKERRQ(ierr);
152992bb6962SLisandro Dalcin   *outis = is;
153092bb6962SLisandro Dalcin 
153192bb6962SLisandro Dalcin   if (!foundpart) {
153292bb6962SLisandro Dalcin 
153392bb6962SLisandro Dalcin     /* Partitioning by contiguous chunks of rows */
153492bb6962SLisandro Dalcin 
153592bb6962SLisandro Dalcin     PetscInt mbs   = (rend-rstart)/bs;
153692bb6962SLisandro Dalcin     PetscInt start = rstart;
153792bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
153892bb6962SLisandro Dalcin       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
153992bb6962SLisandro Dalcin       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
154092bb6962SLisandro Dalcin       start += count;
154192bb6962SLisandro Dalcin     }
154292bb6962SLisandro Dalcin 
154392bb6962SLisandro Dalcin   } else {
154492bb6962SLisandro Dalcin 
154592bb6962SLisandro Dalcin     /* Partitioning by adjacency of diagonal block  */
154692bb6962SLisandro Dalcin 
154792bb6962SLisandro Dalcin     const PetscInt *numbering;
154892bb6962SLisandro Dalcin     PetscInt       *count,nidx,*indices,*newidx,start=0;
154992bb6962SLisandro Dalcin     /* Get node count in each partition */
1550785e854fSJed Brown     ierr = PetscMalloc1(n,&count);CHKERRQ(ierr);
155192bb6962SLisandro Dalcin     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
155292bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
155392bb6962SLisandro Dalcin       for (i=0; i<n; i++) count[i] *= bs;
155492bb6962SLisandro Dalcin     }
155592bb6962SLisandro Dalcin     /* Build indices from node numbering */
155692bb6962SLisandro Dalcin     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1557785e854fSJed Brown     ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr);
155892bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
155992bb6962SLisandro Dalcin     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
156092bb6962SLisandro Dalcin     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
156192bb6962SLisandro Dalcin     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
156292bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1563785e854fSJed Brown       ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr);
15642fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
15652fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
15662fa5cd67SKarl Rupp       }
156792bb6962SLisandro Dalcin       ierr    = PetscFree(indices);CHKERRQ(ierr);
156892bb6962SLisandro Dalcin       nidx   *= bs;
156992bb6962SLisandro Dalcin       indices = newidx;
157092bb6962SLisandro Dalcin     }
157192bb6962SLisandro Dalcin     /* Shift to get global indices */
157292bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] += rstart;
157392bb6962SLisandro Dalcin 
157492bb6962SLisandro Dalcin     /* Build the index sets for each block */
157592bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
157670b3c8c7SBarry Smith       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
157792bb6962SLisandro Dalcin       ierr   = ISSort(is[i]);CHKERRQ(ierr);
157892bb6962SLisandro Dalcin       start += count[i];
157992bb6962SLisandro Dalcin     }
158092bb6962SLisandro Dalcin 
15813bf036e2SBarry Smith     ierr = PetscFree(count);CHKERRQ(ierr);
15823bf036e2SBarry Smith     ierr = PetscFree(indices);CHKERRQ(ierr);
1583fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1584fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
158592bb6962SLisandro Dalcin 
158692bb6962SLisandro Dalcin   }
158792bb6962SLisandro Dalcin   PetscFunctionReturn(0);
158892bb6962SLisandro Dalcin }
158992bb6962SLisandro Dalcin 
159092bb6962SLisandro Dalcin /*@C
159192bb6962SLisandro Dalcin    PCASMDestroySubdomains - Destroys the index sets created with
159292bb6962SLisandro Dalcin    PCASMCreateSubdomains(). Should be called after setting subdomains
159392bb6962SLisandro Dalcin    with PCASMSetLocalSubdomains().
159492bb6962SLisandro Dalcin 
159592bb6962SLisandro Dalcin    Collective
159692bb6962SLisandro Dalcin 
159792bb6962SLisandro Dalcin    Input Parameters:
159892bb6962SLisandro Dalcin +  n - the number of index sets
15992b691e39SMatthew Knepley .  is - the array of index sets
16000298fd71SBarry Smith -  is_local - the array of local index sets, can be NULL
160192bb6962SLisandro Dalcin 
160292bb6962SLisandro Dalcin    Level: advanced
160392bb6962SLisandro Dalcin 
160492bb6962SLisandro Dalcin .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid
160592bb6962SLisandro Dalcin 
160692bb6962SLisandro Dalcin .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains()
160792bb6962SLisandro Dalcin @*/
16087087cfbeSBarry Smith PetscErrorCode  PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
160992bb6962SLisandro Dalcin {
161092bb6962SLisandro Dalcin   PetscInt       i;
161192bb6962SLisandro Dalcin   PetscErrorCode ierr;
16125fd66863SKarl Rupp 
161392bb6962SLisandro Dalcin   PetscFunctionBegin;
1614a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1615a35b7b57SDmitry Karpeev   if (is) {
161692bb6962SLisandro Dalcin     PetscValidPointer(is,2);
1617fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); }
161892bb6962SLisandro Dalcin     ierr = PetscFree(is);CHKERRQ(ierr);
1619a35b7b57SDmitry Karpeev   }
16202b691e39SMatthew Knepley   if (is_local) {
16212b691e39SMatthew Knepley     PetscValidPointer(is_local,3);
1622fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); }
16232b691e39SMatthew Knepley     ierr = PetscFree(is_local);CHKERRQ(ierr);
16242b691e39SMatthew Knepley   }
162592bb6962SLisandro Dalcin   PetscFunctionReturn(0);
162692bb6962SLisandro Dalcin }
162792bb6962SLisandro Dalcin 
16284b9ad928SBarry Smith /*@
16294b9ad928SBarry Smith    PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
16304b9ad928SBarry Smith    preconditioner for a two-dimensional problem on a regular grid.
16314b9ad928SBarry Smith 
16324b9ad928SBarry Smith    Not Collective
16334b9ad928SBarry Smith 
16344b9ad928SBarry Smith    Input Parameters:
16354b9ad928SBarry Smith +  m, n - the number of mesh points in the x and y directions
16364b9ad928SBarry Smith .  M, N - the number of subdomains in the x and y directions
16374b9ad928SBarry Smith .  dof - degrees of freedom per node
16384b9ad928SBarry Smith -  overlap - overlap in mesh lines
16394b9ad928SBarry Smith 
16404b9ad928SBarry Smith    Output Parameters:
16414b9ad928SBarry Smith +  Nsub - the number of subdomains created
16423d03bb48SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
16433d03bb48SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
16444b9ad928SBarry Smith 
16454b9ad928SBarry Smith    Note:
16464b9ad928SBarry Smith    Presently PCAMSCreateSubdomains2d() is valid only for sequential
16474b9ad928SBarry Smith    preconditioners.  More general related routines are
16484b9ad928SBarry Smith    PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().
16494b9ad928SBarry Smith 
16504b9ad928SBarry Smith    Level: advanced
16514b9ad928SBarry Smith 
16524b9ad928SBarry Smith .keywords: PC, ASM, additive Schwarz, create, subdomains, 2D, regular grid
16534b9ad928SBarry Smith 
16544b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
16554b9ad928SBarry Smith           PCASMSetOverlap()
16564b9ad928SBarry Smith @*/
16577087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local)
16584b9ad928SBarry Smith {
16593d03bb48SJed Brown   PetscInt       i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer;
16606849ba73SBarry Smith   PetscErrorCode ierr;
166113f74950SBarry Smith   PetscInt       nidx,*idx,loc,ii,jj,count;
16624b9ad928SBarry Smith 
16634b9ad928SBarry Smith   PetscFunctionBegin;
1664e32f2f54SBarry Smith   if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," ");
16654b9ad928SBarry Smith 
16664b9ad928SBarry Smith   *Nsub     = N*M;
1667854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is);CHKERRQ(ierr);
1668854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is_local);CHKERRQ(ierr);
16694b9ad928SBarry Smith   ystart    = 0;
16703d03bb48SJed Brown   loc_outer = 0;
16714b9ad928SBarry Smith   for (i=0; i<N; i++) {
16724b9ad928SBarry Smith     height = n/N + ((n % N) > i); /* height of subdomain */
1673e32f2f54SBarry Smith     if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n");
16744b9ad928SBarry Smith     yleft  = ystart - overlap; if (yleft < 0) yleft = 0;
16754b9ad928SBarry Smith     yright = ystart + height + overlap; if (yright > n) yright = n;
16764b9ad928SBarry Smith     xstart = 0;
16774b9ad928SBarry Smith     for (j=0; j<M; j++) {
16784b9ad928SBarry Smith       width = m/M + ((m % M) > j); /* width of subdomain */
1679e32f2f54SBarry Smith       if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m");
16804b9ad928SBarry Smith       xleft  = xstart - overlap; if (xleft < 0) xleft = 0;
16814b9ad928SBarry Smith       xright = xstart + width + overlap; if (xright > m) xright = m;
16824b9ad928SBarry Smith       nidx   = (xright - xleft)*(yright - yleft);
1683785e854fSJed Brown       ierr   = PetscMalloc1(nidx,&idx);CHKERRQ(ierr);
16844b9ad928SBarry Smith       loc    = 0;
16854b9ad928SBarry Smith       for (ii=yleft; ii<yright; ii++) {
16864b9ad928SBarry Smith         count = m*ii + xleft;
16872fa5cd67SKarl Rupp         for (jj=xleft; jj<xright; jj++) idx[loc++] = count++;
16884b9ad928SBarry Smith       }
168970b3c8c7SBarry Smith       ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr);
16903d03bb48SJed Brown       if (overlap == 0) {
16913d03bb48SJed Brown         ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr);
16922fa5cd67SKarl Rupp 
16933d03bb48SJed Brown         (*is_local)[loc_outer] = (*is)[loc_outer];
16943d03bb48SJed Brown       } else {
16953d03bb48SJed Brown         for (loc=0,ii=ystart; ii<ystart+height; ii++) {
16963d03bb48SJed Brown           for (jj=xstart; jj<xstart+width; jj++) {
16973d03bb48SJed Brown             idx[loc++] = m*ii + jj;
16983d03bb48SJed Brown           }
16993d03bb48SJed Brown         }
170070b3c8c7SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr);
17013d03bb48SJed Brown       }
17024b9ad928SBarry Smith       ierr    = PetscFree(idx);CHKERRQ(ierr);
17034b9ad928SBarry Smith       xstart += width;
17043d03bb48SJed Brown       loc_outer++;
17054b9ad928SBarry Smith     }
17064b9ad928SBarry Smith     ystart += height;
17074b9ad928SBarry Smith   }
17084b9ad928SBarry Smith   for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); }
17094b9ad928SBarry Smith   PetscFunctionReturn(0);
17104b9ad928SBarry Smith }
17114b9ad928SBarry Smith 
17124b9ad928SBarry Smith /*@C
17134b9ad928SBarry Smith     PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
17144b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
17154b9ad928SBarry Smith 
1716ad4df100SBarry Smith     Not Collective
17174b9ad928SBarry Smith 
17184b9ad928SBarry Smith     Input Parameter:
17194b9ad928SBarry Smith .   pc - the preconditioner context
17204b9ad928SBarry Smith 
17214b9ad928SBarry Smith     Output Parameters:
17224b9ad928SBarry Smith +   n - the number of subdomains for this processor (default value = 1)
17232b691e39SMatthew Knepley .   is - the index sets that define the subdomains for this processor
17240298fd71SBarry Smith -   is_local - the index sets that define the local part of the subdomains for this processor (can be NULL)
17254b9ad928SBarry Smith 
17264b9ad928SBarry Smith 
17274b9ad928SBarry Smith     Notes:
17284b9ad928SBarry Smith     The IS numbering is in the parallel, global numbering of the vector.
17294b9ad928SBarry Smith 
17304b9ad928SBarry Smith     Level: advanced
17314b9ad928SBarry Smith 
17324b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz
17334b9ad928SBarry Smith 
17344b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
17354b9ad928SBarry Smith           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices()
17364b9ad928SBarry Smith @*/
17377087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
17384b9ad928SBarry Smith {
17392a808120SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
174092bb6962SLisandro Dalcin   PetscErrorCode ierr;
1741ace3abfcSBarry Smith   PetscBool      match;
17424b9ad928SBarry Smith 
17434b9ad928SBarry Smith   PetscFunctionBegin;
17440700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
17454482741eSBarry Smith   PetscValidIntPointer(n,2);
174692bb6962SLisandro Dalcin   if (is) PetscValidPointer(is,3);
1747251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
17482a808120SBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM");
17494b9ad928SBarry Smith   if (n) *n = osm->n_local_true;
17504b9ad928SBarry Smith   if (is) *is = osm->is;
17512b691e39SMatthew Knepley   if (is_local) *is_local = osm->is_local;
17524b9ad928SBarry Smith   PetscFunctionReturn(0);
17534b9ad928SBarry Smith }
17544b9ad928SBarry Smith 
17554b9ad928SBarry Smith /*@C
17564b9ad928SBarry Smith     PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
17574b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
17584b9ad928SBarry Smith 
1759ad4df100SBarry Smith     Not Collective
17604b9ad928SBarry Smith 
17614b9ad928SBarry Smith     Input Parameter:
17624b9ad928SBarry Smith .   pc - the preconditioner context
17634b9ad928SBarry Smith 
17644b9ad928SBarry Smith     Output Parameters:
17654b9ad928SBarry Smith +   n - the number of matrices for this processor (default value = 1)
17664b9ad928SBarry Smith -   mat - the matrices
17674b9ad928SBarry Smith 
17684b9ad928SBarry Smith     Level: advanced
17694b9ad928SBarry Smith 
1770cf739d55SBarry Smith     Notes: Call after PCSetUp() (or KSPSetUp()) but before PCApply() (or KSPApply()) and before PCSetUpOnBlocks())
1771cf739d55SBarry Smith 
1772cf739d55SBarry Smith            Usually one would use PCSetModifySubmatrices() to change the submatrices in building the preconditioner.
1773cf739d55SBarry Smith 
17744b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz, block Jacobi
17754b9ad928SBarry Smith 
17764b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1777cf739d55SBarry Smith           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubmatrices()
17784b9ad928SBarry Smith @*/
17797087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
17804b9ad928SBarry Smith {
17814b9ad928SBarry Smith   PC_ASM         *osm;
178292bb6962SLisandro Dalcin   PetscErrorCode ierr;
1783ace3abfcSBarry Smith   PetscBool      match;
17844b9ad928SBarry Smith 
17854b9ad928SBarry Smith   PetscFunctionBegin;
17860700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
178792bb6962SLisandro Dalcin   PetscValidIntPointer(n,2);
178892bb6962SLisandro Dalcin   if (mat) PetscValidPointer(mat,3);
1789ce94432eSBarry Smith   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1790251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
179192bb6962SLisandro Dalcin   if (!match) {
179292bb6962SLisandro Dalcin     if (n) *n = 0;
17930298fd71SBarry Smith     if (mat) *mat = NULL;
179492bb6962SLisandro Dalcin   } else {
17954b9ad928SBarry Smith     osm = (PC_ASM*)pc->data;
17964b9ad928SBarry Smith     if (n) *n = osm->n_local_true;
17974b9ad928SBarry Smith     if (mat) *mat = osm->pmat;
179892bb6962SLisandro Dalcin   }
17994b9ad928SBarry Smith   PetscFunctionReturn(0);
18004b9ad928SBarry Smith }
1801d709ea83SDmitry Karpeev 
1802d709ea83SDmitry Karpeev /*@
1803d709ea83SDmitry Karpeev     PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1804f1ee410cSBarry Smith 
1805d709ea83SDmitry Karpeev     Logically Collective
1806d709ea83SDmitry Karpeev 
1807d709ea83SDmitry Karpeev     Input Parameter:
1808d709ea83SDmitry Karpeev +   pc  - the preconditioner
1809d709ea83SDmitry Karpeev -   flg - boolean indicating whether to use subdomains defined by the DM
1810d709ea83SDmitry Karpeev 
1811d709ea83SDmitry Karpeev     Options Database Key:
1812d709ea83SDmitry Karpeev .   -pc_asm_dm_subdomains
1813d709ea83SDmitry Karpeev 
1814d709ea83SDmitry Karpeev     Level: intermediate
1815d709ea83SDmitry Karpeev 
1816d709ea83SDmitry Karpeev     Notes:
1817d709ea83SDmitry Karpeev     PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(),
1818d709ea83SDmitry Karpeev     so setting either of the first two effectively turns the latter off.
1819d709ea83SDmitry Karpeev 
1820d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1821d709ea83SDmitry Karpeev 
1822d709ea83SDmitry Karpeev .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1823d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1824d709ea83SDmitry Karpeev @*/
1825d709ea83SDmitry Karpeev PetscErrorCode  PCASMSetDMSubdomains(PC pc,PetscBool flg)
1826d709ea83SDmitry Karpeev {
1827d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1828d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1829d709ea83SDmitry Karpeev   PetscBool      match;
1830d709ea83SDmitry Karpeev 
1831d709ea83SDmitry Karpeev   PetscFunctionBegin;
1832d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1833d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
1834d709ea83SDmitry Karpeev   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1835d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1836d709ea83SDmitry Karpeev   if (match) {
1837d709ea83SDmitry Karpeev     osm->dm_subdomains = flg;
1838d709ea83SDmitry Karpeev   }
1839d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1840d709ea83SDmitry Karpeev }
1841d709ea83SDmitry Karpeev 
1842d709ea83SDmitry Karpeev /*@
1843d709ea83SDmitry Karpeev     PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1844d709ea83SDmitry Karpeev     Not Collective
1845d709ea83SDmitry Karpeev 
1846d709ea83SDmitry Karpeev     Input Parameter:
1847d709ea83SDmitry Karpeev .   pc  - the preconditioner
1848d709ea83SDmitry Karpeev 
1849d709ea83SDmitry Karpeev     Output Parameter:
1850d709ea83SDmitry Karpeev .   flg - boolean indicating whether to use subdomains defined by the DM
1851d709ea83SDmitry Karpeev 
1852d709ea83SDmitry Karpeev     Level: intermediate
1853d709ea83SDmitry Karpeev 
1854d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1855d709ea83SDmitry Karpeev 
1856d709ea83SDmitry Karpeev .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1857d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1858d709ea83SDmitry Karpeev @*/
1859d709ea83SDmitry Karpeev PetscErrorCode  PCASMGetDMSubdomains(PC pc,PetscBool* flg)
1860d709ea83SDmitry Karpeev {
1861d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1862d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1863d709ea83SDmitry Karpeev   PetscBool      match;
1864d709ea83SDmitry Karpeev 
1865d709ea83SDmitry Karpeev   PetscFunctionBegin;
1866d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1867d709ea83SDmitry Karpeev   PetscValidPointer(flg,2);
1868d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1869d709ea83SDmitry Karpeev   if (match) {
1870d709ea83SDmitry Karpeev     if (flg) *flg = osm->dm_subdomains;
1871d709ea83SDmitry Karpeev   }
1872d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1873d709ea83SDmitry Karpeev }
187480ec0b7dSPatrick Sanan 
187580ec0b7dSPatrick Sanan /*@
187680ec0b7dSPatrick Sanan      PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string.
187780ec0b7dSPatrick Sanan 
187880ec0b7dSPatrick Sanan    Not Collective
187980ec0b7dSPatrick Sanan 
188080ec0b7dSPatrick Sanan    Input Parameter:
188180ec0b7dSPatrick Sanan .  pc - the PC
188280ec0b7dSPatrick Sanan 
188380ec0b7dSPatrick Sanan    Output Parameter:
1884f1ee410cSBarry Smith .  -pc_asm_sub_mat_type - name of matrix type
188580ec0b7dSPatrick Sanan 
188680ec0b7dSPatrick Sanan    Level: advanced
188780ec0b7dSPatrick Sanan 
188880ec0b7dSPatrick Sanan .keywords: PC, PCASM, MatType, set
188980ec0b7dSPatrick Sanan 
189080ec0b7dSPatrick Sanan .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
189180ec0b7dSPatrick Sanan @*/
189280ec0b7dSPatrick Sanan PetscErrorCode  PCASMGetSubMatType(PC pc,MatType *sub_mat_type){
189380ec0b7dSPatrick Sanan   PetscErrorCode ierr;
189480ec0b7dSPatrick Sanan 
189580ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));CHKERRQ(ierr);
189680ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
189780ec0b7dSPatrick Sanan }
189880ec0b7dSPatrick Sanan 
189980ec0b7dSPatrick Sanan /*@
190080ec0b7dSPatrick Sanan      PCASMSetSubMatType - Set the type of matrix used for ASM subsolves
190180ec0b7dSPatrick Sanan 
190280ec0b7dSPatrick Sanan    Collective on Mat
190380ec0b7dSPatrick Sanan 
190480ec0b7dSPatrick Sanan    Input Parameters:
190580ec0b7dSPatrick Sanan +  pc             - the PC object
190680ec0b7dSPatrick Sanan -  sub_mat_type   - matrix type
190780ec0b7dSPatrick Sanan 
190880ec0b7dSPatrick Sanan    Options Database Key:
190980ec0b7dSPatrick 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.
191080ec0b7dSPatrick Sanan 
191180ec0b7dSPatrick Sanan    Notes:
191280ec0b7dSPatrick Sanan    See "${PETSC_DIR}/include/petscmat.h" for available types
191380ec0b7dSPatrick Sanan 
191480ec0b7dSPatrick Sanan   Level: advanced
191580ec0b7dSPatrick Sanan 
191680ec0b7dSPatrick Sanan .keywords: PC, PCASM, MatType, set
191780ec0b7dSPatrick Sanan 
191880ec0b7dSPatrick Sanan .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
191980ec0b7dSPatrick Sanan @*/
192080ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type)
192180ec0b7dSPatrick Sanan {
192280ec0b7dSPatrick Sanan   PetscErrorCode ierr;
192380ec0b7dSPatrick Sanan 
192480ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));CHKERRQ(ierr);
192580ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
192680ec0b7dSPatrick Sanan }
1927