xref: /petsc/src/ksp/pc/impls/gasm/gasm.c (revision 1575c14d35e6562ef67dce5e5a97f4a0489f95e5)
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);
81*1575c14dSBarry Smith   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
824bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
834bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
84*1575c14dSBarry Smith   ierr = PetscViewerASCIIPopSynchronized(viewer);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);
106*1575c14dSBarry Smith   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
107af538404SDmitry Karpeev   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
108af538404SDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
109*1575c14dSBarry Smith   ierr = PetscViewerASCIIPopSynchronized(viewer);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) {
1493f08860eSBarry Smith         ierr = PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
1504bde246dSDmitry Karpeev         ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr);
1513f08860eSBarry Smith         ierr = PetscViewerRestoreSubViewer(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);
210*1575c14dSBarry Smith     ierr = PetscViewerASCIIPushSynchronized(viewer);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);
213*1575c14dSBarry Smith     ierr = PetscViewerASCIIPopSynchronized(viewer);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);
2303f08860eSBarry Smith           ierr = PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
2316a4f0f73SDmitry Karpeev           ierr = ISGetLocalSize(osm->ois[d],&bsz);CHKERRQ(ierr);
2328f3b4b4dSDmitry 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);
2336a4f0f73SDmitry Karpeev           ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
23406b43e7eSDmitry Karpeev           if (view_subdomains) {
23506b43e7eSDmitry Karpeev             ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr);
23606b43e7eSDmitry Karpeev           }
2376a4f0f73SDmitry Karpeev           if (!pc->setupcalled) {
238*1575c14dSBarry Smith             ierr = PetscViewerASCIIPrintf(sviewer, "Solver not set up yet: PCSetUp() not yet called\n");CHKERRQ(ierr);
2392fa5cd67SKarl Rupp           } else {
24006b43e7eSDmitry Karpeev             ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr);
2416a4f0f73SDmitry Karpeev           }
24206b43e7eSDmitry Karpeev           ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
2437b23a99aSBarry Smith           ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
2443f08860eSBarry Smith           ierr = PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
24506b43e7eSDmitry Karpeev           ++l;
246b1a0a8a3SJed Brown         }
24706b43e7eSDmitry Karpeev       }
248ce94432eSBarry Smith       ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
249b1a0a8a3SJed Brown     }
2508f3b4b4dSDmitry Karpeev     ierr = PetscFree2(numbering,permutation);CHKERRQ(ierr);
251b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
252b1a0a8a3SJed Brown   }
253*1575c14dSBarry Smith   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
254b1a0a8a3SJed Brown   PetscFunctionReturn(0);
255b1a0a8a3SJed Brown }
256b1a0a8a3SJed Brown 
2578f3b4b4dSDmitry Karpeev PETSC_INTERN PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]);
258af538404SDmitry Karpeev 
259b1a0a8a3SJed Brown #undef __FUNCT__
260f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM"
261f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc)
262b1a0a8a3SJed Brown {
263f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
264b1a0a8a3SJed Brown   PetscErrorCode ierr;
265b1a0a8a3SJed Brown   PetscBool      symset,flg;
26688389755SJed Brown   PetscInt       i;
267c10bc1a1SDmitry Karpeev   PetscMPIInt    rank, size;
268b1a0a8a3SJed Brown   MatReuse       scall = MAT_REUSE_MATRIX;
269b1a0a8a3SJed Brown   KSP            ksp;
270b1a0a8a3SJed Brown   PC             subpc;
271b1a0a8a3SJed Brown   const char     *prefix,*pprefix;
272f746d493SDmitry Karpeev   Vec            x,y;
2736a4f0f73SDmitry Karpeev   PetscInt       oni;       /* Number of indices in the i-th local outer subdomain.               */
2746a4f0f73SDmitry Karpeev   const PetscInt *oidxi;    /* Indices from the i-th subdomain local outer subdomain.             */
2756a4f0f73SDmitry Karpeev   PetscInt       on;        /* Number of indices in the disjoint union of local outer subdomains. */
2766a4f0f73SDmitry Karpeev   PetscInt       *oidx;     /* Indices in the disjoint union of local outer subdomains. */
2776a4f0f73SDmitry Karpeev   IS             gois;      /* Disjoint union the global indices of outer subdomains.             */
2786a4f0f73SDmitry Karpeev   IS             goid;      /* Identity IS of the size of the disjoint union of outer subdomains. */
279f746d493SDmitry Karpeev   PetscScalar    *gxarray, *gyarray;
2802fa5cd67SKarl Rupp   PetscInt       gofirst;   /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */
2818f3b4b4dSDmitry Karpeev   PetscInt       num_subdomains    = 0;
2820298fd71SBarry Smith   DM             *subdomain_dm     = NULL;
2838f3b4b4dSDmitry Karpeev   char           **subdomain_names = NULL;
2848f3b4b4dSDmitry Karpeev 
285b1a0a8a3SJed Brown 
286b1a0a8a3SJed Brown   PetscFunctionBegin;
287ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
288ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
289b1a0a8a3SJed Brown   if (!pc->setupcalled) {
290b1a0a8a3SJed Brown 
291b1a0a8a3SJed Brown     if (!osm->type_set) {
292b1a0a8a3SJed Brown       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
2932fa5cd67SKarl Rupp       if (symset && flg) osm->type = PC_GASM_BASIC;
294b1a0a8a3SJed Brown     }
295b1a0a8a3SJed Brown 
2968f3b4b4dSDmitry Karpeev     if (osm->n == PETSC_DETERMINE) {
2978f3b4b4dSDmitry Karpeev       if (osm->N != PETSC_DETERMINE) {
2988f3b4b4dSDmitry Karpeev 	/* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */
2998f3b4b4dSDmitry Karpeev 	ierr = PCGASMCreateSubdomains(pc->pmat,osm->N,&osm->n,&osm->iis);CHKERRQ(ierr);
3008f3b4b4dSDmitry Karpeev       } else if (osm->dm_subdomains && pc->dm) {
3018f3b4b4dSDmitry Karpeev 	/* try pc->dm next, if allowed */
3028f3b4b4dSDmitry Karpeev 	PetscInt  d;
303a35b7b57SDmitry Karpeev 	IS        *inner_subdomain_is, *outer_subdomain_is;
304a35b7b57SDmitry Karpeev 	ierr = DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);CHKERRQ(ierr);
305a35b7b57SDmitry Karpeev 	if (num_subdomains) {
306a35b7b57SDmitry Karpeev 	  ierr = PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);CHKERRQ(ierr);
30769ca1f37SDmitry Karpeev 	}
308a35b7b57SDmitry Karpeev 	for (d = 0; d < num_subdomains; ++d) {
309a35b7b57SDmitry Karpeev 	  if (inner_subdomain_is) {ierr = ISDestroy(&inner_subdomain_is[d]);CHKERRQ(ierr);}
310a35b7b57SDmitry Karpeev 	  if (outer_subdomain_is) {ierr = ISDestroy(&outer_subdomain_is[d]);CHKERRQ(ierr);}
311fdc48646SDmitry Karpeev 	}
312a35b7b57SDmitry Karpeev 	ierr = PetscFree(inner_subdomain_is);CHKERRQ(ierr);
313a35b7b57SDmitry Karpeev 	ierr = PetscFree(outer_subdomain_is);CHKERRQ(ierr);
3148f3b4b4dSDmitry Karpeev       } else {
3158f3b4b4dSDmitry Karpeev 	/* still no subdomains; use one per processor */
31644fe18e5SDmitry Karpeev 	osm->nmax = osm->n = 1;
317ce94432eSBarry Smith 	ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
318f746d493SDmitry Karpeev 	osm->N    = size;
3198f3b4b4dSDmitry Karpeev 	ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr);
320fdc48646SDmitry Karpeev       }
32106b43e7eSDmitry Karpeev     }
322a35b7b57SDmitry Karpeev     if (!osm->iis) {
323a35b7b57SDmitry Karpeev       /*
3248f3b4b4dSDmitry Karpeev        osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied.
3258f3b4b4dSDmitry Karpeev        We create the requisite number of local inner subdomains and then expand them into
3268f3b4b4dSDmitry Karpeev        out subdomains, if necessary.
327a35b7b57SDmitry Karpeev        */
3288f3b4b4dSDmitry Karpeev       ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr);
329a35b7b57SDmitry Karpeev     }
3308f3b4b4dSDmitry Karpeev     if (!osm->ois) {
3318f3b4b4dSDmitry Karpeev       /*
3328f3b4b4dSDmitry Karpeev 	Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap
3338f3b4b4dSDmitry Karpeev 	has been requested, copy the inner subdomains over so they can be modified.
3348f3b4b4dSDmitry Karpeev       */
3358f3b4b4dSDmitry Karpeev       ierr = PetscMalloc1(osm->n,&osm->ois);CHKERRQ(ierr);
3368f3b4b4dSDmitry Karpeev       for (i=0; i<osm->n; ++i) {
3378f3b4b4dSDmitry Karpeev 	if (osm->overlap > 0) { /* With positive overlap, osm->iis[i] will be modified */
3388f3b4b4dSDmitry Karpeev 	  ierr = ISDuplicate(osm->iis[i],(osm->ois)+i);CHKERRQ(ierr);
3398f3b4b4dSDmitry Karpeev 	  ierr = ISCopy(osm->iis[i],osm->ois[i]);CHKERRQ(ierr);
3408f3b4b4dSDmitry Karpeev 	} else {
3418f3b4b4dSDmitry Karpeev 	  ierr      = PetscObjectReference((PetscObject)((osm->iis)[i]));CHKERRQ(ierr);
3428f3b4b4dSDmitry Karpeev 	  osm->ois[i] = osm->iis[i];
3438f3b4b4dSDmitry Karpeev 	}
3448f3b4b4dSDmitry Karpeev       }
3458f3b4b4dSDmitry Karpeev       if (osm->overlap > 0) {
3468f3b4b4dSDmitry Karpeev 	/* Extend the "overlapping" regions by a number of steps */
3478f3b4b4dSDmitry Karpeev 	ierr = MatIncreaseOverlap(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr);
3488f3b4b4dSDmitry Karpeev       }
349b1a0a8a3SJed Brown     }
3506a4f0f73SDmitry Karpeev 
3518f3b4b4dSDmitry Karpeev     /* Now the subdomains are defined.  Determine their global and max local numbers, if necessary. */
3528f3b4b4dSDmitry Karpeev     if (osm->nmax == PETSC_DETERMINE) {
3538f3b4b4dSDmitry Karpeev       PetscMPIInt inwork,outwork;
3548f3b4b4dSDmitry Karpeev       /* determine global number of subdomains and the max number of local subdomains */
3558f3b4b4dSDmitry Karpeev       inwork = osm->n;
3568f3b4b4dSDmitry Karpeev       ierr       = MPI_Allreduce(&inwork,&outwork,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
3578f3b4b4dSDmitry Karpeev       osm->nmax  = outwork;
3588f3b4b4dSDmitry Karpeev     }
3598f3b4b4dSDmitry Karpeev     if (osm->N == PETSC_DETERMINE) {
3608f3b4b4dSDmitry Karpeev       /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */
3618f3b4b4dSDmitry Karpeev       ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,&osm->N,NULL);CHKERRQ(ierr);
3628f3b4b4dSDmitry Karpeev     }
3638f3b4b4dSDmitry Karpeev 
364b1a0a8a3SJed Brown 
365b1a0a8a3SJed Brown     if (osm->sort_indices) {
366f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
3676a4f0f73SDmitry Karpeev         ierr = ISSort(osm->ois[i]);CHKERRQ(ierr);
3686a4f0f73SDmitry Karpeev         ierr = ISSort(osm->iis[i]);CHKERRQ(ierr);
369b1a0a8a3SJed Brown       }
370b1a0a8a3SJed Brown     }
3718f3b4b4dSDmitry Karpeev 
3728f3b4b4dSDmitry Karpeev     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
3738f3b4b4dSDmitry Karpeev     ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr);
3748f3b4b4dSDmitry Karpeev 
3756a4f0f73SDmitry Karpeev     /*
3766a4f0f73SDmitry Karpeev      Merge the ISs, create merged vectors and restrictions.
3776a4f0f73SDmitry Karpeev      */
3786a4f0f73SDmitry Karpeev     /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */
3796a4f0f73SDmitry Karpeev     on = 0;
380f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
3816a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
3826a4f0f73SDmitry Karpeev       on  += oni;
383f746d493SDmitry Karpeev     }
384785e854fSJed Brown     ierr = PetscMalloc1(on, &oidx);CHKERRQ(ierr);
3856a4f0f73SDmitry Karpeev     on   = 0;
386f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
3876a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
3886a4f0f73SDmitry Karpeev       ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr);
3896a4f0f73SDmitry Karpeev       ierr = PetscMemcpy(oidx+on, oidxi, sizeof(PetscInt)*oni);CHKERRQ(ierr);
3906a4f0f73SDmitry Karpeev       ierr = ISRestoreIndices(osm->ois[i], &oidxi);CHKERRQ(ierr);
3916a4f0f73SDmitry Karpeev       on  += oni;
392f746d493SDmitry Karpeev     }
3936a4f0f73SDmitry Karpeev     ierr = ISCreateGeneral(((PetscObject)(pc))->comm, on, oidx, PETSC_OWN_POINTER, &gois);CHKERRQ(ierr);
3942a7a6963SBarry Smith     ierr = MatCreateVecs(pc->pmat,&x,&y);CHKERRQ(ierr);
395ce94432eSBarry Smith     ierr = VecCreateMPI(PetscObjectComm((PetscObject)pc), on, PETSC_DECIDE, &osm->gx);CHKERRQ(ierr);
396f746d493SDmitry Karpeev     ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr);
3970298fd71SBarry Smith     ierr = VecGetOwnershipRange(osm->gx, &gofirst, NULL);CHKERRQ(ierr);
398ce94432eSBarry Smith     ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),on,gofirst,1, &goid);CHKERRQ(ierr);
3996a4f0f73SDmitry Karpeev     ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr);
4006a4f0f73SDmitry Karpeev     ierr = VecDestroy(&x);CHKERRQ(ierr);
4017105d80bSDmitry Karpeev     ierr = ISDestroy(&gois);CHKERRQ(ierr);
4022fa5cd67SKarl Rupp 
4036a4f0f73SDmitry Karpeev     /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */
4042fa5cd67SKarl Rupp     {
4052fa5cd67SKarl Rupp       PetscInt ini;           /* Number of indices the i-th a local inner subdomain. */
4066a4f0f73SDmitry Karpeev       PetscInt in;            /* Number of indices in the disjoint uniont of local inner subdomains. */
4076a4f0f73SDmitry Karpeev       PetscInt *iidx;         /* Global indices in the merged local inner subdomain. */
4086a4f0f73SDmitry Karpeev       PetscInt *ioidx;        /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
4096a4f0f73SDmitry Karpeev       IS       giis;          /* IS for the disjoint union of inner subdomains. */
4106a4f0f73SDmitry Karpeev       IS       giois;         /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
411f746d493SDmitry Karpeev       /**/
4126a4f0f73SDmitry Karpeev       in = 0;
413f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
4146a4f0f73SDmitry Karpeev         ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr);
4156a4f0f73SDmitry Karpeev         in  += ini;
416f746d493SDmitry Karpeev       }
417785e854fSJed Brown       ierr = PetscMalloc1(in, &iidx);CHKERRQ(ierr);
418785e854fSJed Brown       ierr = PetscMalloc1(in, &ioidx);CHKERRQ(ierr);
4190298fd71SBarry Smith       ierr = VecGetOwnershipRange(osm->gx,&gofirst, NULL);CHKERRQ(ierr);
4206a4f0f73SDmitry Karpeev       in   = 0;
4216a4f0f73SDmitry Karpeev       on   = 0;
422f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
4236a4f0f73SDmitry Karpeev         const PetscInt         *iidxi; /* Global indices of the i-th local inner subdomain. */
4246a4f0f73SDmitry Karpeev         ISLocalToGlobalMapping ltogi; /* Map from global to local indices of the i-th outer local subdomain. */
4256a4f0f73SDmitry Karpeev         PetscInt               *ioidxi; /* Local indices of the i-th local inner subdomain within the local outer subdomain. */
4266a4f0f73SDmitry Karpeev         PetscInt               ioni;  /* Number of indices in ioidxi; if ioni != ini the inner subdomain is not a subdomain of the outer subdomain (error). */
4276a4f0f73SDmitry Karpeev         PetscInt               k;
4286a4f0f73SDmitry Karpeev         ierr   = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr);
4296a4f0f73SDmitry Karpeev         ierr   = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
4306a4f0f73SDmitry Karpeev         ierr   = ISGetIndices(osm->iis[i],&iidxi);CHKERRQ(ierr);
4316a4f0f73SDmitry Karpeev         ierr   = PetscMemcpy(iidx+in, iidxi, sizeof(PetscInt)*ini);CHKERRQ(ierr);
4326a4f0f73SDmitry Karpeev         ierr   = ISLocalToGlobalMappingCreateIS(osm->ois[i],&ltogi);CHKERRQ(ierr);
4336a4f0f73SDmitry Karpeev         ioidxi = ioidx+in;
4346a4f0f73SDmitry Karpeev         ierr   = ISGlobalToLocalMappingApply(ltogi,IS_GTOLM_DROP,ini,iidxi,&ioni,ioidxi);CHKERRQ(ierr);
4356a4f0f73SDmitry Karpeev         ierr   = ISLocalToGlobalMappingDestroy(&ltogi);CHKERRQ(ierr);
4366a4f0f73SDmitry Karpeev         ierr   = ISRestoreIndices(osm->iis[i], &iidxi);CHKERRQ(ierr);
4376a4f0f73SDmitry 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);
4382fa5cd67SKarl Rupp         for (k = 0; k < ini; ++k) ioidxi[k] += gofirst+on;
4396a4f0f73SDmitry Karpeev         in += ini;
4406a4f0f73SDmitry Karpeev         on += oni;
441f746d493SDmitry Karpeev       }
442ce94432eSBarry Smith       ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, iidx,  PETSC_OWN_POINTER, &giis);CHKERRQ(ierr);
443ce94432eSBarry Smith       ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, ioidx, PETSC_OWN_POINTER, &giois);CHKERRQ(ierr);
4446a4f0f73SDmitry Karpeev       ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr);
4456a4f0f73SDmitry Karpeev       ierr = VecDestroy(&y);CHKERRQ(ierr);
4466a4f0f73SDmitry Karpeev       ierr = ISDestroy(&giis);CHKERRQ(ierr);
4476a4f0f73SDmitry Karpeev       ierr = ISDestroy(&giois);CHKERRQ(ierr);
448b1a0a8a3SJed Brown     }
4496a4f0f73SDmitry Karpeev     ierr = ISDestroy(&goid);CHKERRQ(ierr);
4502fa5cd67SKarl Rupp 
451f746d493SDmitry Karpeev     /* Create the subdomain work vectors. */
452785e854fSJed Brown     ierr = PetscMalloc1(osm->n,&osm->x);CHKERRQ(ierr);
453785e854fSJed Brown     ierr = PetscMalloc1(osm->n,&osm->y);CHKERRQ(ierr);
454f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr);
455f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr);
4566a4f0f73SDmitry Karpeev     for (i=0, on=0; i<osm->n; ++i, on += oni) {
4576a4f0f73SDmitry Karpeev       PetscInt oNi;
4586a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
4596a4f0f73SDmitry Karpeev       ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr);
4606a4f0f73SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr);
4616a4f0f73SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr);
462b1a0a8a3SJed Brown     }
463f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr);
464f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr);
4658f3b4b4dSDmitry Karpeev     /* Create the subdomain solvers */
466785e854fSJed Brown     ierr = PetscMalloc1(osm->n,&osm->ksp);CHKERRQ(ierr);
4678f3b4b4dSDmitry Karpeev     for (i=0; i<osm->n; i++) {
4688f3b4b4dSDmitry Karpeev       char subprefix[PETSC_MAX_PATH_LEN+1];
4696a4f0f73SDmitry Karpeev       ierr        = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr);
4703bb1ff40SBarry Smith       ierr        = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr);
471b1a0a8a3SJed Brown       ierr        = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
472b1a0a8a3SJed Brown       ierr        = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
4738f3b4b4dSDmitry Karpeev       ierr        = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); /* Why do we need this here? */
4748f3b4b4dSDmitry Karpeev       if (subdomain_dm) {
4758f3b4b4dSDmitry Karpeev 	ierr = KSPSetDM(ksp,subdomain_dm[i]);CHKERRQ(ierr);
4768f3b4b4dSDmitry Karpeev 	ierr = DMDestroy(subdomain_dm+i);CHKERRQ(ierr);
4778f3b4b4dSDmitry Karpeev       }
478b1a0a8a3SJed Brown       ierr        = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
479b1a0a8a3SJed Brown       ierr        = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
4808f3b4b4dSDmitry Karpeev       if (subdomain_names && subdomain_names[i]) {
4818f3b4b4dSDmitry Karpeev 	ierr = PetscSNPrintf(subprefix,PETSC_MAX_PATH_LEN,"sub_%s_",subdomain_names[i]);CHKERRQ(ierr);
4828f3b4b4dSDmitry Karpeev 	ierr = KSPAppendOptionsPrefix(ksp,subprefix);CHKERRQ(ierr);
4838f3b4b4dSDmitry Karpeev 	ierr = PetscFree(subdomain_names[i]);CHKERRQ(ierr);
4848f3b4b4dSDmitry Karpeev       }
485b1a0a8a3SJed Brown       ierr        = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
486b1a0a8a3SJed Brown       osm->ksp[i] = ksp;
487b1a0a8a3SJed Brown     }
4888f3b4b4dSDmitry Karpeev     ierr = PetscFree(subdomain_dm);CHKERRQ(ierr);
4898f3b4b4dSDmitry Karpeev     ierr = PetscFree(subdomain_names);CHKERRQ(ierr);
490b1a0a8a3SJed Brown     scall = MAT_INITIAL_MATRIX;
491b1a0a8a3SJed Brown 
4928f3b4b4dSDmitry Karpeev   } else { /*if (pc->setupcalled)*/
493b1a0a8a3SJed Brown     /*
4948f3b4b4dSDmitry Karpeev        Destroy the submatrices from the previous iteration
495b1a0a8a3SJed Brown     */
496b1a0a8a3SJed Brown     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
497f746d493SDmitry Karpeev       ierr  = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
498b1a0a8a3SJed Brown       scall = MAT_INITIAL_MATRIX;
499b1a0a8a3SJed Brown     }
500b1a0a8a3SJed Brown   }
501b1a0a8a3SJed Brown 
502b1a0a8a3SJed Brown   /*
503f746d493SDmitry Karpeev      Extract out the submatrices.
504b1a0a8a3SJed Brown   */
50581a5c4aaSDmitry Karpeev   if (size > 1) {
5068f3b4b4dSDmitry Karpeev     ierr = MatGetSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr);
5072fa5cd67SKarl Rupp   } else {
5086a4f0f73SDmitry Karpeev     ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr);
50981a5c4aaSDmitry Karpeev   }
510b1a0a8a3SJed Brown   if (scall == MAT_INITIAL_MATRIX) {
511b1a0a8a3SJed Brown     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
512f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
5133bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr);
514b1a0a8a3SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
515b1a0a8a3SJed Brown     }
516b1a0a8a3SJed Brown   }
517b1a0a8a3SJed Brown 
518b1a0a8a3SJed Brown   /* Return control to the user so that the submatrices can be modified (e.g., to apply
519b1a0a8a3SJed Brown      different boundary conditions for the submatrices than for the global problem) */
5206a4f0f73SDmitry Karpeev   ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
521b1a0a8a3SJed Brown 
522b1a0a8a3SJed Brown   /*
5236a4f0f73SDmitry Karpeev      Loop over submatrices putting them into local ksps
524b1a0a8a3SJed Brown   */
525f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
52623ee1639SBarry Smith     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr);
527b1a0a8a3SJed Brown     if (!pc->setupcalled) {
528b1a0a8a3SJed Brown       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
529b1a0a8a3SJed Brown     }
530b1a0a8a3SJed Brown   }
531b1a0a8a3SJed Brown   PetscFunctionReturn(0);
532b1a0a8a3SJed Brown }
533b1a0a8a3SJed Brown 
534b1a0a8a3SJed Brown #undef __FUNCT__
535f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM"
536f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
537b1a0a8a3SJed Brown {
538f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
539b1a0a8a3SJed Brown   PetscErrorCode ierr;
540b1a0a8a3SJed Brown   PetscInt       i;
541b1a0a8a3SJed Brown 
542b1a0a8a3SJed Brown   PetscFunctionBegin;
543f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
544b1a0a8a3SJed Brown     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
545b1a0a8a3SJed Brown   }
546b1a0a8a3SJed Brown   PetscFunctionReturn(0);
547b1a0a8a3SJed Brown }
548b1a0a8a3SJed Brown 
549b1a0a8a3SJed Brown #undef __FUNCT__
550f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM"
551f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y)
552b1a0a8a3SJed Brown {
553f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
554b1a0a8a3SJed Brown   PetscErrorCode ierr;
555f746d493SDmitry Karpeev   PetscInt       i;
556b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
557b1a0a8a3SJed Brown 
558b1a0a8a3SJed Brown   PetscFunctionBegin;
559b1a0a8a3SJed Brown   /*
5606a4f0f73SDmitry Karpeev      Support for limiting the restriction or interpolation only to the inner
561b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
562b1a0a8a3SJed Brown   */
563f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
564b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
565f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
5666a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5672fa5cd67SKarl Rupp   } else {
5686a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
569b1a0a8a3SJed Brown   }
5706a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr);
5716a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
5726a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5732fa5cd67SKarl Rupp   } else {
5746a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5756a4f0f73SDmitry Karpeev   }
57686b96d36SDmitry Karpeev   /* do the subdomain solves */
57786b96d36SDmitry Karpeev   for (i=0; i<osm->n; ++i) {
578b1a0a8a3SJed Brown     ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
579b1a0a8a3SJed Brown   }
5806a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(y);CHKERRQ(ierr);
5816a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
5826a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5836a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);  PetscFunctionReturn(0);
5842fa5cd67SKarl Rupp   } else {
5856a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5866a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);  PetscFunctionReturn(0);
5876a4f0f73SDmitry Karpeev   }
588b1a0a8a3SJed Brown }
589b1a0a8a3SJed Brown 
590b1a0a8a3SJed Brown #undef __FUNCT__
591f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM"
592f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y)
593b1a0a8a3SJed Brown {
594f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
595b1a0a8a3SJed Brown   PetscErrorCode ierr;
596f746d493SDmitry Karpeev   PetscInt       i;
597b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
598b1a0a8a3SJed Brown 
599b1a0a8a3SJed Brown   PetscFunctionBegin;
600b1a0a8a3SJed Brown   /*
601b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
602b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
603b1a0a8a3SJed Brown 
604f746d493SDmitry Karpeev      Note: these are reversed from the PCApply_GASM() because we are applying the
605b1a0a8a3SJed Brown      transpose of the three terms
606b1a0a8a3SJed Brown   */
607f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
608b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
609f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
6106a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
6112fa5cd67SKarl Rupp   } else {
6126a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
613b1a0a8a3SJed Brown   }
6146a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr);
6156a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
6166a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
6172fa5cd67SKarl Rupp   } else {
6186a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
6196a4f0f73SDmitry Karpeev   }
620b1a0a8a3SJed Brown   /* do the local solves */
621ab3e923fSDmitry 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. */
622b1a0a8a3SJed Brown     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
623b1a0a8a3SJed Brown   }
6246a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(y);CHKERRQ(ierr);
6256a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
6266a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6276a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6282fa5cd67SKarl Rupp   } else {
6296a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6306a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6316a4f0f73SDmitry Karpeev   }
632b1a0a8a3SJed Brown   PetscFunctionReturn(0);
633b1a0a8a3SJed Brown }
634b1a0a8a3SJed Brown 
635b1a0a8a3SJed Brown #undef __FUNCT__
636a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM"
637a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc)
638b1a0a8a3SJed Brown {
639f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
640b1a0a8a3SJed Brown   PetscErrorCode ierr;
641b1a0a8a3SJed Brown   PetscInt       i;
642b1a0a8a3SJed Brown 
643b1a0a8a3SJed Brown   PetscFunctionBegin;
644b1a0a8a3SJed Brown   if (osm->ksp) {
645f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
646a06653b4SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
647b1a0a8a3SJed Brown     }
648b1a0a8a3SJed Brown   }
649b1a0a8a3SJed Brown   if (osm->pmat) {
650f746d493SDmitry Karpeev     if (osm->n > 0) {
651ab3e923fSDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
652b1a0a8a3SJed Brown     }
653b1a0a8a3SJed Brown   }
654d34fcf5fSBarry Smith   if (osm->x) {
655f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
656fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
657fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
658b1a0a8a3SJed Brown     }
659d34fcf5fSBarry Smith   }
6606bf464f9SBarry Smith   ierr = VecDestroy(&osm->gx);CHKERRQ(ierr);
6616bf464f9SBarry Smith   ierr = VecDestroy(&osm->gy);CHKERRQ(ierr);
662ab3e923fSDmitry Karpeev 
6636a4f0f73SDmitry Karpeev   ierr     = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr);
6646a4f0f73SDmitry Karpeev   ierr     = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr);
6658f3b4b4dSDmitry Karpeev   if (!osm->user_subdomains) {
6662c112581SDmitry Karpeev     ierr      = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr);
6678f3b4b4dSDmitry Karpeev     osm->N    = PETSC_DETERMINE;
6688f3b4b4dSDmitry Karpeev     osm->nmax = PETSC_DETERMINE;
6698f3b4b4dSDmitry Karpeev   }
670a06653b4SBarry Smith   PetscFunctionReturn(0);
671a06653b4SBarry Smith }
672a06653b4SBarry Smith 
673a06653b4SBarry Smith #undef __FUNCT__
674a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM"
675a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc)
676a06653b4SBarry Smith {
677a06653b4SBarry Smith   PC_GASM        *osm = (PC_GASM*)pc->data;
678a06653b4SBarry Smith   PetscErrorCode ierr;
679a06653b4SBarry Smith   PetscInt       i;
680a06653b4SBarry Smith 
681a06653b4SBarry Smith   PetscFunctionBegin;
682a06653b4SBarry Smith   ierr = PCReset_GASM(pc);CHKERRQ(ierr);
683135757f6SDmitry Karpeev 
684135757f6SDmitry Karpeev   /* PCReset will not destroy subdomains, if user_subdomains is true. */
685135757f6SDmitry Karpeev   ierr = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr);
686135757f6SDmitry Karpeev 
687a06653b4SBarry Smith   if (osm->ksp) {
688a06653b4SBarry Smith     for (i=0; i<osm->n; i++) {
6896bf464f9SBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
690a06653b4SBarry Smith     }
691a06653b4SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
692a06653b4SBarry Smith   }
693a06653b4SBarry Smith   ierr = PetscFree(osm->x);CHKERRQ(ierr);
694a06653b4SBarry Smith   ierr = PetscFree(osm->y);CHKERRQ(ierr);
695c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
696b1a0a8a3SJed Brown   PetscFunctionReturn(0);
697b1a0a8a3SJed Brown }
698b1a0a8a3SJed Brown 
699b1a0a8a3SJed Brown #undef __FUNCT__
700f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM"
7018c34d3f5SBarry Smith static PetscErrorCode PCSetFromOptions_GASM(PetscOptions *PetscOptionsObject,PC pc)
702a6dfd86eSKarl Rupp {
703f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
704b1a0a8a3SJed Brown   PetscErrorCode ierr;
705b1a0a8a3SJed Brown   PetscInt       blocks,ovl;
706b1a0a8a3SJed Brown   PetscBool      symset,flg;
707f746d493SDmitry Karpeev   PCGASMType     gasmtype;
708b1a0a8a3SJed Brown 
709b1a0a8a3SJed Brown   PetscFunctionBegin;
710b1a0a8a3SJed Brown   /* set the type to symmetric if matrix is symmetric */
711b1a0a8a3SJed Brown   if (!osm->type_set && pc->pmat) {
712b1a0a8a3SJed Brown     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
7132fa5cd67SKarl Rupp     if (symset && flg) osm->type = PC_GASM_BASIC;
714b1a0a8a3SJed Brown   }
715e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Generalized additive Schwarz options");CHKERRQ(ierr);
7168f3b4b4dSDmitry 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);
7178f3b4b4dSDmitry Karpeev   ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg);CHKERRQ(ierr);
71865db9045SDmitry Karpeev   if (flg) {
7198f3b4b4dSDmitry Karpeev     ierr = PCGASMSetTotalSubdomains(pc,blocks);CHKERRQ(ierr);
72065db9045SDmitry Karpeev   }
72106b43e7eSDmitry Karpeev   ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
72265db9045SDmitry Karpeev   if (flg) {
72365db9045SDmitry Karpeev     ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr);
724d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
72565db9045SDmitry Karpeev   }
726b1a0a8a3SJed Brown   flg  = PETSC_FALSE;
727f746d493SDmitry Karpeev   ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr);
728f746d493SDmitry Karpeev   if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr);}
729b1a0a8a3SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
730b1a0a8a3SJed Brown   PetscFunctionReturn(0);
731b1a0a8a3SJed Brown }
732b1a0a8a3SJed Brown 
733b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/
734b1a0a8a3SJed Brown 
735b1a0a8a3SJed Brown #undef __FUNCT__
7368f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains"
7378f3b4b4dSDmitry Karpeev /*@
7388f3b4b4dSDmitry Karpeev     PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the
7398f3b4b4dSDmitry Karpeev                                communicator.
7408f3b4b4dSDmitry Karpeev     Logically collective on pc
7418f3b4b4dSDmitry Karpeev 
7428f3b4b4dSDmitry Karpeev     Input Parameters:
7438f3b4b4dSDmitry Karpeev +   pc  - the preconditioner
7448f3b4b4dSDmitry Karpeev -   N   - total number of subdomains
7458f3b4b4dSDmitry Karpeev 
7468f3b4b4dSDmitry Karpeev 
7478f3b4b4dSDmitry Karpeev     Level: beginner
7488f3b4b4dSDmitry Karpeev 
7498f3b4b4dSDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
7508f3b4b4dSDmitry Karpeev 
7518f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap()
7528f3b4b4dSDmitry Karpeev           PCGASMCreateSubdomains2D()
7538f3b4b4dSDmitry Karpeev @*/
7548f3b4b4dSDmitry Karpeev PetscErrorCode  PCGASMSetTotalSubdomains(PC pc,PetscInt N)
7558f3b4b4dSDmitry Karpeev {
7568f3b4b4dSDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
7578f3b4b4dSDmitry Karpeev   PetscMPIInt    size,rank;
7588f3b4b4dSDmitry Karpeev   PetscErrorCode ierr;
7598f3b4b4dSDmitry Karpeev 
7608f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
7618f3b4b4dSDmitry 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);
7628f3b4b4dSDmitry Karpeev   if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before calling PCSetUp().");
7638f3b4b4dSDmitry Karpeev 
7642c112581SDmitry Karpeev   ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr);
7658f3b4b4dSDmitry Karpeev   osm->ois = osm->iis = NULL;
7668f3b4b4dSDmitry Karpeev 
7678f3b4b4dSDmitry Karpeev   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
7688f3b4b4dSDmitry Karpeev   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
7698f3b4b4dSDmitry Karpeev   osm->N    = N;
7708f3b4b4dSDmitry Karpeev   osm->n    = PETSC_DETERMINE;
7718f3b4b4dSDmitry Karpeev   osm->nmax = PETSC_DETERMINE;
7728f3b4b4dSDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
7738f3b4b4dSDmitry Karpeev   PetscFunctionReturn(0);
7748f3b4b4dSDmitry Karpeev }
7758f3b4b4dSDmitry Karpeev 
7768f3b4b4dSDmitry Karpeev #undef __FUNCT__
77706b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM"
778b541e6a4SDmitry Karpeev static PetscErrorCode  PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[])
779b1a0a8a3SJed Brown {
780f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
781b1a0a8a3SJed Brown   PetscErrorCode ierr;
782b1a0a8a3SJed Brown   PetscInt       i;
783b1a0a8a3SJed Brown 
784b1a0a8a3SJed Brown   PetscFunctionBegin;
7858f3b4b4dSDmitry Karpeev   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, got n = %D",n);
7868f3b4b4dSDmitry Karpeev   if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp().");
787b1a0a8a3SJed Brown 
7882c112581SDmitry Karpeev   ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr);
7898f3b4b4dSDmitry Karpeev   osm->iis  = osm->ois = NULL;
7908f3b4b4dSDmitry Karpeev   osm->n    = n;
7918f3b4b4dSDmitry Karpeev   osm->N    = PETSC_DETERMINE;
7928f3b4b4dSDmitry Karpeev   osm->nmax = PETSC_DETERMINE;
793a35b7b57SDmitry Karpeev   if (ois) {
794785e854fSJed Brown     ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr);
7958f3b4b4dSDmitry Karpeev     for (i=0; i<n; i++) {
7968f3b4b4dSDmitry Karpeev       ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);
7978f3b4b4dSDmitry Karpeev       osm->ois[i] = ois[i];
7988f3b4b4dSDmitry Karpeev     }
7998f3b4b4dSDmitry Karpeev     /*
8008f3b4b4dSDmitry Karpeev        Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(),
8018f3b4b4dSDmitry Karpeev        it will be ignored.  To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1.
8028f3b4b4dSDmitry Karpeev     */
803b1a0a8a3SJed Brown     osm->overlap = -1;
804a35b7b57SDmitry Karpeev     if (!iis) {
805785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr);
806a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {
8078f3b4b4dSDmitry Karpeev 	for (i=0; i<n; i++) {
8080e8f3646SDmitry Karpeev 	  ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);
809a35b7b57SDmitry Karpeev 	  osm->iis[i] = ois[i];
810a35b7b57SDmitry Karpeev 	}
811a35b7b57SDmitry Karpeev       }
812a35b7b57SDmitry Karpeev     }
8138f3b4b4dSDmitry Karpeev   }
814a35b7b57SDmitry Karpeev   if (iis) {
815785e854fSJed Brown     ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr);
8168f3b4b4dSDmitry Karpeev     for (i=0; i<n; i++) {
8178f3b4b4dSDmitry Karpeev       ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);
8188f3b4b4dSDmitry Karpeev       osm->iis[i] = iis[i];
8198f3b4b4dSDmitry Karpeev     }
820a35b7b57SDmitry Karpeev     if (!ois) {
821785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr);
822a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {
8238f3b4b4dSDmitry Karpeev 	for (i=0; i<n; i++) {
824a35b7b57SDmitry Karpeev 	  ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);
825a35b7b57SDmitry Karpeev 	  osm->ois[i] = iis[i];
826a35b7b57SDmitry Karpeev 	}
8278f3b4b4dSDmitry Karpeev       }
8288f3b4b4dSDmitry Karpeev       /* If necessary, osm->ois[i] will be expanded using the requested overlap size in PCSetUp_GASM(). */
829a35b7b57SDmitry Karpeev     }
830a35b7b57SDmitry Karpeev   }
8318f3b4b4dSDmitry Karpeev   if (iis)  osm->user_subdomains = PETSC_TRUE;
832b1a0a8a3SJed Brown   PetscFunctionReturn(0);
833b1a0a8a3SJed Brown }
834b1a0a8a3SJed Brown 
835b1a0a8a3SJed Brown 
836b1a0a8a3SJed Brown #undef __FUNCT__
837f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM"
838f7a08781SBarry Smith static PetscErrorCode  PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
839b1a0a8a3SJed Brown {
840f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
841b1a0a8a3SJed Brown 
842b1a0a8a3SJed Brown   PetscFunctionBegin;
843ce94432eSBarry Smith   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
844ce94432eSBarry Smith   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
8452fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
846b1a0a8a3SJed Brown   PetscFunctionReturn(0);
847b1a0a8a3SJed Brown }
848b1a0a8a3SJed Brown 
849b1a0a8a3SJed Brown #undef __FUNCT__
850f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM"
851f7a08781SBarry Smith static PetscErrorCode  PCGASMSetType_GASM(PC pc,PCGASMType type)
852b1a0a8a3SJed Brown {
853f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
854b1a0a8a3SJed Brown 
855b1a0a8a3SJed Brown   PetscFunctionBegin;
856b1a0a8a3SJed Brown   osm->type     = type;
857b1a0a8a3SJed Brown   osm->type_set = PETSC_TRUE;
858b1a0a8a3SJed Brown   PetscFunctionReturn(0);
859b1a0a8a3SJed Brown }
860b1a0a8a3SJed Brown 
861b1a0a8a3SJed Brown #undef __FUNCT__
862f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM"
863f7a08781SBarry Smith static PetscErrorCode  PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort)
864b1a0a8a3SJed Brown {
865f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
866b1a0a8a3SJed Brown 
867b1a0a8a3SJed Brown   PetscFunctionBegin;
868b1a0a8a3SJed Brown   osm->sort_indices = doSort;
869b1a0a8a3SJed Brown   PetscFunctionReturn(0);
870b1a0a8a3SJed Brown }
871b1a0a8a3SJed Brown 
872b1a0a8a3SJed Brown #undef __FUNCT__
873f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM"
87444fe18e5SDmitry Karpeev /*
8758f3b4b4dSDmitry Karpeev    FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed.
87644fe18e5SDmitry Karpeev         In particular, it would upset the global subdomain number calculation.
87744fe18e5SDmitry Karpeev */
878f7a08781SBarry Smith static PetscErrorCode  PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
879b1a0a8a3SJed Brown {
880f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
881b1a0a8a3SJed Brown   PetscErrorCode ierr;
882b1a0a8a3SJed Brown 
883b1a0a8a3SJed Brown   PetscFunctionBegin;
884ce94432eSBarry 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");
885b1a0a8a3SJed Brown 
8862fa5cd67SKarl Rupp   if (n) *n = osm->n;
887ab3e923fSDmitry Karpeev   if (first) {
888ce94432eSBarry Smith     ierr    = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
889ab3e923fSDmitry Karpeev     *first -= osm->n;
890b1a0a8a3SJed Brown   }
891b1a0a8a3SJed Brown   if (ksp) {
892b1a0a8a3SJed Brown     /* Assume that local solves are now different; not necessarily
89306b43e7eSDmitry Karpeev        true, though!  This flag is used only for PCView_GASM() */
894b1a0a8a3SJed Brown     *ksp                        = osm->ksp;
8956a4f0f73SDmitry Karpeev     osm->same_subdomain_solvers = PETSC_FALSE;
896b1a0a8a3SJed Brown   }
897b1a0a8a3SJed Brown   PetscFunctionReturn(0);
898ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */
899b1a0a8a3SJed Brown 
900b1a0a8a3SJed Brown #undef __FUNCT__
90106b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains"
902b1a0a8a3SJed Brown /*@C
90306b43e7eSDmitry Karpeev     PCGASMSetSubdomains - Sets the subdomains for this processor
90406b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
905b1a0a8a3SJed Brown 
9068f3b4b4dSDmitry Karpeev     Collective on pc
907b1a0a8a3SJed Brown 
908b1a0a8a3SJed Brown     Input Parameters:
9098f3b4b4dSDmitry Karpeev +   pc  - the preconditioner object
91006b43e7eSDmitry Karpeev .   n   - the number of subdomains for this processor
9118f3b4b4dSDmitry Karpeev .   iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains)
9128f3b4b4dSDmitry 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)
913b1a0a8a3SJed Brown 
914b1a0a8a3SJed Brown     Notes:
91506b43e7eSDmitry Karpeev     The IS indices use the parallel, global numbering of the vector entries.
9166a4f0f73SDmitry Karpeev     Inner subdomains are those where the correction is applied.
9176a4f0f73SDmitry Karpeev     Outer subdomains are those where the residual necessary to obtain the
9186a4f0f73SDmitry Karpeev     corrections is obtained (see PCGASMType for the use of inner/outer subdomains).
9196a4f0f73SDmitry Karpeev     Both inner and outer subdomains can extend over several processors.
9206a4f0f73SDmitry Karpeev     This processor's portion of a subdomain is known as a local subdomain.
9216a4f0f73SDmitry Karpeev 
9226a4f0f73SDmitry Karpeev     By default the GASM preconditioner uses 1 (local) subdomain per processor.
9236a4f0f73SDmitry Karpeev 
924b1a0a8a3SJed Brown 
925b1a0a8a3SJed Brown     Level: advanced
926b1a0a8a3SJed Brown 
92706b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz
928b1a0a8a3SJed Brown 
9298f3b4b4dSDmitry Karpeev .seealso: PCGASMSetNumSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
93006b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
931b1a0a8a3SJed Brown @*/
9326a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[])
933b1a0a8a3SJed Brown {
9348f3b4b4dSDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
935b1a0a8a3SJed Brown   PetscErrorCode ierr;
936b1a0a8a3SJed Brown 
937b1a0a8a3SJed Brown   PetscFunctionBegin;
938b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9396a4f0f73SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr);
9408f3b4b4dSDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
941b1a0a8a3SJed Brown   PetscFunctionReturn(0);
942b1a0a8a3SJed Brown }
943b1a0a8a3SJed Brown 
944b1a0a8a3SJed Brown 
945b1a0a8a3SJed Brown #undef __FUNCT__
946f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap"
947b1a0a8a3SJed Brown /*@
948f746d493SDmitry Karpeev     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
949b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
9508f3b4b4dSDmitry Karpeev     pc communicator must call this routine.
951b1a0a8a3SJed Brown 
9528f3b4b4dSDmitry Karpeev     Logically Collective on pc
953b1a0a8a3SJed Brown 
954b1a0a8a3SJed Brown     Input Parameters:
955b1a0a8a3SJed Brown +   pc  - the preconditioner context
9568f3b4b4dSDmitry Karpeev -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0)
957b1a0a8a3SJed Brown 
958b1a0a8a3SJed Brown     Options Database Key:
95906b43e7eSDmitry Karpeev .   -pc_gasm_overlap <overlap> - Sets overlap
960b1a0a8a3SJed Brown 
961b1a0a8a3SJed Brown     Notes:
96206b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.  To use
9638f3b4b4dSDmitry Karpeev     multiple subdomain per perocessor or "straddling" subdomains that intersect
9648f3b4b4dSDmitry Karpeev     multiple processors use PCGASMSetSubdomains() (or option -pc_gasm_total_subdomains <n>).
965b1a0a8a3SJed Brown 
9668f3b4b4dSDmitry Karpeev     The overlap defaults to 0, so if one desires that no additional
967b1a0a8a3SJed Brown     overlap be computed beyond what may have been set with a call to
9688f3b4b4dSDmitry Karpeev     PCGASMSetSubdomains(), then ovl must be set to be 0.  In particular, if one does
9698f3b4b4dSDmitry Karpeev     not explicitly set the subdomains in application code, then all overlap would be computed
970f746d493SDmitry Karpeev     internally by PETSc, and using an overlap of 0 would result in an GASM
971b1a0a8a3SJed Brown     variant that is equivalent to the block Jacobi preconditioner.
972b1a0a8a3SJed Brown 
973b1a0a8a3SJed Brown     Note that one can define initial index sets with any overlap via
97406b43e7eSDmitry Karpeev     PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows
97506b43e7eSDmitry Karpeev     PETSc to extend that overlap further, if desired.
976b1a0a8a3SJed Brown 
977b1a0a8a3SJed Brown     Level: intermediate
978b1a0a8a3SJed Brown 
979f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap
980b1a0a8a3SJed Brown 
9818f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
98206b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
983b1a0a8a3SJed Brown @*/
9847087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap(PC pc,PetscInt ovl)
985b1a0a8a3SJed Brown {
986b1a0a8a3SJed Brown   PetscErrorCode ierr;
9878f3b4b4dSDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
988b1a0a8a3SJed Brown 
989b1a0a8a3SJed Brown   PetscFunctionBegin;
990b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
991b1a0a8a3SJed Brown   PetscValidLogicalCollectiveInt(pc,ovl,2);
992f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
9938f3b4b4dSDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
994b1a0a8a3SJed Brown   PetscFunctionReturn(0);
995b1a0a8a3SJed Brown }
996b1a0a8a3SJed Brown 
997b1a0a8a3SJed Brown #undef __FUNCT__
998f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType"
999b1a0a8a3SJed Brown /*@
1000f746d493SDmitry Karpeev     PCGASMSetType - Sets the type of restriction and interpolation used
1001b1a0a8a3SJed Brown     for local problems in the additive Schwarz method.
1002b1a0a8a3SJed Brown 
1003b1a0a8a3SJed Brown     Logically Collective on PC
1004b1a0a8a3SJed Brown 
1005b1a0a8a3SJed Brown     Input Parameters:
1006b1a0a8a3SJed Brown +   pc  - the preconditioner context
1007f746d493SDmitry Karpeev -   type - variant of GASM, one of
1008b1a0a8a3SJed Brown .vb
1009f746d493SDmitry Karpeev       PC_GASM_BASIC       - full interpolation and restriction
1010f746d493SDmitry Karpeev       PC_GASM_RESTRICT    - full restriction, local processor interpolation
1011f746d493SDmitry Karpeev       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
1012f746d493SDmitry Karpeev       PC_GASM_NONE        - local processor restriction and interpolation
1013b1a0a8a3SJed Brown .ve
1014b1a0a8a3SJed Brown 
1015b1a0a8a3SJed Brown     Options Database Key:
1016f746d493SDmitry Karpeev .   -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1017b1a0a8a3SJed Brown 
1018b1a0a8a3SJed Brown     Level: intermediate
1019b1a0a8a3SJed Brown 
1020f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1021b1a0a8a3SJed Brown 
10228f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1023f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1024b1a0a8a3SJed Brown @*/
10257087cfbeSBarry Smith PetscErrorCode  PCGASMSetType(PC pc,PCGASMType type)
1026b1a0a8a3SJed Brown {
1027b1a0a8a3SJed Brown   PetscErrorCode ierr;
1028b1a0a8a3SJed Brown 
1029b1a0a8a3SJed Brown   PetscFunctionBegin;
1030b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1031b1a0a8a3SJed Brown   PetscValidLogicalCollectiveEnum(pc,type,2);
1032f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr);
1033b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1034b1a0a8a3SJed Brown }
1035b1a0a8a3SJed Brown 
1036b1a0a8a3SJed Brown #undef __FUNCT__
1037f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices"
1038b1a0a8a3SJed Brown /*@
1039f746d493SDmitry Karpeev     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
1040b1a0a8a3SJed Brown 
1041b1a0a8a3SJed Brown     Logically Collective on PC
1042b1a0a8a3SJed Brown 
1043b1a0a8a3SJed Brown     Input Parameters:
1044b1a0a8a3SJed Brown +   pc  - the preconditioner context
1045b1a0a8a3SJed Brown -   doSort - sort the subdomain indices
1046b1a0a8a3SJed Brown 
1047b1a0a8a3SJed Brown     Level: intermediate
1048b1a0a8a3SJed Brown 
1049f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1050b1a0a8a3SJed Brown 
10518f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1052f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1053b1a0a8a3SJed Brown @*/
10547087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices(PC pc,PetscBool doSort)
1055b1a0a8a3SJed Brown {
1056b1a0a8a3SJed Brown   PetscErrorCode ierr;
1057b1a0a8a3SJed Brown 
1058b1a0a8a3SJed Brown   PetscFunctionBegin;
1059b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1060b1a0a8a3SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
1061f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
1062b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1063b1a0a8a3SJed Brown }
1064b1a0a8a3SJed Brown 
1065b1a0a8a3SJed Brown #undef __FUNCT__
1066f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP"
1067b1a0a8a3SJed Brown /*@C
1068f746d493SDmitry Karpeev    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
1069b1a0a8a3SJed Brown    this processor.
1070b1a0a8a3SJed Brown 
1071b1a0a8a3SJed Brown    Collective on PC iff first_local is requested
1072b1a0a8a3SJed Brown 
1073b1a0a8a3SJed Brown    Input Parameter:
1074b1a0a8a3SJed Brown .  pc - the preconditioner context
1075b1a0a8a3SJed Brown 
1076b1a0a8a3SJed Brown    Output Parameters:
10770298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
10780298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
10790298fd71SBarry Smith                  all processors must request or all must pass NULL
1080b1a0a8a3SJed Brown -  ksp - the array of KSP contexts
1081b1a0a8a3SJed Brown 
1082b1a0a8a3SJed Brown    Note:
1083f746d493SDmitry Karpeev    After PCGASMGetSubKSP() the array of KSPes is not to be freed
1084b1a0a8a3SJed Brown 
1085b1a0a8a3SJed Brown    Currently for some matrix implementations only 1 block per processor
1086b1a0a8a3SJed Brown    is supported.
1087b1a0a8a3SJed Brown 
1088f746d493SDmitry Karpeev    You must call KSPSetUp() before calling PCGASMGetSubKSP().
1089b1a0a8a3SJed Brown 
1090b1a0a8a3SJed Brown    Level: advanced
1091b1a0a8a3SJed Brown 
1092f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context
1093b1a0a8a3SJed Brown 
10948f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(),
1095f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(),
1096b1a0a8a3SJed Brown @*/
10977087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1098b1a0a8a3SJed Brown {
1099b1a0a8a3SJed Brown   PetscErrorCode ierr;
1100b1a0a8a3SJed Brown 
1101b1a0a8a3SJed Brown   PetscFunctionBegin;
1102b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1103f746d493SDmitry Karpeev   ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
1104b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1105b1a0a8a3SJed Brown }
1106b1a0a8a3SJed Brown 
1107b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/
1108b1a0a8a3SJed Brown /*MC
1109f746d493SDmitry Karpeev    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1110b1a0a8a3SJed Brown            its own KSP object.
1111b1a0a8a3SJed Brown 
1112b1a0a8a3SJed Brown    Options Database Keys:
11138f3b4b4dSDmitry Karpeev +  -pc_gasm_total_subdomains <n>  - Sets total number of local subdomains to be distributed among processors
111406b43e7eSDmitry Karpeev .  -pc_gasm_view_subdomains       - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view
111506b43e7eSDmitry Karpeev .  -pc_gasm_print_subdomains      - activates the printing of subdomain indices in PCSetUp()
111606b43e7eSDmitry Karpeev .  -pc_gasm_overlap <ovl>         - Sets overlap by which to (automatically) extend local subdomains
1117f746d493SDmitry Karpeev -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1118b1a0a8a3SJed Brown 
1119b1a0a8a3SJed Brown      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
1120f746d493SDmitry Karpeev       will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
1121f746d493SDmitry Karpeev       -pc_gasm_type basic to use the standard GASM.
1122b1a0a8a3SJed Brown 
1123b1a0a8a3SJed Brown    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
1124b1a0a8a3SJed Brown      than one processor. Defaults to one block per processor.
1125b1a0a8a3SJed Brown 
1126b1a0a8a3SJed Brown      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1127b1a0a8a3SJed Brown         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
1128b1a0a8a3SJed Brown 
1129f746d493SDmitry Karpeev      To set the options on the solvers separate for each block call PCGASMGetSubKSP()
1130b1a0a8a3SJed Brown          and set the options directly on the resulting KSP object (you can access its PC
1131b1a0a8a3SJed Brown          with KSPGetPC())
1132b1a0a8a3SJed Brown 
1133b1a0a8a3SJed Brown 
1134b1a0a8a3SJed Brown    Level: beginner
1135b1a0a8a3SJed Brown 
1136b1a0a8a3SJed Brown    Concepts: additive Schwarz method
1137b1a0a8a3SJed Brown 
1138b1a0a8a3SJed Brown     References:
1139b1a0a8a3SJed Brown     An additive variant of the Schwarz alternating method for the case of many subregions
1140b1a0a8a3SJed Brown     M Dryja, OB Widlund - Courant Institute, New York University Technical report
1141b1a0a8a3SJed Brown 
1142b1a0a8a3SJed Brown     Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1143b1a0a8a3SJed Brown     Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X.
1144b1a0a8a3SJed Brown 
1145b1a0a8a3SJed Brown .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
114649517cdeSBarry Smith            PCBJACOBI,  PCGASMGetSubKSP(), PCGASMSetSubdomains(),
11478f3b4b4dSDmitry Karpeev            PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType()
1148b1a0a8a3SJed Brown 
1149b1a0a8a3SJed Brown M*/
1150b1a0a8a3SJed Brown 
1151b1a0a8a3SJed Brown #undef __FUNCT__
1152f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM"
11538cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc)
1154b1a0a8a3SJed Brown {
1155b1a0a8a3SJed Brown   PetscErrorCode ierr;
1156f746d493SDmitry Karpeev   PC_GASM        *osm;
1157b1a0a8a3SJed Brown 
1158b1a0a8a3SJed Brown   PetscFunctionBegin;
1159b00a9115SJed Brown   ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr);
11602fa5cd67SKarl Rupp 
11618f3b4b4dSDmitry Karpeev   osm->N                      = PETSC_DETERMINE;
116206b43e7eSDmitry Karpeev   osm->n                      = PETSC_DECIDE;
11638f3b4b4dSDmitry Karpeev   osm->nmax                   = PETSC_DETERMINE;
11648f3b4b4dSDmitry Karpeev   osm->overlap                = 0;
1165b1a0a8a3SJed Brown   osm->ksp                    = 0;
11666a4f0f73SDmitry Karpeev   osm->gorestriction          = 0;
11676a4f0f73SDmitry Karpeev   osm->girestriction          = 0;
1168ab3e923fSDmitry Karpeev   osm->gx                     = 0;
1169ab3e923fSDmitry Karpeev   osm->gy                     = 0;
1170b1a0a8a3SJed Brown   osm->x                      = 0;
1171b1a0a8a3SJed Brown   osm->y                      = 0;
11726a4f0f73SDmitry Karpeev   osm->ois                    = 0;
11736a4f0f73SDmitry Karpeev   osm->iis                    = 0;
1174b1a0a8a3SJed Brown   osm->pmat                   = 0;
1175f746d493SDmitry Karpeev   osm->type                   = PC_GASM_RESTRICT;
11766a4f0f73SDmitry Karpeev   osm->same_subdomain_solvers = PETSC_TRUE;
1177b1a0a8a3SJed Brown   osm->sort_indices           = PETSC_TRUE;
1178d709ea83SDmitry Karpeev   osm->dm_subdomains          = PETSC_FALSE;
1179b1a0a8a3SJed Brown 
1180b1a0a8a3SJed Brown   pc->data                 = (void*)osm;
1181f746d493SDmitry Karpeev   pc->ops->apply           = PCApply_GASM;
1182f746d493SDmitry Karpeev   pc->ops->applytranspose  = PCApplyTranspose_GASM;
1183f746d493SDmitry Karpeev   pc->ops->setup           = PCSetUp_GASM;
1184a06653b4SBarry Smith   pc->ops->reset           = PCReset_GASM;
1185f746d493SDmitry Karpeev   pc->ops->destroy         = PCDestroy_GASM;
1186f746d493SDmitry Karpeev   pc->ops->setfromoptions  = PCSetFromOptions_GASM;
1187f746d493SDmitry Karpeev   pc->ops->setuponblocks   = PCSetUpOnBlocks_GASM;
1188f746d493SDmitry Karpeev   pc->ops->view            = PCView_GASM;
1189b1a0a8a3SJed Brown   pc->ops->applyrichardson = 0;
1190b1a0a8a3SJed Brown 
1191bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM);CHKERRQ(ierr);
1192bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM);CHKERRQ(ierr);
1193bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM);CHKERRQ(ierr);
1194bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM);CHKERRQ(ierr);
1195bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM);CHKERRQ(ierr);
1196b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1197b1a0a8a3SJed Brown }
1198b1a0a8a3SJed Brown 
1199b1a0a8a3SJed Brown 
1200b1a0a8a3SJed Brown #undef __FUNCT__
120106b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains"
12028f3b4b4dSDmitry Karpeev PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[])
1203b1a0a8a3SJed Brown {
1204b1a0a8a3SJed Brown   MatPartitioning mpart;
1205b1a0a8a3SJed Brown   const char      *prefix;
120611bd1e4dSLisandro Dalcin   PetscErrorCode  (*f)(Mat,MatReuse,Mat*);
1207b1a0a8a3SJed Brown   PetscMPIInt     size;
1208b1a0a8a3SJed Brown   PetscInt        i,j,rstart,rend,bs;
120911bd1e4dSLisandro Dalcin   PetscBool       isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
12100298fd71SBarry Smith   Mat             Ad     = NULL, adj;
1211b1a0a8a3SJed Brown   IS              ispart,isnumb,*is;
1212b1a0a8a3SJed Brown   PetscErrorCode  ierr;
1213b1a0a8a3SJed Brown 
1214b1a0a8a3SJed Brown   PetscFunctionBegin;
12158f3b4b4dSDmitry Karpeev   if (nloc < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local subdomains must > 0, got nloc = %D",nloc);
1216b1a0a8a3SJed Brown 
1217b1a0a8a3SJed Brown   /* Get prefix, row distribution, and block size */
1218b1a0a8a3SJed Brown   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1219b1a0a8a3SJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1220b1a0a8a3SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1221b1a0a8a3SJed 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);
1222b1a0a8a3SJed Brown 
1223b1a0a8a3SJed Brown   /* Get diagonal block from matrix if possible */
1224ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
12250005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",&f);CHKERRQ(ierr);
1226b1a0a8a3SJed Brown   if (f) {
122711bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
1228b1a0a8a3SJed Brown   } else if (size == 1) {
122911bd1e4dSLisandro Dalcin     Ad = A;
1230b1a0a8a3SJed Brown   }
1231b1a0a8a3SJed Brown   if (Ad) {
1232251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1233251f4c67SDmitry Karpeev     if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
1234b1a0a8a3SJed Brown   }
12358f3b4b4dSDmitry Karpeev   if (Ad && nloc > 1) {
1236b1a0a8a3SJed Brown     PetscBool  match,done;
1237b1a0a8a3SJed Brown     /* Try to setup a good matrix partitioning if available */
1238b1a0a8a3SJed Brown     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
1239b1a0a8a3SJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
1240b1a0a8a3SJed Brown     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1241251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
1242b1a0a8a3SJed Brown     if (!match) {
1243251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
1244b1a0a8a3SJed Brown     }
1245b1a0a8a3SJed Brown     if (!match) { /* assume a "good" partitioner is available */
12461a83f524SJed Brown       PetscInt       na;
12471a83f524SJed Brown       const PetscInt *ia,*ja;
1248b1a0a8a3SJed Brown       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1249b1a0a8a3SJed Brown       if (done) {
1250b1a0a8a3SJed Brown         /* Build adjacency matrix by hand. Unfortunately a call to
1251b1a0a8a3SJed Brown            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1252b1a0a8a3SJed Brown            remove the block-aij structure and we cannot expect
1253b1a0a8a3SJed Brown            MatPartitioning to split vertices as we need */
12541a83f524SJed Brown         PetscInt       i,j,len,nnz,cnt,*iia=0,*jja=0;
12551a83f524SJed Brown         const PetscInt *row;
1256b1a0a8a3SJed Brown         nnz = 0;
1257b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* count number of nonzeros */
1258b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1259b1a0a8a3SJed Brown           row = ja + ia[i];
1260b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1261b1a0a8a3SJed Brown             if (row[j] == i) { /* don't count diagonal */
1262b1a0a8a3SJed Brown               len--; break;
1263b1a0a8a3SJed Brown             }
1264b1a0a8a3SJed Brown           }
1265b1a0a8a3SJed Brown           nnz += len;
1266b1a0a8a3SJed Brown         }
1267854ce69bSBarry Smith         ierr   = PetscMalloc1(na+1,&iia);CHKERRQ(ierr);
1268854ce69bSBarry Smith         ierr   = PetscMalloc1(nnz,&jja);CHKERRQ(ierr);
1269b1a0a8a3SJed Brown         nnz    = 0;
1270b1a0a8a3SJed Brown         iia[0] = 0;
1271b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* fill adjacency */
1272b1a0a8a3SJed Brown           cnt = 0;
1273b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1274b1a0a8a3SJed Brown           row = ja + ia[i];
1275b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
12762fa5cd67SKarl Rupp             if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */
1277b1a0a8a3SJed Brown           }
1278b1a0a8a3SJed Brown           nnz += cnt;
1279b1a0a8a3SJed Brown           iia[i+1] = nnz;
1280b1a0a8a3SJed Brown         }
1281b1a0a8a3SJed Brown         /* Partitioning of the adjacency matrix */
12820298fd71SBarry Smith         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr);
1283b1a0a8a3SJed Brown         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
12848f3b4b4dSDmitry Karpeev         ierr      = MatPartitioningSetNParts(mpart,nloc);CHKERRQ(ierr);
1285b1a0a8a3SJed Brown         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
1286b1a0a8a3SJed Brown         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
12876bf464f9SBarry Smith         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
1288b1a0a8a3SJed Brown         foundpart = PETSC_TRUE;
1289b1a0a8a3SJed Brown       }
1290b1a0a8a3SJed Brown       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1291b1a0a8a3SJed Brown     }
12926bf464f9SBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1293b1a0a8a3SJed Brown   }
12948f3b4b4dSDmitry Karpeev   ierr = PetscMalloc1(nloc,&is);CHKERRQ(ierr);
1295b1a0a8a3SJed Brown   if (!foundpart) {
1296b1a0a8a3SJed Brown 
1297b1a0a8a3SJed Brown     /* Partitioning by contiguous chunks of rows */
1298b1a0a8a3SJed Brown 
1299b1a0a8a3SJed Brown     PetscInt mbs   = (rend-rstart)/bs;
1300b1a0a8a3SJed Brown     PetscInt start = rstart;
13018f3b4b4dSDmitry Karpeev     for (i=0; i<nloc; i++) {
13028f3b4b4dSDmitry Karpeev       PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs;
1303b1a0a8a3SJed Brown       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
1304b1a0a8a3SJed Brown       start += count;
1305b1a0a8a3SJed Brown     }
1306b1a0a8a3SJed Brown 
1307b1a0a8a3SJed Brown   } else {
1308b1a0a8a3SJed Brown 
1309b1a0a8a3SJed Brown     /* Partitioning by adjacency of diagonal block  */
1310b1a0a8a3SJed Brown 
1311b1a0a8a3SJed Brown     const PetscInt *numbering;
1312b1a0a8a3SJed Brown     PetscInt       *count,nidx,*indices,*newidx,start=0;
1313b1a0a8a3SJed Brown     /* Get node count in each partition */
13148f3b4b4dSDmitry Karpeev     ierr = PetscMalloc1(nloc,&count);CHKERRQ(ierr);
13158f3b4b4dSDmitry Karpeev     ierr = ISPartitioningCount(ispart,nloc,count);CHKERRQ(ierr);
1316b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
13178f3b4b4dSDmitry Karpeev       for (i=0; i<nloc; i++) count[i] *= bs;
1318b1a0a8a3SJed Brown     }
1319b1a0a8a3SJed Brown     /* Build indices from node numbering */
1320b1a0a8a3SJed Brown     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1321785e854fSJed Brown     ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr);
1322b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1323b1a0a8a3SJed Brown     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
1324b1a0a8a3SJed Brown     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
1325b1a0a8a3SJed Brown     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
1326b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1327785e854fSJed Brown       ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr);
13282fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
13292fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
13302fa5cd67SKarl Rupp       }
1331b1a0a8a3SJed Brown       ierr    = PetscFree(indices);CHKERRQ(ierr);
1332b1a0a8a3SJed Brown       nidx   *= bs;
1333b1a0a8a3SJed Brown       indices = newidx;
1334b1a0a8a3SJed Brown     }
1335b1a0a8a3SJed Brown     /* Shift to get global indices */
1336b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] += rstart;
1337b1a0a8a3SJed Brown 
1338b1a0a8a3SJed Brown     /* Build the index sets for each block */
13398f3b4b4dSDmitry Karpeev     for (i=0; i<nloc; i++) {
1340b1a0a8a3SJed Brown       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
1341b1a0a8a3SJed Brown       ierr   = ISSort(is[i]);CHKERRQ(ierr);
1342b1a0a8a3SJed Brown       start += count[i];
1343b1a0a8a3SJed Brown     }
1344b1a0a8a3SJed Brown 
1345302440fdSBarry Smith     ierr = PetscFree(count);CHKERRQ(ierr);
1346302440fdSBarry Smith     ierr = PetscFree(indices);CHKERRQ(ierr);
1347fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1348fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
1349b1a0a8a3SJed Brown   }
13506a4f0f73SDmitry Karpeev   *iis = is;
13518f3b4b4dSDmitry Karpeev   PetscFunctionReturn(0);
13528f3b4b4dSDmitry Karpeev }
13538f3b4b4dSDmitry Karpeev 
13548f3b4b4dSDmitry Karpeev #undef __FUNCT__
13558f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMCreateStraddlingSubdomains"
1356b541e6a4SDmitry Karpeev PETSC_INTERN PetscErrorCode  PCGASMCreateStraddlingSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[])
13578f3b4b4dSDmitry Karpeev {
13588f3b4b4dSDmitry Karpeev   PetscErrorCode  ierr;
13598f3b4b4dSDmitry Karpeev 
13608f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
13618f3b4b4dSDmitry Karpeev   ierr = MatSubdomainsCreateCoalesce(A,N,n,iis);CHKERRQ(ierr);
13628f3b4b4dSDmitry Karpeev   PetscFunctionReturn(0);
13638f3b4b4dSDmitry Karpeev }
13648f3b4b4dSDmitry Karpeev 
13658f3b4b4dSDmitry Karpeev 
13668f3b4b4dSDmitry Karpeev 
13678f3b4b4dSDmitry Karpeev #undef __FUNCT__
13688f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains"
13698f3b4b4dSDmitry Karpeev /*@C
13708f3b4b4dSDmitry Karpeev    PCGASMCreateSubdomains - Creates n index sets defining n nonoverlapping subdomains for the additive
13718f3b4b4dSDmitry Karpeev    Schwarz preconditioner for a any problem based on its matrix.
13728f3b4b4dSDmitry Karpeev 
13738f3b4b4dSDmitry Karpeev    Collective
13748f3b4b4dSDmitry Karpeev 
13758f3b4b4dSDmitry Karpeev    Input Parameters:
13768f3b4b4dSDmitry Karpeev +  A       - The global matrix operator
13778f3b4b4dSDmitry Karpeev -  N       - the number of global subdomains requested
13788f3b4b4dSDmitry Karpeev 
13798f3b4b4dSDmitry Karpeev    Output Parameters:
13808f3b4b4dSDmitry Karpeev +  n   - the number of subdomains created on this processor
13818f3b4b4dSDmitry Karpeev -  iis - the array of index sets defining the local inner subdomains (on which the correction is applied)
13828f3b4b4dSDmitry Karpeev 
13838f3b4b4dSDmitry Karpeev    Level: advanced
13848f3b4b4dSDmitry Karpeev 
13858f3b4b4dSDmitry Karpeev    Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor.
13868f3b4b4dSDmitry Karpeev          When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local.
13878f3b4b4dSDmitry Karpeev          The resulting subdomains can be use in PCGASMSetSubdomains(pc,n,iss,NULL).  The overlapping
13888f3b4b4dSDmitry Karpeev 	 outer subdomains will be automatically generated from these according to the requested amount of
13898f3b4b4dSDmitry Karpeev 	 overlap; this is currently supported only with local subdomains.
13908f3b4b4dSDmitry Karpeev 
13918f3b4b4dSDmitry Karpeev 
13928f3b4b4dSDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
13938f3b4b4dSDmitry Karpeev 
13948f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains()
13958f3b4b4dSDmitry Karpeev @*/
1396b541e6a4SDmitry Karpeev PetscErrorCode  PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[])
13978f3b4b4dSDmitry Karpeev {
13988f3b4b4dSDmitry Karpeev   PetscMPIInt     size;
13998f3b4b4dSDmitry Karpeev   PetscErrorCode  ierr;
14008f3b4b4dSDmitry Karpeev 
14018f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
14028f3b4b4dSDmitry Karpeev   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
14038f3b4b4dSDmitry Karpeev   PetscValidPointer(iis,4);
14048f3b4b4dSDmitry Karpeev 
14058f3b4b4dSDmitry Karpeev   if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %D",N);
14068f3b4b4dSDmitry Karpeev   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
14078f3b4b4dSDmitry Karpeev   if (N >= size) {
14088f3b4b4dSDmitry Karpeev     *n = N/size + (N%size);
14098f3b4b4dSDmitry Karpeev     ierr = PCGASMCreateLocalSubdomains(A,*n,iis);CHKERRQ(ierr);
14106a4f0f73SDmitry Karpeev   } else {
14118f3b4b4dSDmitry Karpeev     ierr = PCGASMCreateStraddlingSubdomains(A,N,n,iis);CHKERRQ(ierr);
14126a4f0f73SDmitry Karpeev   }
1413b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1414b1a0a8a3SJed Brown }
1415b1a0a8a3SJed Brown 
1416b1a0a8a3SJed Brown #undef __FUNCT__
1417f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains"
1418b1a0a8a3SJed Brown /*@C
1419f746d493SDmitry Karpeev    PCGASMDestroySubdomains - Destroys the index sets created with
14208f3b4b4dSDmitry Karpeev    PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be
142106b43e7eSDmitry Karpeev    called after setting subdomains with PCGASMSetSubdomains().
1422b1a0a8a3SJed Brown 
1423b1a0a8a3SJed Brown    Collective
1424b1a0a8a3SJed Brown 
1425b1a0a8a3SJed Brown    Input Parameters:
1426b1a0a8a3SJed Brown +  n   - the number of index sets
14276a4f0f73SDmitry Karpeev .  iis - the array of inner subdomains,
14280298fd71SBarry Smith -  ois - the array of outer subdomains, can be NULL
1429b1a0a8a3SJed Brown 
14306a4f0f73SDmitry Karpeev    Level: intermediate
14316a4f0f73SDmitry Karpeev 
14326a4f0f73SDmitry Karpeev    Notes: this is merely a convenience subroutine that walks each list,
14332c112581SDmitry Karpeev    destroys each IS on the list, and then frees the list. At the end the
14342c112581SDmitry Karpeev    list pointers are set to NULL.
1435b1a0a8a3SJed Brown 
1436f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1437b1a0a8a3SJed Brown 
14388f3b4b4dSDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetSubdomains()
1439b1a0a8a3SJed Brown @*/
14402c112581SDmitry Karpeev PetscErrorCode  PCGASMDestroySubdomains(PetscInt n,IS **iis,IS **ois)
1441b1a0a8a3SJed Brown {
1442b1a0a8a3SJed Brown   PetscInt       i;
1443b1a0a8a3SJed Brown   PetscErrorCode ierr;
14445fd66863SKarl Rupp 
1445b1a0a8a3SJed Brown   PetscFunctionBegin;
1446a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1447a35b7b57SDmitry Karpeev   if (iis) {
14486a4f0f73SDmitry Karpeev     PetscValidPointer(iis,2);
14492c112581SDmitry Karpeev     if (*iis) {
14502c112581SDmitry Karpeev       PetscValidPointer(*iis,2);
1451a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {
14522c112581SDmitry Karpeev         ierr = ISDestroy(&(*iis)[i]);CHKERRQ(ierr);
1453a35b7b57SDmitry Karpeev       }
14542c112581SDmitry Karpeev       ierr = PetscFree((*iis));CHKERRQ(ierr);
14552c112581SDmitry Karpeev     }
1456a35b7b57SDmitry Karpeev   }
14576a4f0f73SDmitry Karpeev   if (ois) {
14582c112581SDmitry Karpeev     PetscValidPointer(ois,3);
14592c112581SDmitry Karpeev     if (*ois) {
14602c112581SDmitry Karpeev       PetscValidPointer(*ois,3);
1461a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {
14622c112581SDmitry Karpeev         ierr = ISDestroy(&(*ois)[i]);CHKERRQ(ierr);
1463a35b7b57SDmitry Karpeev       }
1464135757f6SDmitry Karpeev       ierr = PetscFree((*ois));CHKERRQ(ierr);
14652c112581SDmitry Karpeev     }
1466b1a0a8a3SJed Brown   }
1467b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1468b1a0a8a3SJed Brown }
1469b1a0a8a3SJed Brown 
1470af538404SDmitry Karpeev 
1471af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \
1472af538404SDmitry Karpeev   {                                                                                                       \
1473af538404SDmitry Karpeev     PetscInt first_row = first/M, last_row = last/M+1;                                                     \
1474af538404SDmitry Karpeev     /*                                                                                                    \
1475af538404SDmitry Karpeev      Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners       \
1476af538404SDmitry Karpeev      of the bounding box of the intersection of the subdomain with the local ownership range (local       \
1477af538404SDmitry Karpeev      subdomain).                                                                                          \
1478af538404SDmitry Karpeev      Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows       \
1479af538404SDmitry Karpeev      of the intersection.                                                                                 \
1480af538404SDmitry Karpeev     */                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                             \
1481af538404SDmitry Karpeev     /* ylow_loc is the grid row containing the first element of the local sumbdomain */                   \
1482eec7e3faSDmitry Karpeev     *ylow_loc = PetscMax(first_row,ylow);                                                                    \
1483af538404SDmitry Karpeev     /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1484eec7e3faSDmitry Karpeev     *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft;                                                                            \
1485af538404SDmitry Karpeev     /* yhigh_loc is the grid row above the last local subdomain element */                                                                    \
1486eec7e3faSDmitry Karpeev     *yhigh_loc = PetscMin(last_row,yhigh);                                                                                                     \
1487af538404SDmitry Karpeev     /* xright is the offset of the end of the  local subdomain within its grid row (might actually be outside the local subdomain) */         \
1488eec7e3faSDmitry Karpeev     *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright;                                                                          \
1489af538404SDmitry Karpeev     /* Now compute the size of the local subdomain n. */ \
1490c3518bceSDmitry Karpeev     *n = 0;                                               \
1491eec7e3faSDmitry Karpeev     if (*ylow_loc < *yhigh_loc) {                           \
1492af538404SDmitry Karpeev       PetscInt width = xright-xleft;                     \
1493eec7e3faSDmitry Karpeev       *n += width*(*yhigh_loc-*ylow_loc-1);                 \
1494eec7e3faSDmitry Karpeev       *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \
1495eec7e3faSDmitry Karpeev       *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \
1496af538404SDmitry Karpeev     } \
1497af538404SDmitry Karpeev   }
1498af538404SDmitry Karpeev 
1499c3518bceSDmitry Karpeev 
1500eec7e3faSDmitry Karpeev 
1501eec7e3faSDmitry Karpeev #undef __FUNCT__
1502f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D"
1503b1a0a8a3SJed Brown /*@
1504f746d493SDmitry Karpeev    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1505b1a0a8a3SJed Brown    preconditioner for a two-dimensional problem on a regular grid.
1506b1a0a8a3SJed Brown 
1507af538404SDmitry Karpeev    Collective
1508b1a0a8a3SJed Brown 
1509b1a0a8a3SJed Brown    Input Parameters:
151006b43e7eSDmitry Karpeev +  M, N               - the global number of grid points in the x and y directions
1511af538404SDmitry Karpeev .  Mdomains, Ndomains - the global number of subdomains in the x and y directions
1512b1a0a8a3SJed Brown .  dof                - degrees of freedom per node
1513b1a0a8a3SJed Brown -  overlap            - overlap in mesh lines
1514b1a0a8a3SJed Brown 
1515b1a0a8a3SJed Brown    Output Parameters:
1516af538404SDmitry Karpeev +  Nsub - the number of local subdomains created
15176a4f0f73SDmitry Karpeev .  iis  - array of index sets defining inner (nonoverlapping) subdomains
15186a4f0f73SDmitry Karpeev -  ois  - array of index sets defining outer (overlapping, if overlap > 0) subdomains
1519b1a0a8a3SJed Brown 
1520b1a0a8a3SJed Brown 
1521b1a0a8a3SJed Brown    Level: advanced
1522b1a0a8a3SJed Brown 
1523f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid
1524b1a0a8a3SJed Brown 
15258f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), PCGASMSetOverlap()
1526b1a0a8a3SJed Brown @*/
15276a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMCreateSubdomains2D(PC pc,PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap,PetscInt *nsub,IS **iis,IS **ois)
1528b1a0a8a3SJed Brown {
1529b1a0a8a3SJed Brown   PetscErrorCode ierr;
15302bbb417fSDmitry Karpeev   PetscMPIInt    size, rank;
15312bbb417fSDmitry Karpeev   PetscInt       i, j;
15322bbb417fSDmitry Karpeev   PetscInt       maxheight, maxwidth;
15332bbb417fSDmitry Karpeev   PetscInt       xstart, xleft, xright, xleft_loc, xright_loc;
15342bbb417fSDmitry Karpeev   PetscInt       ystart, ylow,  yhigh,  ylow_loc,  yhigh_loc;
1535eec7e3faSDmitry Karpeev   PetscInt       x[2][2], y[2][2], n[2];
15362bbb417fSDmitry Karpeev   PetscInt       first, last;
15372bbb417fSDmitry Karpeev   PetscInt       nidx, *idx;
15382bbb417fSDmitry Karpeev   PetscInt       ii,jj,s,q,d;
1539af538404SDmitry Karpeev   PetscInt       k,kk;
15402bbb417fSDmitry Karpeev   PetscMPIInt    color;
15412bbb417fSDmitry Karpeev   MPI_Comm       comm, subcomm;
15426a4f0f73SDmitry Karpeev   IS             **xis = 0, **is = ois, **is_local = iis;
1543d34fcf5fSBarry Smith 
1544b1a0a8a3SJed Brown   PetscFunctionBegin;
15452bbb417fSDmitry Karpeev   ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr);
15462bbb417fSDmitry Karpeev   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
15472bbb417fSDmitry Karpeev   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
15482bbb417fSDmitry Karpeev   ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr);
1549d34fcf5fSBarry 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) "
15502bbb417fSDmitry Karpeev                                       "does not respect the number of degrees of freedom per grid point %D", first, last, dof);
1551d34fcf5fSBarry Smith 
1552af538404SDmitry Karpeev   /* Determine the number of domains with nonzero intersections with the local ownership range. */
15532bbb417fSDmitry Karpeev   s      = 0;
1554af538404SDmitry Karpeev   ystart = 0;
1555af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1556af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1557af538404SDmitry 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);
1558eec7e3faSDmitry Karpeev     /* Vertical domain limits with an overlap. */
1559eec7e3faSDmitry Karpeev     ylow   = PetscMax(ystart - overlap,0);
1560eec7e3faSDmitry Karpeev     yhigh  = PetscMin(ystart + maxheight + overlap,N);
1561b1a0a8a3SJed Brown     xstart = 0;
1562af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1563af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1564af538404SDmitry 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);
1565eec7e3faSDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1566eec7e3faSDmitry Karpeev       xleft  = PetscMax(xstart - overlap,0);
1567eec7e3faSDmitry Karpeev       xright = PetscMin(xstart + maxwidth + overlap,M);
1568af538404SDmitry Karpeev       /*
1569af538404SDmitry Karpeev          Determine whether this subdomain intersects this processor's ownership range of pc->pmat.
1570af538404SDmitry Karpeev       */
1571c3518bceSDmitry Karpeev       PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
15722fa5cd67SKarl Rupp       if (nidx) ++s;
1573af538404SDmitry Karpeev       xstart += maxwidth;
1574af538404SDmitry Karpeev     } /* for (i = 0; i < Mdomains; ++i) */
1575af538404SDmitry Karpeev     ystart += maxheight;
1576af538404SDmitry Karpeev   } /* for (j = 0; j < Ndomains; ++j) */
15772fa5cd67SKarl Rupp 
1578af538404SDmitry Karpeev   /* Now we can allocate the necessary number of ISs. */
1579af538404SDmitry Karpeev   *nsub  = s;
1580854ce69bSBarry Smith   ierr   = PetscMalloc1(*nsub,is);CHKERRQ(ierr);
1581854ce69bSBarry Smith   ierr   = PetscMalloc1(*nsub,is_local);CHKERRQ(ierr);
1582af538404SDmitry Karpeev   s      = 0;
1583af538404SDmitry Karpeev   ystart = 0;
1584af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1585af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1586af538404SDmitry 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);
1587af538404SDmitry Karpeev     /* Vertical domain limits with an overlap. */
1588eec7e3faSDmitry Karpeev     y[0][0] = PetscMax(ystart - overlap,0);
1589eec7e3faSDmitry Karpeev     y[0][1] = PetscMin(ystart + maxheight + overlap,N);
1590eec7e3faSDmitry Karpeev     /* Vertical domain limits without an overlap. */
1591eec7e3faSDmitry Karpeev     y[1][0] = ystart;
1592eec7e3faSDmitry Karpeev     y[1][1] = PetscMin(ystart + maxheight,N);
1593eec7e3faSDmitry Karpeev     xstart  = 0;
1594af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1595af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1596af538404SDmitry 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);
1597af538404SDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1598eec7e3faSDmitry Karpeev       x[0][0] = PetscMax(xstart - overlap,0);
1599eec7e3faSDmitry Karpeev       x[0][1] = PetscMin(xstart + maxwidth + overlap,M);
1600eec7e3faSDmitry Karpeev       /* Horizontal domain limits without an overlap. */
1601eec7e3faSDmitry Karpeev       x[1][0] = xstart;
1602eec7e3faSDmitry Karpeev       x[1][1] = PetscMin(xstart+maxwidth,M);
16032bbb417fSDmitry Karpeev       /*
16042bbb417fSDmitry Karpeev          Determine whether this domain intersects this processor's ownership range of pc->pmat.
16052bbb417fSDmitry Karpeev          Do this twice: first for the domains with overlaps, and once without.
16062bbb417fSDmitry Karpeev          During the first pass create the subcommunicators, and use them on the second pass as well.
16072bbb417fSDmitry Karpeev       */
16082bbb417fSDmitry Karpeev       for (q = 0; q < 2; ++q) {
16092bbb417fSDmitry Karpeev         /*
16102bbb417fSDmitry Karpeev           domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
16112bbb417fSDmitry Karpeev           according to whether the domain with an overlap or without is considered.
16122bbb417fSDmitry Karpeev         */
16132bbb417fSDmitry Karpeev         xleft = x[q][0]; xright = x[q][1];
16142bbb417fSDmitry Karpeev         ylow  = y[q][0]; yhigh  = y[q][1];
1615c3518bceSDmitry Karpeev         PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1616af538404SDmitry Karpeev         nidx *= dof;
1617eec7e3faSDmitry Karpeev         n[q]  = nidx;
16182bbb417fSDmitry Karpeev         /*
1619eec7e3faSDmitry Karpeev          Based on the counted number of indices in the local domain *with an overlap*,
1620af538404SDmitry Karpeev          construct a subcommunicator of all the processors supporting this domain.
1621eec7e3faSDmitry Karpeev          Observe that a domain with an overlap might have nontrivial local support,
1622eec7e3faSDmitry Karpeev          while the domain without an overlap might not.  Hence, the decision to participate
1623eec7e3faSDmitry Karpeev          in the subcommunicator must be based on the domain with an overlap.
16242bbb417fSDmitry Karpeev          */
16252bbb417fSDmitry Karpeev         if (q == 0) {
16262fa5cd67SKarl Rupp           if (nidx) color = 1;
16272fa5cd67SKarl Rupp           else color = MPI_UNDEFINED;
16282fa5cd67SKarl Rupp 
16292bbb417fSDmitry Karpeev           ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr);
16302bbb417fSDmitry Karpeev         }
1631af538404SDmitry Karpeev         /*
1632eec7e3faSDmitry Karpeev          Proceed only if the number of local indices *with an overlap* is nonzero.
1633af538404SDmitry Karpeev          */
1634eec7e3faSDmitry Karpeev         if (n[0]) {
16352fa5cd67SKarl Rupp           if (q == 0) xis = is;
1636af538404SDmitry Karpeev           if (q == 1) {
1637af538404SDmitry Karpeev             /*
1638eec7e3faSDmitry Karpeev              The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1639af538404SDmitry Karpeev              Moreover, if the overlap is zero, the two ISs are identical.
1640af538404SDmitry Karpeev              */
1641b1a0a8a3SJed Brown             if (overlap == 0) {
1642eec7e3faSDmitry Karpeev               (*is_local)[s] = (*is)[s];
16432bbb417fSDmitry Karpeev               ierr           = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr);
16442bbb417fSDmitry Karpeev               continue;
1645d34fcf5fSBarry Smith             } else {
16466a4f0f73SDmitry Karpeev               xis     = is_local;
1647eec7e3faSDmitry Karpeev               subcomm = ((PetscObject)(*is)[s])->comm;
16482bbb417fSDmitry Karpeev             }
1649af538404SDmitry Karpeev           } /* if (q == 1) */
16500298fd71SBarry Smith           idx  = NULL;
1651785e854fSJed Brown           ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr);
1652eec7e3faSDmitry Karpeev           if (nidx) {
16532bbb417fSDmitry Karpeev             k = 0;
16542bbb417fSDmitry Karpeev             for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1655af538404SDmitry Karpeev               PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft;
1656af538404SDmitry Karpeev               PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright;
1657af538404SDmitry Karpeev               kk = dof*(M*jj + x0);
16582bbb417fSDmitry Karpeev               for (ii=x0; ii<x1; ++ii) {
16592bbb417fSDmitry Karpeev                 for (d = 0; d < dof; ++d) {
16602bbb417fSDmitry Karpeev                   idx[k++] = kk++;
1661b1a0a8a3SJed Brown                 }
1662b1a0a8a3SJed Brown               }
1663b1a0a8a3SJed Brown             }
1664eec7e3faSDmitry Karpeev           }
16656a4f0f73SDmitry Karpeev           ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr);
1666eec7e3faSDmitry Karpeev         }/* if (n[0]) */
16672bbb417fSDmitry Karpeev       }/* for (q = 0; q < 2; ++q) */
16682fa5cd67SKarl Rupp       if (n[0]) ++s;
1669af538404SDmitry Karpeev       xstart += maxwidth;
1670af538404SDmitry Karpeev     } /* for (i = 0; i < Mdomains; ++i) */
16712bbb417fSDmitry Karpeev     ystart += maxheight;
1672af538404SDmitry Karpeev   } /* for (j = 0; j < Ndomains; ++j) */
1673b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1674b1a0a8a3SJed Brown }
1675b1a0a8a3SJed Brown 
1676b1a0a8a3SJed Brown #undef __FUNCT__
167706b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains"
1678b1a0a8a3SJed Brown /*@C
167906b43e7eSDmitry Karpeev     PCGASMGetSubdomains - Gets the subdomains supported on this processor
168006b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
1681b1a0a8a3SJed Brown 
1682b1a0a8a3SJed Brown     Not Collective
1683b1a0a8a3SJed Brown 
1684b1a0a8a3SJed Brown     Input Parameter:
1685b1a0a8a3SJed Brown .   pc - the preconditioner context
1686b1a0a8a3SJed Brown 
1687b1a0a8a3SJed Brown     Output Parameters:
1688b1a0a8a3SJed Brown +   n   - the number of subdomains for this processor (default value = 1)
16890298fd71SBarry Smith .   iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL)
16900298fd71SBarry Smith -   ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL)
1691b1a0a8a3SJed Brown 
1692b1a0a8a3SJed Brown 
1693b1a0a8a3SJed Brown     Notes:
16946a4f0f73SDmitry Karpeev     The user is responsible for destroying the ISs and freeing the returned arrays.
1695b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
1696b1a0a8a3SJed Brown 
1697b1a0a8a3SJed Brown     Level: advanced
1698b1a0a8a3SJed Brown 
169906b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz
1700b1a0a8a3SJed Brown 
17018f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), PCGASMCreateSubdomains2D(),
17028f3b4b4dSDmitry Karpeev           PCGASMSetSubdomains(), PCGASMGetSubmatrices()
1703b1a0a8a3SJed Brown @*/
17046a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[])
1705b1a0a8a3SJed Brown {
1706f746d493SDmitry Karpeev   PC_GASM        *osm;
1707b1a0a8a3SJed Brown   PetscErrorCode ierr;
1708b1a0a8a3SJed Brown   PetscBool      match;
17096a4f0f73SDmitry Karpeev   PetscInt       i;
17105fd66863SKarl Rupp 
1711b1a0a8a3SJed Brown   PetscFunctionBegin;
1712b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1713251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1714ce94432eSBarry Smith   if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1715f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1716ab3e923fSDmitry Karpeev   if (n) *n = osm->n;
17176a4f0f73SDmitry Karpeev   if (iis) {
1718785e854fSJed Brown     ierr = PetscMalloc1(osm->n, iis);CHKERRQ(ierr);
17196a4f0f73SDmitry Karpeev   }
17206a4f0f73SDmitry Karpeev   if (ois) {
1721785e854fSJed Brown     ierr = PetscMalloc1(osm->n, ois);CHKERRQ(ierr);
17226a4f0f73SDmitry Karpeev   }
17236a4f0f73SDmitry Karpeev   if (iis || ois) {
17246a4f0f73SDmitry Karpeev     for (i = 0; i < osm->n; ++i) {
17256a4f0f73SDmitry Karpeev       if (iis) (*iis)[i] = osm->iis[i];
17266a4f0f73SDmitry Karpeev       if (ois) (*ois)[i] = osm->ois[i];
17276a4f0f73SDmitry Karpeev     }
1728b1a0a8a3SJed Brown   }
1729b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1730b1a0a8a3SJed Brown }
1731b1a0a8a3SJed Brown 
1732b1a0a8a3SJed Brown #undef __FUNCT__
173306b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices"
1734b1a0a8a3SJed Brown /*@C
173506b43e7eSDmitry Karpeev     PCGASMGetSubmatrices - Gets the local submatrices (for this processor
1736b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1737b1a0a8a3SJed Brown 
1738b1a0a8a3SJed Brown     Not Collective
1739b1a0a8a3SJed Brown 
1740b1a0a8a3SJed Brown     Input Parameter:
1741b1a0a8a3SJed Brown .   pc - the preconditioner context
1742b1a0a8a3SJed Brown 
1743b1a0a8a3SJed Brown     Output Parameters:
1744b1a0a8a3SJed Brown +   n   - the number of matrices for this processor (default value = 1)
1745b1a0a8a3SJed Brown -   mat - the matrices
1746b1a0a8a3SJed Brown 
174706b43e7eSDmitry Karpeev     Notes: matrices returned by this routine have the same communicators as the index sets (IS)
17488f3b4b4dSDmitry Karpeev            used to define subdomains in PCGASMSetSubdomains()
1749b1a0a8a3SJed Brown     Level: advanced
1750b1a0a8a3SJed Brown 
1751f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi
1752b1a0a8a3SJed Brown 
17538f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(),
175406b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains()
1755b1a0a8a3SJed Brown @*/
175606b43e7eSDmitry Karpeev PetscErrorCode  PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1757b1a0a8a3SJed Brown {
1758f746d493SDmitry Karpeev   PC_GASM        *osm;
1759b1a0a8a3SJed Brown   PetscErrorCode ierr;
1760b1a0a8a3SJed Brown   PetscBool      match;
1761b1a0a8a3SJed Brown 
1762b1a0a8a3SJed Brown   PetscFunctionBegin;
1763b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1764b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1765b1a0a8a3SJed Brown   if (mat) PetscValidPointer(mat,3);
1766ce94432eSBarry Smith   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1767251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1768ce94432eSBarry Smith   if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1769f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1770ab3e923fSDmitry Karpeev   if (n) *n = osm->n;
1771b1a0a8a3SJed Brown   if (mat) *mat = osm->pmat;
1772b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1773b1a0a8a3SJed Brown }
1774d709ea83SDmitry Karpeev 
1775d709ea83SDmitry Karpeev #undef __FUNCT__
17768f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMSetUseDMSubdomains"
1777d709ea83SDmitry Karpeev /*@
17788f3b4b4dSDmitry Karpeev     PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1779d709ea83SDmitry Karpeev     Logically Collective
1780d709ea83SDmitry Karpeev 
1781d709ea83SDmitry Karpeev     Input Parameter:
1782d709ea83SDmitry Karpeev +   pc  - the preconditioner
1783d709ea83SDmitry Karpeev -   flg - boolean indicating whether to use subdomains defined by the DM
1784d709ea83SDmitry Karpeev 
1785d709ea83SDmitry Karpeev     Options Database Key:
17868f3b4b4dSDmitry Karpeev .   -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains
1787d709ea83SDmitry Karpeev 
1788d709ea83SDmitry Karpeev     Level: intermediate
1789d709ea83SDmitry Karpeev 
1790d709ea83SDmitry Karpeev     Notes:
17918f3b4b4dSDmitry Karpeev     PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(),
17928f3b4b4dSDmitry Karpeev     so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap()
17938f3b4b4dSDmitry Karpeev     automatically turns the latter off.
1794d709ea83SDmitry Karpeev 
1795d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1796d709ea83SDmitry Karpeev 
17978f3b4b4dSDmitry Karpeev .seealso: PCGASMGetUseDMSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap()
1798d709ea83SDmitry Karpeev           PCGASMCreateSubdomains2D()
1799d709ea83SDmitry Karpeev @*/
18008f3b4b4dSDmitry Karpeev PetscErrorCode  PCGASMSetUseDMSubdomains(PC pc,PetscBool flg)
1801d709ea83SDmitry Karpeev {
1802d709ea83SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
1803d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1804d709ea83SDmitry Karpeev   PetscBool      match;
1805d709ea83SDmitry Karpeev 
1806d709ea83SDmitry Karpeev   PetscFunctionBegin;
1807d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1808d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
1809d709ea83SDmitry Karpeev   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1810d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1811d709ea83SDmitry Karpeev   if (match) {
18128f3b4b4dSDmitry Karpeev     if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) {
1813d709ea83SDmitry Karpeev       osm->dm_subdomains = flg;
1814d709ea83SDmitry Karpeev     }
18158f3b4b4dSDmitry Karpeev   }
1816d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1817d709ea83SDmitry Karpeev }
1818d709ea83SDmitry Karpeev 
1819d709ea83SDmitry Karpeev #undef __FUNCT__
18208f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMGetUseDMSubdomains"
1821d709ea83SDmitry Karpeev /*@
18228f3b4b4dSDmitry Karpeev     PCGASMGetUseDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1823d709ea83SDmitry Karpeev     Not Collective
1824d709ea83SDmitry Karpeev 
1825d709ea83SDmitry Karpeev     Input Parameter:
1826d709ea83SDmitry Karpeev .   pc  - the preconditioner
1827d709ea83SDmitry Karpeev 
1828d709ea83SDmitry Karpeev     Output Parameter:
1829d709ea83SDmitry Karpeev .   flg - boolean indicating whether to use subdomains defined by the DM
1830d709ea83SDmitry Karpeev 
1831d709ea83SDmitry Karpeev     Level: intermediate
1832d709ea83SDmitry Karpeev 
1833d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1834d709ea83SDmitry Karpeev 
18358f3b4b4dSDmitry Karpeev .seealso: PCGASMSetUseDMSubdomains(), PCGASMSetOverlap()
1836d709ea83SDmitry Karpeev           PCGASMCreateSubdomains2D()
1837d709ea83SDmitry Karpeev @*/
18388f3b4b4dSDmitry Karpeev PetscErrorCode  PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg)
1839d709ea83SDmitry Karpeev {
1840d709ea83SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
1841d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1842d709ea83SDmitry Karpeev   PetscBool      match;
1843d709ea83SDmitry Karpeev 
1844d709ea83SDmitry Karpeev   PetscFunctionBegin;
1845d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1846d709ea83SDmitry Karpeev   PetscValidPointer(flg,2);
1847d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1848d709ea83SDmitry Karpeev   if (match) {
1849d709ea83SDmitry Karpeev     if (flg) *flg = osm->dm_subdomains;
1850d709ea83SDmitry Karpeev   }
1851d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1852d709ea83SDmitry Karpeev }
1853