xref: /petsc/src/ksp/pc/impls/gasm/gasm.c (revision b9d0fdaa8b8fd04b16ee28e86ef9a83a5d6857d0)
1b1a0a8a3SJed Brown /*
2f746d493SDmitry Karpeev   This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation.
36a4f0f73SDmitry Karpeev   In this version each processor may intersect multiple subdomains and any subdomain may
46a4f0f73SDmitry Karpeev   intersect multiple processors.  Intersections of subdomains with processors are called *local
56a4f0f73SDmitry Karpeev   subdomains*.
6b1a0a8a3SJed Brown 
78f3b4b4dSDmitry Karpeev        N    - total number of distinct global subdomains          (set explicitly in PCGASMSetTotalSubdomains() or implicitly PCGASMSetSubdomains() and then calculated in PCSetUp_GASM())
86a4f0f73SDmitry Karpeev        n    - actual number of local subdomains on this processor (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains())
98f3b4b4dSDmitry Karpeev        nmax - maximum number of local subdomains per processor    (calculated in PCSetUp_GASM())
10b1a0a8a3SJed Brown */
11af0996ceSBarry Smith #include <petsc/private/pcimpl.h>     /*I "petscpc.h" I*/
121e25c274SJed Brown #include <petscdm.h>
13b1a0a8a3SJed Brown 
14b1a0a8a3SJed Brown typedef struct {
15f746d493SDmitry Karpeev   PetscInt    N,n,nmax;
16b1a0a8a3SJed Brown   PetscInt    overlap;                /* overlap requested by user */
178f3b4b4dSDmitry Karpeev   PCGASMType  type;                   /* use reduced interpolation, restriction or both */
188f3b4b4dSDmitry Karpeev   PetscBool   type_set;               /* if user set this value (so won't change it for symmetric problems) */
198f3b4b4dSDmitry Karpeev   PetscBool   same_subdomain_solvers; /* flag indicating whether all local solvers are same */
208f3b4b4dSDmitry Karpeev   PetscBool   sort_indices;           /* flag to sort subdomain indices */
218f3b4b4dSDmitry Karpeev   PetscBool   user_subdomains;        /* whether the user set explicit subdomain index sets -- keep them on PCReset() */
228f3b4b4dSDmitry Karpeev   PetscBool   dm_subdomains;          /* whether DM is allowed to define subdomains */
238f3b4b4dSDmitry Karpeev   IS          *ois;                   /* index sets that define the outer (conceptually, overlapping) subdomains */
248f3b4b4dSDmitry Karpeev   IS          *iis;                   /* index sets that define the inner (conceptually, nonoverlapping) subdomains */
258f3b4b4dSDmitry Karpeev   KSP         *ksp;                /* linear solvers for each subdomain */
268f3b4b4dSDmitry Karpeev   Mat         *pmat;               /* subdomain block matrices */
27f746d493SDmitry Karpeev   Vec         gx,gy;               /* Merged work vectors */
28f746d493SDmitry Karpeev   Vec         *x,*y;               /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */
296a4f0f73SDmitry Karpeev   VecScatter  gorestriction;       /* merged restriction to disjoint union of outer subdomains */
306a4f0f73SDmitry Karpeev   VecScatter  girestriction;       /* merged restriction to disjoint union of inner subdomains */
31f746d493SDmitry Karpeev } PC_GASM;
32b1a0a8a3SJed Brown 
33b1a0a8a3SJed Brown #undef __FUNCT__
348f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMComputeGlobalSubdomainNumbering_Private"
358f3b4b4dSDmitry Karpeev static PetscErrorCode  PCGASMComputeGlobalSubdomainNumbering_Private(PC pc,PetscInt **numbering,PetscInt **permutation)
368f3b4b4dSDmitry Karpeev {
378f3b4b4dSDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
388f3b4b4dSDmitry Karpeev   PetscInt       i;
398f3b4b4dSDmitry Karpeev   PetscErrorCode ierr;
408f3b4b4dSDmitry Karpeev 
418f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
428f3b4b4dSDmitry Karpeev   /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */
438f3b4b4dSDmitry Karpeev   ierr = PetscMalloc2(osm->n,numbering,osm->n,permutation);CHKERRQ(ierr);
448f3b4b4dSDmitry Karpeev   ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->iis,NULL,*numbering);CHKERRQ(ierr);
458f3b4b4dSDmitry Karpeev   for (i = 0; i < osm->n; ++i) (*permutation)[i] = i;
468f3b4b4dSDmitry Karpeev   ierr = PetscSortIntWithPermutation(osm->n,*numbering,*permutation);CHKERRQ(ierr);
478f3b4b4dSDmitry Karpeev   PetscFunctionReturn(0);
488f3b4b4dSDmitry Karpeev }
498f3b4b4dSDmitry Karpeev 
508f3b4b4dSDmitry Karpeev #undef __FUNCT__
5106b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSubdomainView_Private"
5206b43e7eSDmitry Karpeev static PetscErrorCode  PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer)
53af538404SDmitry Karpeev {
54af538404SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
5506b43e7eSDmitry Karpeev   PetscInt       j,nidx;
56af538404SDmitry Karpeev   const PetscInt *idx;
5706b43e7eSDmitry Karpeev   PetscViewer    sviewer;
58af538404SDmitry Karpeev   char           *cidx;
59af538404SDmitry Karpeev   PetscErrorCode ierr;
60af538404SDmitry Karpeev 
61af538404SDmitry Karpeev   PetscFunctionBegin;
62ce94432eSBarry Smith   if (i < 0 || i > osm->n) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Invalid subdomain %D: must nonnegative and less than %D", i, osm->n);
634bde246dSDmitry Karpeev   /* Inner subdomains. */
644bde246dSDmitry Karpeev   ierr = ISGetLocalSize(osm->iis[i], &nidx);CHKERRQ(ierr);
654bde246dSDmitry Karpeev   /*
664bde246dSDmitry Karpeev    No more than 15 characters per index plus a space.
674bde246dSDmitry Karpeev    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
684bde246dSDmitry Karpeev    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
694bde246dSDmitry Karpeev    For nidx == 0, the whole string 16 '\0'.
704bde246dSDmitry Karpeev    */
7189d949e2SBarry Smith   ierr = PetscMalloc1(16*(nidx+1)+1, &cidx);CHKERRQ(ierr);
724bde246dSDmitry Karpeev   ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr);
734bde246dSDmitry Karpeev   ierr = ISGetIndices(osm->iis[i], &idx);CHKERRQ(ierr);
744bde246dSDmitry Karpeev   for (j = 0; j < nidx; ++j) {
754bde246dSDmitry Karpeev     ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr);
764bde246dSDmitry Karpeev   }
774bde246dSDmitry Karpeev   ierr = ISRestoreIndices(osm->iis[i],&idx);CHKERRQ(ierr);
784bde246dSDmitry Karpeev   ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
794bde246dSDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n");CHKERRQ(ierr);
804bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
814bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
824bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
834bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
844bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
854bde246dSDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
864bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
874bde246dSDmitry Karpeev   ierr = PetscFree(cidx);CHKERRQ(ierr);
884bde246dSDmitry Karpeev   /* Outer subdomains. */
896a4f0f73SDmitry Karpeev   ierr = ISGetLocalSize(osm->ois[i], &nidx);CHKERRQ(ierr);
90eec7e3faSDmitry Karpeev   /*
91eec7e3faSDmitry Karpeev    No more than 15 characters per index plus a space.
92eec7e3faSDmitry Karpeev    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
93eec7e3faSDmitry Karpeev    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
94eec7e3faSDmitry Karpeev    For nidx == 0, the whole string 16 '\0'.
95eec7e3faSDmitry Karpeev    */
9689d949e2SBarry Smith   ierr = PetscMalloc1(16*(nidx+1)+1, &cidx);CHKERRQ(ierr);
9706b43e7eSDmitry Karpeev   ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr);
986a4f0f73SDmitry Karpeev   ierr = ISGetIndices(osm->ois[i], &idx);CHKERRQ(ierr);
99af538404SDmitry Karpeev   for (j = 0; j < nidx; ++j) {
100af538404SDmitry Karpeev     ierr = PetscViewerStringSPrintf(sviewer,"%D ", idx[j]);CHKERRQ(ierr);
101af538404SDmitry Karpeev   }
1026bf464f9SBarry Smith   ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
1036a4f0f73SDmitry Karpeev   ierr = ISRestoreIndices(osm->ois[i],&idx);CHKERRQ(ierr);
1046a4f0f73SDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n");CHKERRQ(ierr);
10506b43e7eSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1067b23a99aSBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
107af538404SDmitry Karpeev   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
108af538404SDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1097b23a99aSBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
110af538404SDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
11106b43e7eSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
112af538404SDmitry Karpeev   ierr = PetscFree(cidx);CHKERRQ(ierr);
11306b43e7eSDmitry Karpeev   PetscFunctionReturn(0);
11406b43e7eSDmitry Karpeev }
11506b43e7eSDmitry Karpeev 
11606b43e7eSDmitry Karpeev #undef __FUNCT__
11706b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains"
11806b43e7eSDmitry Karpeev static PetscErrorCode  PCGASMPrintSubdomains(PC pc)
11906b43e7eSDmitry Karpeev {
12006b43e7eSDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
12106b43e7eSDmitry Karpeev   const char     *prefix;
12206b43e7eSDmitry Karpeev   char           fname[PETSC_MAX_PATH_LEN+1];
1238f3b4b4dSDmitry Karpeev   PetscInt       l, d, count;
1248f3b4b4dSDmitry Karpeev   PetscBool      doprint,found;
1250298fd71SBarry Smith   PetscViewer    viewer, sviewer = NULL;
1268f3b4b4dSDmitry Karpeev   PetscInt       *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */
12706b43e7eSDmitry Karpeev   PetscErrorCode ierr;
12806b43e7eSDmitry Karpeev 
12906b43e7eSDmitry Karpeev   PetscFunctionBegin;
13006b43e7eSDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
1318f3b4b4dSDmitry Karpeev   doprint  = PETSC_FALSE;
1328f3b4b4dSDmitry Karpeev   ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&doprint,NULL);CHKERRQ(ierr);
1338f3b4b4dSDmitry Karpeev   if (!doprint) PetscFunctionReturn(0);
13406b43e7eSDmitry Karpeev   ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr);
13506b43e7eSDmitry Karpeev   if (!found) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
136ce94432eSBarry Smith   ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr);
1374bde246dSDmitry Karpeev   /*
1384bde246dSDmitry Karpeev    Make sure the viewer has a name. Otherwise this may cause a deadlock or other weird errors when creating a subcomm viewer:
1394bde246dSDmitry Karpeev    the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed collectively on the comm.
1404bde246dSDmitry Karpeev   */
1414bde246dSDmitry Karpeev   ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr);
1424bde246dSDmitry Karpeev   l    = 0;
1438f3b4b4dSDmitry Karpeev   ierr = PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);CHKERRQ(ierr);
1448f3b4b4dSDmitry Karpeev   for (count = 0; count < osm->N; ++count) {
1454bde246dSDmitry Karpeev     /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
1464bde246dSDmitry Karpeev     if (l<osm->n) {
1474bde246dSDmitry Karpeev       d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
1484bde246dSDmitry Karpeev       if (numbering[d] == count) {
1494bde246dSDmitry Karpeev         ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
1504bde246dSDmitry Karpeev         ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr);
1514bde246dSDmitry Karpeev         ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
1524bde246dSDmitry Karpeev         ++l;
153af538404SDmitry Karpeev       }
1544bde246dSDmitry Karpeev     }
155ce94432eSBarry Smith     ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
1564bde246dSDmitry Karpeev   }
1578f3b4b4dSDmitry Karpeev   ierr = PetscFree2(numbering,permutation);CHKERRQ(ierr);
1584bde246dSDmitry Karpeev   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
159af538404SDmitry Karpeev   PetscFunctionReturn(0);
160af538404SDmitry Karpeev }
161af538404SDmitry Karpeev 
162af538404SDmitry Karpeev 
163af538404SDmitry Karpeev #undef __FUNCT__
164f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM"
165f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer)
166b1a0a8a3SJed Brown {
167f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
168af538404SDmitry Karpeev   const char     *prefix;
169b1a0a8a3SJed Brown   PetscErrorCode ierr;
170af538404SDmitry Karpeev   PetscMPIInt    rank, size;
1718f3b4b4dSDmitry Karpeev   PetscInt       bsz;
17206b43e7eSDmitry Karpeev   PetscBool      iascii,view_subdomains=PETSC_FALSE;
173b1a0a8a3SJed Brown   PetscViewer    sviewer;
1748f3b4b4dSDmitry Karpeev   PetscInt       count, l;
1756a4f0f73SDmitry Karpeev   char           overlap[256]     = "user-defined overlap";
1766a4f0f73SDmitry Karpeev   char           gsubdomains[256] = "unknown total number of subdomains";
17706b43e7eSDmitry Karpeev   char           msubdomains[256] = "unknown max number of local subdomains";
1788f3b4b4dSDmitry Karpeev   PetscInt       *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */
17911aeaf0aSBarry Smith 
180b1a0a8a3SJed Brown   PetscFunctionBegin;
181ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr);
182ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr);
18306b43e7eSDmitry Karpeev 
18406b43e7eSDmitry Karpeev   if (osm->overlap >= 0) {
18506b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %D",osm->overlap);CHKERRQ(ierr);
18606b43e7eSDmitry Karpeev   }
1878f3b4b4dSDmitry Karpeev   if (osm->N != PETSC_DETERMINE) {
1888f3b4b4dSDmitry Karpeev     ierr = PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %D",osm->N);CHKERRQ(ierr);
18906b43e7eSDmitry Karpeev   }
1908f3b4b4dSDmitry Karpeev   if (osm->nmax != PETSC_DETERMINE) {
19106b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr);
19206b43e7eSDmitry Karpeev   }
19306b43e7eSDmitry Karpeev 
194af538404SDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
1950298fd71SBarry Smith   ierr = PetscOptionsGetBool(prefix,"-pc_gasm_view_subdomains",&view_subdomains,NULL);CHKERRQ(ierr);
19606b43e7eSDmitry Karpeev 
19706b43e7eSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
198b1a0a8a3SJed Brown   if (iascii) {
19906b43e7eSDmitry Karpeev     /*
20006b43e7eSDmitry Karpeev      Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer:
20106b43e7eSDmitry Karpeev      the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed
20206b43e7eSDmitry Karpeev      collectively on the comm.
20306b43e7eSDmitry Karpeev      */
20406b43e7eSDmitry Karpeev     ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr);
20506b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Generalized additive Schwarz:\n");CHKERRQ(ierr);
20606b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Restriction/interpolation type: %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr);
20706b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",overlap);CHKERRQ(ierr);
20806b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",gsubdomains);CHKERRQ(ierr);
20906b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",msubdomains);CHKERRQ(ierr);
2107b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
2118f3b4b4dSDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d|%d] number of locally-supported subdomains = %D\n",rank,size,osm->n);CHKERRQ(ierr);
212af538404SDmitry Karpeev     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
2137b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
2146a4f0f73SDmitry Karpeev     /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */
21506b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Subdomain solver info is as follows:\n");CHKERRQ(ierr);
216b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
217b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
21806b43e7eSDmitry Karpeev     /* Make sure that everybody waits for the banner to be printed. */
219ce94432eSBarry Smith     ierr = MPI_Barrier(PetscObjectComm((PetscObject)viewer));CHKERRQ(ierr);
2204bde246dSDmitry Karpeev     /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
2218f3b4b4dSDmitry Karpeev     ierr = PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);CHKERRQ(ierr);
22206b43e7eSDmitry Karpeev     l = 0;
2238f3b4b4dSDmitry Karpeev     for (count = 0; count < osm->N; ++count) {
22406b43e7eSDmitry Karpeev       PetscMPIInt srank, ssize;
22506b43e7eSDmitry Karpeev       if (l<osm->n) {
22606b43e7eSDmitry Karpeev         PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
22706b43e7eSDmitry Karpeev         if (numbering[d] == count) {
2286a4f0f73SDmitry Karpeev           ierr = MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize);CHKERRQ(ierr);
2296a4f0f73SDmitry Karpeev           ierr = MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank);CHKERRQ(ierr);
2306a4f0f73SDmitry Karpeev           ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
2316a4f0f73SDmitry Karpeev           ierr = ISGetLocalSize(osm->ois[d],&bsz);CHKERRQ(ierr);
2327b23a99aSBarry Smith           ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_TRUE);CHKERRQ(ierr);
2338f3b4b4dSDmitry Karpeev           ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d|%d] (subcomm [%d|%d]) local subdomain number %D, local size = %D\n",rank,size,srank,ssize,d,bsz);CHKERRQ(ierr);
2346a4f0f73SDmitry Karpeev           ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
2356a4f0f73SDmitry Karpeev           ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_FALSE);CHKERRQ(ierr);
23606b43e7eSDmitry Karpeev           if (view_subdomains) {
23706b43e7eSDmitry Karpeev             ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr);
23806b43e7eSDmitry Karpeev           }
2396a4f0f73SDmitry Karpeev           if (!pc->setupcalled) {
2406a4f0f73SDmitry Karpeev             PetscViewerASCIIPrintf(sviewer, "Solver not set up yet: PCSetUp() not yet called\n");CHKERRQ(ierr);
2412fa5cd67SKarl Rupp           } else {
24206b43e7eSDmitry Karpeev             ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr);
2436a4f0f73SDmitry Karpeev           }
24406b43e7eSDmitry Karpeev           ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
2457b23a99aSBarry Smith           ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
2466a4f0f73SDmitry Karpeev           ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
24706b43e7eSDmitry Karpeev           ++l;
248b1a0a8a3SJed Brown         }
24906b43e7eSDmitry Karpeev       }
250ce94432eSBarry Smith       ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
251b1a0a8a3SJed Brown     }
2528f3b4b4dSDmitry Karpeev     ierr = PetscFree2(numbering,permutation);CHKERRQ(ierr);
253b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
254b1a0a8a3SJed Brown   }
255b1a0a8a3SJed Brown   PetscFunctionReturn(0);
256b1a0a8a3SJed Brown }
257b1a0a8a3SJed Brown 
2588f3b4b4dSDmitry Karpeev PETSC_INTERN PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]);
259af538404SDmitry Karpeev 
260b1a0a8a3SJed Brown #undef __FUNCT__
261f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM"
262f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc)
263b1a0a8a3SJed Brown {
264f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
265b1a0a8a3SJed Brown   PetscErrorCode ierr;
266b1a0a8a3SJed Brown   PetscBool      symset,flg;
26788389755SJed Brown   PetscInt       i;
268c10bc1a1SDmitry Karpeev   PetscMPIInt    rank, size;
269b1a0a8a3SJed Brown   MatReuse       scall = MAT_REUSE_MATRIX;
270b1a0a8a3SJed Brown   KSP            ksp;
271b1a0a8a3SJed Brown   PC             subpc;
272b1a0a8a3SJed Brown   const char     *prefix,*pprefix;
273f746d493SDmitry Karpeev   Vec            x,y;
2746a4f0f73SDmitry Karpeev   PetscInt       oni;       /* Number of indices in the i-th local outer subdomain.               */
2756a4f0f73SDmitry Karpeev   const PetscInt *oidxi;    /* Indices from the i-th subdomain local outer subdomain.             */
2766a4f0f73SDmitry Karpeev   PetscInt       on;        /* Number of indices in the disjoint union of local outer subdomains. */
2776a4f0f73SDmitry Karpeev   PetscInt       *oidx;     /* Indices in the disjoint union of local outer subdomains. */
2786a4f0f73SDmitry Karpeev   IS             gois;      /* Disjoint union the global indices of outer subdomains.             */
2796a4f0f73SDmitry Karpeev   IS             goid;      /* Identity IS of the size of the disjoint union of outer subdomains. */
280f746d493SDmitry Karpeev   PetscScalar    *gxarray, *gyarray;
281930d09c1SFande Kong   PetscInt       gostart;   /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */
2828f3b4b4dSDmitry Karpeev   PetscInt       num_subdomains    = 0;
2830298fd71SBarry Smith   DM             *subdomain_dm     = NULL;
2848f3b4b4dSDmitry Karpeev   char           **subdomain_names = NULL;
285f771a274SFande Kong   PetscInt       *numbering;
2868f3b4b4dSDmitry Karpeev 
287b1a0a8a3SJed Brown 
288b1a0a8a3SJed Brown   PetscFunctionBegin;
289ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
290ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
291b1a0a8a3SJed Brown   if (!pc->setupcalled) {
292b1a0a8a3SJed Brown 
293b1a0a8a3SJed Brown     if (!osm->type_set) {
294b1a0a8a3SJed Brown       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
2952fa5cd67SKarl Rupp       if (symset && flg) osm->type = PC_GASM_BASIC;
296b1a0a8a3SJed Brown     }
297b1a0a8a3SJed Brown 
2988f3b4b4dSDmitry Karpeev     if (osm->n == PETSC_DETERMINE) {
2998f3b4b4dSDmitry Karpeev       if (osm->N != PETSC_DETERMINE) {
3008f3b4b4dSDmitry Karpeev 	   /* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */
3018f3b4b4dSDmitry Karpeev 	   ierr = PCGASMCreateSubdomains(pc->pmat,osm->N,&osm->n,&osm->iis);CHKERRQ(ierr);
3028f3b4b4dSDmitry Karpeev       } else if (osm->dm_subdomains && pc->dm) {
3038f3b4b4dSDmitry Karpeev 	/* try pc->dm next, if allowed */
3048f3b4b4dSDmitry Karpeev 	PetscInt  d;
305a35b7b57SDmitry Karpeev 	IS       *inner_subdomain_is, *outer_subdomain_is;
306a35b7b57SDmitry Karpeev 	ierr = DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);CHKERRQ(ierr);
307a35b7b57SDmitry Karpeev 	if (num_subdomains) {
308a35b7b57SDmitry Karpeev 	  ierr = PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);CHKERRQ(ierr);
30969ca1f37SDmitry Karpeev 	}
310a35b7b57SDmitry Karpeev 	for (d = 0; d < num_subdomains; ++d) {
311a35b7b57SDmitry Karpeev 	  if (inner_subdomain_is) {ierr = ISDestroy(&inner_subdomain_is[d]);CHKERRQ(ierr);}
312a35b7b57SDmitry Karpeev 	  if (outer_subdomain_is) {ierr = ISDestroy(&outer_subdomain_is[d]);CHKERRQ(ierr);}
313fdc48646SDmitry Karpeev 	}
314a35b7b57SDmitry Karpeev 	ierr = PetscFree(inner_subdomain_is);CHKERRQ(ierr);
315a35b7b57SDmitry Karpeev 	ierr = PetscFree(outer_subdomain_is);CHKERRQ(ierr);
3168f3b4b4dSDmitry Karpeev       } else {
3178f3b4b4dSDmitry Karpeev 	/* still no subdomains; use one per processor */
31844fe18e5SDmitry Karpeev 	osm->nmax = osm->n = 1;
319ce94432eSBarry Smith 	ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
320f746d493SDmitry Karpeev 	osm->N    = size;
3218f3b4b4dSDmitry Karpeev 	ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr);
322fdc48646SDmitry Karpeev       }
32306b43e7eSDmitry Karpeev     }
324a35b7b57SDmitry Karpeev     if (!osm->iis) {
325a35b7b57SDmitry Karpeev       /*
3268f3b4b4dSDmitry Karpeev        osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied.
3278f3b4b4dSDmitry Karpeev        We create the requisite number of local inner subdomains and then expand them into
3288f3b4b4dSDmitry Karpeev        out subdomains, if necessary.
329a35b7b57SDmitry Karpeev        */
3308f3b4b4dSDmitry Karpeev       ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr);
331a35b7b57SDmitry Karpeev     }
3328f3b4b4dSDmitry Karpeev     if (!osm->ois) {
3338f3b4b4dSDmitry Karpeev       /*
3348f3b4b4dSDmitry Karpeev 	Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap
3358f3b4b4dSDmitry Karpeev 	has been requested, copy the inner subdomains over so they can be modified.
3368f3b4b4dSDmitry Karpeev       */
3378f3b4b4dSDmitry Karpeev       ierr = PetscMalloc1(osm->n,&osm->ois);CHKERRQ(ierr);
3388f3b4b4dSDmitry Karpeev       for (i=0; i<osm->n; ++i) {
3398f3b4b4dSDmitry Karpeev 	if (osm->overlap > 0) { /* With positive overlap, osm->iis[i] will be modified */
3408f3b4b4dSDmitry Karpeev 	  ierr = ISDuplicate(osm->iis[i],(osm->ois)+i);CHKERRQ(ierr);
3418f3b4b4dSDmitry Karpeev 	  ierr = ISCopy(osm->iis[i],osm->ois[i]);CHKERRQ(ierr);
3428f3b4b4dSDmitry Karpeev 	} else {
3438f3b4b4dSDmitry Karpeev 	  ierr      = PetscObjectReference((PetscObject)((osm->iis)[i]));CHKERRQ(ierr);
3448f3b4b4dSDmitry Karpeev 	  osm->ois[i] = osm->iis[i];
3458f3b4b4dSDmitry Karpeev 	}
3468f3b4b4dSDmitry Karpeev       }
3478f3b4b4dSDmitry Karpeev       if (osm->overlap > 0) {
3488f3b4b4dSDmitry Karpeev 	   /* Extend the "overlapping" regions by a number of steps */
349d21f2a47SFande Kong 	   ierr = MatIncreaseOverlapSplit(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr);
3508f3b4b4dSDmitry Karpeev       }
351b1a0a8a3SJed Brown     }
3526a4f0f73SDmitry Karpeev 
3538f3b4b4dSDmitry Karpeev     /* Now the subdomains are defined.  Determine their global and max local numbers, if necessary. */
3548f3b4b4dSDmitry Karpeev     if (osm->nmax == PETSC_DETERMINE) {
3558f3b4b4dSDmitry Karpeev       PetscMPIInt inwork,outwork;
3568f3b4b4dSDmitry Karpeev       /* determine global number of subdomains and the max number of local subdomains */
3578f3b4b4dSDmitry Karpeev       inwork = osm->n;
3588f3b4b4dSDmitry Karpeev       ierr       = MPI_Allreduce(&inwork,&outwork,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
3598f3b4b4dSDmitry Karpeev       osm->nmax  = outwork;
3608f3b4b4dSDmitry Karpeev     }
3618f3b4b4dSDmitry Karpeev     if (osm->N == PETSC_DETERMINE) {
3628f3b4b4dSDmitry Karpeev       /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */
3638f3b4b4dSDmitry Karpeev       ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,&osm->N,NULL);CHKERRQ(ierr);
3648f3b4b4dSDmitry Karpeev     }
3658f3b4b4dSDmitry Karpeev 
366b1a0a8a3SJed Brown 
367b1a0a8a3SJed Brown     if (osm->sort_indices) {
368f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
3696a4f0f73SDmitry Karpeev         ierr = ISSort(osm->ois[i]);CHKERRQ(ierr);
3706a4f0f73SDmitry Karpeev         ierr = ISSort(osm->iis[i]);CHKERRQ(ierr);
371b1a0a8a3SJed Brown       }
372b1a0a8a3SJed Brown     }
3738f3b4b4dSDmitry Karpeev     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
3748f3b4b4dSDmitry Karpeev     ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr);
3758f3b4b4dSDmitry Karpeev 
3766a4f0f73SDmitry Karpeev     /*
3776a4f0f73SDmitry Karpeev      Merge the ISs, create merged vectors and restrictions.
3786a4f0f73SDmitry Karpeev      */
3796a4f0f73SDmitry Karpeev     /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */
3806a4f0f73SDmitry Karpeev     on = 0;
381f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
3826a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
3836a4f0f73SDmitry Karpeev       on  += oni;
384f746d493SDmitry Karpeev     }
385785e854fSJed Brown     ierr = PetscMalloc1(on, &oidx);CHKERRQ(ierr);
3866a4f0f73SDmitry Karpeev     on   = 0;
387930d09c1SFande Kong     /*Merge local indices together */
388f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
3896a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
3906a4f0f73SDmitry Karpeev       ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr);
3916a4f0f73SDmitry Karpeev       ierr = PetscMemcpy(oidx+on,oidxi,sizeof(PetscInt)*oni);CHKERRQ(ierr);
3926a4f0f73SDmitry Karpeev       ierr = ISRestoreIndices(osm->ois[i],&oidxi);CHKERRQ(ierr);
3936a4f0f73SDmitry Karpeev       on  += oni;
394f746d493SDmitry Karpeev     }
3956a4f0f73SDmitry Karpeev     ierr = ISCreateGeneral(((PetscObject)(pc))->comm,on,oidx,PETSC_OWN_POINTER,&gois);CHKERRQ(ierr);
3962a7a6963SBarry Smith     ierr = MatCreateVecs(pc->pmat,&x,&y);CHKERRQ(ierr);
397ce94432eSBarry Smith     ierr = VecCreateMPI(PetscObjectComm((PetscObject)pc),on,PETSC_DECIDE,&osm->gx);CHKERRQ(ierr);
398f746d493SDmitry Karpeev     ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr);
399930d09c1SFande Kong     ierr = VecGetOwnershipRange(osm->gx, &gostart, NULL);CHKERRQ(ierr);
400930d09c1SFande Kong     ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),on,gostart,1, &goid);CHKERRQ(ierr);
401930d09c1SFande Kong     /*gois might indices not on local */
4026a4f0f73SDmitry Karpeev     ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr);
403f771a274SFande Kong     ierr = PetscMalloc1(osm->n,&numbering);CHKERRQ(ierr);
404f771a274SFande Kong     ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,NULL,numbering);CHKERRQ(ierr);
4056a4f0f73SDmitry Karpeev     ierr = VecDestroy(&x);CHKERRQ(ierr);
4067105d80bSDmitry Karpeev     ierr = ISDestroy(&gois);CHKERRQ(ierr);
4072fa5cd67SKarl Rupp 
4086a4f0f73SDmitry Karpeev     /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */
4092fa5cd67SKarl Rupp     {
4102fa5cd67SKarl Rupp       PetscInt        ini;           /* Number of indices the i-th a local inner subdomain. */
4116a4f0f73SDmitry Karpeev       PetscInt        in;            /* Number of indices in the disjoint uniont of local inner subdomains. */
4126a4f0f73SDmitry Karpeev       PetscInt       *iidx;          /* Global indices in the merged local inner subdomain. */
4136a4f0f73SDmitry Karpeev       PetscInt       *ioidx;         /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
4146a4f0f73SDmitry Karpeev       IS              giis;          /* IS for the disjoint union of inner subdomains. */
4156a4f0f73SDmitry Karpeev       IS              giois;         /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
416f771a274SFande Kong       PetscScalar    *array;
417f771a274SFande Kong       const PetscInt *indices;
418f771a274SFande Kong       PetscInt        k;
419f746d493SDmitry Karpeev       /**/
4206a4f0f73SDmitry Karpeev       on = 0;
421f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
4226a4f0f73SDmitry Karpeev         ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
4236a4f0f73SDmitry Karpeev         on  += oni;
424f746d493SDmitry Karpeev       }
425f771a274SFande Kong       /*allocate memory */
426f771a274SFande Kong       ierr = PetscMalloc1(on, &iidx);CHKERRQ(ierr);
427f771a274SFande Kong       ierr = PetscMalloc1(on, &ioidx);CHKERRQ(ierr);
428f771a274SFande Kong       ierr = VecGetArray(y,&array);CHKERRQ(ierr);
429f771a274SFande Kong       /*set communicator id */
430f771a274SFande Kong       in   = 0;
431f771a274SFande Kong       for (i=0; i<osm->n; i++) {
432f771a274SFande Kong         ierr   = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr);
433f771a274SFande Kong         for (k = 0; k < ini; ++k){
434f771a274SFande Kong           array[in+k] = numbering[i];
435f771a274SFande Kong         }
436f771a274SFande Kong         in += ini;
437f771a274SFande Kong       }
438f771a274SFande Kong       ierr = VecRestoreArray(y,&array);CHKERRQ(ierr);
439f771a274SFande Kong       ierr = VecScatterBegin(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
440f771a274SFande Kong       ierr = VecScatterEnd(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
441f771a274SFande Kong #if 0
442f771a274SFande Kong       ierr = VecView(osm->gy,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
443f771a274SFande Kong #endif
444f771a274SFande Kong       ierr = VecGetOwnershipRange(osm->gy,&gostart, NULL);CHKERRQ(ierr);
445f771a274SFande Kong       ierr = VecGetArray(osm->gy,&array);CHKERRQ(ierr);
446f771a274SFande Kong       on  = 0;
447f771a274SFande Kong       in  = 0;
448f771a274SFande Kong       for(i=0; i<osm->n;i++){
449f771a274SFande Kong     	ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
450f771a274SFande Kong     	ierr = ISGetIndices(osm->ois[i],&indices);CHKERRQ(ierr);
451f771a274SFande Kong     	for(k=0; k<oni; k++){
452f771a274SFande Kong     	 /*skip overlapping indices */
453f771a274SFande Kong          if(array[on+k] != numbering[i]) continue;
454f771a274SFande Kong          /*record inner indices */
455f771a274SFande Kong          iidx[in]    = indices[k];
456f771a274SFande Kong          ioidx[in++] = gostart+on+k;
457f771a274SFande Kong     	}
458f771a274SFande Kong     	ierr   = ISRestoreIndices(osm->ois[i], &indices);CHKERRQ(ierr);
459f771a274SFande Kong     	on += oni;
460f771a274SFande Kong       }
461f771a274SFande Kong       ierr = VecRestoreArray(osm->gy,&array);CHKERRQ(ierr);
462ce94432eSBarry Smith       ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,iidx,PETSC_OWN_POINTER,&giis);CHKERRQ(ierr);
463ce94432eSBarry Smith       ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,ioidx,PETSC_OWN_POINTER,&giois);CHKERRQ(ierr);
4646a4f0f73SDmitry Karpeev       ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr);
4656a4f0f73SDmitry Karpeev       ierr = VecDestroy(&y);CHKERRQ(ierr);
4666a4f0f73SDmitry Karpeev       ierr = ISDestroy(&giis);CHKERRQ(ierr);
4676a4f0f73SDmitry Karpeev       ierr = ISDestroy(&giois);CHKERRQ(ierr);
468b1a0a8a3SJed Brown     }
4696a4f0f73SDmitry Karpeev     ierr = ISDestroy(&goid);CHKERRQ(ierr);
470f771a274SFande Kong     ierr = PetscFree(numbering);CHKERRQ(ierr);
4712fa5cd67SKarl Rupp 
472f746d493SDmitry Karpeev     /* Create the subdomain work vectors. */
473785e854fSJed Brown     ierr = PetscMalloc1(osm->n,&osm->x);CHKERRQ(ierr);
474785e854fSJed Brown     ierr = PetscMalloc1(osm->n,&osm->y);CHKERRQ(ierr);
475f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr);
476f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr);
4776a4f0f73SDmitry Karpeev     for (i=0, on=0; i<osm->n; ++i, on += oni) {
4786a4f0f73SDmitry Karpeev       PetscInt oNi;
4796a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
480930d09c1SFande Kong       /*on a sub communicator */
4816a4f0f73SDmitry Karpeev       ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr);
4826a4f0f73SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr);
4836a4f0f73SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr);
484b1a0a8a3SJed Brown     }
485f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr);
486f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr);
4878f3b4b4dSDmitry Karpeev     /* Create the subdomain solvers */
488785e854fSJed Brown     ierr = PetscMalloc1(osm->n,&osm->ksp);CHKERRQ(ierr);
4898f3b4b4dSDmitry Karpeev     for (i=0; i<osm->n; i++) {
4908f3b4b4dSDmitry Karpeev       char subprefix[PETSC_MAX_PATH_LEN+1];
4916a4f0f73SDmitry Karpeev       ierr        = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr);
4923bb1ff40SBarry Smith       ierr        = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr);
493b1a0a8a3SJed Brown       ierr        = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
494b1a0a8a3SJed Brown       ierr        = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
4958f3b4b4dSDmitry Karpeev       ierr        = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); /* Why do we need this here? */
4968f3b4b4dSDmitry Karpeev       if (subdomain_dm) {
4978f3b4b4dSDmitry Karpeev 	    ierr = KSPSetDM(ksp,subdomain_dm[i]);CHKERRQ(ierr);
4988f3b4b4dSDmitry Karpeev 	    ierr = DMDestroy(subdomain_dm+i);CHKERRQ(ierr);
4998f3b4b4dSDmitry Karpeev       }
500b1a0a8a3SJed Brown       ierr        = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
501b1a0a8a3SJed Brown       ierr        = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
5028f3b4b4dSDmitry Karpeev       if (subdomain_names && subdomain_names[i]) {
5038f3b4b4dSDmitry Karpeev 	     ierr = PetscSNPrintf(subprefix,PETSC_MAX_PATH_LEN,"sub_%s_",subdomain_names[i]);CHKERRQ(ierr);
5048f3b4b4dSDmitry Karpeev 	     ierr = KSPAppendOptionsPrefix(ksp,subprefix);CHKERRQ(ierr);
5058f3b4b4dSDmitry Karpeev 	     ierr = PetscFree(subdomain_names[i]);CHKERRQ(ierr);
5068f3b4b4dSDmitry Karpeev       }
507b1a0a8a3SJed Brown       ierr        = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
508b1a0a8a3SJed Brown       osm->ksp[i] = ksp;
509b1a0a8a3SJed Brown     }
5108f3b4b4dSDmitry Karpeev     ierr = PetscFree(subdomain_dm);CHKERRQ(ierr);
5118f3b4b4dSDmitry Karpeev     ierr = PetscFree(subdomain_names);CHKERRQ(ierr);
512b1a0a8a3SJed Brown     scall = MAT_INITIAL_MATRIX;
513b1a0a8a3SJed Brown 
5148f3b4b4dSDmitry Karpeev   } else { /*if (pc->setupcalled)*/
515b1a0a8a3SJed Brown     /*
5168f3b4b4dSDmitry Karpeev        Destroy the submatrices from the previous iteration
517b1a0a8a3SJed Brown     */
518b1a0a8a3SJed Brown     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
519f746d493SDmitry Karpeev       ierr  = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
520b1a0a8a3SJed Brown       scall = MAT_INITIAL_MATRIX;
521b1a0a8a3SJed Brown     }
522b1a0a8a3SJed Brown   }
523b1a0a8a3SJed Brown 
524b1a0a8a3SJed Brown   /*
525f746d493SDmitry Karpeev      Extract out the submatrices.
526b1a0a8a3SJed Brown   */
52781a5c4aaSDmitry Karpeev   if (size > 1) {
5288f3b4b4dSDmitry Karpeev     ierr = MatGetSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr);
5292fa5cd67SKarl Rupp   } else {
5306a4f0f73SDmitry Karpeev     ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr);
53181a5c4aaSDmitry Karpeev   }
532b1a0a8a3SJed Brown   if (scall == MAT_INITIAL_MATRIX) {
533b1a0a8a3SJed Brown     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
534f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
5353bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr);
536b1a0a8a3SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
537b1a0a8a3SJed Brown     }
538b1a0a8a3SJed Brown   }
539b1a0a8a3SJed Brown 
540b1a0a8a3SJed Brown   /* Return control to the user so that the submatrices can be modified (e.g., to apply
541b1a0a8a3SJed Brown      different boundary conditions for the submatrices than for the global problem) */
5426a4f0f73SDmitry Karpeev   ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
543b1a0a8a3SJed Brown 
544b1a0a8a3SJed Brown   /*
5456a4f0f73SDmitry Karpeev      Loop over submatrices putting them into local ksps
546b1a0a8a3SJed Brown   */
547f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
54823ee1639SBarry Smith     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr);
549b1a0a8a3SJed Brown     if (!pc->setupcalled) {
550b1a0a8a3SJed Brown       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
551b1a0a8a3SJed Brown     }
552b1a0a8a3SJed Brown   }
553b1a0a8a3SJed Brown   PetscFunctionReturn(0);
554b1a0a8a3SJed Brown }
555b1a0a8a3SJed Brown 
556b1a0a8a3SJed Brown #undef __FUNCT__
557f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM"
558f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
559b1a0a8a3SJed Brown {
560f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
561b1a0a8a3SJed Brown   PetscErrorCode ierr;
562b1a0a8a3SJed Brown   PetscInt       i;
563b1a0a8a3SJed Brown 
564b1a0a8a3SJed Brown   PetscFunctionBegin;
565f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
566b1a0a8a3SJed Brown     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
567b1a0a8a3SJed Brown   }
568b1a0a8a3SJed Brown   PetscFunctionReturn(0);
569b1a0a8a3SJed Brown }
570b1a0a8a3SJed Brown 
571b1a0a8a3SJed Brown #undef __FUNCT__
572f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM"
573f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y)
574b1a0a8a3SJed Brown {
575f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
576b1a0a8a3SJed Brown   PetscErrorCode ierr;
577f746d493SDmitry Karpeev   PetscInt       i;
578b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
579b1a0a8a3SJed Brown 
580b1a0a8a3SJed Brown   PetscFunctionBegin;
581b1a0a8a3SJed Brown   /*
5826a4f0f73SDmitry Karpeev      Support for limiting the restriction or interpolation only to the inner
583b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
584b1a0a8a3SJed Brown   */
585f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
586b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
587f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
5886a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5892fa5cd67SKarl Rupp   } else {
5906a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
591b1a0a8a3SJed Brown   }
5926a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr);
5936a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
5946a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5952fa5cd67SKarl Rupp   } else {
5966a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5976a4f0f73SDmitry Karpeev   }
59886b96d36SDmitry Karpeev   /* do the subdomain solves */
59986b96d36SDmitry Karpeev   for (i=0; i<osm->n; ++i) {
600b1a0a8a3SJed Brown     ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
601b1a0a8a3SJed Brown   }
602930d09c1SFande Kong   /*Do we need to zero y ?? */
6036a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(y);CHKERRQ(ierr);
6046a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
6056a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6066a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);  PetscFunctionReturn(0);
6072fa5cd67SKarl Rupp   } else {
6086a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6096a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);  PetscFunctionReturn(0);
6106a4f0f73SDmitry Karpeev   }
611b1a0a8a3SJed Brown }
612b1a0a8a3SJed Brown 
613b1a0a8a3SJed Brown #undef __FUNCT__
614f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM"
615f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y)
616b1a0a8a3SJed Brown {
617f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
618b1a0a8a3SJed Brown   PetscErrorCode ierr;
619f746d493SDmitry Karpeev   PetscInt       i;
620b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
621b1a0a8a3SJed Brown 
622b1a0a8a3SJed Brown   PetscFunctionBegin;
623b1a0a8a3SJed Brown   /*
624b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
625b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
626b1a0a8a3SJed Brown 
627f746d493SDmitry Karpeev      Note: these are reversed from the PCApply_GASM() because we are applying the
628b1a0a8a3SJed Brown      transpose of the three terms
629b1a0a8a3SJed Brown   */
630f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
631b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
632f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
6336a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
6342fa5cd67SKarl Rupp   } else {
6356a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
636b1a0a8a3SJed Brown   }
6376a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr);
6386a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
6396a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
6402fa5cd67SKarl Rupp   } else {
6416a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
6426a4f0f73SDmitry Karpeev   }
643b1a0a8a3SJed Brown   /* do the local solves */
644ab3e923fSDmitry 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. */
645b1a0a8a3SJed Brown     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
646b1a0a8a3SJed Brown   }
6476a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(y);CHKERRQ(ierr);
6486a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
6496a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6506a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6512fa5cd67SKarl Rupp   } else {
6526a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6536a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6546a4f0f73SDmitry Karpeev   }
655b1a0a8a3SJed Brown   PetscFunctionReturn(0);
656b1a0a8a3SJed Brown }
657b1a0a8a3SJed Brown 
658b1a0a8a3SJed Brown #undef __FUNCT__
659a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM"
660a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc)
661b1a0a8a3SJed Brown {
662f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
663b1a0a8a3SJed Brown   PetscErrorCode ierr;
664b1a0a8a3SJed Brown   PetscInt       i;
665b1a0a8a3SJed Brown 
666b1a0a8a3SJed Brown   PetscFunctionBegin;
667b1a0a8a3SJed Brown   if (osm->ksp) {
668f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
669a06653b4SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
670b1a0a8a3SJed Brown     }
671b1a0a8a3SJed Brown   }
672b1a0a8a3SJed Brown   if (osm->pmat) {
673f746d493SDmitry Karpeev     if (osm->n > 0) {
674ab3e923fSDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
675b1a0a8a3SJed Brown     }
676b1a0a8a3SJed Brown   }
677d34fcf5fSBarry Smith   if (osm->x) {
678f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
679fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
680fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
681b1a0a8a3SJed Brown     }
682d34fcf5fSBarry Smith   }
6836bf464f9SBarry Smith   ierr = VecDestroy(&osm->gx);CHKERRQ(ierr);
6846bf464f9SBarry Smith   ierr = VecDestroy(&osm->gy);CHKERRQ(ierr);
685ab3e923fSDmitry Karpeev 
6866a4f0f73SDmitry Karpeev   ierr     = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr);
6876a4f0f73SDmitry Karpeev   ierr     = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr);
6888f3b4b4dSDmitry Karpeev   if (!osm->user_subdomains) {
6892c112581SDmitry Karpeev     ierr      = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr);
6908f3b4b4dSDmitry Karpeev     osm->N    = PETSC_DETERMINE;
6918f3b4b4dSDmitry Karpeev     osm->nmax = PETSC_DETERMINE;
6928f3b4b4dSDmitry Karpeev   }
693a06653b4SBarry Smith   PetscFunctionReturn(0);
694a06653b4SBarry Smith }
695a06653b4SBarry Smith 
696a06653b4SBarry Smith #undef __FUNCT__
697a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM"
698a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc)
699a06653b4SBarry Smith {
700a06653b4SBarry Smith   PC_GASM        *osm = (PC_GASM*)pc->data;
701a06653b4SBarry Smith   PetscErrorCode ierr;
702a06653b4SBarry Smith   PetscInt       i;
703a06653b4SBarry Smith 
704a06653b4SBarry Smith   PetscFunctionBegin;
705a06653b4SBarry Smith   ierr = PCReset_GASM(pc);CHKERRQ(ierr);
706135757f6SDmitry Karpeev 
707135757f6SDmitry Karpeev   /* PCReset will not destroy subdomains, if user_subdomains is true. */
708135757f6SDmitry Karpeev   ierr = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr);
709135757f6SDmitry Karpeev 
710a06653b4SBarry Smith   if (osm->ksp) {
711a06653b4SBarry Smith     for (i=0; i<osm->n; i++) {
7126bf464f9SBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
713a06653b4SBarry Smith     }
714a06653b4SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
715a06653b4SBarry Smith   }
716a06653b4SBarry Smith   ierr = PetscFree(osm->x);CHKERRQ(ierr);
717a06653b4SBarry Smith   ierr = PetscFree(osm->y);CHKERRQ(ierr);
718c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
719b1a0a8a3SJed Brown   PetscFunctionReturn(0);
720b1a0a8a3SJed Brown }
721b1a0a8a3SJed Brown 
722b1a0a8a3SJed Brown #undef __FUNCT__
723f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM"
7248c34d3f5SBarry Smith static PetscErrorCode PCSetFromOptions_GASM(PetscOptions *PetscOptionsObject,PC pc)
725a6dfd86eSKarl Rupp {
726f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
727b1a0a8a3SJed Brown   PetscErrorCode ierr;
728b1a0a8a3SJed Brown   PetscInt       blocks,ovl;
729b1a0a8a3SJed Brown   PetscBool      symset,flg;
730f746d493SDmitry Karpeev   PCGASMType     gasmtype;
731b1a0a8a3SJed Brown 
732b1a0a8a3SJed Brown   PetscFunctionBegin;
733b1a0a8a3SJed Brown   /* set the type to symmetric if matrix is symmetric */
734b1a0a8a3SJed Brown   if (!osm->type_set && pc->pmat) {
735b1a0a8a3SJed Brown     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
7362fa5cd67SKarl Rupp     if (symset && flg) osm->type = PC_GASM_BASIC;
737b1a0a8a3SJed Brown   }
738e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Generalized additive Schwarz options");CHKERRQ(ierr);
7398f3b4b4dSDmitry 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);
7408f3b4b4dSDmitry Karpeev   ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg);CHKERRQ(ierr);
74165db9045SDmitry Karpeev   if (flg) {
7428f3b4b4dSDmitry Karpeev     ierr = PCGASMSetTotalSubdomains(pc,blocks);CHKERRQ(ierr);
74365db9045SDmitry Karpeev   }
74406b43e7eSDmitry Karpeev   ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
74565db9045SDmitry Karpeev   if (flg) {
74665db9045SDmitry Karpeev     ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr);
747d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
74865db9045SDmitry Karpeev   }
749b1a0a8a3SJed Brown   flg  = PETSC_FALSE;
750f746d493SDmitry Karpeev   ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr);
751f746d493SDmitry Karpeev   if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr);}
752b1a0a8a3SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
753b1a0a8a3SJed Brown   PetscFunctionReturn(0);
754b1a0a8a3SJed Brown }
755b1a0a8a3SJed Brown 
756b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/
757b1a0a8a3SJed Brown 
758b1a0a8a3SJed Brown #undef __FUNCT__
7598f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains"
7608f3b4b4dSDmitry Karpeev /*@
7618f3b4b4dSDmitry Karpeev     PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the
7628f3b4b4dSDmitry Karpeev                                communicator.
7638f3b4b4dSDmitry Karpeev     Logically collective on pc
7648f3b4b4dSDmitry Karpeev 
7658f3b4b4dSDmitry Karpeev     Input Parameters:
7668f3b4b4dSDmitry Karpeev +   pc  - the preconditioner
7678f3b4b4dSDmitry Karpeev -   N   - total number of subdomains
7688f3b4b4dSDmitry Karpeev 
7698f3b4b4dSDmitry Karpeev 
7708f3b4b4dSDmitry Karpeev     Level: beginner
7718f3b4b4dSDmitry Karpeev 
7728f3b4b4dSDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
7738f3b4b4dSDmitry Karpeev 
7748f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap()
7758f3b4b4dSDmitry Karpeev           PCGASMCreateSubdomains2D()
7768f3b4b4dSDmitry Karpeev @*/
7778f3b4b4dSDmitry Karpeev PetscErrorCode  PCGASMSetTotalSubdomains(PC pc,PetscInt N)
7788f3b4b4dSDmitry Karpeev {
7798f3b4b4dSDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
7808f3b4b4dSDmitry Karpeev   PetscMPIInt    size,rank;
7818f3b4b4dSDmitry Karpeev   PetscErrorCode ierr;
7828f3b4b4dSDmitry Karpeev 
7838f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
7848f3b4b4dSDmitry 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);
7858f3b4b4dSDmitry Karpeev   if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before calling PCSetUp().");
7868f3b4b4dSDmitry Karpeev 
7872c112581SDmitry Karpeev   ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr);
7888f3b4b4dSDmitry Karpeev   osm->ois = osm->iis = NULL;
7898f3b4b4dSDmitry Karpeev 
7908f3b4b4dSDmitry Karpeev   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
7918f3b4b4dSDmitry Karpeev   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
7928f3b4b4dSDmitry Karpeev   osm->N    = N;
7938f3b4b4dSDmitry Karpeev   osm->n    = PETSC_DETERMINE;
7948f3b4b4dSDmitry Karpeev   osm->nmax = PETSC_DETERMINE;
7958f3b4b4dSDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
7968f3b4b4dSDmitry Karpeev   PetscFunctionReturn(0);
7978f3b4b4dSDmitry Karpeev }
7988f3b4b4dSDmitry Karpeev 
7998f3b4b4dSDmitry Karpeev #undef __FUNCT__
80006b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM"
801b541e6a4SDmitry Karpeev static PetscErrorCode  PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[])
802b1a0a8a3SJed Brown {
803f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
804b1a0a8a3SJed Brown   PetscErrorCode ierr;
805b1a0a8a3SJed Brown   PetscInt       i;
806b1a0a8a3SJed Brown 
807b1a0a8a3SJed Brown   PetscFunctionBegin;
8088f3b4b4dSDmitry Karpeev   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, got n = %D",n);
8098f3b4b4dSDmitry Karpeev   if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp().");
810b1a0a8a3SJed Brown 
8112c112581SDmitry Karpeev   ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr);
8128f3b4b4dSDmitry Karpeev   osm->iis  = osm->ois = NULL;
8138f3b4b4dSDmitry Karpeev   osm->n    = n;
8148f3b4b4dSDmitry Karpeev   osm->N    = PETSC_DETERMINE;
8158f3b4b4dSDmitry Karpeev   osm->nmax = PETSC_DETERMINE;
816a35b7b57SDmitry Karpeev   if (ois) {
817785e854fSJed Brown     ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr);
8188f3b4b4dSDmitry Karpeev     for (i=0; i<n; i++) {
8198f3b4b4dSDmitry Karpeev       ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);
8208f3b4b4dSDmitry Karpeev       osm->ois[i] = ois[i];
8218f3b4b4dSDmitry Karpeev     }
8228f3b4b4dSDmitry Karpeev     /*
8238f3b4b4dSDmitry Karpeev        Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(),
8248f3b4b4dSDmitry Karpeev        it will be ignored.  To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1.
8258f3b4b4dSDmitry Karpeev     */
826b1a0a8a3SJed Brown     osm->overlap = -1;
827a35b7b57SDmitry Karpeev     if (!iis) {
828785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr);
829a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {
8308f3b4b4dSDmitry Karpeev 	for (i=0; i<n; i++) {
8310e8f3646SDmitry Karpeev 	  ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);
832a35b7b57SDmitry Karpeev 	  osm->iis[i] = ois[i];
833a35b7b57SDmitry Karpeev 	}
834a35b7b57SDmitry Karpeev       }
835a35b7b57SDmitry Karpeev     }
8368f3b4b4dSDmitry Karpeev   }
837a35b7b57SDmitry Karpeev   if (iis) {
838785e854fSJed Brown     ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr);
8398f3b4b4dSDmitry Karpeev     for (i=0; i<n; i++) {
8408f3b4b4dSDmitry Karpeev       ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);
8418f3b4b4dSDmitry Karpeev       osm->iis[i] = iis[i];
8428f3b4b4dSDmitry Karpeev     }
843a35b7b57SDmitry Karpeev     if (!ois) {
844785e854fSJed Brown       ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr);
845a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {
8468f3b4b4dSDmitry Karpeev 	for (i=0; i<n; i++) {
847a35b7b57SDmitry Karpeev 	  ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);
848a35b7b57SDmitry Karpeev 	  osm->ois[i] = iis[i];
849a35b7b57SDmitry Karpeev 	}
8508f3b4b4dSDmitry Karpeev       }
8518f3b4b4dSDmitry Karpeev       /* If necessary, osm->ois[i] will be expanded using the requested overlap size in PCSetUp_GASM(). */
852a35b7b57SDmitry Karpeev     }
853a35b7b57SDmitry Karpeev   }
8548f3b4b4dSDmitry Karpeev   if (iis)  osm->user_subdomains = PETSC_TRUE;
855b1a0a8a3SJed Brown   PetscFunctionReturn(0);
856b1a0a8a3SJed Brown }
857b1a0a8a3SJed Brown 
858b1a0a8a3SJed Brown 
859b1a0a8a3SJed Brown #undef __FUNCT__
860f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM"
861f7a08781SBarry Smith static PetscErrorCode  PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
862b1a0a8a3SJed Brown {
863f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
864b1a0a8a3SJed Brown 
865b1a0a8a3SJed Brown   PetscFunctionBegin;
866ce94432eSBarry Smith   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
867ce94432eSBarry Smith   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
8682fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
869b1a0a8a3SJed Brown   PetscFunctionReturn(0);
870b1a0a8a3SJed Brown }
871b1a0a8a3SJed Brown 
872b1a0a8a3SJed Brown #undef __FUNCT__
873f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM"
874f7a08781SBarry Smith static PetscErrorCode  PCGASMSetType_GASM(PC pc,PCGASMType type)
875b1a0a8a3SJed Brown {
876f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
877b1a0a8a3SJed Brown 
878b1a0a8a3SJed Brown   PetscFunctionBegin;
879b1a0a8a3SJed Brown   osm->type     = type;
880b1a0a8a3SJed Brown   osm->type_set = PETSC_TRUE;
881b1a0a8a3SJed Brown   PetscFunctionReturn(0);
882b1a0a8a3SJed Brown }
883b1a0a8a3SJed Brown 
884b1a0a8a3SJed Brown #undef __FUNCT__
885f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM"
886f7a08781SBarry Smith static PetscErrorCode  PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort)
887b1a0a8a3SJed Brown {
888f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
889b1a0a8a3SJed Brown 
890b1a0a8a3SJed Brown   PetscFunctionBegin;
891b1a0a8a3SJed Brown   osm->sort_indices = doSort;
892b1a0a8a3SJed Brown   PetscFunctionReturn(0);
893b1a0a8a3SJed Brown }
894b1a0a8a3SJed Brown 
895b1a0a8a3SJed Brown #undef __FUNCT__
896f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM"
89744fe18e5SDmitry Karpeev /*
8988f3b4b4dSDmitry Karpeev    FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed.
89944fe18e5SDmitry Karpeev         In particular, it would upset the global subdomain number calculation.
90044fe18e5SDmitry Karpeev */
901f7a08781SBarry Smith static PetscErrorCode  PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
902b1a0a8a3SJed Brown {
903f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
904b1a0a8a3SJed Brown   PetscErrorCode ierr;
905b1a0a8a3SJed Brown 
906b1a0a8a3SJed Brown   PetscFunctionBegin;
907ce94432eSBarry 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");
908b1a0a8a3SJed Brown 
9092fa5cd67SKarl Rupp   if (n) *n = osm->n;
910ab3e923fSDmitry Karpeev   if (first) {
911ce94432eSBarry Smith     ierr    = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
912ab3e923fSDmitry Karpeev     *first -= osm->n;
913b1a0a8a3SJed Brown   }
914b1a0a8a3SJed Brown   if (ksp) {
915b1a0a8a3SJed Brown     /* Assume that local solves are now different; not necessarily
91606b43e7eSDmitry Karpeev        true, though!  This flag is used only for PCView_GASM() */
917b1a0a8a3SJed Brown     *ksp                        = osm->ksp;
9186a4f0f73SDmitry Karpeev     osm->same_subdomain_solvers = PETSC_FALSE;
919b1a0a8a3SJed Brown   }
920b1a0a8a3SJed Brown   PetscFunctionReturn(0);
921ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */
922b1a0a8a3SJed Brown 
923b1a0a8a3SJed Brown #undef __FUNCT__
92406b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains"
925b1a0a8a3SJed Brown /*@C
92606b43e7eSDmitry Karpeev     PCGASMSetSubdomains - Sets the subdomains for this processor
92706b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
928b1a0a8a3SJed Brown 
9298f3b4b4dSDmitry Karpeev     Collective on pc
930b1a0a8a3SJed Brown 
931b1a0a8a3SJed Brown     Input Parameters:
9328f3b4b4dSDmitry Karpeev +   pc  - the preconditioner object
93306b43e7eSDmitry Karpeev .   n   - the number of subdomains for this processor
9348f3b4b4dSDmitry Karpeev .   iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains)
9358f3b4b4dSDmitry 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)
936b1a0a8a3SJed Brown 
937b1a0a8a3SJed Brown     Notes:
93806b43e7eSDmitry Karpeev     The IS indices use the parallel, global numbering of the vector entries.
9396a4f0f73SDmitry Karpeev     Inner subdomains are those where the correction is applied.
9406a4f0f73SDmitry Karpeev     Outer subdomains are those where the residual necessary to obtain the
9416a4f0f73SDmitry Karpeev     corrections is obtained (see PCGASMType for the use of inner/outer subdomains).
9426a4f0f73SDmitry Karpeev     Both inner and outer subdomains can extend over several processors.
9436a4f0f73SDmitry Karpeev     This processor's portion of a subdomain is known as a local subdomain.
9446a4f0f73SDmitry Karpeev 
9456a4f0f73SDmitry Karpeev     By default the GASM preconditioner uses 1 (local) subdomain per processor.
9466a4f0f73SDmitry Karpeev 
947b1a0a8a3SJed Brown 
948b1a0a8a3SJed Brown     Level: advanced
949b1a0a8a3SJed Brown 
95006b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz
951b1a0a8a3SJed Brown 
9528f3b4b4dSDmitry Karpeev .seealso: PCGASMSetNumSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
95306b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
954b1a0a8a3SJed Brown @*/
9556a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[])
956b1a0a8a3SJed Brown {
9578f3b4b4dSDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
958b1a0a8a3SJed Brown   PetscErrorCode ierr;
959b1a0a8a3SJed Brown 
960b1a0a8a3SJed Brown   PetscFunctionBegin;
961b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9626a4f0f73SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr);
9638f3b4b4dSDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
964b1a0a8a3SJed Brown   PetscFunctionReturn(0);
965b1a0a8a3SJed Brown }
966b1a0a8a3SJed Brown 
967b1a0a8a3SJed Brown 
968b1a0a8a3SJed Brown #undef __FUNCT__
969f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap"
970b1a0a8a3SJed Brown /*@
971f746d493SDmitry Karpeev     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
972b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
9738f3b4b4dSDmitry Karpeev     pc communicator must call this routine.
974b1a0a8a3SJed Brown 
9758f3b4b4dSDmitry Karpeev     Logically Collective on pc
976b1a0a8a3SJed Brown 
977b1a0a8a3SJed Brown     Input Parameters:
978b1a0a8a3SJed Brown +   pc  - the preconditioner context
9798f3b4b4dSDmitry Karpeev -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0)
980b1a0a8a3SJed Brown 
981b1a0a8a3SJed Brown     Options Database Key:
98206b43e7eSDmitry Karpeev .   -pc_gasm_overlap <overlap> - Sets overlap
983b1a0a8a3SJed Brown 
984b1a0a8a3SJed Brown     Notes:
98506b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.  To use
9868f3b4b4dSDmitry Karpeev     multiple subdomain per perocessor or "straddling" subdomains that intersect
9878f3b4b4dSDmitry Karpeev     multiple processors use PCGASMSetSubdomains() (or option -pc_gasm_total_subdomains <n>).
988b1a0a8a3SJed Brown 
9898f3b4b4dSDmitry Karpeev     The overlap defaults to 0, so if one desires that no additional
990b1a0a8a3SJed Brown     overlap be computed beyond what may have been set with a call to
9918f3b4b4dSDmitry Karpeev     PCGASMSetSubdomains(), then ovl must be set to be 0.  In particular, if one does
9928f3b4b4dSDmitry Karpeev     not explicitly set the subdomains in application code, then all overlap would be computed
993f746d493SDmitry Karpeev     internally by PETSc, and using an overlap of 0 would result in an GASM
994b1a0a8a3SJed Brown     variant that is equivalent to the block Jacobi preconditioner.
995b1a0a8a3SJed Brown 
996b1a0a8a3SJed Brown     Note that one can define initial index sets with any overlap via
99706b43e7eSDmitry Karpeev     PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows
99806b43e7eSDmitry Karpeev     PETSc to extend that overlap further, if desired.
999b1a0a8a3SJed Brown 
1000b1a0a8a3SJed Brown     Level: intermediate
1001b1a0a8a3SJed Brown 
1002f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap
1003b1a0a8a3SJed Brown 
10048f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
100506b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
1006b1a0a8a3SJed Brown @*/
10077087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap(PC pc,PetscInt ovl)
1008b1a0a8a3SJed Brown {
1009b1a0a8a3SJed Brown   PetscErrorCode ierr;
10108f3b4b4dSDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
1011b1a0a8a3SJed Brown 
1012b1a0a8a3SJed Brown   PetscFunctionBegin;
1013b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1014b1a0a8a3SJed Brown   PetscValidLogicalCollectiveInt(pc,ovl,2);
1015f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
10168f3b4b4dSDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
1017b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1018b1a0a8a3SJed Brown }
1019b1a0a8a3SJed Brown 
1020b1a0a8a3SJed Brown #undef __FUNCT__
1021f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType"
1022b1a0a8a3SJed Brown /*@
1023f746d493SDmitry Karpeev     PCGASMSetType - Sets the type of restriction and interpolation used
1024b1a0a8a3SJed Brown     for local problems in the additive Schwarz method.
1025b1a0a8a3SJed Brown 
1026b1a0a8a3SJed Brown     Logically Collective on PC
1027b1a0a8a3SJed Brown 
1028b1a0a8a3SJed Brown     Input Parameters:
1029b1a0a8a3SJed Brown +   pc  - the preconditioner context
1030f746d493SDmitry Karpeev -   type - variant of GASM, one of
1031b1a0a8a3SJed Brown .vb
1032f746d493SDmitry Karpeev       PC_GASM_BASIC       - full interpolation and restriction
1033f746d493SDmitry Karpeev       PC_GASM_RESTRICT    - full restriction, local processor interpolation
1034f746d493SDmitry Karpeev       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
1035f746d493SDmitry Karpeev       PC_GASM_NONE        - local processor restriction and interpolation
1036b1a0a8a3SJed Brown .ve
1037b1a0a8a3SJed Brown 
1038b1a0a8a3SJed Brown     Options Database Key:
1039f746d493SDmitry Karpeev .   -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1040b1a0a8a3SJed Brown 
1041b1a0a8a3SJed Brown     Level: intermediate
1042b1a0a8a3SJed Brown 
1043f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1044b1a0a8a3SJed Brown 
10458f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1046f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1047b1a0a8a3SJed Brown @*/
10487087cfbeSBarry Smith PetscErrorCode  PCGASMSetType(PC pc,PCGASMType type)
1049b1a0a8a3SJed Brown {
1050b1a0a8a3SJed Brown   PetscErrorCode ierr;
1051b1a0a8a3SJed Brown 
1052b1a0a8a3SJed Brown   PetscFunctionBegin;
1053b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1054b1a0a8a3SJed Brown   PetscValidLogicalCollectiveEnum(pc,type,2);
1055f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr);
1056b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1057b1a0a8a3SJed Brown }
1058b1a0a8a3SJed Brown 
1059b1a0a8a3SJed Brown #undef __FUNCT__
1060f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices"
1061b1a0a8a3SJed Brown /*@
1062f746d493SDmitry Karpeev     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
1063b1a0a8a3SJed Brown 
1064b1a0a8a3SJed Brown     Logically Collective on PC
1065b1a0a8a3SJed Brown 
1066b1a0a8a3SJed Brown     Input Parameters:
1067b1a0a8a3SJed Brown +   pc  - the preconditioner context
1068b1a0a8a3SJed Brown -   doSort - sort the subdomain indices
1069b1a0a8a3SJed Brown 
1070b1a0a8a3SJed Brown     Level: intermediate
1071b1a0a8a3SJed Brown 
1072f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1073b1a0a8a3SJed Brown 
10748f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1075f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1076b1a0a8a3SJed Brown @*/
10777087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices(PC pc,PetscBool doSort)
1078b1a0a8a3SJed Brown {
1079b1a0a8a3SJed Brown   PetscErrorCode ierr;
1080b1a0a8a3SJed Brown 
1081b1a0a8a3SJed Brown   PetscFunctionBegin;
1082b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1083b1a0a8a3SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
1084f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
1085b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1086b1a0a8a3SJed Brown }
1087b1a0a8a3SJed Brown 
1088b1a0a8a3SJed Brown #undef __FUNCT__
1089f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP"
1090b1a0a8a3SJed Brown /*@C
1091f746d493SDmitry Karpeev    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
1092b1a0a8a3SJed Brown    this processor.
1093b1a0a8a3SJed Brown 
1094b1a0a8a3SJed Brown    Collective on PC iff first_local is requested
1095b1a0a8a3SJed Brown 
1096b1a0a8a3SJed Brown    Input Parameter:
1097b1a0a8a3SJed Brown .  pc - the preconditioner context
1098b1a0a8a3SJed Brown 
1099b1a0a8a3SJed Brown    Output Parameters:
11000298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
11010298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
11020298fd71SBarry Smith                  all processors must request or all must pass NULL
1103b1a0a8a3SJed Brown -  ksp - the array of KSP contexts
1104b1a0a8a3SJed Brown 
1105b1a0a8a3SJed Brown    Note:
1106f746d493SDmitry Karpeev    After PCGASMGetSubKSP() the array of KSPes is not to be freed
1107b1a0a8a3SJed Brown 
1108b1a0a8a3SJed Brown    Currently for some matrix implementations only 1 block per processor
1109b1a0a8a3SJed Brown    is supported.
1110b1a0a8a3SJed Brown 
1111f746d493SDmitry Karpeev    You must call KSPSetUp() before calling PCGASMGetSubKSP().
1112b1a0a8a3SJed Brown 
1113b1a0a8a3SJed Brown    Level: advanced
1114b1a0a8a3SJed Brown 
1115f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context
1116b1a0a8a3SJed Brown 
11178f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(),
1118f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(),
1119b1a0a8a3SJed Brown @*/
11207087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1121b1a0a8a3SJed Brown {
1122b1a0a8a3SJed Brown   PetscErrorCode ierr;
1123b1a0a8a3SJed Brown 
1124b1a0a8a3SJed Brown   PetscFunctionBegin;
1125b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1126f746d493SDmitry Karpeev   ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
1127b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1128b1a0a8a3SJed Brown }
1129b1a0a8a3SJed Brown 
1130b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/
1131b1a0a8a3SJed Brown /*MC
1132f746d493SDmitry Karpeev    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1133b1a0a8a3SJed Brown            its own KSP object.
1134b1a0a8a3SJed Brown 
1135b1a0a8a3SJed Brown    Options Database Keys:
11368f3b4b4dSDmitry Karpeev +  -pc_gasm_total_subdomains <n>  - Sets total number of local subdomains to be distributed among processors
113706b43e7eSDmitry Karpeev .  -pc_gasm_view_subdomains       - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view
113806b43e7eSDmitry Karpeev .  -pc_gasm_print_subdomains      - activates the printing of subdomain indices in PCSetUp()
113906b43e7eSDmitry Karpeev .  -pc_gasm_overlap <ovl>         - Sets overlap by which to (automatically) extend local subdomains
1140f746d493SDmitry Karpeev -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1141b1a0a8a3SJed Brown 
1142b1a0a8a3SJed Brown      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
1143f746d493SDmitry Karpeev       will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
1144f746d493SDmitry Karpeev       -pc_gasm_type basic to use the standard GASM.
1145b1a0a8a3SJed Brown 
1146b1a0a8a3SJed Brown    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
1147b1a0a8a3SJed Brown      than one processor. Defaults to one block per processor.
1148b1a0a8a3SJed Brown 
1149b1a0a8a3SJed Brown      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1150b1a0a8a3SJed Brown         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
1151b1a0a8a3SJed Brown 
1152f746d493SDmitry Karpeev      To set the options on the solvers separate for each block call PCGASMGetSubKSP()
1153b1a0a8a3SJed Brown          and set the options directly on the resulting KSP object (you can access its PC
1154b1a0a8a3SJed Brown          with KSPGetPC())
1155b1a0a8a3SJed Brown 
1156b1a0a8a3SJed Brown 
1157b1a0a8a3SJed Brown    Level: beginner
1158b1a0a8a3SJed Brown 
1159b1a0a8a3SJed Brown    Concepts: additive Schwarz method
1160b1a0a8a3SJed Brown 
1161b1a0a8a3SJed Brown     References:
1162b1a0a8a3SJed Brown     An additive variant of the Schwarz alternating method for the case of many subregions
1163b1a0a8a3SJed Brown     M Dryja, OB Widlund - Courant Institute, New York University Technical report
1164b1a0a8a3SJed Brown 
1165b1a0a8a3SJed Brown     Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1166b1a0a8a3SJed Brown     Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X.
1167b1a0a8a3SJed Brown 
1168b1a0a8a3SJed Brown .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
116949517cdeSBarry Smith            PCBJACOBI,  PCGASMGetSubKSP(), PCGASMSetSubdomains(),
11708f3b4b4dSDmitry Karpeev            PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType()
1171b1a0a8a3SJed Brown 
1172b1a0a8a3SJed Brown M*/
1173b1a0a8a3SJed Brown 
1174b1a0a8a3SJed Brown #undef __FUNCT__
1175f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM"
11768cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc)
1177b1a0a8a3SJed Brown {
1178b1a0a8a3SJed Brown   PetscErrorCode ierr;
1179f746d493SDmitry Karpeev   PC_GASM        *osm;
1180b1a0a8a3SJed Brown 
1181b1a0a8a3SJed Brown   PetscFunctionBegin;
1182b00a9115SJed Brown   ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr);
11832fa5cd67SKarl Rupp 
11848f3b4b4dSDmitry Karpeev   osm->N                      = PETSC_DETERMINE;
118506b43e7eSDmitry Karpeev   osm->n                      = PETSC_DECIDE;
11868f3b4b4dSDmitry Karpeev   osm->nmax                   = PETSC_DETERMINE;
11878f3b4b4dSDmitry Karpeev   osm->overlap                = 0;
1188b1a0a8a3SJed Brown   osm->ksp                    = 0;
11896a4f0f73SDmitry Karpeev   osm->gorestriction          = 0;
11906a4f0f73SDmitry Karpeev   osm->girestriction          = 0;
1191ab3e923fSDmitry Karpeev   osm->gx                     = 0;
1192ab3e923fSDmitry Karpeev   osm->gy                     = 0;
1193b1a0a8a3SJed Brown   osm->x                      = 0;
1194b1a0a8a3SJed Brown   osm->y                      = 0;
11956a4f0f73SDmitry Karpeev   osm->ois                    = 0;
11966a4f0f73SDmitry Karpeev   osm->iis                    = 0;
1197b1a0a8a3SJed Brown   osm->pmat                   = 0;
1198f746d493SDmitry Karpeev   osm->type                   = PC_GASM_RESTRICT;
11996a4f0f73SDmitry Karpeev   osm->same_subdomain_solvers = PETSC_TRUE;
1200b1a0a8a3SJed Brown   osm->sort_indices           = PETSC_TRUE;
1201d709ea83SDmitry Karpeev   osm->dm_subdomains          = PETSC_FALSE;
1202b1a0a8a3SJed Brown 
1203b1a0a8a3SJed Brown   pc->data                 = (void*)osm;
1204f746d493SDmitry Karpeev   pc->ops->apply           = PCApply_GASM;
1205f746d493SDmitry Karpeev   pc->ops->applytranspose  = PCApplyTranspose_GASM;
1206f746d493SDmitry Karpeev   pc->ops->setup           = PCSetUp_GASM;
1207a06653b4SBarry Smith   pc->ops->reset           = PCReset_GASM;
1208f746d493SDmitry Karpeev   pc->ops->destroy         = PCDestroy_GASM;
1209f746d493SDmitry Karpeev   pc->ops->setfromoptions  = PCSetFromOptions_GASM;
1210f746d493SDmitry Karpeev   pc->ops->setuponblocks   = PCSetUpOnBlocks_GASM;
1211f746d493SDmitry Karpeev   pc->ops->view            = PCView_GASM;
1212b1a0a8a3SJed Brown   pc->ops->applyrichardson = 0;
1213b1a0a8a3SJed Brown 
1214bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM);CHKERRQ(ierr);
1215bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM);CHKERRQ(ierr);
1216bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM);CHKERRQ(ierr);
1217bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM);CHKERRQ(ierr);
1218bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM);CHKERRQ(ierr);
1219b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1220b1a0a8a3SJed Brown }
1221b1a0a8a3SJed Brown 
1222b1a0a8a3SJed Brown 
1223b1a0a8a3SJed Brown #undef __FUNCT__
122406b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains"
12258f3b4b4dSDmitry Karpeev PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[])
1226b1a0a8a3SJed Brown {
1227b1a0a8a3SJed Brown   MatPartitioning mpart;
1228b1a0a8a3SJed Brown   const char      *prefix;
122911bd1e4dSLisandro Dalcin   PetscErrorCode  (*f)(Mat,MatReuse,Mat*);
1230b1a0a8a3SJed Brown   PetscMPIInt     size;
1231b1a0a8a3SJed Brown   PetscInt        i,j,rstart,rend,bs;
123211bd1e4dSLisandro Dalcin   PetscBool       isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
12330298fd71SBarry Smith   Mat             Ad     = NULL, adj;
1234b1a0a8a3SJed Brown   IS              ispart,isnumb,*is;
1235b1a0a8a3SJed Brown   PetscErrorCode  ierr;
1236b1a0a8a3SJed Brown 
1237b1a0a8a3SJed Brown   PetscFunctionBegin;
12388f3b4b4dSDmitry Karpeev   if (nloc < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local subdomains must > 0, got nloc = %D",nloc);
1239b1a0a8a3SJed Brown 
1240b1a0a8a3SJed Brown   /* Get prefix, row distribution, and block size */
1241b1a0a8a3SJed Brown   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1242b1a0a8a3SJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1243b1a0a8a3SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1244b1a0a8a3SJed 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);
1245b1a0a8a3SJed Brown 
1246b1a0a8a3SJed Brown   /* Get diagonal block from matrix if possible */
1247ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
12480005d66cSJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",&f);CHKERRQ(ierr);
1249b1a0a8a3SJed Brown   if (f) {
125011bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
1251b1a0a8a3SJed Brown   } else if (size == 1) {
125211bd1e4dSLisandro Dalcin     Ad = A;
1253b1a0a8a3SJed Brown   }
1254b1a0a8a3SJed Brown   if (Ad) {
1255251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1256251f4c67SDmitry Karpeev     if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
1257b1a0a8a3SJed Brown   }
12588f3b4b4dSDmitry Karpeev   if (Ad && nloc > 1) {
1259b1a0a8a3SJed Brown     PetscBool  match,done;
1260b1a0a8a3SJed Brown     /* Try to setup a good matrix partitioning if available */
1261b1a0a8a3SJed Brown     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
1262b1a0a8a3SJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
1263b1a0a8a3SJed Brown     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1264251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
1265b1a0a8a3SJed Brown     if (!match) {
1266251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
1267b1a0a8a3SJed Brown     }
1268b1a0a8a3SJed Brown     if (!match) { /* assume a "good" partitioner is available */
12691a83f524SJed Brown       PetscInt       na;
12701a83f524SJed Brown       const PetscInt *ia,*ja;
1271b1a0a8a3SJed Brown       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1272b1a0a8a3SJed Brown       if (done) {
1273b1a0a8a3SJed Brown         /* Build adjacency matrix by hand. Unfortunately a call to
1274b1a0a8a3SJed Brown            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1275b1a0a8a3SJed Brown            remove the block-aij structure and we cannot expect
1276b1a0a8a3SJed Brown            MatPartitioning to split vertices as we need */
12771a83f524SJed Brown         PetscInt       i,j,len,nnz,cnt,*iia=0,*jja=0;
12781a83f524SJed Brown         const PetscInt *row;
1279b1a0a8a3SJed Brown         nnz = 0;
1280b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* count number of nonzeros */
1281b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1282b1a0a8a3SJed Brown           row = ja + ia[i];
1283b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1284b1a0a8a3SJed Brown             if (row[j] == i) { /* don't count diagonal */
1285b1a0a8a3SJed Brown               len--; break;
1286b1a0a8a3SJed Brown             }
1287b1a0a8a3SJed Brown           }
1288b1a0a8a3SJed Brown           nnz += len;
1289b1a0a8a3SJed Brown         }
1290854ce69bSBarry Smith         ierr   = PetscMalloc1(na+1,&iia);CHKERRQ(ierr);
1291854ce69bSBarry Smith         ierr   = PetscMalloc1(nnz,&jja);CHKERRQ(ierr);
1292b1a0a8a3SJed Brown         nnz    = 0;
1293b1a0a8a3SJed Brown         iia[0] = 0;
1294b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* fill adjacency */
1295b1a0a8a3SJed Brown           cnt = 0;
1296b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1297b1a0a8a3SJed Brown           row = ja + ia[i];
1298b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
12992fa5cd67SKarl Rupp             if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */
1300b1a0a8a3SJed Brown           }
1301b1a0a8a3SJed Brown           nnz += cnt;
1302b1a0a8a3SJed Brown           iia[i+1] = nnz;
1303b1a0a8a3SJed Brown         }
1304b1a0a8a3SJed Brown         /* Partitioning of the adjacency matrix */
13050298fd71SBarry Smith         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr);
1306b1a0a8a3SJed Brown         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
13078f3b4b4dSDmitry Karpeev         ierr      = MatPartitioningSetNParts(mpart,nloc);CHKERRQ(ierr);
1308b1a0a8a3SJed Brown         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
1309b1a0a8a3SJed Brown         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
13106bf464f9SBarry Smith         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
1311b1a0a8a3SJed Brown         foundpart = PETSC_TRUE;
1312b1a0a8a3SJed Brown       }
1313b1a0a8a3SJed Brown       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1314b1a0a8a3SJed Brown     }
13156bf464f9SBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1316b1a0a8a3SJed Brown   }
13178f3b4b4dSDmitry Karpeev   ierr = PetscMalloc1(nloc,&is);CHKERRQ(ierr);
1318b1a0a8a3SJed Brown   if (!foundpart) {
1319b1a0a8a3SJed Brown 
1320b1a0a8a3SJed Brown     /* Partitioning by contiguous chunks of rows */
1321b1a0a8a3SJed Brown 
1322b1a0a8a3SJed Brown     PetscInt mbs   = (rend-rstart)/bs;
1323b1a0a8a3SJed Brown     PetscInt start = rstart;
13248f3b4b4dSDmitry Karpeev     for (i=0; i<nloc; i++) {
13258f3b4b4dSDmitry Karpeev       PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs;
1326b1a0a8a3SJed Brown       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
1327b1a0a8a3SJed Brown       start += count;
1328b1a0a8a3SJed Brown     }
1329b1a0a8a3SJed Brown 
1330b1a0a8a3SJed Brown   } else {
1331b1a0a8a3SJed Brown 
1332b1a0a8a3SJed Brown     /* Partitioning by adjacency of diagonal block  */
1333b1a0a8a3SJed Brown 
1334b1a0a8a3SJed Brown     const PetscInt *numbering;
1335b1a0a8a3SJed Brown     PetscInt       *count,nidx,*indices,*newidx,start=0;
1336b1a0a8a3SJed Brown     /* Get node count in each partition */
13378f3b4b4dSDmitry Karpeev     ierr = PetscMalloc1(nloc,&count);CHKERRQ(ierr);
13388f3b4b4dSDmitry Karpeev     ierr = ISPartitioningCount(ispart,nloc,count);CHKERRQ(ierr);
1339b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
13408f3b4b4dSDmitry Karpeev       for (i=0; i<nloc; i++) count[i] *= bs;
1341b1a0a8a3SJed Brown     }
1342b1a0a8a3SJed Brown     /* Build indices from node numbering */
1343b1a0a8a3SJed Brown     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1344785e854fSJed Brown     ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr);
1345b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1346b1a0a8a3SJed Brown     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
1347b1a0a8a3SJed Brown     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
1348b1a0a8a3SJed Brown     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
1349b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1350785e854fSJed Brown       ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr);
13512fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
13522fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
13532fa5cd67SKarl Rupp       }
1354b1a0a8a3SJed Brown       ierr    = PetscFree(indices);CHKERRQ(ierr);
1355b1a0a8a3SJed Brown       nidx   *= bs;
1356b1a0a8a3SJed Brown       indices = newidx;
1357b1a0a8a3SJed Brown     }
1358b1a0a8a3SJed Brown     /* Shift to get global indices */
1359b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] += rstart;
1360b1a0a8a3SJed Brown 
1361b1a0a8a3SJed Brown     /* Build the index sets for each block */
13628f3b4b4dSDmitry Karpeev     for (i=0; i<nloc; i++) {
1363b1a0a8a3SJed Brown       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
1364b1a0a8a3SJed Brown       ierr   = ISSort(is[i]);CHKERRQ(ierr);
1365b1a0a8a3SJed Brown       start += count[i];
1366b1a0a8a3SJed Brown     }
1367b1a0a8a3SJed Brown 
1368302440fdSBarry Smith     ierr = PetscFree(count);CHKERRQ(ierr);
1369302440fdSBarry Smith     ierr = PetscFree(indices);CHKERRQ(ierr);
1370fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1371fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
1372b1a0a8a3SJed Brown   }
13736a4f0f73SDmitry Karpeev   *iis = is;
13748f3b4b4dSDmitry Karpeev   PetscFunctionReturn(0);
13758f3b4b4dSDmitry Karpeev }
13768f3b4b4dSDmitry Karpeev 
13778f3b4b4dSDmitry Karpeev #undef __FUNCT__
13788f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMCreateStraddlingSubdomains"
1379b541e6a4SDmitry Karpeev PETSC_INTERN PetscErrorCode  PCGASMCreateStraddlingSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[])
13808f3b4b4dSDmitry Karpeev {
13818f3b4b4dSDmitry Karpeev   PetscErrorCode  ierr;
13828f3b4b4dSDmitry Karpeev 
13838f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
13848f3b4b4dSDmitry Karpeev   ierr = MatSubdomainsCreateCoalesce(A,N,n,iis);CHKERRQ(ierr);
13858f3b4b4dSDmitry Karpeev   PetscFunctionReturn(0);
13868f3b4b4dSDmitry Karpeev }
13878f3b4b4dSDmitry Karpeev 
13888f3b4b4dSDmitry Karpeev 
13898f3b4b4dSDmitry Karpeev 
13908f3b4b4dSDmitry Karpeev #undef __FUNCT__
13918f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains"
13928f3b4b4dSDmitry Karpeev /*@C
13938f3b4b4dSDmitry Karpeev    PCGASMCreateSubdomains - Creates n index sets defining n nonoverlapping subdomains for the additive
13948f3b4b4dSDmitry Karpeev    Schwarz preconditioner for a any problem based on its matrix.
13958f3b4b4dSDmitry Karpeev 
13968f3b4b4dSDmitry Karpeev    Collective
13978f3b4b4dSDmitry Karpeev 
13988f3b4b4dSDmitry Karpeev    Input Parameters:
13998f3b4b4dSDmitry Karpeev +  A       - The global matrix operator
14008f3b4b4dSDmitry Karpeev -  N       - the number of global subdomains requested
14018f3b4b4dSDmitry Karpeev 
14028f3b4b4dSDmitry Karpeev    Output Parameters:
14038f3b4b4dSDmitry Karpeev +  n   - the number of subdomains created on this processor
14048f3b4b4dSDmitry Karpeev -  iis - the array of index sets defining the local inner subdomains (on which the correction is applied)
14058f3b4b4dSDmitry Karpeev 
14068f3b4b4dSDmitry Karpeev    Level: advanced
14078f3b4b4dSDmitry Karpeev 
14088f3b4b4dSDmitry Karpeev    Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor.
14098f3b4b4dSDmitry Karpeev          When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local.
14108f3b4b4dSDmitry Karpeev          The resulting subdomains can be use in PCGASMSetSubdomains(pc,n,iss,NULL).  The overlapping
14118f3b4b4dSDmitry Karpeev 	 outer subdomains will be automatically generated from these according to the requested amount of
14128f3b4b4dSDmitry Karpeev 	 overlap; this is currently supported only with local subdomains.
14138f3b4b4dSDmitry Karpeev 
14148f3b4b4dSDmitry Karpeev 
14158f3b4b4dSDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
14168f3b4b4dSDmitry Karpeev 
14178f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains()
14188f3b4b4dSDmitry Karpeev @*/
1419b541e6a4SDmitry Karpeev PetscErrorCode  PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[])
14208f3b4b4dSDmitry Karpeev {
14218f3b4b4dSDmitry Karpeev   PetscMPIInt     size;
14228f3b4b4dSDmitry Karpeev   PetscErrorCode  ierr;
14238f3b4b4dSDmitry Karpeev 
14248f3b4b4dSDmitry Karpeev   PetscFunctionBegin;
14258f3b4b4dSDmitry Karpeev   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
14268f3b4b4dSDmitry Karpeev   PetscValidPointer(iis,4);
14278f3b4b4dSDmitry Karpeev 
14288f3b4b4dSDmitry Karpeev   if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %D",N);
14298f3b4b4dSDmitry Karpeev   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
14308f3b4b4dSDmitry Karpeev   if (N >= size) {
14318f3b4b4dSDmitry Karpeev     *n = N/size + (N%size);
14328f3b4b4dSDmitry Karpeev     ierr = PCGASMCreateLocalSubdomains(A,*n,iis);CHKERRQ(ierr);
14336a4f0f73SDmitry Karpeev   } else {
14348f3b4b4dSDmitry Karpeev     ierr = PCGASMCreateStraddlingSubdomains(A,N,n,iis);CHKERRQ(ierr);
14356a4f0f73SDmitry Karpeev   }
1436b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1437b1a0a8a3SJed Brown }
1438b1a0a8a3SJed Brown 
1439b1a0a8a3SJed Brown #undef __FUNCT__
1440f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains"
1441b1a0a8a3SJed Brown /*@C
1442f746d493SDmitry Karpeev    PCGASMDestroySubdomains - Destroys the index sets created with
14438f3b4b4dSDmitry Karpeev    PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be
144406b43e7eSDmitry Karpeev    called after setting subdomains with PCGASMSetSubdomains().
1445b1a0a8a3SJed Brown 
1446b1a0a8a3SJed Brown    Collective
1447b1a0a8a3SJed Brown 
1448b1a0a8a3SJed Brown    Input Parameters:
1449b1a0a8a3SJed Brown +  n   - the number of index sets
14506a4f0f73SDmitry Karpeev .  iis - the array of inner subdomains,
14510298fd71SBarry Smith -  ois - the array of outer subdomains, can be NULL
1452b1a0a8a3SJed Brown 
14536a4f0f73SDmitry Karpeev    Level: intermediate
14546a4f0f73SDmitry Karpeev 
14556a4f0f73SDmitry Karpeev    Notes: this is merely a convenience subroutine that walks each list,
14562c112581SDmitry Karpeev    destroys each IS on the list, and then frees the list. At the end the
14572c112581SDmitry Karpeev    list pointers are set to NULL.
1458b1a0a8a3SJed Brown 
1459f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1460b1a0a8a3SJed Brown 
14618f3b4b4dSDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetSubdomains()
1462b1a0a8a3SJed Brown @*/
14632c112581SDmitry Karpeev PetscErrorCode  PCGASMDestroySubdomains(PetscInt n,IS **iis,IS **ois)
1464b1a0a8a3SJed Brown {
1465b1a0a8a3SJed Brown   PetscInt       i;
1466b1a0a8a3SJed Brown   PetscErrorCode ierr;
14675fd66863SKarl Rupp 
1468b1a0a8a3SJed Brown   PetscFunctionBegin;
1469a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
14706a4f0f73SDmitry Karpeev   if (ois) {
14712c112581SDmitry Karpeev     PetscValidPointer(ois,3);
14722c112581SDmitry Karpeev     if (*ois) {
14732c112581SDmitry Karpeev       PetscValidPointer(*ois,3);
1474a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {
14752c112581SDmitry Karpeev         ierr = ISDestroy(&(*ois)[i]);CHKERRQ(ierr);
1476a35b7b57SDmitry Karpeev       }
1477135757f6SDmitry Karpeev       ierr = PetscFree((*ois));CHKERRQ(ierr);
14782c112581SDmitry Karpeev     }
1479b1a0a8a3SJed Brown   }
1480*b9d0fdaaSFande Kong   if (iis) {
1481*b9d0fdaaSFande Kong     PetscValidPointer(iis,2);
1482*b9d0fdaaSFande Kong     if (*iis) {
1483*b9d0fdaaSFande Kong       PetscValidPointer(*iis,2);
1484*b9d0fdaaSFande Kong       for (i=0; i<n; i++) {
1485*b9d0fdaaSFande Kong         ierr = ISDestroy(&(*iis)[i]);CHKERRQ(ierr);
1486*b9d0fdaaSFande Kong       }
1487*b9d0fdaaSFande Kong       ierr = PetscFree((*iis));CHKERRQ(ierr);
1488*b9d0fdaaSFande Kong     }
1489*b9d0fdaaSFande Kong   }
1490b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1491b1a0a8a3SJed Brown }
1492b1a0a8a3SJed Brown 
1493af538404SDmitry Karpeev 
1494af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \
1495af538404SDmitry Karpeev   {                                                                                                       \
1496af538404SDmitry Karpeev     PetscInt first_row = first/M, last_row = last/M+1;                                                     \
1497af538404SDmitry Karpeev     /*                                                                                                    \
1498af538404SDmitry Karpeev      Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners       \
1499af538404SDmitry Karpeev      of the bounding box of the intersection of the subdomain with the local ownership range (local       \
1500af538404SDmitry Karpeev      subdomain).                                                                                          \
1501af538404SDmitry Karpeev      Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows       \
1502af538404SDmitry Karpeev      of the intersection.                                                                                 \
1503af538404SDmitry Karpeev     */                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                             \
1504af538404SDmitry Karpeev     /* ylow_loc is the grid row containing the first element of the local sumbdomain */                   \
1505eec7e3faSDmitry Karpeev     *ylow_loc = PetscMax(first_row,ylow);                                                                    \
1506af538404SDmitry Karpeev     /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1507eec7e3faSDmitry Karpeev     *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft;                                                                            \
1508af538404SDmitry Karpeev     /* yhigh_loc is the grid row above the last local subdomain element */                                                                    \
1509eec7e3faSDmitry Karpeev     *yhigh_loc = PetscMin(last_row,yhigh);                                                                                                     \
1510af538404SDmitry Karpeev     /* xright is the offset of the end of the  local subdomain within its grid row (might actually be outside the local subdomain) */         \
1511eec7e3faSDmitry Karpeev     *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright;                                                                          \
1512af538404SDmitry Karpeev     /* Now compute the size of the local subdomain n. */ \
1513c3518bceSDmitry Karpeev     *n = 0;                                               \
1514eec7e3faSDmitry Karpeev     if (*ylow_loc < *yhigh_loc) {                           \
1515af538404SDmitry Karpeev       PetscInt width = xright-xleft;                     \
1516eec7e3faSDmitry Karpeev       *n += width*(*yhigh_loc-*ylow_loc-1);                 \
1517eec7e3faSDmitry Karpeev       *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \
1518eec7e3faSDmitry Karpeev       *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \
1519af538404SDmitry Karpeev     } \
1520af538404SDmitry Karpeev   }
1521af538404SDmitry Karpeev 
1522c3518bceSDmitry Karpeev 
1523eec7e3faSDmitry Karpeev 
1524eec7e3faSDmitry Karpeev #undef __FUNCT__
1525f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D"
1526b1a0a8a3SJed Brown /*@
1527f746d493SDmitry Karpeev    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1528b1a0a8a3SJed Brown    preconditioner for a two-dimensional problem on a regular grid.
1529b1a0a8a3SJed Brown 
1530af538404SDmitry Karpeev    Collective
1531b1a0a8a3SJed Brown 
1532b1a0a8a3SJed Brown    Input Parameters:
153306b43e7eSDmitry Karpeev +  M, N               - the global number of grid points in the x and y directions
1534af538404SDmitry Karpeev .  Mdomains, Ndomains - the global number of subdomains in the x and y directions
1535b1a0a8a3SJed Brown .  dof                - degrees of freedom per node
1536b1a0a8a3SJed Brown -  overlap            - overlap in mesh lines
1537b1a0a8a3SJed Brown 
1538b1a0a8a3SJed Brown    Output Parameters:
1539af538404SDmitry Karpeev +  Nsub - the number of local subdomains created
15406a4f0f73SDmitry Karpeev .  iis  - array of index sets defining inner (nonoverlapping) subdomains
15416a4f0f73SDmitry Karpeev -  ois  - array of index sets defining outer (overlapping, if overlap > 0) subdomains
1542b1a0a8a3SJed Brown 
1543b1a0a8a3SJed Brown 
1544b1a0a8a3SJed Brown    Level: advanced
1545b1a0a8a3SJed Brown 
1546f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid
1547b1a0a8a3SJed Brown 
15488f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), PCGASMSetOverlap()
1549b1a0a8a3SJed Brown @*/
15506a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMCreateSubdomains2D(PC pc,PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap,PetscInt *nsub,IS **iis,IS **ois)
1551b1a0a8a3SJed Brown {
1552b1a0a8a3SJed Brown   PetscErrorCode ierr;
15532bbb417fSDmitry Karpeev   PetscMPIInt    size, rank;
15542bbb417fSDmitry Karpeev   PetscInt       i, j;
15552bbb417fSDmitry Karpeev   PetscInt       maxheight, maxwidth;
15562bbb417fSDmitry Karpeev   PetscInt       xstart, xleft, xright, xleft_loc, xright_loc;
15572bbb417fSDmitry Karpeev   PetscInt       ystart, ylow,  yhigh,  ylow_loc,  yhigh_loc;
1558eec7e3faSDmitry Karpeev   PetscInt       x[2][2], y[2][2], n[2];
15592bbb417fSDmitry Karpeev   PetscInt       first, last;
15602bbb417fSDmitry Karpeev   PetscInt       nidx, *idx;
15612bbb417fSDmitry Karpeev   PetscInt       ii,jj,s,q,d;
1562af538404SDmitry Karpeev   PetscInt       k,kk;
15632bbb417fSDmitry Karpeev   PetscMPIInt    color;
15642bbb417fSDmitry Karpeev   MPI_Comm       comm, subcomm;
15656a4f0f73SDmitry Karpeev   IS             **xis = 0, **is = ois, **is_local = iis;
1566d34fcf5fSBarry Smith 
1567b1a0a8a3SJed Brown   PetscFunctionBegin;
15682bbb417fSDmitry Karpeev   ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr);
15692bbb417fSDmitry Karpeev   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
15702bbb417fSDmitry Karpeev   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
15712bbb417fSDmitry Karpeev   ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr);
1572d34fcf5fSBarry 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) "
15732bbb417fSDmitry Karpeev                                       "does not respect the number of degrees of freedom per grid point %D", first, last, dof);
1574d34fcf5fSBarry Smith 
1575af538404SDmitry Karpeev   /* Determine the number of domains with nonzero intersections with the local ownership range. */
15762bbb417fSDmitry Karpeev   s      = 0;
1577af538404SDmitry Karpeev   ystart = 0;
1578af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1579af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1580af538404SDmitry 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);
1581eec7e3faSDmitry Karpeev     /* Vertical domain limits with an overlap. */
1582eec7e3faSDmitry Karpeev     ylow   = PetscMax(ystart - overlap,0);
1583eec7e3faSDmitry Karpeev     yhigh  = PetscMin(ystart + maxheight + overlap,N);
1584b1a0a8a3SJed Brown     xstart = 0;
1585af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1586af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1587af538404SDmitry 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);
1588eec7e3faSDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1589eec7e3faSDmitry Karpeev       xleft  = PetscMax(xstart - overlap,0);
1590eec7e3faSDmitry Karpeev       xright = PetscMin(xstart + maxwidth + overlap,M);
1591af538404SDmitry Karpeev       /*
1592af538404SDmitry Karpeev          Determine whether this subdomain intersects this processor's ownership range of pc->pmat.
1593af538404SDmitry Karpeev       */
1594c3518bceSDmitry Karpeev       PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
15952fa5cd67SKarl Rupp       if (nidx) ++s;
1596af538404SDmitry Karpeev       xstart += maxwidth;
1597af538404SDmitry Karpeev     } /* for (i = 0; i < Mdomains; ++i) */
1598af538404SDmitry Karpeev     ystart += maxheight;
1599af538404SDmitry Karpeev   } /* for (j = 0; j < Ndomains; ++j) */
16002fa5cd67SKarl Rupp 
1601af538404SDmitry Karpeev   /* Now we can allocate the necessary number of ISs. */
1602af538404SDmitry Karpeev   *nsub  = s;
1603854ce69bSBarry Smith   ierr   = PetscMalloc1(*nsub,is);CHKERRQ(ierr);
1604854ce69bSBarry Smith   ierr   = PetscMalloc1(*nsub,is_local);CHKERRQ(ierr);
1605af538404SDmitry Karpeev   s      = 0;
1606af538404SDmitry Karpeev   ystart = 0;
1607af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1608af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1609af538404SDmitry 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);
1610af538404SDmitry Karpeev     /* Vertical domain limits with an overlap. */
1611eec7e3faSDmitry Karpeev     y[0][0] = PetscMax(ystart - overlap,0);
1612eec7e3faSDmitry Karpeev     y[0][1] = PetscMin(ystart + maxheight + overlap,N);
1613eec7e3faSDmitry Karpeev     /* Vertical domain limits without an overlap. */
1614eec7e3faSDmitry Karpeev     y[1][0] = ystart;
1615eec7e3faSDmitry Karpeev     y[1][1] = PetscMin(ystart + maxheight,N);
1616eec7e3faSDmitry Karpeev     xstart  = 0;
1617af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1618af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1619af538404SDmitry 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);
1620af538404SDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1621eec7e3faSDmitry Karpeev       x[0][0] = PetscMax(xstart - overlap,0);
1622eec7e3faSDmitry Karpeev       x[0][1] = PetscMin(xstart + maxwidth + overlap,M);
1623eec7e3faSDmitry Karpeev       /* Horizontal domain limits without an overlap. */
1624eec7e3faSDmitry Karpeev       x[1][0] = xstart;
1625eec7e3faSDmitry Karpeev       x[1][1] = PetscMin(xstart+maxwidth,M);
16262bbb417fSDmitry Karpeev       /*
16272bbb417fSDmitry Karpeev          Determine whether this domain intersects this processor's ownership range of pc->pmat.
16282bbb417fSDmitry Karpeev          Do this twice: first for the domains with overlaps, and once without.
16292bbb417fSDmitry Karpeev          During the first pass create the subcommunicators, and use them on the second pass as well.
16302bbb417fSDmitry Karpeev       */
16312bbb417fSDmitry Karpeev       for (q = 0; q < 2; ++q) {
16322bbb417fSDmitry Karpeev         /*
16332bbb417fSDmitry Karpeev           domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
16342bbb417fSDmitry Karpeev           according to whether the domain with an overlap or without is considered.
16352bbb417fSDmitry Karpeev         */
16362bbb417fSDmitry Karpeev         xleft = x[q][0]; xright = x[q][1];
16372bbb417fSDmitry Karpeev         ylow  = y[q][0]; yhigh  = y[q][1];
1638c3518bceSDmitry Karpeev         PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1639af538404SDmitry Karpeev         nidx *= dof;
1640eec7e3faSDmitry Karpeev         n[q]  = nidx;
16412bbb417fSDmitry Karpeev         /*
1642eec7e3faSDmitry Karpeev          Based on the counted number of indices in the local domain *with an overlap*,
1643af538404SDmitry Karpeev          construct a subcommunicator of all the processors supporting this domain.
1644eec7e3faSDmitry Karpeev          Observe that a domain with an overlap might have nontrivial local support,
1645eec7e3faSDmitry Karpeev          while the domain without an overlap might not.  Hence, the decision to participate
1646eec7e3faSDmitry Karpeev          in the subcommunicator must be based on the domain with an overlap.
16472bbb417fSDmitry Karpeev          */
16482bbb417fSDmitry Karpeev         if (q == 0) {
16492fa5cd67SKarl Rupp           if (nidx) color = 1;
16502fa5cd67SKarl Rupp           else color = MPI_UNDEFINED;
16512fa5cd67SKarl Rupp 
16522bbb417fSDmitry Karpeev           ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr);
16532bbb417fSDmitry Karpeev         }
1654af538404SDmitry Karpeev         /*
1655eec7e3faSDmitry Karpeev          Proceed only if the number of local indices *with an overlap* is nonzero.
1656af538404SDmitry Karpeev          */
1657eec7e3faSDmitry Karpeev         if (n[0]) {
16582fa5cd67SKarl Rupp           if (q == 0) xis = is;
1659af538404SDmitry Karpeev           if (q == 1) {
1660af538404SDmitry Karpeev             /*
1661eec7e3faSDmitry Karpeev              The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1662af538404SDmitry Karpeev              Moreover, if the overlap is zero, the two ISs are identical.
1663af538404SDmitry Karpeev              */
1664b1a0a8a3SJed Brown             if (overlap == 0) {
1665eec7e3faSDmitry Karpeev               (*is_local)[s] = (*is)[s];
16662bbb417fSDmitry Karpeev               ierr           = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr);
16672bbb417fSDmitry Karpeev               continue;
1668d34fcf5fSBarry Smith             } else {
16696a4f0f73SDmitry Karpeev               xis     = is_local;
1670eec7e3faSDmitry Karpeev               subcomm = ((PetscObject)(*is)[s])->comm;
16712bbb417fSDmitry Karpeev             }
1672af538404SDmitry Karpeev           } /* if (q == 1) */
16730298fd71SBarry Smith           idx  = NULL;
1674785e854fSJed Brown           ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr);
1675eec7e3faSDmitry Karpeev           if (nidx) {
16762bbb417fSDmitry Karpeev             k = 0;
16772bbb417fSDmitry Karpeev             for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1678af538404SDmitry Karpeev               PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft;
1679af538404SDmitry Karpeev               PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright;
1680af538404SDmitry Karpeev               kk = dof*(M*jj + x0);
16812bbb417fSDmitry Karpeev               for (ii=x0; ii<x1; ++ii) {
16822bbb417fSDmitry Karpeev                 for (d = 0; d < dof; ++d) {
16832bbb417fSDmitry Karpeev                   idx[k++] = kk++;
1684b1a0a8a3SJed Brown                 }
1685b1a0a8a3SJed Brown               }
1686b1a0a8a3SJed Brown             }
1687eec7e3faSDmitry Karpeev           }
16886a4f0f73SDmitry Karpeev           ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr);
1689eec7e3faSDmitry Karpeev         }/* if (n[0]) */
16902bbb417fSDmitry Karpeev       }/* for (q = 0; q < 2; ++q) */
16912fa5cd67SKarl Rupp       if (n[0]) ++s;
1692af538404SDmitry Karpeev       xstart += maxwidth;
1693af538404SDmitry Karpeev     } /* for (i = 0; i < Mdomains; ++i) */
16942bbb417fSDmitry Karpeev     ystart += maxheight;
1695af538404SDmitry Karpeev   } /* for (j = 0; j < Ndomains; ++j) */
1696b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1697b1a0a8a3SJed Brown }
1698b1a0a8a3SJed Brown 
1699b1a0a8a3SJed Brown #undef __FUNCT__
170006b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains"
1701b1a0a8a3SJed Brown /*@C
170206b43e7eSDmitry Karpeev     PCGASMGetSubdomains - Gets the subdomains supported on this processor
170306b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
1704b1a0a8a3SJed Brown 
1705b1a0a8a3SJed Brown     Not Collective
1706b1a0a8a3SJed Brown 
1707b1a0a8a3SJed Brown     Input Parameter:
1708b1a0a8a3SJed Brown .   pc - the preconditioner context
1709b1a0a8a3SJed Brown 
1710b1a0a8a3SJed Brown     Output Parameters:
1711b1a0a8a3SJed Brown +   n   - the number of subdomains for this processor (default value = 1)
17120298fd71SBarry Smith .   iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL)
17130298fd71SBarry Smith -   ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL)
1714b1a0a8a3SJed Brown 
1715b1a0a8a3SJed Brown 
1716b1a0a8a3SJed Brown     Notes:
17176a4f0f73SDmitry Karpeev     The user is responsible for destroying the ISs and freeing the returned arrays.
1718b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
1719b1a0a8a3SJed Brown 
1720b1a0a8a3SJed Brown     Level: advanced
1721b1a0a8a3SJed Brown 
172206b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz
1723b1a0a8a3SJed Brown 
17248f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), PCGASMCreateSubdomains2D(),
17258f3b4b4dSDmitry Karpeev           PCGASMSetSubdomains(), PCGASMGetSubmatrices()
1726b1a0a8a3SJed Brown @*/
17276a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[])
1728b1a0a8a3SJed Brown {
1729f746d493SDmitry Karpeev   PC_GASM        *osm;
1730b1a0a8a3SJed Brown   PetscErrorCode ierr;
1731b1a0a8a3SJed Brown   PetscBool      match;
17326a4f0f73SDmitry Karpeev   PetscInt       i;
17335fd66863SKarl Rupp 
1734b1a0a8a3SJed Brown   PetscFunctionBegin;
1735b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1736251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1737ce94432eSBarry Smith   if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1738f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1739ab3e923fSDmitry Karpeev   if (n) *n = osm->n;
17406a4f0f73SDmitry Karpeev   if (iis) {
1741785e854fSJed Brown     ierr = PetscMalloc1(osm->n, iis);CHKERRQ(ierr);
17426a4f0f73SDmitry Karpeev   }
17436a4f0f73SDmitry Karpeev   if (ois) {
1744785e854fSJed Brown     ierr = PetscMalloc1(osm->n, ois);CHKERRQ(ierr);
17456a4f0f73SDmitry Karpeev   }
17466a4f0f73SDmitry Karpeev   if (iis || ois) {
17476a4f0f73SDmitry Karpeev     for (i = 0; i < osm->n; ++i) {
17486a4f0f73SDmitry Karpeev       if (iis) (*iis)[i] = osm->iis[i];
17496a4f0f73SDmitry Karpeev       if (ois) (*ois)[i] = osm->ois[i];
17506a4f0f73SDmitry Karpeev     }
1751b1a0a8a3SJed Brown   }
1752b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1753b1a0a8a3SJed Brown }
1754b1a0a8a3SJed Brown 
1755b1a0a8a3SJed Brown #undef __FUNCT__
175606b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices"
1757b1a0a8a3SJed Brown /*@C
175806b43e7eSDmitry Karpeev     PCGASMGetSubmatrices - Gets the local submatrices (for this processor
1759b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1760b1a0a8a3SJed Brown 
1761b1a0a8a3SJed Brown     Not Collective
1762b1a0a8a3SJed Brown 
1763b1a0a8a3SJed Brown     Input Parameter:
1764b1a0a8a3SJed Brown .   pc - the preconditioner context
1765b1a0a8a3SJed Brown 
1766b1a0a8a3SJed Brown     Output Parameters:
1767b1a0a8a3SJed Brown +   n   - the number of matrices for this processor (default value = 1)
1768b1a0a8a3SJed Brown -   mat - the matrices
1769b1a0a8a3SJed Brown 
177006b43e7eSDmitry Karpeev     Notes: matrices returned by this routine have the same communicators as the index sets (IS)
17718f3b4b4dSDmitry Karpeev            used to define subdomains in PCGASMSetSubdomains()
1772b1a0a8a3SJed Brown     Level: advanced
1773b1a0a8a3SJed Brown 
1774f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi
1775b1a0a8a3SJed Brown 
17768f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(),
177706b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains()
1778b1a0a8a3SJed Brown @*/
177906b43e7eSDmitry Karpeev PetscErrorCode  PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1780b1a0a8a3SJed Brown {
1781f746d493SDmitry Karpeev   PC_GASM        *osm;
1782b1a0a8a3SJed Brown   PetscErrorCode ierr;
1783b1a0a8a3SJed Brown   PetscBool      match;
1784b1a0a8a3SJed Brown 
1785b1a0a8a3SJed Brown   PetscFunctionBegin;
1786b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1787b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1788b1a0a8a3SJed Brown   if (mat) PetscValidPointer(mat,3);
1789ce94432eSBarry Smith   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1790251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1791ce94432eSBarry Smith   if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1792f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1793ab3e923fSDmitry Karpeev   if (n) *n = osm->n;
1794b1a0a8a3SJed Brown   if (mat) *mat = osm->pmat;
1795b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1796b1a0a8a3SJed Brown }
1797d709ea83SDmitry Karpeev 
1798d709ea83SDmitry Karpeev #undef __FUNCT__
17998f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMSetUseDMSubdomains"
1800d709ea83SDmitry Karpeev /*@
18018f3b4b4dSDmitry Karpeev     PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1802d709ea83SDmitry Karpeev     Logically Collective
1803d709ea83SDmitry Karpeev 
1804d709ea83SDmitry Karpeev     Input Parameter:
1805d709ea83SDmitry Karpeev +   pc  - the preconditioner
1806d709ea83SDmitry Karpeev -   flg - boolean indicating whether to use subdomains defined by the DM
1807d709ea83SDmitry Karpeev 
1808d709ea83SDmitry Karpeev     Options Database Key:
18098f3b4b4dSDmitry Karpeev .   -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains
1810d709ea83SDmitry Karpeev 
1811d709ea83SDmitry Karpeev     Level: intermediate
1812d709ea83SDmitry Karpeev 
1813d709ea83SDmitry Karpeev     Notes:
18148f3b4b4dSDmitry Karpeev     PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(),
18158f3b4b4dSDmitry Karpeev     so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap()
18168f3b4b4dSDmitry Karpeev     automatically turns the latter off.
1817d709ea83SDmitry Karpeev 
1818d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1819d709ea83SDmitry Karpeev 
18208f3b4b4dSDmitry Karpeev .seealso: PCGASMGetUseDMSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap()
1821d709ea83SDmitry Karpeev           PCGASMCreateSubdomains2D()
1822d709ea83SDmitry Karpeev @*/
18238f3b4b4dSDmitry Karpeev PetscErrorCode  PCGASMSetUseDMSubdomains(PC pc,PetscBool flg)
1824d709ea83SDmitry Karpeev {
1825d709ea83SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
1826d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1827d709ea83SDmitry Karpeev   PetscBool      match;
1828d709ea83SDmitry Karpeev 
1829d709ea83SDmitry Karpeev   PetscFunctionBegin;
1830d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1831d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
1832d709ea83SDmitry Karpeev   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1833d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1834d709ea83SDmitry Karpeev   if (match) {
18358f3b4b4dSDmitry Karpeev     if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) {
1836d709ea83SDmitry Karpeev       osm->dm_subdomains = flg;
1837d709ea83SDmitry Karpeev     }
18388f3b4b4dSDmitry Karpeev   }
1839d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1840d709ea83SDmitry Karpeev }
1841d709ea83SDmitry Karpeev 
1842d709ea83SDmitry Karpeev #undef __FUNCT__
18438f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMGetUseDMSubdomains"
1844d709ea83SDmitry Karpeev /*@
18458f3b4b4dSDmitry Karpeev     PCGASMGetUseDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1846d709ea83SDmitry Karpeev     Not Collective
1847d709ea83SDmitry Karpeev 
1848d709ea83SDmitry Karpeev     Input Parameter:
1849d709ea83SDmitry Karpeev .   pc  - the preconditioner
1850d709ea83SDmitry Karpeev 
1851d709ea83SDmitry Karpeev     Output Parameter:
1852d709ea83SDmitry Karpeev .   flg - boolean indicating whether to use subdomains defined by the DM
1853d709ea83SDmitry Karpeev 
1854d709ea83SDmitry Karpeev     Level: intermediate
1855d709ea83SDmitry Karpeev 
1856d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1857d709ea83SDmitry Karpeev 
18588f3b4b4dSDmitry Karpeev .seealso: PCGASMSetUseDMSubdomains(), PCGASMSetOverlap()
1859d709ea83SDmitry Karpeev           PCGASMCreateSubdomains2D()
1860d709ea83SDmitry Karpeev @*/
18618f3b4b4dSDmitry Karpeev PetscErrorCode  PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg)
1862d709ea83SDmitry Karpeev {
1863d709ea83SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
1864d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1865d709ea83SDmitry Karpeev   PetscBool      match;
1866d709ea83SDmitry Karpeev 
1867d709ea83SDmitry Karpeev   PetscFunctionBegin;
1868d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1869d709ea83SDmitry Karpeev   PetscValidPointer(flg,2);
1870d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1871d709ea83SDmitry Karpeev   if (match) {
1872d709ea83SDmitry Karpeev     if (flg) *flg = osm->dm_subdomains;
1873d709ea83SDmitry Karpeev   }
1874d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1875d709ea83SDmitry Karpeev }
1876