xref: /petsc/src/ksp/pc/impls/gasm/gasm.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
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 */
2965db9045SDmitry Karpeev   PetscBool  use_dm_decomposition;   /* whether DM is allowed to define subdomains */
30f746d493SDmitry Karpeev } PC_GASM;
31b1a0a8a3SJed Brown 
32b1a0a8a3SJed Brown #undef __FUNCT__
3306b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSubdomainView_Private"
3406b43e7eSDmitry Karpeev static PetscErrorCode  PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer)
35af538404SDmitry Karpeev {
36af538404SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
3706b43e7eSDmitry Karpeev   PetscInt       j,nidx;
38af538404SDmitry Karpeev   const PetscInt *idx;
3906b43e7eSDmitry Karpeev   PetscViewer    sviewer;
40af538404SDmitry Karpeev   char           *cidx;
41af538404SDmitry Karpeev   PetscErrorCode ierr;
42af538404SDmitry Karpeev 
43af538404SDmitry Karpeev   PetscFunctionBegin;
44*ce94432eSBarry Smith   if (i < 0 || i > osm->n) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Invalid subdomain %D: must nonnegative and less than %D", i, osm->n);
454bde246dSDmitry Karpeev   /* Inner subdomains. */
464bde246dSDmitry Karpeev   ierr = ISGetLocalSize(osm->iis[i], &nidx);CHKERRQ(ierr);
474bde246dSDmitry Karpeev   /*
484bde246dSDmitry Karpeev    No more than 15 characters per index plus a space.
494bde246dSDmitry Karpeev    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
504bde246dSDmitry Karpeev    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
514bde246dSDmitry Karpeev    For nidx == 0, the whole string 16 '\0'.
524bde246dSDmitry Karpeev    */
534bde246dSDmitry Karpeev   ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr);
544bde246dSDmitry Karpeev   ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr);
554bde246dSDmitry Karpeev   ierr = ISGetIndices(osm->iis[i], &idx);CHKERRQ(ierr);
564bde246dSDmitry Karpeev   for (j = 0; j < nidx; ++j) {
574bde246dSDmitry Karpeev     ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr);
584bde246dSDmitry Karpeev   }
594bde246dSDmitry Karpeev   ierr = ISRestoreIndices(osm->iis[i],&idx);CHKERRQ(ierr);
604bde246dSDmitry Karpeev   ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
614bde246dSDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n");CHKERRQ(ierr);
624bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
634bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
644bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
654bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
664bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
674bde246dSDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
684bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
694bde246dSDmitry Karpeev   ierr = PetscFree(cidx);CHKERRQ(ierr);
704bde246dSDmitry Karpeev   /* Outer subdomains. */
716a4f0f73SDmitry Karpeev   ierr = ISGetLocalSize(osm->ois[i], &nidx);CHKERRQ(ierr);
72eec7e3faSDmitry Karpeev   /*
73eec7e3faSDmitry Karpeev    No more than 15 characters per index plus a space.
74eec7e3faSDmitry Karpeev    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
75eec7e3faSDmitry Karpeev    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
76eec7e3faSDmitry Karpeev    For nidx == 0, the whole string 16 '\0'.
77eec7e3faSDmitry Karpeev    */
78eec7e3faSDmitry Karpeev   ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr);
7906b43e7eSDmitry Karpeev   ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr);
806a4f0f73SDmitry Karpeev   ierr = ISGetIndices(osm->ois[i], &idx);CHKERRQ(ierr);
81af538404SDmitry Karpeev   for (j = 0; j < nidx; ++j) {
82af538404SDmitry Karpeev     ierr = PetscViewerStringSPrintf(sviewer,"%D ", idx[j]);CHKERRQ(ierr);
83af538404SDmitry Karpeev   }
846bf464f9SBarry Smith   ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
856a4f0f73SDmitry Karpeev   ierr = ISRestoreIndices(osm->ois[i],&idx);CHKERRQ(ierr);
866a4f0f73SDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n");CHKERRQ(ierr);
8706b43e7eSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
887b23a99aSBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
89af538404SDmitry Karpeev   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
90af538404SDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
917b23a99aSBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
92af538404SDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
9306b43e7eSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
94af538404SDmitry Karpeev   ierr = PetscFree(cidx);CHKERRQ(ierr);
9506b43e7eSDmitry Karpeev   PetscFunctionReturn(0);
9606b43e7eSDmitry Karpeev }
9706b43e7eSDmitry Karpeev 
9806b43e7eSDmitry Karpeev #undef __FUNCT__
9906b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains"
10006b43e7eSDmitry Karpeev static PetscErrorCode  PCGASMPrintSubdomains(PC pc)
10106b43e7eSDmitry Karpeev {
10206b43e7eSDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
10306b43e7eSDmitry Karpeev   const char     *prefix;
10406b43e7eSDmitry Karpeev   char           fname[PETSC_MAX_PATH_LEN+1];
1054bde246dSDmitry Karpeev   PetscInt       i, l, d, count, gcount, *permutation, *numbering;
10606b43e7eSDmitry Karpeev   PetscBool      found;
1070298fd71SBarry Smith   PetscViewer    viewer, sviewer = NULL;
10806b43e7eSDmitry Karpeev   PetscErrorCode ierr;
10906b43e7eSDmitry Karpeev 
11006b43e7eSDmitry Karpeev   PetscFunctionBegin;
1114bde246dSDmitry Karpeev   ierr = PetscMalloc2(osm->n, PetscInt, &permutation, osm->n, PetscInt, &numbering);CHKERRQ(ierr);
1124bde246dSDmitry Karpeev   for (i = 0; i < osm->n; ++i) permutation[i] = i;
113*ce94432eSBarry Smith   ierr = PetscObjectsGetGlobalNumbering(PetscObjectComm((PetscObject)pc), osm->n, (PetscObject*)osm->ois, &gcount, numbering);CHKERRQ(ierr);
1144bde246dSDmitry Karpeev   ierr = PetscSortIntWithPermutation(osm->n, numbering, permutation);CHKERRQ(ierr);
11506b43e7eSDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
11606b43e7eSDmitry Karpeev   ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr);
11706b43e7eSDmitry Karpeev   if (!found) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
118*ce94432eSBarry Smith   ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr);
1194bde246dSDmitry Karpeev   /*
1204bde246dSDmitry Karpeev    Make sure the viewer has a name. Otherwise this may cause a deadlock or other weird errors when creating a subcomm viewer:
1214bde246dSDmitry Karpeev    the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed collectively on the comm.
1224bde246dSDmitry Karpeev   */
1234bde246dSDmitry Karpeev   ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr);
1244bde246dSDmitry Karpeev   l    = 0;
1254bde246dSDmitry Karpeev   for (count = 0; count < gcount; ++count) {
1264bde246dSDmitry Karpeev     /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
1274bde246dSDmitry Karpeev     if (l<osm->n) {
1284bde246dSDmitry Karpeev       d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
1294bde246dSDmitry Karpeev       if (numbering[d] == count) {
1304bde246dSDmitry Karpeev         ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
1314bde246dSDmitry Karpeev         ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr);
1324bde246dSDmitry Karpeev         ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
1334bde246dSDmitry Karpeev         ++l;
134af538404SDmitry Karpeev       }
1354bde246dSDmitry Karpeev     }
136*ce94432eSBarry Smith     ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
1374bde246dSDmitry Karpeev   }
1384bde246dSDmitry Karpeev   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1397105d80bSDmitry Karpeev   ierr = PetscFree2(permutation,numbering);CHKERRQ(ierr);
140af538404SDmitry Karpeev   PetscFunctionReturn(0);
141af538404SDmitry Karpeev }
142af538404SDmitry Karpeev 
143af538404SDmitry Karpeev 
144af538404SDmitry Karpeev #undef __FUNCT__
145f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM"
146f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer)
147b1a0a8a3SJed Brown {
148f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
149af538404SDmitry Karpeev   const char     *prefix;
150b1a0a8a3SJed Brown   PetscErrorCode ierr;
151af538404SDmitry Karpeev   PetscMPIInt    rank, size;
152b1a0a8a3SJed Brown   PetscInt       i,bsz;
15306b43e7eSDmitry Karpeev   PetscBool      iascii,view_subdomains=PETSC_FALSE;
154b1a0a8a3SJed Brown   PetscViewer    sviewer;
15506b43e7eSDmitry Karpeev   PetscInt       count, l, gcount, *numbering, *permutation;
1566a4f0f73SDmitry Karpeev   char           overlap[256]     = "user-defined overlap";
1576a4f0f73SDmitry Karpeev   char           gsubdomains[256] = "unknown total number of subdomains";
15806b43e7eSDmitry Karpeev   char           lsubdomains[256] = "unknown number of local  subdomains";
15906b43e7eSDmitry Karpeev   char           msubdomains[256] = "unknown max number of local subdomains";
16011aeaf0aSBarry Smith 
161b1a0a8a3SJed Brown   PetscFunctionBegin;
162*ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr);
163*ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr);
16406b43e7eSDmitry Karpeev   ierr = PetscMalloc2(osm->n, PetscInt, &permutation, osm->n, PetscInt, &numbering);CHKERRQ(ierr);
16506b43e7eSDmitry Karpeev   for (i = 0; i < osm->n; ++i) permutation[i] = i;
166*ce94432eSBarry Smith   ierr = PetscObjectsGetGlobalNumbering(PetscObjectComm((PetscObject)pc), osm->n, (PetscObject*)osm->ois, &gcount, numbering);CHKERRQ(ierr);
16706b43e7eSDmitry Karpeev   ierr = PetscSortIntWithPermutation(osm->n, numbering, permutation);CHKERRQ(ierr);
16806b43e7eSDmitry Karpeev 
16906b43e7eSDmitry Karpeev   if (osm->overlap >= 0) {
17006b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %D",osm->overlap);CHKERRQ(ierr);
17106b43e7eSDmitry Karpeev   }
1726a4f0f73SDmitry Karpeev   ierr = PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %D",gcount);CHKERRQ(ierr);
17306b43e7eSDmitry Karpeev   if (osm->N > 0) {
17406b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(lsubdomains, sizeof(gsubdomains), "number of local subdomains = %D",osm->N);CHKERRQ(ierr);
17506b43e7eSDmitry Karpeev   }
17606b43e7eSDmitry Karpeev   if (osm->nmax > 0) {
17706b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr);
17806b43e7eSDmitry Karpeev   }
17906b43e7eSDmitry Karpeev 
180af538404SDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
1810298fd71SBarry Smith   ierr = PetscOptionsGetBool(prefix,"-pc_gasm_view_subdomains",&view_subdomains,NULL);CHKERRQ(ierr);
18206b43e7eSDmitry Karpeev 
18306b43e7eSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
184b1a0a8a3SJed Brown   if (iascii) {
18506b43e7eSDmitry Karpeev     /*
18606b43e7eSDmitry Karpeev      Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer:
18706b43e7eSDmitry Karpeev      the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed
18806b43e7eSDmitry Karpeev      collectively on the comm.
18906b43e7eSDmitry Karpeev      */
19006b43e7eSDmitry Karpeev     ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr);
19106b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Generalized additive Schwarz:\n");CHKERRQ(ierr);
19206b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Restriction/interpolation type: %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr);
19306b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",overlap);CHKERRQ(ierr);
19406b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",gsubdomains);CHKERRQ(ierr);
19506b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",lsubdomains);CHKERRQ(ierr);
19606b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",msubdomains);CHKERRQ(ierr);
1977b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
198eec7e3faSDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d:%d] number of locally-supported subdomains = %D\n",(int)rank,(int)size,osm->n);CHKERRQ(ierr);
199af538404SDmitry Karpeev     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
2007b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
2016a4f0f73SDmitry Karpeev     /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */
20206b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Subdomain solver info is as follows:\n");CHKERRQ(ierr);
203b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
204b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
20506b43e7eSDmitry Karpeev     /* Make sure that everybody waits for the banner to be printed. */
206*ce94432eSBarry Smith     ierr = MPI_Barrier(PetscObjectComm((PetscObject)viewer));CHKERRQ(ierr);
2074bde246dSDmitry Karpeev     /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
20806b43e7eSDmitry Karpeev     l = 0;
20906b43e7eSDmitry Karpeev     for (count = 0; count < gcount; ++count) {
21006b43e7eSDmitry Karpeev       PetscMPIInt srank, ssize;
21106b43e7eSDmitry Karpeev       if (l<osm->n) {
21206b43e7eSDmitry Karpeev         PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
21306b43e7eSDmitry Karpeev         if (numbering[d] == count) {
2146a4f0f73SDmitry Karpeev           ierr = MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize);CHKERRQ(ierr);
2156a4f0f73SDmitry Karpeev           ierr = MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank);CHKERRQ(ierr);
2166a4f0f73SDmitry Karpeev           ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
2176a4f0f73SDmitry Karpeev           ierr = ISGetLocalSize(osm->ois[d],&bsz);CHKERRQ(ierr);
2187b23a99aSBarry Smith           ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_TRUE);CHKERRQ(ierr);
21906b43e7eSDmitry 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);
2206a4f0f73SDmitry Karpeev           ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
2216a4f0f73SDmitry Karpeev           ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_FALSE);CHKERRQ(ierr);
22206b43e7eSDmitry Karpeev           if (view_subdomains) {
22306b43e7eSDmitry Karpeev             ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr);
22406b43e7eSDmitry Karpeev           }
2256a4f0f73SDmitry Karpeev           if (!pc->setupcalled) {
2266a4f0f73SDmitry Karpeev             PetscViewerASCIIPrintf(sviewer, "Solver not set up yet: PCSetUp() not yet called\n");CHKERRQ(ierr);
2272fa5cd67SKarl Rupp           } else {
22806b43e7eSDmitry Karpeev             ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr);
2296a4f0f73SDmitry Karpeev           }
23006b43e7eSDmitry Karpeev           ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
2317b23a99aSBarry Smith           ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
2326a4f0f73SDmitry Karpeev           ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
23306b43e7eSDmitry Karpeev           ++l;
234b1a0a8a3SJed Brown         }
23506b43e7eSDmitry Karpeev       }
236*ce94432eSBarry Smith       ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
237b1a0a8a3SJed Brown     }
238b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
239b1a0a8a3SJed Brown   }
24006b43e7eSDmitry Karpeev   ierr = PetscFree2(permutation,numbering);CHKERRQ(ierr);
241b1a0a8a3SJed Brown   PetscFunctionReturn(0);
242b1a0a8a3SJed Brown }
243b1a0a8a3SJed Brown 
244af538404SDmitry Karpeev 
245af538404SDmitry Karpeev 
246af538404SDmitry Karpeev 
247af538404SDmitry Karpeev 
248b1a0a8a3SJed Brown #undef __FUNCT__
249f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM"
250f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc)
251b1a0a8a3SJed Brown {
252f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
253b1a0a8a3SJed Brown   PetscErrorCode ierr;
254b1a0a8a3SJed Brown   PetscBool      symset,flg;
25588389755SJed Brown   PetscInt       i;
256c10bc1a1SDmitry Karpeev   PetscMPIInt    rank, size;
257b1a0a8a3SJed Brown   MatReuse       scall = MAT_REUSE_MATRIX;
258b1a0a8a3SJed Brown   KSP            ksp;
259b1a0a8a3SJed Brown   PC             subpc;
260b1a0a8a3SJed Brown   const char     *prefix,*pprefix;
261f746d493SDmitry Karpeev   Vec            x,y;
2626a4f0f73SDmitry Karpeev   PetscInt       oni;       /* Number of indices in the i-th local outer subdomain.               */
2636a4f0f73SDmitry Karpeev   const PetscInt *oidxi;    /* Indices from the i-th subdomain local outer subdomain.             */
2646a4f0f73SDmitry Karpeev   PetscInt       on;        /* Number of indices in the disjoint union of local outer subdomains. */
2656a4f0f73SDmitry Karpeev   PetscInt       *oidx;     /* Indices in the disjoint union of local outer subdomains. */
2666a4f0f73SDmitry Karpeev   IS             gois;      /* Disjoint union the global indices of outer subdomains.             */
2676a4f0f73SDmitry Karpeev   IS             goid;      /* Identity IS of the size of the disjoint union of outer subdomains. */
268f746d493SDmitry Karpeev   PetscScalar    *gxarray, *gyarray;
2692fa5cd67SKarl Rupp   PetscInt       gofirst;   /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */
2700298fd71SBarry Smith   DM             *subdomain_dm = NULL;
271b1a0a8a3SJed Brown 
272b1a0a8a3SJed Brown   PetscFunctionBegin;
273*ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
274*ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
275b1a0a8a3SJed Brown   if (!pc->setupcalled) {
276b1a0a8a3SJed Brown 
277b1a0a8a3SJed Brown     if (!osm->type_set) {
278b1a0a8a3SJed Brown       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
2792fa5cd67SKarl Rupp       if (symset && flg) osm->type = PC_GASM_BASIC;
280b1a0a8a3SJed Brown     }
281b1a0a8a3SJed Brown 
28206b43e7eSDmitry Karpeev     /*
28306b43e7eSDmitry Karpeev      If subdomains have been set, then the local number of subdomains, osm->n, is NOT PETSC_DECIDE and is at least 1.
28406b43e7eSDmitry 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.
28506b43e7eSDmitry Karpeev      */
28606b43e7eSDmitry Karpeev     if (osm->n == PETSC_DECIDE) {
28769ca1f37SDmitry Karpeev       /* no subdomains given */
28865db9045SDmitry Karpeev       /* try pc->dm first, if allowed */
28965db9045SDmitry Karpeev       if (osm->use_dm_decomposition && pc->dm) {
290a35b7b57SDmitry Karpeev         PetscInt  num_subdomains, d;
291a35b7b57SDmitry Karpeev         char      **subdomain_names;
292a35b7b57SDmitry Karpeev         IS        *inner_subdomain_is, *outer_subdomain_is;
293a35b7b57SDmitry Karpeev         ierr = DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);CHKERRQ(ierr);
294a35b7b57SDmitry Karpeev         if (num_subdomains) {
295a35b7b57SDmitry Karpeev           ierr = PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);CHKERRQ(ierr);
29669ca1f37SDmitry Karpeev         }
297a35b7b57SDmitry Karpeev         for (d = 0; d < num_subdomains; ++d) {
298a35b7b57SDmitry Karpeev           if (subdomain_names)    {ierr = PetscFree(subdomain_names[d]);CHKERRQ(ierr);}
299a35b7b57SDmitry Karpeev           if (inner_subdomain_is) {ierr = ISDestroy(&inner_subdomain_is[d]);CHKERRQ(ierr);}
300a35b7b57SDmitry Karpeev           if (outer_subdomain_is) {ierr = ISDestroy(&outer_subdomain_is[d]);CHKERRQ(ierr);}
301fdc48646SDmitry Karpeev         }
302a35b7b57SDmitry Karpeev         ierr = PetscFree(subdomain_names);CHKERRQ(ierr);
303a35b7b57SDmitry Karpeev         ierr = PetscFree(inner_subdomain_is);CHKERRQ(ierr);
304a35b7b57SDmitry Karpeev         ierr = PetscFree(outer_subdomain_is);CHKERRQ(ierr);
305fdc48646SDmitry Karpeev       }
30606b43e7eSDmitry Karpeev       if (osm->n == PETSC_DECIDE) { /* still no subdomains; use one per processor */
30744fe18e5SDmitry Karpeev         osm->nmax = osm->n = 1;
308*ce94432eSBarry Smith         ierr      = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
309f746d493SDmitry Karpeev         osm->N    = size;
310fdc48646SDmitry Karpeev       }
31106b43e7eSDmitry Karpeev     }
312a35b7b57SDmitry Karpeev     if (!osm->iis) {
313a35b7b57SDmitry Karpeev       /*
314a35b7b57SDmitry Karpeev        The local number of subdomains was set just above, or in PCGASMSetTotalSubdomains(), or in PCGASMSetSubdomains(),
315a35b7b57SDmitry Karpeev        but the actual subdomains have not been supplied (in PCGASMSetSubdomains()).
316a35b7b57SDmitry Karpeev        We create the requisite number of inner subdomains on PETSC_COMM_SELF (for now).
317a35b7b57SDmitry Karpeev        */
318a35b7b57SDmitry Karpeev       ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->overlap,osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr);
319a35b7b57SDmitry Karpeev     }
32006b43e7eSDmitry Karpeev     if (osm->N == PETSC_DECIDE) {
3216ac3741eSJed Brown       struct {PetscInt max,sum;} inwork,outwork;
322f746d493SDmitry Karpeev       /* determine global number of subdomains and the max number of local subdomains */
3236ac3741eSJed Brown       inwork.max = osm->n;
3246ac3741eSJed Brown       inwork.sum = osm->n;
325*ce94432eSBarry Smith       ierr       = MPI_Allreduce(&inwork,&outwork,1,MPIU_2INT,PetscMaxSum_Op,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
3266ac3741eSJed Brown       osm->nmax  = outwork.max;
3276ac3741eSJed Brown       osm->N     = outwork.sum;
328b1a0a8a3SJed Brown     }
3296a4f0f73SDmitry Karpeev 
330b1a0a8a3SJed Brown     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
331b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
3320298fd71SBarry Smith     ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&flg,NULL);CHKERRQ(ierr);
333f746d493SDmitry Karpeev     if (flg) { ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); }
334b1a0a8a3SJed Brown 
335b1a0a8a3SJed Brown     if (osm->sort_indices) {
336f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
3376a4f0f73SDmitry Karpeev         ierr = ISSort(osm->ois[i]);CHKERRQ(ierr);
3386a4f0f73SDmitry Karpeev         ierr = ISSort(osm->iis[i]);CHKERRQ(ierr);
339b1a0a8a3SJed Brown       }
340b1a0a8a3SJed Brown     }
3416a4f0f73SDmitry Karpeev     /*
3426a4f0f73SDmitry Karpeev      Merge the ISs, create merged vectors and restrictions.
3436a4f0f73SDmitry Karpeev      */
3446a4f0f73SDmitry Karpeev     /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */
3456a4f0f73SDmitry Karpeev     on = 0;
346f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
3476a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
3486a4f0f73SDmitry Karpeev       on  += oni;
349f746d493SDmitry Karpeev     }
3506a4f0f73SDmitry Karpeev     ierr = PetscMalloc(on*sizeof(PetscInt), &oidx);CHKERRQ(ierr);
3516a4f0f73SDmitry Karpeev     on   = 0;
352f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
3536a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
3546a4f0f73SDmitry Karpeev       ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr);
3556a4f0f73SDmitry Karpeev       ierr = PetscMemcpy(oidx+on, oidxi, sizeof(PetscInt)*oni);CHKERRQ(ierr);
3566a4f0f73SDmitry Karpeev       ierr = ISRestoreIndices(osm->ois[i], &oidxi);CHKERRQ(ierr);
3576a4f0f73SDmitry Karpeev       on  += oni;
358f746d493SDmitry Karpeev     }
3596a4f0f73SDmitry Karpeev     ierr = ISCreateGeneral(((PetscObject)(pc))->comm, on, oidx, PETSC_OWN_POINTER, &gois);CHKERRQ(ierr);
360f746d493SDmitry Karpeev     ierr = MatGetVecs(pc->pmat,&x,&y);CHKERRQ(ierr);
361*ce94432eSBarry Smith     ierr = VecCreateMPI(PetscObjectComm((PetscObject)pc), on, PETSC_DECIDE, &osm->gx);CHKERRQ(ierr);
362f746d493SDmitry Karpeev     ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr);
3630298fd71SBarry Smith     ierr = VecGetOwnershipRange(osm->gx, &gofirst, NULL);CHKERRQ(ierr);
364*ce94432eSBarry Smith     ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),on,gofirst,1, &goid);CHKERRQ(ierr);
3656a4f0f73SDmitry Karpeev     ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr);
3666a4f0f73SDmitry Karpeev     ierr = VecDestroy(&x);CHKERRQ(ierr);
3677105d80bSDmitry Karpeev     ierr = ISDestroy(&gois);CHKERRQ(ierr);
3682fa5cd67SKarl Rupp 
3696a4f0f73SDmitry Karpeev     /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */
3702fa5cd67SKarl Rupp     {
3712fa5cd67SKarl Rupp       PetscInt ini;           /* Number of indices the i-th a local inner subdomain. */
3726a4f0f73SDmitry Karpeev       PetscInt in;            /* Number of indices in the disjoint uniont of local inner subdomains. */
3736a4f0f73SDmitry Karpeev       PetscInt *iidx;         /* Global indices in the merged local inner subdomain. */
3746a4f0f73SDmitry Karpeev       PetscInt *ioidx;        /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
3756a4f0f73SDmitry Karpeev       IS       giis;          /* IS for the disjoint union of inner subdomains. */
3766a4f0f73SDmitry Karpeev       IS       giois;         /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
377f746d493SDmitry Karpeev       /**/
3786a4f0f73SDmitry Karpeev       in = 0;
379f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
3806a4f0f73SDmitry Karpeev         ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr);
3816a4f0f73SDmitry Karpeev         in  += ini;
382f746d493SDmitry Karpeev       }
3836a4f0f73SDmitry Karpeev       ierr = PetscMalloc(in*sizeof(PetscInt), &iidx);CHKERRQ(ierr);
3846a4f0f73SDmitry Karpeev       ierr = PetscMalloc(in*sizeof(PetscInt), &ioidx);CHKERRQ(ierr);
3850298fd71SBarry Smith       ierr = VecGetOwnershipRange(osm->gx,&gofirst, NULL);CHKERRQ(ierr);
3866a4f0f73SDmitry Karpeev       in   = 0;
3876a4f0f73SDmitry Karpeev       on   = 0;
388f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
3896a4f0f73SDmitry Karpeev         const PetscInt         *iidxi; /* Global indices of the i-th local inner subdomain. */
3906a4f0f73SDmitry Karpeev         ISLocalToGlobalMapping ltogi; /* Map from global to local indices of the i-th outer local subdomain. */
3916a4f0f73SDmitry Karpeev         PetscInt               *ioidxi; /* Local indices of the i-th local inner subdomain within the local outer subdomain. */
3926a4f0f73SDmitry Karpeev         PetscInt               ioni;  /* Number of indices in ioidxi; if ioni != ini the inner subdomain is not a subdomain of the outer subdomain (error). */
3936a4f0f73SDmitry Karpeev         PetscInt               k;
3946a4f0f73SDmitry Karpeev         ierr   = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr);
3956a4f0f73SDmitry Karpeev         ierr   = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
3966a4f0f73SDmitry Karpeev         ierr   = ISGetIndices(osm->iis[i],&iidxi);CHKERRQ(ierr);
3976a4f0f73SDmitry Karpeev         ierr   = PetscMemcpy(iidx+in, iidxi, sizeof(PetscInt)*ini);CHKERRQ(ierr);
3986a4f0f73SDmitry Karpeev         ierr   = ISLocalToGlobalMappingCreateIS(osm->ois[i],&ltogi);CHKERRQ(ierr);
3996a4f0f73SDmitry Karpeev         ioidxi = ioidx+in;
4006a4f0f73SDmitry Karpeev         ierr   = ISGlobalToLocalMappingApply(ltogi,IS_GTOLM_DROP,ini,iidxi,&ioni,ioidxi);CHKERRQ(ierr);
4016a4f0f73SDmitry Karpeev         ierr   = ISLocalToGlobalMappingDestroy(&ltogi);CHKERRQ(ierr);
4026a4f0f73SDmitry Karpeev         ierr   = ISRestoreIndices(osm->iis[i], &iidxi);CHKERRQ(ierr);
4036a4f0f73SDmitry 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);
4042fa5cd67SKarl Rupp         for (k = 0; k < ini; ++k) ioidxi[k] += gofirst+on;
4056a4f0f73SDmitry Karpeev         in += ini;
4066a4f0f73SDmitry Karpeev         on += oni;
407f746d493SDmitry Karpeev       }
408*ce94432eSBarry Smith       ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, iidx,  PETSC_OWN_POINTER, &giis);CHKERRQ(ierr);
409*ce94432eSBarry Smith       ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, ioidx, PETSC_OWN_POINTER, &giois);CHKERRQ(ierr);
4106a4f0f73SDmitry Karpeev       ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr);
4116a4f0f73SDmitry Karpeev       ierr = VecDestroy(&y);CHKERRQ(ierr);
4126a4f0f73SDmitry Karpeev       ierr = ISDestroy(&giis);CHKERRQ(ierr);
4136a4f0f73SDmitry Karpeev       ierr = ISDestroy(&giois);CHKERRQ(ierr);
414b1a0a8a3SJed Brown     }
4156a4f0f73SDmitry Karpeev     ierr = ISDestroy(&goid);CHKERRQ(ierr);
4162fa5cd67SKarl Rupp 
417f746d493SDmitry Karpeev     /* Create the subdomain work vectors. */
418f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->x);CHKERRQ(ierr);
419f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->y);CHKERRQ(ierr);
420f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr);
421f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr);
4226a4f0f73SDmitry Karpeev     for (i=0, on=0; i<osm->n; ++i, on += oni) {
4236a4f0f73SDmitry Karpeev       PetscInt oNi;
4246a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
4256a4f0f73SDmitry Karpeev       ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr);
4266a4f0f73SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr);
4276a4f0f73SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr);
428b1a0a8a3SJed Brown     }
429f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr);
430f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr);
431b1a0a8a3SJed Brown     /* Create the local solvers */
432f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(KSP*),&osm->ksp);CHKERRQ(ierr);
43344fe18e5SDmitry Karpeev     for (i=0; i<osm->n; i++) { /* KSPs are local */
4346a4f0f73SDmitry Karpeev       ierr        = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr);
435b1a0a8a3SJed Brown       ierr        = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr);
436b1a0a8a3SJed Brown       ierr        = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
437b1a0a8a3SJed Brown       ierr        = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
438b1a0a8a3SJed Brown       ierr        = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
439b1a0a8a3SJed Brown       ierr        = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
440b1a0a8a3SJed Brown       ierr        = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
441b1a0a8a3SJed Brown       ierr        = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
442b1a0a8a3SJed Brown       osm->ksp[i] = ksp;
443b1a0a8a3SJed Brown     }
444b1a0a8a3SJed Brown     scall = MAT_INITIAL_MATRIX;
445b1a0a8a3SJed Brown 
4462fa5cd67SKarl Rupp   } else { /*if (!pc->setupcalled)*/
447b1a0a8a3SJed Brown     /*
448b1a0a8a3SJed Brown        Destroy the blocks from the previous iteration
449b1a0a8a3SJed Brown     */
450b1a0a8a3SJed Brown     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
451f746d493SDmitry Karpeev       ierr  = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
452b1a0a8a3SJed Brown       scall = MAT_INITIAL_MATRIX;
453b1a0a8a3SJed Brown     }
454b1a0a8a3SJed Brown   }
455b1a0a8a3SJed Brown 
456b1a0a8a3SJed Brown   /*
457f746d493SDmitry Karpeev      Extract out the submatrices.
458b1a0a8a3SJed Brown   */
45981a5c4aaSDmitry Karpeev   if (size > 1) {
4606a4f0f73SDmitry Karpeev     ierr = MatGetSubMatricesParallel(pc->pmat,osm->n,osm->ois, osm->ois,scall,&osm->pmat);CHKERRQ(ierr);
4612fa5cd67SKarl Rupp   } else {
4626a4f0f73SDmitry Karpeev     ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->ois, osm->ois,scall,&osm->pmat);CHKERRQ(ierr);
46381a5c4aaSDmitry Karpeev   }
464b1a0a8a3SJed Brown   if (scall == MAT_INITIAL_MATRIX) {
465b1a0a8a3SJed Brown     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
466f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
467b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,osm->pmat[i]);CHKERRQ(ierr);
468b1a0a8a3SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
469b1a0a8a3SJed Brown     }
470b1a0a8a3SJed Brown   }
471b1a0a8a3SJed Brown 
472b1a0a8a3SJed Brown   /* Return control to the user so that the submatrices can be modified (e.g., to apply
473b1a0a8a3SJed Brown      different boundary conditions for the submatrices than for the global problem) */
4746a4f0f73SDmitry Karpeev   ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
475b1a0a8a3SJed Brown 
476b1a0a8a3SJed Brown   /*
4776a4f0f73SDmitry Karpeev      Loop over submatrices putting them into local ksps
478b1a0a8a3SJed Brown   */
479f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
480b1a0a8a3SJed Brown     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr);
481b1a0a8a3SJed Brown     if (!pc->setupcalled) {
482b1a0a8a3SJed Brown       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
483b1a0a8a3SJed Brown     }
484b1a0a8a3SJed Brown   }
485b1a0a8a3SJed Brown   PetscFunctionReturn(0);
486b1a0a8a3SJed Brown }
487b1a0a8a3SJed Brown 
488b1a0a8a3SJed Brown #undef __FUNCT__
489f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM"
490f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
491b1a0a8a3SJed Brown {
492f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
493b1a0a8a3SJed Brown   PetscErrorCode ierr;
494b1a0a8a3SJed Brown   PetscInt       i;
495b1a0a8a3SJed Brown 
496b1a0a8a3SJed Brown   PetscFunctionBegin;
497f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
498b1a0a8a3SJed Brown     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
499b1a0a8a3SJed Brown   }
500b1a0a8a3SJed Brown   PetscFunctionReturn(0);
501b1a0a8a3SJed Brown }
502b1a0a8a3SJed Brown 
503b1a0a8a3SJed Brown #undef __FUNCT__
504f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM"
505f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y)
506b1a0a8a3SJed Brown {
507f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
508b1a0a8a3SJed Brown   PetscErrorCode ierr;
509f746d493SDmitry Karpeev   PetscInt       i;
510b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
511b1a0a8a3SJed Brown 
512b1a0a8a3SJed Brown   PetscFunctionBegin;
513b1a0a8a3SJed Brown   /*
5146a4f0f73SDmitry Karpeev      Support for limiting the restriction or interpolation only to the inner
515b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
516b1a0a8a3SJed Brown   */
517f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
518b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
519f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
5206a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5212fa5cd67SKarl Rupp   } else {
5226a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
523b1a0a8a3SJed Brown   }
5246a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr);
5256a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
5266a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5272fa5cd67SKarl Rupp   } else {
5286a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5296a4f0f73SDmitry Karpeev   }
53086b96d36SDmitry Karpeev   /* do the subdomain solves */
53186b96d36SDmitry Karpeev   for (i=0; i<osm->n; ++i) {
532b1a0a8a3SJed Brown     ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
533b1a0a8a3SJed Brown   }
5346a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(y);CHKERRQ(ierr);
5356a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
5366a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5376a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);  PetscFunctionReturn(0);
5382fa5cd67SKarl Rupp   } else {
5396a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5406a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);  PetscFunctionReturn(0);
5416a4f0f73SDmitry Karpeev   }
542b1a0a8a3SJed Brown }
543b1a0a8a3SJed Brown 
544b1a0a8a3SJed Brown #undef __FUNCT__
545f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM"
546f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y)
547b1a0a8a3SJed Brown {
548f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
549b1a0a8a3SJed Brown   PetscErrorCode ierr;
550f746d493SDmitry Karpeev   PetscInt       i;
551b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
552b1a0a8a3SJed Brown 
553b1a0a8a3SJed Brown   PetscFunctionBegin;
554b1a0a8a3SJed Brown   /*
555b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
556b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
557b1a0a8a3SJed Brown 
558f746d493SDmitry Karpeev      Note: these are reversed from the PCApply_GASM() because we are applying the
559b1a0a8a3SJed Brown      transpose of the three terms
560b1a0a8a3SJed Brown   */
561f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
562b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
563f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
5646a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5652fa5cd67SKarl Rupp   } else {
5666a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
567b1a0a8a3SJed Brown   }
5686a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr);
5696a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
5706a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5712fa5cd67SKarl Rupp   } else {
5726a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5736a4f0f73SDmitry Karpeev   }
574b1a0a8a3SJed Brown   /* do the local solves */
575ab3e923fSDmitry 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. */
576b1a0a8a3SJed Brown     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
577b1a0a8a3SJed Brown   }
5786a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(y);CHKERRQ(ierr);
5796a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
5806a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5816a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5822fa5cd67SKarl Rupp   } else {
5836a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5846a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5856a4f0f73SDmitry Karpeev   }
586b1a0a8a3SJed Brown   PetscFunctionReturn(0);
587b1a0a8a3SJed Brown }
588b1a0a8a3SJed Brown 
589b1a0a8a3SJed Brown #undef __FUNCT__
590a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM"
591a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc)
592b1a0a8a3SJed Brown {
593f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
594b1a0a8a3SJed Brown   PetscErrorCode ierr;
595b1a0a8a3SJed Brown   PetscInt       i;
596b1a0a8a3SJed Brown 
597b1a0a8a3SJed Brown   PetscFunctionBegin;
598b1a0a8a3SJed Brown   if (osm->ksp) {
599f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
600a06653b4SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
601b1a0a8a3SJed Brown     }
602b1a0a8a3SJed Brown   }
603b1a0a8a3SJed Brown   if (osm->pmat) {
604f746d493SDmitry Karpeev     if (osm->n > 0) {
605ab3e923fSDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
606b1a0a8a3SJed Brown     }
607b1a0a8a3SJed Brown   }
608d34fcf5fSBarry Smith   if (osm->x) {
609f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
610fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
611fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
612b1a0a8a3SJed Brown     }
613d34fcf5fSBarry Smith   }
6146bf464f9SBarry Smith   ierr = VecDestroy(&osm->gx);CHKERRQ(ierr);
6156bf464f9SBarry Smith   ierr = VecDestroy(&osm->gy);CHKERRQ(ierr);
616ab3e923fSDmitry Karpeev 
6176a4f0f73SDmitry Karpeev   ierr     = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr);
6186a4f0f73SDmitry Karpeev   ierr     = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr);
6196a4f0f73SDmitry Karpeev   ierr     = PCGASMDestroySubdomains(osm->n,osm->ois,osm->iis);CHKERRQ(ierr);
6206a4f0f73SDmitry Karpeev   osm->ois = 0;
6216a4f0f73SDmitry Karpeev   osm->iis = 0;
622a06653b4SBarry Smith   PetscFunctionReturn(0);
623a06653b4SBarry Smith }
624a06653b4SBarry Smith 
625a06653b4SBarry Smith #undef __FUNCT__
626a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM"
627a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc)
628a06653b4SBarry Smith {
629a06653b4SBarry Smith   PC_GASM        *osm = (PC_GASM*)pc->data;
630a06653b4SBarry Smith   PetscErrorCode ierr;
631a06653b4SBarry Smith   PetscInt       i;
632a06653b4SBarry Smith 
633a06653b4SBarry Smith   PetscFunctionBegin;
634a06653b4SBarry Smith   ierr = PCReset_GASM(pc);CHKERRQ(ierr);
635a06653b4SBarry Smith   if (osm->ksp) {
636a06653b4SBarry Smith     for (i=0; i<osm->n; i++) {
6376bf464f9SBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
638a06653b4SBarry Smith     }
639a06653b4SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
640a06653b4SBarry Smith   }
641a06653b4SBarry Smith   ierr = PetscFree(osm->x);CHKERRQ(ierr);
642a06653b4SBarry Smith   ierr = PetscFree(osm->y);CHKERRQ(ierr);
643c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
644b1a0a8a3SJed Brown   PetscFunctionReturn(0);
645b1a0a8a3SJed Brown }
646b1a0a8a3SJed Brown 
647b1a0a8a3SJed Brown #undef __FUNCT__
648f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM"
649a6dfd86eSKarl Rupp static PetscErrorCode PCSetFromOptions_GASM(PC pc)
650a6dfd86eSKarl Rupp {
651f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
652b1a0a8a3SJed Brown   PetscErrorCode ierr;
653b1a0a8a3SJed Brown   PetscInt       blocks,ovl;
654b1a0a8a3SJed Brown   PetscBool      symset,flg;
655f746d493SDmitry Karpeev   PCGASMType     gasmtype;
656b1a0a8a3SJed Brown 
657b1a0a8a3SJed Brown   PetscFunctionBegin;
658b1a0a8a3SJed Brown   /* set the type to symmetric if matrix is symmetric */
659b1a0a8a3SJed Brown   if (!osm->type_set && pc->pmat) {
660b1a0a8a3SJed Brown     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
6612fa5cd67SKarl Rupp     if (symset && flg) osm->type = PC_GASM_BASIC;
662b1a0a8a3SJed Brown   }
66344fe18e5SDmitry Karpeev   ierr = PetscOptionsHead("Generalized additive Schwarz options");CHKERRQ(ierr);
6646a4f0f73SDmitry Karpeev   ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
66565db9045SDmitry Karpeev   if (flg) {
6666a4f0f73SDmitry Karpeev     osm->create_local = PETSC_TRUE;
66765db9045SDmitry Karpeev     ierr = PetscOptionsBool("-pc_gasm_subdomains_create_local","Whether to make autocreated subdomains local (true by default)","PCGASMSetTotalSubdomains",osm->create_local,&osm->create_local,NULL);CHKERRQ(ierr);
66865db9045SDmitry Karpeev     ierr = PCGASMSetTotalSubdomains(pc,blocks,osm->create_local);CHKERRQ(ierr);
66965db9045SDmitry Karpeev     osm->use_dm_decomposition = PETSC_FALSE;
67065db9045SDmitry Karpeev   }
67106b43e7eSDmitry Karpeev   ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
67265db9045SDmitry Karpeev   if (flg) {
67365db9045SDmitry Karpeev     ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr);
67465db9045SDmitry Karpeev     osm->use_dm_decomposition = PETSC_FALSE;
67565db9045SDmitry Karpeev   }
676b1a0a8a3SJed Brown   flg  = PETSC_FALSE;
677f746d493SDmitry Karpeev   ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr);
678f746d493SDmitry Karpeev   if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr); }
679b1a0a8a3SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
680b1a0a8a3SJed Brown   PetscFunctionReturn(0);
681b1a0a8a3SJed Brown }
682b1a0a8a3SJed Brown 
683b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/
684b1a0a8a3SJed Brown 
685b1a0a8a3SJed Brown EXTERN_C_BEGIN
686b1a0a8a3SJed Brown #undef __FUNCT__
68706b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM"
688a35b7b57SDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[])
689b1a0a8a3SJed Brown {
690f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
691b1a0a8a3SJed Brown   PetscErrorCode ierr;
692b1a0a8a3SJed Brown   PetscInt       i;
693b1a0a8a3SJed Brown 
694b1a0a8a3SJed Brown   PetscFunctionBegin;
69569ca1f37SDmitry Karpeev   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, n = %D",n);
696*ce94432eSBarry Smith   if (pc->setupcalled && (n != osm->n || iis || ois)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp().");
697b1a0a8a3SJed Brown 
698b1a0a8a3SJed Brown   if (!pc->setupcalled) {
69993c1ec32SDmitry Karpeev     osm->n   = n;
7006a4f0f73SDmitry Karpeev     osm->ois = 0;
7016a4f0f73SDmitry Karpeev     osm->iis = 0;
702a35b7b57SDmitry Karpeev     if (ois) {
703a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);}
704b1a0a8a3SJed Brown     }
705a35b7b57SDmitry Karpeev     if (iis) {
706a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);}
707b1a0a8a3SJed Brown     }
7086a4f0f73SDmitry Karpeev     ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr);
709a35b7b57SDmitry Karpeev     if (ois) {
7106a4f0f73SDmitry Karpeev       ierr = PetscMalloc(n*sizeof(IS),&osm->ois);CHKERRQ(ierr);
7112fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->ois[i] = ois[i];
7126a4f0f73SDmitry Karpeev       /* Flag indicating that the user has set outer subdomains, so PCGASM should not increase their size. */
713b1a0a8a3SJed Brown       osm->overlap = -1;
714a35b7b57SDmitry Karpeev       if (!iis) {
7156a4f0f73SDmitry Karpeev         ierr = PetscMalloc(n*sizeof(IS),&osm->iis);CHKERRQ(ierr);
716a35b7b57SDmitry Karpeev         for (i=0; i<n; i++) {
717a35b7b57SDmitry Karpeev           for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);}
718a35b7b57SDmitry Karpeev           osm->iis[i] = ois[i];
719a35b7b57SDmitry Karpeev         }
720a35b7b57SDmitry Karpeev       }
721a35b7b57SDmitry Karpeev     }
722a35b7b57SDmitry Karpeev     if (iis) {
723a35b7b57SDmitry Karpeev       ierr = PetscMalloc(n*sizeof(IS),&osm->iis);CHKERRQ(ierr);
7242fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->iis[i] = iis[i];
725a35b7b57SDmitry Karpeev       if (!ois) {
726a35b7b57SDmitry Karpeev         ierr = PetscMalloc(n*sizeof(IS),&osm->ois);CHKERRQ(ierr);
727a35b7b57SDmitry Karpeev         for (i=0; i<n; i++) {
728a35b7b57SDmitry Karpeev           for (i=0; i<n; i++) {
729a35b7b57SDmitry Karpeev             ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);
730a35b7b57SDmitry Karpeev             osm->ois[i] = iis[i];
731a35b7b57SDmitry Karpeev           }
732a35b7b57SDmitry Karpeev         }
733a35b7b57SDmitry Karpeev         if (osm->overlap > 0) {
734a35b7b57SDmitry Karpeev           /* Extend the "overlapping" regions by a number of steps */
735a35b7b57SDmitry Karpeev           ierr = MatIncreaseOverlap(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr);
736a35b7b57SDmitry Karpeev         }
737a35b7b57SDmitry Karpeev       }
738b1a0a8a3SJed Brown     }
739b1a0a8a3SJed Brown   }
740b1a0a8a3SJed Brown   PetscFunctionReturn(0);
741b1a0a8a3SJed Brown }
742b1a0a8a3SJed Brown EXTERN_C_END
743b1a0a8a3SJed Brown 
744b1a0a8a3SJed Brown EXTERN_C_BEGIN
745b1a0a8a3SJed Brown #undef __FUNCT__
7466a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains_GASM"
7470adebc6cSBarry Smith PetscErrorCode  PCGASMSetTotalSubdomains_GASM(PC pc,PetscInt N, PetscBool create_local)
7480adebc6cSBarry Smith {
749f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
750b1a0a8a3SJed Brown   PetscErrorCode ierr;
751b1a0a8a3SJed Brown   PetscMPIInt    rank,size;
752b1a0a8a3SJed Brown   PetscInt       n;
7536a4f0f73SDmitry Karpeev   PetscInt       Nmin, Nmax;
7545fd66863SKarl Rupp 
755b1a0a8a3SJed Brown   PetscFunctionBegin;
756*ce94432eSBarry Smith   if (!create_local) SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "No suppor for autocreation of nonlocal subdomains.");
757*ce94432eSBarry Smith   if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be > 0, N = %D",N);
758*ce94432eSBarry Smith   ierr = MPI_Allreduce(&N,&Nmin,1,MPIU_INT,MPIU_MIN,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
759*ce94432eSBarry Smith   ierr = MPI_Allreduce(&N,&Nmax,1,MPIU_INT,MPIU_MAX,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
760*ce94432eSBarry Smith   if (Nmin != Nmax) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "All processors must use the same number of subdomains.  min(N) = %D != %D = max(N)", Nmin, Nmax);
761b1a0a8a3SJed Brown 
7626a4f0f73SDmitry Karpeev   osm->create_local = create_local;
763b1a0a8a3SJed Brown   /*
764b1a0a8a3SJed Brown      Split the subdomains equally among all processors
765b1a0a8a3SJed Brown   */
766*ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
767*ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
768b1a0a8a3SJed Brown   n    = N/size + ((N % size) > rank);
7696a4f0f73SDmitry 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);
770f746d493SDmitry Karpeev   if (pc->setupcalled && n != osm->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before PCSetUp().");
771b1a0a8a3SJed Brown   if (!pc->setupcalled) {
7726a4f0f73SDmitry Karpeev     ierr      = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr);
77393c1ec32SDmitry Karpeev     osm->N    = N;
774f746d493SDmitry Karpeev     osm->n    = n;
77506b43e7eSDmitry Karpeev     osm->nmax = N/size + ((N%size) ? 1 : 0);
7766a4f0f73SDmitry Karpeev     osm->ois  = 0;
7776a4f0f73SDmitry Karpeev     osm->iis  = 0;
778b1a0a8a3SJed Brown   }
779b1a0a8a3SJed Brown   PetscFunctionReturn(0);
78006b43e7eSDmitry Karpeev }
781b1a0a8a3SJed Brown EXTERN_C_END
782b1a0a8a3SJed Brown 
783b1a0a8a3SJed Brown EXTERN_C_BEGIN
784b1a0a8a3SJed Brown #undef __FUNCT__
785f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM"
7867087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
787b1a0a8a3SJed Brown {
788f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
789b1a0a8a3SJed Brown 
790b1a0a8a3SJed Brown   PetscFunctionBegin;
791*ce94432eSBarry Smith   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
792*ce94432eSBarry Smith   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
7932fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
794b1a0a8a3SJed Brown   PetscFunctionReturn(0);
795b1a0a8a3SJed Brown }
796b1a0a8a3SJed Brown EXTERN_C_END
797b1a0a8a3SJed Brown 
798b1a0a8a3SJed Brown EXTERN_C_BEGIN
799b1a0a8a3SJed Brown #undef __FUNCT__
800f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM"
8017087cfbeSBarry Smith PetscErrorCode  PCGASMSetType_GASM(PC pc,PCGASMType type)
802b1a0a8a3SJed Brown {
803f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
804b1a0a8a3SJed Brown 
805b1a0a8a3SJed Brown   PetscFunctionBegin;
806b1a0a8a3SJed Brown   osm->type     = type;
807b1a0a8a3SJed Brown   osm->type_set = PETSC_TRUE;
808b1a0a8a3SJed Brown   PetscFunctionReturn(0);
809b1a0a8a3SJed Brown }
810b1a0a8a3SJed Brown EXTERN_C_END
811b1a0a8a3SJed Brown 
812b1a0a8a3SJed Brown EXTERN_C_BEGIN
813b1a0a8a3SJed Brown #undef __FUNCT__
814f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM"
8157087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort)
816b1a0a8a3SJed Brown {
817f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
818b1a0a8a3SJed Brown 
819b1a0a8a3SJed Brown   PetscFunctionBegin;
820b1a0a8a3SJed Brown   osm->sort_indices = doSort;
821b1a0a8a3SJed Brown   PetscFunctionReturn(0);
822b1a0a8a3SJed Brown }
823b1a0a8a3SJed Brown EXTERN_C_END
824b1a0a8a3SJed Brown 
825b1a0a8a3SJed Brown EXTERN_C_BEGIN
826b1a0a8a3SJed Brown #undef __FUNCT__
827f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM"
82844fe18e5SDmitry Karpeev /*
82944fe18e5SDmitry Karpeev    FIX: This routine might need to be modified once multiple ranks per subdomain are allowed.
83044fe18e5SDmitry Karpeev         In particular, it would upset the global subdomain number calculation.
83144fe18e5SDmitry Karpeev */
8327087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
833b1a0a8a3SJed Brown {
834f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
835b1a0a8a3SJed Brown   PetscErrorCode ierr;
836b1a0a8a3SJed Brown 
837b1a0a8a3SJed Brown   PetscFunctionBegin;
838*ce94432eSBarry Smith   if (osm->n < 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here");
839b1a0a8a3SJed Brown 
8402fa5cd67SKarl Rupp   if (n) *n = osm->n;
841ab3e923fSDmitry Karpeev   if (first) {
842*ce94432eSBarry Smith     ierr    = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
843ab3e923fSDmitry Karpeev     *first -= osm->n;
844b1a0a8a3SJed Brown   }
845b1a0a8a3SJed Brown   if (ksp) {
846b1a0a8a3SJed Brown     /* Assume that local solves are now different; not necessarily
84706b43e7eSDmitry Karpeev        true, though!  This flag is used only for PCView_GASM() */
848b1a0a8a3SJed Brown     *ksp                        = osm->ksp;
8496a4f0f73SDmitry Karpeev     osm->same_subdomain_solvers = PETSC_FALSE;
850b1a0a8a3SJed Brown   }
851b1a0a8a3SJed Brown   PetscFunctionReturn(0);
852ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */
853b1a0a8a3SJed Brown EXTERN_C_END
854b1a0a8a3SJed Brown 
855b1a0a8a3SJed Brown 
856b1a0a8a3SJed Brown #undef __FUNCT__
85706b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains"
858b1a0a8a3SJed Brown /*@C
85906b43e7eSDmitry Karpeev     PCGASMSetSubdomains - Sets the subdomains for this processor
86006b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
861b1a0a8a3SJed Brown 
862b1a0a8a3SJed Brown     Collective on PC
863b1a0a8a3SJed Brown 
864b1a0a8a3SJed Brown     Input Parameters:
865b1a0a8a3SJed Brown +   pc  - the preconditioner context
86606b43e7eSDmitry Karpeev .   n   - the number of subdomains for this processor
8676a4f0f73SDmitry Karpeev .   iis - the index sets that define this processor's local inner subdomains
8680298fd71SBarry Smith          (or NULL for PETSc to determine subdomains)
8696a4f0f73SDmitry Karpeev -   ois- the index sets that define this processor's local outer subdomains
8700298fd71SBarry Smith          (or NULL to use the same as iis)
871b1a0a8a3SJed Brown 
872b1a0a8a3SJed Brown     Notes:
87306b43e7eSDmitry Karpeev     The IS indices use the parallel, global numbering of the vector entries.
8746a4f0f73SDmitry Karpeev     Inner subdomains are those where the correction is applied.
8756a4f0f73SDmitry Karpeev     Outer subdomains are those where the residual necessary to obtain the
8766a4f0f73SDmitry Karpeev     corrections is obtained (see PCGASMType for the use of inner/outer subdomains).
8776a4f0f73SDmitry Karpeev     Both inner and outer subdomains can extend over several processors.
8786a4f0f73SDmitry Karpeev     This processor's portion of a subdomain is known as a local subdomain.
8796a4f0f73SDmitry Karpeev 
8806a4f0f73SDmitry Karpeev     By default the GASM preconditioner uses 1 (local) subdomain per processor.
8816a4f0f73SDmitry Karpeev     Use PCGASMSetTotalSubdomains() to set the total number of subdomains across
8826a4f0f73SDmitry Karpeev     all processors that PCGASM will create automatically, and to specify whether
8836a4f0f73SDmitry Karpeev     they should be local or not.
8846a4f0f73SDmitry Karpeev 
885b1a0a8a3SJed Brown 
886b1a0a8a3SJed Brown     Level: advanced
887b1a0a8a3SJed Brown 
88806b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz
889b1a0a8a3SJed Brown 
8906a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
89106b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
892b1a0a8a3SJed Brown @*/
8936a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[])
894b1a0a8a3SJed Brown {
895b1a0a8a3SJed Brown   PetscErrorCode ierr;
896b1a0a8a3SJed Brown 
897b1a0a8a3SJed Brown   PetscFunctionBegin;
898b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8996a4f0f73SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr);
900b1a0a8a3SJed Brown   PetscFunctionReturn(0);
901b1a0a8a3SJed Brown }
902b1a0a8a3SJed Brown 
903b1a0a8a3SJed Brown #undef __FUNCT__
9046a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains"
905b1a0a8a3SJed Brown /*@C
9066a4f0f73SDmitry Karpeev     PCGASMSetTotalSubdomains - Sets the total number of subdomains to use in the generalized additive
9076a4f0f73SDmitry Karpeev     Schwarz preconditioner.  The number of subdomains is cumulative across all processors in pc's
9086a4f0f73SDmitry Karpeev     communicator. Either all or no processors in the PC communicator must call this routine with
9096a4f0f73SDmitry Karpeev     the same N.  The subdomains will be created automatically during PCSetUp().
910b1a0a8a3SJed Brown 
911b1a0a8a3SJed Brown     Collective on PC
912b1a0a8a3SJed Brown 
913b1a0a8a3SJed Brown     Input Parameters:
914b1a0a8a3SJed Brown +   pc           - the preconditioner context
9156a4f0f73SDmitry Karpeev .   N            - the total number of subdomains cumulative across all processors
9166a4f0f73SDmitry Karpeev -   create_local - whether the subdomains to be created are to be local
917b1a0a8a3SJed Brown 
918b1a0a8a3SJed Brown     Options Database Key:
9196a4f0f73SDmitry Karpeev     To set the total number of subdomains and let PCGASM autocreate them, rather than specify the index sets, use the following options:
9206a4f0f73SDmitry Karpeev +    -pc_gasm_total_subdomains <n>                  - sets the total number of subdomains to be autocreated by PCGASM
9216a4f0f73SDmitry Karpeev -    -pc_gasm_subdomains_create_local <true|false>  - whether autocreated subdomains should be local or not (default is true)
922b1a0a8a3SJed Brown 
92306b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.
924b1a0a8a3SJed Brown 
9256a4f0f73SDmitry Karpeev 
9266a4f0f73SDmitry Karpeev     Use PCGASMSetSubdomains() to set subdomains explicitly or to set different numbers
9276a4f0f73SDmitry Karpeev     of subdomains per processor.
928b1a0a8a3SJed Brown 
929b1a0a8a3SJed Brown     Level: advanced
930b1a0a8a3SJed Brown 
931f746d493SDmitry Karpeev .keywords: PC, GASM, set, total, global, subdomains, additive Schwarz
932b1a0a8a3SJed Brown 
93306b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
934f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
935b1a0a8a3SJed Brown @*/
9366a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetTotalSubdomains(PC pc,PetscInt N, PetscBool create_local)
937b1a0a8a3SJed Brown {
938b1a0a8a3SJed Brown   PetscErrorCode ierr;
939b1a0a8a3SJed Brown 
940b1a0a8a3SJed Brown   PetscFunctionBegin;
941b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9426a4f0f73SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetTotalSubdomains_C",(PC,PetscInt,PetscBool),(pc,N,create_local));CHKERRQ(ierr);
943b1a0a8a3SJed Brown   PetscFunctionReturn(0);
944b1a0a8a3SJed Brown }
945b1a0a8a3SJed Brown 
946b1a0a8a3SJed Brown #undef __FUNCT__
947f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap"
948b1a0a8a3SJed Brown /*@
949f746d493SDmitry Karpeev     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
950b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
951b1a0a8a3SJed Brown     PC communicator must call this routine.
952b1a0a8a3SJed Brown 
953b1a0a8a3SJed Brown     Logically Collective on PC
954b1a0a8a3SJed Brown 
955b1a0a8a3SJed Brown     Input Parameters:
956b1a0a8a3SJed Brown +   pc  - the preconditioner context
957b1a0a8a3SJed Brown -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
958b1a0a8a3SJed Brown 
959b1a0a8a3SJed Brown     Options Database Key:
96006b43e7eSDmitry Karpeev .   -pc_gasm_overlap <overlap> - Sets overlap
961b1a0a8a3SJed Brown 
962b1a0a8a3SJed Brown     Notes:
96306b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.  To use
9646a4f0f73SDmitry Karpeev     multiple subdomain per perocessor, see PCGASMSetTotalSubdomains() or
9656a4f0f73SDmitry Karpeev     PCGASMSetSubdomains() (and the option -pc_gasm_total_subdomains <n>).
966b1a0a8a3SJed Brown 
967b1a0a8a3SJed Brown     The overlap defaults to 1, so if one desires that no additional
968b1a0a8a3SJed Brown     overlap be computed beyond what may have been set with a call to
9696a4f0f73SDmitry Karpeev     PCGASMSetTotalSubdomains() or PCGASMSetSubdomains(), then ovl
970b1a0a8a3SJed Brown     must be set to be 0.  In particular, if one does not explicitly set
97106b43e7eSDmitry Karpeev     the subdomains in application code, then all overlap would be computed
972f746d493SDmitry Karpeev     internally by PETSc, and using an overlap of 0 would result in an GASM
973b1a0a8a3SJed Brown     variant that is equivalent to the block Jacobi preconditioner.
974b1a0a8a3SJed Brown 
975b1a0a8a3SJed Brown     Note that one can define initial index sets with any overlap via
97606b43e7eSDmitry Karpeev     PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows
97706b43e7eSDmitry Karpeev     PETSc to extend that overlap further, if desired.
978b1a0a8a3SJed Brown 
979b1a0a8a3SJed Brown     Level: intermediate
980b1a0a8a3SJed Brown 
981f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap
982b1a0a8a3SJed Brown 
9836a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
98406b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
985b1a0a8a3SJed Brown @*/
9867087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap(PC pc,PetscInt ovl)
987b1a0a8a3SJed Brown {
988b1a0a8a3SJed Brown   PetscErrorCode ierr;
989b1a0a8a3SJed Brown 
990b1a0a8a3SJed Brown   PetscFunctionBegin;
991b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
992b1a0a8a3SJed Brown   PetscValidLogicalCollectiveInt(pc,ovl,2);
993f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
994b1a0a8a3SJed Brown   PetscFunctionReturn(0);
995b1a0a8a3SJed Brown }
996b1a0a8a3SJed Brown 
997b1a0a8a3SJed Brown #undef __FUNCT__
998f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType"
999b1a0a8a3SJed Brown /*@
1000f746d493SDmitry Karpeev     PCGASMSetType - Sets the type of restriction and interpolation used
1001b1a0a8a3SJed Brown     for local problems in the additive Schwarz method.
1002b1a0a8a3SJed Brown 
1003b1a0a8a3SJed Brown     Logically Collective on PC
1004b1a0a8a3SJed Brown 
1005b1a0a8a3SJed Brown     Input Parameters:
1006b1a0a8a3SJed Brown +   pc  - the preconditioner context
1007f746d493SDmitry Karpeev -   type - variant of GASM, one of
1008b1a0a8a3SJed Brown .vb
1009f746d493SDmitry Karpeev       PC_GASM_BASIC       - full interpolation and restriction
1010f746d493SDmitry Karpeev       PC_GASM_RESTRICT    - full restriction, local processor interpolation
1011f746d493SDmitry Karpeev       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
1012f746d493SDmitry Karpeev       PC_GASM_NONE        - local processor restriction and interpolation
1013b1a0a8a3SJed Brown .ve
1014b1a0a8a3SJed Brown 
1015b1a0a8a3SJed Brown     Options Database Key:
1016f746d493SDmitry Karpeev .   -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1017b1a0a8a3SJed Brown 
1018b1a0a8a3SJed Brown     Level: intermediate
1019b1a0a8a3SJed Brown 
1020f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1021b1a0a8a3SJed Brown 
10226a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1023f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1024b1a0a8a3SJed Brown @*/
10257087cfbeSBarry Smith PetscErrorCode  PCGASMSetType(PC pc,PCGASMType type)
1026b1a0a8a3SJed Brown {
1027b1a0a8a3SJed Brown   PetscErrorCode ierr;
1028b1a0a8a3SJed Brown 
1029b1a0a8a3SJed Brown   PetscFunctionBegin;
1030b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1031b1a0a8a3SJed Brown   PetscValidLogicalCollectiveEnum(pc,type,2);
1032f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr);
1033b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1034b1a0a8a3SJed Brown }
1035b1a0a8a3SJed Brown 
1036b1a0a8a3SJed Brown #undef __FUNCT__
1037f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices"
1038b1a0a8a3SJed Brown /*@
1039f746d493SDmitry Karpeev     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
1040b1a0a8a3SJed Brown 
1041b1a0a8a3SJed Brown     Logically Collective on PC
1042b1a0a8a3SJed Brown 
1043b1a0a8a3SJed Brown     Input Parameters:
1044b1a0a8a3SJed Brown +   pc  - the preconditioner context
1045b1a0a8a3SJed Brown -   doSort - sort the subdomain indices
1046b1a0a8a3SJed Brown 
1047b1a0a8a3SJed Brown     Level: intermediate
1048b1a0a8a3SJed Brown 
1049f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1050b1a0a8a3SJed Brown 
10516a4f0f73SDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(),
1052f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1053b1a0a8a3SJed Brown @*/
10547087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices(PC pc,PetscBool doSort)
1055b1a0a8a3SJed Brown {
1056b1a0a8a3SJed Brown   PetscErrorCode ierr;
1057b1a0a8a3SJed Brown 
1058b1a0a8a3SJed Brown   PetscFunctionBegin;
1059b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1060b1a0a8a3SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
1061f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
1062b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1063b1a0a8a3SJed Brown }
1064b1a0a8a3SJed Brown 
1065b1a0a8a3SJed Brown #undef __FUNCT__
1066f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP"
1067b1a0a8a3SJed Brown /*@C
1068f746d493SDmitry Karpeev    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
1069b1a0a8a3SJed Brown    this processor.
1070b1a0a8a3SJed Brown 
1071b1a0a8a3SJed Brown    Collective on PC iff first_local is requested
1072b1a0a8a3SJed Brown 
1073b1a0a8a3SJed Brown    Input Parameter:
1074b1a0a8a3SJed Brown .  pc - the preconditioner context
1075b1a0a8a3SJed Brown 
1076b1a0a8a3SJed Brown    Output Parameters:
10770298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
10780298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
10790298fd71SBarry Smith                  all processors must request or all must pass NULL
1080b1a0a8a3SJed Brown -  ksp - the array of KSP contexts
1081b1a0a8a3SJed Brown 
1082b1a0a8a3SJed Brown    Note:
1083f746d493SDmitry Karpeev    After PCGASMGetSubKSP() the array of KSPes is not to be freed
1084b1a0a8a3SJed Brown 
1085b1a0a8a3SJed Brown    Currently for some matrix implementations only 1 block per processor
1086b1a0a8a3SJed Brown    is supported.
1087b1a0a8a3SJed Brown 
1088f746d493SDmitry Karpeev    You must call KSPSetUp() before calling PCGASMGetSubKSP().
1089b1a0a8a3SJed Brown 
1090b1a0a8a3SJed Brown    Level: advanced
1091b1a0a8a3SJed Brown 
1092f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context
1093b1a0a8a3SJed Brown 
10946a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap(),
1095f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(),
1096b1a0a8a3SJed Brown @*/
10977087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1098b1a0a8a3SJed Brown {
1099b1a0a8a3SJed Brown   PetscErrorCode ierr;
1100b1a0a8a3SJed Brown 
1101b1a0a8a3SJed Brown   PetscFunctionBegin;
1102b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1103f746d493SDmitry Karpeev   ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
1104b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1105b1a0a8a3SJed Brown }
1106b1a0a8a3SJed Brown 
1107b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/
1108b1a0a8a3SJed Brown /*MC
1109f746d493SDmitry Karpeev    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1110b1a0a8a3SJed Brown            its own KSP object.
1111b1a0a8a3SJed Brown 
1112b1a0a8a3SJed Brown    Options Database Keys:
111306b43e7eSDmitry Karpeev +  -pc_gasm_total_block_count <n> - Sets total number of local subdomains (known as blocks) to be distributed among processors
111406b43e7eSDmitry Karpeev .  -pc_gasm_view_subdomains       - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view
111506b43e7eSDmitry Karpeev .  -pc_gasm_print_subdomains      - activates the printing of subdomain indices in PCSetUp()
111606b43e7eSDmitry Karpeev .  -pc_gasm_overlap <ovl>         - Sets overlap by which to (automatically) extend local subdomains
1117f746d493SDmitry Karpeev -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1118b1a0a8a3SJed Brown 
1119b1a0a8a3SJed Brown      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
1120f746d493SDmitry Karpeev       will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
1121f746d493SDmitry Karpeev       -pc_gasm_type basic to use the standard GASM.
1122b1a0a8a3SJed Brown 
1123b1a0a8a3SJed Brown    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
1124b1a0a8a3SJed Brown      than one processor. Defaults to one block per processor.
1125b1a0a8a3SJed Brown 
1126b1a0a8a3SJed Brown      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1127b1a0a8a3SJed Brown         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
1128b1a0a8a3SJed Brown 
1129f746d493SDmitry Karpeev      To set the options on the solvers separate for each block call PCGASMGetSubKSP()
1130b1a0a8a3SJed Brown          and set the options directly on the resulting KSP object (you can access its PC
1131b1a0a8a3SJed Brown          with KSPGetPC())
1132b1a0a8a3SJed Brown 
1133b1a0a8a3SJed Brown 
1134b1a0a8a3SJed Brown    Level: beginner
1135b1a0a8a3SJed Brown 
1136b1a0a8a3SJed Brown    Concepts: additive Schwarz method
1137b1a0a8a3SJed Brown 
1138b1a0a8a3SJed Brown     References:
1139b1a0a8a3SJed Brown     An additive variant of the Schwarz alternating method for the case of many subregions
1140b1a0a8a3SJed Brown     M Dryja, OB Widlund - Courant Institute, New York University Technical report
1141b1a0a8a3SJed Brown 
1142b1a0a8a3SJed Brown     Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1143b1a0a8a3SJed Brown     Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X.
1144b1a0a8a3SJed Brown 
1145b1a0a8a3SJed Brown .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
114606b43e7eSDmitry Karpeev            PCBJACOBI, PCGASMSetUseTrueLocal(), PCGASMGetSubKSP(), PCGASMSetSubdomains(),
11476a4f0f73SDmitry Karpeev            PCGASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType()
1148b1a0a8a3SJed Brown 
1149b1a0a8a3SJed Brown M*/
1150b1a0a8a3SJed Brown 
1151b1a0a8a3SJed Brown EXTERN_C_BEGIN
1152b1a0a8a3SJed Brown #undef __FUNCT__
1153f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM"
11547087cfbeSBarry Smith PetscErrorCode  PCCreate_GASM(PC pc)
1155b1a0a8a3SJed Brown {
1156b1a0a8a3SJed Brown   PetscErrorCode ierr;
1157f746d493SDmitry Karpeev   PC_GASM        *osm;
1158b1a0a8a3SJed Brown 
1159b1a0a8a3SJed Brown   PetscFunctionBegin;
1160f746d493SDmitry Karpeev   ierr = PetscNewLog(pc,PC_GASM,&osm);CHKERRQ(ierr);
11612fa5cd67SKarl Rupp 
1162ab3e923fSDmitry Karpeev   osm->N                      = PETSC_DECIDE;
116306b43e7eSDmitry Karpeev   osm->n                      = PETSC_DECIDE;
1164ab3e923fSDmitry Karpeev   osm->nmax                   = 0;
1165b1a0a8a3SJed Brown   osm->overlap                = 1;
1166b1a0a8a3SJed Brown   osm->ksp                    = 0;
11676a4f0f73SDmitry Karpeev   osm->gorestriction          = 0;
11686a4f0f73SDmitry Karpeev   osm->girestriction          = 0;
1169ab3e923fSDmitry Karpeev   osm->gx                     = 0;
1170ab3e923fSDmitry Karpeev   osm->gy                     = 0;
1171b1a0a8a3SJed Brown   osm->x                      = 0;
1172b1a0a8a3SJed Brown   osm->y                      = 0;
11736a4f0f73SDmitry Karpeev   osm->ois                    = 0;
11746a4f0f73SDmitry Karpeev   osm->iis                    = 0;
1175b1a0a8a3SJed Brown   osm->pmat                   = 0;
1176f746d493SDmitry Karpeev   osm->type                   = PC_GASM_RESTRICT;
11776a4f0f73SDmitry Karpeev   osm->same_subdomain_solvers = PETSC_TRUE;
1178b1a0a8a3SJed Brown   osm->sort_indices           = PETSC_TRUE;
117965db9045SDmitry Karpeev   osm->use_dm_decomposition   = PETSC_TRUE;
1180b1a0a8a3SJed Brown 
1181b1a0a8a3SJed Brown   pc->data                 = (void*)osm;
1182f746d493SDmitry Karpeev   pc->ops->apply           = PCApply_GASM;
1183f746d493SDmitry Karpeev   pc->ops->applytranspose  = PCApplyTranspose_GASM;
1184f746d493SDmitry Karpeev   pc->ops->setup           = PCSetUp_GASM;
1185a06653b4SBarry Smith   pc->ops->reset           = PCReset_GASM;
1186f746d493SDmitry Karpeev   pc->ops->destroy         = PCDestroy_GASM;
1187f746d493SDmitry Karpeev   pc->ops->setfromoptions  = PCSetFromOptions_GASM;
1188f746d493SDmitry Karpeev   pc->ops->setuponblocks   = PCSetUpOnBlocks_GASM;
1189f746d493SDmitry Karpeev   pc->ops->view            = PCView_GASM;
1190b1a0a8a3SJed Brown   pc->ops->applyrichardson = 0;
1191b1a0a8a3SJed Brown 
119206b43e7eSDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSubdomains_C","PCGASMSetSubdomains_GASM",
119306b43e7eSDmitry Karpeev                                            PCGASMSetSubdomains_GASM);CHKERRQ(ierr);
11946a4f0f73SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetTotalSubdomains_C","PCGASMSetTotalSubdomains_GASM",
11956a4f0f73SDmitry Karpeev                                            PCGASMSetTotalSubdomains_GASM);CHKERRQ(ierr);
1196f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetOverlap_C","PCGASMSetOverlap_GASM",
1197f746d493SDmitry Karpeev                                            PCGASMSetOverlap_GASM);CHKERRQ(ierr);
1198f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetType_C","PCGASMSetType_GASM",
1199f746d493SDmitry Karpeev                                            PCGASMSetType_GASM);CHKERRQ(ierr);
1200f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSortIndices_C","PCGASMSetSortIndices_GASM",
1201f746d493SDmitry Karpeev                                            PCGASMSetSortIndices_GASM);CHKERRQ(ierr);
1202f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMGetSubKSP_C","PCGASMGetSubKSP_GASM",
1203f746d493SDmitry Karpeev                                            PCGASMGetSubKSP_GASM);CHKERRQ(ierr);
1204b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1205b1a0a8a3SJed Brown }
1206b1a0a8a3SJed Brown EXTERN_C_END
1207b1a0a8a3SJed Brown 
1208b1a0a8a3SJed Brown 
1209b1a0a8a3SJed Brown #undef __FUNCT__
121006b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains"
1211b1a0a8a3SJed Brown /*@C
121206b43e7eSDmitry Karpeev    PCGASMCreateLocalSubdomains - Creates n local index sets for the overlapping
121306b43e7eSDmitry Karpeev    Schwarz preconditioner for a any problem based on its matrix.
1214b1a0a8a3SJed Brown 
1215b1a0a8a3SJed Brown    Collective
1216b1a0a8a3SJed Brown 
1217b1a0a8a3SJed Brown    Input Parameters:
1218b1a0a8a3SJed Brown +  A       - The global matrix operator
12196a4f0f73SDmitry Karpeev .  overlap - amount of overlap in outer subdomains
122006b43e7eSDmitry Karpeev -  n       - the number of local subdomains
1221b1a0a8a3SJed Brown 
1222b1a0a8a3SJed Brown    Output Parameters:
12236a4f0f73SDmitry Karpeev +  iis - the array of index sets defining the local inner subdomains (on which the correction is applied)
12246a4f0f73SDmitry Karpeev -  ois - the array of index sets defining the local outer subdomains (on which the residual is computed)
1225b1a0a8a3SJed Brown 
1226b1a0a8a3SJed Brown    Level: advanced
1227b1a0a8a3SJed Brown 
12286a4f0f73SDmitry Karpeev    Note: this generates n nonoverlapping local inner subdomains on PETSC_COMM_SELF;
12296a4f0f73SDmitry Karpeev          PCGASM will generate the overlap from these if you use them in PCGASMSetSubdomains() and set a
12306a4f0f73SDmitry Karpeev          nonzero overlap with PCGASMSetOverlap()
1231b1a0a8a3SJed Brown 
1232b1a0a8a3SJed Brown     In the Fortran version you must provide the array outis[] already allocated of length n.
1233b1a0a8a3SJed Brown 
1234f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1235b1a0a8a3SJed Brown 
123606b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains()
1237b1a0a8a3SJed Brown @*/
12386a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt overlap, PetscInt n, IS *iis[], IS *ois[])
1239b1a0a8a3SJed Brown {
1240b1a0a8a3SJed Brown   MatPartitioning mpart;
1241b1a0a8a3SJed Brown   const char      *prefix;
124211bd1e4dSLisandro Dalcin   PetscErrorCode  (*f)(Mat,MatReuse,Mat*);
1243b1a0a8a3SJed Brown   PetscMPIInt     size;
1244b1a0a8a3SJed Brown   PetscInt        i,j,rstart,rend,bs;
124511bd1e4dSLisandro Dalcin   PetscBool       isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
12460298fd71SBarry Smith   Mat             Ad     = NULL, adj;
1247b1a0a8a3SJed Brown   IS              ispart,isnumb,*is;
1248b1a0a8a3SJed Brown   PetscErrorCode  ierr;
1249b1a0a8a3SJed Brown 
1250b1a0a8a3SJed Brown   PetscFunctionBegin;
1251b1a0a8a3SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
12526a4f0f73SDmitry Karpeev   PetscValidPointer(iis,4);
1253b1a0a8a3SJed Brown   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
1254b1a0a8a3SJed Brown 
1255b1a0a8a3SJed Brown   /* Get prefix, row distribution, and block size */
1256b1a0a8a3SJed Brown   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1257b1a0a8a3SJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1258b1a0a8a3SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1259b1a0a8a3SJed 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);
1260b1a0a8a3SJed Brown 
1261b1a0a8a3SJed Brown   /* Get diagonal block from matrix if possible */
1262*ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1263b1a0a8a3SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr);
1264b1a0a8a3SJed Brown   if (f) {
126511bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
1266b1a0a8a3SJed Brown   } else if (size == 1) {
126711bd1e4dSLisandro Dalcin     Ad = A;
1268b1a0a8a3SJed Brown   }
1269b1a0a8a3SJed Brown   if (Ad) {
1270251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1271251f4c67SDmitry Karpeev     if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
1272b1a0a8a3SJed Brown   }
1273b1a0a8a3SJed Brown   if (Ad && n > 1) {
1274b1a0a8a3SJed Brown     PetscBool  match,done;
1275b1a0a8a3SJed Brown     /* Try to setup a good matrix partitioning if available */
1276b1a0a8a3SJed Brown     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
1277b1a0a8a3SJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
1278b1a0a8a3SJed Brown     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1279251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
1280b1a0a8a3SJed Brown     if (!match) {
1281251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
1282b1a0a8a3SJed Brown     }
1283b1a0a8a3SJed Brown     if (!match) { /* assume a "good" partitioner is available */
12841a83f524SJed Brown       PetscInt       na;
12851a83f524SJed Brown       const PetscInt *ia,*ja;
1286b1a0a8a3SJed Brown       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1287b1a0a8a3SJed Brown       if (done) {
1288b1a0a8a3SJed Brown         /* Build adjacency matrix by hand. Unfortunately a call to
1289b1a0a8a3SJed Brown            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1290b1a0a8a3SJed Brown            remove the block-aij structure and we cannot expect
1291b1a0a8a3SJed Brown            MatPartitioning to split vertices as we need */
12921a83f524SJed Brown         PetscInt       i,j,len,nnz,cnt,*iia=0,*jja=0;
12931a83f524SJed Brown         const PetscInt *row;
1294b1a0a8a3SJed Brown         nnz = 0;
1295b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* count number of nonzeros */
1296b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1297b1a0a8a3SJed Brown           row = ja + ia[i];
1298b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1299b1a0a8a3SJed Brown             if (row[j] == i) { /* don't count diagonal */
1300b1a0a8a3SJed Brown               len--; break;
1301b1a0a8a3SJed Brown             }
1302b1a0a8a3SJed Brown           }
1303b1a0a8a3SJed Brown           nnz += len;
1304b1a0a8a3SJed Brown         }
1305b1a0a8a3SJed Brown         ierr   = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr);
1306b1a0a8a3SJed Brown         ierr   = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr);
1307b1a0a8a3SJed Brown         nnz    = 0;
1308b1a0a8a3SJed Brown         iia[0] = 0;
1309b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* fill adjacency */
1310b1a0a8a3SJed Brown           cnt = 0;
1311b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1312b1a0a8a3SJed Brown           row = ja + ia[i];
1313b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
13142fa5cd67SKarl Rupp             if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */
1315b1a0a8a3SJed Brown           }
1316b1a0a8a3SJed Brown           nnz += cnt;
1317b1a0a8a3SJed Brown           iia[i+1] = nnz;
1318b1a0a8a3SJed Brown         }
1319b1a0a8a3SJed Brown         /* Partitioning of the adjacency matrix */
13200298fd71SBarry Smith         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr);
1321b1a0a8a3SJed Brown         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
1322b1a0a8a3SJed Brown         ierr      = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
1323b1a0a8a3SJed Brown         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
1324b1a0a8a3SJed Brown         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
13256bf464f9SBarry Smith         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
1326b1a0a8a3SJed Brown         foundpart = PETSC_TRUE;
1327b1a0a8a3SJed Brown       }
1328b1a0a8a3SJed Brown       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1329b1a0a8a3SJed Brown     }
13306bf464f9SBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1331b1a0a8a3SJed Brown   }
1332b1a0a8a3SJed Brown   ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr);
1333b1a0a8a3SJed Brown   if (!foundpart) {
1334b1a0a8a3SJed Brown 
1335b1a0a8a3SJed Brown     /* Partitioning by contiguous chunks of rows */
1336b1a0a8a3SJed Brown 
1337b1a0a8a3SJed Brown     PetscInt mbs   = (rend-rstart)/bs;
1338b1a0a8a3SJed Brown     PetscInt start = rstart;
1339b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1340b1a0a8a3SJed Brown       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
1341b1a0a8a3SJed Brown       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
1342b1a0a8a3SJed Brown       start += count;
1343b1a0a8a3SJed Brown     }
1344b1a0a8a3SJed Brown 
1345b1a0a8a3SJed Brown   } else {
1346b1a0a8a3SJed Brown 
1347b1a0a8a3SJed Brown     /* Partitioning by adjacency of diagonal block  */
1348b1a0a8a3SJed Brown 
1349b1a0a8a3SJed Brown     const PetscInt *numbering;
1350b1a0a8a3SJed Brown     PetscInt       *count,nidx,*indices,*newidx,start=0;
1351b1a0a8a3SJed Brown     /* Get node count in each partition */
1352b1a0a8a3SJed Brown     ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr);
1353b1a0a8a3SJed Brown     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
1354b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1355b1a0a8a3SJed Brown       for (i=0; i<n; i++) count[i] *= bs;
1356b1a0a8a3SJed Brown     }
1357b1a0a8a3SJed Brown     /* Build indices from node numbering */
1358b1a0a8a3SJed Brown     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1359b1a0a8a3SJed Brown     ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr);
1360b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1361b1a0a8a3SJed Brown     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
1362b1a0a8a3SJed Brown     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
1363b1a0a8a3SJed Brown     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
1364b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1365b1a0a8a3SJed Brown       ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr);
13662fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
13672fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
13682fa5cd67SKarl Rupp       }
1369b1a0a8a3SJed Brown       ierr    = PetscFree(indices);CHKERRQ(ierr);
1370b1a0a8a3SJed Brown       nidx   *= bs;
1371b1a0a8a3SJed Brown       indices = newidx;
1372b1a0a8a3SJed Brown     }
1373b1a0a8a3SJed Brown     /* Shift to get global indices */
1374b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] += rstart;
1375b1a0a8a3SJed Brown 
1376b1a0a8a3SJed Brown     /* Build the index sets for each block */
1377b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1378b1a0a8a3SJed Brown       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
1379b1a0a8a3SJed Brown       ierr   = ISSort(is[i]);CHKERRQ(ierr);
1380b1a0a8a3SJed Brown       start += count[i];
1381b1a0a8a3SJed Brown     }
1382b1a0a8a3SJed Brown 
1383b1a0a8a3SJed Brown     ierr = PetscFree(count);
1384b1a0a8a3SJed Brown     ierr = PetscFree(indices);
1385fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1386fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
1387b1a0a8a3SJed Brown   }
13886a4f0f73SDmitry Karpeev   *iis = is;
13896a4f0f73SDmitry Karpeev   if (!ois) PetscFunctionReturn(0);
13906a4f0f73SDmitry Karpeev   /*
13916a4f0f73SDmitry Karpeev    Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap
13926a4f0f73SDmitry Karpeev    has been requested, copy the inner subdomains over so they can be modified.
13936a4f0f73SDmitry Karpeev    */
13946a4f0f73SDmitry Karpeev   ierr = PetscMalloc(n*sizeof(IS),ois);CHKERRQ(ierr);
13956a4f0f73SDmitry Karpeev   for (i=0; i<n; ++i) {
13966a4f0f73SDmitry Karpeev     if (overlap > 0) { /* With positive overlap, (*iis)[i] will be modified */
13976a4f0f73SDmitry Karpeev       ierr = ISDuplicate((*iis)[i],(*ois)+i);CHKERRQ(ierr);
13986a4f0f73SDmitry Karpeev       ierr = ISCopy((*iis)[i],(*ois)[i]);CHKERRQ(ierr);
13996a4f0f73SDmitry Karpeev     } else {
14006a4f0f73SDmitry Karpeev       ierr      = PetscObjectReference((PetscObject)(*iis)[i]);CHKERRQ(ierr);
14016a4f0f73SDmitry Karpeev       (*ois)[i] = (*iis)[i];
14026a4f0f73SDmitry Karpeev     }
14036a4f0f73SDmitry Karpeev   }
14046a4f0f73SDmitry Karpeev   if (overlap > 0) {
14056a4f0f73SDmitry Karpeev     /* Extend the "overlapping" regions by a number of steps */
14066a4f0f73SDmitry Karpeev     ierr = MatIncreaseOverlap(A,n,*ois,overlap);CHKERRQ(ierr);
14076a4f0f73SDmitry Karpeev   }
1408b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1409b1a0a8a3SJed Brown }
1410b1a0a8a3SJed Brown 
1411b1a0a8a3SJed Brown #undef __FUNCT__
1412f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains"
1413b1a0a8a3SJed Brown /*@C
1414f746d493SDmitry Karpeev    PCGASMDestroySubdomains - Destroys the index sets created with
141506b43e7eSDmitry Karpeev    PCGASMCreateLocalSubdomains() or PCGASMCreateSubdomains2D. Should be
141606b43e7eSDmitry Karpeev    called after setting subdomains with PCGASMSetSubdomains().
1417b1a0a8a3SJed Brown 
1418b1a0a8a3SJed Brown    Collective
1419b1a0a8a3SJed Brown 
1420b1a0a8a3SJed Brown    Input Parameters:
1421b1a0a8a3SJed Brown +  n   - the number of index sets
14226a4f0f73SDmitry Karpeev .  iis - the array of inner subdomains,
14230298fd71SBarry Smith -  ois - the array of outer subdomains, can be NULL
1424b1a0a8a3SJed Brown 
14256a4f0f73SDmitry Karpeev    Level: intermediate
14266a4f0f73SDmitry Karpeev 
14276a4f0f73SDmitry Karpeev    Notes: this is merely a convenience subroutine that walks each list,
14286a4f0f73SDmitry Karpeev    destroys each IS on the list, and then frees the list.
1429b1a0a8a3SJed Brown 
1430f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1431b1a0a8a3SJed Brown 
143206b43e7eSDmitry Karpeev .seealso: PCGASMCreateLocalSubdomains(), PCGASMSetSubdomains()
1433b1a0a8a3SJed Brown @*/
14346a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMDestroySubdomains(PetscInt n, IS iis[], IS ois[])
1435b1a0a8a3SJed Brown {
1436b1a0a8a3SJed Brown   PetscInt       i;
1437b1a0a8a3SJed Brown   PetscErrorCode ierr;
14385fd66863SKarl Rupp 
1439b1a0a8a3SJed Brown   PetscFunctionBegin;
1440a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1441a35b7b57SDmitry Karpeev   if (iis) {
14426a4f0f73SDmitry Karpeev     PetscValidPointer(iis,2);
1443a35b7b57SDmitry Karpeev     for (i=0; i<n; i++) {
1444a35b7b57SDmitry Karpeev       ierr = ISDestroy(&iis[i]);CHKERRQ(ierr);
1445a35b7b57SDmitry Karpeev     }
14466a4f0f73SDmitry Karpeev     ierr = PetscFree(iis);CHKERRQ(ierr);
1447a35b7b57SDmitry Karpeev   }
14486a4f0f73SDmitry Karpeev   if (ois) {
1449a35b7b57SDmitry Karpeev     for (i=0; i<n; i++) {
1450a35b7b57SDmitry Karpeev       ierr = ISDestroy(&ois[i]);CHKERRQ(ierr);
1451a35b7b57SDmitry Karpeev     }
14526a4f0f73SDmitry Karpeev     ierr = PetscFree(ois);CHKERRQ(ierr);
1453b1a0a8a3SJed Brown   }
1454b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1455b1a0a8a3SJed Brown }
1456b1a0a8a3SJed Brown 
1457af538404SDmitry Karpeev 
1458af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \
1459af538404SDmitry Karpeev   {                                                                                                       \
1460af538404SDmitry Karpeev     PetscInt first_row = first/M, last_row = last/M+1;                                                     \
1461af538404SDmitry Karpeev     /*                                                                                                    \
1462af538404SDmitry Karpeev      Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners       \
1463af538404SDmitry Karpeev      of the bounding box of the intersection of the subdomain with the local ownership range (local       \
1464af538404SDmitry Karpeev      subdomain).                                                                                          \
1465af538404SDmitry Karpeev      Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows       \
1466af538404SDmitry Karpeev      of the intersection.                                                                                 \
1467af538404SDmitry Karpeev     */                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                             \
1468af538404SDmitry Karpeev     /* ylow_loc is the grid row containing the first element of the local sumbdomain */                   \
1469eec7e3faSDmitry Karpeev     *ylow_loc = PetscMax(first_row,ylow);                                                                    \
1470af538404SDmitry Karpeev     /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1471eec7e3faSDmitry Karpeev     *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft;                                                                            \
1472af538404SDmitry Karpeev     /* yhigh_loc is the grid row above the last local subdomain element */                                                                    \
1473eec7e3faSDmitry Karpeev     *yhigh_loc = PetscMin(last_row,yhigh);                                                                                                     \
1474af538404SDmitry Karpeev     /* xright is the offset of the end of the  local subdomain within its grid row (might actually be outside the local subdomain) */         \
1475eec7e3faSDmitry Karpeev     *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright;                                                                          \
1476af538404SDmitry Karpeev     /* Now compute the size of the local subdomain n. */ \
1477c3518bceSDmitry Karpeev     *n = 0;                                               \
1478eec7e3faSDmitry Karpeev     if (*ylow_loc < *yhigh_loc) {                           \
1479af538404SDmitry Karpeev       PetscInt width = xright-xleft;                     \
1480eec7e3faSDmitry Karpeev       *n += width*(*yhigh_loc-*ylow_loc-1);                 \
1481eec7e3faSDmitry Karpeev       *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \
1482eec7e3faSDmitry Karpeev       *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \
1483af538404SDmitry Karpeev     } \
1484af538404SDmitry Karpeev   }
1485af538404SDmitry Karpeev 
1486c3518bceSDmitry Karpeev 
1487eec7e3faSDmitry Karpeev 
1488eec7e3faSDmitry Karpeev #undef __FUNCT__
1489f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D"
1490b1a0a8a3SJed Brown /*@
1491f746d493SDmitry Karpeev    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1492b1a0a8a3SJed Brown    preconditioner for a two-dimensional problem on a regular grid.
1493b1a0a8a3SJed Brown 
1494af538404SDmitry Karpeev    Collective
1495b1a0a8a3SJed Brown 
1496b1a0a8a3SJed Brown    Input Parameters:
149706b43e7eSDmitry Karpeev +  M, N               - the global number of grid points in the x and y directions
1498af538404SDmitry Karpeev .  Mdomains, Ndomains - the global number of subdomains in the x and y directions
1499b1a0a8a3SJed Brown .  dof                - degrees of freedom per node
1500b1a0a8a3SJed Brown -  overlap            - overlap in mesh lines
1501b1a0a8a3SJed Brown 
1502b1a0a8a3SJed Brown    Output Parameters:
1503af538404SDmitry Karpeev +  Nsub - the number of local subdomains created
15046a4f0f73SDmitry Karpeev .  iis  - array of index sets defining inner (nonoverlapping) subdomains
15056a4f0f73SDmitry Karpeev -  ois  - array of index sets defining outer (overlapping, if overlap > 0) subdomains
1506b1a0a8a3SJed Brown 
1507b1a0a8a3SJed Brown 
1508b1a0a8a3SJed Brown    Level: advanced
1509b1a0a8a3SJed Brown 
1510f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid
1511b1a0a8a3SJed Brown 
15126a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1513f746d493SDmitry Karpeev           PCGASMSetOverlap()
1514b1a0a8a3SJed Brown @*/
15156a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMCreateSubdomains2D(PC pc, PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap, PetscInt *nsub,IS **iis,IS **ois)
1516b1a0a8a3SJed Brown {
1517b1a0a8a3SJed Brown   PetscErrorCode ierr;
15182bbb417fSDmitry Karpeev   PetscMPIInt    size, rank;
15192bbb417fSDmitry Karpeev   PetscInt       i, j;
15202bbb417fSDmitry Karpeev   PetscInt       maxheight, maxwidth;
15212bbb417fSDmitry Karpeev   PetscInt       xstart, xleft, xright, xleft_loc, xright_loc;
15222bbb417fSDmitry Karpeev   PetscInt       ystart, ylow,  yhigh,  ylow_loc,  yhigh_loc;
1523eec7e3faSDmitry Karpeev   PetscInt       x[2][2], y[2][2], n[2];
15242bbb417fSDmitry Karpeev   PetscInt       first, last;
15252bbb417fSDmitry Karpeev   PetscInt       nidx, *idx;
15262bbb417fSDmitry Karpeev   PetscInt       ii,jj,s,q,d;
1527af538404SDmitry Karpeev   PetscInt       k,kk;
15282bbb417fSDmitry Karpeev   PetscMPIInt    color;
15292bbb417fSDmitry Karpeev   MPI_Comm       comm, subcomm;
15306a4f0f73SDmitry Karpeev   IS             **xis = 0, **is = ois, **is_local = iis;
1531d34fcf5fSBarry Smith 
1532b1a0a8a3SJed Brown   PetscFunctionBegin;
15332bbb417fSDmitry Karpeev   ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr);
15342bbb417fSDmitry Karpeev   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
15352bbb417fSDmitry Karpeev   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
15362bbb417fSDmitry Karpeev   ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr);
1537d34fcf5fSBarry 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) "
15382bbb417fSDmitry Karpeev                                       "does not respect the number of degrees of freedom per grid point %D", first, last, dof);
1539d34fcf5fSBarry Smith 
1540af538404SDmitry Karpeev   /* Determine the number of domains with nonzero intersections with the local ownership range. */
15412bbb417fSDmitry Karpeev   s      = 0;
1542af538404SDmitry Karpeev   ystart = 0;
1543af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1544af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1545af538404SDmitry 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);
1546eec7e3faSDmitry Karpeev     /* Vertical domain limits with an overlap. */
1547eec7e3faSDmitry Karpeev     ylow   = PetscMax(ystart - overlap,0);
1548eec7e3faSDmitry Karpeev     yhigh  = PetscMin(ystart + maxheight + overlap,N);
1549b1a0a8a3SJed Brown     xstart = 0;
1550af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1551af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1552af538404SDmitry 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);
1553eec7e3faSDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1554eec7e3faSDmitry Karpeev       xleft  = PetscMax(xstart - overlap,0);
1555eec7e3faSDmitry Karpeev       xright = PetscMin(xstart + maxwidth + overlap,M);
1556af538404SDmitry Karpeev       /*
1557af538404SDmitry Karpeev          Determine whether this subdomain intersects this processor's ownership range of pc->pmat.
1558af538404SDmitry Karpeev       */
1559c3518bceSDmitry Karpeev       PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
15602fa5cd67SKarl Rupp       if (nidx) ++s;
1561af538404SDmitry Karpeev       xstart += maxwidth;
1562af538404SDmitry Karpeev     } /* for (i = 0; i < Mdomains; ++i) */
1563af538404SDmitry Karpeev     ystart += maxheight;
1564af538404SDmitry Karpeev   } /* for (j = 0; j < Ndomains; ++j) */
15652fa5cd67SKarl Rupp 
1566af538404SDmitry Karpeev   /* Now we can allocate the necessary number of ISs. */
1567af538404SDmitry Karpeev   *nsub  = s;
1568af538404SDmitry Karpeev   ierr   = PetscMalloc((*nsub)*sizeof(IS*),is);CHKERRQ(ierr);
1569af538404SDmitry Karpeev   ierr   = PetscMalloc((*nsub)*sizeof(IS*),is_local);CHKERRQ(ierr);
1570af538404SDmitry Karpeev   s      = 0;
1571af538404SDmitry Karpeev   ystart = 0;
1572af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1573af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1574af538404SDmitry 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);
1575af538404SDmitry Karpeev     /* Vertical domain limits with an overlap. */
1576eec7e3faSDmitry Karpeev     y[0][0] = PetscMax(ystart - overlap,0);
1577eec7e3faSDmitry Karpeev     y[0][1] = PetscMin(ystart + maxheight + overlap,N);
1578eec7e3faSDmitry Karpeev     /* Vertical domain limits without an overlap. */
1579eec7e3faSDmitry Karpeev     y[1][0] = ystart;
1580eec7e3faSDmitry Karpeev     y[1][1] = PetscMin(ystart + maxheight,N);
1581eec7e3faSDmitry Karpeev     xstart  = 0;
1582af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1583af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1584af538404SDmitry 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);
1585af538404SDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1586eec7e3faSDmitry Karpeev       x[0][0] = PetscMax(xstart - overlap,0);
1587eec7e3faSDmitry Karpeev       x[0][1] = PetscMin(xstart + maxwidth + overlap,M);
1588eec7e3faSDmitry Karpeev       /* Horizontal domain limits without an overlap. */
1589eec7e3faSDmitry Karpeev       x[1][0] = xstart;
1590eec7e3faSDmitry Karpeev       x[1][1] = PetscMin(xstart+maxwidth,M);
15912bbb417fSDmitry Karpeev       /*
15922bbb417fSDmitry Karpeev          Determine whether this domain intersects this processor's ownership range of pc->pmat.
15932bbb417fSDmitry Karpeev          Do this twice: first for the domains with overlaps, and once without.
15942bbb417fSDmitry Karpeev          During the first pass create the subcommunicators, and use them on the second pass as well.
15952bbb417fSDmitry Karpeev       */
15962bbb417fSDmitry Karpeev       for (q = 0; q < 2; ++q) {
15972bbb417fSDmitry Karpeev         /*
15982bbb417fSDmitry Karpeev           domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
15992bbb417fSDmitry Karpeev           according to whether the domain with an overlap or without is considered.
16002bbb417fSDmitry Karpeev         */
16012bbb417fSDmitry Karpeev         xleft = x[q][0]; xright = x[q][1];
16022bbb417fSDmitry Karpeev         ylow  = y[q][0]; yhigh  = y[q][1];
1603c3518bceSDmitry Karpeev         PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1604af538404SDmitry Karpeev         nidx *= dof;
1605eec7e3faSDmitry Karpeev         n[q]  = nidx;
16062bbb417fSDmitry Karpeev         /*
1607eec7e3faSDmitry Karpeev          Based on the counted number of indices in the local domain *with an overlap*,
1608af538404SDmitry Karpeev          construct a subcommunicator of all the processors supporting this domain.
1609eec7e3faSDmitry Karpeev          Observe that a domain with an overlap might have nontrivial local support,
1610eec7e3faSDmitry Karpeev          while the domain without an overlap might not.  Hence, the decision to participate
1611eec7e3faSDmitry Karpeev          in the subcommunicator must be based on the domain with an overlap.
16122bbb417fSDmitry Karpeev          */
16132bbb417fSDmitry Karpeev         if (q == 0) {
16142fa5cd67SKarl Rupp           if (nidx) color = 1;
16152fa5cd67SKarl Rupp           else color = MPI_UNDEFINED;
16162fa5cd67SKarl Rupp 
16172bbb417fSDmitry Karpeev           ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr);
16182bbb417fSDmitry Karpeev         }
1619af538404SDmitry Karpeev         /*
1620eec7e3faSDmitry Karpeev          Proceed only if the number of local indices *with an overlap* is nonzero.
1621af538404SDmitry Karpeev          */
1622eec7e3faSDmitry Karpeev         if (n[0]) {
16232fa5cd67SKarl Rupp           if (q == 0) xis = is;
1624af538404SDmitry Karpeev           if (q == 1) {
1625af538404SDmitry Karpeev             /*
1626eec7e3faSDmitry Karpeev              The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1627af538404SDmitry Karpeev              Moreover, if the overlap is zero, the two ISs are identical.
1628af538404SDmitry Karpeev              */
1629b1a0a8a3SJed Brown             if (overlap == 0) {
1630eec7e3faSDmitry Karpeev               (*is_local)[s] = (*is)[s];
16312bbb417fSDmitry Karpeev               ierr           = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr);
16322bbb417fSDmitry Karpeev               continue;
1633d34fcf5fSBarry Smith             } else {
16346a4f0f73SDmitry Karpeev               xis     = is_local;
1635eec7e3faSDmitry Karpeev               subcomm = ((PetscObject)(*is)[s])->comm;
16362bbb417fSDmitry Karpeev             }
1637af538404SDmitry Karpeev           } /* if (q == 1) */
16380298fd71SBarry Smith           idx  = NULL;
16392bbb417fSDmitry Karpeev           ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1640eec7e3faSDmitry Karpeev           if (nidx) {
16412bbb417fSDmitry Karpeev             k = 0;
16422bbb417fSDmitry Karpeev             for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1643af538404SDmitry Karpeev               PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft;
1644af538404SDmitry Karpeev               PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright;
1645af538404SDmitry Karpeev               kk = dof*(M*jj + x0);
16462bbb417fSDmitry Karpeev               for (ii=x0; ii<x1; ++ii) {
16472bbb417fSDmitry Karpeev                 for (d = 0; d < dof; ++d) {
16482bbb417fSDmitry Karpeev                   idx[k++] = kk++;
1649b1a0a8a3SJed Brown                 }
1650b1a0a8a3SJed Brown               }
1651b1a0a8a3SJed Brown             }
1652eec7e3faSDmitry Karpeev           }
16536a4f0f73SDmitry Karpeev           ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr);
1654eec7e3faSDmitry Karpeev         }/* if (n[0]) */
16552bbb417fSDmitry Karpeev       }/* for (q = 0; q < 2; ++q) */
16562fa5cd67SKarl Rupp       if (n[0]) ++s;
1657af538404SDmitry Karpeev       xstart += maxwidth;
1658af538404SDmitry Karpeev     } /* for (i = 0; i < Mdomains; ++i) */
16592bbb417fSDmitry Karpeev     ystart += maxheight;
1660af538404SDmitry Karpeev   } /* for (j = 0; j < Ndomains; ++j) */
1661b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1662b1a0a8a3SJed Brown }
1663b1a0a8a3SJed Brown 
1664b1a0a8a3SJed Brown #undef __FUNCT__
166506b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains"
1666b1a0a8a3SJed Brown /*@C
166706b43e7eSDmitry Karpeev     PCGASMGetSubdomains - Gets the subdomains supported on this processor
166806b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
1669b1a0a8a3SJed Brown 
1670b1a0a8a3SJed Brown     Not Collective
1671b1a0a8a3SJed Brown 
1672b1a0a8a3SJed Brown     Input Parameter:
1673b1a0a8a3SJed Brown .   pc - the preconditioner context
1674b1a0a8a3SJed Brown 
1675b1a0a8a3SJed Brown     Output Parameters:
1676b1a0a8a3SJed Brown +   n   - the number of subdomains for this processor (default value = 1)
16770298fd71SBarry Smith .   iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL)
16780298fd71SBarry Smith -   ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL)
1679b1a0a8a3SJed Brown 
1680b1a0a8a3SJed Brown 
1681b1a0a8a3SJed Brown     Notes:
16826a4f0f73SDmitry Karpeev     The user is responsible for destroying the ISs and freeing the returned arrays.
1683b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
1684b1a0a8a3SJed Brown 
1685b1a0a8a3SJed Brown     Level: advanced
1686b1a0a8a3SJed Brown 
168706b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz
1688b1a0a8a3SJed Brown 
16896a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
169006b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubmatrices()
1691b1a0a8a3SJed Brown @*/
16926a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[])
1693b1a0a8a3SJed Brown {
1694f746d493SDmitry Karpeev   PC_GASM        *osm;
1695b1a0a8a3SJed Brown   PetscErrorCode ierr;
1696b1a0a8a3SJed Brown   PetscBool      match;
16976a4f0f73SDmitry Karpeev   PetscInt       i;
16985fd66863SKarl Rupp 
1699b1a0a8a3SJed Brown   PetscFunctionBegin;
1700b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1701251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1702*ce94432eSBarry Smith   if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1703f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1704ab3e923fSDmitry Karpeev   if (n) *n = osm->n;
17056a4f0f73SDmitry Karpeev   if (iis) {
17066a4f0f73SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(IS), iis);CHKERRQ(ierr);
17076a4f0f73SDmitry Karpeev   }
17086a4f0f73SDmitry Karpeev   if (ois) {
17096a4f0f73SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(IS), ois);CHKERRQ(ierr);
17106a4f0f73SDmitry Karpeev   }
17116a4f0f73SDmitry Karpeev   if (iis || ois) {
17126a4f0f73SDmitry Karpeev     for (i = 0; i < osm->n; ++i) {
17136a4f0f73SDmitry Karpeev       if (iis) (*iis)[i] = osm->iis[i];
17146a4f0f73SDmitry Karpeev       if (ois) (*ois)[i] = osm->ois[i];
17156a4f0f73SDmitry Karpeev     }
1716b1a0a8a3SJed Brown   }
1717b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1718b1a0a8a3SJed Brown }
1719b1a0a8a3SJed Brown 
1720b1a0a8a3SJed Brown #undef __FUNCT__
172106b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices"
1722b1a0a8a3SJed Brown /*@C
172306b43e7eSDmitry Karpeev     PCGASMGetSubmatrices - Gets the local submatrices (for this processor
1724b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1725b1a0a8a3SJed Brown 
1726b1a0a8a3SJed Brown     Not Collective
1727b1a0a8a3SJed Brown 
1728b1a0a8a3SJed Brown     Input Parameter:
1729b1a0a8a3SJed Brown .   pc - the preconditioner context
1730b1a0a8a3SJed Brown 
1731b1a0a8a3SJed Brown     Output Parameters:
1732b1a0a8a3SJed Brown +   n   - the number of matrices for this processor (default value = 1)
1733b1a0a8a3SJed Brown -   mat - the matrices
1734b1a0a8a3SJed Brown 
173506b43e7eSDmitry Karpeev     Notes: matrices returned by this routine have the same communicators as the index sets (IS)
173606b43e7eSDmitry Karpeev            used to define subdomains in PCGASMSetSubdomains(), or PETSC_COMM_SELF, if the
17376a4f0f73SDmitry Karpeev            subdomains were defined using PCGASMSetTotalSubdomains().
1738b1a0a8a3SJed Brown     Level: advanced
1739b1a0a8a3SJed Brown 
1740f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi
1741b1a0a8a3SJed Brown 
17426a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomain(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
174306b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains()
1744b1a0a8a3SJed Brown @*/
174506b43e7eSDmitry Karpeev PetscErrorCode  PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1746b1a0a8a3SJed Brown {
1747f746d493SDmitry Karpeev   PC_GASM        *osm;
1748b1a0a8a3SJed Brown   PetscErrorCode ierr;
1749b1a0a8a3SJed Brown   PetscBool      match;
1750b1a0a8a3SJed Brown 
1751b1a0a8a3SJed Brown   PetscFunctionBegin;
1752b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1753b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1754b1a0a8a3SJed Brown   if (mat) PetscValidPointer(mat,3);
1755*ce94432eSBarry Smith   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1756251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1757*ce94432eSBarry Smith   if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1758f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1759ab3e923fSDmitry Karpeev   if (n) *n = osm->n;
1760b1a0a8a3SJed Brown   if (mat) *mat = osm->pmat;
1761b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1762b1a0a8a3SJed Brown }
1763