xref: /petsc/src/ksp/pc/impls/gasm/gasm.c (revision 930d09c13c6be24c34efcb20ee858969724c0dd9)
1b1a0a8a3SJed Brown /*
2f746d493SDmitry Karpeev   This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation.
36a4f0f73SDmitry Karpeev   In this version each processor may intersect multiple subdomains and any subdomain may
46a4f0f73SDmitry Karpeev   intersect multiple processors.  Intersections of subdomains with processors are called *local
56a4f0f73SDmitry Karpeev   subdomains*.
6b1a0a8a3SJed Brown 
78f3b4b4dSDmitry Karpeev        N    - total number of distinct global subdomains          (set explicitly in PCGASMSetTotalSubdomains() or implicitly PCGASMSetSubdomains() and then calculated in PCSetUp_GASM())
86a4f0f73SDmitry Karpeev        n    - actual number of local subdomains on this processor (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains())
98f3b4b4dSDmitry Karpeev        nmax - maximum number of local subdomains per processor    (calculated in PCSetUp_GASM())
10b1a0a8a3SJed Brown */
11af0996ceSBarry Smith #include <petsc/private/pcimpl.h>     /*I "petscpc.h" I*/
121e25c274SJed Brown #include <petscdm.h>
13b1a0a8a3SJed Brown 
14b1a0a8a3SJed Brown typedef struct {
15f746d493SDmitry Karpeev   PetscInt    N,n,nmax;
16b1a0a8a3SJed Brown   PetscInt    overlap;                /* overlap requested by user */
178f3b4b4dSDmitry Karpeev   PCGASMType  type;                   /* use reduced interpolation, restriction or both */
188f3b4b4dSDmitry Karpeev   PetscBool   type_set;               /* if user set this value (so won't change it for symmetric problems) */
198f3b4b4dSDmitry Karpeev   PetscBool   same_subdomain_solvers; /* flag indicating whether all local solvers are same */
208f3b4b4dSDmitry Karpeev   PetscBool   sort_indices;           /* flag to sort subdomain indices */
218f3b4b4dSDmitry Karpeev   PetscBool   user_subdomains;        /* whether the user set explicit subdomain index sets -- keep them on PCReset() */
228f3b4b4dSDmitry Karpeev   PetscBool   dm_subdomains;          /* whether DM is allowed to define subdomains */
238f3b4b4dSDmitry Karpeev   IS          *ois;                   /* index sets that define the outer (conceptually, overlapping) subdomains */
248f3b4b4dSDmitry Karpeev   IS          *iis;                   /* index sets that define the inner (conceptually, nonoverlapping) subdomains */
258f3b4b4dSDmitry Karpeev   KSP         *ksp;                /* linear solvers for each subdomain */
268f3b4b4dSDmitry Karpeev   Mat         *pmat;               /* subdomain block matrices */
27f746d493SDmitry Karpeev   Vec         gx,gy;               /* Merged work vectors */
28f746d493SDmitry Karpeev   Vec         *x,*y;               /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */
296a4f0f73SDmitry Karpeev   VecScatter  gorestriction;       /* merged restriction to disjoint union of outer subdomains */
306a4f0f73SDmitry Karpeev   VecScatter  girestriction;       /* merged restriction to disjoint union of inner subdomains */
31f746d493SDmitry Karpeev } PC_GASM;
32b1a0a8a3SJed Brown 
33b1a0a8a3SJed Brown #undef __FUNCT__
348f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMComputeGlobalSubdomainNumbering_Private"
358f3b4b4dSDmitry Karpeev static PetscErrorCode  PCGASMComputeGlobalSubdomainNumbering_Private(PC pc,PetscInt **numbering,PetscInt **permutation)
368f3b4b4dSDmitry Karpeev {
378f3b4b4dSDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
388f3b4b4dSDmitry Karpeev   PetscInt       i;
398f3b4b4dSDmitry Karpeev   PetscErrorCode ierr;
408f3b4b4dSDmitry Karpeev 
418f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
428f3b4b4dSDmitry Karpeev   /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */
438f3b4b4dSDmitry Karpeev   ierr = PetscMalloc2(osm->n,numbering,osm->n,permutation);CHKERRQ(ierr);
448f3b4b4dSDmitry Karpeev   ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->iis,NULL,*numbering);CHKERRQ(ierr);
458f3b4b4dSDmitry Karpeev   for (i = 0; i < osm->n; ++i) (*permutation)[i] = i;
468f3b4b4dSDmitry Karpeev   ierr = PetscSortIntWithPermutation(osm->n,*numbering,*permutation);CHKERRQ(ierr);
478f3b4b4dSDmitry Karpeev   PetscFunctionReturn(0);
488f3b4b4dSDmitry Karpeev }
498f3b4b4dSDmitry Karpeev 
508f3b4b4dSDmitry Karpeev #undef __FUNCT__
5106b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSubdomainView_Private"
5206b43e7eSDmitry Karpeev static PetscErrorCode  PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer)
53af538404SDmitry Karpeev {
54af538404SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
5506b43e7eSDmitry Karpeev   PetscInt       j,nidx;
56af538404SDmitry Karpeev   const PetscInt *idx;
5706b43e7eSDmitry Karpeev   PetscViewer    sviewer;
58af538404SDmitry Karpeev   char           *cidx;
59af538404SDmitry Karpeev   PetscErrorCode ierr;
60af538404SDmitry Karpeev 
61af538404SDmitry Karpeev   PetscFunctionBegin;
62ce94432eSBarry Smith   if (i < 0 || i > osm->n) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Invalid subdomain %D: must nonnegative and less than %D", i, osm->n);
634bde246dSDmitry Karpeev   /* Inner subdomains. */
644bde246dSDmitry Karpeev   ierr = ISGetLocalSize(osm->iis[i], &nidx);CHKERRQ(ierr);
654bde246dSDmitry Karpeev   /*
664bde246dSDmitry Karpeev    No more than 15 characters per index plus a space.
674bde246dSDmitry Karpeev    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
684bde246dSDmitry Karpeev    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
694bde246dSDmitry Karpeev    For nidx == 0, the whole string 16 '\0'.
704bde246dSDmitry Karpeev    */
7189d949e2SBarry Smith   ierr = PetscMalloc1(16*(nidx+1)+1, &cidx);CHKERRQ(ierr);
724bde246dSDmitry Karpeev   ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr);
734bde246dSDmitry Karpeev   ierr = ISGetIndices(osm->iis[i], &idx);CHKERRQ(ierr);
744bde246dSDmitry Karpeev   for (j = 0; j < nidx; ++j) {
754bde246dSDmitry Karpeev     ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr);
764bde246dSDmitry Karpeev   }
774bde246dSDmitry Karpeev   ierr = ISRestoreIndices(osm->iis[i],&idx);CHKERRQ(ierr);
784bde246dSDmitry Karpeev   ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
794bde246dSDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n");CHKERRQ(ierr);
804bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
814bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
824bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
834bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
844bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
854bde246dSDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
864bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
874bde246dSDmitry Karpeev   ierr = PetscFree(cidx);CHKERRQ(ierr);
884bde246dSDmitry Karpeev   /* Outer subdomains. */
896a4f0f73SDmitry Karpeev   ierr = ISGetLocalSize(osm->ois[i], &nidx);CHKERRQ(ierr);
90eec7e3faSDmitry Karpeev   /*
91eec7e3faSDmitry Karpeev    No more than 15 characters per index plus a space.
92eec7e3faSDmitry Karpeev    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
93eec7e3faSDmitry Karpeev    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
94eec7e3faSDmitry Karpeev    For nidx == 0, the whole string 16 '\0'.
95eec7e3faSDmitry Karpeev    */
9689d949e2SBarry Smith   ierr = PetscMalloc1(16*(nidx+1)+1, &cidx);CHKERRQ(ierr);
9706b43e7eSDmitry Karpeev   ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr);
986a4f0f73SDmitry Karpeev   ierr = ISGetIndices(osm->ois[i], &idx);CHKERRQ(ierr);
99af538404SDmitry Karpeev   for (j = 0; j < nidx; ++j) {
100af538404SDmitry Karpeev     ierr = PetscViewerStringSPrintf(sviewer,"%D ", idx[j]);CHKERRQ(ierr);
101af538404SDmitry Karpeev   }
1026bf464f9SBarry Smith   ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
1036a4f0f73SDmitry Karpeev   ierr = ISRestoreIndices(osm->ois[i],&idx);CHKERRQ(ierr);
1046a4f0f73SDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n");CHKERRQ(ierr);
10506b43e7eSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1067b23a99aSBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
107af538404SDmitry Karpeev   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
108af538404SDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1097b23a99aSBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
110af538404SDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
11106b43e7eSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
112af538404SDmitry Karpeev   ierr = PetscFree(cidx);CHKERRQ(ierr);
11306b43e7eSDmitry Karpeev   PetscFunctionReturn(0);
11406b43e7eSDmitry Karpeev }
11506b43e7eSDmitry Karpeev 
11606b43e7eSDmitry Karpeev #undef __FUNCT__
11706b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains"
11806b43e7eSDmitry Karpeev static PetscErrorCode  PCGASMPrintSubdomains(PC pc)
11906b43e7eSDmitry Karpeev {
12006b43e7eSDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
12106b43e7eSDmitry Karpeev   const char     *prefix;
12206b43e7eSDmitry Karpeev   char           fname[PETSC_MAX_PATH_LEN+1];
1238f3b4b4dSDmitry Karpeev   PetscInt       l, d, count;
1248f3b4b4dSDmitry Karpeev   PetscBool      doprint,found;
1250298fd71SBarry Smith   PetscViewer    viewer, sviewer = NULL;
1268f3b4b4dSDmitry Karpeev   PetscInt       *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */
12706b43e7eSDmitry Karpeev   PetscErrorCode ierr;
12806b43e7eSDmitry Karpeev 
12906b43e7eSDmitry Karpeev   PetscFunctionBegin;
13006b43e7eSDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
1318f3b4b4dSDmitry Karpeev   doprint  = PETSC_FALSE;
1328f3b4b4dSDmitry Karpeev   ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&doprint,NULL);CHKERRQ(ierr);
1338f3b4b4dSDmitry Karpeev   if (!doprint) PetscFunctionReturn(0);
13406b43e7eSDmitry Karpeev   ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr);
13506b43e7eSDmitry Karpeev   if (!found) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
136ce94432eSBarry Smith   ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr);
1374bde246dSDmitry Karpeev   /*
1384bde246dSDmitry Karpeev    Make sure the viewer has a name. Otherwise this may cause a deadlock or other weird errors when creating a subcomm viewer:
1394bde246dSDmitry Karpeev    the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed collectively on the comm.
1404bde246dSDmitry Karpeev   */
1414bde246dSDmitry Karpeev   ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr);
1424bde246dSDmitry Karpeev   l    = 0;
1438f3b4b4dSDmitry Karpeev   ierr = PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);CHKERRQ(ierr);
1448f3b4b4dSDmitry Karpeev   for (count = 0; count < osm->N; ++count) {
1454bde246dSDmitry Karpeev     /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
1464bde246dSDmitry Karpeev     if (l<osm->n) {
1474bde246dSDmitry Karpeev       d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
1484bde246dSDmitry Karpeev       if (numbering[d] == count) {
1494bde246dSDmitry Karpeev         ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
1504bde246dSDmitry Karpeev         ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr);
1514bde246dSDmitry Karpeev         ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
1524bde246dSDmitry Karpeev         ++l;
153af538404SDmitry Karpeev       }
1544bde246dSDmitry Karpeev     }
155ce94432eSBarry Smith     ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
1564bde246dSDmitry Karpeev   }
1578f3b4b4dSDmitry Karpeev   ierr = PetscFree2(numbering,permutation);CHKERRQ(ierr);
1584bde246dSDmitry Karpeev   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
159af538404SDmitry Karpeev   PetscFunctionReturn(0);
160af538404SDmitry Karpeev }
161af538404SDmitry Karpeev 
162af538404SDmitry Karpeev 
163af538404SDmitry Karpeev #undef __FUNCT__
164f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM"
165f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer)
166b1a0a8a3SJed Brown {
167f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
168af538404SDmitry Karpeev   const char     *prefix;
169b1a0a8a3SJed Brown   PetscErrorCode ierr;
170af538404SDmitry Karpeev   PetscMPIInt    rank, size;
1718f3b4b4dSDmitry Karpeev   PetscInt       bsz;
17206b43e7eSDmitry Karpeev   PetscBool      iascii,view_subdomains=PETSC_FALSE;
173b1a0a8a3SJed Brown   PetscViewer    sviewer;
1748f3b4b4dSDmitry Karpeev   PetscInt       count, l;
1756a4f0f73SDmitry Karpeev   char           overlap[256]     = "user-defined overlap";
1766a4f0f73SDmitry Karpeev   char           gsubdomains[256] = "unknown total number of subdomains";
17706b43e7eSDmitry Karpeev   char           msubdomains[256] = "unknown max number of local subdomains";
1788f3b4b4dSDmitry Karpeev   PetscInt       *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */
17911aeaf0aSBarry Smith 
180b1a0a8a3SJed Brown   PetscFunctionBegin;
181ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr);
182ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr);
18306b43e7eSDmitry Karpeev 
18406b43e7eSDmitry Karpeev   if (osm->overlap >= 0) {
18506b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %D",osm->overlap);CHKERRQ(ierr);
18606b43e7eSDmitry Karpeev   }
1878f3b4b4dSDmitry Karpeev   if (osm->N != PETSC_DETERMINE) {
1888f3b4b4dSDmitry Karpeev     ierr = PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %D",osm->N);CHKERRQ(ierr);
18906b43e7eSDmitry Karpeev   }
1908f3b4b4dSDmitry Karpeev   if (osm->nmax != PETSC_DETERMINE) {
19106b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr);
19206b43e7eSDmitry Karpeev   }
19306b43e7eSDmitry Karpeev 
194af538404SDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
1950298fd71SBarry Smith   ierr = PetscOptionsGetBool(prefix,"-pc_gasm_view_subdomains",&view_subdomains,NULL);CHKERRQ(ierr);
19606b43e7eSDmitry Karpeev 
19706b43e7eSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
198b1a0a8a3SJed Brown   if (iascii) {
19906b43e7eSDmitry Karpeev     /*
20006b43e7eSDmitry Karpeev      Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer:
20106b43e7eSDmitry Karpeev      the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed
20206b43e7eSDmitry Karpeev      collectively on the comm.
20306b43e7eSDmitry Karpeev      */
20406b43e7eSDmitry Karpeev     ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr);
20506b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Generalized additive Schwarz:\n");CHKERRQ(ierr);
20606b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Restriction/interpolation type: %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr);
20706b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",overlap);CHKERRQ(ierr);
20806b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",gsubdomains);CHKERRQ(ierr);
20906b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",msubdomains);CHKERRQ(ierr);
2107b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
2118f3b4b4dSDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d|%d] number of locally-supported subdomains = %D\n",rank,size,osm->n);CHKERRQ(ierr);
212af538404SDmitry Karpeev     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
2137b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
2146a4f0f73SDmitry Karpeev     /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */
21506b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Subdomain solver info is as follows:\n");CHKERRQ(ierr);
216b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
217b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
21806b43e7eSDmitry Karpeev     /* Make sure that everybody waits for the banner to be printed. */
219ce94432eSBarry Smith     ierr = MPI_Barrier(PetscObjectComm((PetscObject)viewer));CHKERRQ(ierr);
2204bde246dSDmitry Karpeev     /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
2218f3b4b4dSDmitry Karpeev     ierr = PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);CHKERRQ(ierr);
22206b43e7eSDmitry Karpeev     l = 0;
2238f3b4b4dSDmitry Karpeev     for (count = 0; count < osm->N; ++count) {
22406b43e7eSDmitry Karpeev       PetscMPIInt srank, ssize;
22506b43e7eSDmitry Karpeev       if (l<osm->n) {
22606b43e7eSDmitry Karpeev         PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
22706b43e7eSDmitry Karpeev         if (numbering[d] == count) {
2286a4f0f73SDmitry Karpeev           ierr = MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize);CHKERRQ(ierr);
2296a4f0f73SDmitry Karpeev           ierr = MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank);CHKERRQ(ierr);
2306a4f0f73SDmitry Karpeev           ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
2316a4f0f73SDmitry Karpeev           ierr = ISGetLocalSize(osm->ois[d],&bsz);CHKERRQ(ierr);
2327b23a99aSBarry Smith           ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_TRUE);CHKERRQ(ierr);
2338f3b4b4dSDmitry Karpeev           ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d|%d] (subcomm [%d|%d]) local subdomain number %D, local size = %D\n",rank,size,srank,ssize,d,bsz);CHKERRQ(ierr);
2346a4f0f73SDmitry Karpeev           ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
2356a4f0f73SDmitry Karpeev           ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_FALSE);CHKERRQ(ierr);
23606b43e7eSDmitry Karpeev           if (view_subdomains) {
23706b43e7eSDmitry Karpeev             ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr);
23806b43e7eSDmitry Karpeev           }
2396a4f0f73SDmitry Karpeev           if (!pc->setupcalled) {
2406a4f0f73SDmitry Karpeev             PetscViewerASCIIPrintf(sviewer, "Solver not set up yet: PCSetUp() not yet called\n");CHKERRQ(ierr);
2412fa5cd67SKarl Rupp           } else {
24206b43e7eSDmitry Karpeev             ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr);
2436a4f0f73SDmitry Karpeev           }
24406b43e7eSDmitry Karpeev           ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
2457b23a99aSBarry Smith           ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
2466a4f0f73SDmitry Karpeev           ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
24706b43e7eSDmitry Karpeev           ++l;
248b1a0a8a3SJed Brown         }
24906b43e7eSDmitry Karpeev       }
250ce94432eSBarry Smith       ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
251b1a0a8a3SJed Brown     }
2528f3b4b4dSDmitry Karpeev     ierr = PetscFree2(numbering,permutation);CHKERRQ(ierr);
253b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
254b1a0a8a3SJed Brown   }
255b1a0a8a3SJed Brown   PetscFunctionReturn(0);
256b1a0a8a3SJed Brown }
257b1a0a8a3SJed Brown 
2588f3b4b4dSDmitry Karpeev PETSC_INTERN PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]);
259af538404SDmitry Karpeev 
260b1a0a8a3SJed Brown #undef __FUNCT__
261f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM"
262f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc)
263b1a0a8a3SJed Brown {
264f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
265b1a0a8a3SJed Brown   PetscErrorCode ierr;
266b1a0a8a3SJed Brown   PetscBool      symset,flg;
26788389755SJed Brown   PetscInt       i;
268c10bc1a1SDmitry Karpeev   PetscMPIInt    rank, size;
269b1a0a8a3SJed Brown   MatReuse       scall = MAT_REUSE_MATRIX;
270b1a0a8a3SJed Brown   KSP            ksp;
271b1a0a8a3SJed Brown   PC             subpc;
272b1a0a8a3SJed Brown   const char     *prefix,*pprefix;
273f746d493SDmitry Karpeev   Vec            x,y;
2746a4f0f73SDmitry Karpeev   PetscInt       oni;       /* Number of indices in the i-th local outer subdomain.               */
2756a4f0f73SDmitry Karpeev   const PetscInt *oidxi;    /* Indices from the i-th subdomain local outer subdomain.             */
2766a4f0f73SDmitry Karpeev   PetscInt       on;        /* Number of indices in the disjoint union of local outer subdomains. */
2776a4f0f73SDmitry Karpeev   PetscInt       *oidx;     /* Indices in the disjoint union of local outer subdomains. */
2786a4f0f73SDmitry Karpeev   IS             gois;      /* Disjoint union the global indices of outer subdomains.             */
2796a4f0f73SDmitry Karpeev   IS             goid;      /* Identity IS of the size of the disjoint union of outer subdomains. */
280f746d493SDmitry Karpeev   PetscScalar    *gxarray, *gyarray;
281*930d09c1SFande Kong   PetscInt       gostart;   /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */
2828f3b4b4dSDmitry Karpeev   PetscInt       num_subdomains    = 0;
2830298fd71SBarry Smith   DM             *subdomain_dm     = NULL;
2848f3b4b4dSDmitry Karpeev   char           **subdomain_names = NULL;
2858f3b4b4dSDmitry Karpeev 
286b1a0a8a3SJed Brown 
287b1a0a8a3SJed Brown   PetscFunctionBegin;
288ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
289ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
290b1a0a8a3SJed Brown   if (!pc->setupcalled) {
291b1a0a8a3SJed Brown 
292b1a0a8a3SJed Brown     if (!osm->type_set) {
293b1a0a8a3SJed Brown       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
2942fa5cd67SKarl Rupp       if (symset && flg) osm->type = PC_GASM_BASIC;
295b1a0a8a3SJed Brown     }
296b1a0a8a3SJed Brown 
2978f3b4b4dSDmitry Karpeev     if (osm->n == PETSC_DETERMINE) {
2988f3b4b4dSDmitry Karpeev       if (osm->N != PETSC_DETERMINE) {
2998f3b4b4dSDmitry Karpeev 	   /* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */
3008f3b4b4dSDmitry Karpeev 	   ierr = PCGASMCreateSubdomains(pc->pmat,osm->N,&osm->n,&osm->iis);CHKERRQ(ierr);
3018f3b4b4dSDmitry Karpeev       } else if (osm->dm_subdomains && pc->dm) {
3028f3b4b4dSDmitry Karpeev 	/* try pc->dm next, if allowed */
3038f3b4b4dSDmitry Karpeev 	PetscInt  d;
304a35b7b57SDmitry Karpeev 	IS       *inner_subdomain_is, *outer_subdomain_is;
305a35b7b57SDmitry Karpeev 	ierr = DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);CHKERRQ(ierr);
306a35b7b57SDmitry Karpeev 	if (num_subdomains) {
307a35b7b57SDmitry Karpeev 	  ierr = PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);CHKERRQ(ierr);
30869ca1f37SDmitry Karpeev 	}
309a35b7b57SDmitry Karpeev 	for (d = 0; d < num_subdomains; ++d) {
310a35b7b57SDmitry Karpeev 	  if (inner_subdomain_is) {ierr = ISDestroy(&inner_subdomain_is[d]);CHKERRQ(ierr);}
311a35b7b57SDmitry Karpeev 	  if (outer_subdomain_is) {ierr = ISDestroy(&outer_subdomain_is[d]);CHKERRQ(ierr);}
312fdc48646SDmitry Karpeev 	}
313a35b7b57SDmitry Karpeev 	ierr = PetscFree(inner_subdomain_is);CHKERRQ(ierr);
314a35b7b57SDmitry Karpeev 	ierr = PetscFree(outer_subdomain_is);CHKERRQ(ierr);
3158f3b4b4dSDmitry Karpeev       } else {
3168f3b4b4dSDmitry Karpeev 	/* still no subdomains; use one per processor */
31744fe18e5SDmitry Karpeev 	osm->nmax = osm->n = 1;
318ce94432eSBarry Smith 	ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
319f746d493SDmitry Karpeev 	osm->N    = size;
3208f3b4b4dSDmitry Karpeev 	ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr);
321fdc48646SDmitry Karpeev       }
32206b43e7eSDmitry Karpeev     }
323a35b7b57SDmitry Karpeev     if (!osm->iis) {
324a35b7b57SDmitry Karpeev       /*
3258f3b4b4dSDmitry Karpeev        osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied.
3268f3b4b4dSDmitry Karpeev        We create the requisite number of local inner subdomains and then expand them into
3278f3b4b4dSDmitry Karpeev        out subdomains, if necessary.
328a35b7b57SDmitry Karpeev        */
3298f3b4b4dSDmitry Karpeev       ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr);
330a35b7b57SDmitry Karpeev     }
3318f3b4b4dSDmitry Karpeev     if (!osm->ois) {
3328f3b4b4dSDmitry Karpeev       /*
3338f3b4b4dSDmitry Karpeev 	Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap
3348f3b4b4dSDmitry Karpeev 	has been requested, copy the inner subdomains over so they can be modified.
3358f3b4b4dSDmitry Karpeev       */
3368f3b4b4dSDmitry Karpeev       ierr = PetscMalloc1(osm->n,&osm->ois);CHKERRQ(ierr);
3378f3b4b4dSDmitry Karpeev       for (i=0; i<osm->n; ++i) {
3388f3b4b4dSDmitry Karpeev 	if (osm->overlap > 0) { /* With positive overlap, osm->iis[i] will be modified */
3398f3b4b4dSDmitry Karpeev 	  ierr = ISDuplicate(osm->iis[i],(osm->ois)+i);CHKERRQ(ierr);
3408f3b4b4dSDmitry Karpeev 	  ierr = ISCopy(osm->iis[i],osm->ois[i]);CHKERRQ(ierr);
3418f3b4b4dSDmitry Karpeev 	} else {
3428f3b4b4dSDmitry Karpeev 	  ierr      = PetscObjectReference((PetscObject)((osm->iis)[i]));CHKERRQ(ierr);
3438f3b4b4dSDmitry Karpeev 	  osm->ois[i] = osm->iis[i];
3448f3b4b4dSDmitry Karpeev 	}
3458f3b4b4dSDmitry Karpeev       }
3468f3b4b4dSDmitry Karpeev       if (osm->overlap > 0) {
3478f3b4b4dSDmitry Karpeev 	   /* Extend the "overlapping" regions by a number of steps */
348d21f2a47SFande Kong 	   ierr = MatIncreaseOverlapSplit(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr);
3498f3b4b4dSDmitry Karpeev       }
350b1a0a8a3SJed Brown     }
3516a4f0f73SDmitry Karpeev 
3528f3b4b4dSDmitry Karpeev     /* Now the subdomains are defined.  Determine their global and max local numbers, if necessary. */
3538f3b4b4dSDmitry Karpeev     if (osm->nmax == PETSC_DETERMINE) {
3548f3b4b4dSDmitry Karpeev       PetscMPIInt inwork,outwork;
3558f3b4b4dSDmitry Karpeev       /* determine global number of subdomains and the max number of local subdomains */
3568f3b4b4dSDmitry Karpeev       inwork = osm->n;
3578f3b4b4dSDmitry Karpeev       ierr       = MPI_Allreduce(&inwork,&outwork,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
3588f3b4b4dSDmitry Karpeev       osm->nmax  = outwork;
3598f3b4b4dSDmitry Karpeev     }
3608f3b4b4dSDmitry Karpeev     if (osm->N == PETSC_DETERMINE) {
3618f3b4b4dSDmitry Karpeev       /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */
3628f3b4b4dSDmitry Karpeev       ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,&osm->N,NULL);CHKERRQ(ierr);
3638f3b4b4dSDmitry Karpeev     }
3648f3b4b4dSDmitry Karpeev 
365b1a0a8a3SJed Brown 
366b1a0a8a3SJed Brown     if (osm->sort_indices) {
367f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
3686a4f0f73SDmitry Karpeev         ierr = ISSort(osm->ois[i]);CHKERRQ(ierr);
3696a4f0f73SDmitry Karpeev         ierr = ISSort(osm->iis[i]);CHKERRQ(ierr);
370b1a0a8a3SJed Brown       }
371b1a0a8a3SJed Brown     }
372*930d09c1SFande Kong #if 1
373*930d09c1SFande Kong     for (i=0; i<osm->n; i++){
374*930d09c1SFande Kong        ierr = ISView(osm->ois[i],NULL);CHKERRQ(ierr);
375*930d09c1SFande Kong        ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRQ(ierr);
376*930d09c1SFande Kong        ierr = ISView(osm->iis[i],NULL);CHKERRQ(ierr);
377*930d09c1SFande Kong        ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRQ(ierr);
378*930d09c1SFande Kong     }
379*930d09c1SFande Kong #endif
3808f3b4b4dSDmitry Karpeev     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
3818f3b4b4dSDmitry Karpeev     ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr);
3828f3b4b4dSDmitry Karpeev 
3836a4f0f73SDmitry Karpeev     /*
3846a4f0f73SDmitry Karpeev      Merge the ISs, create merged vectors and restrictions.
3856a4f0f73SDmitry Karpeev      */
3866a4f0f73SDmitry Karpeev     /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */
3876a4f0f73SDmitry Karpeev     on = 0;
388f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
3896a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
3906a4f0f73SDmitry Karpeev       on  += oni;
391f746d493SDmitry Karpeev     }
392785e854fSJed Brown     ierr = PetscMalloc1(on, &oidx);CHKERRQ(ierr);
3936a4f0f73SDmitry Karpeev     on   = 0;
394*930d09c1SFande Kong     /*Merge local indices together */
395f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
3966a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
3976a4f0f73SDmitry Karpeev       ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr);
3986a4f0f73SDmitry Karpeev       ierr = PetscMemcpy(oidx+on, oidxi, sizeof(PetscInt)*oni);CHKERRQ(ierr);
3996a4f0f73SDmitry Karpeev       ierr = ISRestoreIndices(osm->ois[i], &oidxi);CHKERRQ(ierr);
4006a4f0f73SDmitry Karpeev       on  += oni;
401f746d493SDmitry Karpeev     }
4026a4f0f73SDmitry Karpeev     ierr = ISCreateGeneral(((PetscObject)(pc))->comm,on,oidx,PETSC_OWN_POINTER,&gois);CHKERRQ(ierr);
4032a7a6963SBarry Smith     ierr = MatCreateVecs(pc->pmat,&x,&y);CHKERRQ(ierr);
404ce94432eSBarry Smith     ierr = VecCreateMPI(PetscObjectComm((PetscObject)pc),on,PETSC_DECIDE,&osm->gx);CHKERRQ(ierr);
405f746d493SDmitry Karpeev     ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr);
406*930d09c1SFande Kong     ierr = VecGetOwnershipRange(osm->gx, &gostart, NULL);CHKERRQ(ierr);
407*930d09c1SFande Kong     ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),on,gostart,1, &goid);CHKERRQ(ierr);
408*930d09c1SFande Kong     /*gois might indices not on local */
4096a4f0f73SDmitry Karpeev     ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr);
4106a4f0f73SDmitry Karpeev     ierr = VecDestroy(&x);CHKERRQ(ierr);
4117105d80bSDmitry Karpeev     ierr = ISDestroy(&gois);CHKERRQ(ierr);
4122fa5cd67SKarl Rupp 
4136a4f0f73SDmitry Karpeev     /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */
4142fa5cd67SKarl Rupp     {
4152fa5cd67SKarl Rupp       PetscInt  ini;           /* Number of indices the i-th a local inner subdomain. */
4166a4f0f73SDmitry Karpeev       PetscInt  in;            /* Number of indices in the disjoint uniont of local inner subdomains. */
4176a4f0f73SDmitry Karpeev       PetscInt *iidx;         /* Global indices in the merged local inner subdomain. */
4186a4f0f73SDmitry Karpeev       PetscInt *ioidx;        /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
4196a4f0f73SDmitry Karpeev       IS        giis;          /* IS for the disjoint union of inner subdomains. */
4206a4f0f73SDmitry Karpeev       IS        giois;         /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
421f746d493SDmitry Karpeev       /**/
4226a4f0f73SDmitry Karpeev       in = 0;
423f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
4246a4f0f73SDmitry Karpeev         ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr);
4256a4f0f73SDmitry Karpeev         in  += ini;
426f746d493SDmitry Karpeev       }
427785e854fSJed Brown       ierr = PetscMalloc1(in, &iidx);CHKERRQ(ierr);
428785e854fSJed Brown       ierr = PetscMalloc1(in, &ioidx);CHKERRQ(ierr);
429*930d09c1SFande Kong       ierr = VecGetOwnershipRange(osm->gx,&gostart, NULL);CHKERRQ(ierr);
4306a4f0f73SDmitry Karpeev       in   = 0;
4316a4f0f73SDmitry Karpeev       on   = 0;
432f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
4336a4f0f73SDmitry Karpeev         const PetscInt         *iidxi; /* Global indices of the i-th local inner subdomain. */
4346a4f0f73SDmitry Karpeev         ISLocalToGlobalMapping ltogi; /* Map from global to local indices of the i-th outer local subdomain. */
4356a4f0f73SDmitry Karpeev         PetscInt               *ioidxi; /* Local indices of the i-th local inner subdomain within the local outer subdomain. */
4366a4f0f73SDmitry Karpeev         PetscInt               ioni;  /* Number of indices in ioidxi; if ioni != ini the inner subdomain is not a subdomain of the outer subdomain (error). */
4376a4f0f73SDmitry Karpeev         PetscInt               k;
4386a4f0f73SDmitry Karpeev         ierr   = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr);
4396a4f0f73SDmitry Karpeev         ierr   = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
4406a4f0f73SDmitry Karpeev         ierr   = ISGetIndices(osm->iis[i],&iidxi);CHKERRQ(ierr);
4416a4f0f73SDmitry Karpeev         ierr   = PetscMemcpy(iidx+in, iidxi, sizeof(PetscInt)*ini);CHKERRQ(ierr);
4426a4f0f73SDmitry Karpeev         ierr   = ISLocalToGlobalMappingCreateIS(osm->ois[i],&ltogi);CHKERRQ(ierr);
4436a4f0f73SDmitry Karpeev         ioidxi = ioidx+in;
444*930d09c1SFande Kong         /*This function is not scalable */
4456a4f0f73SDmitry Karpeev         ierr   = ISGlobalToLocalMappingApply(ltogi,IS_GTOLM_DROP,ini,iidxi,&ioni,ioidxi);CHKERRQ(ierr);
4466a4f0f73SDmitry Karpeev         ierr   = ISLocalToGlobalMappingDestroy(&ltogi);CHKERRQ(ierr);
4476a4f0f73SDmitry Karpeev         ierr   = ISRestoreIndices(osm->iis[i], &iidxi);CHKERRQ(ierr);
4486a4f0f73SDmitry Karpeev         if (ioni != ini) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inner subdomain %D contains %D indices outside of its outer subdomain", i, ini - ioni);
449*930d09c1SFande Kong         for (k = 0; k < ini; ++k) ioidxi[k] += gostart+on;
4506a4f0f73SDmitry Karpeev         in += ini;
4516a4f0f73SDmitry Karpeev         on += oni;
452f746d493SDmitry Karpeev       }
453ce94432eSBarry Smith       ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, iidx,  PETSC_OWN_POINTER, &giis);CHKERRQ(ierr);
454ce94432eSBarry Smith       ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, ioidx, PETSC_OWN_POINTER, &giois);CHKERRQ(ierr);
4556a4f0f73SDmitry Karpeev       ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr);
4566a4f0f73SDmitry Karpeev       ierr = VecDestroy(&y);CHKERRQ(ierr);
4576a4f0f73SDmitry Karpeev       ierr = ISDestroy(&giis);CHKERRQ(ierr);
4586a4f0f73SDmitry Karpeev       ierr = ISDestroy(&giois);CHKERRQ(ierr);
459b1a0a8a3SJed Brown     }
4606a4f0f73SDmitry Karpeev     ierr = ISDestroy(&goid);CHKERRQ(ierr);
4612fa5cd67SKarl Rupp 
462f746d493SDmitry Karpeev     /* Create the subdomain work vectors. */
463785e854fSJed Brown     ierr = PetscMalloc1(osm->n,&osm->x);CHKERRQ(ierr);
464785e854fSJed Brown     ierr = PetscMalloc1(osm->n,&osm->y);CHKERRQ(ierr);
465f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr);
466f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr);
4676a4f0f73SDmitry Karpeev     for (i=0, on=0; i<osm->n; ++i, on += oni) {
4686a4f0f73SDmitry Karpeev       PetscInt oNi;
4696a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
470*930d09c1SFande Kong       /*on a sub communicator */
4716a4f0f73SDmitry Karpeev       ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr);
4726a4f0f73SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr);
4736a4f0f73SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr);
474b1a0a8a3SJed Brown     }
475f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr);
476f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr);
4778f3b4b4dSDmitry Karpeev     /* Create the subdomain solvers */
478785e854fSJed Brown     ierr = PetscMalloc1(osm->n,&osm->ksp);CHKERRQ(ierr);
4798f3b4b4dSDmitry Karpeev     for (i=0; i<osm->n; i++) {
4808f3b4b4dSDmitry Karpeev       char subprefix[PETSC_MAX_PATH_LEN+1];
4816a4f0f73SDmitry Karpeev       ierr        = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr);
4823bb1ff40SBarry Smith       ierr        = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr);
483b1a0a8a3SJed Brown       ierr        = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
484b1a0a8a3SJed Brown       ierr        = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
4858f3b4b4dSDmitry Karpeev       ierr        = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); /* Why do we need this here? */
4868f3b4b4dSDmitry Karpeev       if (subdomain_dm) {
4878f3b4b4dSDmitry Karpeev 	    ierr = KSPSetDM(ksp,subdomain_dm[i]);CHKERRQ(ierr);
4888f3b4b4dSDmitry Karpeev 	    ierr = DMDestroy(subdomain_dm+i);CHKERRQ(ierr);
4898f3b4b4dSDmitry Karpeev       }
490b1a0a8a3SJed Brown       ierr        = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
491b1a0a8a3SJed Brown       ierr        = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
4928f3b4b4dSDmitry Karpeev       if (subdomain_names && subdomain_names[i]) {
4938f3b4b4dSDmitry Karpeev 	     ierr = PetscSNPrintf(subprefix,PETSC_MAX_PATH_LEN,"sub_%s_",subdomain_names[i]);CHKERRQ(ierr);
4948f3b4b4dSDmitry Karpeev 	     ierr = KSPAppendOptionsPrefix(ksp,subprefix);CHKERRQ(ierr);
4958f3b4b4dSDmitry Karpeev 	     ierr = PetscFree(subdomain_names[i]);CHKERRQ(ierr);
4968f3b4b4dSDmitry Karpeev       }
497b1a0a8a3SJed Brown       ierr        = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
498b1a0a8a3SJed Brown       osm->ksp[i] = ksp;
499b1a0a8a3SJed Brown     }
5008f3b4b4dSDmitry Karpeev     ierr = PetscFree(subdomain_dm);CHKERRQ(ierr);
5018f3b4b4dSDmitry Karpeev     ierr = PetscFree(subdomain_names);CHKERRQ(ierr);
502b1a0a8a3SJed Brown     scall = MAT_INITIAL_MATRIX;
503b1a0a8a3SJed Brown 
5048f3b4b4dSDmitry Karpeev   } else { /*if (pc->setupcalled)*/
505b1a0a8a3SJed Brown     /*
5068f3b4b4dSDmitry Karpeev        Destroy the submatrices from the previous iteration
507b1a0a8a3SJed Brown     */
508b1a0a8a3SJed Brown     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
509f746d493SDmitry Karpeev       ierr  = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
510b1a0a8a3SJed Brown       scall = MAT_INITIAL_MATRIX;
511b1a0a8a3SJed Brown     }
512b1a0a8a3SJed Brown   }
513b1a0a8a3SJed Brown 
514b1a0a8a3SJed Brown   /*
515f746d493SDmitry Karpeev      Extract out the submatrices.
516b1a0a8a3SJed Brown   */
51781a5c4aaSDmitry Karpeev   if (size > 1) {
5188f3b4b4dSDmitry Karpeev     ierr = MatGetSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr);
5192fa5cd67SKarl Rupp   } else {
5206a4f0f73SDmitry Karpeev     ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr);
52181a5c4aaSDmitry Karpeev   }
522b1a0a8a3SJed Brown   if (scall == MAT_INITIAL_MATRIX) {
523b1a0a8a3SJed Brown     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
524f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
5253bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr);
526b1a0a8a3SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
527b1a0a8a3SJed Brown     }
528b1a0a8a3SJed Brown   }
529b1a0a8a3SJed Brown 
530b1a0a8a3SJed Brown   /* Return control to the user so that the submatrices can be modified (e.g., to apply
531b1a0a8a3SJed Brown      different boundary conditions for the submatrices than for the global problem) */
5326a4f0f73SDmitry Karpeev   ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
533b1a0a8a3SJed Brown 
534b1a0a8a3SJed Brown   /*
5356a4f0f73SDmitry Karpeev      Loop over submatrices putting them into local ksps
536b1a0a8a3SJed Brown   */
537f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
53823ee1639SBarry Smith     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr);
539b1a0a8a3SJed Brown     if (!pc->setupcalled) {
540b1a0a8a3SJed Brown       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
541b1a0a8a3SJed Brown     }
542b1a0a8a3SJed Brown   }
543b1a0a8a3SJed Brown   PetscFunctionReturn(0);
544b1a0a8a3SJed Brown }
545b1a0a8a3SJed Brown 
546b1a0a8a3SJed Brown #undef __FUNCT__
547f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM"
548f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
549b1a0a8a3SJed Brown {
550f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
551b1a0a8a3SJed Brown   PetscErrorCode ierr;
552b1a0a8a3SJed Brown   PetscInt       i;
553b1a0a8a3SJed Brown 
554b1a0a8a3SJed Brown   PetscFunctionBegin;
555f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
556b1a0a8a3SJed Brown     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
557b1a0a8a3SJed Brown   }
558b1a0a8a3SJed Brown   PetscFunctionReturn(0);
559b1a0a8a3SJed Brown }
560b1a0a8a3SJed Brown 
561b1a0a8a3SJed Brown #undef __FUNCT__
562f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM"
563f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y)
564b1a0a8a3SJed Brown {
565f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
566b1a0a8a3SJed Brown   PetscErrorCode ierr;
567f746d493SDmitry Karpeev   PetscInt       i;
568b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
569b1a0a8a3SJed Brown 
570b1a0a8a3SJed Brown   PetscFunctionBegin;
571b1a0a8a3SJed Brown   /*
5726a4f0f73SDmitry Karpeev      Support for limiting the restriction or interpolation only to the inner
573b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
574b1a0a8a3SJed Brown   */
575f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
576b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
577f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
5786a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5792fa5cd67SKarl Rupp   } else {
5806a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
581b1a0a8a3SJed Brown   }
5826a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr);
5836a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
5846a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5852fa5cd67SKarl Rupp   } else {
5866a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5876a4f0f73SDmitry Karpeev   }
58886b96d36SDmitry Karpeev   /* do the subdomain solves */
58986b96d36SDmitry Karpeev   for (i=0; i<osm->n; ++i) {
590b1a0a8a3SJed Brown     ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
591b1a0a8a3SJed Brown   }
592*930d09c1SFande Kong   /*Do we need to zero y ?? */
5936a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(y);CHKERRQ(ierr);
5946a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
5956a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5966a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);  PetscFunctionReturn(0);
5972fa5cd67SKarl Rupp   } else {
5986a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5996a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);  PetscFunctionReturn(0);
6006a4f0f73SDmitry Karpeev   }
601b1a0a8a3SJed Brown }
602b1a0a8a3SJed Brown 
603b1a0a8a3SJed Brown #undef __FUNCT__
604f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM"
605f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y)
606b1a0a8a3SJed Brown {
607f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
608b1a0a8a3SJed Brown   PetscErrorCode ierr;
609f746d493SDmitry Karpeev   PetscInt       i;
610b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
611b1a0a8a3SJed Brown 
612b1a0a8a3SJed Brown   PetscFunctionBegin;
613b1a0a8a3SJed Brown   /*
614b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
615b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
616b1a0a8a3SJed Brown 
617f746d493SDmitry Karpeev      Note: these are reversed from the PCApply_GASM() because we are applying the
618b1a0a8a3SJed Brown      transpose of the three terms
619b1a0a8a3SJed Brown   */
620f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
621b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
622f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
6236a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
6242fa5cd67SKarl Rupp   } else {
6256a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
626b1a0a8a3SJed Brown   }
6276a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr);
6286a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
6296a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
6302fa5cd67SKarl Rupp   } else {
6316a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
6326a4f0f73SDmitry Karpeev   }
633b1a0a8a3SJed Brown   /* do the local solves */
634ab3e923fSDmitry Karpeev   for (i=0; i<osm->n; ++i) { /* Note that the solves are local, so we can go to osm->n, rather than osm->nmax. */
635b1a0a8a3SJed Brown     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
636b1a0a8a3SJed Brown   }
6376a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(y);CHKERRQ(ierr);
6386a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
6396a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6406a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6412fa5cd67SKarl Rupp   } else {
6426a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6436a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6446a4f0f73SDmitry Karpeev   }
645b1a0a8a3SJed Brown   PetscFunctionReturn(0);
646b1a0a8a3SJed Brown }
647b1a0a8a3SJed Brown 
648b1a0a8a3SJed Brown #undef __FUNCT__
649a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM"
650a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc)
651b1a0a8a3SJed Brown {
652f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
653b1a0a8a3SJed Brown   PetscErrorCode ierr;
654b1a0a8a3SJed Brown   PetscInt       i;
655b1a0a8a3SJed Brown 
656b1a0a8a3SJed Brown   PetscFunctionBegin;
657b1a0a8a3SJed Brown   if (osm->ksp) {
658f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
659a06653b4SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
660b1a0a8a3SJed Brown     }
661b1a0a8a3SJed Brown   }
662b1a0a8a3SJed Brown   if (osm->pmat) {
663f746d493SDmitry Karpeev     if (osm->n > 0) {
664ab3e923fSDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
665b1a0a8a3SJed Brown     }
666b1a0a8a3SJed Brown   }
667d34fcf5fSBarry Smith   if (osm->x) {
668f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
669fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
670fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
671b1a0a8a3SJed Brown     }
672d34fcf5fSBarry Smith   }
6736bf464f9SBarry Smith   ierr = VecDestroy(&osm->gx);CHKERRQ(ierr);
6746bf464f9SBarry Smith   ierr = VecDestroy(&osm->gy);CHKERRQ(ierr);
675ab3e923fSDmitry Karpeev 
6766a4f0f73SDmitry Karpeev   ierr     = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr);
6776a4f0f73SDmitry Karpeev   ierr     = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr);
6788f3b4b4dSDmitry Karpeev   if (!osm->user_subdomains) {
6792c112581SDmitry Karpeev     ierr      = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr);
6808f3b4b4dSDmitry Karpeev     osm->N    = PETSC_DETERMINE;
6818f3b4b4dSDmitry Karpeev     osm->nmax = PETSC_DETERMINE;
6828f3b4b4dSDmitry Karpeev   }
683a06653b4SBarry Smith   PetscFunctionReturn(0);
684a06653b4SBarry Smith }
685a06653b4SBarry Smith 
686a06653b4SBarry Smith #undef __FUNCT__
687a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM"
688a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc)
689a06653b4SBarry Smith {
690a06653b4SBarry Smith   PC_GASM        *osm = (PC_GASM*)pc->data;
691a06653b4SBarry Smith   PetscErrorCode ierr;
692a06653b4SBarry Smith   PetscInt       i;
693a06653b4SBarry Smith 
694a06653b4SBarry Smith   PetscFunctionBegin;
695a06653b4SBarry Smith   ierr = PCReset_GASM(pc);CHKERRQ(ierr);
696135757f6SDmitry Karpeev 
697135757f6SDmitry Karpeev   /* PCReset will not destroy subdomains, if user_subdomains is true. */
698135757f6SDmitry Karpeev   ierr = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr);
699135757f6SDmitry Karpeev 
700a06653b4SBarry Smith   if (osm->ksp) {
701a06653b4SBarry Smith     for (i=0; i<osm->n; i++) {
7026bf464f9SBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
703a06653b4SBarry Smith     }
704a06653b4SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
705a06653b4SBarry Smith   }
706a06653b4SBarry Smith   ierr = PetscFree(osm->x);CHKERRQ(ierr);
707a06653b4SBarry Smith   ierr = PetscFree(osm->y);CHKERRQ(ierr);
708c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
709b1a0a8a3SJed Brown   PetscFunctionReturn(0);
710b1a0a8a3SJed Brown }
711b1a0a8a3SJed Brown 
712b1a0a8a3SJed Brown #undef __FUNCT__
713f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM"
7148c34d3f5SBarry Smith static PetscErrorCode PCSetFromOptions_GASM(PetscOptions *PetscOptionsObject,PC pc)
715a6dfd86eSKarl Rupp {
716f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
717b1a0a8a3SJed Brown   PetscErrorCode ierr;
718b1a0a8a3SJed Brown   PetscInt       blocks,ovl;
719b1a0a8a3SJed Brown   PetscBool      symset,flg;
720f746d493SDmitry Karpeev   PCGASMType     gasmtype;
721b1a0a8a3SJed Brown 
722b1a0a8a3SJed Brown   PetscFunctionBegin;
723b1a0a8a3SJed Brown   /* set the type to symmetric if matrix is symmetric */
724b1a0a8a3SJed Brown   if (!osm->type_set && pc->pmat) {
725b1a0a8a3SJed Brown     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
7262fa5cd67SKarl Rupp     if (symset && flg) osm->type = PC_GASM_BASIC;
727b1a0a8a3SJed Brown   }
728e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Generalized additive Schwarz options");CHKERRQ(ierr);
7298f3b4b4dSDmitry Karpeev   ierr = PetscOptionsBool("-pc_gasm_use_dm_subdomains","If subdomains aren't set, use DMCreateDomainDecomposition() to define subdomains.","PCGASMSetUseDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr);
7308f3b4b4dSDmitry Karpeev   ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg);CHKERRQ(ierr);
73165db9045SDmitry Karpeev   if (flg) {
7328f3b4b4dSDmitry Karpeev     ierr = PCGASMSetTotalSubdomains(pc,blocks);CHKERRQ(ierr);
73365db9045SDmitry Karpeev   }
73406b43e7eSDmitry Karpeev   ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
73565db9045SDmitry Karpeev   if (flg) {
73665db9045SDmitry Karpeev     ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr);
737d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
73865db9045SDmitry Karpeev   }
739b1a0a8a3SJed Brown   flg  = PETSC_FALSE;
740f746d493SDmitry Karpeev   ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr);
741f746d493SDmitry Karpeev   if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr);}
742b1a0a8a3SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
743b1a0a8a3SJed Brown   PetscFunctionReturn(0);
744b1a0a8a3SJed Brown }
745b1a0a8a3SJed Brown 
746b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/
747b1a0a8a3SJed Brown 
748b1a0a8a3SJed Brown #undef __FUNCT__
7498f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains"
7508f3b4b4dSDmitry Karpeev /*@
7518f3b4b4dSDmitry Karpeev     PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the
7528f3b4b4dSDmitry Karpeev                                communicator.
7538f3b4b4dSDmitry Karpeev     Logically collective on pc
7548f3b4b4dSDmitry Karpeev 
7558f3b4b4dSDmitry Karpeev     Input Parameters:
7568f3b4b4dSDmitry Karpeev +   pc  - the preconditioner
7578f3b4b4dSDmitry Karpeev -   N   - total number of subdomains
7588f3b4b4dSDmitry Karpeev 
7598f3b4b4dSDmitry Karpeev 
7608f3b4b4dSDmitry Karpeev     Level: beginner
7618f3b4b4dSDmitry Karpeev 
7628f3b4b4dSDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
7638f3b4b4dSDmitry Karpeev 
7648f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap()
7658f3b4b4dSDmitry Karpeev           PCGASMCreateSubdomains2D()
7668f3b4b4dSDmitry Karpeev @*/
7678f3b4b4dSDmitry Karpeev PetscErrorCode  PCGASMSetTotalSubdomains(PC pc,PetscInt N)
7688f3b4b4dSDmitry Karpeev {
7698f3b4b4dSDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
7708f3b4b4dSDmitry Karpeev   PetscMPIInt    size,rank;
7718f3b4b4dSDmitry Karpeev   PetscErrorCode ierr;
7728f3b4b4dSDmitry Karpeev 
7738f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
7748f3b4b4dSDmitry Karpeev   if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be 1 or more, got N = %D",N);
7758f3b4b4dSDmitry Karpeev   if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before calling PCSetUp().");
7768f3b4b4dSDmitry Karpeev 
7772c112581SDmitry Karpeev   ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr);
7788f3b4b4dSDmitry Karpeev   osm->ois = osm->iis = NULL;
7798f3b4b4dSDmitry Karpeev 
7808f3b4b4dSDmitry Karpeev   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
7818f3b4b4dSDmitry Karpeev   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
7828f3b4b4dSDmitry Karpeev   osm->N    = N;
7838f3b4b4dSDmitry Karpeev   osm->n    = PETSC_DETERMINE;
7848f3b4b4dSDmitry Karpeev   osm->nmax = PETSC_DETERMINE;
7858f3b4b4dSDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
7868f3b4b4dSDmitry Karpeev   PetscFunctionReturn(0);
7878f3b4b4dSDmitry Karpeev }
7888f3b4b4dSDmitry Karpeev 
7898f3b4b4dSDmitry Karpeev #undef __FUNCT__
79006b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM"
791b541e6a4SDmitry Karpeev static PetscErrorCode  PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[])
792b1a0a8a3SJed Brown {
793f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
794b1a0a8a3SJed Brown   PetscErrorCode ierr;
795b1a0a8a3SJed Brown   PetscInt       i;
796b1a0a8a3SJed Brown 
797b1a0a8a3SJed Brown   PetscFunctionBegin;
7988f3b4b4dSDmitry Karpeev   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, got n = %D",n);
7998f3b4b4dSDmitry Karpeev   if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp().");
800b1a0a8a3SJed Brown 
8012c112581SDmitry Karpeev   ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr);
8028f3b4b4dSDmitry Karpeev   osm->iis  = osm->ois = NULL;
8038f3b4b4dSDmitry Karpeev   osm->n    = n;
8048f3b4b4dSDmitry Karpeev   osm->N    = PETSC_DETERMINE;
8058f3b4b4dSDmitry Karpeev   osm->nmax = PETSC_DETERMINE;
806a35b7b57SDmitry Karpeev   if (ois) {
807785e854fSJed Brown     ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr);
8088f3b4b4dSDmitry Karpeev     for (i=0; i<n; i++) {
8098f3b4b4dSDmitry Karpeev       ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);
8108f3b4b4dSDmitry Karpeev       osm->ois[i] = ois[i];
8118f3b4b4dSDmitry Karpeev     }
8128f3b4b4dSDmitry Karpeev     /*
8138f3b4b4dSDmitry Karpeev        Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(),
8148f3b4b4dSDmitry Karpeev        it will be ignored.  To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1.
8158f3b4b4dSDmitry Karpeev     */
816b1a0a8a3SJed Brown     osm->overlap = -1;
817a35b7b57SDmitry Karpeev     if (!iis) {
818785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr);
819a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {
8208f3b4b4dSDmitry Karpeev 	for (i=0; i<n; i++) {
8210e8f3646SDmitry Karpeev 	  ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);
822a35b7b57SDmitry Karpeev 	  osm->iis[i] = ois[i];
823a35b7b57SDmitry Karpeev 	}
824a35b7b57SDmitry Karpeev       }
825a35b7b57SDmitry Karpeev     }
8268f3b4b4dSDmitry Karpeev   }
827a35b7b57SDmitry Karpeev   if (iis) {
828785e854fSJed Brown     ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr);
8298f3b4b4dSDmitry Karpeev     for (i=0; i<n; i++) {
8308f3b4b4dSDmitry Karpeev       ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);
8318f3b4b4dSDmitry Karpeev       osm->iis[i] = iis[i];
8328f3b4b4dSDmitry Karpeev     }
833a35b7b57SDmitry Karpeev     if (!ois) {
834785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr);
835a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {
8368f3b4b4dSDmitry Karpeev 	for (i=0; i<n; i++) {
837a35b7b57SDmitry Karpeev 	  ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);
838a35b7b57SDmitry Karpeev 	  osm->ois[i] = iis[i];
839a35b7b57SDmitry Karpeev 	}
8408f3b4b4dSDmitry Karpeev       }
8418f3b4b4dSDmitry Karpeev       /* If necessary, osm->ois[i] will be expanded using the requested overlap size in PCSetUp_GASM(). */
842a35b7b57SDmitry Karpeev     }
843a35b7b57SDmitry Karpeev   }
8448f3b4b4dSDmitry Karpeev   if (iis)  osm->user_subdomains = PETSC_TRUE;
845b1a0a8a3SJed Brown   PetscFunctionReturn(0);
846b1a0a8a3SJed Brown }
847b1a0a8a3SJed Brown 
848b1a0a8a3SJed Brown 
849b1a0a8a3SJed Brown #undef __FUNCT__
850f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM"
851f7a08781SBarry Smith static PetscErrorCode  PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
852b1a0a8a3SJed Brown {
853f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
854b1a0a8a3SJed Brown 
855b1a0a8a3SJed Brown   PetscFunctionBegin;
856ce94432eSBarry Smith   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
857ce94432eSBarry Smith   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
8582fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
859b1a0a8a3SJed Brown   PetscFunctionReturn(0);
860b1a0a8a3SJed Brown }
861b1a0a8a3SJed Brown 
862b1a0a8a3SJed Brown #undef __FUNCT__
863f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM"
864f7a08781SBarry Smith static PetscErrorCode  PCGASMSetType_GASM(PC pc,PCGASMType type)
865b1a0a8a3SJed Brown {
866f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
867b1a0a8a3SJed Brown 
868b1a0a8a3SJed Brown   PetscFunctionBegin;
869b1a0a8a3SJed Brown   osm->type     = type;
870b1a0a8a3SJed Brown   osm->type_set = PETSC_TRUE;
871b1a0a8a3SJed Brown   PetscFunctionReturn(0);
872b1a0a8a3SJed Brown }
873b1a0a8a3SJed Brown 
874b1a0a8a3SJed Brown #undef __FUNCT__
875f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM"
876f7a08781SBarry Smith static PetscErrorCode  PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort)
877b1a0a8a3SJed Brown {
878f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
879b1a0a8a3SJed Brown 
880b1a0a8a3SJed Brown   PetscFunctionBegin;
881b1a0a8a3SJed Brown   osm->sort_indices = doSort;
882b1a0a8a3SJed Brown   PetscFunctionReturn(0);
883b1a0a8a3SJed Brown }
884b1a0a8a3SJed Brown 
885b1a0a8a3SJed Brown #undef __FUNCT__
886f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM"
88744fe18e5SDmitry Karpeev /*
8888f3b4b4dSDmitry Karpeev    FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed.
88944fe18e5SDmitry Karpeev         In particular, it would upset the global subdomain number calculation.
89044fe18e5SDmitry Karpeev */
891f7a08781SBarry Smith static PetscErrorCode  PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
892b1a0a8a3SJed Brown {
893f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
894b1a0a8a3SJed Brown   PetscErrorCode ierr;
895b1a0a8a3SJed Brown 
896b1a0a8a3SJed Brown   PetscFunctionBegin;
897ce94432eSBarry Smith   if (osm->n < 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here");
898b1a0a8a3SJed Brown 
8992fa5cd67SKarl Rupp   if (n) *n = osm->n;
900ab3e923fSDmitry Karpeev   if (first) {
901ce94432eSBarry Smith     ierr    = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
902ab3e923fSDmitry Karpeev     *first -= osm->n;
903b1a0a8a3SJed Brown   }
904b1a0a8a3SJed Brown   if (ksp) {
905b1a0a8a3SJed Brown     /* Assume that local solves are now different; not necessarily
90606b43e7eSDmitry Karpeev        true, though!  This flag is used only for PCView_GASM() */
907b1a0a8a3SJed Brown     *ksp                        = osm->ksp;
9086a4f0f73SDmitry Karpeev     osm->same_subdomain_solvers = PETSC_FALSE;
909b1a0a8a3SJed Brown   }
910b1a0a8a3SJed Brown   PetscFunctionReturn(0);
911ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */
912b1a0a8a3SJed Brown 
913b1a0a8a3SJed Brown #undef __FUNCT__
91406b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains"
915b1a0a8a3SJed Brown /*@C
91606b43e7eSDmitry Karpeev     PCGASMSetSubdomains - Sets the subdomains for this processor
91706b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
918b1a0a8a3SJed Brown 
9198f3b4b4dSDmitry Karpeev     Collective on pc
920b1a0a8a3SJed Brown 
921b1a0a8a3SJed Brown     Input Parameters:
9228f3b4b4dSDmitry Karpeev +   pc  - the preconditioner object
92306b43e7eSDmitry Karpeev .   n   - the number of subdomains for this processor
9248f3b4b4dSDmitry Karpeev .   iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains)
9258f3b4b4dSDmitry Karpeev -   ois - the index sets that define the outer subdomains (or NULL to use the same as iis, or to construct by expanding iis by the requested overlap)
926b1a0a8a3SJed Brown 
927b1a0a8a3SJed Brown     Notes:
92806b43e7eSDmitry Karpeev     The IS indices use the parallel, global numbering of the vector entries.
9296a4f0f73SDmitry Karpeev     Inner subdomains are those where the correction is applied.
9306a4f0f73SDmitry Karpeev     Outer subdomains are those where the residual necessary to obtain the
9316a4f0f73SDmitry Karpeev     corrections is obtained (see PCGASMType for the use of inner/outer subdomains).
9326a4f0f73SDmitry Karpeev     Both inner and outer subdomains can extend over several processors.
9336a4f0f73SDmitry Karpeev     This processor's portion of a subdomain is known as a local subdomain.
9346a4f0f73SDmitry Karpeev 
9356a4f0f73SDmitry Karpeev     By default the GASM preconditioner uses 1 (local) subdomain per processor.
9366a4f0f73SDmitry Karpeev 
937b1a0a8a3SJed Brown 
938b1a0a8a3SJed Brown     Level: advanced
939b1a0a8a3SJed Brown 
94006b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz
941b1a0a8a3SJed Brown 
9428f3b4b4dSDmitry Karpeev .seealso: PCGASMSetNumSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
94306b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
944b1a0a8a3SJed Brown @*/
9456a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[])
946b1a0a8a3SJed Brown {
9478f3b4b4dSDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
948b1a0a8a3SJed Brown   PetscErrorCode ierr;
949b1a0a8a3SJed Brown 
950b1a0a8a3SJed Brown   PetscFunctionBegin;
951b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9526a4f0f73SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr);
9538f3b4b4dSDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
954b1a0a8a3SJed Brown   PetscFunctionReturn(0);
955b1a0a8a3SJed Brown }
956b1a0a8a3SJed Brown 
957b1a0a8a3SJed Brown 
958b1a0a8a3SJed Brown #undef __FUNCT__
959f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap"
960b1a0a8a3SJed Brown /*@
961f746d493SDmitry Karpeev     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
962b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
9638f3b4b4dSDmitry Karpeev     pc communicator must call this routine.
964b1a0a8a3SJed Brown 
9658f3b4b4dSDmitry Karpeev     Logically Collective on pc
966b1a0a8a3SJed Brown 
967b1a0a8a3SJed Brown     Input Parameters:
968b1a0a8a3SJed Brown +   pc  - the preconditioner context
9698f3b4b4dSDmitry Karpeev -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0)
970b1a0a8a3SJed Brown 
971b1a0a8a3SJed Brown     Options Database Key:
97206b43e7eSDmitry Karpeev .   -pc_gasm_overlap <overlap> - Sets overlap
973b1a0a8a3SJed Brown 
974b1a0a8a3SJed Brown     Notes:
97506b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.  To use
9768f3b4b4dSDmitry Karpeev     multiple subdomain per perocessor or "straddling" subdomains that intersect
9778f3b4b4dSDmitry Karpeev     multiple processors use PCGASMSetSubdomains() (or option -pc_gasm_total_subdomains <n>).
978b1a0a8a3SJed Brown 
9798f3b4b4dSDmitry Karpeev     The overlap defaults to 0, so if one desires that no additional
980b1a0a8a3SJed Brown     overlap be computed beyond what may have been set with a call to
9818f3b4b4dSDmitry Karpeev     PCGASMSetSubdomains(), then ovl must be set to be 0.  In particular, if one does
9828f3b4b4dSDmitry Karpeev     not explicitly set the subdomains in application code, then all overlap would be computed
983f746d493SDmitry Karpeev     internally by PETSc, and using an overlap of 0 would result in an GASM
984b1a0a8a3SJed Brown     variant that is equivalent to the block Jacobi preconditioner.
985b1a0a8a3SJed Brown 
986b1a0a8a3SJed Brown     Note that one can define initial index sets with any overlap via
98706b43e7eSDmitry Karpeev     PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows
98806b43e7eSDmitry Karpeev     PETSc to extend that overlap further, if desired.
989b1a0a8a3SJed Brown 
990b1a0a8a3SJed Brown     Level: intermediate
991b1a0a8a3SJed Brown 
992f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap
993b1a0a8a3SJed Brown 
9948f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
99506b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
996b1a0a8a3SJed Brown @*/
9977087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap(PC pc,PetscInt ovl)
998b1a0a8a3SJed Brown {
999b1a0a8a3SJed Brown   PetscErrorCode ierr;
10008f3b4b4dSDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
1001b1a0a8a3SJed Brown 
1002b1a0a8a3SJed Brown   PetscFunctionBegin;
1003b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1004b1a0a8a3SJed Brown   PetscValidLogicalCollectiveInt(pc,ovl,2);
1005f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
10068f3b4b4dSDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
1007b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1008b1a0a8a3SJed Brown }
1009b1a0a8a3SJed Brown 
1010b1a0a8a3SJed Brown #undef __FUNCT__
1011f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType"
1012b1a0a8a3SJed Brown /*@
1013f746d493SDmitry Karpeev     PCGASMSetType - Sets the type of restriction and interpolation used
1014b1a0a8a3SJed Brown     for local problems in the additive Schwarz method.
1015b1a0a8a3SJed Brown 
1016b1a0a8a3SJed Brown     Logically Collective on PC
1017b1a0a8a3SJed Brown 
1018b1a0a8a3SJed Brown     Input Parameters:
1019b1a0a8a3SJed Brown +   pc  - the preconditioner context
1020f746d493SDmitry Karpeev -   type - variant of GASM, one of
1021b1a0a8a3SJed Brown .vb
1022f746d493SDmitry Karpeev       PC_GASM_BASIC       - full interpolation and restriction
1023f746d493SDmitry Karpeev       PC_GASM_RESTRICT    - full restriction, local processor interpolation
1024f746d493SDmitry Karpeev       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
1025f746d493SDmitry Karpeev       PC_GASM_NONE        - local processor restriction and interpolation
1026b1a0a8a3SJed Brown .ve
1027b1a0a8a3SJed Brown 
1028b1a0a8a3SJed Brown     Options Database Key:
1029f746d493SDmitry Karpeev .   -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1030b1a0a8a3SJed Brown 
1031b1a0a8a3SJed Brown     Level: intermediate
1032b1a0a8a3SJed Brown 
1033f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1034b1a0a8a3SJed Brown 
10358f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1036f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1037b1a0a8a3SJed Brown @*/
10387087cfbeSBarry Smith PetscErrorCode  PCGASMSetType(PC pc,PCGASMType type)
1039b1a0a8a3SJed Brown {
1040b1a0a8a3SJed Brown   PetscErrorCode ierr;
1041b1a0a8a3SJed Brown 
1042b1a0a8a3SJed Brown   PetscFunctionBegin;
1043b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1044b1a0a8a3SJed Brown   PetscValidLogicalCollectiveEnum(pc,type,2);
1045f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr);
1046b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1047b1a0a8a3SJed Brown }
1048b1a0a8a3SJed Brown 
1049b1a0a8a3SJed Brown #undef __FUNCT__
1050f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices"
1051b1a0a8a3SJed Brown /*@
1052f746d493SDmitry Karpeev     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
1053b1a0a8a3SJed Brown 
1054b1a0a8a3SJed Brown     Logically Collective on PC
1055b1a0a8a3SJed Brown 
1056b1a0a8a3SJed Brown     Input Parameters:
1057b1a0a8a3SJed Brown +   pc  - the preconditioner context
1058b1a0a8a3SJed Brown -   doSort - sort the subdomain indices
1059b1a0a8a3SJed Brown 
1060b1a0a8a3SJed Brown     Level: intermediate
1061b1a0a8a3SJed Brown 
1062f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1063b1a0a8a3SJed Brown 
10648f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1065f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1066b1a0a8a3SJed Brown @*/
10677087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices(PC pc,PetscBool doSort)
1068b1a0a8a3SJed Brown {
1069b1a0a8a3SJed Brown   PetscErrorCode ierr;
1070b1a0a8a3SJed Brown 
1071b1a0a8a3SJed Brown   PetscFunctionBegin;
1072b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1073b1a0a8a3SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
1074f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
1075b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1076b1a0a8a3SJed Brown }
1077b1a0a8a3SJed Brown 
1078b1a0a8a3SJed Brown #undef __FUNCT__
1079f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP"
1080b1a0a8a3SJed Brown /*@C
1081f746d493SDmitry Karpeev    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
1082b1a0a8a3SJed Brown    this processor.
1083b1a0a8a3SJed Brown 
1084b1a0a8a3SJed Brown    Collective on PC iff first_local is requested
1085b1a0a8a3SJed Brown 
1086b1a0a8a3SJed Brown    Input Parameter:
1087b1a0a8a3SJed Brown .  pc - the preconditioner context
1088b1a0a8a3SJed Brown 
1089b1a0a8a3SJed Brown    Output Parameters:
10900298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
10910298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
10920298fd71SBarry Smith                  all processors must request or all must pass NULL
1093b1a0a8a3SJed Brown -  ksp - the array of KSP contexts
1094b1a0a8a3SJed Brown 
1095b1a0a8a3SJed Brown    Note:
1096f746d493SDmitry Karpeev    After PCGASMGetSubKSP() the array of KSPes is not to be freed
1097b1a0a8a3SJed Brown 
1098b1a0a8a3SJed Brown    Currently for some matrix implementations only 1 block per processor
1099b1a0a8a3SJed Brown    is supported.
1100b1a0a8a3SJed Brown 
1101f746d493SDmitry Karpeev    You must call KSPSetUp() before calling PCGASMGetSubKSP().
1102b1a0a8a3SJed Brown 
1103b1a0a8a3SJed Brown    Level: advanced
1104b1a0a8a3SJed Brown 
1105f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context
1106b1a0a8a3SJed Brown 
11078f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(),
1108f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(),
1109b1a0a8a3SJed Brown @*/
11107087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1111b1a0a8a3SJed Brown {
1112b1a0a8a3SJed Brown   PetscErrorCode ierr;
1113b1a0a8a3SJed Brown 
1114b1a0a8a3SJed Brown   PetscFunctionBegin;
1115b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1116f746d493SDmitry Karpeev   ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
1117b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1118b1a0a8a3SJed Brown }
1119b1a0a8a3SJed Brown 
1120b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/
1121b1a0a8a3SJed Brown /*MC
1122f746d493SDmitry Karpeev    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1123b1a0a8a3SJed Brown            its own KSP object.
1124b1a0a8a3SJed Brown 
1125b1a0a8a3SJed Brown    Options Database Keys:
11268f3b4b4dSDmitry Karpeev +  -pc_gasm_total_subdomains <n>  - Sets total number of local subdomains to be distributed among processors
112706b43e7eSDmitry Karpeev .  -pc_gasm_view_subdomains       - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view
112806b43e7eSDmitry Karpeev .  -pc_gasm_print_subdomains      - activates the printing of subdomain indices in PCSetUp()
112906b43e7eSDmitry Karpeev .  -pc_gasm_overlap <ovl>         - Sets overlap by which to (automatically) extend local subdomains
1130f746d493SDmitry Karpeev -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1131b1a0a8a3SJed Brown 
1132b1a0a8a3SJed Brown      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
1133f746d493SDmitry Karpeev       will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
1134f746d493SDmitry Karpeev       -pc_gasm_type basic to use the standard GASM.
1135b1a0a8a3SJed Brown 
1136b1a0a8a3SJed Brown    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
1137b1a0a8a3SJed Brown      than one processor. Defaults to one block per processor.
1138b1a0a8a3SJed Brown 
1139b1a0a8a3SJed Brown      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1140b1a0a8a3SJed Brown         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
1141b1a0a8a3SJed Brown 
1142f746d493SDmitry Karpeev      To set the options on the solvers separate for each block call PCGASMGetSubKSP()
1143b1a0a8a3SJed Brown          and set the options directly on the resulting KSP object (you can access its PC
1144b1a0a8a3SJed Brown          with KSPGetPC())
1145b1a0a8a3SJed Brown 
1146b1a0a8a3SJed Brown 
1147b1a0a8a3SJed Brown    Level: beginner
1148b1a0a8a3SJed Brown 
1149b1a0a8a3SJed Brown    Concepts: additive Schwarz method
1150b1a0a8a3SJed Brown 
1151b1a0a8a3SJed Brown     References:
1152b1a0a8a3SJed Brown     An additive variant of the Schwarz alternating method for the case of many subregions
1153b1a0a8a3SJed Brown     M Dryja, OB Widlund - Courant Institute, New York University Technical report
1154b1a0a8a3SJed Brown 
1155b1a0a8a3SJed Brown     Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1156b1a0a8a3SJed Brown     Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X.
1157b1a0a8a3SJed Brown 
1158b1a0a8a3SJed Brown .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
115949517cdeSBarry Smith            PCBJACOBI,  PCGASMGetSubKSP(), PCGASMSetSubdomains(),
11608f3b4b4dSDmitry Karpeev            PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType()
1161b1a0a8a3SJed Brown 
1162b1a0a8a3SJed Brown M*/
1163b1a0a8a3SJed Brown 
1164b1a0a8a3SJed Brown #undef __FUNCT__
1165f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM"
11668cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc)
1167b1a0a8a3SJed Brown {
1168b1a0a8a3SJed Brown   PetscErrorCode ierr;
1169f746d493SDmitry Karpeev   PC_GASM        *osm;
1170b1a0a8a3SJed Brown 
1171b1a0a8a3SJed Brown   PetscFunctionBegin;
1172b00a9115SJed Brown   ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr);
11732fa5cd67SKarl Rupp 
11748f3b4b4dSDmitry Karpeev   osm->N                      = PETSC_DETERMINE;
117506b43e7eSDmitry Karpeev   osm->n                      = PETSC_DECIDE;
11768f3b4b4dSDmitry Karpeev   osm->nmax                   = PETSC_DETERMINE;
11778f3b4b4dSDmitry Karpeev   osm->overlap                = 0;
1178b1a0a8a3SJed Brown   osm->ksp                    = 0;
11796a4f0f73SDmitry Karpeev   osm->gorestriction          = 0;
11806a4f0f73SDmitry Karpeev   osm->girestriction          = 0;
1181ab3e923fSDmitry Karpeev   osm->gx                     = 0;
1182ab3e923fSDmitry Karpeev   osm->gy                     = 0;
1183b1a0a8a3SJed Brown   osm->x                      = 0;
1184b1a0a8a3SJed Brown   osm->y                      = 0;
11856a4f0f73SDmitry Karpeev   osm->ois                    = 0;
11866a4f0f73SDmitry Karpeev   osm->iis                    = 0;
1187b1a0a8a3SJed Brown   osm->pmat                   = 0;
1188f746d493SDmitry Karpeev   osm->type                   = PC_GASM_RESTRICT;
11896a4f0f73SDmitry Karpeev   osm->same_subdomain_solvers = PETSC_TRUE;
1190b1a0a8a3SJed Brown   osm->sort_indices           = PETSC_TRUE;
1191d709ea83SDmitry Karpeev   osm->dm_subdomains          = PETSC_FALSE;
1192b1a0a8a3SJed Brown 
1193b1a0a8a3SJed Brown   pc->data                 = (void*)osm;
1194f746d493SDmitry Karpeev   pc->ops->apply           = PCApply_GASM;
1195f746d493SDmitry Karpeev   pc->ops->applytranspose  = PCApplyTranspose_GASM;
1196f746d493SDmitry Karpeev   pc->ops->setup           = PCSetUp_GASM;
1197a06653b4SBarry Smith   pc->ops->reset           = PCReset_GASM;
1198f746d493SDmitry Karpeev   pc->ops->destroy         = PCDestroy_GASM;
1199f746d493SDmitry Karpeev   pc->ops->setfromoptions  = PCSetFromOptions_GASM;
1200f746d493SDmitry Karpeev   pc->ops->setuponblocks   = PCSetUpOnBlocks_GASM;
1201f746d493SDmitry Karpeev   pc->ops->view            = PCView_GASM;
1202b1a0a8a3SJed Brown   pc->ops->applyrichardson = 0;
1203b1a0a8a3SJed Brown 
1204bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM);CHKERRQ(ierr);
1205bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM);CHKERRQ(ierr);
1206bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM);CHKERRQ(ierr);
1207bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM);CHKERRQ(ierr);
1208bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM);CHKERRQ(ierr);
1209b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1210b1a0a8a3SJed Brown }
1211b1a0a8a3SJed Brown 
1212b1a0a8a3SJed Brown 
1213b1a0a8a3SJed Brown #undef __FUNCT__
121406b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains"
12158f3b4b4dSDmitry Karpeev PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[])
1216b1a0a8a3SJed Brown {
1217b1a0a8a3SJed Brown   MatPartitioning mpart;
1218b1a0a8a3SJed Brown   const char      *prefix;
121911bd1e4dSLisandro Dalcin   PetscErrorCode  (*f)(Mat,MatReuse,Mat*);
1220b1a0a8a3SJed Brown   PetscMPIInt     size;
1221b1a0a8a3SJed Brown   PetscInt        i,j,rstart,rend,bs;
122211bd1e4dSLisandro Dalcin   PetscBool       isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
12230298fd71SBarry Smith   Mat             Ad     = NULL, adj;
1224b1a0a8a3SJed Brown   IS              ispart,isnumb,*is;
1225b1a0a8a3SJed Brown   PetscErrorCode  ierr;
1226b1a0a8a3SJed Brown 
1227b1a0a8a3SJed Brown   PetscFunctionBegin;
12288f3b4b4dSDmitry Karpeev   if (nloc < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local subdomains must > 0, got nloc = %D",nloc);
1229b1a0a8a3SJed Brown 
1230b1a0a8a3SJed Brown   /* Get prefix, row distribution, and block size */
1231b1a0a8a3SJed Brown   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1232b1a0a8a3SJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1233b1a0a8a3SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1234b1a0a8a3SJed Brown   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);
1235b1a0a8a3SJed Brown 
1236b1a0a8a3SJed Brown   /* Get diagonal block from matrix if possible */
1237ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
12380005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",&f);CHKERRQ(ierr);
1239b1a0a8a3SJed Brown   if (f) {
124011bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
1241b1a0a8a3SJed Brown   } else if (size == 1) {
124211bd1e4dSLisandro Dalcin     Ad = A;
1243b1a0a8a3SJed Brown   }
1244b1a0a8a3SJed Brown   if (Ad) {
1245251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1246251f4c67SDmitry Karpeev     if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
1247b1a0a8a3SJed Brown   }
12488f3b4b4dSDmitry Karpeev   if (Ad && nloc > 1) {
1249b1a0a8a3SJed Brown     PetscBool  match,done;
1250b1a0a8a3SJed Brown     /* Try to setup a good matrix partitioning if available */
1251b1a0a8a3SJed Brown     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
1252b1a0a8a3SJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
1253b1a0a8a3SJed Brown     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1254251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
1255b1a0a8a3SJed Brown     if (!match) {
1256251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
1257b1a0a8a3SJed Brown     }
1258b1a0a8a3SJed Brown     if (!match) { /* assume a "good" partitioner is available */
12591a83f524SJed Brown       PetscInt       na;
12601a83f524SJed Brown       const PetscInt *ia,*ja;
1261b1a0a8a3SJed Brown       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1262b1a0a8a3SJed Brown       if (done) {
1263b1a0a8a3SJed Brown         /* Build adjacency matrix by hand. Unfortunately a call to
1264b1a0a8a3SJed Brown            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1265b1a0a8a3SJed Brown            remove the block-aij structure and we cannot expect
1266b1a0a8a3SJed Brown            MatPartitioning to split vertices as we need */
12671a83f524SJed Brown         PetscInt       i,j,len,nnz,cnt,*iia=0,*jja=0;
12681a83f524SJed Brown         const PetscInt *row;
1269b1a0a8a3SJed Brown         nnz = 0;
1270b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* count number of nonzeros */
1271b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1272b1a0a8a3SJed Brown           row = ja + ia[i];
1273b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1274b1a0a8a3SJed Brown             if (row[j] == i) { /* don't count diagonal */
1275b1a0a8a3SJed Brown               len--; break;
1276b1a0a8a3SJed Brown             }
1277b1a0a8a3SJed Brown           }
1278b1a0a8a3SJed Brown           nnz += len;
1279b1a0a8a3SJed Brown         }
1280854ce69bSBarry Smith         ierr   = PetscMalloc1(na+1,&iia);CHKERRQ(ierr);
1281854ce69bSBarry Smith         ierr   = PetscMalloc1(nnz,&jja);CHKERRQ(ierr);
1282b1a0a8a3SJed Brown         nnz    = 0;
1283b1a0a8a3SJed Brown         iia[0] = 0;
1284b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* fill adjacency */
1285b1a0a8a3SJed Brown           cnt = 0;
1286b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1287b1a0a8a3SJed Brown           row = ja + ia[i];
1288b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
12892fa5cd67SKarl Rupp             if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */
1290b1a0a8a3SJed Brown           }
1291b1a0a8a3SJed Brown           nnz += cnt;
1292b1a0a8a3SJed Brown           iia[i+1] = nnz;
1293b1a0a8a3SJed Brown         }
1294b1a0a8a3SJed Brown         /* Partitioning of the adjacency matrix */
12950298fd71SBarry Smith         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr);
1296b1a0a8a3SJed Brown         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
12978f3b4b4dSDmitry Karpeev         ierr      = MatPartitioningSetNParts(mpart,nloc);CHKERRQ(ierr);
1298b1a0a8a3SJed Brown         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
1299b1a0a8a3SJed Brown         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
13006bf464f9SBarry Smith         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
1301b1a0a8a3SJed Brown         foundpart = PETSC_TRUE;
1302b1a0a8a3SJed Brown       }
1303b1a0a8a3SJed Brown       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1304b1a0a8a3SJed Brown     }
13056bf464f9SBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1306b1a0a8a3SJed Brown   }
13078f3b4b4dSDmitry Karpeev   ierr = PetscMalloc1(nloc,&is);CHKERRQ(ierr);
1308b1a0a8a3SJed Brown   if (!foundpart) {
1309b1a0a8a3SJed Brown 
1310b1a0a8a3SJed Brown     /* Partitioning by contiguous chunks of rows */
1311b1a0a8a3SJed Brown 
1312b1a0a8a3SJed Brown     PetscInt mbs   = (rend-rstart)/bs;
1313b1a0a8a3SJed Brown     PetscInt start = rstart;
13148f3b4b4dSDmitry Karpeev     for (i=0; i<nloc; i++) {
13158f3b4b4dSDmitry Karpeev       PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs;
1316b1a0a8a3SJed Brown       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
1317b1a0a8a3SJed Brown       start += count;
1318b1a0a8a3SJed Brown     }
1319b1a0a8a3SJed Brown 
1320b1a0a8a3SJed Brown   } else {
1321b1a0a8a3SJed Brown 
1322b1a0a8a3SJed Brown     /* Partitioning by adjacency of diagonal block  */
1323b1a0a8a3SJed Brown 
1324b1a0a8a3SJed Brown     const PetscInt *numbering;
1325b1a0a8a3SJed Brown     PetscInt       *count,nidx,*indices,*newidx,start=0;
1326b1a0a8a3SJed Brown     /* Get node count in each partition */
13278f3b4b4dSDmitry Karpeev     ierr = PetscMalloc1(nloc,&count);CHKERRQ(ierr);
13288f3b4b4dSDmitry Karpeev     ierr = ISPartitioningCount(ispart,nloc,count);CHKERRQ(ierr);
1329b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
13308f3b4b4dSDmitry Karpeev       for (i=0; i<nloc; i++) count[i] *= bs;
1331b1a0a8a3SJed Brown     }
1332b1a0a8a3SJed Brown     /* Build indices from node numbering */
1333b1a0a8a3SJed Brown     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1334785e854fSJed Brown     ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr);
1335b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1336b1a0a8a3SJed Brown     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
1337b1a0a8a3SJed Brown     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
1338b1a0a8a3SJed Brown     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
1339b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1340785e854fSJed Brown       ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr);
13412fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
13422fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
13432fa5cd67SKarl Rupp       }
1344b1a0a8a3SJed Brown       ierr    = PetscFree(indices);CHKERRQ(ierr);
1345b1a0a8a3SJed Brown       nidx   *= bs;
1346b1a0a8a3SJed Brown       indices = newidx;
1347b1a0a8a3SJed Brown     }
1348b1a0a8a3SJed Brown     /* Shift to get global indices */
1349b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] += rstart;
1350b1a0a8a3SJed Brown 
1351b1a0a8a3SJed Brown     /* Build the index sets for each block */
13528f3b4b4dSDmitry Karpeev     for (i=0; i<nloc; i++) {
1353b1a0a8a3SJed Brown       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
1354b1a0a8a3SJed Brown       ierr   = ISSort(is[i]);CHKERRQ(ierr);
1355b1a0a8a3SJed Brown       start += count[i];
1356b1a0a8a3SJed Brown     }
1357b1a0a8a3SJed Brown 
1358302440fdSBarry Smith     ierr = PetscFree(count);CHKERRQ(ierr);
1359302440fdSBarry Smith     ierr = PetscFree(indices);CHKERRQ(ierr);
1360fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1361fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
1362b1a0a8a3SJed Brown   }
13636a4f0f73SDmitry Karpeev   *iis = is;
13648f3b4b4dSDmitry Karpeev   PetscFunctionReturn(0);
13658f3b4b4dSDmitry Karpeev }
13668f3b4b4dSDmitry Karpeev 
13678f3b4b4dSDmitry Karpeev #undef __FUNCT__
13688f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMCreateStraddlingSubdomains"
1369b541e6a4SDmitry Karpeev PETSC_INTERN PetscErrorCode  PCGASMCreateStraddlingSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[])
13708f3b4b4dSDmitry Karpeev {
13718f3b4b4dSDmitry Karpeev   PetscErrorCode  ierr;
13728f3b4b4dSDmitry Karpeev 
13738f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
13748f3b4b4dSDmitry Karpeev   ierr = MatSubdomainsCreateCoalesce(A,N,n,iis);CHKERRQ(ierr);
13758f3b4b4dSDmitry Karpeev   PetscFunctionReturn(0);
13768f3b4b4dSDmitry Karpeev }
13778f3b4b4dSDmitry Karpeev 
13788f3b4b4dSDmitry Karpeev 
13798f3b4b4dSDmitry Karpeev 
13808f3b4b4dSDmitry Karpeev #undef __FUNCT__
13818f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains"
13828f3b4b4dSDmitry Karpeev /*@C
13838f3b4b4dSDmitry Karpeev    PCGASMCreateSubdomains - Creates n index sets defining n nonoverlapping subdomains for the additive
13848f3b4b4dSDmitry Karpeev    Schwarz preconditioner for a any problem based on its matrix.
13858f3b4b4dSDmitry Karpeev 
13868f3b4b4dSDmitry Karpeev    Collective
13878f3b4b4dSDmitry Karpeev 
13888f3b4b4dSDmitry Karpeev    Input Parameters:
13898f3b4b4dSDmitry Karpeev +  A       - The global matrix operator
13908f3b4b4dSDmitry Karpeev -  N       - the number of global subdomains requested
13918f3b4b4dSDmitry Karpeev 
13928f3b4b4dSDmitry Karpeev    Output Parameters:
13938f3b4b4dSDmitry Karpeev +  n   - the number of subdomains created on this processor
13948f3b4b4dSDmitry Karpeev -  iis - the array of index sets defining the local inner subdomains (on which the correction is applied)
13958f3b4b4dSDmitry Karpeev 
13968f3b4b4dSDmitry Karpeev    Level: advanced
13978f3b4b4dSDmitry Karpeev 
13988f3b4b4dSDmitry Karpeev    Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor.
13998f3b4b4dSDmitry Karpeev          When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local.
14008f3b4b4dSDmitry Karpeev          The resulting subdomains can be use in PCGASMSetSubdomains(pc,n,iss,NULL).  The overlapping
14018f3b4b4dSDmitry Karpeev 	 outer subdomains will be automatically generated from these according to the requested amount of
14028f3b4b4dSDmitry Karpeev 	 overlap; this is currently supported only with local subdomains.
14038f3b4b4dSDmitry Karpeev 
14048f3b4b4dSDmitry Karpeev 
14058f3b4b4dSDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
14068f3b4b4dSDmitry Karpeev 
14078f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains()
14088f3b4b4dSDmitry Karpeev @*/
1409b541e6a4SDmitry Karpeev PetscErrorCode  PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[])
14108f3b4b4dSDmitry Karpeev {
14118f3b4b4dSDmitry Karpeev   PetscMPIInt     size;
14128f3b4b4dSDmitry Karpeev   PetscErrorCode  ierr;
14138f3b4b4dSDmitry Karpeev 
14148f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
14158f3b4b4dSDmitry Karpeev   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
14168f3b4b4dSDmitry Karpeev   PetscValidPointer(iis,4);
14178f3b4b4dSDmitry Karpeev 
14188f3b4b4dSDmitry Karpeev   if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %D",N);
14198f3b4b4dSDmitry Karpeev   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
14208f3b4b4dSDmitry Karpeev   if (N >= size) {
14218f3b4b4dSDmitry Karpeev     *n = N/size + (N%size);
14228f3b4b4dSDmitry Karpeev     ierr = PCGASMCreateLocalSubdomains(A,*n,iis);CHKERRQ(ierr);
14236a4f0f73SDmitry Karpeev   } else {
14248f3b4b4dSDmitry Karpeev     ierr = PCGASMCreateStraddlingSubdomains(A,N,n,iis);CHKERRQ(ierr);
14256a4f0f73SDmitry Karpeev   }
1426b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1427b1a0a8a3SJed Brown }
1428b1a0a8a3SJed Brown 
1429b1a0a8a3SJed Brown #undef __FUNCT__
1430f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains"
1431b1a0a8a3SJed Brown /*@C
1432f746d493SDmitry Karpeev    PCGASMDestroySubdomains - Destroys the index sets created with
14338f3b4b4dSDmitry Karpeev    PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be
143406b43e7eSDmitry Karpeev    called after setting subdomains with PCGASMSetSubdomains().
1435b1a0a8a3SJed Brown 
1436b1a0a8a3SJed Brown    Collective
1437b1a0a8a3SJed Brown 
1438b1a0a8a3SJed Brown    Input Parameters:
1439b1a0a8a3SJed Brown +  n   - the number of index sets
14406a4f0f73SDmitry Karpeev .  iis - the array of inner subdomains,
14410298fd71SBarry Smith -  ois - the array of outer subdomains, can be NULL
1442b1a0a8a3SJed Brown 
14436a4f0f73SDmitry Karpeev    Level: intermediate
14446a4f0f73SDmitry Karpeev 
14456a4f0f73SDmitry Karpeev    Notes: this is merely a convenience subroutine that walks each list,
14462c112581SDmitry Karpeev    destroys each IS on the list, and then frees the list. At the end the
14472c112581SDmitry Karpeev    list pointers are set to NULL.
1448b1a0a8a3SJed Brown 
1449f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1450b1a0a8a3SJed Brown 
14518f3b4b4dSDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetSubdomains()
1452b1a0a8a3SJed Brown @*/
14532c112581SDmitry Karpeev PetscErrorCode  PCGASMDestroySubdomains(PetscInt n,IS **iis,IS **ois)
1454b1a0a8a3SJed Brown {
1455b1a0a8a3SJed Brown   PetscInt       i;
1456b1a0a8a3SJed Brown   PetscErrorCode ierr;
14575fd66863SKarl Rupp 
1458b1a0a8a3SJed Brown   PetscFunctionBegin;
1459a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1460a35b7b57SDmitry Karpeev   if (iis) {
14616a4f0f73SDmitry Karpeev     PetscValidPointer(iis,2);
14622c112581SDmitry Karpeev     if (*iis) {
14632c112581SDmitry Karpeev       PetscValidPointer(*iis,2);
1464a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {
14652c112581SDmitry Karpeev         ierr = ISDestroy(&(*iis)[i]);CHKERRQ(ierr);
1466a35b7b57SDmitry Karpeev       }
14672c112581SDmitry Karpeev       ierr = PetscFree((*iis));CHKERRQ(ierr);
14682c112581SDmitry Karpeev     }
1469a35b7b57SDmitry Karpeev   }
14706a4f0f73SDmitry Karpeev   if (ois) {
14712c112581SDmitry Karpeev     PetscValidPointer(ois,3);
14722c112581SDmitry Karpeev     if (*ois) {
14732c112581SDmitry Karpeev       PetscValidPointer(*ois,3);
1474a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {
14752c112581SDmitry Karpeev         ierr = ISDestroy(&(*ois)[i]);CHKERRQ(ierr);
1476a35b7b57SDmitry Karpeev       }
1477135757f6SDmitry Karpeev       ierr = PetscFree((*ois));CHKERRQ(ierr);
14782c112581SDmitry Karpeev     }
1479b1a0a8a3SJed Brown   }
1480b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1481b1a0a8a3SJed Brown }
1482b1a0a8a3SJed Brown 
1483af538404SDmitry Karpeev 
1484af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \
1485af538404SDmitry Karpeev   {                                                                                                       \
1486af538404SDmitry Karpeev     PetscInt first_row = first/M, last_row = last/M+1;                                                     \
1487af538404SDmitry Karpeev     /*                                                                                                    \
1488af538404SDmitry Karpeev      Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners       \
1489af538404SDmitry Karpeev      of the bounding box of the intersection of the subdomain with the local ownership range (local       \
1490af538404SDmitry Karpeev      subdomain).                                                                                          \
1491af538404SDmitry Karpeev      Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows       \
1492af538404SDmitry Karpeev      of the intersection.                                                                                 \
1493af538404SDmitry Karpeev     */                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                             \
1494af538404SDmitry Karpeev     /* ylow_loc is the grid row containing the first element of the local sumbdomain */                   \
1495eec7e3faSDmitry Karpeev     *ylow_loc = PetscMax(first_row,ylow);                                                                    \
1496af538404SDmitry Karpeev     /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1497eec7e3faSDmitry Karpeev     *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft;                                                                            \
1498af538404SDmitry Karpeev     /* yhigh_loc is the grid row above the last local subdomain element */                                                                    \
1499eec7e3faSDmitry Karpeev     *yhigh_loc = PetscMin(last_row,yhigh);                                                                                                     \
1500af538404SDmitry Karpeev     /* xright is the offset of the end of the  local subdomain within its grid row (might actually be outside the local subdomain) */         \
1501eec7e3faSDmitry Karpeev     *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright;                                                                          \
1502af538404SDmitry Karpeev     /* Now compute the size of the local subdomain n. */ \
1503c3518bceSDmitry Karpeev     *n = 0;                                               \
1504eec7e3faSDmitry Karpeev     if (*ylow_loc < *yhigh_loc) {                           \
1505af538404SDmitry Karpeev       PetscInt width = xright-xleft;                     \
1506eec7e3faSDmitry Karpeev       *n += width*(*yhigh_loc-*ylow_loc-1);                 \
1507eec7e3faSDmitry Karpeev       *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \
1508eec7e3faSDmitry Karpeev       *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \
1509af538404SDmitry Karpeev     } \
1510af538404SDmitry Karpeev   }
1511af538404SDmitry Karpeev 
1512c3518bceSDmitry Karpeev 
1513eec7e3faSDmitry Karpeev 
1514eec7e3faSDmitry Karpeev #undef __FUNCT__
1515f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D"
1516b1a0a8a3SJed Brown /*@
1517f746d493SDmitry Karpeev    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1518b1a0a8a3SJed Brown    preconditioner for a two-dimensional problem on a regular grid.
1519b1a0a8a3SJed Brown 
1520af538404SDmitry Karpeev    Collective
1521b1a0a8a3SJed Brown 
1522b1a0a8a3SJed Brown    Input Parameters:
152306b43e7eSDmitry Karpeev +  M, N               - the global number of grid points in the x and y directions
1524af538404SDmitry Karpeev .  Mdomains, Ndomains - the global number of subdomains in the x and y directions
1525b1a0a8a3SJed Brown .  dof                - degrees of freedom per node
1526b1a0a8a3SJed Brown -  overlap            - overlap in mesh lines
1527b1a0a8a3SJed Brown 
1528b1a0a8a3SJed Brown    Output Parameters:
1529af538404SDmitry Karpeev +  Nsub - the number of local subdomains created
15306a4f0f73SDmitry Karpeev .  iis  - array of index sets defining inner (nonoverlapping) subdomains
15316a4f0f73SDmitry Karpeev -  ois  - array of index sets defining outer (overlapping, if overlap > 0) subdomains
1532b1a0a8a3SJed Brown 
1533b1a0a8a3SJed Brown 
1534b1a0a8a3SJed Brown    Level: advanced
1535b1a0a8a3SJed Brown 
1536f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid
1537b1a0a8a3SJed Brown 
15388f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), PCGASMSetOverlap()
1539b1a0a8a3SJed Brown @*/
15406a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMCreateSubdomains2D(PC pc,PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap,PetscInt *nsub,IS **iis,IS **ois)
1541b1a0a8a3SJed Brown {
1542b1a0a8a3SJed Brown   PetscErrorCode ierr;
15432bbb417fSDmitry Karpeev   PetscMPIInt    size, rank;
15442bbb417fSDmitry Karpeev   PetscInt       i, j;
15452bbb417fSDmitry Karpeev   PetscInt       maxheight, maxwidth;
15462bbb417fSDmitry Karpeev   PetscInt       xstart, xleft, xright, xleft_loc, xright_loc;
15472bbb417fSDmitry Karpeev   PetscInt       ystart, ylow,  yhigh,  ylow_loc,  yhigh_loc;
1548eec7e3faSDmitry Karpeev   PetscInt       x[2][2], y[2][2], n[2];
15492bbb417fSDmitry Karpeev   PetscInt       first, last;
15502bbb417fSDmitry Karpeev   PetscInt       nidx, *idx;
15512bbb417fSDmitry Karpeev   PetscInt       ii,jj,s,q,d;
1552af538404SDmitry Karpeev   PetscInt       k,kk;
15532bbb417fSDmitry Karpeev   PetscMPIInt    color;
15542bbb417fSDmitry Karpeev   MPI_Comm       comm, subcomm;
15556a4f0f73SDmitry Karpeev   IS             **xis = 0, **is = ois, **is_local = iis;
1556d34fcf5fSBarry Smith 
1557b1a0a8a3SJed Brown   PetscFunctionBegin;
15582bbb417fSDmitry Karpeev   ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr);
15592bbb417fSDmitry Karpeev   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
15602bbb417fSDmitry Karpeev   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
15612bbb417fSDmitry Karpeev   ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr);
1562d34fcf5fSBarry Smith   if (first%dof || last%dof) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Matrix row partitioning unsuitable for domain decomposition: local row range (%D,%D) "
15632bbb417fSDmitry Karpeev                                       "does not respect the number of degrees of freedom per grid point %D", first, last, dof);
1564d34fcf5fSBarry Smith 
1565af538404SDmitry Karpeev   /* Determine the number of domains with nonzero intersections with the local ownership range. */
15662bbb417fSDmitry Karpeev   s      = 0;
1567af538404SDmitry Karpeev   ystart = 0;
1568af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1569af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1570af538404SDmitry Karpeev     if (maxheight < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the vertical directon for mesh height %D", Ndomains, N);
1571eec7e3faSDmitry Karpeev     /* Vertical domain limits with an overlap. */
1572eec7e3faSDmitry Karpeev     ylow   = PetscMax(ystart - overlap,0);
1573eec7e3faSDmitry Karpeev     yhigh  = PetscMin(ystart + maxheight + overlap,N);
1574b1a0a8a3SJed Brown     xstart = 0;
1575af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1576af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1577af538404SDmitry Karpeev       if (maxwidth < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the horizontal direction for mesh width %D", Mdomains, M);
1578eec7e3faSDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1579eec7e3faSDmitry Karpeev       xleft  = PetscMax(xstart - overlap,0);
1580eec7e3faSDmitry Karpeev       xright = PetscMin(xstart + maxwidth + overlap,M);
1581af538404SDmitry Karpeev       /*
1582af538404SDmitry Karpeev          Determine whether this subdomain intersects this processor's ownership range of pc->pmat.
1583af538404SDmitry Karpeev       */
1584c3518bceSDmitry Karpeev       PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
15852fa5cd67SKarl Rupp       if (nidx) ++s;
1586af538404SDmitry Karpeev       xstart += maxwidth;
1587af538404SDmitry Karpeev     } /* for (i = 0; i < Mdomains; ++i) */
1588af538404SDmitry Karpeev     ystart += maxheight;
1589af538404SDmitry Karpeev   } /* for (j = 0; j < Ndomains; ++j) */
15902fa5cd67SKarl Rupp 
1591af538404SDmitry Karpeev   /* Now we can allocate the necessary number of ISs. */
1592af538404SDmitry Karpeev   *nsub  = s;
1593854ce69bSBarry Smith   ierr   = PetscMalloc1(*nsub,is);CHKERRQ(ierr);
1594854ce69bSBarry Smith   ierr   = PetscMalloc1(*nsub,is_local);CHKERRQ(ierr);
1595af538404SDmitry Karpeev   s      = 0;
1596af538404SDmitry Karpeev   ystart = 0;
1597af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1598af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1599af538404SDmitry Karpeev     if (maxheight < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the vertical directon for mesh height %D", Ndomains, N);
1600af538404SDmitry Karpeev     /* Vertical domain limits with an overlap. */
1601eec7e3faSDmitry Karpeev     y[0][0] = PetscMax(ystart - overlap,0);
1602eec7e3faSDmitry Karpeev     y[0][1] = PetscMin(ystart + maxheight + overlap,N);
1603eec7e3faSDmitry Karpeev     /* Vertical domain limits without an overlap. */
1604eec7e3faSDmitry Karpeev     y[1][0] = ystart;
1605eec7e3faSDmitry Karpeev     y[1][1] = PetscMin(ystart + maxheight,N);
1606eec7e3faSDmitry Karpeev     xstart  = 0;
1607af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1608af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1609af538404SDmitry Karpeev       if (maxwidth < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the horizontal direction for mesh width %D", Mdomains, M);
1610af538404SDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1611eec7e3faSDmitry Karpeev       x[0][0] = PetscMax(xstart - overlap,0);
1612eec7e3faSDmitry Karpeev       x[0][1] = PetscMin(xstart + maxwidth + overlap,M);
1613eec7e3faSDmitry Karpeev       /* Horizontal domain limits without an overlap. */
1614eec7e3faSDmitry Karpeev       x[1][0] = xstart;
1615eec7e3faSDmitry Karpeev       x[1][1] = PetscMin(xstart+maxwidth,M);
16162bbb417fSDmitry Karpeev       /*
16172bbb417fSDmitry Karpeev          Determine whether this domain intersects this processor's ownership range of pc->pmat.
16182bbb417fSDmitry Karpeev          Do this twice: first for the domains with overlaps, and once without.
16192bbb417fSDmitry Karpeev          During the first pass create the subcommunicators, and use them on the second pass as well.
16202bbb417fSDmitry Karpeev       */
16212bbb417fSDmitry Karpeev       for (q = 0; q < 2; ++q) {
16222bbb417fSDmitry Karpeev         /*
16232bbb417fSDmitry Karpeev           domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
16242bbb417fSDmitry Karpeev           according to whether the domain with an overlap or without is considered.
16252bbb417fSDmitry Karpeev         */
16262bbb417fSDmitry Karpeev         xleft = x[q][0]; xright = x[q][1];
16272bbb417fSDmitry Karpeev         ylow  = y[q][0]; yhigh  = y[q][1];
1628c3518bceSDmitry Karpeev         PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1629af538404SDmitry Karpeev         nidx *= dof;
1630eec7e3faSDmitry Karpeev         n[q]  = nidx;
16312bbb417fSDmitry Karpeev         /*
1632eec7e3faSDmitry Karpeev          Based on the counted number of indices in the local domain *with an overlap*,
1633af538404SDmitry Karpeev          construct a subcommunicator of all the processors supporting this domain.
1634eec7e3faSDmitry Karpeev          Observe that a domain with an overlap might have nontrivial local support,
1635eec7e3faSDmitry Karpeev          while the domain without an overlap might not.  Hence, the decision to participate
1636eec7e3faSDmitry Karpeev          in the subcommunicator must be based on the domain with an overlap.
16372bbb417fSDmitry Karpeev          */
16382bbb417fSDmitry Karpeev         if (q == 0) {
16392fa5cd67SKarl Rupp           if (nidx) color = 1;
16402fa5cd67SKarl Rupp           else color = MPI_UNDEFINED;
16412fa5cd67SKarl Rupp 
16422bbb417fSDmitry Karpeev           ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr);
16432bbb417fSDmitry Karpeev         }
1644af538404SDmitry Karpeev         /*
1645eec7e3faSDmitry Karpeev          Proceed only if the number of local indices *with an overlap* is nonzero.
1646af538404SDmitry Karpeev          */
1647eec7e3faSDmitry Karpeev         if (n[0]) {
16482fa5cd67SKarl Rupp           if (q == 0) xis = is;
1649af538404SDmitry Karpeev           if (q == 1) {
1650af538404SDmitry Karpeev             /*
1651eec7e3faSDmitry Karpeev              The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1652af538404SDmitry Karpeev              Moreover, if the overlap is zero, the two ISs are identical.
1653af538404SDmitry Karpeev              */
1654b1a0a8a3SJed Brown             if (overlap == 0) {
1655eec7e3faSDmitry Karpeev               (*is_local)[s] = (*is)[s];
16562bbb417fSDmitry Karpeev               ierr           = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr);
16572bbb417fSDmitry Karpeev               continue;
1658d34fcf5fSBarry Smith             } else {
16596a4f0f73SDmitry Karpeev               xis     = is_local;
1660eec7e3faSDmitry Karpeev               subcomm = ((PetscObject)(*is)[s])->comm;
16612bbb417fSDmitry Karpeev             }
1662af538404SDmitry Karpeev           } /* if (q == 1) */
16630298fd71SBarry Smith           idx  = NULL;
1664785e854fSJed Brown           ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr);
1665eec7e3faSDmitry Karpeev           if (nidx) {
16662bbb417fSDmitry Karpeev             k = 0;
16672bbb417fSDmitry Karpeev             for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1668af538404SDmitry Karpeev               PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft;
1669af538404SDmitry Karpeev               PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright;
1670af538404SDmitry Karpeev               kk = dof*(M*jj + x0);
16712bbb417fSDmitry Karpeev               for (ii=x0; ii<x1; ++ii) {
16722bbb417fSDmitry Karpeev                 for (d = 0; d < dof; ++d) {
16732bbb417fSDmitry Karpeev                   idx[k++] = kk++;
1674b1a0a8a3SJed Brown                 }
1675b1a0a8a3SJed Brown               }
1676b1a0a8a3SJed Brown             }
1677eec7e3faSDmitry Karpeev           }
16786a4f0f73SDmitry Karpeev           ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr);
1679eec7e3faSDmitry Karpeev         }/* if (n[0]) */
16802bbb417fSDmitry Karpeev       }/* for (q = 0; q < 2; ++q) */
16812fa5cd67SKarl Rupp       if (n[0]) ++s;
1682af538404SDmitry Karpeev       xstart += maxwidth;
1683af538404SDmitry Karpeev     } /* for (i = 0; i < Mdomains; ++i) */
16842bbb417fSDmitry Karpeev     ystart += maxheight;
1685af538404SDmitry Karpeev   } /* for (j = 0; j < Ndomains; ++j) */
1686b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1687b1a0a8a3SJed Brown }
1688b1a0a8a3SJed Brown 
1689b1a0a8a3SJed Brown #undef __FUNCT__
169006b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains"
1691b1a0a8a3SJed Brown /*@C
169206b43e7eSDmitry Karpeev     PCGASMGetSubdomains - Gets the subdomains supported on this processor
169306b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
1694b1a0a8a3SJed Brown 
1695b1a0a8a3SJed Brown     Not Collective
1696b1a0a8a3SJed Brown 
1697b1a0a8a3SJed Brown     Input Parameter:
1698b1a0a8a3SJed Brown .   pc - the preconditioner context
1699b1a0a8a3SJed Brown 
1700b1a0a8a3SJed Brown     Output Parameters:
1701b1a0a8a3SJed Brown +   n   - the number of subdomains for this processor (default value = 1)
17020298fd71SBarry Smith .   iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL)
17030298fd71SBarry Smith -   ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL)
1704b1a0a8a3SJed Brown 
1705b1a0a8a3SJed Brown 
1706b1a0a8a3SJed Brown     Notes:
17076a4f0f73SDmitry Karpeev     The user is responsible for destroying the ISs and freeing the returned arrays.
1708b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
1709b1a0a8a3SJed Brown 
1710b1a0a8a3SJed Brown     Level: advanced
1711b1a0a8a3SJed Brown 
171206b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz
1713b1a0a8a3SJed Brown 
17148f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), PCGASMCreateSubdomains2D(),
17158f3b4b4dSDmitry Karpeev           PCGASMSetSubdomains(), PCGASMGetSubmatrices()
1716b1a0a8a3SJed Brown @*/
17176a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[])
1718b1a0a8a3SJed Brown {
1719f746d493SDmitry Karpeev   PC_GASM        *osm;
1720b1a0a8a3SJed Brown   PetscErrorCode ierr;
1721b1a0a8a3SJed Brown   PetscBool      match;
17226a4f0f73SDmitry Karpeev   PetscInt       i;
17235fd66863SKarl Rupp 
1724b1a0a8a3SJed Brown   PetscFunctionBegin;
1725b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1726251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1727ce94432eSBarry Smith   if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1728f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1729ab3e923fSDmitry Karpeev   if (n) *n = osm->n;
17306a4f0f73SDmitry Karpeev   if (iis) {
1731785e854fSJed Brown     ierr = PetscMalloc1(osm->n, iis);CHKERRQ(ierr);
17326a4f0f73SDmitry Karpeev   }
17336a4f0f73SDmitry Karpeev   if (ois) {
1734785e854fSJed Brown     ierr = PetscMalloc1(osm->n, ois);CHKERRQ(ierr);
17356a4f0f73SDmitry Karpeev   }
17366a4f0f73SDmitry Karpeev   if (iis || ois) {
17376a4f0f73SDmitry Karpeev     for (i = 0; i < osm->n; ++i) {
17386a4f0f73SDmitry Karpeev       if (iis) (*iis)[i] = osm->iis[i];
17396a4f0f73SDmitry Karpeev       if (ois) (*ois)[i] = osm->ois[i];
17406a4f0f73SDmitry Karpeev     }
1741b1a0a8a3SJed Brown   }
1742b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1743b1a0a8a3SJed Brown }
1744b1a0a8a3SJed Brown 
1745b1a0a8a3SJed Brown #undef __FUNCT__
174606b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices"
1747b1a0a8a3SJed Brown /*@C
174806b43e7eSDmitry Karpeev     PCGASMGetSubmatrices - Gets the local submatrices (for this processor
1749b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1750b1a0a8a3SJed Brown 
1751b1a0a8a3SJed Brown     Not Collective
1752b1a0a8a3SJed Brown 
1753b1a0a8a3SJed Brown     Input Parameter:
1754b1a0a8a3SJed Brown .   pc - the preconditioner context
1755b1a0a8a3SJed Brown 
1756b1a0a8a3SJed Brown     Output Parameters:
1757b1a0a8a3SJed Brown +   n   - the number of matrices for this processor (default value = 1)
1758b1a0a8a3SJed Brown -   mat - the matrices
1759b1a0a8a3SJed Brown 
176006b43e7eSDmitry Karpeev     Notes: matrices returned by this routine have the same communicators as the index sets (IS)
17618f3b4b4dSDmitry Karpeev            used to define subdomains in PCGASMSetSubdomains()
1762b1a0a8a3SJed Brown     Level: advanced
1763b1a0a8a3SJed Brown 
1764f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi
1765b1a0a8a3SJed Brown 
17668f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(),
176706b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains()
1768b1a0a8a3SJed Brown @*/
176906b43e7eSDmitry Karpeev PetscErrorCode  PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1770b1a0a8a3SJed Brown {
1771f746d493SDmitry Karpeev   PC_GASM        *osm;
1772b1a0a8a3SJed Brown   PetscErrorCode ierr;
1773b1a0a8a3SJed Brown   PetscBool      match;
1774b1a0a8a3SJed Brown 
1775b1a0a8a3SJed Brown   PetscFunctionBegin;
1776b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1777b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1778b1a0a8a3SJed Brown   if (mat) PetscValidPointer(mat,3);
1779ce94432eSBarry Smith   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1780251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1781ce94432eSBarry Smith   if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1782f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1783ab3e923fSDmitry Karpeev   if (n) *n = osm->n;
1784b1a0a8a3SJed Brown   if (mat) *mat = osm->pmat;
1785b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1786b1a0a8a3SJed Brown }
1787d709ea83SDmitry Karpeev 
1788d709ea83SDmitry Karpeev #undef __FUNCT__
17898f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMSetUseDMSubdomains"
1790d709ea83SDmitry Karpeev /*@
17918f3b4b4dSDmitry Karpeev     PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1792d709ea83SDmitry Karpeev     Logically Collective
1793d709ea83SDmitry Karpeev 
1794d709ea83SDmitry Karpeev     Input Parameter:
1795d709ea83SDmitry Karpeev +   pc  - the preconditioner
1796d709ea83SDmitry Karpeev -   flg - boolean indicating whether to use subdomains defined by the DM
1797d709ea83SDmitry Karpeev 
1798d709ea83SDmitry Karpeev     Options Database Key:
17998f3b4b4dSDmitry Karpeev .   -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains
1800d709ea83SDmitry Karpeev 
1801d709ea83SDmitry Karpeev     Level: intermediate
1802d709ea83SDmitry Karpeev 
1803d709ea83SDmitry Karpeev     Notes:
18048f3b4b4dSDmitry Karpeev     PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(),
18058f3b4b4dSDmitry Karpeev     so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap()
18068f3b4b4dSDmitry Karpeev     automatically turns the latter off.
1807d709ea83SDmitry Karpeev 
1808d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1809d709ea83SDmitry Karpeev 
18108f3b4b4dSDmitry Karpeev .seealso: PCGASMGetUseDMSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap()
1811d709ea83SDmitry Karpeev           PCGASMCreateSubdomains2D()
1812d709ea83SDmitry Karpeev @*/
18138f3b4b4dSDmitry Karpeev PetscErrorCode  PCGASMSetUseDMSubdomains(PC pc,PetscBool flg)
1814d709ea83SDmitry Karpeev {
1815d709ea83SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
1816d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1817d709ea83SDmitry Karpeev   PetscBool      match;
1818d709ea83SDmitry Karpeev 
1819d709ea83SDmitry Karpeev   PetscFunctionBegin;
1820d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1821d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
1822d709ea83SDmitry Karpeev   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1823d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1824d709ea83SDmitry Karpeev   if (match) {
18258f3b4b4dSDmitry Karpeev     if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) {
1826d709ea83SDmitry Karpeev       osm->dm_subdomains = flg;
1827d709ea83SDmitry Karpeev     }
18288f3b4b4dSDmitry Karpeev   }
1829d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1830d709ea83SDmitry Karpeev }
1831d709ea83SDmitry Karpeev 
1832d709ea83SDmitry Karpeev #undef __FUNCT__
18338f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMGetUseDMSubdomains"
1834d709ea83SDmitry Karpeev /*@
18358f3b4b4dSDmitry Karpeev     PCGASMGetUseDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1836d709ea83SDmitry Karpeev     Not Collective
1837d709ea83SDmitry Karpeev 
1838d709ea83SDmitry Karpeev     Input Parameter:
1839d709ea83SDmitry Karpeev .   pc  - the preconditioner
1840d709ea83SDmitry Karpeev 
1841d709ea83SDmitry Karpeev     Output Parameter:
1842d709ea83SDmitry Karpeev .   flg - boolean indicating whether to use subdomains defined by the DM
1843d709ea83SDmitry Karpeev 
1844d709ea83SDmitry Karpeev     Level: intermediate
1845d709ea83SDmitry Karpeev 
1846d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1847d709ea83SDmitry Karpeev 
18488f3b4b4dSDmitry Karpeev .seealso: PCGASMSetUseDMSubdomains(), PCGASMSetOverlap()
1849d709ea83SDmitry Karpeev           PCGASMCreateSubdomains2D()
1850d709ea83SDmitry Karpeev @*/
18518f3b4b4dSDmitry Karpeev PetscErrorCode  PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg)
1852d709ea83SDmitry Karpeev {
1853d709ea83SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
1854d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1855d709ea83SDmitry Karpeev   PetscBool      match;
1856d709ea83SDmitry Karpeev 
1857d709ea83SDmitry Karpeev   PetscFunctionBegin;
1858d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1859d709ea83SDmitry Karpeev   PetscValidPointer(flg,2);
1860d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1861d709ea83SDmitry Karpeev   if (match) {
1862d709ea83SDmitry Karpeev     if (flg) *flg = osm->dm_subdomains;
1863d709ea83SDmitry Karpeev   }
1864d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1865d709ea83SDmitry Karpeev }
1866