xref: /petsc/src/ksp/pc/impls/asm/asm.c (revision 534a8f05a7a8aff70dd8cfd53d9cd834400a8dbf)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith   This file defines an additive Schwarz preconditioner for any Mat implementation.
34b9ad928SBarry Smith 
44b9ad928SBarry Smith   Note that each processor may have any number of subdomains. But in order to
54b9ad928SBarry Smith   deal easily with the VecScatter(), we treat each processor as if it has the
64b9ad928SBarry Smith   same number of subdomains.
74b9ad928SBarry Smith 
84b9ad928SBarry Smith        n - total number of true subdomains on all processors
94b9ad928SBarry Smith        n_local_true - actual number of subdomains on this processor
104b9ad928SBarry Smith        n_local = maximum over all processors of n_local_true
114b9ad928SBarry Smith */
12af0996ceSBarry Smith #include <petsc/private/pcimpl.h>     /*I "petscpc.h" I*/
131e25c274SJed Brown #include <petscdm.h>
144b9ad928SBarry Smith 
154b9ad928SBarry Smith typedef struct {
1613f74950SBarry Smith   PetscInt   n, n_local, n_local_true;
1713f74950SBarry Smith   PetscInt   overlap;             /* overlap requested by user */
184b9ad928SBarry Smith   KSP        *ksp;                /* linear solvers for each block */
19b0de9f23SBarry Smith   VecScatter restriction;         /* mapping from global to overlapping (process) subdomain*/
20b0de9f23SBarry Smith   VecScatter *lrestriction;       /* mapping from subregion to overlapping (process) subdomain */
21b0de9f23SBarry Smith   VecScatter *lprolongation;      /* mapping from non-overlapping subregion to overlapping (process) subdomain; used for restrict additive version of algorithms */
221dd8081eSeaulisa   Vec        lx, ly;              /* work vectors */
231dd8081eSeaulisa   Vec        *x,*y;               /* work vectors */
241dd8081eSeaulisa   IS         lis;                 /* index set that defines each overlapping multiplicative (process) subdomain */
252b691e39SMatthew Knepley   IS         *is;                 /* index set that defines each overlapping subdomain */
262b691e39SMatthew Knepley   IS         *is_local;           /* index set that defines each non-overlapping subdomain, may be NULL */
274b9ad928SBarry Smith   Mat        *mat,*pmat;          /* mat is not currently used */
284b9ad928SBarry Smith   PCASMType  type;                /* use reduced interpolation, restriction or both */
29ace3abfcSBarry Smith   PetscBool  type_set;            /* if user set this value (so won't change it for symmetric problems) */
30ace3abfcSBarry Smith   PetscBool  same_local_solves;   /* flag indicating whether all local solvers are same */
31ace3abfcSBarry Smith   PetscBool  sort_indices;        /* flag to sort subdomain indices */
32d709ea83SDmitry Karpeev   PetscBool  dm_subdomains;       /* whether DM is allowed to define subdomains */
3312cd4985SMatthew G. Knepley   PCCompositeType loctype;        /* the type of composition for local solves */
3480ec0b7dSPatrick Sanan   MatType    sub_mat_type;        /* the type of Mat used for subdomain solves (can be MATSAME or NULL) */
35fb745f2cSMatthew G. Knepley   /* For multiplicative solve */
36235cc4ceSMatthew G. Knepley   Mat       *lmats;               /* submatrices for overlapping multiplicative (process) subdomain */
374b9ad928SBarry Smith } PC_ASM;
384b9ad928SBarry Smith 
396849ba73SBarry Smith static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer)
404b9ad928SBarry Smith {
4192bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
426849ba73SBarry Smith   PetscErrorCode ierr;
4313f74950SBarry Smith   PetscMPIInt    rank;
444d219a6aSLisandro Dalcin   PetscInt       i,bsz;
45ace3abfcSBarry Smith   PetscBool      iascii,isstring;
464b9ad928SBarry Smith   PetscViewer    sviewer;
474b9ad928SBarry Smith 
484b9ad928SBarry Smith   PetscFunctionBegin;
49251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
50251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
5132077d6dSBarry Smith   if (iascii) {
523d03bb48SJed Brown     char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set";
538caf3d72SBarry Smith     if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof(overlaps),"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);}
548caf3d72SBarry Smith     if (osm->n > 0) {ierr = PetscSNPrintf(blocks,sizeof(blocks),"total subdomain blocks = %D",osm->n);CHKERRQ(ierr);}
55efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  %s, %s\n",blocks,overlaps);CHKERRQ(ierr);
56efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  restriction/interpolation type - %s\n",PCASMTypes[osm->type]);CHKERRQ(ierr);
57e64f0791SPatrick Sanan     if (osm->dm_subdomains) {ierr = PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: using DM to define subdomains\n");CHKERRQ(ierr);}
5812cd4985SMatthew G. Knepley     if (osm->loctype != PC_COMPOSITE_ADDITIVE) {ierr = PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: local solve composition type - %s\n",PCCompositeTypes[osm->loctype]);CHKERRQ(ierr);}
59ce94432eSBarry Smith     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
6092bb6962SLisandro Dalcin     if (osm->same_local_solves) {
614d219a6aSLisandro Dalcin       if (osm->ksp) {
624b9ad928SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Local solve is same for all blocks, in the following KSP and PC objects:\n");CHKERRQ(ierr);
63c5e4d11fSDmitry Karpeev         ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
644d219a6aSLisandro Dalcin         if (!rank) {
654b9ad928SBarry Smith           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
6692bb6962SLisandro Dalcin           ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);
674b9ad928SBarry Smith           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
684b9ad928SBarry Smith         }
69c5e4d11fSDmitry Karpeev         ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
704d219a6aSLisandro Dalcin       }
714b9ad928SBarry Smith     } else {
72c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
734d219a6aSLisandro Dalcin       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,osm->n_local_true);CHKERRQ(ierr);
744d219a6aSLisandro Dalcin       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
754b9ad928SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr);
764b9ad928SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
774d219a6aSLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
78c5e4d11fSDmitry Karpeev       ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
79dd2fa690SBarry Smith       for (i=0; i<osm->n_local_true; i++) {
804d219a6aSLisandro Dalcin         ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr);
814d219a6aSLisandro Dalcin         ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr);
8292bb6962SLisandro Dalcin         ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr);
834d219a6aSLisandro Dalcin         ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
844b9ad928SBarry Smith       }
85c5e4d11fSDmitry Karpeev       ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
864b9ad928SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
874b9ad928SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
88c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
894b9ad928SBarry Smith     }
904b9ad928SBarry Smith   } else if (isstring) {
914d219a6aSLisandro Dalcin     ierr = PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCASMTypes[osm->type]);CHKERRQ(ierr);
92c5e4d11fSDmitry Karpeev     ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
9392bb6962SLisandro Dalcin     if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);}
94c5e4d11fSDmitry Karpeev     ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
954b9ad928SBarry Smith   }
964b9ad928SBarry Smith   PetscFunctionReturn(0);
974b9ad928SBarry Smith }
984b9ad928SBarry Smith 
9992bb6962SLisandro Dalcin static PetscErrorCode PCASMPrintSubdomains(PC pc)
10092bb6962SLisandro Dalcin {
10192bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
10292bb6962SLisandro Dalcin   const char     *prefix;
10392bb6962SLisandro Dalcin   char           fname[PETSC_MAX_PATH_LEN+1];
104643535aeSDmitry Karpeev   PetscViewer    viewer, sviewer;
105643535aeSDmitry Karpeev   char           *s;
10692bb6962SLisandro Dalcin   PetscInt       i,j,nidx;
10792bb6962SLisandro Dalcin   const PetscInt *idx;
108643535aeSDmitry Karpeev   PetscMPIInt    rank, size;
10992bb6962SLisandro Dalcin   PetscErrorCode ierr;
110846783a0SBarry Smith 
11192bb6962SLisandro Dalcin   PetscFunctionBegin;
112ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr);
113ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr);
11492bb6962SLisandro Dalcin   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
115c5929fdfSBarry Smith   ierr = PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr);
11692bb6962SLisandro Dalcin   if (fname[0] == 0) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
117ce94432eSBarry Smith   ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr);
118643535aeSDmitry Karpeev   for (i=0; i<osm->n_local; i++) {
119643535aeSDmitry Karpeev     if (i < osm->n_local_true) {
12092bb6962SLisandro Dalcin       ierr = ISGetLocalSize(osm->is[i],&nidx);CHKERRQ(ierr);
12192bb6962SLisandro Dalcin       ierr = ISGetIndices(osm->is[i],&idx);CHKERRQ(ierr);
122643535aeSDmitry Karpeev       /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
12336a9e3b9SBarry Smith #define len  16*(nidx+1)+512
12436a9e3b9SBarry Smith       ierr = PetscMalloc1(len,&s);CHKERRQ(ierr);
12536a9e3b9SBarry Smith       ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer);CHKERRQ(ierr);
12636a9e3b9SBarry Smith #undef len
127643535aeSDmitry Karpeev       ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D with overlap:\n", rank, size, i);CHKERRQ(ierr);
12892bb6962SLisandro Dalcin       for (j=0; j<nidx; j++) {
129643535aeSDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr);
13092bb6962SLisandro Dalcin       }
13192bb6962SLisandro Dalcin       ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr);
132643535aeSDmitry Karpeev       ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr);
133643535aeSDmitry Karpeev       ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
134c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
135643535aeSDmitry Karpeev       ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr);
136643535aeSDmitry Karpeev       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
137c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
138643535aeSDmitry Karpeev       ierr = PetscFree(s);CHKERRQ(ierr);
1392b691e39SMatthew Knepley       if (osm->is_local) {
140643535aeSDmitry Karpeev         /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
14136a9e3b9SBarry Smith #define len  16*(nidx+1)+512
14236a9e3b9SBarry Smith         ierr = PetscMalloc1(len, &s);CHKERRQ(ierr);
14336a9e3b9SBarry Smith         ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer);CHKERRQ(ierr);
14436a9e3b9SBarry Smith #undef len
14509d011a0SDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D without overlap:\n", rank, size, i);CHKERRQ(ierr);
1462b691e39SMatthew Knepley         ierr = ISGetLocalSize(osm->is_local[i],&nidx);CHKERRQ(ierr);
1472b691e39SMatthew Knepley         ierr = ISGetIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
1482b691e39SMatthew Knepley         for (j=0; j<nidx; j++) {
14909d011a0SDmitry Karpeev           ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr);
1502b691e39SMatthew Knepley         }
1512b691e39SMatthew Knepley         ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
15209d011a0SDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr);
153643535aeSDmitry Karpeev         ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
154c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
155643535aeSDmitry Karpeev         ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr);
156643535aeSDmitry Karpeev         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
157c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
158643535aeSDmitry Karpeev         ierr = PetscFree(s);CHKERRQ(ierr);
159643535aeSDmitry Karpeev       }
1602fa5cd67SKarl Rupp     } else {
161643535aeSDmitry Karpeev       /* Participate in collective viewer calls. */
162c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
163643535aeSDmitry Karpeev       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
164c5e4d11fSDmitry Karpeev       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
165643535aeSDmitry Karpeev       /* Assume either all ranks have is_local or none do. */
166643535aeSDmitry Karpeev       if (osm->is_local) {
167c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
168643535aeSDmitry Karpeev         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
169c5e4d11fSDmitry Karpeev         ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
170643535aeSDmitry Karpeev       }
171fdc48646SDmitry Karpeev     }
17292bb6962SLisandro Dalcin   }
17392bb6962SLisandro Dalcin   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
174fcfd50ebSBarry Smith   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
17592bb6962SLisandro Dalcin   PetscFunctionReturn(0);
17692bb6962SLisandro Dalcin }
17792bb6962SLisandro Dalcin 
1786849ba73SBarry Smith static PetscErrorCode PCSetUp_ASM(PC pc)
1794b9ad928SBarry Smith {
1804b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
1816849ba73SBarry Smith   PetscErrorCode ierr;
182ace3abfcSBarry Smith   PetscBool      symset,flg;
18387e86712SBarry Smith   PetscInt       i,m,m_local;
1844b9ad928SBarry Smith   MatReuse       scall = MAT_REUSE_MATRIX;
1854b9ad928SBarry Smith   IS             isl;
1864b9ad928SBarry Smith   KSP            ksp;
1874b9ad928SBarry Smith   PC             subpc;
1882dcb1b2aSMatthew Knepley   const char     *prefix,*pprefix;
18923ce1328SBarry Smith   Vec            vec;
1900298fd71SBarry Smith   DM             *domain_dm = NULL;
1914b9ad928SBarry Smith 
1924b9ad928SBarry Smith   PetscFunctionBegin;
1934b9ad928SBarry Smith   if (!pc->setupcalled) {
194265a2120SBarry Smith     PetscInt m;
19592bb6962SLisandro Dalcin 
19692bb6962SLisandro Dalcin     if (!osm->type_set) {
19792bb6962SLisandro Dalcin       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
1982fa5cd67SKarl Rupp       if (symset && flg) osm->type = PC_ASM_BASIC;
19992bb6962SLisandro Dalcin     }
20092bb6962SLisandro Dalcin 
2012b837212SDmitry Karpeev     /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */
2022b837212SDmitry Karpeev     if (osm->n_local_true == PETSC_DECIDE) {
20369ca1f37SDmitry Karpeev       /* no subdomains given */
20465db9045SDmitry Karpeev       /* try pc->dm first, if allowed */
205d709ea83SDmitry Karpeev       if (osm->dm_subdomains && pc->dm) {
206feb221f8SDmitry Karpeev         PetscInt  num_domains, d;
207feb221f8SDmitry Karpeev         char      **domain_names;
2088d4ac253SDmitry Karpeev         IS        *inner_domain_is, *outer_domain_is;
2098d4ac253SDmitry Karpeev         ierr = DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm);CHKERRQ(ierr);
210704f0395SPatrick Sanan         osm->overlap = -1; /* We do not want to increase the overlap of the IS.
211704f0395SPatrick Sanan                               A future improvement of this code might allow one to use
212704f0395SPatrick Sanan                               DM-defined subdomains and also increase the overlap,
213704f0395SPatrick Sanan                               but that is not currently supported */
21469ca1f37SDmitry Karpeev         if (num_domains) {
2158d4ac253SDmitry Karpeev           ierr = PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is);CHKERRQ(ierr);
21669ca1f37SDmitry Karpeev         }
217feb221f8SDmitry Karpeev         for (d = 0; d < num_domains; ++d) {
218a35b7b57SDmitry Karpeev           if (domain_names)    {ierr = PetscFree(domain_names[d]);CHKERRQ(ierr);}
219a35b7b57SDmitry Karpeev           if (inner_domain_is) {ierr = ISDestroy(&inner_domain_is[d]);CHKERRQ(ierr);}
220a35b7b57SDmitry Karpeev           if (outer_domain_is) {ierr = ISDestroy(&outer_domain_is[d]);CHKERRQ(ierr);}
221feb221f8SDmitry Karpeev         }
222feb221f8SDmitry Karpeev         ierr = PetscFree(domain_names);CHKERRQ(ierr);
2238d4ac253SDmitry Karpeev         ierr = PetscFree(inner_domain_is);CHKERRQ(ierr);
2248d4ac253SDmitry Karpeev         ierr = PetscFree(outer_domain_is);CHKERRQ(ierr);
225feb221f8SDmitry Karpeev       }
2262b837212SDmitry Karpeev       if (osm->n_local_true == PETSC_DECIDE) {
22769ca1f37SDmitry Karpeev         /* still no subdomains; use one subdomain per processor */
2282b837212SDmitry Karpeev         osm->n_local_true = 1;
229feb221f8SDmitry Karpeev       }
2302b837212SDmitry Karpeev     }
2312b837212SDmitry Karpeev     { /* determine the global and max number of subdomains */
2326ac3741eSJed Brown       struct {PetscInt max,sum;} inwork,outwork;
233c10200c1SHong Zhang       PetscMPIInt size;
234c10200c1SHong Zhang 
2356ac3741eSJed Brown       inwork.max   = osm->n_local_true;
2366ac3741eSJed Brown       inwork.sum   = osm->n_local_true;
237367daffbSBarry Smith       ierr         = MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
2386ac3741eSJed Brown       osm->n_local = outwork.max;
2396ac3741eSJed Brown       osm->n       = outwork.sum;
240c10200c1SHong Zhang 
241c10200c1SHong Zhang       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
242c10200c1SHong Zhang       if (outwork.max == 1 && outwork.sum == size) {
2437dae84e0SHong Zhang         /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */
244c10200c1SHong Zhang         ierr = MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE);CHKERRQ(ierr);
245c10200c1SHong Zhang       }
2464b9ad928SBarry Smith     }
24778904715SLisandro Dalcin     if (!osm->is) { /* create the index sets */
24878904715SLisandro Dalcin       ierr = PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);CHKERRQ(ierr);
2494b9ad928SBarry Smith     }
250f5234e35SJed Brown     if (osm->n_local_true > 1 && !osm->is_local) {
251785e854fSJed Brown       ierr = PetscMalloc1(osm->n_local_true,&osm->is_local);CHKERRQ(ierr);
252f5234e35SJed Brown       for (i=0; i<osm->n_local_true; i++) {
253f5234e35SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
254f5234e35SJed Brown           ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr);
255f5234e35SJed Brown           ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr);
256f5234e35SJed Brown         } else {
257f5234e35SJed Brown           ierr             = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr);
258f5234e35SJed Brown           osm->is_local[i] = osm->is[i];
259f5234e35SJed Brown         }
260f5234e35SJed Brown       }
261f5234e35SJed Brown     }
26292bb6962SLisandro Dalcin     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
26390d69ab7SBarry Smith     flg  = PETSC_FALSE;
264c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,prefix,"-pc_asm_print_subdomains",&flg,NULL);CHKERRQ(ierr);
26592bb6962SLisandro Dalcin     if (flg) { ierr = PCASMPrintSubdomains(pc);CHKERRQ(ierr); }
2664b9ad928SBarry Smith 
2673d03bb48SJed Brown     if (osm->overlap > 0) {
2684b9ad928SBarry Smith       /* Extend the "overlapping" regions by a number of steps */
26992bb6962SLisandro Dalcin       ierr = MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);CHKERRQ(ierr);
2703d03bb48SJed Brown     }
2716ed231c7SMatthew Knepley     if (osm->sort_indices) {
27292bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2734b9ad928SBarry Smith         ierr = ISSort(osm->is[i]);CHKERRQ(ierr);
2742b691e39SMatthew Knepley         if (osm->is_local) {
2752b691e39SMatthew Knepley           ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr);
2762b691e39SMatthew Knepley         }
2774b9ad928SBarry Smith       }
2786ed231c7SMatthew Knepley     }
2794b9ad928SBarry Smith 
280f6991133SBarry Smith     if (!osm->ksp) {
28178904715SLisandro Dalcin       /* Create the local solvers */
282785e854fSJed Brown       ierr = PetscMalloc1(osm->n_local_true,&osm->ksp);CHKERRQ(ierr);
283feb221f8SDmitry Karpeev       if (domain_dm) {
284feb221f8SDmitry Karpeev         ierr = PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n");CHKERRQ(ierr);
285feb221f8SDmitry Karpeev       }
28692bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2874b9ad928SBarry Smith         ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
288422a814eSBarry Smith         ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr);
2893bb1ff40SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr);
29092bb6962SLisandro Dalcin         ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
2914b9ad928SBarry Smith         ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
2924b9ad928SBarry Smith         ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
2934b9ad928SBarry Smith         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
2944b9ad928SBarry Smith         ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
2954b9ad928SBarry Smith         ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
296feb221f8SDmitry Karpeev         if (domain_dm) {
297feb221f8SDmitry Karpeev           ierr = KSPSetDM(ksp, domain_dm[i]);CHKERRQ(ierr);
298feb221f8SDmitry Karpeev           ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr);
299feb221f8SDmitry Karpeev           ierr = DMDestroy(&domain_dm[i]);CHKERRQ(ierr);
300feb221f8SDmitry Karpeev         }
3014b9ad928SBarry Smith         osm->ksp[i] = ksp;
3024b9ad928SBarry Smith       }
303feb221f8SDmitry Karpeev       if (domain_dm) {
304feb221f8SDmitry Karpeev         ierr = PetscFree(domain_dm);CHKERRQ(ierr);
305feb221f8SDmitry Karpeev       }
306f6991133SBarry Smith     }
3071dd8081eSeaulisa 
308fb745f2cSMatthew G. Knepley     ierr = ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis);CHKERRQ(ierr);
309fb745f2cSMatthew G. Knepley     ierr = ISSortRemoveDups(osm->lis);CHKERRQ(ierr);
310fb745f2cSMatthew G. Knepley     ierr = ISGetLocalSize(osm->lis, &m);CHKERRQ(ierr);
311fb745f2cSMatthew G. Knepley     ierr = VecCreateSeq(PETSC_COMM_SELF, m, &osm->lx);CHKERRQ(ierr);
312fb745f2cSMatthew G. Knepley     ierr = VecDuplicate(osm->lx, &osm->ly);CHKERRQ(ierr);
3131dd8081eSeaulisa 
3144b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
3154b9ad928SBarry Smith   } else {
3164b9ad928SBarry Smith     /*
3174b9ad928SBarry Smith        Destroy the blocks from the previous iteration
3184b9ad928SBarry Smith     */
319e09e08ccSBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
3204b9ad928SBarry Smith       ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
3214b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
3224b9ad928SBarry Smith     }
3234b9ad928SBarry Smith   }
3244b9ad928SBarry Smith 
32592bb6962SLisandro Dalcin   /*
32692bb6962SLisandro Dalcin      Extract out the submatrices
32792bb6962SLisandro Dalcin   */
3287dae84e0SHong Zhang   ierr = MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr);
32992bb6962SLisandro Dalcin   if (scall == MAT_INITIAL_MATRIX) {
33092bb6962SLisandro Dalcin     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
33192bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
3323bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr);
33378904715SLisandro Dalcin       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
33492bb6962SLisandro Dalcin     }
33592bb6962SLisandro Dalcin   }
33680ec0b7dSPatrick Sanan 
33780ec0b7dSPatrick Sanan   /* Convert the types of the submatrices (if needbe) */
33880ec0b7dSPatrick Sanan   if (osm->sub_mat_type) {
33980ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; i++) {
34080ec0b7dSPatrick Sanan       ierr = MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]));CHKERRQ(ierr);
34180ec0b7dSPatrick Sanan     }
34280ec0b7dSPatrick Sanan   }
34380ec0b7dSPatrick Sanan 
34480ec0b7dSPatrick Sanan   if(!pc->setupcalled){
34580ec0b7dSPatrick Sanan     /* Create the local work vectors (from the local matrices) and scatter contexts */
34680ec0b7dSPatrick Sanan     ierr = MatCreateVecs(pc->pmat,&vec,0);CHKERRQ(ierr);
3475b3c0d42Seaulisa 
3481dd8081eSeaulisa     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()");
3491dd8081eSeaulisa     if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) {
3501dd8081eSeaulisa       ierr = PetscMalloc1(osm->n_local_true,&osm->lprolongation);CHKERRQ(ierr);
3511dd8081eSeaulisa     }
3521dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->lrestriction);CHKERRQ(ierr);
3531dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->x);CHKERRQ(ierr);
3541dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->y);CHKERRQ(ierr);
355347574c9Seaulisa 
356347574c9Seaulisa     ierr = ISGetLocalSize(osm->lis,&m);CHKERRQ(ierr);
357347574c9Seaulisa     ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
3589448b7f1SJunchao Zhang     ierr = VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->restriction);CHKERRQ(ierr);
359347574c9Seaulisa     ierr = ISDestroy(&isl);CHKERRQ(ierr);
360347574c9Seaulisa 
361347574c9Seaulisa 
36280ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; ++i) {
3635b3c0d42Seaulisa       ISLocalToGlobalMapping ltog;
3645b3c0d42Seaulisa       IS                     isll;
3655b3c0d42Seaulisa       const PetscInt         *idx_is;
3665b3c0d42Seaulisa       PetscInt               *idx_lis,nout;
3675b3c0d42Seaulisa 
36855817e79SBarry Smith       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
36955817e79SBarry Smith       ierr = MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);CHKERRQ(ierr);
37055817e79SBarry Smith       ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr);
37155817e79SBarry Smith 
372b0de9f23SBarry Smith       /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */
3735b3c0d42Seaulisa       ierr = ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);CHKERRQ(ierr);
3745b3c0d42Seaulisa       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
3755b3c0d42Seaulisa       ierr = ISGetIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
3765b3c0d42Seaulisa       ierr = PetscMalloc1(m,&idx_lis);CHKERRQ(ierr);
3775b3c0d42Seaulisa       ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis);CHKERRQ(ierr);
3785b3c0d42Seaulisa       if (nout != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis");
3795b3c0d42Seaulisa       ierr = ISRestoreIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
3805b3c0d42Seaulisa       ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
381d8b3b5e3Seaulisa       ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
3825b3c0d42Seaulisa       ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
3839448b7f1SJunchao Zhang       ierr = VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lrestriction[i]);CHKERRQ(ierr);
3845b3c0d42Seaulisa       ierr = ISDestroy(&isll);CHKERRQ(ierr);
3855b3c0d42Seaulisa       ierr = ISDestroy(&isl);CHKERRQ(ierr);
386b0de9f23SBarry Smith       if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overalapping is_local[i] entries */
387d8b3b5e3Seaulisa         ISLocalToGlobalMapping ltog;
388d8b3b5e3Seaulisa         IS                     isll,isll_local;
389d8b3b5e3Seaulisa         const PetscInt         *idx_local;
390d8b3b5e3Seaulisa         PetscInt               *idx1, *idx2, nout;
391d8b3b5e3Seaulisa 
392d8b3b5e3Seaulisa         ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr);
393d8b3b5e3Seaulisa         ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
394d8b3b5e3Seaulisa 
395d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],&ltog);CHKERRQ(ierr);
396d8b3b5e3Seaulisa         ierr = PetscMalloc1(m_local,&idx1);CHKERRQ(ierr);
397d8b3b5e3Seaulisa         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1);CHKERRQ(ierr);
398d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
399d8b3b5e3Seaulisa         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is");
400d8b3b5e3Seaulisa         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
401d8b3b5e3Seaulisa 
402d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);CHKERRQ(ierr);
403d8b3b5e3Seaulisa         ierr = PetscMalloc1(m_local,&idx2);CHKERRQ(ierr);
404d8b3b5e3Seaulisa         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2);CHKERRQ(ierr);
405d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
406d8b3b5e3Seaulisa         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis");
407d8b3b5e3Seaulisa         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local);CHKERRQ(ierr);
408d8b3b5e3Seaulisa 
409d8b3b5e3Seaulisa         ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
4109448b7f1SJunchao Zhang         ierr = VecScatterCreate(osm->y[i],isll,osm->ly,isll_local,&osm->lprolongation[i]);CHKERRQ(ierr);
411d8b3b5e3Seaulisa 
412d8b3b5e3Seaulisa         ierr = ISDestroy(&isll);CHKERRQ(ierr);
413d8b3b5e3Seaulisa         ierr = ISDestroy(&isll_local);CHKERRQ(ierr);
414d8b3b5e3Seaulisa       }
41580ec0b7dSPatrick Sanan     }
41680ec0b7dSPatrick Sanan     ierr = VecDestroy(&vec);CHKERRQ(ierr);
41780ec0b7dSPatrick Sanan   }
41880ec0b7dSPatrick Sanan 
419fb745f2cSMatthew G. Knepley   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
420235cc4ceSMatthew G. Knepley     IS      *cis;
421235cc4ceSMatthew G. Knepley     PetscInt c;
422235cc4ceSMatthew G. Knepley 
423235cc4ceSMatthew G. Knepley     ierr = PetscMalloc1(osm->n_local_true, &cis);CHKERRQ(ierr);
424235cc4ceSMatthew G. Knepley     for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis;
4257dae84e0SHong Zhang     ierr = MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);CHKERRQ(ierr);
426235cc4ceSMatthew G. Knepley     ierr = PetscFree(cis);CHKERRQ(ierr);
427fb745f2cSMatthew G. Knepley   }
4284b9ad928SBarry Smith 
4294b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
4304b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
43192bb6962SLisandro Dalcin   ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
4324b9ad928SBarry Smith 
43392bb6962SLisandro Dalcin   /*
43492bb6962SLisandro Dalcin      Loop over subdomains putting them into local ksp
43592bb6962SLisandro Dalcin   */
43692bb6962SLisandro Dalcin   for (i=0; i<osm->n_local_true; i++) {
43723ee1639SBarry Smith     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr);
43892bb6962SLisandro Dalcin     if (!pc->setupcalled) {
439bf108f30SBarry Smith       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
4404b9ad928SBarry Smith     }
44192bb6962SLisandro Dalcin   }
4424b9ad928SBarry Smith   PetscFunctionReturn(0);
4434b9ad928SBarry Smith }
4444b9ad928SBarry Smith 
4456849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc)
4464b9ad928SBarry Smith {
4474b9ad928SBarry Smith   PC_ASM             *osm = (PC_ASM*)pc->data;
4486849ba73SBarry Smith   PetscErrorCode     ierr;
44913f74950SBarry Smith   PetscInt           i;
450539c167fSBarry Smith   KSPConvergedReason reason;
4514b9ad928SBarry Smith 
4524b9ad928SBarry Smith   PetscFunctionBegin;
4534b9ad928SBarry Smith   for (i=0; i<osm->n_local_true; i++) {
4544b9ad928SBarry Smith     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
455539c167fSBarry Smith     ierr = KSPGetConvergedReason(osm->ksp[i],&reason);CHKERRQ(ierr);
456c0decd05SBarry Smith     if (reason == KSP_DIVERGED_PC_FAILED) {
457261222b2SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
458e0eafd54SHong Zhang     }
4594b9ad928SBarry Smith   }
4604b9ad928SBarry Smith   PetscFunctionReturn(0);
4614b9ad928SBarry Smith }
4624b9ad928SBarry Smith 
4636849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y)
4644b9ad928SBarry Smith {
4654b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
4666849ba73SBarry Smith   PetscErrorCode ierr;
4671dd8081eSeaulisa   PetscInt       i,n_local_true = osm->n_local_true;
4684b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
4694b9ad928SBarry Smith 
4704b9ad928SBarry Smith   PetscFunctionBegin;
4714b9ad928SBarry Smith   /*
4724b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
4734b9ad928SBarry Smith      subdomain values (leaving the other values 0).
4744b9ad928SBarry Smith   */
4754b9ad928SBarry Smith   if (!(osm->type & PC_ASM_RESTRICT)) {
4764b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
4774b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
4781dd8081eSeaulisa     ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr);
4794b9ad928SBarry Smith   }
480347574c9Seaulisa   if (!(osm->type & PC_ASM_INTERPOLATE)) {
481347574c9Seaulisa     reverse = SCATTER_REVERSE_LOCAL;
482347574c9Seaulisa   }
4834b9ad928SBarry Smith 
484347574c9Seaulisa   if(osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE){
485b0de9f23SBarry Smith     /* zero the global and the local solutions */
48612cd4985SMatthew G. Knepley     ierr = VecZeroEntries(y);CHKERRQ(ierr);
4875b3c0d42Seaulisa     ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
488347574c9Seaulisa 
489b0de9f23SBarry Smith     /* Copy the global RHS to local RHS including the ghost nodes */
4901dd8081eSeaulisa     ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
4911dd8081eSeaulisa     ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
492347574c9Seaulisa 
493b0de9f23SBarry Smith     /* Restrict local RHS to the overlapping 0-block RHS */
4941dd8081eSeaulisa     ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
4951dd8081eSeaulisa     ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0],  INSERT_VALUES, forward);CHKERRQ(ierr);
496d8b3b5e3Seaulisa 
49712cd4985SMatthew G. Knepley     /* do the local solves */
49812cd4985SMatthew G. Knepley     for (i = 0; i < n_local_true; ++i) {
499347574c9Seaulisa 
500b0de9f23SBarry Smith       /* solve the overlapping i-block */
501fa2bb9feSLisandro Dalcin       ierr = PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
50212cd4985SMatthew G. Knepley       ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
503c0decd05SBarry Smith       ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr);
504fa2bb9feSLisandro Dalcin       ierr = PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
505d8b3b5e3Seaulisa 
506b0de9f23SBarry Smith       if (osm->lprolongation) { /* interpolate the non-overalapping i-block solution to the local solution (only for restrictive additive) */
5071dd8081eSeaulisa         ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
5081dd8081eSeaulisa         ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
509d8b3b5e3Seaulisa       }
510b0de9f23SBarry Smith       else{ /* interpolate the overalapping i-block solution to the local solution */
5111dd8081eSeaulisa         ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
5121dd8081eSeaulisa         ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
513d8b3b5e3Seaulisa       }
514347574c9Seaulisa 
515347574c9Seaulisa       if (i < n_local_true-1) {
516b0de9f23SBarry Smith         /* Restrict local RHS to the overlapping (i+1)-block RHS */
5171dd8081eSeaulisa         ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
5181dd8081eSeaulisa         ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
519347574c9Seaulisa 
520347574c9Seaulisa         if ( osm->loctype == PC_COMPOSITE_MULTIPLICATIVE){
521b0de9f23SBarry Smith           /* udpdate the overlapping (i+1)-block RHS using the current local solution */
522347574c9Seaulisa           ierr = MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]);CHKERRQ(ierr);
523347574c9Seaulisa           ierr = VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1]); CHKERRQ(ierr);
5247c3d802fSMatthew G. Knepley         }
52512cd4985SMatthew G. Knepley       }
52612cd4985SMatthew G. Knepley     }
527b0de9f23SBarry Smith     /* Add the local solution to the global solution including the ghost nodes */
5281dd8081eSeaulisa     ierr = VecScatterBegin(osm->restriction, osm->ly, y,  ADD_VALUES, reverse);CHKERRQ(ierr);
5291dd8081eSeaulisa     ierr = VecScatterEnd(osm->restriction,  osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
530347574c9Seaulisa   }else{
531347574c9Seaulisa     SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
53212cd4985SMatthew G. Knepley   }
5334b9ad928SBarry Smith   PetscFunctionReturn(0);
5344b9ad928SBarry Smith }
5354b9ad928SBarry Smith 
5366849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y)
5374b9ad928SBarry Smith {
5384b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
5396849ba73SBarry Smith   PetscErrorCode ierr;
5401dd8081eSeaulisa   PetscInt       i,n_local_true = osm->n_local_true;
5414b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
5424b9ad928SBarry Smith 
5434b9ad928SBarry Smith   PetscFunctionBegin;
5444b9ad928SBarry Smith   /*
5454b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
5464b9ad928SBarry Smith      subdomain values (leaving the other values 0).
5474b9ad928SBarry Smith 
5484b9ad928SBarry Smith      Note: these are reversed from the PCApply_ASM() because we are applying the
5494b9ad928SBarry Smith      transpose of the three terms
5504b9ad928SBarry Smith   */
551d8b3b5e3Seaulisa 
5524b9ad928SBarry Smith   if (!(osm->type & PC_ASM_INTERPOLATE)) {
5534b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
5544b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
5551dd8081eSeaulisa     ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr);
5564b9ad928SBarry Smith   }
5572fa5cd67SKarl Rupp   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;
5584b9ad928SBarry Smith 
559b0de9f23SBarry Smith   /* zero the global and the local solutions */
56010bd88b9SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
561d8b3b5e3Seaulisa   ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
562d8b3b5e3Seaulisa 
563b0de9f23SBarry Smith   /* Copy the global RHS to local RHS including the ghost nodes */
5641dd8081eSeaulisa   ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
5651dd8081eSeaulisa   ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
566d8b3b5e3Seaulisa 
567b0de9f23SBarry Smith   /* Restrict local RHS to the overlapping 0-block RHS */
5681dd8081eSeaulisa   ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
5691dd8081eSeaulisa   ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0],  INSERT_VALUES, forward);CHKERRQ(ierr);
570d8b3b5e3Seaulisa 
5714b9ad928SBarry Smith   /* do the local solves */
572d8b3b5e3Seaulisa   for (i = 0; i < n_local_true; ++i) {
573d8b3b5e3Seaulisa 
574b0de9f23SBarry Smith     /* solve the overlapping i-block */
575fa2bb9feSLisandro Dalcin     ierr = PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
576e1bcd54cSBarry Smith     ierr = KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
577c0decd05SBarry Smith     ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr);
578fa2bb9feSLisandro Dalcin     ierr = PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
579d8b3b5e3Seaulisa 
580b0de9f23SBarry Smith     if (osm->lprolongation) { /* interpolate the non-overalapping i-block solution to the local solution */
5811dd8081eSeaulisa      ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
5821dd8081eSeaulisa      ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
583b9845e0eSMatthew Knepley     }
584b0de9f23SBarry Smith     else{ /* interpolate the overalapping i-block solution to the local solution */
5851dd8081eSeaulisa       ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
5861dd8081eSeaulisa       ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
5874b9ad928SBarry Smith     }
588d8b3b5e3Seaulisa 
589d8b3b5e3Seaulisa     if (i < n_local_true-1) {
590b0de9f23SBarry Smith       /* Restrict local RHS to the overlapping (i+1)-block RHS */
5911dd8081eSeaulisa       ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
5921dd8081eSeaulisa       ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
5934b9ad928SBarry Smith     }
5944b9ad928SBarry Smith   }
595b0de9f23SBarry Smith   /* Add the local solution to the global solution including the ghost nodes */
5961dd8081eSeaulisa   ierr = VecScatterBegin(osm->restriction, osm->ly, y,  ADD_VALUES, reverse);CHKERRQ(ierr);
5971dd8081eSeaulisa   ierr = VecScatterEnd(osm->restriction,  osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
598d8b3b5e3Seaulisa 
5994b9ad928SBarry Smith   PetscFunctionReturn(0);
600d8b3b5e3Seaulisa 
6014b9ad928SBarry Smith }
6024b9ad928SBarry Smith 
603e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc)
6044b9ad928SBarry Smith {
6054b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
6066849ba73SBarry Smith   PetscErrorCode ierr;
60713f74950SBarry Smith   PetscInt       i;
6084b9ad928SBarry Smith 
6094b9ad928SBarry Smith   PetscFunctionBegin;
61092bb6962SLisandro Dalcin   if (osm->ksp) {
61192bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
612e91c6855SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
61392bb6962SLisandro Dalcin     }
61492bb6962SLisandro Dalcin   }
615e09e08ccSBarry Smith   if (osm->pmat) {
61692bb6962SLisandro Dalcin     if (osm->n_local_true > 0) {
61730a70a9aSHong Zhang       ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
61892bb6962SLisandro Dalcin     }
61992bb6962SLisandro Dalcin   }
6201dd8081eSeaulisa   if (osm->lrestriction) {
6211dd8081eSeaulisa     ierr = VecScatterDestroy(&osm->restriction);CHKERRQ(ierr);
6221dd8081eSeaulisa     for (i=0; i<osm->n_local_true; i++) {
6231dd8081eSeaulisa       ierr = VecScatterDestroy(&osm->lrestriction[i]);CHKERRQ(ierr);
6241dd8081eSeaulisa       if (osm->lprolongation) {ierr = VecScatterDestroy(&osm->lprolongation[i]);CHKERRQ(ierr);}
625fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
626fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
6274b9ad928SBarry Smith     }
6281dd8081eSeaulisa     ierr = PetscFree(osm->lrestriction);CHKERRQ(ierr);
6291dd8081eSeaulisa     if (osm->lprolongation) {ierr = PetscFree(osm->lprolongation);CHKERRQ(ierr);}
63005b42c5fSBarry Smith     ierr = PetscFree(osm->x);CHKERRQ(ierr);
63178904715SLisandro Dalcin     ierr = PetscFree(osm->y);CHKERRQ(ierr);
6321dd8081eSeaulisa 
63392bb6962SLisandro Dalcin   }
6342b691e39SMatthew Knepley   ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
635fb745f2cSMatthew G. Knepley   ierr = ISDestroy(&osm->lis);CHKERRQ(ierr);
636fb745f2cSMatthew G. Knepley   ierr = VecDestroy(&osm->lx);CHKERRQ(ierr);
637fb745f2cSMatthew G. Knepley   ierr = VecDestroy(&osm->ly);CHKERRQ(ierr);
638347574c9Seaulisa   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
639347574c9Seaulisa     ierr = MatDestroyMatrices(osm->n_local_true, &osm->lmats);CHKERRQ(ierr);
640fb745f2cSMatthew G. Knepley   }
6412fa5cd67SKarl Rupp 
64280ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
64380ec0b7dSPatrick Sanan 
644e91c6855SBarry Smith   osm->is       = 0;
645e91c6855SBarry Smith   osm->is_local = 0;
646e91c6855SBarry Smith   PetscFunctionReturn(0);
647e91c6855SBarry Smith }
648e91c6855SBarry Smith 
649e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc)
650e91c6855SBarry Smith {
651e91c6855SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
652e91c6855SBarry Smith   PetscErrorCode ierr;
653e91c6855SBarry Smith   PetscInt       i;
654e91c6855SBarry Smith 
655e91c6855SBarry Smith   PetscFunctionBegin;
656e91c6855SBarry Smith   ierr = PCReset_ASM(pc);CHKERRQ(ierr);
657e91c6855SBarry Smith   if (osm->ksp) {
658e91c6855SBarry Smith     for (i=0; i<osm->n_local_true; i++) {
659fcfd50ebSBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
660e91c6855SBarry Smith     }
661e91c6855SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
662e91c6855SBarry Smith   }
663e91c6855SBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
6644b9ad928SBarry Smith   PetscFunctionReturn(0);
6654b9ad928SBarry Smith }
6664b9ad928SBarry Smith 
6674416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc)
6684b9ad928SBarry Smith {
6694b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
6706849ba73SBarry Smith   PetscErrorCode ierr;
6719dcbbd2bSBarry Smith   PetscInt       blocks,ovl;
672ace3abfcSBarry Smith   PetscBool      symset,flg;
67392bb6962SLisandro Dalcin   PCASMType      asmtype;
67412cd4985SMatthew G. Knepley   PCCompositeType loctype;
67580ec0b7dSPatrick Sanan   char           sub_mat_type[256];
6764b9ad928SBarry Smith 
6774b9ad928SBarry Smith   PetscFunctionBegin;
678bf108f30SBarry Smith   /* set the type to symmetric if matrix is symmetric */
67992bb6962SLisandro Dalcin   if (!osm->type_set && pc->pmat) {
68092bb6962SLisandro Dalcin     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
6812fa5cd67SKarl Rupp     if (symset && flg) osm->type = PC_ASM_BASIC;
682bf108f30SBarry Smith   }
683e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");CHKERRQ(ierr);
684d709ea83SDmitry Karpeev   ierr = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr);
68590d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
68665db9045SDmitry Karpeev   if (flg) {
687f77a5242SKarl Rupp     ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr);
688d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
68965db9045SDmitry Karpeev   }
690342c94f9SMatthew G. Knepley   ierr = PetscOptionsInt("-pc_asm_local_blocks","Number of local subdomains","PCASMSetLocalSubdomains",osm->n_local_true,&blocks,&flg);CHKERRQ(ierr);
691342c94f9SMatthew G. Knepley   if (flg) {
692342c94f9SMatthew G. Knepley     ierr = PCASMSetLocalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr);
693342c94f9SMatthew G. Knepley     osm->dm_subdomains = PETSC_FALSE;
694342c94f9SMatthew G. Knepley   }
69590d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
69665db9045SDmitry Karpeev   if (flg) {
69765db9045SDmitry Karpeev     ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr);
698d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
69965db9045SDmitry Karpeev   }
70090d69ab7SBarry Smith   flg  = PETSC_FALSE;
70190d69ab7SBarry Smith   ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
70292bb6962SLisandro Dalcin   if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); }
70312cd4985SMatthew G. Knepley   flg  = PETSC_FALSE;
70412cd4985SMatthew G. Knepley   ierr = PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr);
70512cd4985SMatthew G. Knepley   if (flg) {ierr = PCASMSetLocalType(pc,loctype);CHKERRQ(ierr); }
70680ec0b7dSPatrick Sanan   ierr = PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);CHKERRQ(ierr);
70780ec0b7dSPatrick Sanan   if(flg){
708459726d8SSatish Balay     ierr = PCASMSetSubMatType(pc,sub_mat_type);CHKERRQ(ierr);
70980ec0b7dSPatrick Sanan   }
7104b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7114b9ad928SBarry Smith   PetscFunctionReturn(0);
7124b9ad928SBarry Smith }
7134b9ad928SBarry Smith 
7144b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/
7154b9ad928SBarry Smith 
7161e6b0712SBarry Smith static PetscErrorCode  PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[])
7174b9ad928SBarry Smith {
7184b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
71992bb6962SLisandro Dalcin   PetscErrorCode ierr;
72092bb6962SLisandro Dalcin   PetscInt       i;
7214b9ad928SBarry Smith 
7224b9ad928SBarry Smith   PetscFunctionBegin;
723e32f2f54SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
724ce94432eSBarry 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().");
725e7e72b3dSBarry Smith 
7264b9ad928SBarry Smith   if (!pc->setupcalled) {
72792bb6962SLisandro Dalcin     if (is) {
72892bb6962SLisandro Dalcin       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
72992bb6962SLisandro Dalcin     }
730832fc9a5SMatthew Knepley     if (is_local) {
731832fc9a5SMatthew Knepley       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
732832fc9a5SMatthew Knepley     }
7332b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
7342fa5cd67SKarl Rupp 
7354b9ad928SBarry Smith     osm->n_local_true = n;
73692bb6962SLisandro Dalcin     osm->is           = 0;
7372b691e39SMatthew Knepley     osm->is_local     = 0;
73892bb6962SLisandro Dalcin     if (is) {
739785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is);CHKERRQ(ierr);
7402fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is[i] = is[i];
7413d03bb48SJed Brown       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
7423d03bb48SJed Brown       osm->overlap = -1;
74392bb6962SLisandro Dalcin     }
7442b691e39SMatthew Knepley     if (is_local) {
745785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is_local);CHKERRQ(ierr);
7462fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is_local[i] = is_local[i];
747a35b7b57SDmitry Karpeev       if (!is) {
748785e854fSJed Brown         ierr = PetscMalloc1(osm->n_local_true,&osm->is);CHKERRQ(ierr);
749a35b7b57SDmitry Karpeev         for (i=0; i<osm->n_local_true; i++) {
750a35b7b57SDmitry Karpeev           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
751a35b7b57SDmitry Karpeev             ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr);
752a35b7b57SDmitry Karpeev             ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr);
753a35b7b57SDmitry Karpeev           } else {
754a35b7b57SDmitry Karpeev             ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr);
755a35b7b57SDmitry Karpeev             osm->is[i] = osm->is_local[i];
756a35b7b57SDmitry Karpeev           }
757a35b7b57SDmitry Karpeev         }
758a35b7b57SDmitry Karpeev       }
7592b691e39SMatthew Knepley     }
7604b9ad928SBarry Smith   }
7614b9ad928SBarry Smith   PetscFunctionReturn(0);
7624b9ad928SBarry Smith }
7634b9ad928SBarry Smith 
7641e6b0712SBarry Smith static PetscErrorCode  PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local)
7654b9ad928SBarry Smith {
7664b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
7676849ba73SBarry Smith   PetscErrorCode ierr;
76813f74950SBarry Smith   PetscMPIInt    rank,size;
76978904715SLisandro Dalcin   PetscInt       n;
7704b9ad928SBarry Smith 
7714b9ad928SBarry Smith   PetscFunctionBegin;
772ce94432eSBarry Smith   if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
773ce94432eSBarry 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.");
7744b9ad928SBarry Smith 
7754b9ad928SBarry Smith   /*
776880770e9SJed Brown      Split the subdomains equally among all processors
7774b9ad928SBarry Smith   */
778ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
779ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
7804b9ad928SBarry Smith   n    = N/size + ((N % size) > rank);
781e32f2f54SBarry 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);
782e32f2f54SBarry Smith   if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp().");
7834b9ad928SBarry Smith   if (!pc->setupcalled) {
7842b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
7852fa5cd67SKarl Rupp 
7864b9ad928SBarry Smith     osm->n_local_true = n;
7874b9ad928SBarry Smith     osm->is           = 0;
7882b691e39SMatthew Knepley     osm->is_local     = 0;
7894b9ad928SBarry Smith   }
7904b9ad928SBarry Smith   PetscFunctionReturn(0);
7914b9ad928SBarry Smith }
7924b9ad928SBarry Smith 
7931e6b0712SBarry Smith static PetscErrorCode  PCASMSetOverlap_ASM(PC pc,PetscInt ovl)
7944b9ad928SBarry Smith {
79592bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
7964b9ad928SBarry Smith 
7974b9ad928SBarry Smith   PetscFunctionBegin;
798ce94432eSBarry Smith   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
799ce94432eSBarry Smith   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp().");
8002fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
8014b9ad928SBarry Smith   PetscFunctionReturn(0);
8024b9ad928SBarry Smith }
8034b9ad928SBarry Smith 
8041e6b0712SBarry Smith static PetscErrorCode  PCASMSetType_ASM(PC pc,PCASMType type)
8054b9ad928SBarry Smith {
80692bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
8074b9ad928SBarry Smith 
8084b9ad928SBarry Smith   PetscFunctionBegin;
8094b9ad928SBarry Smith   osm->type     = type;
810bf108f30SBarry Smith   osm->type_set = PETSC_TRUE;
8114b9ad928SBarry Smith   PetscFunctionReturn(0);
8124b9ad928SBarry Smith }
8134b9ad928SBarry Smith 
814c60c7ad4SBarry Smith static PetscErrorCode  PCASMGetType_ASM(PC pc,PCASMType *type)
815c60c7ad4SBarry Smith {
816c60c7ad4SBarry Smith   PC_ASM *osm = (PC_ASM*)pc->data;
817c60c7ad4SBarry Smith 
818c60c7ad4SBarry Smith   PetscFunctionBegin;
819c60c7ad4SBarry Smith   *type = osm->type;
820c60c7ad4SBarry Smith   PetscFunctionReturn(0);
821c60c7ad4SBarry Smith }
822c60c7ad4SBarry Smith 
82312cd4985SMatthew G. Knepley static PetscErrorCode  PCASMSetLocalType_ASM(PC pc, PCCompositeType type)
82412cd4985SMatthew G. Knepley {
82512cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
82612cd4985SMatthew G. Knepley 
82712cd4985SMatthew G. Knepley   PetscFunctionBegin;
828d0ecd4dfSBarry 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");
82912cd4985SMatthew G. Knepley   osm->loctype = type;
83012cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
83112cd4985SMatthew G. Knepley }
83212cd4985SMatthew G. Knepley 
83312cd4985SMatthew G. Knepley static PetscErrorCode  PCASMGetLocalType_ASM(PC pc, PCCompositeType *type)
83412cd4985SMatthew G. Knepley {
83512cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
83612cd4985SMatthew G. Knepley 
83712cd4985SMatthew G. Knepley   PetscFunctionBegin;
83812cd4985SMatthew G. Knepley   *type = osm->loctype;
83912cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
84012cd4985SMatthew G. Knepley }
84112cd4985SMatthew G. Knepley 
8421e6b0712SBarry Smith static PetscErrorCode  PCASMSetSortIndices_ASM(PC pc,PetscBool  doSort)
8436ed231c7SMatthew Knepley {
8446ed231c7SMatthew Knepley   PC_ASM *osm = (PC_ASM*)pc->data;
8456ed231c7SMatthew Knepley 
8466ed231c7SMatthew Knepley   PetscFunctionBegin;
8476ed231c7SMatthew Knepley   osm->sort_indices = doSort;
8486ed231c7SMatthew Knepley   PetscFunctionReturn(0);
8496ed231c7SMatthew Knepley }
8506ed231c7SMatthew Knepley 
8511e6b0712SBarry Smith static PetscErrorCode  PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
8524b9ad928SBarry Smith {
85392bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
854dfbe8321SBarry Smith   PetscErrorCode ierr;
8554b9ad928SBarry Smith 
8564b9ad928SBarry Smith   PetscFunctionBegin;
857ce94432eSBarry 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");
8584b9ad928SBarry Smith 
8592fa5cd67SKarl Rupp   if (n_local) *n_local = osm->n_local_true;
86092bb6962SLisandro Dalcin   if (first_local) {
861ce94432eSBarry Smith     ierr          = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
86292bb6962SLisandro Dalcin     *first_local -= osm->n_local_true;
86392bb6962SLisandro Dalcin   }
86492bb6962SLisandro Dalcin   if (ksp) {
86592bb6962SLisandro Dalcin     /* Assume that local solves are now different; not necessarily
86692bb6962SLisandro Dalcin        true though!  This flag is used only for PCView_ASM() */
86792bb6962SLisandro Dalcin     *ksp                   = osm->ksp;
86892bb6962SLisandro Dalcin     osm->same_local_solves = PETSC_FALSE;
86992bb6962SLisandro Dalcin   }
8704b9ad928SBarry Smith   PetscFunctionReturn(0);
8714b9ad928SBarry Smith }
8724b9ad928SBarry Smith 
87380ec0b7dSPatrick Sanan static PetscErrorCode  PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type)
87480ec0b7dSPatrick Sanan {
87580ec0b7dSPatrick Sanan   PC_ASM         *osm = (PC_ASM*)pc->data;
87680ec0b7dSPatrick Sanan 
87780ec0b7dSPatrick Sanan   PetscFunctionBegin;
87880ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
87980ec0b7dSPatrick Sanan   PetscValidPointer(sub_mat_type,2);
88080ec0b7dSPatrick Sanan   *sub_mat_type = osm->sub_mat_type;
88180ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
88280ec0b7dSPatrick Sanan }
88380ec0b7dSPatrick Sanan 
88480ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type)
88580ec0b7dSPatrick Sanan {
88680ec0b7dSPatrick Sanan   PetscErrorCode    ierr;
88780ec0b7dSPatrick Sanan   PC_ASM            *osm = (PC_ASM*)pc->data;
88880ec0b7dSPatrick Sanan 
88980ec0b7dSPatrick Sanan   PetscFunctionBegin;
89080ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
89180ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
89280ec0b7dSPatrick Sanan   ierr = PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);CHKERRQ(ierr);
89380ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
89480ec0b7dSPatrick Sanan }
89580ec0b7dSPatrick Sanan 
8964b9ad928SBarry Smith /*@C
8971093a601SBarry Smith     PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner.
8984b9ad928SBarry Smith 
899d083f849SBarry Smith     Collective on pc
9004b9ad928SBarry Smith 
9014b9ad928SBarry Smith     Input Parameters:
9024b9ad928SBarry Smith +   pc - the preconditioner context
9034b9ad928SBarry Smith .   n - the number of subdomains for this processor (default value = 1)
9048c03b21aSDmitry Karpeev .   is - the index set that defines the subdomains for this processor
9050298fd71SBarry Smith          (or NULL for PETSc to determine subdomains)
906f1ee410cSBarry 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
907f1ee410cSBarry Smith          (or NULL to not provide these)
9084b9ad928SBarry Smith 
909342c94f9SMatthew G. Knepley     Options Database Key:
910342c94f9SMatthew G. Knepley     To set the total number of subdomain blocks rather than specify the
911342c94f9SMatthew G. Knepley     index sets, use the option
912342c94f9SMatthew G. Knepley .    -pc_asm_local_blocks <blks> - Sets local blocks
913342c94f9SMatthew G. Knepley 
9144b9ad928SBarry Smith     Notes:
9151093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
9164b9ad928SBarry Smith 
9174b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
9184b9ad928SBarry Smith 
9194b9ad928SBarry Smith     Use PCASMSetTotalSubdomains() to set the subdomains for all processors.
9204b9ad928SBarry Smith 
921f1ee410cSBarry Smith     If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated
922f1ee410cSBarry Smith     back to form the global solution (this is the standard restricted additive Schwarz method)
923f1ee410cSBarry Smith 
924f1ee410cSBarry 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
925f1ee410cSBarry Smith     no code to handle that case.
926f1ee410cSBarry Smith 
9274b9ad928SBarry Smith     Level: advanced
9284b9ad928SBarry Smith 
9294b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
930f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType()
9314b9ad928SBarry Smith @*/
9327087cfbeSBarry Smith PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
9334b9ad928SBarry Smith {
9347bb14e67SBarry Smith   PetscErrorCode ierr;
9354b9ad928SBarry Smith 
9364b9ad928SBarry Smith   PetscFunctionBegin;
9370700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9387bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr);
9394b9ad928SBarry Smith   PetscFunctionReturn(0);
9404b9ad928SBarry Smith }
9414b9ad928SBarry Smith 
9424b9ad928SBarry Smith /*@C
943feb221f8SDmitry Karpeev     PCASMSetTotalSubdomains - Sets the subdomains for all processors for the
9444b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
9454b9ad928SBarry Smith     PC communicator must call this routine, with the same index sets.
9464b9ad928SBarry Smith 
947d083f849SBarry Smith     Collective on pc
9484b9ad928SBarry Smith 
9494b9ad928SBarry Smith     Input Parameters:
9504b9ad928SBarry Smith +   pc - the preconditioner context
951feb221f8SDmitry Karpeev .   N  - the number of subdomains for all processors
952feb221f8SDmitry Karpeev .   is - the index sets that define the subdomains for all processors
953dfaa0c5fSBarry Smith          (or NULL to ask PETSc to determine the subdomains)
9542b691e39SMatthew Knepley -   is_local - the index sets that define the local part of the subdomains for this processor
955f1ee410cSBarry Smith          (or NULL to not provide this information)
9564b9ad928SBarry Smith 
9574b9ad928SBarry Smith     Options Database Key:
9584b9ad928SBarry Smith     To set the total number of subdomain blocks rather than specify the
9594b9ad928SBarry Smith     index sets, use the option
9604b9ad928SBarry Smith .    -pc_asm_blocks <blks> - Sets total blocks
9614b9ad928SBarry Smith 
9624b9ad928SBarry Smith     Notes:
963f1ee410cSBarry Smith     Currently you cannot use this to set the actual subdomains with the argument is or is_local.
9644b9ad928SBarry Smith 
9654b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
9664b9ad928SBarry Smith 
9674b9ad928SBarry Smith     These index sets cannot be destroyed until after completion of the
9684b9ad928SBarry Smith     linear solves for which the ASM preconditioner is being used.
9694b9ad928SBarry Smith 
9704b9ad928SBarry Smith     Use PCASMSetLocalSubdomains() to set local subdomains.
9714b9ad928SBarry Smith 
9721093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
9731093a601SBarry Smith 
9744b9ad928SBarry Smith     Level: advanced
9754b9ad928SBarry Smith 
9764b9ad928SBarry Smith .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
9774b9ad928SBarry Smith           PCASMCreateSubdomains2D()
9784b9ad928SBarry Smith @*/
9797087cfbeSBarry Smith PetscErrorCode  PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[])
9804b9ad928SBarry Smith {
9817bb14e67SBarry Smith   PetscErrorCode ierr;
9824b9ad928SBarry Smith 
9834b9ad928SBarry Smith   PetscFunctionBegin;
9840700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9857bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr);
9864b9ad928SBarry Smith   PetscFunctionReturn(0);
9874b9ad928SBarry Smith }
9884b9ad928SBarry Smith 
9894b9ad928SBarry Smith /*@
9904b9ad928SBarry Smith     PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
9914b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
992f1ee410cSBarry Smith     PC communicator must call this routine.
9934b9ad928SBarry Smith 
994d083f849SBarry Smith     Logically Collective on pc
9954b9ad928SBarry Smith 
9964b9ad928SBarry Smith     Input Parameters:
9974b9ad928SBarry Smith +   pc  - the preconditioner context
9984b9ad928SBarry Smith -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
9994b9ad928SBarry Smith 
10004b9ad928SBarry Smith     Options Database Key:
10014b9ad928SBarry Smith .   -pc_asm_overlap <ovl> - Sets overlap
10024b9ad928SBarry Smith 
10034b9ad928SBarry Smith     Notes:
10044b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.  To use
10054b9ad928SBarry Smith     multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
10064b9ad928SBarry Smith     PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).
10074b9ad928SBarry Smith 
10084b9ad928SBarry Smith     The overlap defaults to 1, so if one desires that no additional
10094b9ad928SBarry Smith     overlap be computed beyond what may have been set with a call to
10104b9ad928SBarry Smith     PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
10114b9ad928SBarry Smith     must be set to be 0.  In particular, if one does not explicitly set
10124b9ad928SBarry Smith     the subdomains an application code, then all overlap would be computed
10134b9ad928SBarry Smith     internally by PETSc, and using an overlap of 0 would result in an ASM
10144b9ad928SBarry Smith     variant that is equivalent to the block Jacobi preconditioner.
10154b9ad928SBarry Smith 
1016f1ee410cSBarry Smith     The default algorithm used by PETSc to increase overlap is fast, but not scalable,
1017f1ee410cSBarry Smith     use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
1018f1ee410cSBarry Smith 
10194b9ad928SBarry Smith     Note that one can define initial index sets with any overlap via
1020f1ee410cSBarry Smith     PCASMSetLocalSubdomains(); the routine
10214b9ad928SBarry Smith     PCASMSetOverlap() merely allows PETSc to extend that overlap further
10224b9ad928SBarry Smith     if desired.
10234b9ad928SBarry Smith 
10244b9ad928SBarry Smith     Level: intermediate
10254b9ad928SBarry Smith 
10264b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
1027f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap()
10284b9ad928SBarry Smith @*/
10297087cfbeSBarry Smith PetscErrorCode  PCASMSetOverlap(PC pc,PetscInt ovl)
10304b9ad928SBarry Smith {
10317bb14e67SBarry Smith   PetscErrorCode ierr;
10324b9ad928SBarry Smith 
10334b9ad928SBarry Smith   PetscFunctionBegin;
10340700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1035c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,ovl,2);
10367bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
10374b9ad928SBarry Smith   PetscFunctionReturn(0);
10384b9ad928SBarry Smith }
10394b9ad928SBarry Smith 
10404b9ad928SBarry Smith /*@
10414b9ad928SBarry Smith     PCASMSetType - Sets the type of restriction and interpolation used
10424b9ad928SBarry Smith     for local problems in the additive Schwarz method.
10434b9ad928SBarry Smith 
1044d083f849SBarry Smith     Logically Collective on pc
10454b9ad928SBarry Smith 
10464b9ad928SBarry Smith     Input Parameters:
10474b9ad928SBarry Smith +   pc  - the preconditioner context
10484b9ad928SBarry Smith -   type - variant of ASM, one of
10494b9ad928SBarry Smith .vb
10504b9ad928SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
10514b9ad928SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
10524b9ad928SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
10534b9ad928SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
10544b9ad928SBarry Smith .ve
10554b9ad928SBarry Smith 
10564b9ad928SBarry Smith     Options Database Key:
10574b9ad928SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
10584b9ad928SBarry Smith 
105995452b02SPatrick Sanan     Notes:
106095452b02SPatrick Sanan     if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected
1061f1ee410cSBarry Smith     to limit the local processor interpolation
1062f1ee410cSBarry Smith 
10634b9ad928SBarry Smith     Level: intermediate
10644b9ad928SBarry Smith 
10654b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1066f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType()
10674b9ad928SBarry Smith @*/
10687087cfbeSBarry Smith PetscErrorCode  PCASMSetType(PC pc,PCASMType type)
10694b9ad928SBarry Smith {
10707bb14e67SBarry Smith   PetscErrorCode ierr;
10714b9ad928SBarry Smith 
10724b9ad928SBarry Smith   PetscFunctionBegin;
10730700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1074c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
10757bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr);
10764b9ad928SBarry Smith   PetscFunctionReturn(0);
10774b9ad928SBarry Smith }
10784b9ad928SBarry Smith 
1079c60c7ad4SBarry Smith /*@
1080c60c7ad4SBarry Smith     PCASMGetType - Gets the type of restriction and interpolation used
1081c60c7ad4SBarry Smith     for local problems in the additive Schwarz method.
1082c60c7ad4SBarry Smith 
1083d083f849SBarry Smith     Logically Collective on pc
1084c60c7ad4SBarry Smith 
1085c60c7ad4SBarry Smith     Input Parameter:
1086c60c7ad4SBarry Smith .   pc  - the preconditioner context
1087c60c7ad4SBarry Smith 
1088c60c7ad4SBarry Smith     Output Parameter:
1089c60c7ad4SBarry Smith .   type - variant of ASM, one of
1090c60c7ad4SBarry Smith 
1091c60c7ad4SBarry Smith .vb
1092c60c7ad4SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
1093c60c7ad4SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1094c60c7ad4SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1095c60c7ad4SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
1096c60c7ad4SBarry Smith .ve
1097c60c7ad4SBarry Smith 
1098c60c7ad4SBarry Smith     Options Database Key:
1099c60c7ad4SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
1100c60c7ad4SBarry Smith 
1101c60c7ad4SBarry Smith     Level: intermediate
1102c60c7ad4SBarry Smith 
1103c60c7ad4SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1104f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType()
1105c60c7ad4SBarry Smith @*/
1106c60c7ad4SBarry Smith PetscErrorCode  PCASMGetType(PC pc,PCASMType *type)
1107c60c7ad4SBarry Smith {
1108c60c7ad4SBarry Smith   PetscErrorCode ierr;
1109c60c7ad4SBarry Smith 
1110c60c7ad4SBarry Smith   PetscFunctionBegin;
1111c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1112c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));CHKERRQ(ierr);
1113c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1114c60c7ad4SBarry Smith }
1115c60c7ad4SBarry Smith 
111612cd4985SMatthew G. Knepley /*@
111712cd4985SMatthew G. Knepley   PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method.
111812cd4985SMatthew G. Knepley 
1119d083f849SBarry Smith   Logically Collective on pc
112012cd4985SMatthew G. Knepley 
112112cd4985SMatthew G. Knepley   Input Parameters:
112212cd4985SMatthew G. Knepley + pc  - the preconditioner context
112312cd4985SMatthew G. Knepley - type - type of composition, one of
112412cd4985SMatthew G. Knepley .vb
112512cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
112612cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
112712cd4985SMatthew G. Knepley .ve
112812cd4985SMatthew G. Knepley 
112912cd4985SMatthew G. Knepley   Options Database Key:
113012cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
113112cd4985SMatthew G. Knepley 
113212cd4985SMatthew G. Knepley   Level: intermediate
113312cd4985SMatthew G. Knepley 
1134f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
113512cd4985SMatthew G. Knepley @*/
113612cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type)
113712cd4985SMatthew G. Knepley {
113812cd4985SMatthew G. Knepley   PetscErrorCode ierr;
113912cd4985SMatthew G. Knepley 
114012cd4985SMatthew G. Knepley   PetscFunctionBegin;
114112cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
114212cd4985SMatthew G. Knepley   PetscValidLogicalCollectiveEnum(pc, type, 2);
114312cd4985SMatthew G. Knepley   ierr = PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));CHKERRQ(ierr);
114412cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
114512cd4985SMatthew G. Knepley }
114612cd4985SMatthew G. Knepley 
114712cd4985SMatthew G. Knepley /*@
114812cd4985SMatthew G. Knepley   PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method.
114912cd4985SMatthew G. Knepley 
1150d083f849SBarry Smith   Logically Collective on pc
115112cd4985SMatthew G. Knepley 
115212cd4985SMatthew G. Knepley   Input Parameter:
115312cd4985SMatthew G. Knepley . pc  - the preconditioner context
115412cd4985SMatthew G. Knepley 
115512cd4985SMatthew G. Knepley   Output Parameter:
115612cd4985SMatthew G. Knepley . type - type of composition, one of
115712cd4985SMatthew G. Knepley .vb
115812cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
115912cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
116012cd4985SMatthew G. Knepley .ve
116112cd4985SMatthew G. Knepley 
116212cd4985SMatthew G. Knepley   Options Database Key:
116312cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
116412cd4985SMatthew G. Knepley 
116512cd4985SMatthew G. Knepley   Level: intermediate
116612cd4985SMatthew G. Knepley 
1167f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
116812cd4985SMatthew G. Knepley @*/
116912cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type)
117012cd4985SMatthew G. Knepley {
117112cd4985SMatthew G. Knepley   PetscErrorCode ierr;
117212cd4985SMatthew G. Knepley 
117312cd4985SMatthew G. Knepley   PetscFunctionBegin;
117412cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
117512cd4985SMatthew G. Knepley   PetscValidPointer(type, 2);
117612cd4985SMatthew G. Knepley   ierr = PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));CHKERRQ(ierr);
117712cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
117812cd4985SMatthew G. Knepley }
117912cd4985SMatthew G. Knepley 
11806ed231c7SMatthew Knepley /*@
11816ed231c7SMatthew Knepley     PCASMSetSortIndices - Determines whether subdomain indices are sorted.
11826ed231c7SMatthew Knepley 
1183d083f849SBarry Smith     Logically Collective on pc
11846ed231c7SMatthew Knepley 
11856ed231c7SMatthew Knepley     Input Parameters:
11866ed231c7SMatthew Knepley +   pc  - the preconditioner context
11876ed231c7SMatthew Knepley -   doSort - sort the subdomain indices
11886ed231c7SMatthew Knepley 
11896ed231c7SMatthew Knepley     Level: intermediate
11906ed231c7SMatthew Knepley 
11916ed231c7SMatthew Knepley .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
11926ed231c7SMatthew Knepley           PCASMCreateSubdomains2D()
11936ed231c7SMatthew Knepley @*/
11947087cfbeSBarry Smith PetscErrorCode  PCASMSetSortIndices(PC pc,PetscBool doSort)
11956ed231c7SMatthew Knepley {
11967bb14e67SBarry Smith   PetscErrorCode ierr;
11976ed231c7SMatthew Knepley 
11986ed231c7SMatthew Knepley   PetscFunctionBegin;
11990700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1200acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
12017bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
12026ed231c7SMatthew Knepley   PetscFunctionReturn(0);
12036ed231c7SMatthew Knepley }
12046ed231c7SMatthew Knepley 
12054b9ad928SBarry Smith /*@C
12064b9ad928SBarry Smith    PCASMGetSubKSP - Gets the local KSP contexts for all blocks on
12074b9ad928SBarry Smith    this processor.
12084b9ad928SBarry Smith 
1209d083f849SBarry Smith    Collective on pc iff first_local is requested
12104b9ad928SBarry Smith 
12114b9ad928SBarry Smith    Input Parameter:
12124b9ad928SBarry Smith .  pc - the preconditioner context
12134b9ad928SBarry Smith 
12144b9ad928SBarry Smith    Output Parameters:
12150298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
12160298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
12170298fd71SBarry Smith                  all processors must request or all must pass NULL
12184b9ad928SBarry Smith -  ksp - the array of KSP contexts
12194b9ad928SBarry Smith 
12204b9ad928SBarry Smith    Note:
1221d29017ddSJed Brown    After PCASMGetSubKSP() the array of KSPes is not to be freed.
12224b9ad928SBarry Smith 
12234b9ad928SBarry Smith    You must call KSPSetUp() before calling PCASMGetSubKSP().
12244b9ad928SBarry Smith 
1225d29017ddSJed Brown    Fortran note:
12262bf68e3eSBarry 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.
1227d29017ddSJed Brown 
12284b9ad928SBarry Smith    Level: advanced
12294b9ad928SBarry Smith 
12304b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
12314b9ad928SBarry Smith           PCASMCreateSubdomains2D(),
12324b9ad928SBarry Smith @*/
12337087cfbeSBarry Smith PetscErrorCode  PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
12344b9ad928SBarry Smith {
12357bb14e67SBarry Smith   PetscErrorCode ierr;
12364b9ad928SBarry Smith 
12374b9ad928SBarry Smith   PetscFunctionBegin;
12380700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12397bb14e67SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
12404b9ad928SBarry Smith   PetscFunctionReturn(0);
12414b9ad928SBarry Smith }
12424b9ad928SBarry Smith 
12434b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
12444b9ad928SBarry Smith /*MC
12453b09bd56SBarry Smith    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
12464b9ad928SBarry Smith            its own KSP object.
12474b9ad928SBarry Smith 
12484b9ad928SBarry Smith    Options Database Keys:
124949517cdeSBarry Smith +  -pc_asm_blocks <blks> - Sets total blocks
12504b9ad928SBarry Smith .  -pc_asm_overlap <ovl> - Sets overlap
1251f1ee410cSBarry Smith .  -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict
1252f1ee410cSBarry Smith -  -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive
12534b9ad928SBarry Smith 
12543b09bd56SBarry Smith      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
12553b09bd56SBarry Smith       will get a different convergence rate due to the default option of -pc_asm_type restrict. Use
12563b09bd56SBarry Smith       -pc_asm_type basic to use the standard ASM.
12573b09bd56SBarry Smith 
125895452b02SPatrick Sanan    Notes:
125995452b02SPatrick Sanan     Each processor can have one or more blocks, but a block cannot be shared by more
1260f1ee410cSBarry Smith      than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor.
12614b9ad928SBarry Smith 
12623b09bd56SBarry Smith      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1263d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
12644b9ad928SBarry Smith 
1265a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCASMGetSubKSP()
12664b9ad928SBarry Smith          and set the options directly on the resulting KSP object (you can access its PC
12674b9ad928SBarry Smith          with KSPGetPC())
12684b9ad928SBarry Smith 
12694b9ad928SBarry Smith    Level: beginner
12704b9ad928SBarry Smith 
1271c582cd25SBarry Smith     References:
127296a0c994SBarry Smith +   1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
127396a0c994SBarry Smith      Courant Institute, New York University Technical report
127496a0c994SBarry Smith -   1. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
127596a0c994SBarry Smith     Cambridge University Press.
1276c582cd25SBarry Smith 
12774b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1278f1ee410cSBarry Smith            PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType()
1279f1ee410cSBarry Smith            PCASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType
1280e09e08ccSBarry Smith 
12814b9ad928SBarry Smith M*/
12824b9ad928SBarry Smith 
12838cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc)
12844b9ad928SBarry Smith {
1285dfbe8321SBarry Smith   PetscErrorCode ierr;
12864b9ad928SBarry Smith   PC_ASM         *osm;
12874b9ad928SBarry Smith 
12884b9ad928SBarry Smith   PetscFunctionBegin;
1289b00a9115SJed Brown   ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr);
12902fa5cd67SKarl Rupp 
12914b9ad928SBarry Smith   osm->n                 = PETSC_DECIDE;
12924b9ad928SBarry Smith   osm->n_local           = 0;
12932b837212SDmitry Karpeev   osm->n_local_true      = PETSC_DECIDE;
12944b9ad928SBarry Smith   osm->overlap           = 1;
12954b9ad928SBarry Smith   osm->ksp               = 0;
12962b691e39SMatthew Knepley   osm->restriction       = 0;
12975b3c0d42Seaulisa   osm->lprolongation     = 0;
12981dd8081eSeaulisa   osm->lrestriction      = 0;
129992bb6962SLisandro Dalcin   osm->x                 = 0;
130092bb6962SLisandro Dalcin   osm->y                 = 0;
13014b9ad928SBarry Smith   osm->is                = 0;
13022b691e39SMatthew Knepley   osm->is_local          = 0;
13034b9ad928SBarry Smith   osm->mat               = 0;
13044b9ad928SBarry Smith   osm->pmat              = 0;
13054b9ad928SBarry Smith   osm->type              = PC_ASM_RESTRICT;
130612cd4985SMatthew G. Knepley   osm->loctype           = PC_COMPOSITE_ADDITIVE;
13074b9ad928SBarry Smith   osm->same_local_solves = PETSC_TRUE;
13086ed231c7SMatthew Knepley   osm->sort_indices      = PETSC_TRUE;
1309d709ea83SDmitry Karpeev   osm->dm_subdomains     = PETSC_FALSE;
131080ec0b7dSPatrick Sanan   osm->sub_mat_type      = NULL;
13114b9ad928SBarry Smith 
131292bb6962SLisandro Dalcin   pc->data                 = (void*)osm;
13134b9ad928SBarry Smith   pc->ops->apply           = PCApply_ASM;
13144b9ad928SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_ASM;
13154b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_ASM;
1316e91c6855SBarry Smith   pc->ops->reset           = PCReset_ASM;
13174b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_ASM;
13184b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_ASM;
13194b9ad928SBarry Smith   pc->ops->setuponblocks   = PCSetUpOnBlocks_ASM;
13204b9ad928SBarry Smith   pc->ops->view            = PCView_ASM;
13214b9ad928SBarry Smith   pc->ops->applyrichardson = 0;
13224b9ad928SBarry Smith 
1323bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr);
1324bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr);
1325bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);CHKERRQ(ierr);
1326bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);CHKERRQ(ierr);
1327c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);CHKERRQ(ierr);
132812cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);CHKERRQ(ierr);
132912cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);CHKERRQ(ierr);
1330bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);CHKERRQ(ierr);
1331bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);CHKERRQ(ierr);
133280ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);CHKERRQ(ierr);
133380ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);CHKERRQ(ierr);
13344b9ad928SBarry Smith   PetscFunctionReturn(0);
13354b9ad928SBarry Smith }
13364b9ad928SBarry Smith 
133792bb6962SLisandro Dalcin /*@C
133892bb6962SLisandro Dalcin    PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
133992bb6962SLisandro Dalcin    preconditioner for a any problem on a general grid.
134092bb6962SLisandro Dalcin 
134192bb6962SLisandro Dalcin    Collective
134292bb6962SLisandro Dalcin 
134392bb6962SLisandro Dalcin    Input Parameters:
134492bb6962SLisandro Dalcin +  A - The global matrix operator
134592bb6962SLisandro Dalcin -  n - the number of local blocks
134692bb6962SLisandro Dalcin 
134792bb6962SLisandro Dalcin    Output Parameters:
134892bb6962SLisandro Dalcin .  outis - the array of index sets defining the subdomains
134992bb6962SLisandro Dalcin 
135092bb6962SLisandro Dalcin    Level: advanced
135192bb6962SLisandro Dalcin 
13527d6bfa3bSBarry Smith    Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap
13537d6bfa3bSBarry Smith     from these if you use PCASMSetLocalSubdomains()
13547d6bfa3bSBarry Smith 
13557d6bfa3bSBarry Smith     In the Fortran version you must provide the array outis[] already allocated of length n.
13567d6bfa3bSBarry Smith 
135792bb6962SLisandro Dalcin .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains()
135892bb6962SLisandro Dalcin @*/
13597087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
136092bb6962SLisandro Dalcin {
136192bb6962SLisandro Dalcin   MatPartitioning mpart;
136292bb6962SLisandro Dalcin   const char      *prefix;
136392bb6962SLisandro Dalcin   PetscInt        i,j,rstart,rend,bs;
1364976e8c5aSLisandro Dalcin   PetscBool       hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
13650298fd71SBarry Smith   Mat             Ad     = NULL, adj;
136692bb6962SLisandro Dalcin   IS              ispart,isnumb,*is;
136792bb6962SLisandro Dalcin   PetscErrorCode  ierr;
136892bb6962SLisandro Dalcin 
136992bb6962SLisandro Dalcin   PetscFunctionBegin;
13700700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
137192bb6962SLisandro Dalcin   PetscValidPointer(outis,3);
1372c1235816SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
137392bb6962SLisandro Dalcin 
137492bb6962SLisandro Dalcin   /* Get prefix, row distribution, and block size */
137592bb6962SLisandro Dalcin   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
137692bb6962SLisandro Dalcin   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
137792bb6962SLisandro Dalcin   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
137865e19b50SBarry 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);
137965e19b50SBarry Smith 
138092bb6962SLisandro Dalcin   /* Get diagonal block from matrix if possible */
1381976e8c5aSLisandro Dalcin   ierr = MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr);
1382976e8c5aSLisandro Dalcin   if (hasop) {
138311bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
138492bb6962SLisandro Dalcin   }
138592bb6962SLisandro Dalcin   if (Ad) {
1386b9e7e5c1SBarry Smith     ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1387b9e7e5c1SBarry Smith     if (!isbaij) {ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
138892bb6962SLisandro Dalcin   }
138992bb6962SLisandro Dalcin   if (Ad && n > 1) {
1390ace3abfcSBarry Smith     PetscBool match,done;
139192bb6962SLisandro Dalcin     /* Try to setup a good matrix partitioning if available */
139292bb6962SLisandro Dalcin     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
139392bb6962SLisandro Dalcin     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
139492bb6962SLisandro Dalcin     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1395251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
139692bb6962SLisandro Dalcin     if (!match) {
1397251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
139892bb6962SLisandro Dalcin     }
139992bb6962SLisandro Dalcin     if (!match) { /* assume a "good" partitioner is available */
14001a83f524SJed Brown       PetscInt       na;
14011a83f524SJed Brown       const PetscInt *ia,*ja;
140292bb6962SLisandro Dalcin       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
140392bb6962SLisandro Dalcin       if (done) {
140492bb6962SLisandro Dalcin         /* Build adjacency matrix by hand. Unfortunately a call to
140592bb6962SLisandro Dalcin            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
140692bb6962SLisandro Dalcin            remove the block-aij structure and we cannot expect
140792bb6962SLisandro Dalcin            MatPartitioning to split vertices as we need */
14081a83f524SJed Brown         PetscInt       i,j,len,nnz,cnt,*iia=0,*jja=0;
14091a83f524SJed Brown         const PetscInt *row;
141092bb6962SLisandro Dalcin         nnz = 0;
141192bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* count number of nonzeros */
141292bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
141392bb6962SLisandro Dalcin           row = ja + ia[i];
141492bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
141592bb6962SLisandro Dalcin             if (row[j] == i) { /* don't count diagonal */
141692bb6962SLisandro Dalcin               len--; break;
141792bb6962SLisandro Dalcin             }
141892bb6962SLisandro Dalcin           }
141992bb6962SLisandro Dalcin           nnz += len;
142092bb6962SLisandro Dalcin         }
1421854ce69bSBarry Smith         ierr   = PetscMalloc1(na+1,&iia);CHKERRQ(ierr);
1422854ce69bSBarry Smith         ierr   = PetscMalloc1(nnz,&jja);CHKERRQ(ierr);
142392bb6962SLisandro Dalcin         nnz    = 0;
142492bb6962SLisandro Dalcin         iia[0] = 0;
142592bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* fill adjacency */
142692bb6962SLisandro Dalcin           cnt = 0;
142792bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
142892bb6962SLisandro Dalcin           row = ja + ia[i];
142992bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
143092bb6962SLisandro Dalcin             if (row[j] != i) { /* if not diagonal */
143192bb6962SLisandro Dalcin               jja[nnz+cnt++] = row[j];
143292bb6962SLisandro Dalcin             }
143392bb6962SLisandro Dalcin           }
143492bb6962SLisandro Dalcin           nnz     += cnt;
143592bb6962SLisandro Dalcin           iia[i+1] = nnz;
143692bb6962SLisandro Dalcin         }
143792bb6962SLisandro Dalcin         /* Partitioning of the adjacency matrix */
14380298fd71SBarry Smith         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr);
143992bb6962SLisandro Dalcin         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
144092bb6962SLisandro Dalcin         ierr      = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
144192bb6962SLisandro Dalcin         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
144292bb6962SLisandro Dalcin         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
1443fcfd50ebSBarry Smith         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
144492bb6962SLisandro Dalcin         foundpart = PETSC_TRUE;
144592bb6962SLisandro Dalcin       }
144692bb6962SLisandro Dalcin       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
144792bb6962SLisandro Dalcin     }
1448fcfd50ebSBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
144992bb6962SLisandro Dalcin   }
145092bb6962SLisandro Dalcin 
1451785e854fSJed Brown   ierr   = PetscMalloc1(n,&is);CHKERRQ(ierr);
145292bb6962SLisandro Dalcin   *outis = is;
145392bb6962SLisandro Dalcin 
145492bb6962SLisandro Dalcin   if (!foundpart) {
145592bb6962SLisandro Dalcin 
145692bb6962SLisandro Dalcin     /* Partitioning by contiguous chunks of rows */
145792bb6962SLisandro Dalcin 
145892bb6962SLisandro Dalcin     PetscInt mbs   = (rend-rstart)/bs;
145992bb6962SLisandro Dalcin     PetscInt start = rstart;
146092bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
146192bb6962SLisandro Dalcin       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
146292bb6962SLisandro Dalcin       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
146392bb6962SLisandro Dalcin       start += count;
146492bb6962SLisandro Dalcin     }
146592bb6962SLisandro Dalcin 
146692bb6962SLisandro Dalcin   } else {
146792bb6962SLisandro Dalcin 
146892bb6962SLisandro Dalcin     /* Partitioning by adjacency of diagonal block  */
146992bb6962SLisandro Dalcin 
147092bb6962SLisandro Dalcin     const PetscInt *numbering;
147192bb6962SLisandro Dalcin     PetscInt       *count,nidx,*indices,*newidx,start=0;
147292bb6962SLisandro Dalcin     /* Get node count in each partition */
1473785e854fSJed Brown     ierr = PetscMalloc1(n,&count);CHKERRQ(ierr);
147492bb6962SLisandro Dalcin     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
147592bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
147692bb6962SLisandro Dalcin       for (i=0; i<n; i++) count[i] *= bs;
147792bb6962SLisandro Dalcin     }
147892bb6962SLisandro Dalcin     /* Build indices from node numbering */
147992bb6962SLisandro Dalcin     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1480785e854fSJed Brown     ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr);
148192bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
148292bb6962SLisandro Dalcin     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
148392bb6962SLisandro Dalcin     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
148492bb6962SLisandro Dalcin     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
148592bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1486785e854fSJed Brown       ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr);
14872fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
14882fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
14892fa5cd67SKarl Rupp       }
149092bb6962SLisandro Dalcin       ierr    = PetscFree(indices);CHKERRQ(ierr);
149192bb6962SLisandro Dalcin       nidx   *= bs;
149292bb6962SLisandro Dalcin       indices = newidx;
149392bb6962SLisandro Dalcin     }
149492bb6962SLisandro Dalcin     /* Shift to get global indices */
149592bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] += rstart;
149692bb6962SLisandro Dalcin 
149792bb6962SLisandro Dalcin     /* Build the index sets for each block */
149892bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
149970b3c8c7SBarry Smith       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
150092bb6962SLisandro Dalcin       ierr   = ISSort(is[i]);CHKERRQ(ierr);
150192bb6962SLisandro Dalcin       start += count[i];
150292bb6962SLisandro Dalcin     }
150392bb6962SLisandro Dalcin 
15043bf036e2SBarry Smith     ierr = PetscFree(count);CHKERRQ(ierr);
15053bf036e2SBarry Smith     ierr = PetscFree(indices);CHKERRQ(ierr);
1506fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1507fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
150892bb6962SLisandro Dalcin 
150992bb6962SLisandro Dalcin   }
151092bb6962SLisandro Dalcin   PetscFunctionReturn(0);
151192bb6962SLisandro Dalcin }
151292bb6962SLisandro Dalcin 
151392bb6962SLisandro Dalcin /*@C
151492bb6962SLisandro Dalcin    PCASMDestroySubdomains - Destroys the index sets created with
151592bb6962SLisandro Dalcin    PCASMCreateSubdomains(). Should be called after setting subdomains
151692bb6962SLisandro Dalcin    with PCASMSetLocalSubdomains().
151792bb6962SLisandro Dalcin 
151892bb6962SLisandro Dalcin    Collective
151992bb6962SLisandro Dalcin 
152092bb6962SLisandro Dalcin    Input Parameters:
152192bb6962SLisandro Dalcin +  n - the number of index sets
15222b691e39SMatthew Knepley .  is - the array of index sets
15230298fd71SBarry Smith -  is_local - the array of local index sets, can be NULL
152492bb6962SLisandro Dalcin 
152592bb6962SLisandro Dalcin    Level: advanced
152692bb6962SLisandro Dalcin 
152792bb6962SLisandro Dalcin .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains()
152892bb6962SLisandro Dalcin @*/
15297087cfbeSBarry Smith PetscErrorCode  PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
153092bb6962SLisandro Dalcin {
153192bb6962SLisandro Dalcin   PetscInt       i;
153292bb6962SLisandro Dalcin   PetscErrorCode ierr;
15335fd66863SKarl Rupp 
153492bb6962SLisandro Dalcin   PetscFunctionBegin;
1535a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1536a35b7b57SDmitry Karpeev   if (is) {
153792bb6962SLisandro Dalcin     PetscValidPointer(is,2);
1538fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); }
153992bb6962SLisandro Dalcin     ierr = PetscFree(is);CHKERRQ(ierr);
1540a35b7b57SDmitry Karpeev   }
15412b691e39SMatthew Knepley   if (is_local) {
15422b691e39SMatthew Knepley     PetscValidPointer(is_local,3);
1543fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); }
15442b691e39SMatthew Knepley     ierr = PetscFree(is_local);CHKERRQ(ierr);
15452b691e39SMatthew Knepley   }
154692bb6962SLisandro Dalcin   PetscFunctionReturn(0);
154792bb6962SLisandro Dalcin }
154892bb6962SLisandro Dalcin 
15494b9ad928SBarry Smith /*@
15504b9ad928SBarry Smith    PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
15514b9ad928SBarry Smith    preconditioner for a two-dimensional problem on a regular grid.
15524b9ad928SBarry Smith 
15534b9ad928SBarry Smith    Not Collective
15544b9ad928SBarry Smith 
15554b9ad928SBarry Smith    Input Parameters:
15564b9ad928SBarry Smith +  m, n - the number of mesh points in the x and y directions
15574b9ad928SBarry Smith .  M, N - the number of subdomains in the x and y directions
15584b9ad928SBarry Smith .  dof - degrees of freedom per node
15594b9ad928SBarry Smith -  overlap - overlap in mesh lines
15604b9ad928SBarry Smith 
15614b9ad928SBarry Smith    Output Parameters:
15624b9ad928SBarry Smith +  Nsub - the number of subdomains created
15633d03bb48SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
15643d03bb48SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
15654b9ad928SBarry Smith 
15664b9ad928SBarry Smith    Note:
15674b9ad928SBarry Smith    Presently PCAMSCreateSubdomains2d() is valid only for sequential
15684b9ad928SBarry Smith    preconditioners.  More general related routines are
15694b9ad928SBarry Smith    PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().
15704b9ad928SBarry Smith 
15714b9ad928SBarry Smith    Level: advanced
15724b9ad928SBarry Smith 
15734b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
15744b9ad928SBarry Smith           PCASMSetOverlap()
15754b9ad928SBarry Smith @*/
15767087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local)
15774b9ad928SBarry Smith {
15783d03bb48SJed Brown   PetscInt       i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer;
15796849ba73SBarry Smith   PetscErrorCode ierr;
158013f74950SBarry Smith   PetscInt       nidx,*idx,loc,ii,jj,count;
15814b9ad928SBarry Smith 
15824b9ad928SBarry Smith   PetscFunctionBegin;
1583e32f2f54SBarry Smith   if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," ");
15844b9ad928SBarry Smith 
15854b9ad928SBarry Smith   *Nsub     = N*M;
1586854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is);CHKERRQ(ierr);
1587854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is_local);CHKERRQ(ierr);
15884b9ad928SBarry Smith   ystart    = 0;
15893d03bb48SJed Brown   loc_outer = 0;
15904b9ad928SBarry Smith   for (i=0; i<N; i++) {
15914b9ad928SBarry Smith     height = n/N + ((n % N) > i); /* height of subdomain */
1592e32f2f54SBarry Smith     if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n");
15934b9ad928SBarry Smith     yleft  = ystart - overlap; if (yleft < 0) yleft = 0;
15944b9ad928SBarry Smith     yright = ystart + height + overlap; if (yright > n) yright = n;
15954b9ad928SBarry Smith     xstart = 0;
15964b9ad928SBarry Smith     for (j=0; j<M; j++) {
15974b9ad928SBarry Smith       width = m/M + ((m % M) > j); /* width of subdomain */
1598e32f2f54SBarry Smith       if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m");
15994b9ad928SBarry Smith       xleft  = xstart - overlap; if (xleft < 0) xleft = 0;
16004b9ad928SBarry Smith       xright = xstart + width + overlap; if (xright > m) xright = m;
16014b9ad928SBarry Smith       nidx   = (xright - xleft)*(yright - yleft);
1602785e854fSJed Brown       ierr   = PetscMalloc1(nidx,&idx);CHKERRQ(ierr);
16034b9ad928SBarry Smith       loc    = 0;
16044b9ad928SBarry Smith       for (ii=yleft; ii<yright; ii++) {
16054b9ad928SBarry Smith         count = m*ii + xleft;
16062fa5cd67SKarl Rupp         for (jj=xleft; jj<xright; jj++) idx[loc++] = count++;
16074b9ad928SBarry Smith       }
160870b3c8c7SBarry Smith       ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr);
16093d03bb48SJed Brown       if (overlap == 0) {
16103d03bb48SJed Brown         ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr);
16112fa5cd67SKarl Rupp 
16123d03bb48SJed Brown         (*is_local)[loc_outer] = (*is)[loc_outer];
16133d03bb48SJed Brown       } else {
16143d03bb48SJed Brown         for (loc=0,ii=ystart; ii<ystart+height; ii++) {
16153d03bb48SJed Brown           for (jj=xstart; jj<xstart+width; jj++) {
16163d03bb48SJed Brown             idx[loc++] = m*ii + jj;
16173d03bb48SJed Brown           }
16183d03bb48SJed Brown         }
161970b3c8c7SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr);
16203d03bb48SJed Brown       }
16214b9ad928SBarry Smith       ierr    = PetscFree(idx);CHKERRQ(ierr);
16224b9ad928SBarry Smith       xstart += width;
16233d03bb48SJed Brown       loc_outer++;
16244b9ad928SBarry Smith     }
16254b9ad928SBarry Smith     ystart += height;
16264b9ad928SBarry Smith   }
16274b9ad928SBarry Smith   for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); }
16284b9ad928SBarry Smith   PetscFunctionReturn(0);
16294b9ad928SBarry Smith }
16304b9ad928SBarry Smith 
16314b9ad928SBarry Smith /*@C
16324b9ad928SBarry Smith     PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
16334b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
16344b9ad928SBarry Smith 
1635ad4df100SBarry Smith     Not Collective
16364b9ad928SBarry Smith 
16374b9ad928SBarry Smith     Input Parameter:
16384b9ad928SBarry Smith .   pc - the preconditioner context
16394b9ad928SBarry Smith 
16404b9ad928SBarry Smith     Output Parameters:
16414b9ad928SBarry Smith +   n - the number of subdomains for this processor (default value = 1)
16422b691e39SMatthew Knepley .   is - the index sets that define the subdomains for this processor
16430298fd71SBarry Smith -   is_local - the index sets that define the local part of the subdomains for this processor (can be NULL)
16444b9ad928SBarry Smith 
16454b9ad928SBarry Smith 
16464b9ad928SBarry Smith     Notes:
16474b9ad928SBarry Smith     The IS numbering is in the parallel, global numbering of the vector.
16484b9ad928SBarry Smith 
16494b9ad928SBarry Smith     Level: advanced
16504b9ad928SBarry Smith 
16514b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
16524b9ad928SBarry Smith           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices()
16534b9ad928SBarry Smith @*/
16547087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
16554b9ad928SBarry Smith {
16562a808120SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
165792bb6962SLisandro Dalcin   PetscErrorCode ierr;
1658ace3abfcSBarry Smith   PetscBool      match;
16594b9ad928SBarry Smith 
16604b9ad928SBarry Smith   PetscFunctionBegin;
16610700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
16624482741eSBarry Smith   PetscValidIntPointer(n,2);
166392bb6962SLisandro Dalcin   if (is) PetscValidPointer(is,3);
1664251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
16652a808120SBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM");
16664b9ad928SBarry Smith   if (n) *n = osm->n_local_true;
16674b9ad928SBarry Smith   if (is) *is = osm->is;
16682b691e39SMatthew Knepley   if (is_local) *is_local = osm->is_local;
16694b9ad928SBarry Smith   PetscFunctionReturn(0);
16704b9ad928SBarry Smith }
16714b9ad928SBarry Smith 
16724b9ad928SBarry Smith /*@C
16734b9ad928SBarry Smith     PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
16744b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
16754b9ad928SBarry Smith 
1676ad4df100SBarry Smith     Not Collective
16774b9ad928SBarry Smith 
16784b9ad928SBarry Smith     Input Parameter:
16794b9ad928SBarry Smith .   pc - the preconditioner context
16804b9ad928SBarry Smith 
16814b9ad928SBarry Smith     Output Parameters:
16824b9ad928SBarry Smith +   n - the number of matrices for this processor (default value = 1)
16834b9ad928SBarry Smith -   mat - the matrices
16844b9ad928SBarry Smith 
16854b9ad928SBarry Smith     Level: advanced
16864b9ad928SBarry Smith 
168795452b02SPatrick Sanan     Notes:
168895452b02SPatrick Sanan     Call after PCSetUp() (or KSPSetUp()) but before PCApply() (or KSPApply()) and before PCSetUpOnBlocks())
1689cf739d55SBarry Smith 
1690cf739d55SBarry Smith            Usually one would use PCSetModifySubmatrices() to change the submatrices in building the preconditioner.
1691cf739d55SBarry Smith 
16924b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1693cf739d55SBarry Smith           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubmatrices()
16944b9ad928SBarry Smith @*/
16957087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
16964b9ad928SBarry Smith {
16974b9ad928SBarry Smith   PC_ASM         *osm;
169892bb6962SLisandro Dalcin   PetscErrorCode ierr;
1699ace3abfcSBarry Smith   PetscBool      match;
17004b9ad928SBarry Smith 
17014b9ad928SBarry Smith   PetscFunctionBegin;
17020700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
170392bb6962SLisandro Dalcin   PetscValidIntPointer(n,2);
170492bb6962SLisandro Dalcin   if (mat) PetscValidPointer(mat,3);
1705ce94432eSBarry Smith   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1706251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
170792bb6962SLisandro Dalcin   if (!match) {
170892bb6962SLisandro Dalcin     if (n) *n = 0;
17090298fd71SBarry Smith     if (mat) *mat = NULL;
171092bb6962SLisandro Dalcin   } else {
17114b9ad928SBarry Smith     osm = (PC_ASM*)pc->data;
17124b9ad928SBarry Smith     if (n) *n = osm->n_local_true;
17134b9ad928SBarry Smith     if (mat) *mat = osm->pmat;
171492bb6962SLisandro Dalcin   }
17154b9ad928SBarry Smith   PetscFunctionReturn(0);
17164b9ad928SBarry Smith }
1717d709ea83SDmitry Karpeev 
1718d709ea83SDmitry Karpeev /*@
1719d709ea83SDmitry Karpeev     PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1720f1ee410cSBarry Smith 
1721d709ea83SDmitry Karpeev     Logically Collective
1722d709ea83SDmitry Karpeev 
1723d709ea83SDmitry Karpeev     Input Parameter:
1724d709ea83SDmitry Karpeev +   pc  - the preconditioner
1725d709ea83SDmitry Karpeev -   flg - boolean indicating whether to use subdomains defined by the DM
1726d709ea83SDmitry Karpeev 
1727d709ea83SDmitry Karpeev     Options Database Key:
1728d709ea83SDmitry Karpeev .   -pc_asm_dm_subdomains
1729d709ea83SDmitry Karpeev 
1730d709ea83SDmitry Karpeev     Level: intermediate
1731d709ea83SDmitry Karpeev 
1732d709ea83SDmitry Karpeev     Notes:
1733d709ea83SDmitry Karpeev     PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(),
1734d709ea83SDmitry Karpeev     so setting either of the first two effectively turns the latter off.
1735d709ea83SDmitry Karpeev 
1736d709ea83SDmitry Karpeev .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1737d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1738d709ea83SDmitry Karpeev @*/
1739d709ea83SDmitry Karpeev PetscErrorCode  PCASMSetDMSubdomains(PC pc,PetscBool flg)
1740d709ea83SDmitry Karpeev {
1741d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1742d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1743d709ea83SDmitry Karpeev   PetscBool      match;
1744d709ea83SDmitry Karpeev 
1745d709ea83SDmitry Karpeev   PetscFunctionBegin;
1746d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1747d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
1748d709ea83SDmitry Karpeev   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1749d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1750d709ea83SDmitry Karpeev   if (match) {
1751d709ea83SDmitry Karpeev     osm->dm_subdomains = flg;
1752d709ea83SDmitry Karpeev   }
1753d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1754d709ea83SDmitry Karpeev }
1755d709ea83SDmitry Karpeev 
1756d709ea83SDmitry Karpeev /*@
1757d709ea83SDmitry Karpeev     PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1758d709ea83SDmitry Karpeev     Not Collective
1759d709ea83SDmitry Karpeev 
1760d709ea83SDmitry Karpeev     Input Parameter:
1761d709ea83SDmitry Karpeev .   pc  - the preconditioner
1762d709ea83SDmitry Karpeev 
1763d709ea83SDmitry Karpeev     Output Parameter:
1764d709ea83SDmitry Karpeev .   flg - boolean indicating whether to use subdomains defined by the DM
1765d709ea83SDmitry Karpeev 
1766d709ea83SDmitry Karpeev     Level: intermediate
1767d709ea83SDmitry Karpeev 
1768d709ea83SDmitry Karpeev .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1769d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1770d709ea83SDmitry Karpeev @*/
1771d709ea83SDmitry Karpeev PetscErrorCode  PCASMGetDMSubdomains(PC pc,PetscBool* flg)
1772d709ea83SDmitry Karpeev {
1773d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1774d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1775d709ea83SDmitry Karpeev   PetscBool      match;
1776d709ea83SDmitry Karpeev 
1777d709ea83SDmitry Karpeev   PetscFunctionBegin;
1778d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1779*534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg,2);
1780d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1781d709ea83SDmitry Karpeev   if (match) {
1782d709ea83SDmitry Karpeev     if (flg) *flg = osm->dm_subdomains;
1783d709ea83SDmitry Karpeev   }
1784d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1785d709ea83SDmitry Karpeev }
178680ec0b7dSPatrick Sanan 
178780ec0b7dSPatrick Sanan /*@
178880ec0b7dSPatrick Sanan      PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string.
178980ec0b7dSPatrick Sanan 
179080ec0b7dSPatrick Sanan    Not Collective
179180ec0b7dSPatrick Sanan 
179280ec0b7dSPatrick Sanan    Input Parameter:
179380ec0b7dSPatrick Sanan .  pc - the PC
179480ec0b7dSPatrick Sanan 
179580ec0b7dSPatrick Sanan    Output Parameter:
1796f1ee410cSBarry Smith .  -pc_asm_sub_mat_type - name of matrix type
179780ec0b7dSPatrick Sanan 
179880ec0b7dSPatrick Sanan    Level: advanced
179980ec0b7dSPatrick Sanan 
180080ec0b7dSPatrick Sanan .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
180180ec0b7dSPatrick Sanan @*/
180280ec0b7dSPatrick Sanan PetscErrorCode  PCASMGetSubMatType(PC pc,MatType *sub_mat_type){
180380ec0b7dSPatrick Sanan   PetscErrorCode ierr;
180480ec0b7dSPatrick Sanan 
180580ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));CHKERRQ(ierr);
180680ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
180780ec0b7dSPatrick Sanan }
180880ec0b7dSPatrick Sanan 
180980ec0b7dSPatrick Sanan /*@
181080ec0b7dSPatrick Sanan      PCASMSetSubMatType - Set the type of matrix used for ASM subsolves
181180ec0b7dSPatrick Sanan 
181280ec0b7dSPatrick Sanan    Collective on Mat
181380ec0b7dSPatrick Sanan 
181480ec0b7dSPatrick Sanan    Input Parameters:
181580ec0b7dSPatrick Sanan +  pc             - the PC object
181680ec0b7dSPatrick Sanan -  sub_mat_type   - matrix type
181780ec0b7dSPatrick Sanan 
181880ec0b7dSPatrick Sanan    Options Database Key:
181980ec0b7dSPatrick 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.
182080ec0b7dSPatrick Sanan 
182180ec0b7dSPatrick Sanan    Notes:
182280ec0b7dSPatrick Sanan    See "${PETSC_DIR}/include/petscmat.h" for available types
182380ec0b7dSPatrick Sanan 
182480ec0b7dSPatrick Sanan   Level: advanced
182580ec0b7dSPatrick Sanan 
182680ec0b7dSPatrick Sanan .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
182780ec0b7dSPatrick Sanan @*/
182880ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type)
182980ec0b7dSPatrick Sanan {
183080ec0b7dSPatrick Sanan   PetscErrorCode ierr;
183180ec0b7dSPatrick Sanan 
183280ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));CHKERRQ(ierr);
183380ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
183480ec0b7dSPatrick Sanan }
1835