xref: /petsc/src/ksp/pc/impls/gasm/gasm.c (revision 2fa5cd679192b9b390e47ae2d0650965e6b1d9fa)
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 
76a4f0f73SDmitry Karpeev        N    - total number of local subdomains on all processors  (set in PCGASMSetTotalSubdomains() or calculated in PCSetUp_GASM())
86a4f0f73SDmitry Karpeev        n    - actual number of local subdomains on this processor (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains())
96a4f0f73SDmitry Karpeev        nmax - maximum number of local subdomains per processor    (calculated in PCGASMSetTotalSubdomains() or in PCSetUp_GASM())
10b1a0a8a3SJed Brown */
11b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>     /*I "petscpc.h" I*/
12b1a0a8a3SJed Brown 
13b1a0a8a3SJed Brown typedef struct {
14f746d493SDmitry Karpeev   PetscInt   N,n,nmax;
15b1a0a8a3SJed Brown   PetscInt   overlap;             /* overlap requested by user */
16b1a0a8a3SJed Brown   KSP        *ksp;                /* linear solvers for each block */
17f746d493SDmitry Karpeev   Vec        gx,gy;               /* Merged work vectors */
18f746d493SDmitry Karpeev   Vec        *x,*y;               /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */
196a4f0f73SDmitry Karpeev   VecScatter gorestriction;       /* merged restriction to disjoint union of outer subdomains */
206a4f0f73SDmitry Karpeev   VecScatter girestriction;       /* merged restriction to disjoint union of inner subdomains */
216a4f0f73SDmitry Karpeev   IS         *ois;                /* index sets that define the outer (conceptually, overlapping) subdomains */
226a4f0f73SDmitry Karpeev   IS         *iis;                /* index sets that define the inner (conceptually, nonoverlapping) subdomains */
23f746d493SDmitry Karpeev   Mat        *pmat;               /* subdomain block matrices */
24f746d493SDmitry Karpeev   PCGASMType type;                /* use reduced interpolation, restriction or both */
256a4f0f73SDmitry Karpeev   PetscBool  create_local;           /* whether the autocreated subdomains are local or not. */
26b1a0a8a3SJed Brown   PetscBool  type_set;               /* if user set this value (so won't change it for symmetric problems) */
276a4f0f73SDmitry Karpeev   PetscBool  same_subdomain_solvers; /* flag indicating whether all local solvers are same */
28b1a0a8a3SJed Brown   PetscBool  sort_indices;           /* flag to sort subdomain indices */
29f746d493SDmitry Karpeev } PC_GASM;
30b1a0a8a3SJed Brown 
31b1a0a8a3SJed Brown #undef __FUNCT__
3206b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSubdomainView_Private"
3306b43e7eSDmitry Karpeev static PetscErrorCode  PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer)
34af538404SDmitry Karpeev {
35af538404SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
3606b43e7eSDmitry Karpeev   PetscInt       j,nidx;
37af538404SDmitry Karpeev   const PetscInt *idx;
3806b43e7eSDmitry Karpeev   PetscViewer    sviewer;
39af538404SDmitry Karpeev   char           *cidx;
40af538404SDmitry Karpeev   PetscErrorCode ierr;
41af538404SDmitry Karpeev 
42af538404SDmitry Karpeev   PetscFunctionBegin;
4306b43e7eSDmitry Karpeev   if (i < 0 || i > osm->n) SETERRQ2(((PetscObject)viewer)->comm, PETSC_ERR_ARG_WRONG, "Invalid subdomain %D: must nonnegative and less than %D", i, osm->n);
444bde246dSDmitry Karpeev   /* Inner subdomains. */
454bde246dSDmitry Karpeev   ierr = ISGetLocalSize(osm->iis[i], &nidx);CHKERRQ(ierr);
464bde246dSDmitry Karpeev   /*
474bde246dSDmitry Karpeev    No more than 15 characters per index plus a space.
484bde246dSDmitry Karpeev    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
494bde246dSDmitry Karpeev    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
504bde246dSDmitry Karpeev    For nidx == 0, the whole string 16 '\0'.
514bde246dSDmitry Karpeev    */
524bde246dSDmitry Karpeev   ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr);
534bde246dSDmitry Karpeev   ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr);
544bde246dSDmitry Karpeev   ierr = ISGetIndices(osm->iis[i], &idx);CHKERRQ(ierr);
554bde246dSDmitry Karpeev   for (j = 0; j < nidx; ++j) {
564bde246dSDmitry Karpeev     ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr);
574bde246dSDmitry Karpeev   }
584bde246dSDmitry Karpeev   ierr = ISRestoreIndices(osm->iis[i],&idx);CHKERRQ(ierr);
594bde246dSDmitry Karpeev   ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
604bde246dSDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n");CHKERRQ(ierr);
614bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
624bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
634bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
644bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
654bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
664bde246dSDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
674bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
684bde246dSDmitry Karpeev   ierr = PetscFree(cidx);CHKERRQ(ierr);
694bde246dSDmitry Karpeev   /* Outer subdomains. */
706a4f0f73SDmitry Karpeev   ierr = ISGetLocalSize(osm->ois[i], &nidx);CHKERRQ(ierr);
71eec7e3faSDmitry Karpeev   /*
72eec7e3faSDmitry Karpeev    No more than 15 characters per index plus a space.
73eec7e3faSDmitry Karpeev    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
74eec7e3faSDmitry Karpeev    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
75eec7e3faSDmitry Karpeev    For nidx == 0, the whole string 16 '\0'.
76eec7e3faSDmitry Karpeev    */
77eec7e3faSDmitry Karpeev   ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr);
7806b43e7eSDmitry Karpeev   ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr);
796a4f0f73SDmitry Karpeev   ierr = ISGetIndices(osm->ois[i], &idx);CHKERRQ(ierr);
80af538404SDmitry Karpeev   for (j = 0; j < nidx; ++j) {
81af538404SDmitry Karpeev     ierr = PetscViewerStringSPrintf(sviewer,"%D ", idx[j]);CHKERRQ(ierr);
82af538404SDmitry Karpeev   }
836bf464f9SBarry Smith   ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
846a4f0f73SDmitry Karpeev   ierr = ISRestoreIndices(osm->ois[i],&idx);CHKERRQ(ierr);
856a4f0f73SDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n");CHKERRQ(ierr);
8606b43e7eSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
877b23a99aSBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
88af538404SDmitry Karpeev   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
89af538404SDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
907b23a99aSBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
91af538404SDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
9206b43e7eSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
93af538404SDmitry Karpeev   ierr = PetscFree(cidx);CHKERRQ(ierr);
9406b43e7eSDmitry Karpeev   PetscFunctionReturn(0);
9506b43e7eSDmitry Karpeev }
9606b43e7eSDmitry Karpeev 
9706b43e7eSDmitry Karpeev #undef __FUNCT__
9806b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains"
9906b43e7eSDmitry Karpeev static PetscErrorCode  PCGASMPrintSubdomains(PC pc)
10006b43e7eSDmitry Karpeev {
10106b43e7eSDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
10206b43e7eSDmitry Karpeev   const char     *prefix;
10306b43e7eSDmitry Karpeev   char           fname[PETSC_MAX_PATH_LEN+1];
1044bde246dSDmitry Karpeev   PetscInt       i, l, d, count, gcount, *permutation, *numbering;
10506b43e7eSDmitry Karpeev   PetscBool      found;
1064bde246dSDmitry Karpeev   PetscViewer    viewer, sviewer = PETSC_NULL;
10706b43e7eSDmitry Karpeev   PetscErrorCode ierr;
10806b43e7eSDmitry Karpeev 
10906b43e7eSDmitry Karpeev   PetscFunctionBegin;
1104bde246dSDmitry Karpeev   ierr = PetscMalloc2(osm->n, PetscInt, &permutation, osm->n, PetscInt, &numbering);CHKERRQ(ierr);
1114bde246dSDmitry Karpeev   for (i = 0; i < osm->n; ++i) permutation[i] = i;
1124bde246dSDmitry Karpeev   ierr = PetscObjectsGetGlobalNumbering(((PetscObject)pc)->comm, osm->n, (PetscObject*)osm->ois, &gcount, numbering);CHKERRQ(ierr);
1134bde246dSDmitry Karpeev   ierr = PetscSortIntWithPermutation(osm->n, numbering, permutation);CHKERRQ(ierr);
11406b43e7eSDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
11506b43e7eSDmitry Karpeev   ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr);
11606b43e7eSDmitry Karpeev   if (!found) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
1174bde246dSDmitry Karpeev   ierr = PetscViewerASCIIOpen(((PetscObject)pc)->comm,fname,&viewer);CHKERRQ(ierr);
1184bde246dSDmitry Karpeev   /*
1194bde246dSDmitry Karpeev    Make sure the viewer has a name. Otherwise this may cause a deadlock or other weird errors when creating a subcomm viewer:
1204bde246dSDmitry Karpeev    the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed collectively on the comm.
1214bde246dSDmitry Karpeev   */
1224bde246dSDmitry Karpeev   ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr);
1234bde246dSDmitry Karpeev   l    = 0;
1244bde246dSDmitry Karpeev   for (count = 0; count < gcount; ++count) {
1254bde246dSDmitry Karpeev     /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
1264bde246dSDmitry Karpeev     if (l<osm->n) {
1274bde246dSDmitry Karpeev       d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
1284bde246dSDmitry Karpeev       if (numbering[d] == count) {
1294bde246dSDmitry Karpeev         ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
1304bde246dSDmitry Karpeev         ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr);
1314bde246dSDmitry Karpeev         ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
1324bde246dSDmitry Karpeev         ++l;
133af538404SDmitry Karpeev       }
1344bde246dSDmitry Karpeev     }
1354bde246dSDmitry Karpeev     ierr = MPI_Barrier(((PetscObject)pc)->comm);CHKERRQ(ierr);
1364bde246dSDmitry Karpeev   }
1374bde246dSDmitry Karpeev   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1387105d80bSDmitry Karpeev   ierr = PetscFree2(permutation,numbering);CHKERRQ(ierr);
139af538404SDmitry Karpeev   PetscFunctionReturn(0);
140af538404SDmitry Karpeev }
141af538404SDmitry Karpeev 
142af538404SDmitry Karpeev 
143af538404SDmitry Karpeev #undef __FUNCT__
144f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM"
145f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer)
146b1a0a8a3SJed Brown {
147f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
148af538404SDmitry Karpeev   const char     *prefix;
149b1a0a8a3SJed Brown   PetscErrorCode ierr;
150af538404SDmitry Karpeev   PetscMPIInt    rank, size;
151b1a0a8a3SJed Brown   PetscInt       i,bsz;
15206b43e7eSDmitry Karpeev   PetscBool      iascii,view_subdomains=PETSC_FALSE;
153b1a0a8a3SJed Brown   PetscViewer    sviewer;
15406b43e7eSDmitry Karpeev   PetscInt       count, l, gcount, *numbering, *permutation;
1556a4f0f73SDmitry Karpeev   char           overlap[256]     = "user-defined overlap";
1566a4f0f73SDmitry Karpeev   char           gsubdomains[256] = "unknown total number of subdomains";
15706b43e7eSDmitry Karpeev   char           lsubdomains[256] = "unknown number of local  subdomains";
15806b43e7eSDmitry Karpeev   char           msubdomains[256] = "unknown max number of local subdomains";
15911aeaf0aSBarry Smith 
160b1a0a8a3SJed Brown   PetscFunctionBegin;
161af538404SDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm, &size);CHKERRQ(ierr);
162af538404SDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm, &rank);CHKERRQ(ierr);
16306b43e7eSDmitry Karpeev   ierr = PetscMalloc2(osm->n, PetscInt, &permutation, osm->n, PetscInt, &numbering);CHKERRQ(ierr);
16406b43e7eSDmitry Karpeev   for (i = 0; i < osm->n; ++i) permutation[i] = i;
1656a4f0f73SDmitry Karpeev   ierr = PetscObjectsGetGlobalNumbering(((PetscObject)pc)->comm, osm->n, (PetscObject*)osm->ois, &gcount, numbering);CHKERRQ(ierr);
16606b43e7eSDmitry Karpeev   ierr = PetscSortIntWithPermutation(osm->n, numbering, permutation);CHKERRQ(ierr);
16706b43e7eSDmitry Karpeev 
16806b43e7eSDmitry Karpeev   if (osm->overlap >= 0) {
16906b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %D",osm->overlap);CHKERRQ(ierr);
17006b43e7eSDmitry Karpeev   }
1716a4f0f73SDmitry Karpeev   ierr = PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %D",gcount);CHKERRQ(ierr);
17206b43e7eSDmitry Karpeev   if (osm->N > 0) {
17306b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(lsubdomains, sizeof(gsubdomains), "number of local subdomains = %D",osm->N);CHKERRQ(ierr);
17406b43e7eSDmitry Karpeev   }
17506b43e7eSDmitry Karpeev   if (osm->nmax > 0) {
17606b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr);
17706b43e7eSDmitry Karpeev   }
17806b43e7eSDmitry Karpeev 
179af538404SDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
18006b43e7eSDmitry Karpeev   ierr = PetscOptionsGetBool(prefix,"-pc_gasm_view_subdomains",&view_subdomains,PETSC_NULL);CHKERRQ(ierr);
18106b43e7eSDmitry Karpeev 
18206b43e7eSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
183b1a0a8a3SJed Brown   if (iascii) {
18406b43e7eSDmitry Karpeev     /*
18506b43e7eSDmitry Karpeev      Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer:
18606b43e7eSDmitry Karpeev      the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed
18706b43e7eSDmitry Karpeev      collectively on the comm.
18806b43e7eSDmitry Karpeev      */
18906b43e7eSDmitry Karpeev     ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr);
19006b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Generalized additive Schwarz:\n");CHKERRQ(ierr);
19106b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Restriction/interpolation type: %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr);
19206b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",overlap);CHKERRQ(ierr);
19306b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",gsubdomains);CHKERRQ(ierr);
19406b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",lsubdomains);CHKERRQ(ierr);
19506b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",msubdomains);CHKERRQ(ierr);
1967b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
197eec7e3faSDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d:%d] number of locally-supported subdomains = %D\n",(int)rank,(int)size,osm->n);CHKERRQ(ierr);
198af538404SDmitry Karpeev     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1997b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
2006a4f0f73SDmitry Karpeev     /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */
20106b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Subdomain solver info is as follows:\n");CHKERRQ(ierr);
202b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
203b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
20406b43e7eSDmitry Karpeev     /* Make sure that everybody waits for the banner to be printed. */
20506b43e7eSDmitry Karpeev     ierr = MPI_Barrier(((PetscObject)viewer)->comm);CHKERRQ(ierr);
2064bde246dSDmitry Karpeev     /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
20706b43e7eSDmitry Karpeev     l = 0;
20806b43e7eSDmitry Karpeev     for (count = 0; count < gcount; ++count) {
20906b43e7eSDmitry Karpeev       PetscMPIInt srank, ssize;
21006b43e7eSDmitry Karpeev       if (l<osm->n) {
21106b43e7eSDmitry Karpeev         PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
21206b43e7eSDmitry Karpeev         if (numbering[d] == count) {
2136a4f0f73SDmitry Karpeev           ierr = MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize);CHKERRQ(ierr);
2146a4f0f73SDmitry Karpeev           ierr = MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank);CHKERRQ(ierr);
2156a4f0f73SDmitry Karpeev           ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
2166a4f0f73SDmitry Karpeev           ierr = ISGetLocalSize(osm->ois[d],&bsz);CHKERRQ(ierr);
2177b23a99aSBarry Smith           ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_TRUE);CHKERRQ(ierr);
21806b43e7eSDmitry Karpeev           ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%D:%D] (subcomm [%D:%D]) local subdomain number %D, local size = %D\n",(int)rank,(int)size,(int)srank,(int)ssize,d,bsz);CHKERRQ(ierr);
2196a4f0f73SDmitry Karpeev           ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
2206a4f0f73SDmitry Karpeev           ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_FALSE);CHKERRQ(ierr);
22106b43e7eSDmitry Karpeev           if (view_subdomains) {
22206b43e7eSDmitry Karpeev             ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr);
22306b43e7eSDmitry Karpeev           }
2246a4f0f73SDmitry Karpeev           if (!pc->setupcalled) {
2256a4f0f73SDmitry Karpeev             PetscViewerASCIIPrintf(sviewer, "Solver not set up yet: PCSetUp() not yet called\n");CHKERRQ(ierr);
226*2fa5cd67SKarl Rupp           } else {
22706b43e7eSDmitry Karpeev             ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr);
2286a4f0f73SDmitry Karpeev           }
22906b43e7eSDmitry Karpeev           ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
2307b23a99aSBarry Smith           ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
2316a4f0f73SDmitry Karpeev           ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
23206b43e7eSDmitry Karpeev           ++l;
233b1a0a8a3SJed Brown         }
23406b43e7eSDmitry Karpeev       }
23506b43e7eSDmitry Karpeev       ierr = MPI_Barrier(((PetscObject)pc)->comm);CHKERRQ(ierr);
236b1a0a8a3SJed Brown     }
237b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
238b1a0a8a3SJed Brown   }
23906b43e7eSDmitry Karpeev   ierr = PetscFree2(permutation,numbering);CHKERRQ(ierr);
240b1a0a8a3SJed Brown   PetscFunctionReturn(0);
241b1a0a8a3SJed Brown }
242b1a0a8a3SJed Brown 
243af538404SDmitry Karpeev 
244af538404SDmitry Karpeev 
245af538404SDmitry Karpeev 
246af538404SDmitry Karpeev 
247b1a0a8a3SJed Brown #undef __FUNCT__
248f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM"
249f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc)
250b1a0a8a3SJed Brown {
251f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
252b1a0a8a3SJed Brown   PetscErrorCode ierr;
253b1a0a8a3SJed Brown   PetscBool      symset,flg;
25488389755SJed Brown   PetscInt       i;
255c10bc1a1SDmitry Karpeev   PetscMPIInt    rank, size;
256b1a0a8a3SJed Brown   MatReuse       scall = MAT_REUSE_MATRIX;
257b1a0a8a3SJed Brown   KSP            ksp;
258b1a0a8a3SJed Brown   PC             subpc;
259b1a0a8a3SJed Brown   const char     *prefix,*pprefix;
260f746d493SDmitry Karpeev   Vec            x,y;
2616a4f0f73SDmitry Karpeev   PetscInt       oni;       /* Number of indices in the i-th local outer subdomain.               */
2626a4f0f73SDmitry Karpeev   const PetscInt *oidxi;    /* Indices from the i-th subdomain local outer subdomain.             */
2636a4f0f73SDmitry Karpeev   PetscInt       on;        /* Number of indices in the disjoint union of local outer subdomains. */
2646a4f0f73SDmitry Karpeev   PetscInt       *oidx;     /* Indices in the disjoint union of local outer subdomains. */
2656a4f0f73SDmitry Karpeev   IS             gois;      /* Disjoint union the global indices of outer subdomains.             */
2666a4f0f73SDmitry Karpeev   IS             goid;      /* Identity IS of the size of the disjoint union of outer subdomains. */
267f746d493SDmitry Karpeev   PetscScalar    *gxarray, *gyarray;
268*2fa5cd67SKarl Rupp   PetscInt       gofirst;   /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */
269a35b7b57SDmitry Karpeev   DM             *subdomain_dm = PETSC_NULL;
270b1a0a8a3SJed Brown 
271b1a0a8a3SJed Brown   PetscFunctionBegin;
272c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
273c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
274b1a0a8a3SJed Brown   if (!pc->setupcalled) {
275b1a0a8a3SJed Brown 
276b1a0a8a3SJed Brown     if (!osm->type_set) {
277b1a0a8a3SJed Brown       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
278*2fa5cd67SKarl Rupp       if (symset && flg) osm->type = PC_GASM_BASIC;
279b1a0a8a3SJed Brown     }
280b1a0a8a3SJed Brown 
28106b43e7eSDmitry Karpeev     /*
28206b43e7eSDmitry Karpeev      If subdomains have been set, then the local number of subdomains, osm->n, is NOT PETSC_DECIDE and is at least 1.
28306b43e7eSDmitry Karpeev      The total number of subdomains, osm->N is not necessarily set, might be PETSC_DECIDE, and then will have to be calculated from osm->n.
28406b43e7eSDmitry Karpeev      */
28506b43e7eSDmitry Karpeev     if (osm->n == PETSC_DECIDE) {
28669ca1f37SDmitry Karpeev       /* no subdomains given */
28769ca1f37SDmitry Karpeev       /* try pc->dm first */
288fdc48646SDmitry Karpeev       if (pc->dm) {
289fdc48646SDmitry Karpeev         char      ddm_name[1024];
290fdc48646SDmitry Karpeev         DM        ddm;
291fdc48646SDmitry Karpeev         PetscBool flg;
292a35b7b57SDmitry Karpeev         PetscInt  num_subdomains, d;
293a35b7b57SDmitry Karpeev         char      **subdomain_names;
294a35b7b57SDmitry Karpeev         IS        *inner_subdomain_is, *outer_subdomain_is;
295fdc48646SDmitry Karpeev         /* Allow the user to request a decomposition DM by name */
296fdc48646SDmitry Karpeev         ierr = PetscStrncpy(ddm_name, "", 1024);CHKERRQ(ierr);
297a35b7b57SDmitry Karpeev         ierr = PetscOptionsString("-pc_gasm_decomposition","Name of the DM defining the decomposition", "PCSetDM",ddm_name,ddm_name,1024,&flg);CHKERRQ(ierr);
298fdc48646SDmitry Karpeev         if (flg) {
29916621825SDmitry Karpeev           ierr = DMCreateDomainDecompositionDM(pc->dm, ddm_name, &ddm);CHKERRQ(ierr);
300f23aa3ddSBarry Smith           if (!ddm) SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown DM decomposition name %s", ddm_name);
301fdc48646SDmitry Karpeev           ierr = PetscInfo(pc,"Using decomposition DM defined using options database\n");CHKERRQ(ierr);
302fdc48646SDmitry Karpeev           ierr = PCSetDM(pc,ddm);CHKERRQ(ierr);
303fdc48646SDmitry Karpeev         }
304a35b7b57SDmitry Karpeev         ierr = DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);CHKERRQ(ierr);
305a35b7b57SDmitry Karpeev         if (num_subdomains) {
306a35b7b57SDmitry Karpeev           ierr = PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);CHKERRQ(ierr);
30769ca1f37SDmitry Karpeev         }
308a35b7b57SDmitry Karpeev         for (d = 0; d < num_subdomains; ++d) {
309a35b7b57SDmitry Karpeev           if (subdomain_names)    {ierr = PetscFree(subdomain_names[d]);CHKERRQ(ierr);}
310a35b7b57SDmitry Karpeev           if (inner_subdomain_is) {ierr = ISDestroy(&inner_subdomain_is[d]);CHKERRQ(ierr);}
311a35b7b57SDmitry Karpeev           if (outer_subdomain_is) {ierr = ISDestroy(&outer_subdomain_is[d]);CHKERRQ(ierr);}
312fdc48646SDmitry Karpeev         }
313a35b7b57SDmitry Karpeev         ierr = PetscFree(subdomain_names);CHKERRQ(ierr);
314a35b7b57SDmitry Karpeev         ierr = PetscFree(inner_subdomain_is);CHKERRQ(ierr);
315a35b7b57SDmitry Karpeev         ierr = PetscFree(outer_subdomain_is);CHKERRQ(ierr);
316fdc48646SDmitry Karpeev       }
31706b43e7eSDmitry Karpeev       if (osm->n == PETSC_DECIDE) { /* still no subdomains; use one per processor */
31844fe18e5SDmitry Karpeev         osm->nmax = osm->n = 1;
319b1a0a8a3SJed Brown         ierr      = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
320f746d493SDmitry Karpeev         osm->N    = size;
321fdc48646SDmitry Karpeev       }
32206b43e7eSDmitry Karpeev     }
323a35b7b57SDmitry Karpeev     if (!osm->iis) {
324a35b7b57SDmitry Karpeev       /*
325a35b7b57SDmitry Karpeev        The local number of subdomains was set just above, or in PCGASMSetTotalSubdomains(), or in PCGASMSetSubdomains(),
326a35b7b57SDmitry Karpeev        but the actual subdomains have not been supplied (in PCGASMSetSubdomains()).
327a35b7b57SDmitry Karpeev        We create the requisite number of inner subdomains on PETSC_COMM_SELF (for now).
328a35b7b57SDmitry Karpeev        */
329a35b7b57SDmitry Karpeev       ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->overlap,osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr);
330a35b7b57SDmitry Karpeev     }
33106b43e7eSDmitry Karpeev     if (osm->N == PETSC_DECIDE) {
3326ac3741eSJed Brown       struct {PetscInt max,sum;} inwork,outwork;
333f746d493SDmitry Karpeev       /* determine global number of subdomains and the max number of local subdomains */
3346ac3741eSJed Brown       inwork.max = osm->n;
3356ac3741eSJed Brown       inwork.sum = osm->n;
3366ac3741eSJed Brown       ierr       = MPI_Allreduce(&inwork,&outwork,1,MPIU_2INT,PetscMaxSum_Op,((PetscObject)pc)->comm);CHKERRQ(ierr);
3376ac3741eSJed Brown       osm->nmax  = outwork.max;
3386ac3741eSJed Brown       osm->N     = outwork.sum;
339b1a0a8a3SJed Brown     }
3406a4f0f73SDmitry Karpeev 
341b1a0a8a3SJed Brown     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
342b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
343f746d493SDmitry Karpeev     ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&flg,PETSC_NULL);CHKERRQ(ierr);
344f746d493SDmitry Karpeev     if (flg) { ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); }
345b1a0a8a3SJed Brown 
346b1a0a8a3SJed Brown     if (osm->sort_indices) {
347f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
3486a4f0f73SDmitry Karpeev         ierr = ISSort(osm->ois[i]);CHKERRQ(ierr);
3496a4f0f73SDmitry Karpeev         ierr = ISSort(osm->iis[i]);CHKERRQ(ierr);
350b1a0a8a3SJed Brown       }
351b1a0a8a3SJed Brown     }
3526a4f0f73SDmitry Karpeev     /*
3536a4f0f73SDmitry Karpeev      Merge the ISs, create merged vectors and restrictions.
3546a4f0f73SDmitry Karpeev      */
3556a4f0f73SDmitry Karpeev     /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */
3566a4f0f73SDmitry Karpeev     on = 0;
357f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
3586a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
3596a4f0f73SDmitry Karpeev       on  += oni;
360f746d493SDmitry Karpeev     }
3616a4f0f73SDmitry Karpeev     ierr = PetscMalloc(on*sizeof(PetscInt), &oidx);CHKERRQ(ierr);
3626a4f0f73SDmitry Karpeev     on   = 0;
363f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
3646a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
3656a4f0f73SDmitry Karpeev       ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr);
3666a4f0f73SDmitry Karpeev       ierr = PetscMemcpy(oidx+on, oidxi, sizeof(PetscInt)*oni);CHKERRQ(ierr);
3676a4f0f73SDmitry Karpeev       ierr = ISRestoreIndices(osm->ois[i], &oidxi);CHKERRQ(ierr);
3686a4f0f73SDmitry Karpeev       on  += oni;
369f746d493SDmitry Karpeev     }
3706a4f0f73SDmitry Karpeev     ierr = ISCreateGeneral(((PetscObject)(pc))->comm, on, oidx, PETSC_OWN_POINTER, &gois);CHKERRQ(ierr);
371f746d493SDmitry Karpeev     ierr = MatGetVecs(pc->pmat,&x,&y);CHKERRQ(ierr);
3726a4f0f73SDmitry Karpeev     ierr = VecCreateMPI(((PetscObject)pc)->comm, on, PETSC_DECIDE, &osm->gx);CHKERRQ(ierr);
373f746d493SDmitry Karpeev     ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr);
3746a4f0f73SDmitry Karpeev     ierr = VecGetOwnershipRange(osm->gx, &gofirst, PETSC_NULL);CHKERRQ(ierr);
3756a4f0f73SDmitry Karpeev     ierr = ISCreateStride(((PetscObject)pc)->comm,on,gofirst,1, &goid);CHKERRQ(ierr);
3766a4f0f73SDmitry Karpeev     ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr);
3776a4f0f73SDmitry Karpeev     ierr = VecDestroy(&x);CHKERRQ(ierr);
3787105d80bSDmitry Karpeev     ierr = ISDestroy(&gois);CHKERRQ(ierr);
379*2fa5cd67SKarl Rupp 
3806a4f0f73SDmitry Karpeev     /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */
381*2fa5cd67SKarl Rupp     {
382*2fa5cd67SKarl Rupp       PetscInt ini;           /* Number of indices the i-th a local inner subdomain. */
3836a4f0f73SDmitry Karpeev       PetscInt in;            /* Number of indices in the disjoint uniont of local inner subdomains. */
3846a4f0f73SDmitry Karpeev       PetscInt *iidx;         /* Global indices in the merged local inner subdomain. */
3856a4f0f73SDmitry Karpeev       PetscInt *ioidx;        /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
3866a4f0f73SDmitry Karpeev       IS       giis;          /* IS for the disjoint union of inner subdomains. */
3876a4f0f73SDmitry Karpeev       IS       giois;         /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
388f746d493SDmitry Karpeev       /**/
3896a4f0f73SDmitry Karpeev       in = 0;
390f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
3916a4f0f73SDmitry Karpeev         ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr);
3926a4f0f73SDmitry Karpeev         in  += ini;
393f746d493SDmitry Karpeev       }
3946a4f0f73SDmitry Karpeev       ierr = PetscMalloc(in*sizeof(PetscInt), &iidx);CHKERRQ(ierr);
3956a4f0f73SDmitry Karpeev       ierr = PetscMalloc(in*sizeof(PetscInt), &ioidx);CHKERRQ(ierr);
3966a4f0f73SDmitry Karpeev       ierr = VecGetOwnershipRange(osm->gx,&gofirst, PETSC_NULL);CHKERRQ(ierr);
3976a4f0f73SDmitry Karpeev       in   = 0;
3986a4f0f73SDmitry Karpeev       on   = 0;
399f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
4006a4f0f73SDmitry Karpeev         const PetscInt         *iidxi; /* Global indices of the i-th local inner subdomain. */
4016a4f0f73SDmitry Karpeev         ISLocalToGlobalMapping ltogi; /* Map from global to local indices of the i-th outer local subdomain. */
4026a4f0f73SDmitry Karpeev         PetscInt               *ioidxi; /* Local indices of the i-th local inner subdomain within the local outer subdomain. */
4036a4f0f73SDmitry Karpeev         PetscInt               ioni;  /* Number of indices in ioidxi; if ioni != ini the inner subdomain is not a subdomain of the outer subdomain (error). */
4046a4f0f73SDmitry Karpeev         PetscInt               k;
4056a4f0f73SDmitry Karpeev         ierr   = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr);
4066a4f0f73SDmitry Karpeev         ierr   = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
4076a4f0f73SDmitry Karpeev         ierr   = ISGetIndices(osm->iis[i],&iidxi);CHKERRQ(ierr);
4086a4f0f73SDmitry Karpeev         ierr   = PetscMemcpy(iidx+in, iidxi, sizeof(PetscInt)*ini);CHKERRQ(ierr);
4096a4f0f73SDmitry Karpeev         ierr   = ISLocalToGlobalMappingCreateIS(osm->ois[i],&ltogi);CHKERRQ(ierr);
4106a4f0f73SDmitry Karpeev         ioidxi = ioidx+in;
4116a4f0f73SDmitry Karpeev         ierr   = ISGlobalToLocalMappingApply(ltogi,IS_GTOLM_DROP,ini,iidxi,&ioni,ioidxi);CHKERRQ(ierr);
4126a4f0f73SDmitry Karpeev         ierr   = ISLocalToGlobalMappingDestroy(&ltogi);CHKERRQ(ierr);
4136a4f0f73SDmitry Karpeev         ierr   = ISRestoreIndices(osm->iis[i], &iidxi);CHKERRQ(ierr);
4146a4f0f73SDmitry Karpeev         if (ioni != ini) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inner subdomain %D contains %D indices outside of its outer subdomain", i, ini - ioni);
415*2fa5cd67SKarl Rupp         for (k = 0; k < ini; ++k) ioidxi[k] += gofirst+on;
4166a4f0f73SDmitry Karpeev         in += ini;
4176a4f0f73SDmitry Karpeev         on += oni;
418f746d493SDmitry Karpeev       }
4196a4f0f73SDmitry Karpeev       ierr = ISCreateGeneral(((PetscObject)pc)->comm, in, iidx,  PETSC_OWN_POINTER, &giis);CHKERRQ(ierr);
4206a4f0f73SDmitry Karpeev       ierr = ISCreateGeneral(((PetscObject)pc)->comm, in, ioidx, PETSC_OWN_POINTER, &giois);CHKERRQ(ierr);
4216a4f0f73SDmitry Karpeev       ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr);
4226a4f0f73SDmitry Karpeev       ierr = VecDestroy(&y);CHKERRQ(ierr);
4236a4f0f73SDmitry Karpeev       ierr = ISDestroy(&giis);CHKERRQ(ierr);
4246a4f0f73SDmitry Karpeev       ierr = ISDestroy(&giois);CHKERRQ(ierr);
425b1a0a8a3SJed Brown     }
4266a4f0f73SDmitry Karpeev     ierr = ISDestroy(&goid);CHKERRQ(ierr);
427*2fa5cd67SKarl Rupp 
428f746d493SDmitry Karpeev     /* Create the subdomain work vectors. */
429f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->x);CHKERRQ(ierr);
430f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->y);CHKERRQ(ierr);
431f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr);
432f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr);
4336a4f0f73SDmitry Karpeev     for (i=0, on=0; i<osm->n; ++i, on += oni) {
4346a4f0f73SDmitry Karpeev       PetscInt oNi;
4356a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
4366a4f0f73SDmitry Karpeev       ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr);
4376a4f0f73SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr);
4386a4f0f73SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr);
439b1a0a8a3SJed Brown     }
440f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr);
441f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr);
442b1a0a8a3SJed Brown     /* Create the local solvers */
443f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(KSP*),&osm->ksp);CHKERRQ(ierr);
44444fe18e5SDmitry Karpeev     for (i=0; i<osm->n; i++) { /* KSPs are local */
4456a4f0f73SDmitry Karpeev       ierr        = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr);
446b1a0a8a3SJed Brown       ierr        = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr);
447b1a0a8a3SJed Brown       ierr        = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
448b1a0a8a3SJed Brown       ierr        = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
449b1a0a8a3SJed Brown       ierr        = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
450b1a0a8a3SJed Brown       ierr        = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
451b1a0a8a3SJed Brown       ierr        = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
452b1a0a8a3SJed Brown       ierr        = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
453b1a0a8a3SJed Brown       osm->ksp[i] = ksp;
454b1a0a8a3SJed Brown     }
455b1a0a8a3SJed Brown     scall = MAT_INITIAL_MATRIX;
456b1a0a8a3SJed Brown 
457*2fa5cd67SKarl Rupp   } else { /*if (!pc->setupcalled)*/
458b1a0a8a3SJed Brown     /*
459b1a0a8a3SJed Brown        Destroy the blocks from the previous iteration
460b1a0a8a3SJed Brown     */
461b1a0a8a3SJed Brown     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
462f746d493SDmitry Karpeev       ierr  = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
463b1a0a8a3SJed Brown       scall = MAT_INITIAL_MATRIX;
464b1a0a8a3SJed Brown     }
465b1a0a8a3SJed Brown   }
466b1a0a8a3SJed Brown 
467b1a0a8a3SJed Brown   /*
468f746d493SDmitry Karpeev      Extract out the submatrices.
469b1a0a8a3SJed Brown   */
47081a5c4aaSDmitry Karpeev   if (size > 1) {
4716a4f0f73SDmitry Karpeev     ierr = MatGetSubMatricesParallel(pc->pmat,osm->n,osm->ois, osm->ois,scall,&osm->pmat);CHKERRQ(ierr);
472*2fa5cd67SKarl Rupp   } else {
4736a4f0f73SDmitry Karpeev     ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->ois, osm->ois,scall,&osm->pmat);CHKERRQ(ierr);
47481a5c4aaSDmitry Karpeev   }
475b1a0a8a3SJed Brown   if (scall == MAT_INITIAL_MATRIX) {
476b1a0a8a3SJed Brown     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
477f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
478b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,osm->pmat[i]);CHKERRQ(ierr);
479b1a0a8a3SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
480b1a0a8a3SJed Brown     }
481b1a0a8a3SJed Brown   }
482b1a0a8a3SJed Brown 
483b1a0a8a3SJed Brown   /* Return control to the user so that the submatrices can be modified (e.g., to apply
484b1a0a8a3SJed Brown      different boundary conditions for the submatrices than for the global problem) */
4856a4f0f73SDmitry Karpeev   ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
486b1a0a8a3SJed Brown 
487b1a0a8a3SJed Brown   /*
4886a4f0f73SDmitry Karpeev      Loop over submatrices putting them into local ksps
489b1a0a8a3SJed Brown   */
490f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
491b1a0a8a3SJed Brown     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr);
492b1a0a8a3SJed Brown     if (!pc->setupcalled) {
493b1a0a8a3SJed Brown       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
494b1a0a8a3SJed Brown     }
495b1a0a8a3SJed Brown   }
496b1a0a8a3SJed Brown   PetscFunctionReturn(0);
497b1a0a8a3SJed Brown }
498b1a0a8a3SJed Brown 
499b1a0a8a3SJed Brown #undef __FUNCT__
500f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM"
501f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
502b1a0a8a3SJed Brown {
503f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
504b1a0a8a3SJed Brown   PetscErrorCode ierr;
505b1a0a8a3SJed Brown   PetscInt       i;
506b1a0a8a3SJed Brown 
507b1a0a8a3SJed Brown   PetscFunctionBegin;
508f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
509b1a0a8a3SJed Brown     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
510b1a0a8a3SJed Brown   }
511b1a0a8a3SJed Brown   PetscFunctionReturn(0);
512b1a0a8a3SJed Brown }
513b1a0a8a3SJed Brown 
514b1a0a8a3SJed Brown #undef __FUNCT__
515f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM"
516f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y)
517b1a0a8a3SJed Brown {
518f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
519b1a0a8a3SJed Brown   PetscErrorCode ierr;
520f746d493SDmitry Karpeev   PetscInt       i;
521b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
522b1a0a8a3SJed Brown 
523b1a0a8a3SJed Brown   PetscFunctionBegin;
524b1a0a8a3SJed Brown   /*
5256a4f0f73SDmitry Karpeev      Support for limiting the restriction or interpolation only to the inner
526b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
527b1a0a8a3SJed Brown   */
528f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
529b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
530f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
5316a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
532*2fa5cd67SKarl Rupp   } else {
5336a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
534b1a0a8a3SJed Brown   }
5356a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr);
5366a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
5376a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
538*2fa5cd67SKarl Rupp   } else {
5396a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5406a4f0f73SDmitry Karpeev   }
54186b96d36SDmitry Karpeev   /* do the subdomain solves */
54286b96d36SDmitry Karpeev   for (i=0; i<osm->n; ++i) {
543b1a0a8a3SJed Brown     ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
544b1a0a8a3SJed Brown   }
5456a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(y);CHKERRQ(ierr);
5466a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
5476a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5486a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);  PetscFunctionReturn(0);
549*2fa5cd67SKarl Rupp   } else {
5506a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5516a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);  PetscFunctionReturn(0);
5526a4f0f73SDmitry Karpeev   }
553b1a0a8a3SJed Brown }
554b1a0a8a3SJed Brown 
555b1a0a8a3SJed Brown #undef __FUNCT__
556f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM"
557f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y)
558b1a0a8a3SJed Brown {
559f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
560b1a0a8a3SJed Brown   PetscErrorCode ierr;
561f746d493SDmitry Karpeev   PetscInt       i;
562b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
563b1a0a8a3SJed Brown 
564b1a0a8a3SJed Brown   PetscFunctionBegin;
565b1a0a8a3SJed Brown   /*
566b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
567b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
568b1a0a8a3SJed Brown 
569f746d493SDmitry Karpeev      Note: these are reversed from the PCApply_GASM() because we are applying the
570b1a0a8a3SJed Brown      transpose of the three terms
571b1a0a8a3SJed Brown   */
572f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
573b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
574f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
5756a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
576*2fa5cd67SKarl Rupp   } else {
5776a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
578b1a0a8a3SJed Brown   }
5796a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr);
5806a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
5816a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
582*2fa5cd67SKarl Rupp   } else {
5836a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5846a4f0f73SDmitry Karpeev   }
585b1a0a8a3SJed Brown   /* do the local solves */
586ab3e923fSDmitry 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. */
587b1a0a8a3SJed Brown     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
588b1a0a8a3SJed Brown   }
5896a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(y);CHKERRQ(ierr);
5906a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
5916a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5926a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
593*2fa5cd67SKarl Rupp   } else {
5946a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5956a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5966a4f0f73SDmitry Karpeev   }
597b1a0a8a3SJed Brown   PetscFunctionReturn(0);
598b1a0a8a3SJed Brown }
599b1a0a8a3SJed Brown 
600b1a0a8a3SJed Brown #undef __FUNCT__
601a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM"
602a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc)
603b1a0a8a3SJed Brown {
604f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
605b1a0a8a3SJed Brown   PetscErrorCode ierr;
606b1a0a8a3SJed Brown   PetscInt       i;
607b1a0a8a3SJed Brown 
608b1a0a8a3SJed Brown   PetscFunctionBegin;
609b1a0a8a3SJed Brown   if (osm->ksp) {
610f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
611a06653b4SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
612b1a0a8a3SJed Brown     }
613b1a0a8a3SJed Brown   }
614b1a0a8a3SJed Brown   if (osm->pmat) {
615f746d493SDmitry Karpeev     if (osm->n > 0) {
616ab3e923fSDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
617b1a0a8a3SJed Brown     }
618b1a0a8a3SJed Brown   }
619d34fcf5fSBarry Smith   if (osm->x) {
620f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
621fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
622fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
623b1a0a8a3SJed Brown     }
624d34fcf5fSBarry Smith   }
6256bf464f9SBarry Smith   ierr = VecDestroy(&osm->gx);CHKERRQ(ierr);
6266bf464f9SBarry Smith   ierr = VecDestroy(&osm->gy);CHKERRQ(ierr);
627ab3e923fSDmitry Karpeev 
6286a4f0f73SDmitry Karpeev   ierr     = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr);
6296a4f0f73SDmitry Karpeev   ierr     = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr);
6306a4f0f73SDmitry Karpeev   ierr     = PCGASMDestroySubdomains(osm->n,osm->ois,osm->iis);CHKERRQ(ierr);
6316a4f0f73SDmitry Karpeev   osm->ois = 0;
6326a4f0f73SDmitry Karpeev   osm->iis = 0;
633a06653b4SBarry Smith   PetscFunctionReturn(0);
634a06653b4SBarry Smith }
635a06653b4SBarry Smith 
636a06653b4SBarry Smith #undef __FUNCT__
637a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM"
638a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc)
639a06653b4SBarry Smith {
640a06653b4SBarry Smith   PC_GASM        *osm = (PC_GASM*)pc->data;
641a06653b4SBarry Smith   PetscErrorCode ierr;
642a06653b4SBarry Smith   PetscInt       i;
643a06653b4SBarry Smith 
644a06653b4SBarry Smith   PetscFunctionBegin;
645a06653b4SBarry Smith   ierr = PCReset_GASM(pc);CHKERRQ(ierr);
646a06653b4SBarry Smith   if (osm->ksp) {
647a06653b4SBarry Smith     for (i=0; i<osm->n; i++) {
6486bf464f9SBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
649a06653b4SBarry Smith     }
650a06653b4SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
651a06653b4SBarry Smith   }
652a06653b4SBarry Smith   ierr = PetscFree(osm->x);CHKERRQ(ierr);
653a06653b4SBarry Smith   ierr = PetscFree(osm->y);CHKERRQ(ierr);
654c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
655b1a0a8a3SJed Brown   PetscFunctionReturn(0);
656b1a0a8a3SJed Brown }
657b1a0a8a3SJed Brown 
658b1a0a8a3SJed Brown #undef __FUNCT__
659f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM"
660a6dfd86eSKarl Rupp static PetscErrorCode PCSetFromOptions_GASM(PC pc)
661a6dfd86eSKarl Rupp {
662f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
663b1a0a8a3SJed Brown   PetscErrorCode ierr;
664b1a0a8a3SJed Brown   PetscInt       blocks,ovl;
665b1a0a8a3SJed Brown   PetscBool      symset,flg;
666f746d493SDmitry Karpeev   PCGASMType     gasmtype;
667b1a0a8a3SJed Brown 
668b1a0a8a3SJed Brown   PetscFunctionBegin;
669b1a0a8a3SJed Brown   /* set the type to symmetric if matrix is symmetric */
670b1a0a8a3SJed Brown   if (!osm->type_set && pc->pmat) {
671b1a0a8a3SJed Brown     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
672*2fa5cd67SKarl Rupp     if (symset && flg) osm->type = PC_GASM_BASIC;
673b1a0a8a3SJed Brown   }
67444fe18e5SDmitry Karpeev   ierr = PetscOptionsHead("Generalized additive Schwarz options");CHKERRQ(ierr);
6756a4f0f73SDmitry Karpeev   ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
676*2fa5cd67SKarl Rupp 
6776a4f0f73SDmitry Karpeev   osm->create_local = PETSC_TRUE;
678*2fa5cd67SKarl Rupp 
6796a4f0f73SDmitry Karpeev   ierr = PetscOptionsBool("-pc_gasm_subdomains_create_local","Whether to make autocreated subdomains local (true by default)","PCGASMSetTotalSubdomains",osm->create_local,&osm->create_local,&flg);CHKERRQ(ierr);
6806a4f0f73SDmitry Karpeev   if (!osm->create_local) SETERRQ(((PetscObject)pc)->comm, PETSC_ERR_SUP, "No support for autocreation of nonlocal subdomains yet.");
6816a4f0f73SDmitry Karpeev 
6826a4f0f73SDmitry Karpeev   if (flg) {ierr = PCGASMSetTotalSubdomains(pc,blocks,osm->create_local);CHKERRQ(ierr); }
68306b43e7eSDmitry Karpeev   ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
684f746d493SDmitry Karpeev   if (flg) {ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); }
685b1a0a8a3SJed Brown   flg  = PETSC_FALSE;
686f746d493SDmitry Karpeev   ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr);
687f746d493SDmitry Karpeev   if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr); }
688b1a0a8a3SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
689b1a0a8a3SJed Brown   PetscFunctionReturn(0);
690b1a0a8a3SJed Brown }
691b1a0a8a3SJed Brown 
692b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/
693b1a0a8a3SJed Brown 
694b1a0a8a3SJed Brown EXTERN_C_BEGIN
695b1a0a8a3SJed Brown #undef __FUNCT__
69606b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM"
697a35b7b57SDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[])
698b1a0a8a3SJed Brown {
699f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
700b1a0a8a3SJed Brown   PetscErrorCode ierr;
701b1a0a8a3SJed Brown   PetscInt       i;
702b1a0a8a3SJed Brown 
703b1a0a8a3SJed Brown   PetscFunctionBegin;
70469ca1f37SDmitry Karpeev   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, n = %D",n);
705a35b7b57SDmitry Karpeev   if (pc->setupcalled && (n != osm->n || iis || ois)) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp().");
706b1a0a8a3SJed Brown 
707b1a0a8a3SJed Brown   if (!pc->setupcalled) {
70893c1ec32SDmitry Karpeev     osm->n   = n;
7096a4f0f73SDmitry Karpeev     osm->ois = 0;
7106a4f0f73SDmitry Karpeev     osm->iis = 0;
711a35b7b57SDmitry Karpeev     if (ois) {
712a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);}
713b1a0a8a3SJed Brown     }
714a35b7b57SDmitry Karpeev     if (iis) {
715a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);}
716b1a0a8a3SJed Brown     }
7176a4f0f73SDmitry Karpeev     ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr);
718a35b7b57SDmitry Karpeev     if (ois) {
7196a4f0f73SDmitry Karpeev       ierr = PetscMalloc(n*sizeof(IS),&osm->ois);CHKERRQ(ierr);
720*2fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->ois[i] = ois[i];
7216a4f0f73SDmitry Karpeev       /* Flag indicating that the user has set outer subdomains, so PCGASM should not increase their size. */
722b1a0a8a3SJed Brown       osm->overlap = -1;
723a35b7b57SDmitry Karpeev       if (!iis) {
7246a4f0f73SDmitry Karpeev         ierr = PetscMalloc(n*sizeof(IS),&osm->iis);CHKERRQ(ierr);
725a35b7b57SDmitry Karpeev         for (i=0; i<n; i++) {
726a35b7b57SDmitry Karpeev           for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);}
727a35b7b57SDmitry Karpeev           osm->iis[i] = ois[i];
728a35b7b57SDmitry Karpeev         }
729a35b7b57SDmitry Karpeev       }
730a35b7b57SDmitry Karpeev     }
731a35b7b57SDmitry Karpeev     if (iis) {
732a35b7b57SDmitry Karpeev       ierr = PetscMalloc(n*sizeof(IS),&osm->iis);CHKERRQ(ierr);
733*2fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->iis[i] = iis[i];
734a35b7b57SDmitry Karpeev       if (!ois) {
735a35b7b57SDmitry Karpeev         ierr = PetscMalloc(n*sizeof(IS),&osm->ois);CHKERRQ(ierr);
736a35b7b57SDmitry Karpeev         for (i=0; i<n; i++) {
737a35b7b57SDmitry Karpeev           for (i=0; i<n; i++) {
738a35b7b57SDmitry Karpeev             ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);
739a35b7b57SDmitry Karpeev             osm->ois[i] = iis[i];
740a35b7b57SDmitry Karpeev           }
741a35b7b57SDmitry Karpeev         }
742a35b7b57SDmitry Karpeev         if (osm->overlap > 0) {
743a35b7b57SDmitry Karpeev           /* Extend the "overlapping" regions by a number of steps */
744a35b7b57SDmitry Karpeev           ierr = MatIncreaseOverlap(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr);
745a35b7b57SDmitry Karpeev         }
746a35b7b57SDmitry Karpeev       }
747b1a0a8a3SJed Brown     }
748b1a0a8a3SJed Brown   }
749b1a0a8a3SJed Brown   PetscFunctionReturn(0);
750b1a0a8a3SJed Brown }
751b1a0a8a3SJed Brown EXTERN_C_END
752b1a0a8a3SJed Brown 
753b1a0a8a3SJed Brown EXTERN_C_BEGIN
754b1a0a8a3SJed Brown #undef __FUNCT__
7556a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains_GASM"
7560adebc6cSBarry Smith PetscErrorCode  PCGASMSetTotalSubdomains_GASM(PC pc,PetscInt N, PetscBool create_local)
7570adebc6cSBarry Smith {
758f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
759b1a0a8a3SJed Brown   PetscErrorCode ierr;
760b1a0a8a3SJed Brown   PetscMPIInt    rank,size;
761b1a0a8a3SJed Brown   PetscInt       n;
7626a4f0f73SDmitry Karpeev   PetscInt       Nmin, Nmax;
7635fd66863SKarl Rupp 
764b1a0a8a3SJed Brown   PetscFunctionBegin;
7656a4f0f73SDmitry Karpeev   if (!create_local) SETERRQ(((PetscObject)pc)->comm, PETSC_ERR_SUP, "No suppor for autocreation of nonlocal subdomains.");
7666a4f0f73SDmitry Karpeev   if (N < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be > 0, N = %D",N);
767adcec1e5SJed Brown   ierr = MPI_Allreduce(&N,&Nmin,1,MPIU_INT,MPIU_MIN,((PetscObject)pc)->comm);CHKERRQ(ierr);
768adcec1e5SJed Brown   ierr = MPI_Allreduce(&N,&Nmax,1,MPIU_INT,MPIU_MAX,((PetscObject)pc)->comm);CHKERRQ(ierr);
769f23aa3ddSBarry Smith   if (Nmin != Nmax) SETERRQ2(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONG, "All processors must use the same number of subdomains.  min(N) = %D != %D = max(N)", Nmin, Nmax);
770b1a0a8a3SJed Brown 
7716a4f0f73SDmitry Karpeev   osm->create_local = create_local;
772b1a0a8a3SJed Brown   /*
773b1a0a8a3SJed Brown      Split the subdomains equally among all processors
774b1a0a8a3SJed Brown   */
775b1a0a8a3SJed Brown   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
776b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
777b1a0a8a3SJed Brown   n    = N/size + ((N % size) > rank);
7786a4f0f73SDmitry Karpeev   if (!n) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Process %d must have at least one subdomain: total processors %d total blocks %D",(int)rank,(int)size,N);
779f746d493SDmitry Karpeev   if (pc->setupcalled && n != osm->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before PCSetUp().");
780b1a0a8a3SJed Brown   if (!pc->setupcalled) {
7816a4f0f73SDmitry Karpeev     ierr      = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr);
78293c1ec32SDmitry Karpeev     osm->N    = N;
783f746d493SDmitry Karpeev     osm->n    = n;
78406b43e7eSDmitry Karpeev     osm->nmax = N/size + ((N%size) ? 1 : 0);
7856a4f0f73SDmitry Karpeev     osm->ois  = 0;
7866a4f0f73SDmitry Karpeev     osm->iis  = 0;
787b1a0a8a3SJed Brown   }
788b1a0a8a3SJed Brown   PetscFunctionReturn(0);
78906b43e7eSDmitry Karpeev }
790b1a0a8a3SJed Brown EXTERN_C_END
791b1a0a8a3SJed Brown 
792b1a0a8a3SJed Brown EXTERN_C_BEGIN
793b1a0a8a3SJed Brown #undef __FUNCT__
794f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM"
7957087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
796b1a0a8a3SJed Brown {
797f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
798b1a0a8a3SJed Brown 
799b1a0a8a3SJed Brown   PetscFunctionBegin;
800b1a0a8a3SJed Brown   if (ovl < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
801f746d493SDmitry Karpeev   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
802*2fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
803b1a0a8a3SJed Brown   PetscFunctionReturn(0);
804b1a0a8a3SJed Brown }
805b1a0a8a3SJed Brown EXTERN_C_END
806b1a0a8a3SJed Brown 
807b1a0a8a3SJed Brown EXTERN_C_BEGIN
808b1a0a8a3SJed Brown #undef __FUNCT__
809f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM"
8107087cfbeSBarry Smith PetscErrorCode  PCGASMSetType_GASM(PC pc,PCGASMType type)
811b1a0a8a3SJed Brown {
812f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
813b1a0a8a3SJed Brown 
814b1a0a8a3SJed Brown   PetscFunctionBegin;
815b1a0a8a3SJed Brown   osm->type     = type;
816b1a0a8a3SJed Brown   osm->type_set = PETSC_TRUE;
817b1a0a8a3SJed Brown   PetscFunctionReturn(0);
818b1a0a8a3SJed Brown }
819b1a0a8a3SJed Brown EXTERN_C_END
820b1a0a8a3SJed Brown 
821b1a0a8a3SJed Brown EXTERN_C_BEGIN
822b1a0a8a3SJed Brown #undef __FUNCT__
823f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM"
8247087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort)
825b1a0a8a3SJed Brown {
826f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
827b1a0a8a3SJed Brown 
828b1a0a8a3SJed Brown   PetscFunctionBegin;
829b1a0a8a3SJed Brown   osm->sort_indices = doSort;
830b1a0a8a3SJed Brown   PetscFunctionReturn(0);
831b1a0a8a3SJed Brown }
832b1a0a8a3SJed Brown EXTERN_C_END
833b1a0a8a3SJed Brown 
834b1a0a8a3SJed Brown EXTERN_C_BEGIN
835b1a0a8a3SJed Brown #undef __FUNCT__
836f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM"
83744fe18e5SDmitry Karpeev /*
83844fe18e5SDmitry Karpeev    FIX: This routine might need to be modified once multiple ranks per subdomain are allowed.
83944fe18e5SDmitry Karpeev         In particular, it would upset the global subdomain number calculation.
84044fe18e5SDmitry Karpeev */
8417087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
842b1a0a8a3SJed Brown {
843f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
844b1a0a8a3SJed Brown   PetscErrorCode ierr;
845b1a0a8a3SJed Brown 
846b1a0a8a3SJed Brown   PetscFunctionBegin;
847ab3e923fSDmitry Karpeev   if (osm->n < 1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here");
848b1a0a8a3SJed Brown 
849*2fa5cd67SKarl Rupp   if (n) *n = osm->n;
850ab3e923fSDmitry Karpeev   if (first) {
851ab3e923fSDmitry Karpeev     ierr    = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr);
852ab3e923fSDmitry Karpeev     *first -= osm->n;
853b1a0a8a3SJed Brown   }
854b1a0a8a3SJed Brown   if (ksp) {
855b1a0a8a3SJed Brown     /* Assume that local solves are now different; not necessarily
85606b43e7eSDmitry Karpeev        true, though!  This flag is used only for PCView_GASM() */
857b1a0a8a3SJed Brown     *ksp                        = osm->ksp;
8586a4f0f73SDmitry Karpeev     osm->same_subdomain_solvers = PETSC_FALSE;
859b1a0a8a3SJed Brown   }
860b1a0a8a3SJed Brown   PetscFunctionReturn(0);
861ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */
862b1a0a8a3SJed Brown EXTERN_C_END
863b1a0a8a3SJed Brown 
864b1a0a8a3SJed Brown 
865b1a0a8a3SJed Brown #undef __FUNCT__
86606b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains"
867b1a0a8a3SJed Brown /*@C
86806b43e7eSDmitry Karpeev     PCGASMSetSubdomains - Sets the subdomains for this processor
86906b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
870b1a0a8a3SJed Brown 
871b1a0a8a3SJed Brown     Collective on PC
872b1a0a8a3SJed Brown 
873b1a0a8a3SJed Brown     Input Parameters:
874b1a0a8a3SJed Brown +   pc  - the preconditioner context
87506b43e7eSDmitry Karpeev .   n   - the number of subdomains for this processor
8766a4f0f73SDmitry Karpeev .   iis - the index sets that define this processor's local inner subdomains
877b1a0a8a3SJed Brown          (or PETSC_NULL for PETSc to determine subdomains)
8786a4f0f73SDmitry Karpeev -   ois- the index sets that define this processor's local outer subdomains
8796a4f0f73SDmitry Karpeev          (or PETSC_NULL to use the same as iis)
880b1a0a8a3SJed Brown 
881b1a0a8a3SJed Brown     Notes:
88206b43e7eSDmitry Karpeev     The IS indices use the parallel, global numbering of the vector entries.
8836a4f0f73SDmitry Karpeev     Inner subdomains are those where the correction is applied.
8846a4f0f73SDmitry Karpeev     Outer subdomains are those where the residual necessary to obtain the
8856a4f0f73SDmitry Karpeev     corrections is obtained (see PCGASMType for the use of inner/outer subdomains).
8866a4f0f73SDmitry Karpeev     Both inner and outer subdomains can extend over several processors.
8876a4f0f73SDmitry Karpeev     This processor's portion of a subdomain is known as a local subdomain.
8886a4f0f73SDmitry Karpeev 
8896a4f0f73SDmitry Karpeev     By default the GASM preconditioner uses 1 (local) subdomain per processor.
8906a4f0f73SDmitry Karpeev     Use PCGASMSetTotalSubdomains() to set the total number of subdomains across
8916a4f0f73SDmitry Karpeev     all processors that PCGASM will create automatically, and to specify whether
8926a4f0f73SDmitry Karpeev     they should be local or not.
8936a4f0f73SDmitry Karpeev 
894b1a0a8a3SJed Brown 
895b1a0a8a3SJed Brown     Level: advanced
896b1a0a8a3SJed Brown 
89706b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz
898b1a0a8a3SJed Brown 
8996a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
90006b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
901b1a0a8a3SJed Brown @*/
9026a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[])
903b1a0a8a3SJed Brown {
904b1a0a8a3SJed Brown   PetscErrorCode ierr;
905b1a0a8a3SJed Brown 
906b1a0a8a3SJed Brown   PetscFunctionBegin;
907b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9086a4f0f73SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr);
909b1a0a8a3SJed Brown   PetscFunctionReturn(0);
910b1a0a8a3SJed Brown }
911b1a0a8a3SJed Brown 
912b1a0a8a3SJed Brown #undef __FUNCT__
9136a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains"
914b1a0a8a3SJed Brown /*@C
9156a4f0f73SDmitry Karpeev     PCGASMSetTotalSubdomains - Sets the total number of subdomains to use in the generalized additive
9166a4f0f73SDmitry Karpeev     Schwarz preconditioner.  The number of subdomains is cumulative across all processors in pc's
9176a4f0f73SDmitry Karpeev     communicator. Either all or no processors in the PC communicator must call this routine with
9186a4f0f73SDmitry Karpeev     the same N.  The subdomains will be created automatically during PCSetUp().
919b1a0a8a3SJed Brown 
920b1a0a8a3SJed Brown     Collective on PC
921b1a0a8a3SJed Brown 
922b1a0a8a3SJed Brown     Input Parameters:
923b1a0a8a3SJed Brown +   pc           - the preconditioner context
9246a4f0f73SDmitry Karpeev .   N            - the total number of subdomains cumulative across all processors
9256a4f0f73SDmitry Karpeev -   create_local - whether the subdomains to be created are to be local
926b1a0a8a3SJed Brown 
927b1a0a8a3SJed Brown     Options Database Key:
9286a4f0f73SDmitry Karpeev     To set the total number of subdomains and let PCGASM autocreate them, rather than specify the index sets, use the following options:
9296a4f0f73SDmitry Karpeev +    -pc_gasm_total_subdomains <n>                  - sets the total number of subdomains to be autocreated by PCGASM
9306a4f0f73SDmitry Karpeev -    -pc_gasm_subdomains_create_local <true|false>  - whether autocreated subdomains should be local or not (default is true)
931b1a0a8a3SJed Brown 
93206b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.
933b1a0a8a3SJed Brown 
9346a4f0f73SDmitry Karpeev 
9356a4f0f73SDmitry Karpeev     Use PCGASMSetSubdomains() to set subdomains explicitly or to set different numbers
9366a4f0f73SDmitry Karpeev     of subdomains per processor.
937b1a0a8a3SJed Brown 
938b1a0a8a3SJed Brown     Level: advanced
939b1a0a8a3SJed Brown 
940f746d493SDmitry Karpeev .keywords: PC, GASM, set, total, global, subdomains, additive Schwarz
941b1a0a8a3SJed Brown 
94206b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
943f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
944b1a0a8a3SJed Brown @*/
9456a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetTotalSubdomains(PC pc,PetscInt N, PetscBool create_local)
946b1a0a8a3SJed Brown {
947b1a0a8a3SJed Brown   PetscErrorCode ierr;
948b1a0a8a3SJed Brown 
949b1a0a8a3SJed Brown   PetscFunctionBegin;
950b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9516a4f0f73SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetTotalSubdomains_C",(PC,PetscInt,PetscBool),(pc,N,create_local));CHKERRQ(ierr);
952b1a0a8a3SJed Brown   PetscFunctionReturn(0);
953b1a0a8a3SJed Brown }
954b1a0a8a3SJed Brown 
955b1a0a8a3SJed Brown #undef __FUNCT__
956f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap"
957b1a0a8a3SJed Brown /*@
958f746d493SDmitry Karpeev     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
959b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
960b1a0a8a3SJed Brown     PC communicator must call this routine.
961b1a0a8a3SJed Brown 
962b1a0a8a3SJed Brown     Logically Collective on PC
963b1a0a8a3SJed Brown 
964b1a0a8a3SJed Brown     Input Parameters:
965b1a0a8a3SJed Brown +   pc  - the preconditioner context
966b1a0a8a3SJed Brown -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
967b1a0a8a3SJed Brown 
968b1a0a8a3SJed Brown     Options Database Key:
96906b43e7eSDmitry Karpeev .   -pc_gasm_overlap <overlap> - Sets overlap
970b1a0a8a3SJed Brown 
971b1a0a8a3SJed Brown     Notes:
97206b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.  To use
9736a4f0f73SDmitry Karpeev     multiple subdomain per perocessor, see PCGASMSetTotalSubdomains() or
9746a4f0f73SDmitry Karpeev     PCGASMSetSubdomains() (and the option -pc_gasm_total_subdomains <n>).
975b1a0a8a3SJed Brown 
976b1a0a8a3SJed Brown     The overlap defaults to 1, so if one desires that no additional
977b1a0a8a3SJed Brown     overlap be computed beyond what may have been set with a call to
9786a4f0f73SDmitry Karpeev     PCGASMSetTotalSubdomains() or PCGASMSetSubdomains(), then ovl
979b1a0a8a3SJed Brown     must be set to be 0.  In particular, if one does not explicitly set
98006b43e7eSDmitry Karpeev     the subdomains in application code, then all overlap would be computed
981f746d493SDmitry Karpeev     internally by PETSc, and using an overlap of 0 would result in an GASM
982b1a0a8a3SJed Brown     variant that is equivalent to the block Jacobi preconditioner.
983b1a0a8a3SJed Brown 
984b1a0a8a3SJed Brown     Note that one can define initial index sets with any overlap via
98506b43e7eSDmitry Karpeev     PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows
98606b43e7eSDmitry Karpeev     PETSc to extend that overlap further, if desired.
987b1a0a8a3SJed Brown 
988b1a0a8a3SJed Brown     Level: intermediate
989b1a0a8a3SJed Brown 
990f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap
991b1a0a8a3SJed Brown 
9926a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
99306b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
994b1a0a8a3SJed Brown @*/
9957087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap(PC pc,PetscInt ovl)
996b1a0a8a3SJed Brown {
997b1a0a8a3SJed Brown   PetscErrorCode ierr;
998b1a0a8a3SJed Brown 
999b1a0a8a3SJed Brown   PetscFunctionBegin;
1000b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1001b1a0a8a3SJed Brown   PetscValidLogicalCollectiveInt(pc,ovl,2);
1002f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
1003b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1004b1a0a8a3SJed Brown }
1005b1a0a8a3SJed Brown 
1006b1a0a8a3SJed Brown #undef __FUNCT__
1007f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType"
1008b1a0a8a3SJed Brown /*@
1009f746d493SDmitry Karpeev     PCGASMSetType - Sets the type of restriction and interpolation used
1010b1a0a8a3SJed Brown     for local problems in the additive Schwarz method.
1011b1a0a8a3SJed Brown 
1012b1a0a8a3SJed Brown     Logically Collective on PC
1013b1a0a8a3SJed Brown 
1014b1a0a8a3SJed Brown     Input Parameters:
1015b1a0a8a3SJed Brown +   pc  - the preconditioner context
1016f746d493SDmitry Karpeev -   type - variant of GASM, one of
1017b1a0a8a3SJed Brown .vb
1018f746d493SDmitry Karpeev       PC_GASM_BASIC       - full interpolation and restriction
1019f746d493SDmitry Karpeev       PC_GASM_RESTRICT    - full restriction, local processor interpolation
1020f746d493SDmitry Karpeev       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
1021f746d493SDmitry Karpeev       PC_GASM_NONE        - local processor restriction and interpolation
1022b1a0a8a3SJed Brown .ve
1023b1a0a8a3SJed Brown 
1024b1a0a8a3SJed Brown     Options Database Key:
1025f746d493SDmitry Karpeev .   -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1026b1a0a8a3SJed Brown 
1027b1a0a8a3SJed Brown     Level: intermediate
1028b1a0a8a3SJed Brown 
1029f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1030b1a0a8a3SJed Brown 
10316a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1032f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1033b1a0a8a3SJed Brown @*/
10347087cfbeSBarry Smith PetscErrorCode  PCGASMSetType(PC pc,PCGASMType type)
1035b1a0a8a3SJed Brown {
1036b1a0a8a3SJed Brown   PetscErrorCode ierr;
1037b1a0a8a3SJed Brown 
1038b1a0a8a3SJed Brown   PetscFunctionBegin;
1039b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1040b1a0a8a3SJed Brown   PetscValidLogicalCollectiveEnum(pc,type,2);
1041f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr);
1042b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1043b1a0a8a3SJed Brown }
1044b1a0a8a3SJed Brown 
1045b1a0a8a3SJed Brown #undef __FUNCT__
1046f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices"
1047b1a0a8a3SJed Brown /*@
1048f746d493SDmitry Karpeev     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
1049b1a0a8a3SJed Brown 
1050b1a0a8a3SJed Brown     Logically Collective on PC
1051b1a0a8a3SJed Brown 
1052b1a0a8a3SJed Brown     Input Parameters:
1053b1a0a8a3SJed Brown +   pc  - the preconditioner context
1054b1a0a8a3SJed Brown -   doSort - sort the subdomain indices
1055b1a0a8a3SJed Brown 
1056b1a0a8a3SJed Brown     Level: intermediate
1057b1a0a8a3SJed Brown 
1058f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1059b1a0a8a3SJed Brown 
10606a4f0f73SDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(),
1061f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1062b1a0a8a3SJed Brown @*/
10637087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices(PC pc,PetscBool doSort)
1064b1a0a8a3SJed Brown {
1065b1a0a8a3SJed Brown   PetscErrorCode ierr;
1066b1a0a8a3SJed Brown 
1067b1a0a8a3SJed Brown   PetscFunctionBegin;
1068b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1069b1a0a8a3SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
1070f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
1071b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1072b1a0a8a3SJed Brown }
1073b1a0a8a3SJed Brown 
1074b1a0a8a3SJed Brown #undef __FUNCT__
1075f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP"
1076b1a0a8a3SJed Brown /*@C
1077f746d493SDmitry Karpeev    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
1078b1a0a8a3SJed Brown    this processor.
1079b1a0a8a3SJed Brown 
1080b1a0a8a3SJed Brown    Collective on PC iff first_local is requested
1081b1a0a8a3SJed Brown 
1082b1a0a8a3SJed Brown    Input Parameter:
1083b1a0a8a3SJed Brown .  pc - the preconditioner context
1084b1a0a8a3SJed Brown 
1085b1a0a8a3SJed Brown    Output Parameters:
1086b1a0a8a3SJed Brown +  n_local - the number of blocks on this processor or PETSC_NULL
1087b1a0a8a3SJed Brown .  first_local - the global number of the first block on this processor or PETSC_NULL,
1088b1a0a8a3SJed Brown                  all processors must request or all must pass PETSC_NULL
1089b1a0a8a3SJed Brown -  ksp - the array of KSP contexts
1090b1a0a8a3SJed Brown 
1091b1a0a8a3SJed Brown    Note:
1092f746d493SDmitry Karpeev    After PCGASMGetSubKSP() the array of KSPes is not to be freed
1093b1a0a8a3SJed Brown 
1094b1a0a8a3SJed Brown    Currently for some matrix implementations only 1 block per processor
1095b1a0a8a3SJed Brown    is supported.
1096b1a0a8a3SJed Brown 
1097f746d493SDmitry Karpeev    You must call KSPSetUp() before calling PCGASMGetSubKSP().
1098b1a0a8a3SJed Brown 
1099b1a0a8a3SJed Brown    Level: advanced
1100b1a0a8a3SJed Brown 
1101f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context
1102b1a0a8a3SJed Brown 
11036a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap(),
1104f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(),
1105b1a0a8a3SJed Brown @*/
11067087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1107b1a0a8a3SJed Brown {
1108b1a0a8a3SJed Brown   PetscErrorCode ierr;
1109b1a0a8a3SJed Brown 
1110b1a0a8a3SJed Brown   PetscFunctionBegin;
1111b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1112f746d493SDmitry Karpeev   ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
1113b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1114b1a0a8a3SJed Brown }
1115b1a0a8a3SJed Brown 
1116b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/
1117b1a0a8a3SJed Brown /*MC
1118f746d493SDmitry Karpeev    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1119b1a0a8a3SJed Brown            its own KSP object.
1120b1a0a8a3SJed Brown 
1121b1a0a8a3SJed Brown    Options Database Keys:
112206b43e7eSDmitry Karpeev +  -pc_gasm_total_block_count <n> - Sets total number of local subdomains (known as blocks) to be distributed among processors
112306b43e7eSDmitry Karpeev .  -pc_gasm_view_subdomains       - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view
112406b43e7eSDmitry Karpeev .  -pc_gasm_print_subdomains      - activates the printing of subdomain indices in PCSetUp()
112506b43e7eSDmitry Karpeev .  -pc_gasm_overlap <ovl>         - Sets overlap by which to (automatically) extend local subdomains
1126f746d493SDmitry Karpeev -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1127b1a0a8a3SJed Brown 
1128b1a0a8a3SJed Brown      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
1129f746d493SDmitry Karpeev       will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
1130f746d493SDmitry Karpeev       -pc_gasm_type basic to use the standard GASM.
1131b1a0a8a3SJed Brown 
1132b1a0a8a3SJed Brown    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
1133b1a0a8a3SJed Brown      than one processor. Defaults to one block per processor.
1134b1a0a8a3SJed Brown 
1135b1a0a8a3SJed Brown      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1136b1a0a8a3SJed Brown         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
1137b1a0a8a3SJed Brown 
1138f746d493SDmitry Karpeev      To set the options on the solvers separate for each block call PCGASMGetSubKSP()
1139b1a0a8a3SJed Brown          and set the options directly on the resulting KSP object (you can access its PC
1140b1a0a8a3SJed Brown          with KSPGetPC())
1141b1a0a8a3SJed Brown 
1142b1a0a8a3SJed Brown 
1143b1a0a8a3SJed Brown    Level: beginner
1144b1a0a8a3SJed Brown 
1145b1a0a8a3SJed Brown    Concepts: additive Schwarz method
1146b1a0a8a3SJed Brown 
1147b1a0a8a3SJed Brown     References:
1148b1a0a8a3SJed Brown     An additive variant of the Schwarz alternating method for the case of many subregions
1149b1a0a8a3SJed Brown     M Dryja, OB Widlund - Courant Institute, New York University Technical report
1150b1a0a8a3SJed Brown 
1151b1a0a8a3SJed Brown     Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1152b1a0a8a3SJed Brown     Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X.
1153b1a0a8a3SJed Brown 
1154b1a0a8a3SJed Brown .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
115506b43e7eSDmitry Karpeev            PCBJACOBI, PCGASMSetUseTrueLocal(), PCGASMGetSubKSP(), PCGASMSetSubdomains(),
11566a4f0f73SDmitry Karpeev            PCGASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType()
1157b1a0a8a3SJed Brown 
1158b1a0a8a3SJed Brown M*/
1159b1a0a8a3SJed Brown 
1160b1a0a8a3SJed Brown EXTERN_C_BEGIN
1161b1a0a8a3SJed Brown #undef __FUNCT__
1162f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM"
11637087cfbeSBarry Smith PetscErrorCode  PCCreate_GASM(PC pc)
1164b1a0a8a3SJed Brown {
1165b1a0a8a3SJed Brown   PetscErrorCode ierr;
1166f746d493SDmitry Karpeev   PC_GASM        *osm;
1167b1a0a8a3SJed Brown 
1168b1a0a8a3SJed Brown   PetscFunctionBegin;
1169f746d493SDmitry Karpeev   ierr = PetscNewLog(pc,PC_GASM,&osm);CHKERRQ(ierr);
1170*2fa5cd67SKarl Rupp 
1171ab3e923fSDmitry Karpeev   osm->N                      = PETSC_DECIDE;
117206b43e7eSDmitry Karpeev   osm->n                      = PETSC_DECIDE;
1173ab3e923fSDmitry Karpeev   osm->nmax                   = 0;
1174b1a0a8a3SJed Brown   osm->overlap                = 1;
1175b1a0a8a3SJed Brown   osm->ksp                    = 0;
11766a4f0f73SDmitry Karpeev   osm->gorestriction          = 0;
11776a4f0f73SDmitry Karpeev   osm->girestriction          = 0;
1178ab3e923fSDmitry Karpeev   osm->gx                     = 0;
1179ab3e923fSDmitry Karpeev   osm->gy                     = 0;
1180b1a0a8a3SJed Brown   osm->x                      = 0;
1181b1a0a8a3SJed Brown   osm->y                      = 0;
11826a4f0f73SDmitry Karpeev   osm->ois                    = 0;
11836a4f0f73SDmitry Karpeev   osm->iis                    = 0;
1184b1a0a8a3SJed Brown   osm->pmat                   = 0;
1185f746d493SDmitry Karpeev   osm->type                   = PC_GASM_RESTRICT;
11866a4f0f73SDmitry Karpeev   osm->same_subdomain_solvers = PETSC_TRUE;
1187b1a0a8a3SJed Brown   osm->sort_indices           = PETSC_TRUE;
1188b1a0a8a3SJed Brown 
1189b1a0a8a3SJed Brown   pc->data                 = (void*)osm;
1190f746d493SDmitry Karpeev   pc->ops->apply           = PCApply_GASM;
1191f746d493SDmitry Karpeev   pc->ops->applytranspose  = PCApplyTranspose_GASM;
1192f746d493SDmitry Karpeev   pc->ops->setup           = PCSetUp_GASM;
1193a06653b4SBarry Smith   pc->ops->reset           = PCReset_GASM;
1194f746d493SDmitry Karpeev   pc->ops->destroy         = PCDestroy_GASM;
1195f746d493SDmitry Karpeev   pc->ops->setfromoptions  = PCSetFromOptions_GASM;
1196f746d493SDmitry Karpeev   pc->ops->setuponblocks   = PCSetUpOnBlocks_GASM;
1197f746d493SDmitry Karpeev   pc->ops->view            = PCView_GASM;
1198b1a0a8a3SJed Brown   pc->ops->applyrichardson = 0;
1199b1a0a8a3SJed Brown 
120006b43e7eSDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSubdomains_C","PCGASMSetSubdomains_GASM",
120106b43e7eSDmitry Karpeev                                            PCGASMSetSubdomains_GASM);CHKERRQ(ierr);
12026a4f0f73SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetTotalSubdomains_C","PCGASMSetTotalSubdomains_GASM",
12036a4f0f73SDmitry Karpeev                                            PCGASMSetTotalSubdomains_GASM);CHKERRQ(ierr);
1204f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetOverlap_C","PCGASMSetOverlap_GASM",
1205f746d493SDmitry Karpeev                                            PCGASMSetOverlap_GASM);CHKERRQ(ierr);
1206f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetType_C","PCGASMSetType_GASM",
1207f746d493SDmitry Karpeev                                            PCGASMSetType_GASM);CHKERRQ(ierr);
1208f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSortIndices_C","PCGASMSetSortIndices_GASM",
1209f746d493SDmitry Karpeev                                            PCGASMSetSortIndices_GASM);CHKERRQ(ierr);
1210f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMGetSubKSP_C","PCGASMGetSubKSP_GASM",
1211f746d493SDmitry Karpeev                                            PCGASMGetSubKSP_GASM);CHKERRQ(ierr);
1212b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1213b1a0a8a3SJed Brown }
1214b1a0a8a3SJed Brown EXTERN_C_END
1215b1a0a8a3SJed Brown 
1216b1a0a8a3SJed Brown 
1217b1a0a8a3SJed Brown #undef __FUNCT__
121806b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains"
1219b1a0a8a3SJed Brown /*@C
122006b43e7eSDmitry Karpeev    PCGASMCreateLocalSubdomains - Creates n local index sets for the overlapping
122106b43e7eSDmitry Karpeev    Schwarz preconditioner for a any problem based on its matrix.
1222b1a0a8a3SJed Brown 
1223b1a0a8a3SJed Brown    Collective
1224b1a0a8a3SJed Brown 
1225b1a0a8a3SJed Brown    Input Parameters:
1226b1a0a8a3SJed Brown +  A       - The global matrix operator
12276a4f0f73SDmitry Karpeev .  overlap - amount of overlap in outer subdomains
122806b43e7eSDmitry Karpeev -  n       - the number of local subdomains
1229b1a0a8a3SJed Brown 
1230b1a0a8a3SJed Brown    Output Parameters:
12316a4f0f73SDmitry Karpeev +  iis - the array of index sets defining the local inner subdomains (on which the correction is applied)
12326a4f0f73SDmitry Karpeev -  ois - the array of index sets defining the local outer subdomains (on which the residual is computed)
1233b1a0a8a3SJed Brown 
1234b1a0a8a3SJed Brown    Level: advanced
1235b1a0a8a3SJed Brown 
12366a4f0f73SDmitry Karpeev    Note: this generates n nonoverlapping local inner subdomains on PETSC_COMM_SELF;
12376a4f0f73SDmitry Karpeev          PCGASM will generate the overlap from these if you use them in PCGASMSetSubdomains() and set a
12386a4f0f73SDmitry Karpeev          nonzero overlap with PCGASMSetOverlap()
1239b1a0a8a3SJed Brown 
1240b1a0a8a3SJed Brown     In the Fortran version you must provide the array outis[] already allocated of length n.
1241b1a0a8a3SJed Brown 
1242f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1243b1a0a8a3SJed Brown 
124406b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains()
1245b1a0a8a3SJed Brown @*/
12466a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt overlap, PetscInt n, IS *iis[], IS *ois[])
1247b1a0a8a3SJed Brown {
1248b1a0a8a3SJed Brown   MatPartitioning mpart;
1249b1a0a8a3SJed Brown   const char      *prefix;
125011bd1e4dSLisandro Dalcin   PetscErrorCode  (*f)(Mat,MatReuse,Mat*);
1251b1a0a8a3SJed Brown   PetscMPIInt     size;
1252b1a0a8a3SJed Brown   PetscInt        i,j,rstart,rend,bs;
125311bd1e4dSLisandro Dalcin   PetscBool       isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
1254b1a0a8a3SJed Brown   Mat             Ad     = PETSC_NULL, adj;
1255b1a0a8a3SJed Brown   IS              ispart,isnumb,*is;
1256b1a0a8a3SJed Brown   PetscErrorCode  ierr;
1257b1a0a8a3SJed Brown 
1258b1a0a8a3SJed Brown   PetscFunctionBegin;
1259b1a0a8a3SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
12606a4f0f73SDmitry Karpeev   PetscValidPointer(iis,4);
1261b1a0a8a3SJed Brown   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
1262b1a0a8a3SJed Brown 
1263b1a0a8a3SJed Brown   /* Get prefix, row distribution, and block size */
1264b1a0a8a3SJed Brown   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1265b1a0a8a3SJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1266b1a0a8a3SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1267b1a0a8a3SJed 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);
1268b1a0a8a3SJed Brown 
1269b1a0a8a3SJed Brown   /* Get diagonal block from matrix if possible */
1270b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1271b1a0a8a3SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr);
1272b1a0a8a3SJed Brown   if (f) {
127311bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
1274b1a0a8a3SJed Brown   } else if (size == 1) {
127511bd1e4dSLisandro Dalcin     Ad = A;
1276b1a0a8a3SJed Brown   }
1277b1a0a8a3SJed Brown   if (Ad) {
1278251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1279251f4c67SDmitry Karpeev     if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
1280b1a0a8a3SJed Brown   }
1281b1a0a8a3SJed Brown   if (Ad && n > 1) {
1282b1a0a8a3SJed Brown     PetscBool  match,done;
1283b1a0a8a3SJed Brown     /* Try to setup a good matrix partitioning if available */
1284b1a0a8a3SJed Brown     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
1285b1a0a8a3SJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
1286b1a0a8a3SJed Brown     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1287251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
1288b1a0a8a3SJed Brown     if (!match) {
1289251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
1290b1a0a8a3SJed Brown     }
1291b1a0a8a3SJed Brown     if (!match) { /* assume a "good" partitioner is available */
12921a83f524SJed Brown       PetscInt       na;
12931a83f524SJed Brown       const PetscInt *ia,*ja;
1294b1a0a8a3SJed Brown       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1295b1a0a8a3SJed Brown       if (done) {
1296b1a0a8a3SJed Brown         /* Build adjacency matrix by hand. Unfortunately a call to
1297b1a0a8a3SJed Brown            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1298b1a0a8a3SJed Brown            remove the block-aij structure and we cannot expect
1299b1a0a8a3SJed Brown            MatPartitioning to split vertices as we need */
13001a83f524SJed Brown         PetscInt       i,j,len,nnz,cnt,*iia=0,*jja=0;
13011a83f524SJed Brown         const PetscInt *row;
1302b1a0a8a3SJed Brown         nnz = 0;
1303b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* count number of nonzeros */
1304b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1305b1a0a8a3SJed Brown           row = ja + ia[i];
1306b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1307b1a0a8a3SJed Brown             if (row[j] == i) { /* don't count diagonal */
1308b1a0a8a3SJed Brown               len--; break;
1309b1a0a8a3SJed Brown             }
1310b1a0a8a3SJed Brown           }
1311b1a0a8a3SJed Brown           nnz += len;
1312b1a0a8a3SJed Brown         }
1313b1a0a8a3SJed Brown         ierr   = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr);
1314b1a0a8a3SJed Brown         ierr   = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr);
1315b1a0a8a3SJed Brown         nnz    = 0;
1316b1a0a8a3SJed Brown         iia[0] = 0;
1317b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* fill adjacency */
1318b1a0a8a3SJed Brown           cnt = 0;
1319b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1320b1a0a8a3SJed Brown           row = ja + ia[i];
1321b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1322*2fa5cd67SKarl Rupp             if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */
1323b1a0a8a3SJed Brown           }
1324b1a0a8a3SJed Brown           nnz += cnt;
1325b1a0a8a3SJed Brown           iia[i+1] = nnz;
1326b1a0a8a3SJed Brown         }
1327b1a0a8a3SJed Brown         /* Partitioning of the adjacency matrix */
1328b1a0a8a3SJed Brown         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);CHKERRQ(ierr);
1329b1a0a8a3SJed Brown         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
1330b1a0a8a3SJed Brown         ierr      = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
1331b1a0a8a3SJed Brown         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
1332b1a0a8a3SJed Brown         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
13336bf464f9SBarry Smith         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
1334b1a0a8a3SJed Brown         foundpart = PETSC_TRUE;
1335b1a0a8a3SJed Brown       }
1336b1a0a8a3SJed Brown       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1337b1a0a8a3SJed Brown     }
13386bf464f9SBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1339b1a0a8a3SJed Brown   }
1340b1a0a8a3SJed Brown   ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr);
1341b1a0a8a3SJed Brown   if (!foundpart) {
1342b1a0a8a3SJed Brown 
1343b1a0a8a3SJed Brown     /* Partitioning by contiguous chunks of rows */
1344b1a0a8a3SJed Brown 
1345b1a0a8a3SJed Brown     PetscInt mbs   = (rend-rstart)/bs;
1346b1a0a8a3SJed Brown     PetscInt start = rstart;
1347b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1348b1a0a8a3SJed Brown       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
1349b1a0a8a3SJed Brown       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
1350b1a0a8a3SJed Brown       start += count;
1351b1a0a8a3SJed Brown     }
1352b1a0a8a3SJed Brown 
1353b1a0a8a3SJed Brown   } else {
1354b1a0a8a3SJed Brown 
1355b1a0a8a3SJed Brown     /* Partitioning by adjacency of diagonal block  */
1356b1a0a8a3SJed Brown 
1357b1a0a8a3SJed Brown     const PetscInt *numbering;
1358b1a0a8a3SJed Brown     PetscInt       *count,nidx,*indices,*newidx,start=0;
1359b1a0a8a3SJed Brown     /* Get node count in each partition */
1360b1a0a8a3SJed Brown     ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr);
1361b1a0a8a3SJed Brown     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
1362b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1363b1a0a8a3SJed Brown       for (i=0; i<n; i++) count[i] *= bs;
1364b1a0a8a3SJed Brown     }
1365b1a0a8a3SJed Brown     /* Build indices from node numbering */
1366b1a0a8a3SJed Brown     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1367b1a0a8a3SJed Brown     ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr);
1368b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1369b1a0a8a3SJed Brown     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
1370b1a0a8a3SJed Brown     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
1371b1a0a8a3SJed Brown     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
1372b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1373b1a0a8a3SJed Brown       ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr);
1374*2fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
1375*2fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
1376*2fa5cd67SKarl Rupp       }
1377b1a0a8a3SJed Brown       ierr    = PetscFree(indices);CHKERRQ(ierr);
1378b1a0a8a3SJed Brown       nidx   *= bs;
1379b1a0a8a3SJed Brown       indices = newidx;
1380b1a0a8a3SJed Brown     }
1381b1a0a8a3SJed Brown     /* Shift to get global indices */
1382b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] += rstart;
1383b1a0a8a3SJed Brown 
1384b1a0a8a3SJed Brown     /* Build the index sets for each block */
1385b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1386b1a0a8a3SJed Brown       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
1387b1a0a8a3SJed Brown       ierr   = ISSort(is[i]);CHKERRQ(ierr);
1388b1a0a8a3SJed Brown       start += count[i];
1389b1a0a8a3SJed Brown     }
1390b1a0a8a3SJed Brown 
1391b1a0a8a3SJed Brown     ierr = PetscFree(count);
1392b1a0a8a3SJed Brown     ierr = PetscFree(indices);
1393fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1394fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
1395b1a0a8a3SJed Brown   }
13966a4f0f73SDmitry Karpeev   *iis = is;
13976a4f0f73SDmitry Karpeev   if (!ois) PetscFunctionReturn(0);
13986a4f0f73SDmitry Karpeev   /*
13996a4f0f73SDmitry Karpeev    Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap
14006a4f0f73SDmitry Karpeev    has been requested, copy the inner subdomains over so they can be modified.
14016a4f0f73SDmitry Karpeev    */
14026a4f0f73SDmitry Karpeev   ierr = PetscMalloc(n*sizeof(IS),ois);CHKERRQ(ierr);
14036a4f0f73SDmitry Karpeev   for (i=0; i<n; ++i) {
14046a4f0f73SDmitry Karpeev     if (overlap > 0) { /* With positive overlap, (*iis)[i] will be modified */
14056a4f0f73SDmitry Karpeev       ierr = ISDuplicate((*iis)[i],(*ois)+i);CHKERRQ(ierr);
14066a4f0f73SDmitry Karpeev       ierr = ISCopy((*iis)[i],(*ois)[i]);CHKERRQ(ierr);
14076a4f0f73SDmitry Karpeev     } else {
14086a4f0f73SDmitry Karpeev       ierr      = PetscObjectReference((PetscObject)(*iis)[i]);CHKERRQ(ierr);
14096a4f0f73SDmitry Karpeev       (*ois)[i] = (*iis)[i];
14106a4f0f73SDmitry Karpeev     }
14116a4f0f73SDmitry Karpeev   }
14126a4f0f73SDmitry Karpeev   if (overlap > 0) {
14136a4f0f73SDmitry Karpeev     /* Extend the "overlapping" regions by a number of steps */
14146a4f0f73SDmitry Karpeev     ierr = MatIncreaseOverlap(A,n,*ois,overlap);CHKERRQ(ierr);
14156a4f0f73SDmitry Karpeev   }
1416b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1417b1a0a8a3SJed Brown }
1418b1a0a8a3SJed Brown 
1419b1a0a8a3SJed Brown #undef __FUNCT__
1420f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains"
1421b1a0a8a3SJed Brown /*@C
1422f746d493SDmitry Karpeev    PCGASMDestroySubdomains - Destroys the index sets created with
142306b43e7eSDmitry Karpeev    PCGASMCreateLocalSubdomains() or PCGASMCreateSubdomains2D. Should be
142406b43e7eSDmitry Karpeev    called after setting subdomains with PCGASMSetSubdomains().
1425b1a0a8a3SJed Brown 
1426b1a0a8a3SJed Brown    Collective
1427b1a0a8a3SJed Brown 
1428b1a0a8a3SJed Brown    Input Parameters:
1429b1a0a8a3SJed Brown +  n   - the number of index sets
14306a4f0f73SDmitry Karpeev .  iis - the array of inner subdomains,
14316a4f0f73SDmitry Karpeev -  ois - the array of outer subdomains, can be PETSC_NULL
1432b1a0a8a3SJed Brown 
14336a4f0f73SDmitry Karpeev    Level: intermediate
14346a4f0f73SDmitry Karpeev 
14356a4f0f73SDmitry Karpeev    Notes: this is merely a convenience subroutine that walks each list,
14366a4f0f73SDmitry Karpeev    destroys each IS on the list, and then frees the list.
1437b1a0a8a3SJed Brown 
1438f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1439b1a0a8a3SJed Brown 
144006b43e7eSDmitry Karpeev .seealso: PCGASMCreateLocalSubdomains(), PCGASMSetSubdomains()
1441b1a0a8a3SJed Brown @*/
14426a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMDestroySubdomains(PetscInt n, IS iis[], IS ois[])
1443b1a0a8a3SJed Brown {
1444b1a0a8a3SJed Brown   PetscInt       i;
1445b1a0a8a3SJed Brown   PetscErrorCode ierr;
14465fd66863SKarl Rupp 
1447b1a0a8a3SJed Brown   PetscFunctionBegin;
1448a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1449a35b7b57SDmitry Karpeev   if (iis) {
14506a4f0f73SDmitry Karpeev     PetscValidPointer(iis,2);
1451a35b7b57SDmitry Karpeev     for (i=0; i<n; i++) {
1452a35b7b57SDmitry Karpeev       ierr = ISDestroy(&iis[i]);CHKERRQ(ierr);
1453a35b7b57SDmitry Karpeev     }
14546a4f0f73SDmitry Karpeev     ierr = PetscFree(iis);CHKERRQ(ierr);
1455a35b7b57SDmitry Karpeev   }
14566a4f0f73SDmitry Karpeev   if (ois) {
1457a35b7b57SDmitry Karpeev     for (i=0; i<n; i++) {
1458a35b7b57SDmitry Karpeev       ierr = ISDestroy(&ois[i]);CHKERRQ(ierr);
1459a35b7b57SDmitry Karpeev     }
14606a4f0f73SDmitry Karpeev     ierr = PetscFree(ois);CHKERRQ(ierr);
1461b1a0a8a3SJed Brown   }
1462b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1463b1a0a8a3SJed Brown }
1464b1a0a8a3SJed Brown 
1465af538404SDmitry Karpeev 
1466af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \
1467af538404SDmitry Karpeev   {                                                                                                       \
1468af538404SDmitry Karpeev     PetscInt first_row = first/M, last_row = last/M+1;                                                     \
1469af538404SDmitry Karpeev     /*                                                                                                    \
1470af538404SDmitry Karpeev      Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners       \
1471af538404SDmitry Karpeev      of the bounding box of the intersection of the subdomain with the local ownership range (local       \
1472af538404SDmitry Karpeev      subdomain).                                                                                          \
1473af538404SDmitry Karpeev      Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows       \
1474af538404SDmitry Karpeev      of the intersection.                                                                                 \
1475af538404SDmitry Karpeev     */                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                             \
1476af538404SDmitry Karpeev     /* ylow_loc is the grid row containing the first element of the local sumbdomain */                   \
1477eec7e3faSDmitry Karpeev     *ylow_loc = PetscMax(first_row,ylow);                                                                    \
1478af538404SDmitry Karpeev     /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1479eec7e3faSDmitry Karpeev     *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft;                                                                            \
1480af538404SDmitry Karpeev     /* yhigh_loc is the grid row above the last local subdomain element */                                                                    \
1481eec7e3faSDmitry Karpeev     *yhigh_loc = PetscMin(last_row,yhigh);                                                                                                     \
1482af538404SDmitry Karpeev     /* xright is the offset of the end of the  local subdomain within its grid row (might actually be outside the local subdomain) */         \
1483eec7e3faSDmitry Karpeev     *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright;                                                                          \
1484af538404SDmitry Karpeev     /* Now compute the size of the local subdomain n. */ \
1485c3518bceSDmitry Karpeev     *n = 0;                                               \
1486eec7e3faSDmitry Karpeev     if (*ylow_loc < *yhigh_loc) {                           \
1487af538404SDmitry Karpeev       PetscInt width = xright-xleft;                     \
1488eec7e3faSDmitry Karpeev       *n += width*(*yhigh_loc-*ylow_loc-1);                 \
1489eec7e3faSDmitry Karpeev       *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \
1490eec7e3faSDmitry Karpeev       *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \
1491af538404SDmitry Karpeev     } \
1492af538404SDmitry Karpeev   }
1493af538404SDmitry Karpeev 
1494c3518bceSDmitry Karpeev 
1495eec7e3faSDmitry Karpeev 
1496eec7e3faSDmitry Karpeev #undef __FUNCT__
1497f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D"
1498b1a0a8a3SJed Brown /*@
1499f746d493SDmitry Karpeev    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1500b1a0a8a3SJed Brown    preconditioner for a two-dimensional problem on a regular grid.
1501b1a0a8a3SJed Brown 
1502af538404SDmitry Karpeev    Collective
1503b1a0a8a3SJed Brown 
1504b1a0a8a3SJed Brown    Input Parameters:
150506b43e7eSDmitry Karpeev +  M, N               - the global number of grid points in the x and y directions
1506af538404SDmitry Karpeev .  Mdomains, Ndomains - the global number of subdomains in the x and y directions
1507b1a0a8a3SJed Brown .  dof                - degrees of freedom per node
1508b1a0a8a3SJed Brown -  overlap            - overlap in mesh lines
1509b1a0a8a3SJed Brown 
1510b1a0a8a3SJed Brown    Output Parameters:
1511af538404SDmitry Karpeev +  Nsub - the number of local subdomains created
15126a4f0f73SDmitry Karpeev .  iis  - array of index sets defining inner (nonoverlapping) subdomains
15136a4f0f73SDmitry Karpeev -  ois  - array of index sets defining outer (overlapping, if overlap > 0) subdomains
1514b1a0a8a3SJed Brown 
1515b1a0a8a3SJed Brown 
1516b1a0a8a3SJed Brown    Level: advanced
1517b1a0a8a3SJed Brown 
1518f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid
1519b1a0a8a3SJed Brown 
15206a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1521f746d493SDmitry Karpeev           PCGASMSetOverlap()
1522b1a0a8a3SJed Brown @*/
15236a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMCreateSubdomains2D(PC pc, PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap, PetscInt *nsub,IS **iis,IS **ois)
1524b1a0a8a3SJed Brown {
1525b1a0a8a3SJed Brown   PetscErrorCode ierr;
15262bbb417fSDmitry Karpeev   PetscMPIInt    size, rank;
15272bbb417fSDmitry Karpeev   PetscInt       i, j;
15282bbb417fSDmitry Karpeev   PetscInt       maxheight, maxwidth;
15292bbb417fSDmitry Karpeev   PetscInt       xstart, xleft, xright, xleft_loc, xright_loc;
15302bbb417fSDmitry Karpeev   PetscInt       ystart, ylow,  yhigh,  ylow_loc,  yhigh_loc;
1531eec7e3faSDmitry Karpeev   PetscInt       x[2][2], y[2][2], n[2];
15322bbb417fSDmitry Karpeev   PetscInt       first, last;
15332bbb417fSDmitry Karpeev   PetscInt       nidx, *idx;
15342bbb417fSDmitry Karpeev   PetscInt       ii,jj,s,q,d;
1535af538404SDmitry Karpeev   PetscInt       k,kk;
15362bbb417fSDmitry Karpeev   PetscMPIInt    color;
15372bbb417fSDmitry Karpeev   MPI_Comm       comm, subcomm;
15386a4f0f73SDmitry Karpeev   IS             **xis = 0, **is = ois, **is_local = iis;
1539d34fcf5fSBarry Smith 
1540b1a0a8a3SJed Brown   PetscFunctionBegin;
15412bbb417fSDmitry Karpeev   ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr);
15422bbb417fSDmitry Karpeev   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
15432bbb417fSDmitry Karpeev   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
15442bbb417fSDmitry Karpeev   ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr);
1545d34fcf5fSBarry 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) "
15462bbb417fSDmitry Karpeev                                       "does not respect the number of degrees of freedom per grid point %D", first, last, dof);
1547d34fcf5fSBarry Smith 
1548af538404SDmitry Karpeev   /* Determine the number of domains with nonzero intersections with the local ownership range. */
15492bbb417fSDmitry Karpeev   s      = 0;
1550af538404SDmitry Karpeev   ystart = 0;
1551af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1552af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1553af538404SDmitry 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);
1554eec7e3faSDmitry Karpeev     /* Vertical domain limits with an overlap. */
1555eec7e3faSDmitry Karpeev     ylow   = PetscMax(ystart - overlap,0);
1556eec7e3faSDmitry Karpeev     yhigh  = PetscMin(ystart + maxheight + overlap,N);
1557b1a0a8a3SJed Brown     xstart = 0;
1558af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1559af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1560af538404SDmitry 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);
1561eec7e3faSDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1562eec7e3faSDmitry Karpeev       xleft  = PetscMax(xstart - overlap,0);
1563eec7e3faSDmitry Karpeev       xright = PetscMin(xstart + maxwidth + overlap,M);
1564af538404SDmitry Karpeev       /*
1565af538404SDmitry Karpeev          Determine whether this subdomain intersects this processor's ownership range of pc->pmat.
1566af538404SDmitry Karpeev       */
1567c3518bceSDmitry Karpeev       PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1568*2fa5cd67SKarl Rupp       if (nidx) ++s;
1569af538404SDmitry Karpeev       xstart += maxwidth;
1570af538404SDmitry Karpeev     } /* for (i = 0; i < Mdomains; ++i) */
1571af538404SDmitry Karpeev     ystart += maxheight;
1572af538404SDmitry Karpeev   } /* for (j = 0; j < Ndomains; ++j) */
1573*2fa5cd67SKarl Rupp 
1574af538404SDmitry Karpeev   /* Now we can allocate the necessary number of ISs. */
1575af538404SDmitry Karpeev   *nsub  = s;
1576af538404SDmitry Karpeev   ierr   = PetscMalloc((*nsub)*sizeof(IS*),is);CHKERRQ(ierr);
1577af538404SDmitry Karpeev   ierr   = PetscMalloc((*nsub)*sizeof(IS*),is_local);CHKERRQ(ierr);
1578af538404SDmitry Karpeev   s      = 0;
1579af538404SDmitry Karpeev   ystart = 0;
1580af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1581af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1582af538404SDmitry 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);
1583af538404SDmitry Karpeev     /* Vertical domain limits with an overlap. */
1584eec7e3faSDmitry Karpeev     y[0][0] = PetscMax(ystart - overlap,0);
1585eec7e3faSDmitry Karpeev     y[0][1] = PetscMin(ystart + maxheight + overlap,N);
1586eec7e3faSDmitry Karpeev     /* Vertical domain limits without an overlap. */
1587eec7e3faSDmitry Karpeev     y[1][0] = ystart;
1588eec7e3faSDmitry Karpeev     y[1][1] = PetscMin(ystart + maxheight,N);
1589eec7e3faSDmitry Karpeev     xstart  = 0;
1590af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1591af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1592af538404SDmitry 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);
1593af538404SDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1594eec7e3faSDmitry Karpeev       x[0][0] = PetscMax(xstart - overlap,0);
1595eec7e3faSDmitry Karpeev       x[0][1] = PetscMin(xstart + maxwidth + overlap,M);
1596eec7e3faSDmitry Karpeev       /* Horizontal domain limits without an overlap. */
1597eec7e3faSDmitry Karpeev       x[1][0] = xstart;
1598eec7e3faSDmitry Karpeev       x[1][1] = PetscMin(xstart+maxwidth,M);
15992bbb417fSDmitry Karpeev       /*
16002bbb417fSDmitry Karpeev          Determine whether this domain intersects this processor's ownership range of pc->pmat.
16012bbb417fSDmitry Karpeev          Do this twice: first for the domains with overlaps, and once without.
16022bbb417fSDmitry Karpeev          During the first pass create the subcommunicators, and use them on the second pass as well.
16032bbb417fSDmitry Karpeev       */
16042bbb417fSDmitry Karpeev       for (q = 0; q < 2; ++q) {
16052bbb417fSDmitry Karpeev         /*
16062bbb417fSDmitry Karpeev           domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
16072bbb417fSDmitry Karpeev           according to whether the domain with an overlap or without is considered.
16082bbb417fSDmitry Karpeev         */
16092bbb417fSDmitry Karpeev         xleft = x[q][0]; xright = x[q][1];
16102bbb417fSDmitry Karpeev         ylow  = y[q][0]; yhigh  = y[q][1];
1611c3518bceSDmitry Karpeev         PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1612af538404SDmitry Karpeev         nidx *= dof;
1613eec7e3faSDmitry Karpeev         n[q]  = nidx;
16142bbb417fSDmitry Karpeev         /*
1615eec7e3faSDmitry Karpeev          Based on the counted number of indices in the local domain *with an overlap*,
1616af538404SDmitry Karpeev          construct a subcommunicator of all the processors supporting this domain.
1617eec7e3faSDmitry Karpeev          Observe that a domain with an overlap might have nontrivial local support,
1618eec7e3faSDmitry Karpeev          while the domain without an overlap might not.  Hence, the decision to participate
1619eec7e3faSDmitry Karpeev          in the subcommunicator must be based on the domain with an overlap.
16202bbb417fSDmitry Karpeev          */
16212bbb417fSDmitry Karpeev         if (q == 0) {
1622*2fa5cd67SKarl Rupp           if (nidx) color = 1;
1623*2fa5cd67SKarl Rupp           else color = MPI_UNDEFINED;
1624*2fa5cd67SKarl Rupp 
16252bbb417fSDmitry Karpeev           ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr);
16262bbb417fSDmitry Karpeev         }
1627af538404SDmitry Karpeev         /*
1628eec7e3faSDmitry Karpeev          Proceed only if the number of local indices *with an overlap* is nonzero.
1629af538404SDmitry Karpeev          */
1630eec7e3faSDmitry Karpeev         if (n[0]) {
1631*2fa5cd67SKarl Rupp           if (q == 0) xis = is;
1632af538404SDmitry Karpeev           if (q == 1) {
1633af538404SDmitry Karpeev             /*
1634eec7e3faSDmitry Karpeev              The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1635af538404SDmitry Karpeev              Moreover, if the overlap is zero, the two ISs are identical.
1636af538404SDmitry Karpeev              */
1637b1a0a8a3SJed Brown             if (overlap == 0) {
1638eec7e3faSDmitry Karpeev               (*is_local)[s] = (*is)[s];
16392bbb417fSDmitry Karpeev               ierr           = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr);
16402bbb417fSDmitry Karpeev               continue;
1641d34fcf5fSBarry Smith             } else {
16426a4f0f73SDmitry Karpeev               xis     = is_local;
1643eec7e3faSDmitry Karpeev               subcomm = ((PetscObject)(*is)[s])->comm;
16442bbb417fSDmitry Karpeev             }
1645af538404SDmitry Karpeev           } /* if (q == 1) */
1646eec7e3faSDmitry Karpeev           idx  = PETSC_NULL;
16472bbb417fSDmitry Karpeev           ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1648eec7e3faSDmitry Karpeev           if (nidx) {
16492bbb417fSDmitry Karpeev             k = 0;
16502bbb417fSDmitry Karpeev             for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1651af538404SDmitry Karpeev               PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft;
1652af538404SDmitry Karpeev               PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright;
1653af538404SDmitry Karpeev               kk = dof*(M*jj + x0);
16542bbb417fSDmitry Karpeev               for (ii=x0; ii<x1; ++ii) {
16552bbb417fSDmitry Karpeev                 for (d = 0; d < dof; ++d) {
16562bbb417fSDmitry Karpeev                   idx[k++] = kk++;
1657b1a0a8a3SJed Brown                 }
1658b1a0a8a3SJed Brown               }
1659b1a0a8a3SJed Brown             }
1660eec7e3faSDmitry Karpeev           }
16616a4f0f73SDmitry Karpeev           ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr);
1662eec7e3faSDmitry Karpeev         }/* if (n[0]) */
16632bbb417fSDmitry Karpeev       }/* for (q = 0; q < 2; ++q) */
1664*2fa5cd67SKarl Rupp       if (n[0]) ++s;
1665af538404SDmitry Karpeev       xstart += maxwidth;
1666af538404SDmitry Karpeev     } /* for (i = 0; i < Mdomains; ++i) */
16672bbb417fSDmitry Karpeev     ystart += maxheight;
1668af538404SDmitry Karpeev   } /* for (j = 0; j < Ndomains; ++j) */
1669b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1670b1a0a8a3SJed Brown }
1671b1a0a8a3SJed Brown 
1672b1a0a8a3SJed Brown #undef __FUNCT__
167306b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains"
1674b1a0a8a3SJed Brown /*@C
167506b43e7eSDmitry Karpeev     PCGASMGetSubdomains - Gets the subdomains supported on this processor
167606b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
1677b1a0a8a3SJed Brown 
1678b1a0a8a3SJed Brown     Not Collective
1679b1a0a8a3SJed Brown 
1680b1a0a8a3SJed Brown     Input Parameter:
1681b1a0a8a3SJed Brown .   pc - the preconditioner context
1682b1a0a8a3SJed Brown 
1683b1a0a8a3SJed Brown     Output Parameters:
1684b1a0a8a3SJed Brown +   n   - the number of subdomains for this processor (default value = 1)
16856a4f0f73SDmitry Karpeev .   iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be PETSC_NULL)
16866a4f0f73SDmitry Karpeev -   ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be PETSC_NULL)
1687b1a0a8a3SJed Brown 
1688b1a0a8a3SJed Brown 
1689b1a0a8a3SJed Brown     Notes:
16906a4f0f73SDmitry Karpeev     The user is responsible for destroying the ISs and freeing the returned arrays.
1691b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
1692b1a0a8a3SJed Brown 
1693b1a0a8a3SJed Brown     Level: advanced
1694b1a0a8a3SJed Brown 
169506b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz
1696b1a0a8a3SJed Brown 
16976a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
169806b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubmatrices()
1699b1a0a8a3SJed Brown @*/
17006a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[])
1701b1a0a8a3SJed Brown {
1702f746d493SDmitry Karpeev   PC_GASM        *osm;
1703b1a0a8a3SJed Brown   PetscErrorCode ierr;
1704b1a0a8a3SJed Brown   PetscBool      match;
17056a4f0f73SDmitry Karpeev   PetscInt       i;
17065fd66863SKarl Rupp 
1707b1a0a8a3SJed Brown   PetscFunctionBegin;
1708b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1709251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1710f23aa3ddSBarry Smith   if (!match) SETERRQ2(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1711f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1712ab3e923fSDmitry Karpeev   if (n) *n = osm->n;
17136a4f0f73SDmitry Karpeev   if (iis) {
17146a4f0f73SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(IS), iis);CHKERRQ(ierr);
17156a4f0f73SDmitry Karpeev   }
17166a4f0f73SDmitry Karpeev   if (ois) {
17176a4f0f73SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(IS), ois);CHKERRQ(ierr);
17186a4f0f73SDmitry Karpeev   }
17196a4f0f73SDmitry Karpeev   if (iis || ois) {
17206a4f0f73SDmitry Karpeev     for (i = 0; i < osm->n; ++i) {
17216a4f0f73SDmitry Karpeev       if (iis) (*iis)[i] = osm->iis[i];
17226a4f0f73SDmitry Karpeev       if (ois) (*ois)[i] = osm->ois[i];
17236a4f0f73SDmitry Karpeev     }
1724b1a0a8a3SJed Brown   }
1725b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1726b1a0a8a3SJed Brown }
1727b1a0a8a3SJed Brown 
1728b1a0a8a3SJed Brown #undef __FUNCT__
172906b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices"
1730b1a0a8a3SJed Brown /*@C
173106b43e7eSDmitry Karpeev     PCGASMGetSubmatrices - Gets the local submatrices (for this processor
1732b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1733b1a0a8a3SJed Brown 
1734b1a0a8a3SJed Brown     Not Collective
1735b1a0a8a3SJed Brown 
1736b1a0a8a3SJed Brown     Input Parameter:
1737b1a0a8a3SJed Brown .   pc - the preconditioner context
1738b1a0a8a3SJed Brown 
1739b1a0a8a3SJed Brown     Output Parameters:
1740b1a0a8a3SJed Brown +   n   - the number of matrices for this processor (default value = 1)
1741b1a0a8a3SJed Brown -   mat - the matrices
1742b1a0a8a3SJed Brown 
174306b43e7eSDmitry Karpeev     Notes: matrices returned by this routine have the same communicators as the index sets (IS)
174406b43e7eSDmitry Karpeev            used to define subdomains in PCGASMSetSubdomains(), or PETSC_COMM_SELF, if the
17456a4f0f73SDmitry Karpeev            subdomains were defined using PCGASMSetTotalSubdomains().
1746b1a0a8a3SJed Brown     Level: advanced
1747b1a0a8a3SJed Brown 
1748f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi
1749b1a0a8a3SJed Brown 
17506a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomain(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
175106b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains()
1752b1a0a8a3SJed Brown @*/
175306b43e7eSDmitry Karpeev PetscErrorCode  PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1754b1a0a8a3SJed Brown {
1755f746d493SDmitry Karpeev   PC_GASM        *osm;
1756b1a0a8a3SJed Brown   PetscErrorCode ierr;
1757b1a0a8a3SJed Brown   PetscBool      match;
1758b1a0a8a3SJed Brown 
1759b1a0a8a3SJed Brown   PetscFunctionBegin;
1760b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1761b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1762b1a0a8a3SJed Brown   if (mat) PetscValidPointer(mat,3);
1763b1a0a8a3SJed Brown   if (!pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1764251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
176506b43e7eSDmitry Karpeev   if (!match) SETERRQ2(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1766f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1767ab3e923fSDmitry Karpeev   if (n) *n = osm->n;
1768b1a0a8a3SJed Brown   if (mat) *mat = osm->pmat;
1769b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1770b1a0a8a3SJed Brown }
1771