xref: /petsc/src/ksp/pc/impls/asm/asm.c (revision 8d76b567012869fc44a72bfbee2ffc0c796c6d4b)
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) {
62*8d76b567SPatrick Sanan         ierr = PetscViewerASCIIPrintf(viewer,"  Local solver is the same for all blocks, as in the following KSP and PC objects on rank 0:\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;
18257501b6eSBarry Smith   PetscBool      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 
1962b837212SDmitry Karpeev     /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */
1972b837212SDmitry Karpeev     if (osm->n_local_true == PETSC_DECIDE) {
19869ca1f37SDmitry Karpeev       /* no subdomains given */
19965db9045SDmitry Karpeev       /* try pc->dm first, if allowed */
200d709ea83SDmitry Karpeev       if (osm->dm_subdomains && pc->dm) {
201feb221f8SDmitry Karpeev         PetscInt  num_domains, d;
202feb221f8SDmitry Karpeev         char      **domain_names;
2038d4ac253SDmitry Karpeev         IS        *inner_domain_is, *outer_domain_is;
2048d4ac253SDmitry Karpeev         ierr = DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm);CHKERRQ(ierr);
205704f0395SPatrick Sanan         osm->overlap = -1; /* We do not want to increase the overlap of the IS.
206704f0395SPatrick Sanan                               A future improvement of this code might allow one to use
207704f0395SPatrick Sanan                               DM-defined subdomains and also increase the overlap,
208704f0395SPatrick Sanan                               but that is not currently supported */
20969ca1f37SDmitry Karpeev         if (num_domains) {
2108d4ac253SDmitry Karpeev           ierr = PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is);CHKERRQ(ierr);
21169ca1f37SDmitry Karpeev         }
212feb221f8SDmitry Karpeev         for (d = 0; d < num_domains; ++d) {
213a35b7b57SDmitry Karpeev           if (domain_names)    {ierr = PetscFree(domain_names[d]);CHKERRQ(ierr);}
214a35b7b57SDmitry Karpeev           if (inner_domain_is) {ierr = ISDestroy(&inner_domain_is[d]);CHKERRQ(ierr);}
215a35b7b57SDmitry Karpeev           if (outer_domain_is) {ierr = ISDestroy(&outer_domain_is[d]);CHKERRQ(ierr);}
216feb221f8SDmitry Karpeev         }
217feb221f8SDmitry Karpeev         ierr = PetscFree(domain_names);CHKERRQ(ierr);
2188d4ac253SDmitry Karpeev         ierr = PetscFree(inner_domain_is);CHKERRQ(ierr);
2198d4ac253SDmitry Karpeev         ierr = PetscFree(outer_domain_is);CHKERRQ(ierr);
220feb221f8SDmitry Karpeev       }
2212b837212SDmitry Karpeev       if (osm->n_local_true == PETSC_DECIDE) {
22269ca1f37SDmitry Karpeev         /* still no subdomains; use one subdomain per processor */
2232b837212SDmitry Karpeev         osm->n_local_true = 1;
224feb221f8SDmitry Karpeev       }
2252b837212SDmitry Karpeev     }
2262b837212SDmitry Karpeev     { /* determine the global and max number of subdomains */
2276ac3741eSJed Brown       struct {PetscInt max,sum;} inwork,outwork;
228c10200c1SHong Zhang       PetscMPIInt size;
229c10200c1SHong Zhang 
2306ac3741eSJed Brown       inwork.max   = osm->n_local_true;
2316ac3741eSJed Brown       inwork.sum   = osm->n_local_true;
232367daffbSBarry Smith       ierr         = MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
2336ac3741eSJed Brown       osm->n_local = outwork.max;
2346ac3741eSJed Brown       osm->n       = outwork.sum;
235c10200c1SHong Zhang 
236c10200c1SHong Zhang       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
237c10200c1SHong Zhang       if (outwork.max == 1 && outwork.sum == size) {
2387dae84e0SHong Zhang         /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */
239c10200c1SHong Zhang         ierr = MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE);CHKERRQ(ierr);
240c10200c1SHong Zhang       }
2414b9ad928SBarry Smith     }
24278904715SLisandro Dalcin     if (!osm->is) { /* create the index sets */
24378904715SLisandro Dalcin       ierr = PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);CHKERRQ(ierr);
2444b9ad928SBarry Smith     }
245f5234e35SJed Brown     if (osm->n_local_true > 1 && !osm->is_local) {
246785e854fSJed Brown       ierr = PetscMalloc1(osm->n_local_true,&osm->is_local);CHKERRQ(ierr);
247f5234e35SJed Brown       for (i=0; i<osm->n_local_true; i++) {
248f5234e35SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
249f5234e35SJed Brown           ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr);
250f5234e35SJed Brown           ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr);
251f5234e35SJed Brown         } else {
252f5234e35SJed Brown           ierr             = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr);
253f5234e35SJed Brown           osm->is_local[i] = osm->is[i];
254f5234e35SJed Brown         }
255f5234e35SJed Brown       }
256f5234e35SJed Brown     }
25792bb6962SLisandro Dalcin     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
25890d69ab7SBarry Smith     flg  = PETSC_FALSE;
259c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,prefix,"-pc_asm_print_subdomains",&flg,NULL);CHKERRQ(ierr);
26092bb6962SLisandro Dalcin     if (flg) { ierr = PCASMPrintSubdomains(pc);CHKERRQ(ierr); }
2614b9ad928SBarry Smith 
2623d03bb48SJed Brown     if (osm->overlap > 0) {
2634b9ad928SBarry Smith       /* Extend the "overlapping" regions by a number of steps */
26492bb6962SLisandro Dalcin       ierr = MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);CHKERRQ(ierr);
2653d03bb48SJed Brown     }
2666ed231c7SMatthew Knepley     if (osm->sort_indices) {
26792bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2684b9ad928SBarry Smith         ierr = ISSort(osm->is[i]);CHKERRQ(ierr);
2692b691e39SMatthew Knepley         if (osm->is_local) {
2702b691e39SMatthew Knepley           ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr);
2712b691e39SMatthew Knepley         }
2724b9ad928SBarry Smith       }
2736ed231c7SMatthew Knepley     }
2744b9ad928SBarry Smith 
275f6991133SBarry Smith     if (!osm->ksp) {
27678904715SLisandro Dalcin       /* Create the local solvers */
277785e854fSJed Brown       ierr = PetscMalloc1(osm->n_local_true,&osm->ksp);CHKERRQ(ierr);
278feb221f8SDmitry Karpeev       if (domain_dm) {
279feb221f8SDmitry Karpeev         ierr = PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n");CHKERRQ(ierr);
280feb221f8SDmitry Karpeev       }
28192bb6962SLisandro Dalcin       for (i=0; i<osm->n_local_true; i++) {
2824b9ad928SBarry Smith         ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
283422a814eSBarry Smith         ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr);
2843bb1ff40SBarry Smith         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr);
28592bb6962SLisandro Dalcin         ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
2864b9ad928SBarry Smith         ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
2874b9ad928SBarry Smith         ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
2884b9ad928SBarry Smith         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
2894b9ad928SBarry Smith         ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
2904b9ad928SBarry Smith         ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
291feb221f8SDmitry Karpeev         if (domain_dm) {
292feb221f8SDmitry Karpeev           ierr = KSPSetDM(ksp, domain_dm[i]);CHKERRQ(ierr);
293feb221f8SDmitry Karpeev           ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr);
294feb221f8SDmitry Karpeev           ierr = DMDestroy(&domain_dm[i]);CHKERRQ(ierr);
295feb221f8SDmitry Karpeev         }
2964b9ad928SBarry Smith         osm->ksp[i] = ksp;
2974b9ad928SBarry Smith       }
298feb221f8SDmitry Karpeev       if (domain_dm) {
299feb221f8SDmitry Karpeev         ierr = PetscFree(domain_dm);CHKERRQ(ierr);
300feb221f8SDmitry Karpeev       }
301f6991133SBarry Smith     }
3021dd8081eSeaulisa 
303fb745f2cSMatthew G. Knepley     ierr = ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis);CHKERRQ(ierr);
304fb745f2cSMatthew G. Knepley     ierr = ISSortRemoveDups(osm->lis);CHKERRQ(ierr);
305fb745f2cSMatthew G. Knepley     ierr = ISGetLocalSize(osm->lis, &m);CHKERRQ(ierr);
306fb745f2cSMatthew G. Knepley     ierr = VecCreateSeq(PETSC_COMM_SELF, m, &osm->lx);CHKERRQ(ierr);
307fb745f2cSMatthew G. Knepley     ierr = VecDuplicate(osm->lx, &osm->ly);CHKERRQ(ierr);
3081dd8081eSeaulisa 
3094b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
3104b9ad928SBarry Smith   } else {
3114b9ad928SBarry Smith     /*
3124b9ad928SBarry Smith        Destroy the blocks from the previous iteration
3134b9ad928SBarry Smith     */
314e09e08ccSBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
3154b9ad928SBarry Smith       ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
3164b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
3174b9ad928SBarry Smith     }
3184b9ad928SBarry Smith   }
3194b9ad928SBarry Smith 
32092bb6962SLisandro Dalcin   /*
32192bb6962SLisandro Dalcin      Extract out the submatrices
32292bb6962SLisandro Dalcin   */
3237dae84e0SHong Zhang   ierr = MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr);
32492bb6962SLisandro Dalcin   if (scall == MAT_INITIAL_MATRIX) {
32592bb6962SLisandro Dalcin     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
32692bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
3273bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr);
32878904715SLisandro Dalcin       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
32992bb6962SLisandro Dalcin     }
33092bb6962SLisandro Dalcin   }
33180ec0b7dSPatrick Sanan 
33280ec0b7dSPatrick Sanan   /* Convert the types of the submatrices (if needbe) */
33380ec0b7dSPatrick Sanan   if (osm->sub_mat_type) {
33480ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; i++) {
33580ec0b7dSPatrick Sanan       ierr = MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]));CHKERRQ(ierr);
33680ec0b7dSPatrick Sanan     }
33780ec0b7dSPatrick Sanan   }
33880ec0b7dSPatrick Sanan 
33980ec0b7dSPatrick Sanan   if(!pc->setupcalled){
34080ec0b7dSPatrick Sanan     /* Create the local work vectors (from the local matrices) and scatter contexts */
34180ec0b7dSPatrick Sanan     ierr = MatCreateVecs(pc->pmat,&vec,0);CHKERRQ(ierr);
3425b3c0d42Seaulisa 
3431dd8081eSeaulisa     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()");
3441dd8081eSeaulisa     if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) {
3451dd8081eSeaulisa       ierr = PetscMalloc1(osm->n_local_true,&osm->lprolongation);CHKERRQ(ierr);
3461dd8081eSeaulisa     }
3471dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->lrestriction);CHKERRQ(ierr);
3481dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->x);CHKERRQ(ierr);
3491dd8081eSeaulisa     ierr = PetscMalloc1(osm->n_local_true,&osm->y);CHKERRQ(ierr);
350347574c9Seaulisa 
351347574c9Seaulisa     ierr = ISGetLocalSize(osm->lis,&m);CHKERRQ(ierr);
352347574c9Seaulisa     ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
3539448b7f1SJunchao Zhang     ierr = VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->restriction);CHKERRQ(ierr);
354347574c9Seaulisa     ierr = ISDestroy(&isl);CHKERRQ(ierr);
355347574c9Seaulisa 
356347574c9Seaulisa 
35780ec0b7dSPatrick Sanan     for (i=0; i<osm->n_local_true; ++i) {
3585b3c0d42Seaulisa       ISLocalToGlobalMapping ltog;
3595b3c0d42Seaulisa       IS                     isll;
3605b3c0d42Seaulisa       const PetscInt         *idx_is;
3615b3c0d42Seaulisa       PetscInt               *idx_lis,nout;
3625b3c0d42Seaulisa 
36355817e79SBarry Smith       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
36455817e79SBarry Smith       ierr = MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);CHKERRQ(ierr);
36555817e79SBarry Smith       ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr);
36655817e79SBarry Smith 
367b0de9f23SBarry Smith       /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */
3685b3c0d42Seaulisa       ierr = ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);CHKERRQ(ierr);
3695b3c0d42Seaulisa       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
3705b3c0d42Seaulisa       ierr = ISGetIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
3715b3c0d42Seaulisa       ierr = PetscMalloc1(m,&idx_lis);CHKERRQ(ierr);
3725b3c0d42Seaulisa       ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis);CHKERRQ(ierr);
3735b3c0d42Seaulisa       if (nout != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis");
3745b3c0d42Seaulisa       ierr = ISRestoreIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
3755b3c0d42Seaulisa       ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
376d8b3b5e3Seaulisa       ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
3775b3c0d42Seaulisa       ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
3789448b7f1SJunchao Zhang       ierr = VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lrestriction[i]);CHKERRQ(ierr);
3795b3c0d42Seaulisa       ierr = ISDestroy(&isll);CHKERRQ(ierr);
3805b3c0d42Seaulisa       ierr = ISDestroy(&isl);CHKERRQ(ierr);
381b0de9f23SBarry Smith       if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overalapping is_local[i] entries */
382d8b3b5e3Seaulisa         ISLocalToGlobalMapping ltog;
383d8b3b5e3Seaulisa         IS                     isll,isll_local;
384d8b3b5e3Seaulisa         const PetscInt         *idx_local;
385d8b3b5e3Seaulisa         PetscInt               *idx1, *idx2, nout;
386d8b3b5e3Seaulisa 
387d8b3b5e3Seaulisa         ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr);
388d8b3b5e3Seaulisa         ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
389d8b3b5e3Seaulisa 
390d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],&ltog);CHKERRQ(ierr);
391d8b3b5e3Seaulisa         ierr = PetscMalloc1(m_local,&idx1);CHKERRQ(ierr);
392d8b3b5e3Seaulisa         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1);CHKERRQ(ierr);
393d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
394d8b3b5e3Seaulisa         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is");
395d8b3b5e3Seaulisa         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
396d8b3b5e3Seaulisa 
397d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);CHKERRQ(ierr);
398d8b3b5e3Seaulisa         ierr = PetscMalloc1(m_local,&idx2);CHKERRQ(ierr);
399d8b3b5e3Seaulisa         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2);CHKERRQ(ierr);
400d8b3b5e3Seaulisa         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
401d8b3b5e3Seaulisa         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis");
402d8b3b5e3Seaulisa         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local);CHKERRQ(ierr);
403d8b3b5e3Seaulisa 
404d8b3b5e3Seaulisa         ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
4059448b7f1SJunchao Zhang         ierr = VecScatterCreate(osm->y[i],isll,osm->ly,isll_local,&osm->lprolongation[i]);CHKERRQ(ierr);
406d8b3b5e3Seaulisa 
407d8b3b5e3Seaulisa         ierr = ISDestroy(&isll);CHKERRQ(ierr);
408d8b3b5e3Seaulisa         ierr = ISDestroy(&isll_local);CHKERRQ(ierr);
409d8b3b5e3Seaulisa       }
41080ec0b7dSPatrick Sanan     }
41180ec0b7dSPatrick Sanan     ierr = VecDestroy(&vec);CHKERRQ(ierr);
41280ec0b7dSPatrick Sanan   }
41380ec0b7dSPatrick Sanan 
414fb745f2cSMatthew G. Knepley   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
415235cc4ceSMatthew G. Knepley     IS      *cis;
416235cc4ceSMatthew G. Knepley     PetscInt c;
417235cc4ceSMatthew G. Knepley 
418235cc4ceSMatthew G. Knepley     ierr = PetscMalloc1(osm->n_local_true, &cis);CHKERRQ(ierr);
419235cc4ceSMatthew G. Knepley     for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis;
4207dae84e0SHong Zhang     ierr = MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);CHKERRQ(ierr);
421235cc4ceSMatthew G. Knepley     ierr = PetscFree(cis);CHKERRQ(ierr);
422fb745f2cSMatthew G. Knepley   }
4234b9ad928SBarry Smith 
4244b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
4254b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
42692bb6962SLisandro Dalcin   ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
4274b9ad928SBarry Smith 
42892bb6962SLisandro Dalcin   /*
42992bb6962SLisandro Dalcin      Loop over subdomains putting them into local ksp
43092bb6962SLisandro Dalcin   */
43192bb6962SLisandro Dalcin   for (i=0; i<osm->n_local_true; i++) {
43223ee1639SBarry Smith     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr);
43392bb6962SLisandro Dalcin     if (!pc->setupcalled) {
434bf108f30SBarry Smith       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
4354b9ad928SBarry Smith     }
43692bb6962SLisandro Dalcin   }
4374b9ad928SBarry Smith   PetscFunctionReturn(0);
4384b9ad928SBarry Smith }
4394b9ad928SBarry Smith 
4406849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc)
4414b9ad928SBarry Smith {
4424b9ad928SBarry Smith   PC_ASM             *osm = (PC_ASM*)pc->data;
4436849ba73SBarry Smith   PetscErrorCode     ierr;
44413f74950SBarry Smith   PetscInt           i;
445539c167fSBarry Smith   KSPConvergedReason reason;
4464b9ad928SBarry Smith 
4474b9ad928SBarry Smith   PetscFunctionBegin;
4484b9ad928SBarry Smith   for (i=0; i<osm->n_local_true; i++) {
4494b9ad928SBarry Smith     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
450539c167fSBarry Smith     ierr = KSPGetConvergedReason(osm->ksp[i],&reason);CHKERRQ(ierr);
451c0decd05SBarry Smith     if (reason == KSP_DIVERGED_PC_FAILED) {
452261222b2SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
453e0eafd54SHong Zhang     }
4544b9ad928SBarry Smith   }
4554b9ad928SBarry Smith   PetscFunctionReturn(0);
4564b9ad928SBarry Smith }
4574b9ad928SBarry Smith 
4586849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y)
4594b9ad928SBarry Smith {
4604b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
4616849ba73SBarry Smith   PetscErrorCode ierr;
4621dd8081eSeaulisa   PetscInt       i,n_local_true = osm->n_local_true;
4634b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
4644b9ad928SBarry Smith 
4654b9ad928SBarry Smith   PetscFunctionBegin;
4664b9ad928SBarry Smith   /*
4674b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
4684b9ad928SBarry Smith      subdomain values (leaving the other values 0).
4694b9ad928SBarry Smith   */
4704b9ad928SBarry Smith   if (!(osm->type & PC_ASM_RESTRICT)) {
4714b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
4724b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
4731dd8081eSeaulisa     ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr);
4744b9ad928SBarry Smith   }
475347574c9Seaulisa   if (!(osm->type & PC_ASM_INTERPOLATE)) {
476347574c9Seaulisa     reverse = SCATTER_REVERSE_LOCAL;
477347574c9Seaulisa   }
4784b9ad928SBarry Smith 
479347574c9Seaulisa   if(osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE){
480b0de9f23SBarry Smith     /* zero the global and the local solutions */
48112cd4985SMatthew G. Knepley     ierr = VecZeroEntries(y);CHKERRQ(ierr);
4825b3c0d42Seaulisa     ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
483347574c9Seaulisa 
484b0de9f23SBarry Smith     /* Copy the global RHS to local RHS including the ghost nodes */
4851dd8081eSeaulisa     ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
4861dd8081eSeaulisa     ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
487347574c9Seaulisa 
488b0de9f23SBarry Smith     /* Restrict local RHS to the overlapping 0-block RHS */
4891dd8081eSeaulisa     ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
4901dd8081eSeaulisa     ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0],  INSERT_VALUES, forward);CHKERRQ(ierr);
491d8b3b5e3Seaulisa 
49212cd4985SMatthew G. Knepley     /* do the local solves */
49312cd4985SMatthew G. Knepley     for (i = 0; i < n_local_true; ++i) {
494347574c9Seaulisa 
495b0de9f23SBarry Smith       /* solve the overlapping i-block */
496fa2bb9feSLisandro Dalcin       ierr = PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
49712cd4985SMatthew G. Knepley       ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
498c0decd05SBarry Smith       ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr);
499fa2bb9feSLisandro Dalcin       ierr = PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
500d8b3b5e3Seaulisa 
501b0de9f23SBarry Smith       if (osm->lprolongation) { /* interpolate the non-overalapping i-block solution to the local solution (only for restrictive additive) */
5021dd8081eSeaulisa         ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
5031dd8081eSeaulisa         ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
504d8b3b5e3Seaulisa       }
505b0de9f23SBarry Smith       else{ /* interpolate the overalapping i-block solution to the local solution */
5061dd8081eSeaulisa         ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
5071dd8081eSeaulisa         ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
508d8b3b5e3Seaulisa       }
509347574c9Seaulisa 
510347574c9Seaulisa       if (i < n_local_true-1) {
511b0de9f23SBarry Smith         /* Restrict local RHS to the overlapping (i+1)-block RHS */
5121dd8081eSeaulisa         ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
5131dd8081eSeaulisa         ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
514347574c9Seaulisa 
515347574c9Seaulisa         if ( osm->loctype == PC_COMPOSITE_MULTIPLICATIVE){
5168966356dSPierre Jolivet           /* update the overlapping (i+1)-block RHS using the current local solution */
517347574c9Seaulisa           ierr = MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]);CHKERRQ(ierr);
518347574c9Seaulisa           ierr = VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1]); CHKERRQ(ierr);
5197c3d802fSMatthew G. Knepley         }
52012cd4985SMatthew G. Knepley       }
52112cd4985SMatthew G. Knepley     }
522b0de9f23SBarry Smith     /* Add the local solution to the global solution including the ghost nodes */
5231dd8081eSeaulisa     ierr = VecScatterBegin(osm->restriction, osm->ly, y,  ADD_VALUES, reverse);CHKERRQ(ierr);
5241dd8081eSeaulisa     ierr = VecScatterEnd(osm->restriction,  osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
525347574c9Seaulisa   }else{
526347574c9Seaulisa     SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
52712cd4985SMatthew G. Knepley   }
5284b9ad928SBarry Smith   PetscFunctionReturn(0);
5294b9ad928SBarry Smith }
5304b9ad928SBarry Smith 
5316849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y)
5324b9ad928SBarry Smith {
5334b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
5346849ba73SBarry Smith   PetscErrorCode ierr;
5351dd8081eSeaulisa   PetscInt       i,n_local_true = osm->n_local_true;
5364b9ad928SBarry Smith   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
5374b9ad928SBarry Smith 
5384b9ad928SBarry Smith   PetscFunctionBegin;
5394b9ad928SBarry Smith   /*
5404b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
5414b9ad928SBarry Smith      subdomain values (leaving the other values 0).
5424b9ad928SBarry Smith 
5434b9ad928SBarry Smith      Note: these are reversed from the PCApply_ASM() because we are applying the
5444b9ad928SBarry Smith      transpose of the three terms
5454b9ad928SBarry Smith   */
546d8b3b5e3Seaulisa 
5474b9ad928SBarry Smith   if (!(osm->type & PC_ASM_INTERPOLATE)) {
5484b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
5494b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
5501dd8081eSeaulisa     ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr);
5514b9ad928SBarry Smith   }
5522fa5cd67SKarl Rupp   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;
5534b9ad928SBarry Smith 
554b0de9f23SBarry Smith   /* zero the global and the local solutions */
55510bd88b9SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
556d8b3b5e3Seaulisa   ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
557d8b3b5e3Seaulisa 
558b0de9f23SBarry Smith   /* Copy the global RHS to local RHS including the ghost nodes */
5591dd8081eSeaulisa   ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
5601dd8081eSeaulisa   ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
561d8b3b5e3Seaulisa 
562b0de9f23SBarry Smith   /* Restrict local RHS to the overlapping 0-block RHS */
5631dd8081eSeaulisa   ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
5641dd8081eSeaulisa   ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0],  INSERT_VALUES, forward);CHKERRQ(ierr);
565d8b3b5e3Seaulisa 
5664b9ad928SBarry Smith   /* do the local solves */
567d8b3b5e3Seaulisa   for (i = 0; i < n_local_true; ++i) {
568d8b3b5e3Seaulisa 
569b0de9f23SBarry Smith     /* solve the overlapping i-block */
570fa2bb9feSLisandro Dalcin     ierr = PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
571e1bcd54cSBarry Smith     ierr = KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
572c0decd05SBarry Smith     ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr);
573fa2bb9feSLisandro Dalcin     ierr = PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
574d8b3b5e3Seaulisa 
575b0de9f23SBarry Smith     if (osm->lprolongation) { /* interpolate the non-overalapping i-block solution to the local solution */
5761dd8081eSeaulisa      ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
5771dd8081eSeaulisa      ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
578b9845e0eSMatthew Knepley     }
579b0de9f23SBarry Smith     else{ /* interpolate the overalapping i-block solution to the local solution */
5801dd8081eSeaulisa       ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
5811dd8081eSeaulisa       ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
5824b9ad928SBarry Smith     }
583d8b3b5e3Seaulisa 
584d8b3b5e3Seaulisa     if (i < n_local_true-1) {
585b0de9f23SBarry Smith       /* Restrict local RHS to the overlapping (i+1)-block RHS */
5861dd8081eSeaulisa       ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
5871dd8081eSeaulisa       ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
5884b9ad928SBarry Smith     }
5894b9ad928SBarry Smith   }
590b0de9f23SBarry Smith   /* Add the local solution to the global solution including the ghost nodes */
5911dd8081eSeaulisa   ierr = VecScatterBegin(osm->restriction, osm->ly, y,  ADD_VALUES, reverse);CHKERRQ(ierr);
5921dd8081eSeaulisa   ierr = VecScatterEnd(osm->restriction,  osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
593d8b3b5e3Seaulisa 
5944b9ad928SBarry Smith   PetscFunctionReturn(0);
595d8b3b5e3Seaulisa 
5964b9ad928SBarry Smith }
5974b9ad928SBarry Smith 
598e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc)
5994b9ad928SBarry Smith {
6004b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
6016849ba73SBarry Smith   PetscErrorCode ierr;
60213f74950SBarry Smith   PetscInt       i;
6034b9ad928SBarry Smith 
6044b9ad928SBarry Smith   PetscFunctionBegin;
60592bb6962SLisandro Dalcin   if (osm->ksp) {
60692bb6962SLisandro Dalcin     for (i=0; i<osm->n_local_true; i++) {
607e91c6855SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
60892bb6962SLisandro Dalcin     }
60992bb6962SLisandro Dalcin   }
610e09e08ccSBarry Smith   if (osm->pmat) {
61192bb6962SLisandro Dalcin     if (osm->n_local_true > 0) {
61230a70a9aSHong Zhang       ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
61392bb6962SLisandro Dalcin     }
61492bb6962SLisandro Dalcin   }
6151dd8081eSeaulisa   if (osm->lrestriction) {
6161dd8081eSeaulisa     ierr = VecScatterDestroy(&osm->restriction);CHKERRQ(ierr);
6171dd8081eSeaulisa     for (i=0; i<osm->n_local_true; i++) {
6181dd8081eSeaulisa       ierr = VecScatterDestroy(&osm->lrestriction[i]);CHKERRQ(ierr);
6191dd8081eSeaulisa       if (osm->lprolongation) {ierr = VecScatterDestroy(&osm->lprolongation[i]);CHKERRQ(ierr);}
620fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
621fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
6224b9ad928SBarry Smith     }
6231dd8081eSeaulisa     ierr = PetscFree(osm->lrestriction);CHKERRQ(ierr);
6241dd8081eSeaulisa     if (osm->lprolongation) {ierr = PetscFree(osm->lprolongation);CHKERRQ(ierr);}
62505b42c5fSBarry Smith     ierr = PetscFree(osm->x);CHKERRQ(ierr);
62678904715SLisandro Dalcin     ierr = PetscFree(osm->y);CHKERRQ(ierr);
6271dd8081eSeaulisa 
62892bb6962SLisandro Dalcin   }
6292b691e39SMatthew Knepley   ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
630fb745f2cSMatthew G. Knepley   ierr = ISDestroy(&osm->lis);CHKERRQ(ierr);
631fb745f2cSMatthew G. Knepley   ierr = VecDestroy(&osm->lx);CHKERRQ(ierr);
632fb745f2cSMatthew G. Knepley   ierr = VecDestroy(&osm->ly);CHKERRQ(ierr);
633347574c9Seaulisa   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
634347574c9Seaulisa     ierr = MatDestroyMatrices(osm->n_local_true, &osm->lmats);CHKERRQ(ierr);
635fb745f2cSMatthew G. Knepley   }
6362fa5cd67SKarl Rupp 
63780ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
63880ec0b7dSPatrick Sanan 
639e91c6855SBarry Smith   osm->is       = 0;
640e91c6855SBarry Smith   osm->is_local = 0;
641e91c6855SBarry Smith   PetscFunctionReturn(0);
642e91c6855SBarry Smith }
643e91c6855SBarry Smith 
644e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc)
645e91c6855SBarry Smith {
646e91c6855SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
647e91c6855SBarry Smith   PetscErrorCode ierr;
648e91c6855SBarry Smith   PetscInt       i;
649e91c6855SBarry Smith 
650e91c6855SBarry Smith   PetscFunctionBegin;
651e91c6855SBarry Smith   ierr = PCReset_ASM(pc);CHKERRQ(ierr);
652e91c6855SBarry Smith   if (osm->ksp) {
653e91c6855SBarry Smith     for (i=0; i<osm->n_local_true; i++) {
654fcfd50ebSBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
655e91c6855SBarry Smith     }
656e91c6855SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
657e91c6855SBarry Smith   }
658e91c6855SBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
6594b9ad928SBarry Smith   PetscFunctionReturn(0);
6604b9ad928SBarry Smith }
6614b9ad928SBarry Smith 
6624416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc)
6634b9ad928SBarry Smith {
6644b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
6656849ba73SBarry Smith   PetscErrorCode ierr;
6669dcbbd2bSBarry Smith   PetscInt       blocks,ovl;
66757501b6eSBarry Smith   PetscBool      flg;
66892bb6962SLisandro Dalcin   PCASMType      asmtype;
66912cd4985SMatthew G. Knepley   PCCompositeType loctype;
67080ec0b7dSPatrick Sanan   char           sub_mat_type[256];
6714b9ad928SBarry Smith 
6724b9ad928SBarry Smith   PetscFunctionBegin;
673e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");CHKERRQ(ierr);
674d709ea83SDmitry Karpeev   ierr = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr);
67590d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
67665db9045SDmitry Karpeev   if (flg) {
677f77a5242SKarl Rupp     ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr);
678d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
67965db9045SDmitry Karpeev   }
680342c94f9SMatthew G. Knepley   ierr = PetscOptionsInt("-pc_asm_local_blocks","Number of local subdomains","PCASMSetLocalSubdomains",osm->n_local_true,&blocks,&flg);CHKERRQ(ierr);
681342c94f9SMatthew G. Knepley   if (flg) {
682342c94f9SMatthew G. Knepley     ierr = PCASMSetLocalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr);
683342c94f9SMatthew G. Knepley     osm->dm_subdomains = PETSC_FALSE;
684342c94f9SMatthew G. Knepley   }
68590d69ab7SBarry Smith   ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
68665db9045SDmitry Karpeev   if (flg) {
68765db9045SDmitry Karpeev     ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr);
688d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
68965db9045SDmitry Karpeev   }
69090d69ab7SBarry Smith   flg  = PETSC_FALSE;
69190d69ab7SBarry Smith   ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
69292bb6962SLisandro Dalcin   if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); }
69312cd4985SMatthew G. Knepley   flg  = PETSC_FALSE;
69412cd4985SMatthew G. Knepley   ierr = PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr);
69512cd4985SMatthew G. Knepley   if (flg) {ierr = PCASMSetLocalType(pc,loctype);CHKERRQ(ierr); }
69680ec0b7dSPatrick Sanan   ierr = PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);CHKERRQ(ierr);
69780ec0b7dSPatrick Sanan   if(flg){
698459726d8SSatish Balay     ierr = PCASMSetSubMatType(pc,sub_mat_type);CHKERRQ(ierr);
69980ec0b7dSPatrick Sanan   }
7004b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7014b9ad928SBarry Smith   PetscFunctionReturn(0);
7024b9ad928SBarry Smith }
7034b9ad928SBarry Smith 
7044b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/
7054b9ad928SBarry Smith 
7061e6b0712SBarry Smith static PetscErrorCode  PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[])
7074b9ad928SBarry Smith {
7084b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
70992bb6962SLisandro Dalcin   PetscErrorCode ierr;
71092bb6962SLisandro Dalcin   PetscInt       i;
7114b9ad928SBarry Smith 
7124b9ad928SBarry Smith   PetscFunctionBegin;
713e32f2f54SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
714ce94432eSBarry 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().");
715e7e72b3dSBarry Smith 
7164b9ad928SBarry Smith   if (!pc->setupcalled) {
71792bb6962SLisandro Dalcin     if (is) {
71892bb6962SLisandro Dalcin       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
71992bb6962SLisandro Dalcin     }
720832fc9a5SMatthew Knepley     if (is_local) {
721832fc9a5SMatthew Knepley       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
722832fc9a5SMatthew Knepley     }
7232b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
7242fa5cd67SKarl Rupp 
7254b9ad928SBarry Smith     osm->n_local_true = n;
72692bb6962SLisandro Dalcin     osm->is           = 0;
7272b691e39SMatthew Knepley     osm->is_local     = 0;
72892bb6962SLisandro Dalcin     if (is) {
729785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is);CHKERRQ(ierr);
7302fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is[i] = is[i];
7313d03bb48SJed Brown       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
7323d03bb48SJed Brown       osm->overlap = -1;
73392bb6962SLisandro Dalcin     }
7342b691e39SMatthew Knepley     if (is_local) {
735785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->is_local);CHKERRQ(ierr);
7362fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->is_local[i] = is_local[i];
737a35b7b57SDmitry Karpeev       if (!is) {
738785e854fSJed Brown         ierr = PetscMalloc1(osm->n_local_true,&osm->is);CHKERRQ(ierr);
739a35b7b57SDmitry Karpeev         for (i=0; i<osm->n_local_true; i++) {
740a35b7b57SDmitry Karpeev           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
741a35b7b57SDmitry Karpeev             ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr);
742a35b7b57SDmitry Karpeev             ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr);
743a35b7b57SDmitry Karpeev           } else {
744a35b7b57SDmitry Karpeev             ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr);
745a35b7b57SDmitry Karpeev             osm->is[i] = osm->is_local[i];
746a35b7b57SDmitry Karpeev           }
747a35b7b57SDmitry Karpeev         }
748a35b7b57SDmitry Karpeev       }
7492b691e39SMatthew Knepley     }
7504b9ad928SBarry Smith   }
7514b9ad928SBarry Smith   PetscFunctionReturn(0);
7524b9ad928SBarry Smith }
7534b9ad928SBarry Smith 
7541e6b0712SBarry Smith static PetscErrorCode  PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local)
7554b9ad928SBarry Smith {
7564b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
7576849ba73SBarry Smith   PetscErrorCode ierr;
75813f74950SBarry Smith   PetscMPIInt    rank,size;
75978904715SLisandro Dalcin   PetscInt       n;
7604b9ad928SBarry Smith 
7614b9ad928SBarry Smith   PetscFunctionBegin;
762ce94432eSBarry Smith   if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
763ce94432eSBarry 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.");
7644b9ad928SBarry Smith 
7654b9ad928SBarry Smith   /*
766880770e9SJed Brown      Split the subdomains equally among all processors
7674b9ad928SBarry Smith   */
768ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
769ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
7704b9ad928SBarry Smith   n    = N/size + ((N % size) > rank);
771e32f2f54SBarry 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);
772e32f2f54SBarry Smith   if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp().");
7734b9ad928SBarry Smith   if (!pc->setupcalled) {
7742b691e39SMatthew Knepley     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
7752fa5cd67SKarl Rupp 
7764b9ad928SBarry Smith     osm->n_local_true = n;
7774b9ad928SBarry Smith     osm->is           = 0;
7782b691e39SMatthew Knepley     osm->is_local     = 0;
7794b9ad928SBarry Smith   }
7804b9ad928SBarry Smith   PetscFunctionReturn(0);
7814b9ad928SBarry Smith }
7824b9ad928SBarry Smith 
7831e6b0712SBarry Smith static PetscErrorCode  PCASMSetOverlap_ASM(PC pc,PetscInt ovl)
7844b9ad928SBarry Smith {
78592bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
7864b9ad928SBarry Smith 
7874b9ad928SBarry Smith   PetscFunctionBegin;
788ce94432eSBarry Smith   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
789ce94432eSBarry Smith   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp().");
7902fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
7914b9ad928SBarry Smith   PetscFunctionReturn(0);
7924b9ad928SBarry Smith }
7934b9ad928SBarry Smith 
7941e6b0712SBarry Smith static PetscErrorCode  PCASMSetType_ASM(PC pc,PCASMType type)
7954b9ad928SBarry Smith {
79692bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM*)pc->data;
7974b9ad928SBarry Smith 
7984b9ad928SBarry Smith   PetscFunctionBegin;
7994b9ad928SBarry Smith   osm->type     = type;
800bf108f30SBarry Smith   osm->type_set = PETSC_TRUE;
8014b9ad928SBarry Smith   PetscFunctionReturn(0);
8024b9ad928SBarry Smith }
8034b9ad928SBarry Smith 
804c60c7ad4SBarry Smith static PetscErrorCode  PCASMGetType_ASM(PC pc,PCASMType *type)
805c60c7ad4SBarry Smith {
806c60c7ad4SBarry Smith   PC_ASM *osm = (PC_ASM*)pc->data;
807c60c7ad4SBarry Smith 
808c60c7ad4SBarry Smith   PetscFunctionBegin;
809c60c7ad4SBarry Smith   *type = osm->type;
810c60c7ad4SBarry Smith   PetscFunctionReturn(0);
811c60c7ad4SBarry Smith }
812c60c7ad4SBarry Smith 
81312cd4985SMatthew G. Knepley static PetscErrorCode  PCASMSetLocalType_ASM(PC pc, PCCompositeType type)
81412cd4985SMatthew G. Knepley {
81512cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *) pc->data;
81612cd4985SMatthew G. Knepley 
81712cd4985SMatthew G. Knepley   PetscFunctionBegin;
818d0ecd4dfSBarry 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");
81912cd4985SMatthew G. Knepley   osm->loctype = type;
82012cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
82112cd4985SMatthew G. Knepley }
82212cd4985SMatthew G. Knepley 
82312cd4985SMatthew G. Knepley static PetscErrorCode  PCASMGetLocalType_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;
82812cd4985SMatthew G. Knepley   *type = osm->loctype;
82912cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
83012cd4985SMatthew G. Knepley }
83112cd4985SMatthew G. Knepley 
8321e6b0712SBarry Smith static PetscErrorCode  PCASMSetSortIndices_ASM(PC pc,PetscBool  doSort)
8336ed231c7SMatthew Knepley {
8346ed231c7SMatthew Knepley   PC_ASM *osm = (PC_ASM*)pc->data;
8356ed231c7SMatthew Knepley 
8366ed231c7SMatthew Knepley   PetscFunctionBegin;
8376ed231c7SMatthew Knepley   osm->sort_indices = doSort;
8386ed231c7SMatthew Knepley   PetscFunctionReturn(0);
8396ed231c7SMatthew Knepley }
8406ed231c7SMatthew Knepley 
8411e6b0712SBarry Smith static PetscErrorCode  PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
8424b9ad928SBarry Smith {
84392bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM*)pc->data;
844dfbe8321SBarry Smith   PetscErrorCode ierr;
8454b9ad928SBarry Smith 
8464b9ad928SBarry Smith   PetscFunctionBegin;
84734a84908Sprj-   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");
8484b9ad928SBarry Smith 
8492fa5cd67SKarl Rupp   if (n_local) *n_local = osm->n_local_true;
85092bb6962SLisandro Dalcin   if (first_local) {
851ce94432eSBarry Smith     ierr          = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
85292bb6962SLisandro Dalcin     *first_local -= osm->n_local_true;
85392bb6962SLisandro Dalcin   }
85492bb6962SLisandro Dalcin   if (ksp) {
85592bb6962SLisandro Dalcin     /* Assume that local solves are now different; not necessarily
85692bb6962SLisandro Dalcin        true though!  This flag is used only for PCView_ASM() */
85792bb6962SLisandro Dalcin     *ksp                   = osm->ksp;
85892bb6962SLisandro Dalcin     osm->same_local_solves = PETSC_FALSE;
85992bb6962SLisandro Dalcin   }
8604b9ad928SBarry Smith   PetscFunctionReturn(0);
8614b9ad928SBarry Smith }
8624b9ad928SBarry Smith 
86380ec0b7dSPatrick Sanan static PetscErrorCode  PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type)
86480ec0b7dSPatrick Sanan {
86580ec0b7dSPatrick Sanan   PC_ASM         *osm = (PC_ASM*)pc->data;
86680ec0b7dSPatrick Sanan 
86780ec0b7dSPatrick Sanan   PetscFunctionBegin;
86880ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
86980ec0b7dSPatrick Sanan   PetscValidPointer(sub_mat_type,2);
87080ec0b7dSPatrick Sanan   *sub_mat_type = osm->sub_mat_type;
87180ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
87280ec0b7dSPatrick Sanan }
87380ec0b7dSPatrick Sanan 
87480ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type)
87580ec0b7dSPatrick Sanan {
87680ec0b7dSPatrick Sanan   PetscErrorCode    ierr;
87780ec0b7dSPatrick Sanan   PC_ASM            *osm = (PC_ASM*)pc->data;
87880ec0b7dSPatrick Sanan 
87980ec0b7dSPatrick Sanan   PetscFunctionBegin;
88080ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
88180ec0b7dSPatrick Sanan   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
88280ec0b7dSPatrick Sanan   ierr = PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);CHKERRQ(ierr);
88380ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
88480ec0b7dSPatrick Sanan }
88580ec0b7dSPatrick Sanan 
8864b9ad928SBarry Smith /*@C
8871093a601SBarry Smith     PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner.
8884b9ad928SBarry Smith 
889d083f849SBarry Smith     Collective on pc
8904b9ad928SBarry Smith 
8914b9ad928SBarry Smith     Input Parameters:
8924b9ad928SBarry Smith +   pc - the preconditioner context
8934b9ad928SBarry Smith .   n - the number of subdomains for this processor (default value = 1)
8948c03b21aSDmitry Karpeev .   is - the index set that defines the subdomains for this processor
8950298fd71SBarry Smith          (or NULL for PETSc to determine subdomains)
896f1ee410cSBarry 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
897f1ee410cSBarry Smith          (or NULL to not provide these)
8984b9ad928SBarry Smith 
899342c94f9SMatthew G. Knepley     Options Database Key:
900342c94f9SMatthew G. Knepley     To set the total number of subdomain blocks rather than specify the
901342c94f9SMatthew G. Knepley     index sets, use the option
902342c94f9SMatthew G. Knepley .    -pc_asm_local_blocks <blks> - Sets local blocks
903342c94f9SMatthew G. Knepley 
9044b9ad928SBarry Smith     Notes:
9051093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
9064b9ad928SBarry Smith 
9074b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
9084b9ad928SBarry Smith 
9094b9ad928SBarry Smith     Use PCASMSetTotalSubdomains() to set the subdomains for all processors.
9104b9ad928SBarry Smith 
911f1ee410cSBarry Smith     If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated
912f1ee410cSBarry Smith     back to form the global solution (this is the standard restricted additive Schwarz method)
913f1ee410cSBarry Smith 
914f1ee410cSBarry 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
915f1ee410cSBarry Smith     no code to handle that case.
916f1ee410cSBarry Smith 
9174b9ad928SBarry Smith     Level: advanced
9184b9ad928SBarry Smith 
9194b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
920f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType()
9214b9ad928SBarry Smith @*/
9227087cfbeSBarry Smith PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
9234b9ad928SBarry Smith {
9247bb14e67SBarry Smith   PetscErrorCode ierr;
9254b9ad928SBarry Smith 
9264b9ad928SBarry Smith   PetscFunctionBegin;
9270700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9287bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr);
9294b9ad928SBarry Smith   PetscFunctionReturn(0);
9304b9ad928SBarry Smith }
9314b9ad928SBarry Smith 
9324b9ad928SBarry Smith /*@C
933feb221f8SDmitry Karpeev     PCASMSetTotalSubdomains - Sets the subdomains for all processors for the
9344b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
9354b9ad928SBarry Smith     PC communicator must call this routine, with the same index sets.
9364b9ad928SBarry Smith 
937d083f849SBarry Smith     Collective on pc
9384b9ad928SBarry Smith 
9394b9ad928SBarry Smith     Input Parameters:
9404b9ad928SBarry Smith +   pc - the preconditioner context
941feb221f8SDmitry Karpeev .   N  - the number of subdomains for all processors
942feb221f8SDmitry Karpeev .   is - the index sets that define the subdomains for all processors
943dfaa0c5fSBarry Smith          (or NULL to ask PETSc to determine the subdomains)
9442b691e39SMatthew Knepley -   is_local - the index sets that define the local part of the subdomains for this processor
945f1ee410cSBarry Smith          (or NULL to not provide this information)
9464b9ad928SBarry Smith 
9474b9ad928SBarry Smith     Options Database Key:
9484b9ad928SBarry Smith     To set the total number of subdomain blocks rather than specify the
9494b9ad928SBarry Smith     index sets, use the option
9504b9ad928SBarry Smith .    -pc_asm_blocks <blks> - Sets total blocks
9514b9ad928SBarry Smith 
9524b9ad928SBarry Smith     Notes:
953f1ee410cSBarry Smith     Currently you cannot use this to set the actual subdomains with the argument is or is_local.
9544b9ad928SBarry Smith 
9554b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
9564b9ad928SBarry Smith 
9574b9ad928SBarry Smith     These index sets cannot be destroyed until after completion of the
9584b9ad928SBarry Smith     linear solves for which the ASM preconditioner is being used.
9594b9ad928SBarry Smith 
9604b9ad928SBarry Smith     Use PCASMSetLocalSubdomains() to set local subdomains.
9614b9ad928SBarry Smith 
9621093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
9631093a601SBarry Smith 
9644b9ad928SBarry Smith     Level: advanced
9654b9ad928SBarry Smith 
9664b9ad928SBarry Smith .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
9674b9ad928SBarry Smith           PCASMCreateSubdomains2D()
9684b9ad928SBarry Smith @*/
9697087cfbeSBarry Smith PetscErrorCode  PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[])
9704b9ad928SBarry Smith {
9717bb14e67SBarry Smith   PetscErrorCode ierr;
9724b9ad928SBarry Smith 
9734b9ad928SBarry Smith   PetscFunctionBegin;
9740700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9757bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr);
9764b9ad928SBarry Smith   PetscFunctionReturn(0);
9774b9ad928SBarry Smith }
9784b9ad928SBarry Smith 
9794b9ad928SBarry Smith /*@
9804b9ad928SBarry Smith     PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
9814b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
982f1ee410cSBarry Smith     PC communicator must call this routine.
9834b9ad928SBarry Smith 
984d083f849SBarry Smith     Logically Collective on pc
9854b9ad928SBarry Smith 
9864b9ad928SBarry Smith     Input Parameters:
9874b9ad928SBarry Smith +   pc  - the preconditioner context
9884b9ad928SBarry Smith -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
9894b9ad928SBarry Smith 
9904b9ad928SBarry Smith     Options Database Key:
9914b9ad928SBarry Smith .   -pc_asm_overlap <ovl> - Sets overlap
9924b9ad928SBarry Smith 
9934b9ad928SBarry Smith     Notes:
9944b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.  To use
9954b9ad928SBarry Smith     multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
9964b9ad928SBarry Smith     PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).
9974b9ad928SBarry Smith 
9984b9ad928SBarry Smith     The overlap defaults to 1, so if one desires that no additional
9994b9ad928SBarry Smith     overlap be computed beyond what may have been set with a call to
10004b9ad928SBarry Smith     PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
10014b9ad928SBarry Smith     must be set to be 0.  In particular, if one does not explicitly set
10024b9ad928SBarry Smith     the subdomains an application code, then all overlap would be computed
10034b9ad928SBarry Smith     internally by PETSc, and using an overlap of 0 would result in an ASM
10044b9ad928SBarry Smith     variant that is equivalent to the block Jacobi preconditioner.
10054b9ad928SBarry Smith 
1006f1ee410cSBarry Smith     The default algorithm used by PETSc to increase overlap is fast, but not scalable,
1007f1ee410cSBarry Smith     use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
1008f1ee410cSBarry Smith 
10094b9ad928SBarry Smith     Note that one can define initial index sets with any overlap via
1010f1ee410cSBarry Smith     PCASMSetLocalSubdomains(); the routine
10114b9ad928SBarry Smith     PCASMSetOverlap() merely allows PETSc to extend that overlap further
10124b9ad928SBarry Smith     if desired.
10134b9ad928SBarry Smith 
10144b9ad928SBarry Smith     Level: intermediate
10154b9ad928SBarry Smith 
10164b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
1017f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap()
10184b9ad928SBarry Smith @*/
10197087cfbeSBarry Smith PetscErrorCode  PCASMSetOverlap(PC pc,PetscInt ovl)
10204b9ad928SBarry Smith {
10217bb14e67SBarry Smith   PetscErrorCode ierr;
10224b9ad928SBarry Smith 
10234b9ad928SBarry Smith   PetscFunctionBegin;
10240700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1025c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,ovl,2);
10267bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
10274b9ad928SBarry Smith   PetscFunctionReturn(0);
10284b9ad928SBarry Smith }
10294b9ad928SBarry Smith 
10304b9ad928SBarry Smith /*@
10314b9ad928SBarry Smith     PCASMSetType - Sets the type of restriction and interpolation used
10324b9ad928SBarry Smith     for local problems in the additive Schwarz method.
10334b9ad928SBarry Smith 
1034d083f849SBarry Smith     Logically Collective on pc
10354b9ad928SBarry Smith 
10364b9ad928SBarry Smith     Input Parameters:
10374b9ad928SBarry Smith +   pc  - the preconditioner context
10384b9ad928SBarry Smith -   type - variant of ASM, one of
10394b9ad928SBarry Smith .vb
10404b9ad928SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
10414b9ad928SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
10424b9ad928SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
10434b9ad928SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
10444b9ad928SBarry Smith .ve
10454b9ad928SBarry Smith 
10464b9ad928SBarry Smith     Options Database Key:
10474b9ad928SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
10484b9ad928SBarry Smith 
104995452b02SPatrick Sanan     Notes:
105095452b02SPatrick Sanan     if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected
1051f1ee410cSBarry Smith     to limit the local processor interpolation
1052f1ee410cSBarry Smith 
10534b9ad928SBarry Smith     Level: intermediate
10544b9ad928SBarry Smith 
10554b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1056f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType()
10574b9ad928SBarry Smith @*/
10587087cfbeSBarry Smith PetscErrorCode  PCASMSetType(PC pc,PCASMType type)
10594b9ad928SBarry Smith {
10607bb14e67SBarry Smith   PetscErrorCode ierr;
10614b9ad928SBarry Smith 
10624b9ad928SBarry Smith   PetscFunctionBegin;
10630700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1064c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
10657bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr);
10664b9ad928SBarry Smith   PetscFunctionReturn(0);
10674b9ad928SBarry Smith }
10684b9ad928SBarry Smith 
1069c60c7ad4SBarry Smith /*@
1070c60c7ad4SBarry Smith     PCASMGetType - Gets the type of restriction and interpolation used
1071c60c7ad4SBarry Smith     for local problems in the additive Schwarz method.
1072c60c7ad4SBarry Smith 
1073d083f849SBarry Smith     Logically Collective on pc
1074c60c7ad4SBarry Smith 
1075c60c7ad4SBarry Smith     Input Parameter:
1076c60c7ad4SBarry Smith .   pc  - the preconditioner context
1077c60c7ad4SBarry Smith 
1078c60c7ad4SBarry Smith     Output Parameter:
1079c60c7ad4SBarry Smith .   type - variant of ASM, one of
1080c60c7ad4SBarry Smith 
1081c60c7ad4SBarry Smith .vb
1082c60c7ad4SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
1083c60c7ad4SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1084c60c7ad4SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1085c60c7ad4SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
1086c60c7ad4SBarry Smith .ve
1087c60c7ad4SBarry Smith 
1088c60c7ad4SBarry Smith     Options Database Key:
1089c60c7ad4SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
1090c60c7ad4SBarry Smith 
1091c60c7ad4SBarry Smith     Level: intermediate
1092c60c7ad4SBarry Smith 
1093c60c7ad4SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1094f1ee410cSBarry Smith           PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType()
1095c60c7ad4SBarry Smith @*/
1096c60c7ad4SBarry Smith PetscErrorCode  PCASMGetType(PC pc,PCASMType *type)
1097c60c7ad4SBarry Smith {
1098c60c7ad4SBarry Smith   PetscErrorCode ierr;
1099c60c7ad4SBarry Smith 
1100c60c7ad4SBarry Smith   PetscFunctionBegin;
1101c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1102c60c7ad4SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));CHKERRQ(ierr);
1103c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1104c60c7ad4SBarry Smith }
1105c60c7ad4SBarry Smith 
110612cd4985SMatthew G. Knepley /*@
110712cd4985SMatthew G. Knepley   PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method.
110812cd4985SMatthew G. Knepley 
1109d083f849SBarry Smith   Logically Collective on pc
111012cd4985SMatthew G. Knepley 
111112cd4985SMatthew G. Knepley   Input Parameters:
111212cd4985SMatthew G. Knepley + pc  - the preconditioner context
111312cd4985SMatthew G. Knepley - type - type of composition, one of
111412cd4985SMatthew G. Knepley .vb
111512cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
111612cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
111712cd4985SMatthew G. Knepley .ve
111812cd4985SMatthew G. Knepley 
111912cd4985SMatthew G. Knepley   Options Database Key:
112012cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
112112cd4985SMatthew G. Knepley 
112212cd4985SMatthew G. Knepley   Level: intermediate
112312cd4985SMatthew G. Knepley 
1124f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
112512cd4985SMatthew G. Knepley @*/
112612cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type)
112712cd4985SMatthew G. Knepley {
112812cd4985SMatthew G. Knepley   PetscErrorCode ierr;
112912cd4985SMatthew G. Knepley 
113012cd4985SMatthew G. Knepley   PetscFunctionBegin;
113112cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
113212cd4985SMatthew G. Knepley   PetscValidLogicalCollectiveEnum(pc, type, 2);
113312cd4985SMatthew G. Knepley   ierr = PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));CHKERRQ(ierr);
113412cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
113512cd4985SMatthew G. Knepley }
113612cd4985SMatthew G. Knepley 
113712cd4985SMatthew G. Knepley /*@
113812cd4985SMatthew G. Knepley   PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method.
113912cd4985SMatthew G. Knepley 
1140d083f849SBarry Smith   Logically Collective on pc
114112cd4985SMatthew G. Knepley 
114212cd4985SMatthew G. Knepley   Input Parameter:
114312cd4985SMatthew G. Knepley . pc  - the preconditioner context
114412cd4985SMatthew G. Knepley 
114512cd4985SMatthew G. Knepley   Output Parameter:
114612cd4985SMatthew G. Knepley . type - type of composition, one of
114712cd4985SMatthew G. Knepley .vb
114812cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
114912cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
115012cd4985SMatthew G. Knepley .ve
115112cd4985SMatthew G. Knepley 
115212cd4985SMatthew G. Knepley   Options Database Key:
115312cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
115412cd4985SMatthew G. Knepley 
115512cd4985SMatthew G. Knepley   Level: intermediate
115612cd4985SMatthew G. Knepley 
1157f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
115812cd4985SMatthew G. Knepley @*/
115912cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type)
116012cd4985SMatthew G. Knepley {
116112cd4985SMatthew G. Knepley   PetscErrorCode ierr;
116212cd4985SMatthew G. Knepley 
116312cd4985SMatthew G. Knepley   PetscFunctionBegin;
116412cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
116512cd4985SMatthew G. Knepley   PetscValidPointer(type, 2);
116612cd4985SMatthew G. Knepley   ierr = PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));CHKERRQ(ierr);
116712cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
116812cd4985SMatthew G. Knepley }
116912cd4985SMatthew G. Knepley 
11706ed231c7SMatthew Knepley /*@
11716ed231c7SMatthew Knepley     PCASMSetSortIndices - Determines whether subdomain indices are sorted.
11726ed231c7SMatthew Knepley 
1173d083f849SBarry Smith     Logically Collective on pc
11746ed231c7SMatthew Knepley 
11756ed231c7SMatthew Knepley     Input Parameters:
11766ed231c7SMatthew Knepley +   pc  - the preconditioner context
11776ed231c7SMatthew Knepley -   doSort - sort the subdomain indices
11786ed231c7SMatthew Knepley 
11796ed231c7SMatthew Knepley     Level: intermediate
11806ed231c7SMatthew Knepley 
11816ed231c7SMatthew Knepley .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
11826ed231c7SMatthew Knepley           PCASMCreateSubdomains2D()
11836ed231c7SMatthew Knepley @*/
11847087cfbeSBarry Smith PetscErrorCode  PCASMSetSortIndices(PC pc,PetscBool doSort)
11856ed231c7SMatthew Knepley {
11867bb14e67SBarry Smith   PetscErrorCode ierr;
11876ed231c7SMatthew Knepley 
11886ed231c7SMatthew Knepley   PetscFunctionBegin;
11890700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1190acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
11917bb14e67SBarry Smith   ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
11926ed231c7SMatthew Knepley   PetscFunctionReturn(0);
11936ed231c7SMatthew Knepley }
11946ed231c7SMatthew Knepley 
11954b9ad928SBarry Smith /*@C
11964b9ad928SBarry Smith    PCASMGetSubKSP - Gets the local KSP contexts for all blocks on
11974b9ad928SBarry Smith    this processor.
11984b9ad928SBarry Smith 
1199d083f849SBarry Smith    Collective on pc iff first_local is requested
12004b9ad928SBarry Smith 
12014b9ad928SBarry Smith    Input Parameter:
12024b9ad928SBarry Smith .  pc - the preconditioner context
12034b9ad928SBarry Smith 
12044b9ad928SBarry Smith    Output Parameters:
12050298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
12060298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
12070298fd71SBarry Smith                  all processors must request or all must pass NULL
12084b9ad928SBarry Smith -  ksp - the array of KSP contexts
12094b9ad928SBarry Smith 
12104b9ad928SBarry Smith    Note:
1211d29017ddSJed Brown    After PCASMGetSubKSP() the array of KSPes is not to be freed.
12124b9ad928SBarry Smith 
12134b9ad928SBarry Smith    You must call KSPSetUp() before calling PCASMGetSubKSP().
12144b9ad928SBarry Smith 
1215d29017ddSJed Brown    Fortran note:
12162bf68e3eSBarry 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.
1217d29017ddSJed Brown 
12184b9ad928SBarry Smith    Level: advanced
12194b9ad928SBarry Smith 
12204b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
12214b9ad928SBarry Smith           PCASMCreateSubdomains2D(),
12224b9ad928SBarry Smith @*/
12237087cfbeSBarry Smith PetscErrorCode  PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
12244b9ad928SBarry Smith {
12257bb14e67SBarry Smith   PetscErrorCode ierr;
12264b9ad928SBarry Smith 
12274b9ad928SBarry Smith   PetscFunctionBegin;
12280700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12297bb14e67SBarry Smith   ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
12304b9ad928SBarry Smith   PetscFunctionReturn(0);
12314b9ad928SBarry Smith }
12324b9ad928SBarry Smith 
12334b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
12344b9ad928SBarry Smith /*MC
12353b09bd56SBarry Smith    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
12364b9ad928SBarry Smith            its own KSP object.
12374b9ad928SBarry Smith 
12384b9ad928SBarry Smith    Options Database Keys:
123949517cdeSBarry Smith +  -pc_asm_blocks <blks> - Sets total blocks
12404b9ad928SBarry Smith .  -pc_asm_overlap <ovl> - Sets overlap
1241f1ee410cSBarry Smith .  -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict
1242f1ee410cSBarry Smith -  -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive
12434b9ad928SBarry Smith 
12443b09bd56SBarry Smith      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
12453b09bd56SBarry Smith       will get a different convergence rate due to the default option of -pc_asm_type restrict. Use
12463b09bd56SBarry Smith       -pc_asm_type basic to use the standard ASM.
12473b09bd56SBarry Smith 
124895452b02SPatrick Sanan    Notes:
124995452b02SPatrick Sanan     Each processor can have one or more blocks, but a block cannot be shared by more
1250f1ee410cSBarry Smith      than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor.
12514b9ad928SBarry Smith 
12523b09bd56SBarry Smith      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1253d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
12544b9ad928SBarry Smith 
1255a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCASMGetSubKSP()
12564b9ad928SBarry Smith          and set the options directly on the resulting KSP object (you can access its PC
12574b9ad928SBarry Smith          with KSPGetPC())
12584b9ad928SBarry Smith 
12594b9ad928SBarry Smith    Level: beginner
12604b9ad928SBarry Smith 
1261c582cd25SBarry Smith     References:
126296a0c994SBarry Smith +   1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
126396a0c994SBarry Smith      Courant Institute, New York University Technical report
12646d33885cSprj- -   2. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
126596a0c994SBarry Smith     Cambridge University Press.
1266c582cd25SBarry Smith 
12674b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1268f1ee410cSBarry Smith            PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType()
126934a84908Sprj-            PCASMSetTotalSubdomains(), PCSetModifySubMatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType
1270e09e08ccSBarry Smith 
12714b9ad928SBarry Smith M*/
12724b9ad928SBarry Smith 
12738cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc)
12744b9ad928SBarry Smith {
1275dfbe8321SBarry Smith   PetscErrorCode ierr;
12764b9ad928SBarry Smith   PC_ASM         *osm;
12774b9ad928SBarry Smith 
12784b9ad928SBarry Smith   PetscFunctionBegin;
1279b00a9115SJed Brown   ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr);
12802fa5cd67SKarl Rupp 
12814b9ad928SBarry Smith   osm->n                 = PETSC_DECIDE;
12824b9ad928SBarry Smith   osm->n_local           = 0;
12832b837212SDmitry Karpeev   osm->n_local_true      = PETSC_DECIDE;
12844b9ad928SBarry Smith   osm->overlap           = 1;
12854b9ad928SBarry Smith   osm->ksp               = 0;
12862b691e39SMatthew Knepley   osm->restriction       = 0;
12875b3c0d42Seaulisa   osm->lprolongation     = 0;
12881dd8081eSeaulisa   osm->lrestriction      = 0;
128992bb6962SLisandro Dalcin   osm->x                 = 0;
129092bb6962SLisandro Dalcin   osm->y                 = 0;
12914b9ad928SBarry Smith   osm->is                = 0;
12922b691e39SMatthew Knepley   osm->is_local          = 0;
12934b9ad928SBarry Smith   osm->mat               = 0;
12944b9ad928SBarry Smith   osm->pmat              = 0;
12954b9ad928SBarry Smith   osm->type              = PC_ASM_RESTRICT;
129612cd4985SMatthew G. Knepley   osm->loctype           = PC_COMPOSITE_ADDITIVE;
12974b9ad928SBarry Smith   osm->same_local_solves = PETSC_TRUE;
12986ed231c7SMatthew Knepley   osm->sort_indices      = PETSC_TRUE;
1299d709ea83SDmitry Karpeev   osm->dm_subdomains     = PETSC_FALSE;
130080ec0b7dSPatrick Sanan   osm->sub_mat_type      = NULL;
13014b9ad928SBarry Smith 
130292bb6962SLisandro Dalcin   pc->data                 = (void*)osm;
13034b9ad928SBarry Smith   pc->ops->apply           = PCApply_ASM;
13044b9ad928SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_ASM;
13054b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_ASM;
1306e91c6855SBarry Smith   pc->ops->reset           = PCReset_ASM;
13074b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_ASM;
13084b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_ASM;
13094b9ad928SBarry Smith   pc->ops->setuponblocks   = PCSetUpOnBlocks_ASM;
13104b9ad928SBarry Smith   pc->ops->view            = PCView_ASM;
13114b9ad928SBarry Smith   pc->ops->applyrichardson = 0;
13124b9ad928SBarry Smith 
1313bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr);
1314bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr);
1315bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);CHKERRQ(ierr);
1316bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);CHKERRQ(ierr);
1317c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);CHKERRQ(ierr);
131812cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);CHKERRQ(ierr);
131912cd4985SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);CHKERRQ(ierr);
1320bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);CHKERRQ(ierr);
1321bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);CHKERRQ(ierr);
132280ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);CHKERRQ(ierr);
132380ec0b7dSPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);CHKERRQ(ierr);
13244b9ad928SBarry Smith   PetscFunctionReturn(0);
13254b9ad928SBarry Smith }
13264b9ad928SBarry Smith 
132792bb6962SLisandro Dalcin /*@C
132892bb6962SLisandro Dalcin    PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
132992bb6962SLisandro Dalcin    preconditioner for a any problem on a general grid.
133092bb6962SLisandro Dalcin 
133192bb6962SLisandro Dalcin    Collective
133292bb6962SLisandro Dalcin 
133392bb6962SLisandro Dalcin    Input Parameters:
133492bb6962SLisandro Dalcin +  A - The global matrix operator
133592bb6962SLisandro Dalcin -  n - the number of local blocks
133692bb6962SLisandro Dalcin 
133792bb6962SLisandro Dalcin    Output Parameters:
133892bb6962SLisandro Dalcin .  outis - the array of index sets defining the subdomains
133992bb6962SLisandro Dalcin 
134092bb6962SLisandro Dalcin    Level: advanced
134192bb6962SLisandro Dalcin 
13427d6bfa3bSBarry Smith    Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap
13437d6bfa3bSBarry Smith     from these if you use PCASMSetLocalSubdomains()
13447d6bfa3bSBarry Smith 
13457d6bfa3bSBarry Smith     In the Fortran version you must provide the array outis[] already allocated of length n.
13467d6bfa3bSBarry Smith 
134792bb6962SLisandro Dalcin .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains()
134892bb6962SLisandro Dalcin @*/
13497087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
135092bb6962SLisandro Dalcin {
135192bb6962SLisandro Dalcin   MatPartitioning mpart;
135292bb6962SLisandro Dalcin   const char      *prefix;
135392bb6962SLisandro Dalcin   PetscInt        i,j,rstart,rend,bs;
1354976e8c5aSLisandro Dalcin   PetscBool       hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
13550298fd71SBarry Smith   Mat             Ad     = NULL, adj;
135692bb6962SLisandro Dalcin   IS              ispart,isnumb,*is;
135792bb6962SLisandro Dalcin   PetscErrorCode  ierr;
135892bb6962SLisandro Dalcin 
135992bb6962SLisandro Dalcin   PetscFunctionBegin;
13600700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
136192bb6962SLisandro Dalcin   PetscValidPointer(outis,3);
1362c1235816SBarry Smith   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
136392bb6962SLisandro Dalcin 
136492bb6962SLisandro Dalcin   /* Get prefix, row distribution, and block size */
136592bb6962SLisandro Dalcin   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
136692bb6962SLisandro Dalcin   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
136792bb6962SLisandro Dalcin   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
136865e19b50SBarry 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);
136965e19b50SBarry Smith 
137092bb6962SLisandro Dalcin   /* Get diagonal block from matrix if possible */
1371976e8c5aSLisandro Dalcin   ierr = MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr);
1372976e8c5aSLisandro Dalcin   if (hasop) {
137311bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
137492bb6962SLisandro Dalcin   }
137592bb6962SLisandro Dalcin   if (Ad) {
1376b9e7e5c1SBarry Smith     ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1377b9e7e5c1SBarry Smith     if (!isbaij) {ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
137892bb6962SLisandro Dalcin   }
137992bb6962SLisandro Dalcin   if (Ad && n > 1) {
1380ace3abfcSBarry Smith     PetscBool match,done;
138192bb6962SLisandro Dalcin     /* Try to setup a good matrix partitioning if available */
138292bb6962SLisandro Dalcin     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
138392bb6962SLisandro Dalcin     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
138492bb6962SLisandro Dalcin     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1385251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
138692bb6962SLisandro Dalcin     if (!match) {
1387251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
138892bb6962SLisandro Dalcin     }
138992bb6962SLisandro Dalcin     if (!match) { /* assume a "good" partitioner is available */
13901a83f524SJed Brown       PetscInt       na;
13911a83f524SJed Brown       const PetscInt *ia,*ja;
139292bb6962SLisandro Dalcin       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
139392bb6962SLisandro Dalcin       if (done) {
139492bb6962SLisandro Dalcin         /* Build adjacency matrix by hand. Unfortunately a call to
139592bb6962SLisandro Dalcin            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
139692bb6962SLisandro Dalcin            remove the block-aij structure and we cannot expect
139792bb6962SLisandro Dalcin            MatPartitioning to split vertices as we need */
13981a83f524SJed Brown         PetscInt       i,j,len,nnz,cnt,*iia=0,*jja=0;
13991a83f524SJed Brown         const PetscInt *row;
140092bb6962SLisandro Dalcin         nnz = 0;
140192bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* count number of nonzeros */
140292bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
140392bb6962SLisandro Dalcin           row = ja + ia[i];
140492bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
140592bb6962SLisandro Dalcin             if (row[j] == i) { /* don't count diagonal */
140692bb6962SLisandro Dalcin               len--; break;
140792bb6962SLisandro Dalcin             }
140892bb6962SLisandro Dalcin           }
140992bb6962SLisandro Dalcin           nnz += len;
141092bb6962SLisandro Dalcin         }
1411854ce69bSBarry Smith         ierr   = PetscMalloc1(na+1,&iia);CHKERRQ(ierr);
1412854ce69bSBarry Smith         ierr   = PetscMalloc1(nnz,&jja);CHKERRQ(ierr);
141392bb6962SLisandro Dalcin         nnz    = 0;
141492bb6962SLisandro Dalcin         iia[0] = 0;
141592bb6962SLisandro Dalcin         for (i=0; i<na; i++) { /* fill adjacency */
141692bb6962SLisandro Dalcin           cnt = 0;
141792bb6962SLisandro Dalcin           len = ia[i+1] - ia[i];
141892bb6962SLisandro Dalcin           row = ja + ia[i];
141992bb6962SLisandro Dalcin           for (j=0; j<len; j++) {
142092bb6962SLisandro Dalcin             if (row[j] != i) { /* if not diagonal */
142192bb6962SLisandro Dalcin               jja[nnz+cnt++] = row[j];
142292bb6962SLisandro Dalcin             }
142392bb6962SLisandro Dalcin           }
142492bb6962SLisandro Dalcin           nnz     += cnt;
142592bb6962SLisandro Dalcin           iia[i+1] = nnz;
142692bb6962SLisandro Dalcin         }
142792bb6962SLisandro Dalcin         /* Partitioning of the adjacency matrix */
14280298fd71SBarry Smith         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr);
142992bb6962SLisandro Dalcin         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
143092bb6962SLisandro Dalcin         ierr      = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
143192bb6962SLisandro Dalcin         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
143292bb6962SLisandro Dalcin         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
1433fcfd50ebSBarry Smith         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
143492bb6962SLisandro Dalcin         foundpart = PETSC_TRUE;
143592bb6962SLisandro Dalcin       }
143692bb6962SLisandro Dalcin       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
143792bb6962SLisandro Dalcin     }
1438fcfd50ebSBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
143992bb6962SLisandro Dalcin   }
144092bb6962SLisandro Dalcin 
1441785e854fSJed Brown   ierr   = PetscMalloc1(n,&is);CHKERRQ(ierr);
144292bb6962SLisandro Dalcin   *outis = is;
144392bb6962SLisandro Dalcin 
144492bb6962SLisandro Dalcin   if (!foundpart) {
144592bb6962SLisandro Dalcin 
144692bb6962SLisandro Dalcin     /* Partitioning by contiguous chunks of rows */
144792bb6962SLisandro Dalcin 
144892bb6962SLisandro Dalcin     PetscInt mbs   = (rend-rstart)/bs;
144992bb6962SLisandro Dalcin     PetscInt start = rstart;
145092bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
145192bb6962SLisandro Dalcin       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
145292bb6962SLisandro Dalcin       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
145392bb6962SLisandro Dalcin       start += count;
145492bb6962SLisandro Dalcin     }
145592bb6962SLisandro Dalcin 
145692bb6962SLisandro Dalcin   } else {
145792bb6962SLisandro Dalcin 
145892bb6962SLisandro Dalcin     /* Partitioning by adjacency of diagonal block  */
145992bb6962SLisandro Dalcin 
146092bb6962SLisandro Dalcin     const PetscInt *numbering;
146192bb6962SLisandro Dalcin     PetscInt       *count,nidx,*indices,*newidx,start=0;
146292bb6962SLisandro Dalcin     /* Get node count in each partition */
1463785e854fSJed Brown     ierr = PetscMalloc1(n,&count);CHKERRQ(ierr);
146492bb6962SLisandro Dalcin     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
146592bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
146692bb6962SLisandro Dalcin       for (i=0; i<n; i++) count[i] *= bs;
146792bb6962SLisandro Dalcin     }
146892bb6962SLisandro Dalcin     /* Build indices from node numbering */
146992bb6962SLisandro Dalcin     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1470785e854fSJed Brown     ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr);
147192bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
147292bb6962SLisandro Dalcin     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
147392bb6962SLisandro Dalcin     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
147492bb6962SLisandro Dalcin     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
147592bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1476785e854fSJed Brown       ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr);
14772fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
14782fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
14792fa5cd67SKarl Rupp       }
148092bb6962SLisandro Dalcin       ierr    = PetscFree(indices);CHKERRQ(ierr);
148192bb6962SLisandro Dalcin       nidx   *= bs;
148292bb6962SLisandro Dalcin       indices = newidx;
148392bb6962SLisandro Dalcin     }
148492bb6962SLisandro Dalcin     /* Shift to get global indices */
148592bb6962SLisandro Dalcin     for (i=0; i<nidx; i++) indices[i] += rstart;
148692bb6962SLisandro Dalcin 
148792bb6962SLisandro Dalcin     /* Build the index sets for each block */
148892bb6962SLisandro Dalcin     for (i=0; i<n; i++) {
148970b3c8c7SBarry Smith       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
149092bb6962SLisandro Dalcin       ierr   = ISSort(is[i]);CHKERRQ(ierr);
149192bb6962SLisandro Dalcin       start += count[i];
149292bb6962SLisandro Dalcin     }
149392bb6962SLisandro Dalcin 
14943bf036e2SBarry Smith     ierr = PetscFree(count);CHKERRQ(ierr);
14953bf036e2SBarry Smith     ierr = PetscFree(indices);CHKERRQ(ierr);
1496fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1497fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
149892bb6962SLisandro Dalcin 
149992bb6962SLisandro Dalcin   }
150092bb6962SLisandro Dalcin   PetscFunctionReturn(0);
150192bb6962SLisandro Dalcin }
150292bb6962SLisandro Dalcin 
150392bb6962SLisandro Dalcin /*@C
150492bb6962SLisandro Dalcin    PCASMDestroySubdomains - Destroys the index sets created with
150592bb6962SLisandro Dalcin    PCASMCreateSubdomains(). Should be called after setting subdomains
150692bb6962SLisandro Dalcin    with PCASMSetLocalSubdomains().
150792bb6962SLisandro Dalcin 
150892bb6962SLisandro Dalcin    Collective
150992bb6962SLisandro Dalcin 
151092bb6962SLisandro Dalcin    Input Parameters:
151192bb6962SLisandro Dalcin +  n - the number of index sets
15122b691e39SMatthew Knepley .  is - the array of index sets
15130298fd71SBarry Smith -  is_local - the array of local index sets, can be NULL
151492bb6962SLisandro Dalcin 
151592bb6962SLisandro Dalcin    Level: advanced
151692bb6962SLisandro Dalcin 
151792bb6962SLisandro Dalcin .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains()
151892bb6962SLisandro Dalcin @*/
15197087cfbeSBarry Smith PetscErrorCode  PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
152092bb6962SLisandro Dalcin {
152192bb6962SLisandro Dalcin   PetscInt       i;
152292bb6962SLisandro Dalcin   PetscErrorCode ierr;
15235fd66863SKarl Rupp 
152492bb6962SLisandro Dalcin   PetscFunctionBegin;
1525a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1526a35b7b57SDmitry Karpeev   if (is) {
152792bb6962SLisandro Dalcin     PetscValidPointer(is,2);
1528fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); }
152992bb6962SLisandro Dalcin     ierr = PetscFree(is);CHKERRQ(ierr);
1530a35b7b57SDmitry Karpeev   }
15312b691e39SMatthew Knepley   if (is_local) {
15322b691e39SMatthew Knepley     PetscValidPointer(is_local,3);
1533fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); }
15342b691e39SMatthew Knepley     ierr = PetscFree(is_local);CHKERRQ(ierr);
15352b691e39SMatthew Knepley   }
153692bb6962SLisandro Dalcin   PetscFunctionReturn(0);
153792bb6962SLisandro Dalcin }
153892bb6962SLisandro Dalcin 
15394b9ad928SBarry Smith /*@
15404b9ad928SBarry Smith    PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
15414b9ad928SBarry Smith    preconditioner for a two-dimensional problem on a regular grid.
15424b9ad928SBarry Smith 
15434b9ad928SBarry Smith    Not Collective
15444b9ad928SBarry Smith 
15454b9ad928SBarry Smith    Input Parameters:
15464b9ad928SBarry Smith +  m, n - the number of mesh points in the x and y directions
15474b9ad928SBarry Smith .  M, N - the number of subdomains in the x and y directions
15484b9ad928SBarry Smith .  dof - degrees of freedom per node
15494b9ad928SBarry Smith -  overlap - overlap in mesh lines
15504b9ad928SBarry Smith 
15514b9ad928SBarry Smith    Output Parameters:
15524b9ad928SBarry Smith +  Nsub - the number of subdomains created
15533d03bb48SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
15543d03bb48SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
15554b9ad928SBarry Smith 
15564b9ad928SBarry Smith    Note:
15574b9ad928SBarry Smith    Presently PCAMSCreateSubdomains2d() is valid only for sequential
15584b9ad928SBarry Smith    preconditioners.  More general related routines are
15594b9ad928SBarry Smith    PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().
15604b9ad928SBarry Smith 
15614b9ad928SBarry Smith    Level: advanced
15624b9ad928SBarry Smith 
15634b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
15644b9ad928SBarry Smith           PCASMSetOverlap()
15654b9ad928SBarry Smith @*/
15667087cfbeSBarry Smith PetscErrorCode  PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local)
15674b9ad928SBarry Smith {
15683d03bb48SJed Brown   PetscInt       i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer;
15696849ba73SBarry Smith   PetscErrorCode ierr;
157013f74950SBarry Smith   PetscInt       nidx,*idx,loc,ii,jj,count;
15714b9ad928SBarry Smith 
15724b9ad928SBarry Smith   PetscFunctionBegin;
1573e32f2f54SBarry Smith   if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," ");
15744b9ad928SBarry Smith 
15754b9ad928SBarry Smith   *Nsub     = N*M;
1576854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is);CHKERRQ(ierr);
1577854ce69bSBarry Smith   ierr      = PetscMalloc1(*Nsub,is_local);CHKERRQ(ierr);
15784b9ad928SBarry Smith   ystart    = 0;
15793d03bb48SJed Brown   loc_outer = 0;
15804b9ad928SBarry Smith   for (i=0; i<N; i++) {
15814b9ad928SBarry Smith     height = n/N + ((n % N) > i); /* height of subdomain */
1582e32f2f54SBarry Smith     if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n");
15834b9ad928SBarry Smith     yleft  = ystart - overlap; if (yleft < 0) yleft = 0;
15844b9ad928SBarry Smith     yright = ystart + height + overlap; if (yright > n) yright = n;
15854b9ad928SBarry Smith     xstart = 0;
15864b9ad928SBarry Smith     for (j=0; j<M; j++) {
15874b9ad928SBarry Smith       width = m/M + ((m % M) > j); /* width of subdomain */
1588e32f2f54SBarry Smith       if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m");
15894b9ad928SBarry Smith       xleft  = xstart - overlap; if (xleft < 0) xleft = 0;
15904b9ad928SBarry Smith       xright = xstart + width + overlap; if (xright > m) xright = m;
15914b9ad928SBarry Smith       nidx   = (xright - xleft)*(yright - yleft);
1592785e854fSJed Brown       ierr   = PetscMalloc1(nidx,&idx);CHKERRQ(ierr);
15934b9ad928SBarry Smith       loc    = 0;
15944b9ad928SBarry Smith       for (ii=yleft; ii<yright; ii++) {
15954b9ad928SBarry Smith         count = m*ii + xleft;
15962fa5cd67SKarl Rupp         for (jj=xleft; jj<xright; jj++) idx[loc++] = count++;
15974b9ad928SBarry Smith       }
159870b3c8c7SBarry Smith       ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr);
15993d03bb48SJed Brown       if (overlap == 0) {
16003d03bb48SJed Brown         ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr);
16012fa5cd67SKarl Rupp 
16023d03bb48SJed Brown         (*is_local)[loc_outer] = (*is)[loc_outer];
16033d03bb48SJed Brown       } else {
16043d03bb48SJed Brown         for (loc=0,ii=ystart; ii<ystart+height; ii++) {
16053d03bb48SJed Brown           for (jj=xstart; jj<xstart+width; jj++) {
16063d03bb48SJed Brown             idx[loc++] = m*ii + jj;
16073d03bb48SJed Brown           }
16083d03bb48SJed Brown         }
160970b3c8c7SBarry Smith         ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr);
16103d03bb48SJed Brown       }
16114b9ad928SBarry Smith       ierr    = PetscFree(idx);CHKERRQ(ierr);
16124b9ad928SBarry Smith       xstart += width;
16133d03bb48SJed Brown       loc_outer++;
16144b9ad928SBarry Smith     }
16154b9ad928SBarry Smith     ystart += height;
16164b9ad928SBarry Smith   }
16174b9ad928SBarry Smith   for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); }
16184b9ad928SBarry Smith   PetscFunctionReturn(0);
16194b9ad928SBarry Smith }
16204b9ad928SBarry Smith 
16214b9ad928SBarry Smith /*@C
16224b9ad928SBarry Smith     PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
16234b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
16244b9ad928SBarry Smith 
1625ad4df100SBarry Smith     Not Collective
16264b9ad928SBarry Smith 
16274b9ad928SBarry Smith     Input Parameter:
16284b9ad928SBarry Smith .   pc - the preconditioner context
16294b9ad928SBarry Smith 
16304b9ad928SBarry Smith     Output Parameters:
16314b9ad928SBarry Smith +   n - the number of subdomains for this processor (default value = 1)
16322b691e39SMatthew Knepley .   is - the index sets that define the subdomains for this processor
16330298fd71SBarry Smith -   is_local - the index sets that define the local part of the subdomains for this processor (can be NULL)
16344b9ad928SBarry Smith 
16354b9ad928SBarry Smith 
16364b9ad928SBarry Smith     Notes:
16374b9ad928SBarry Smith     The IS numbering is in the parallel, global numbering of the vector.
16384b9ad928SBarry Smith 
16394b9ad928SBarry Smith     Level: advanced
16404b9ad928SBarry Smith 
16414b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
16424b9ad928SBarry Smith           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices()
16434b9ad928SBarry Smith @*/
16447087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
16454b9ad928SBarry Smith {
16462a808120SBarry Smith   PC_ASM         *osm = (PC_ASM*)pc->data;
164792bb6962SLisandro Dalcin   PetscErrorCode ierr;
1648ace3abfcSBarry Smith   PetscBool      match;
16494b9ad928SBarry Smith 
16504b9ad928SBarry Smith   PetscFunctionBegin;
16510700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
16524482741eSBarry Smith   PetscValidIntPointer(n,2);
165392bb6962SLisandro Dalcin   if (is) PetscValidPointer(is,3);
1654251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
16552a808120SBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM");
16564b9ad928SBarry Smith   if (n) *n = osm->n_local_true;
16574b9ad928SBarry Smith   if (is) *is = osm->is;
16582b691e39SMatthew Knepley   if (is_local) *is_local = osm->is_local;
16594b9ad928SBarry Smith   PetscFunctionReturn(0);
16604b9ad928SBarry Smith }
16614b9ad928SBarry Smith 
16624b9ad928SBarry Smith /*@C
16634b9ad928SBarry Smith     PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
16644b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
16654b9ad928SBarry Smith 
1666ad4df100SBarry Smith     Not Collective
16674b9ad928SBarry Smith 
16684b9ad928SBarry Smith     Input Parameter:
16694b9ad928SBarry Smith .   pc - the preconditioner context
16704b9ad928SBarry Smith 
16714b9ad928SBarry Smith     Output Parameters:
16724b9ad928SBarry Smith +   n - the number of matrices for this processor (default value = 1)
16734b9ad928SBarry Smith -   mat - the matrices
16744b9ad928SBarry Smith 
16754b9ad928SBarry Smith     Level: advanced
16764b9ad928SBarry Smith 
167795452b02SPatrick Sanan     Notes:
167834a84908Sprj-     Call after PCSetUp() (or KSPSetUp()) but before PCApply() and before PCSetUpOnBlocks())
1679cf739d55SBarry Smith 
168034a84908Sprj-            Usually one would use PCSetModifySubMatrices() to change the submatrices in building the preconditioner.
1681cf739d55SBarry Smith 
16824b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
168334a84908Sprj-           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubMatrices()
16844b9ad928SBarry Smith @*/
16857087cfbeSBarry Smith PetscErrorCode  PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
16864b9ad928SBarry Smith {
16874b9ad928SBarry Smith   PC_ASM         *osm;
168892bb6962SLisandro Dalcin   PetscErrorCode ierr;
1689ace3abfcSBarry Smith   PetscBool      match;
16904b9ad928SBarry Smith 
16914b9ad928SBarry Smith   PetscFunctionBegin;
16920700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
169392bb6962SLisandro Dalcin   PetscValidIntPointer(n,2);
169492bb6962SLisandro Dalcin   if (mat) PetscValidPointer(mat,3);
169534a84908Sprj-   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp().");
1696251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
169792bb6962SLisandro Dalcin   if (!match) {
169892bb6962SLisandro Dalcin     if (n) *n = 0;
16990298fd71SBarry Smith     if (mat) *mat = NULL;
170092bb6962SLisandro Dalcin   } else {
17014b9ad928SBarry Smith     osm = (PC_ASM*)pc->data;
17024b9ad928SBarry Smith     if (n) *n = osm->n_local_true;
17034b9ad928SBarry Smith     if (mat) *mat = osm->pmat;
170492bb6962SLisandro Dalcin   }
17054b9ad928SBarry Smith   PetscFunctionReturn(0);
17064b9ad928SBarry Smith }
1707d709ea83SDmitry Karpeev 
1708d709ea83SDmitry Karpeev /*@
1709d709ea83SDmitry Karpeev     PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1710f1ee410cSBarry Smith 
1711d709ea83SDmitry Karpeev     Logically Collective
1712d709ea83SDmitry Karpeev 
1713d709ea83SDmitry Karpeev     Input Parameter:
1714d709ea83SDmitry Karpeev +   pc  - the preconditioner
1715d709ea83SDmitry Karpeev -   flg - boolean indicating whether to use subdomains defined by the DM
1716d709ea83SDmitry Karpeev 
1717d709ea83SDmitry Karpeev     Options Database Key:
1718d709ea83SDmitry Karpeev .   -pc_asm_dm_subdomains
1719d709ea83SDmitry Karpeev 
1720d709ea83SDmitry Karpeev     Level: intermediate
1721d709ea83SDmitry Karpeev 
1722d709ea83SDmitry Karpeev     Notes:
1723d709ea83SDmitry Karpeev     PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(),
1724d709ea83SDmitry Karpeev     so setting either of the first two effectively turns the latter off.
1725d709ea83SDmitry Karpeev 
1726d709ea83SDmitry Karpeev .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1727d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1728d709ea83SDmitry Karpeev @*/
1729d709ea83SDmitry Karpeev PetscErrorCode  PCASMSetDMSubdomains(PC pc,PetscBool flg)
1730d709ea83SDmitry Karpeev {
1731d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1732d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1733d709ea83SDmitry Karpeev   PetscBool      match;
1734d709ea83SDmitry Karpeev 
1735d709ea83SDmitry Karpeev   PetscFunctionBegin;
1736d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1737d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
1738d709ea83SDmitry Karpeev   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1739d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1740d709ea83SDmitry Karpeev   if (match) {
1741d709ea83SDmitry Karpeev     osm->dm_subdomains = flg;
1742d709ea83SDmitry Karpeev   }
1743d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1744d709ea83SDmitry Karpeev }
1745d709ea83SDmitry Karpeev 
1746d709ea83SDmitry Karpeev /*@
1747d709ea83SDmitry Karpeev     PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1748d709ea83SDmitry Karpeev     Not Collective
1749d709ea83SDmitry Karpeev 
1750d709ea83SDmitry Karpeev     Input Parameter:
1751d709ea83SDmitry Karpeev .   pc  - the preconditioner
1752d709ea83SDmitry Karpeev 
1753d709ea83SDmitry Karpeev     Output Parameter:
1754d709ea83SDmitry Karpeev .   flg - boolean indicating whether to use subdomains defined by the DM
1755d709ea83SDmitry Karpeev 
1756d709ea83SDmitry Karpeev     Level: intermediate
1757d709ea83SDmitry Karpeev 
1758d709ea83SDmitry Karpeev .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1759d709ea83SDmitry Karpeev           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1760d709ea83SDmitry Karpeev @*/
1761d709ea83SDmitry Karpeev PetscErrorCode  PCASMGetDMSubdomains(PC pc,PetscBool* flg)
1762d709ea83SDmitry Karpeev {
1763d709ea83SDmitry Karpeev   PC_ASM         *osm = (PC_ASM*)pc->data;
1764d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1765d709ea83SDmitry Karpeev   PetscBool      match;
1766d709ea83SDmitry Karpeev 
1767d709ea83SDmitry Karpeev   PetscFunctionBegin;
1768d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1769534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg,2);
1770d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1771d709ea83SDmitry Karpeev   if (match) {
1772d709ea83SDmitry Karpeev     if (flg) *flg = osm->dm_subdomains;
1773d709ea83SDmitry Karpeev   }
1774d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1775d709ea83SDmitry Karpeev }
177680ec0b7dSPatrick Sanan 
177780ec0b7dSPatrick Sanan /*@
177880ec0b7dSPatrick Sanan      PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string.
177980ec0b7dSPatrick Sanan 
178080ec0b7dSPatrick Sanan    Not Collective
178180ec0b7dSPatrick Sanan 
178280ec0b7dSPatrick Sanan    Input Parameter:
178380ec0b7dSPatrick Sanan .  pc - the PC
178480ec0b7dSPatrick Sanan 
178580ec0b7dSPatrick Sanan    Output Parameter:
1786f1ee410cSBarry Smith .  -pc_asm_sub_mat_type - name of matrix type
178780ec0b7dSPatrick Sanan 
178880ec0b7dSPatrick Sanan    Level: advanced
178980ec0b7dSPatrick Sanan 
179080ec0b7dSPatrick Sanan .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
179180ec0b7dSPatrick Sanan @*/
179280ec0b7dSPatrick Sanan PetscErrorCode  PCASMGetSubMatType(PC pc,MatType *sub_mat_type){
179380ec0b7dSPatrick Sanan   PetscErrorCode ierr;
179480ec0b7dSPatrick Sanan 
179580ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));CHKERRQ(ierr);
179680ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
179780ec0b7dSPatrick Sanan }
179880ec0b7dSPatrick Sanan 
179980ec0b7dSPatrick Sanan /*@
180080ec0b7dSPatrick Sanan      PCASMSetSubMatType - Set the type of matrix used for ASM subsolves
180180ec0b7dSPatrick Sanan 
180280ec0b7dSPatrick Sanan    Collective on Mat
180380ec0b7dSPatrick Sanan 
180480ec0b7dSPatrick Sanan    Input Parameters:
180580ec0b7dSPatrick Sanan +  pc             - the PC object
180680ec0b7dSPatrick Sanan -  sub_mat_type   - matrix type
180780ec0b7dSPatrick Sanan 
180880ec0b7dSPatrick Sanan    Options Database Key:
180980ec0b7dSPatrick 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.
181080ec0b7dSPatrick Sanan 
181180ec0b7dSPatrick Sanan    Notes:
181280ec0b7dSPatrick Sanan    See "${PETSC_DIR}/include/petscmat.h" for available types
181380ec0b7dSPatrick Sanan 
181480ec0b7dSPatrick Sanan   Level: advanced
181580ec0b7dSPatrick Sanan 
181680ec0b7dSPatrick Sanan .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
181780ec0b7dSPatrick Sanan @*/
181880ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type)
181980ec0b7dSPatrick Sanan {
182080ec0b7dSPatrick Sanan   PetscErrorCode ierr;
182180ec0b7dSPatrick Sanan 
182280ec0b7dSPatrick Sanan   ierr = PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));CHKERRQ(ierr);
182380ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
182480ec0b7dSPatrick Sanan }
1825