xref: /petsc/src/ksp/pc/impls/gasm/gasm.c (revision d709ea831a5fc3427233a16b23cda3f1f7499f01)
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 */
29*d709ea83SDmitry Karpeev   PetscBool  dm_subdomains;          /* 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;
44ce94432eSBarry 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;
113ce94432eSBarry 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); };
118ce94432eSBarry 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     }
136ce94432eSBarry 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;
162ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr);
163ce94432eSBarry 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;
166ce94432eSBarry 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. */
206ce94432eSBarry 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       }
236ce94432eSBarry 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;
273ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
274ce94432eSBarry 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 */
289*d709ea83SDmitry Karpeev       if (osm->dm_subdomains && 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;
308ce94432eSBarry 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;
325ce94432eSBarry 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);
361ce94432eSBarry 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);
364ce94432eSBarry 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       }
408ce94432eSBarry Smith       ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, iidx,  PETSC_OWN_POINTER, &giis);CHKERRQ(ierr);
409ce94432eSBarry 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);
664*d709ea83SDmitry Karpeev   ierr = PetscOptionsBool("-pc_gasm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCGASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr);
6656a4f0f73SDmitry Karpeev   ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
66665db9045SDmitry Karpeev   if (flg) {
6676a4f0f73SDmitry Karpeev     osm->create_local = PETSC_TRUE;
66865db9045SDmitry 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);
66965db9045SDmitry Karpeev     ierr = PCGASMSetTotalSubdomains(pc,blocks,osm->create_local);CHKERRQ(ierr);
670*d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
67165db9045SDmitry Karpeev   }
67206b43e7eSDmitry Karpeev   ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
67365db9045SDmitry Karpeev   if (flg) {
67465db9045SDmitry Karpeev     ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr);
675*d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
67665db9045SDmitry Karpeev   }
677b1a0a8a3SJed Brown   flg  = PETSC_FALSE;
678f746d493SDmitry Karpeev   ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr);
679f746d493SDmitry Karpeev   if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr); }
680b1a0a8a3SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
681b1a0a8a3SJed Brown   PetscFunctionReturn(0);
682b1a0a8a3SJed Brown }
683b1a0a8a3SJed Brown 
684b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/
685b1a0a8a3SJed Brown 
686b1a0a8a3SJed Brown EXTERN_C_BEGIN
687b1a0a8a3SJed Brown #undef __FUNCT__
68806b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM"
689a35b7b57SDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[])
690b1a0a8a3SJed Brown {
691f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
692b1a0a8a3SJed Brown   PetscErrorCode ierr;
693b1a0a8a3SJed Brown   PetscInt       i;
694b1a0a8a3SJed Brown 
695b1a0a8a3SJed Brown   PetscFunctionBegin;
69669ca1f37SDmitry Karpeev   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, n = %D",n);
697ce94432eSBarry Smith   if (pc->setupcalled && (n != osm->n || iis || ois)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp().");
698b1a0a8a3SJed Brown 
699b1a0a8a3SJed Brown   if (!pc->setupcalled) {
70093c1ec32SDmitry Karpeev     osm->n   = n;
7016a4f0f73SDmitry Karpeev     osm->ois = 0;
7026a4f0f73SDmitry Karpeev     osm->iis = 0;
703a35b7b57SDmitry Karpeev     if (ois) {
704a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);}
705b1a0a8a3SJed Brown     }
706a35b7b57SDmitry Karpeev     if (iis) {
707a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);}
708b1a0a8a3SJed Brown     }
7096a4f0f73SDmitry Karpeev     ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr);
710a35b7b57SDmitry Karpeev     if (ois) {
7116a4f0f73SDmitry Karpeev       ierr = PetscMalloc(n*sizeof(IS),&osm->ois);CHKERRQ(ierr);
7122fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->ois[i] = ois[i];
7136a4f0f73SDmitry Karpeev       /* Flag indicating that the user has set outer subdomains, so PCGASM should not increase their size. */
714b1a0a8a3SJed Brown       osm->overlap = -1;
715a35b7b57SDmitry Karpeev       if (!iis) {
7166a4f0f73SDmitry Karpeev         ierr = PetscMalloc(n*sizeof(IS),&osm->iis);CHKERRQ(ierr);
717a35b7b57SDmitry Karpeev         for (i=0; i<n; i++) {
718a35b7b57SDmitry Karpeev           for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);}
719a35b7b57SDmitry Karpeev           osm->iis[i] = ois[i];
720a35b7b57SDmitry Karpeev         }
721a35b7b57SDmitry Karpeev       }
722a35b7b57SDmitry Karpeev     }
723a35b7b57SDmitry Karpeev     if (iis) {
724a35b7b57SDmitry Karpeev       ierr = PetscMalloc(n*sizeof(IS),&osm->iis);CHKERRQ(ierr);
7252fa5cd67SKarl Rupp       for (i=0; i<n; i++) osm->iis[i] = iis[i];
726a35b7b57SDmitry Karpeev       if (!ois) {
727a35b7b57SDmitry Karpeev         ierr = PetscMalloc(n*sizeof(IS),&osm->ois);CHKERRQ(ierr);
728a35b7b57SDmitry Karpeev         for (i=0; i<n; i++) {
729a35b7b57SDmitry Karpeev           for (i=0; i<n; i++) {
730a35b7b57SDmitry Karpeev             ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);
731a35b7b57SDmitry Karpeev             osm->ois[i] = iis[i];
732a35b7b57SDmitry Karpeev           }
733a35b7b57SDmitry Karpeev         }
734a35b7b57SDmitry Karpeev         if (osm->overlap > 0) {
735a35b7b57SDmitry Karpeev           /* Extend the "overlapping" regions by a number of steps */
736a35b7b57SDmitry Karpeev           ierr = MatIncreaseOverlap(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr);
737a35b7b57SDmitry Karpeev         }
738a35b7b57SDmitry Karpeev       }
739b1a0a8a3SJed Brown     }
740b1a0a8a3SJed Brown   }
741b1a0a8a3SJed Brown   PetscFunctionReturn(0);
742b1a0a8a3SJed Brown }
743b1a0a8a3SJed Brown EXTERN_C_END
744b1a0a8a3SJed Brown 
745b1a0a8a3SJed Brown EXTERN_C_BEGIN
746b1a0a8a3SJed Brown #undef __FUNCT__
7476a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains_GASM"
7480adebc6cSBarry Smith PetscErrorCode  PCGASMSetTotalSubdomains_GASM(PC pc,PetscInt N, PetscBool create_local)
7490adebc6cSBarry Smith {
750f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
751b1a0a8a3SJed Brown   PetscErrorCode ierr;
752b1a0a8a3SJed Brown   PetscMPIInt    rank,size;
753b1a0a8a3SJed Brown   PetscInt       n;
7546a4f0f73SDmitry Karpeev   PetscInt       Nmin, Nmax;
7555fd66863SKarl Rupp 
756b1a0a8a3SJed Brown   PetscFunctionBegin;
757ce94432eSBarry Smith   if (!create_local) SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "No suppor for autocreation of nonlocal subdomains.");
758ce94432eSBarry Smith   if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be > 0, N = %D",N);
759ce94432eSBarry Smith   ierr = MPI_Allreduce(&N,&Nmin,1,MPIU_INT,MPIU_MIN,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
760ce94432eSBarry Smith   ierr = MPI_Allreduce(&N,&Nmax,1,MPIU_INT,MPIU_MAX,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
761ce94432eSBarry 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);
762b1a0a8a3SJed Brown 
7636a4f0f73SDmitry Karpeev   osm->create_local = create_local;
764b1a0a8a3SJed Brown   /*
765b1a0a8a3SJed Brown      Split the subdomains equally among all processors
766b1a0a8a3SJed Brown   */
767ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr);
768ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
769b1a0a8a3SJed Brown   n    = N/size + ((N % size) > rank);
7706a4f0f73SDmitry 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);
771f746d493SDmitry Karpeev   if (pc->setupcalled && n != osm->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before PCSetUp().");
772b1a0a8a3SJed Brown   if (!pc->setupcalled) {
7736a4f0f73SDmitry Karpeev     ierr      = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr);
77493c1ec32SDmitry Karpeev     osm->N    = N;
775f746d493SDmitry Karpeev     osm->n    = n;
77606b43e7eSDmitry Karpeev     osm->nmax = N/size + ((N%size) ? 1 : 0);
7776a4f0f73SDmitry Karpeev     osm->ois  = 0;
7786a4f0f73SDmitry Karpeev     osm->iis  = 0;
779b1a0a8a3SJed Brown   }
780b1a0a8a3SJed Brown   PetscFunctionReturn(0);
78106b43e7eSDmitry Karpeev }
782b1a0a8a3SJed Brown EXTERN_C_END
783b1a0a8a3SJed Brown 
784b1a0a8a3SJed Brown EXTERN_C_BEGIN
785b1a0a8a3SJed Brown #undef __FUNCT__
786f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM"
7877087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
788b1a0a8a3SJed Brown {
789f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
790b1a0a8a3SJed Brown 
791b1a0a8a3SJed Brown   PetscFunctionBegin;
792ce94432eSBarry Smith   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
793ce94432eSBarry Smith   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
7942fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
795b1a0a8a3SJed Brown   PetscFunctionReturn(0);
796b1a0a8a3SJed Brown }
797b1a0a8a3SJed Brown EXTERN_C_END
798b1a0a8a3SJed Brown 
799b1a0a8a3SJed Brown EXTERN_C_BEGIN
800b1a0a8a3SJed Brown #undef __FUNCT__
801f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM"
8027087cfbeSBarry Smith PetscErrorCode  PCGASMSetType_GASM(PC pc,PCGASMType type)
803b1a0a8a3SJed Brown {
804f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
805b1a0a8a3SJed Brown 
806b1a0a8a3SJed Brown   PetscFunctionBegin;
807b1a0a8a3SJed Brown   osm->type     = type;
808b1a0a8a3SJed Brown   osm->type_set = PETSC_TRUE;
809b1a0a8a3SJed Brown   PetscFunctionReturn(0);
810b1a0a8a3SJed Brown }
811b1a0a8a3SJed Brown EXTERN_C_END
812b1a0a8a3SJed Brown 
813b1a0a8a3SJed Brown EXTERN_C_BEGIN
814b1a0a8a3SJed Brown #undef __FUNCT__
815f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM"
8167087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort)
817b1a0a8a3SJed Brown {
818f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
819b1a0a8a3SJed Brown 
820b1a0a8a3SJed Brown   PetscFunctionBegin;
821b1a0a8a3SJed Brown   osm->sort_indices = doSort;
822b1a0a8a3SJed Brown   PetscFunctionReturn(0);
823b1a0a8a3SJed Brown }
824b1a0a8a3SJed Brown EXTERN_C_END
825b1a0a8a3SJed Brown 
826b1a0a8a3SJed Brown EXTERN_C_BEGIN
827b1a0a8a3SJed Brown #undef __FUNCT__
828f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM"
82944fe18e5SDmitry Karpeev /*
83044fe18e5SDmitry Karpeev    FIX: This routine might need to be modified once multiple ranks per subdomain are allowed.
83144fe18e5SDmitry Karpeev         In particular, it would upset the global subdomain number calculation.
83244fe18e5SDmitry Karpeev */
8337087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
834b1a0a8a3SJed Brown {
835f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
836b1a0a8a3SJed Brown   PetscErrorCode ierr;
837b1a0a8a3SJed Brown 
838b1a0a8a3SJed Brown   PetscFunctionBegin;
839ce94432eSBarry 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");
840b1a0a8a3SJed Brown 
8412fa5cd67SKarl Rupp   if (n) *n = osm->n;
842ab3e923fSDmitry Karpeev   if (first) {
843ce94432eSBarry Smith     ierr    = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr);
844ab3e923fSDmitry Karpeev     *first -= osm->n;
845b1a0a8a3SJed Brown   }
846b1a0a8a3SJed Brown   if (ksp) {
847b1a0a8a3SJed Brown     /* Assume that local solves are now different; not necessarily
84806b43e7eSDmitry Karpeev        true, though!  This flag is used only for PCView_GASM() */
849b1a0a8a3SJed Brown     *ksp                        = osm->ksp;
8506a4f0f73SDmitry Karpeev     osm->same_subdomain_solvers = PETSC_FALSE;
851b1a0a8a3SJed Brown   }
852b1a0a8a3SJed Brown   PetscFunctionReturn(0);
853ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */
854b1a0a8a3SJed Brown EXTERN_C_END
855b1a0a8a3SJed Brown 
856b1a0a8a3SJed Brown 
857b1a0a8a3SJed Brown #undef __FUNCT__
85806b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains"
859b1a0a8a3SJed Brown /*@C
86006b43e7eSDmitry Karpeev     PCGASMSetSubdomains - Sets the subdomains for this processor
86106b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
862b1a0a8a3SJed Brown 
863b1a0a8a3SJed Brown     Collective on PC
864b1a0a8a3SJed Brown 
865b1a0a8a3SJed Brown     Input Parameters:
866b1a0a8a3SJed Brown +   pc  - the preconditioner context
86706b43e7eSDmitry Karpeev .   n   - the number of subdomains for this processor
8686a4f0f73SDmitry Karpeev .   iis - the index sets that define this processor's local inner subdomains
8690298fd71SBarry Smith          (or NULL for PETSc to determine subdomains)
8706a4f0f73SDmitry Karpeev -   ois- the index sets that define this processor's local outer subdomains
8710298fd71SBarry Smith          (or NULL to use the same as iis)
872b1a0a8a3SJed Brown 
873b1a0a8a3SJed Brown     Notes:
87406b43e7eSDmitry Karpeev     The IS indices use the parallel, global numbering of the vector entries.
8756a4f0f73SDmitry Karpeev     Inner subdomains are those where the correction is applied.
8766a4f0f73SDmitry Karpeev     Outer subdomains are those where the residual necessary to obtain the
8776a4f0f73SDmitry Karpeev     corrections is obtained (see PCGASMType for the use of inner/outer subdomains).
8786a4f0f73SDmitry Karpeev     Both inner and outer subdomains can extend over several processors.
8796a4f0f73SDmitry Karpeev     This processor's portion of a subdomain is known as a local subdomain.
8806a4f0f73SDmitry Karpeev 
8816a4f0f73SDmitry Karpeev     By default the GASM preconditioner uses 1 (local) subdomain per processor.
8826a4f0f73SDmitry Karpeev     Use PCGASMSetTotalSubdomains() to set the total number of subdomains across
8836a4f0f73SDmitry Karpeev     all processors that PCGASM will create automatically, and to specify whether
8846a4f0f73SDmitry Karpeev     they should be local or not.
8856a4f0f73SDmitry Karpeev 
886b1a0a8a3SJed Brown 
887b1a0a8a3SJed Brown     Level: advanced
888b1a0a8a3SJed Brown 
88906b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz
890b1a0a8a3SJed Brown 
8916a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
89206b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
893b1a0a8a3SJed Brown @*/
8946a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[])
895b1a0a8a3SJed Brown {
896b1a0a8a3SJed Brown   PetscErrorCode ierr;
897b1a0a8a3SJed Brown 
898b1a0a8a3SJed Brown   PetscFunctionBegin;
899b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9006a4f0f73SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr);
901b1a0a8a3SJed Brown   PetscFunctionReturn(0);
902b1a0a8a3SJed Brown }
903b1a0a8a3SJed Brown 
904b1a0a8a3SJed Brown #undef __FUNCT__
9056a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains"
906b1a0a8a3SJed Brown /*@C
9076a4f0f73SDmitry Karpeev     PCGASMSetTotalSubdomains - Sets the total number of subdomains to use in the generalized additive
9086a4f0f73SDmitry Karpeev     Schwarz preconditioner.  The number of subdomains is cumulative across all processors in pc's
9096a4f0f73SDmitry Karpeev     communicator. Either all or no processors in the PC communicator must call this routine with
9106a4f0f73SDmitry Karpeev     the same N.  The subdomains will be created automatically during PCSetUp().
911b1a0a8a3SJed Brown 
912b1a0a8a3SJed Brown     Collective on PC
913b1a0a8a3SJed Brown 
914b1a0a8a3SJed Brown     Input Parameters:
915b1a0a8a3SJed Brown +   pc           - the preconditioner context
9166a4f0f73SDmitry Karpeev .   N            - the total number of subdomains cumulative across all processors
9176a4f0f73SDmitry Karpeev -   create_local - whether the subdomains to be created are to be local
918b1a0a8a3SJed Brown 
919b1a0a8a3SJed Brown     Options Database Key:
9206a4f0f73SDmitry Karpeev     To set the total number of subdomains and let PCGASM autocreate them, rather than specify the index sets, use the following options:
9216a4f0f73SDmitry Karpeev +    -pc_gasm_total_subdomains <n>                  - sets the total number of subdomains to be autocreated by PCGASM
9226a4f0f73SDmitry Karpeev -    -pc_gasm_subdomains_create_local <true|false>  - whether autocreated subdomains should be local or not (default is true)
923b1a0a8a3SJed Brown 
92406b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.
925b1a0a8a3SJed Brown 
9266a4f0f73SDmitry Karpeev 
9276a4f0f73SDmitry Karpeev     Use PCGASMSetSubdomains() to set subdomains explicitly or to set different numbers
9286a4f0f73SDmitry Karpeev     of subdomains per processor.
929b1a0a8a3SJed Brown 
930b1a0a8a3SJed Brown     Level: advanced
931b1a0a8a3SJed Brown 
932f746d493SDmitry Karpeev .keywords: PC, GASM, set, total, global, subdomains, additive Schwarz
933b1a0a8a3SJed Brown 
93406b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
935f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
936b1a0a8a3SJed Brown @*/
9376a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetTotalSubdomains(PC pc,PetscInt N, PetscBool create_local)
938b1a0a8a3SJed Brown {
939b1a0a8a3SJed Brown   PetscErrorCode ierr;
940b1a0a8a3SJed Brown 
941b1a0a8a3SJed Brown   PetscFunctionBegin;
942b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9436a4f0f73SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetTotalSubdomains_C",(PC,PetscInt,PetscBool),(pc,N,create_local));CHKERRQ(ierr);
944b1a0a8a3SJed Brown   PetscFunctionReturn(0);
945b1a0a8a3SJed Brown }
946b1a0a8a3SJed Brown 
947b1a0a8a3SJed Brown #undef __FUNCT__
948f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap"
949b1a0a8a3SJed Brown /*@
950f746d493SDmitry Karpeev     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
951b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
952b1a0a8a3SJed Brown     PC communicator must call this routine.
953b1a0a8a3SJed Brown 
954b1a0a8a3SJed Brown     Logically Collective on PC
955b1a0a8a3SJed Brown 
956b1a0a8a3SJed Brown     Input Parameters:
957b1a0a8a3SJed Brown +   pc  - the preconditioner context
958b1a0a8a3SJed Brown -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
959b1a0a8a3SJed Brown 
960b1a0a8a3SJed Brown     Options Database Key:
96106b43e7eSDmitry Karpeev .   -pc_gasm_overlap <overlap> - Sets overlap
962b1a0a8a3SJed Brown 
963b1a0a8a3SJed Brown     Notes:
96406b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.  To use
9656a4f0f73SDmitry Karpeev     multiple subdomain per perocessor, see PCGASMSetTotalSubdomains() or
9666a4f0f73SDmitry Karpeev     PCGASMSetSubdomains() (and the option -pc_gasm_total_subdomains <n>).
967b1a0a8a3SJed Brown 
968b1a0a8a3SJed Brown     The overlap defaults to 1, so if one desires that no additional
969b1a0a8a3SJed Brown     overlap be computed beyond what may have been set with a call to
9706a4f0f73SDmitry Karpeev     PCGASMSetTotalSubdomains() or PCGASMSetSubdomains(), then ovl
971b1a0a8a3SJed Brown     must be set to be 0.  In particular, if one does not explicitly set
97206b43e7eSDmitry Karpeev     the subdomains in application code, then all overlap would be computed
973f746d493SDmitry Karpeev     internally by PETSc, and using an overlap of 0 would result in an GASM
974b1a0a8a3SJed Brown     variant that is equivalent to the block Jacobi preconditioner.
975b1a0a8a3SJed Brown 
976b1a0a8a3SJed Brown     Note that one can define initial index sets with any overlap via
97706b43e7eSDmitry Karpeev     PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows
97806b43e7eSDmitry Karpeev     PETSc to extend that overlap further, if desired.
979b1a0a8a3SJed Brown 
980b1a0a8a3SJed Brown     Level: intermediate
981b1a0a8a3SJed Brown 
982f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap
983b1a0a8a3SJed Brown 
9846a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
98506b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
986b1a0a8a3SJed Brown @*/
9877087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap(PC pc,PetscInt ovl)
988b1a0a8a3SJed Brown {
989b1a0a8a3SJed Brown   PetscErrorCode ierr;
990b1a0a8a3SJed Brown 
991b1a0a8a3SJed Brown   PetscFunctionBegin;
992b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
993b1a0a8a3SJed Brown   PetscValidLogicalCollectiveInt(pc,ovl,2);
994f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
995b1a0a8a3SJed Brown   PetscFunctionReturn(0);
996b1a0a8a3SJed Brown }
997b1a0a8a3SJed Brown 
998b1a0a8a3SJed Brown #undef __FUNCT__
999f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType"
1000b1a0a8a3SJed Brown /*@
1001f746d493SDmitry Karpeev     PCGASMSetType - Sets the type of restriction and interpolation used
1002b1a0a8a3SJed Brown     for local problems in the additive Schwarz method.
1003b1a0a8a3SJed Brown 
1004b1a0a8a3SJed Brown     Logically Collective on PC
1005b1a0a8a3SJed Brown 
1006b1a0a8a3SJed Brown     Input Parameters:
1007b1a0a8a3SJed Brown +   pc  - the preconditioner context
1008f746d493SDmitry Karpeev -   type - variant of GASM, one of
1009b1a0a8a3SJed Brown .vb
1010f746d493SDmitry Karpeev       PC_GASM_BASIC       - full interpolation and restriction
1011f746d493SDmitry Karpeev       PC_GASM_RESTRICT    - full restriction, local processor interpolation
1012f746d493SDmitry Karpeev       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
1013f746d493SDmitry Karpeev       PC_GASM_NONE        - local processor restriction and interpolation
1014b1a0a8a3SJed Brown .ve
1015b1a0a8a3SJed Brown 
1016b1a0a8a3SJed Brown     Options Database Key:
1017f746d493SDmitry Karpeev .   -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1018b1a0a8a3SJed Brown 
1019b1a0a8a3SJed Brown     Level: intermediate
1020b1a0a8a3SJed Brown 
1021f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1022b1a0a8a3SJed Brown 
10236a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1024f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1025b1a0a8a3SJed Brown @*/
10267087cfbeSBarry Smith PetscErrorCode  PCGASMSetType(PC pc,PCGASMType type)
1027b1a0a8a3SJed Brown {
1028b1a0a8a3SJed Brown   PetscErrorCode ierr;
1029b1a0a8a3SJed Brown 
1030b1a0a8a3SJed Brown   PetscFunctionBegin;
1031b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1032b1a0a8a3SJed Brown   PetscValidLogicalCollectiveEnum(pc,type,2);
1033f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr);
1034b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1035b1a0a8a3SJed Brown }
1036b1a0a8a3SJed Brown 
1037b1a0a8a3SJed Brown #undef __FUNCT__
1038f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices"
1039b1a0a8a3SJed Brown /*@
1040f746d493SDmitry Karpeev     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
1041b1a0a8a3SJed Brown 
1042b1a0a8a3SJed Brown     Logically Collective on PC
1043b1a0a8a3SJed Brown 
1044b1a0a8a3SJed Brown     Input Parameters:
1045b1a0a8a3SJed Brown +   pc  - the preconditioner context
1046b1a0a8a3SJed Brown -   doSort - sort the subdomain indices
1047b1a0a8a3SJed Brown 
1048b1a0a8a3SJed Brown     Level: intermediate
1049b1a0a8a3SJed Brown 
1050f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1051b1a0a8a3SJed Brown 
10526a4f0f73SDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(),
1053f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1054b1a0a8a3SJed Brown @*/
10557087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices(PC pc,PetscBool doSort)
1056b1a0a8a3SJed Brown {
1057b1a0a8a3SJed Brown   PetscErrorCode ierr;
1058b1a0a8a3SJed Brown 
1059b1a0a8a3SJed Brown   PetscFunctionBegin;
1060b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1061b1a0a8a3SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
1062f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
1063b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1064b1a0a8a3SJed Brown }
1065b1a0a8a3SJed Brown 
1066b1a0a8a3SJed Brown #undef __FUNCT__
1067f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP"
1068b1a0a8a3SJed Brown /*@C
1069f746d493SDmitry Karpeev    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
1070b1a0a8a3SJed Brown    this processor.
1071b1a0a8a3SJed Brown 
1072b1a0a8a3SJed Brown    Collective on PC iff first_local is requested
1073b1a0a8a3SJed Brown 
1074b1a0a8a3SJed Brown    Input Parameter:
1075b1a0a8a3SJed Brown .  pc - the preconditioner context
1076b1a0a8a3SJed Brown 
1077b1a0a8a3SJed Brown    Output Parameters:
10780298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
10790298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
10800298fd71SBarry Smith                  all processors must request or all must pass NULL
1081b1a0a8a3SJed Brown -  ksp - the array of KSP contexts
1082b1a0a8a3SJed Brown 
1083b1a0a8a3SJed Brown    Note:
1084f746d493SDmitry Karpeev    After PCGASMGetSubKSP() the array of KSPes is not to be freed
1085b1a0a8a3SJed Brown 
1086b1a0a8a3SJed Brown    Currently for some matrix implementations only 1 block per processor
1087b1a0a8a3SJed Brown    is supported.
1088b1a0a8a3SJed Brown 
1089f746d493SDmitry Karpeev    You must call KSPSetUp() before calling PCGASMGetSubKSP().
1090b1a0a8a3SJed Brown 
1091b1a0a8a3SJed Brown    Level: advanced
1092b1a0a8a3SJed Brown 
1093f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context
1094b1a0a8a3SJed Brown 
10956a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap(),
1096f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(),
1097b1a0a8a3SJed Brown @*/
10987087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1099b1a0a8a3SJed Brown {
1100b1a0a8a3SJed Brown   PetscErrorCode ierr;
1101b1a0a8a3SJed Brown 
1102b1a0a8a3SJed Brown   PetscFunctionBegin;
1103b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1104f746d493SDmitry Karpeev   ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
1105b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1106b1a0a8a3SJed Brown }
1107b1a0a8a3SJed Brown 
1108b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/
1109b1a0a8a3SJed Brown /*MC
1110f746d493SDmitry Karpeev    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1111b1a0a8a3SJed Brown            its own KSP object.
1112b1a0a8a3SJed Brown 
1113b1a0a8a3SJed Brown    Options Database Keys:
111406b43e7eSDmitry Karpeev +  -pc_gasm_total_block_count <n> - Sets total number of local subdomains (known as blocks) to be distributed among processors
111506b43e7eSDmitry Karpeev .  -pc_gasm_view_subdomains       - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view
111606b43e7eSDmitry Karpeev .  -pc_gasm_print_subdomains      - activates the printing of subdomain indices in PCSetUp()
111706b43e7eSDmitry Karpeev .  -pc_gasm_overlap <ovl>         - Sets overlap by which to (automatically) extend local subdomains
1118f746d493SDmitry Karpeev -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1119b1a0a8a3SJed Brown 
1120b1a0a8a3SJed Brown      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
1121f746d493SDmitry Karpeev       will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
1122f746d493SDmitry Karpeev       -pc_gasm_type basic to use the standard GASM.
1123b1a0a8a3SJed Brown 
1124b1a0a8a3SJed Brown    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
1125b1a0a8a3SJed Brown      than one processor. Defaults to one block per processor.
1126b1a0a8a3SJed Brown 
1127b1a0a8a3SJed Brown      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1128b1a0a8a3SJed Brown         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
1129b1a0a8a3SJed Brown 
1130f746d493SDmitry Karpeev      To set the options on the solvers separate for each block call PCGASMGetSubKSP()
1131b1a0a8a3SJed Brown          and set the options directly on the resulting KSP object (you can access its PC
1132b1a0a8a3SJed Brown          with KSPGetPC())
1133b1a0a8a3SJed Brown 
1134b1a0a8a3SJed Brown 
1135b1a0a8a3SJed Brown    Level: beginner
1136b1a0a8a3SJed Brown 
1137b1a0a8a3SJed Brown    Concepts: additive Schwarz method
1138b1a0a8a3SJed Brown 
1139b1a0a8a3SJed Brown     References:
1140b1a0a8a3SJed Brown     An additive variant of the Schwarz alternating method for the case of many subregions
1141b1a0a8a3SJed Brown     M Dryja, OB Widlund - Courant Institute, New York University Technical report
1142b1a0a8a3SJed Brown 
1143b1a0a8a3SJed Brown     Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1144b1a0a8a3SJed Brown     Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X.
1145b1a0a8a3SJed Brown 
1146b1a0a8a3SJed Brown .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
114706b43e7eSDmitry Karpeev            PCBJACOBI, PCGASMSetUseTrueLocal(), PCGASMGetSubKSP(), PCGASMSetSubdomains(),
11486a4f0f73SDmitry Karpeev            PCGASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType()
1149b1a0a8a3SJed Brown 
1150b1a0a8a3SJed Brown M*/
1151b1a0a8a3SJed Brown 
1152b1a0a8a3SJed Brown EXTERN_C_BEGIN
1153b1a0a8a3SJed Brown #undef __FUNCT__
1154f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM"
11557087cfbeSBarry Smith PetscErrorCode  PCCreate_GASM(PC pc)
1156b1a0a8a3SJed Brown {
1157b1a0a8a3SJed Brown   PetscErrorCode ierr;
1158f746d493SDmitry Karpeev   PC_GASM        *osm;
1159b1a0a8a3SJed Brown 
1160b1a0a8a3SJed Brown   PetscFunctionBegin;
1161f746d493SDmitry Karpeev   ierr = PetscNewLog(pc,PC_GASM,&osm);CHKERRQ(ierr);
11622fa5cd67SKarl Rupp 
1163ab3e923fSDmitry Karpeev   osm->N                      = PETSC_DECIDE;
116406b43e7eSDmitry Karpeev   osm->n                      = PETSC_DECIDE;
1165ab3e923fSDmitry Karpeev   osm->nmax                   = 0;
1166b1a0a8a3SJed Brown   osm->overlap                = 1;
1167b1a0a8a3SJed Brown   osm->ksp                    = 0;
11686a4f0f73SDmitry Karpeev   osm->gorestriction          = 0;
11696a4f0f73SDmitry Karpeev   osm->girestriction          = 0;
1170ab3e923fSDmitry Karpeev   osm->gx                     = 0;
1171ab3e923fSDmitry Karpeev   osm->gy                     = 0;
1172b1a0a8a3SJed Brown   osm->x                      = 0;
1173b1a0a8a3SJed Brown   osm->y                      = 0;
11746a4f0f73SDmitry Karpeev   osm->ois                    = 0;
11756a4f0f73SDmitry Karpeev   osm->iis                    = 0;
1176b1a0a8a3SJed Brown   osm->pmat                   = 0;
1177f746d493SDmitry Karpeev   osm->type                   = PC_GASM_RESTRICT;
11786a4f0f73SDmitry Karpeev   osm->same_subdomain_solvers = PETSC_TRUE;
1179b1a0a8a3SJed Brown   osm->sort_indices           = PETSC_TRUE;
1180*d709ea83SDmitry Karpeev   osm->dm_subdomains          = PETSC_FALSE;
1181b1a0a8a3SJed Brown 
1182b1a0a8a3SJed Brown   pc->data                 = (void*)osm;
1183f746d493SDmitry Karpeev   pc->ops->apply           = PCApply_GASM;
1184f746d493SDmitry Karpeev   pc->ops->applytranspose  = PCApplyTranspose_GASM;
1185f746d493SDmitry Karpeev   pc->ops->setup           = PCSetUp_GASM;
1186a06653b4SBarry Smith   pc->ops->reset           = PCReset_GASM;
1187f746d493SDmitry Karpeev   pc->ops->destroy         = PCDestroy_GASM;
1188f746d493SDmitry Karpeev   pc->ops->setfromoptions  = PCSetFromOptions_GASM;
1189f746d493SDmitry Karpeev   pc->ops->setuponblocks   = PCSetUpOnBlocks_GASM;
1190f746d493SDmitry Karpeev   pc->ops->view            = PCView_GASM;
1191b1a0a8a3SJed Brown   pc->ops->applyrichardson = 0;
1192b1a0a8a3SJed Brown 
119306b43e7eSDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSubdomains_C","PCGASMSetSubdomains_GASM",
119406b43e7eSDmitry Karpeev                                            PCGASMSetSubdomains_GASM);CHKERRQ(ierr);
11956a4f0f73SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetTotalSubdomains_C","PCGASMSetTotalSubdomains_GASM",
11966a4f0f73SDmitry Karpeev                                            PCGASMSetTotalSubdomains_GASM);CHKERRQ(ierr);
1197f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetOverlap_C","PCGASMSetOverlap_GASM",
1198f746d493SDmitry Karpeev                                            PCGASMSetOverlap_GASM);CHKERRQ(ierr);
1199f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetType_C","PCGASMSetType_GASM",
1200f746d493SDmitry Karpeev                                            PCGASMSetType_GASM);CHKERRQ(ierr);
1201f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSortIndices_C","PCGASMSetSortIndices_GASM",
1202f746d493SDmitry Karpeev                                            PCGASMSetSortIndices_GASM);CHKERRQ(ierr);
1203f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMGetSubKSP_C","PCGASMGetSubKSP_GASM",
1204f746d493SDmitry Karpeev                                            PCGASMGetSubKSP_GASM);CHKERRQ(ierr);
1205b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1206b1a0a8a3SJed Brown }
1207b1a0a8a3SJed Brown EXTERN_C_END
1208b1a0a8a3SJed Brown 
1209b1a0a8a3SJed Brown 
1210b1a0a8a3SJed Brown #undef __FUNCT__
121106b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains"
1212b1a0a8a3SJed Brown /*@C
121306b43e7eSDmitry Karpeev    PCGASMCreateLocalSubdomains - Creates n local index sets for the overlapping
121406b43e7eSDmitry Karpeev    Schwarz preconditioner for a any problem based on its matrix.
1215b1a0a8a3SJed Brown 
1216b1a0a8a3SJed Brown    Collective
1217b1a0a8a3SJed Brown 
1218b1a0a8a3SJed Brown    Input Parameters:
1219b1a0a8a3SJed Brown +  A       - The global matrix operator
12206a4f0f73SDmitry Karpeev .  overlap - amount of overlap in outer subdomains
122106b43e7eSDmitry Karpeev -  n       - the number of local subdomains
1222b1a0a8a3SJed Brown 
1223b1a0a8a3SJed Brown    Output Parameters:
12246a4f0f73SDmitry Karpeev +  iis - the array of index sets defining the local inner subdomains (on which the correction is applied)
12256a4f0f73SDmitry Karpeev -  ois - the array of index sets defining the local outer subdomains (on which the residual is computed)
1226b1a0a8a3SJed Brown 
1227b1a0a8a3SJed Brown    Level: advanced
1228b1a0a8a3SJed Brown 
12296a4f0f73SDmitry Karpeev    Note: this generates n nonoverlapping local inner subdomains on PETSC_COMM_SELF;
12306a4f0f73SDmitry Karpeev          PCGASM will generate the overlap from these if you use them in PCGASMSetSubdomains() and set a
12316a4f0f73SDmitry Karpeev          nonzero overlap with PCGASMSetOverlap()
1232b1a0a8a3SJed Brown 
1233b1a0a8a3SJed Brown     In the Fortran version you must provide the array outis[] already allocated of length n.
1234b1a0a8a3SJed Brown 
1235f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1236b1a0a8a3SJed Brown 
123706b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains()
1238b1a0a8a3SJed Brown @*/
12396a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt overlap, PetscInt n, IS *iis[], IS *ois[])
1240b1a0a8a3SJed Brown {
1241b1a0a8a3SJed Brown   MatPartitioning mpart;
1242b1a0a8a3SJed Brown   const char      *prefix;
124311bd1e4dSLisandro Dalcin   PetscErrorCode  (*f)(Mat,MatReuse,Mat*);
1244b1a0a8a3SJed Brown   PetscMPIInt     size;
1245b1a0a8a3SJed Brown   PetscInt        i,j,rstart,rend,bs;
124611bd1e4dSLisandro Dalcin   PetscBool       isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
12470298fd71SBarry Smith   Mat             Ad     = NULL, adj;
1248b1a0a8a3SJed Brown   IS              ispart,isnumb,*is;
1249b1a0a8a3SJed Brown   PetscErrorCode  ierr;
1250b1a0a8a3SJed Brown 
1251b1a0a8a3SJed Brown   PetscFunctionBegin;
1252b1a0a8a3SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
12536a4f0f73SDmitry Karpeev   PetscValidPointer(iis,4);
1254b1a0a8a3SJed Brown   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
1255b1a0a8a3SJed Brown 
1256b1a0a8a3SJed Brown   /* Get prefix, row distribution, and block size */
1257b1a0a8a3SJed Brown   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1258b1a0a8a3SJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1259b1a0a8a3SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1260b1a0a8a3SJed 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);
1261b1a0a8a3SJed Brown 
1262b1a0a8a3SJed Brown   /* Get diagonal block from matrix if possible */
1263ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1264b1a0a8a3SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr);
1265b1a0a8a3SJed Brown   if (f) {
126611bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
1267b1a0a8a3SJed Brown   } else if (size == 1) {
126811bd1e4dSLisandro Dalcin     Ad = A;
1269b1a0a8a3SJed Brown   }
1270b1a0a8a3SJed Brown   if (Ad) {
1271251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1272251f4c67SDmitry Karpeev     if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
1273b1a0a8a3SJed Brown   }
1274b1a0a8a3SJed Brown   if (Ad && n > 1) {
1275b1a0a8a3SJed Brown     PetscBool  match,done;
1276b1a0a8a3SJed Brown     /* Try to setup a good matrix partitioning if available */
1277b1a0a8a3SJed Brown     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
1278b1a0a8a3SJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
1279b1a0a8a3SJed Brown     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1280251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
1281b1a0a8a3SJed Brown     if (!match) {
1282251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
1283b1a0a8a3SJed Brown     }
1284b1a0a8a3SJed Brown     if (!match) { /* assume a "good" partitioner is available */
12851a83f524SJed Brown       PetscInt       na;
12861a83f524SJed Brown       const PetscInt *ia,*ja;
1287b1a0a8a3SJed Brown       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1288b1a0a8a3SJed Brown       if (done) {
1289b1a0a8a3SJed Brown         /* Build adjacency matrix by hand. Unfortunately a call to
1290b1a0a8a3SJed Brown            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1291b1a0a8a3SJed Brown            remove the block-aij structure and we cannot expect
1292b1a0a8a3SJed Brown            MatPartitioning to split vertices as we need */
12931a83f524SJed Brown         PetscInt       i,j,len,nnz,cnt,*iia=0,*jja=0;
12941a83f524SJed Brown         const PetscInt *row;
1295b1a0a8a3SJed Brown         nnz = 0;
1296b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* count number of nonzeros */
1297b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1298b1a0a8a3SJed Brown           row = ja + ia[i];
1299b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1300b1a0a8a3SJed Brown             if (row[j] == i) { /* don't count diagonal */
1301b1a0a8a3SJed Brown               len--; break;
1302b1a0a8a3SJed Brown             }
1303b1a0a8a3SJed Brown           }
1304b1a0a8a3SJed Brown           nnz += len;
1305b1a0a8a3SJed Brown         }
1306b1a0a8a3SJed Brown         ierr   = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr);
1307b1a0a8a3SJed Brown         ierr   = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr);
1308b1a0a8a3SJed Brown         nnz    = 0;
1309b1a0a8a3SJed Brown         iia[0] = 0;
1310b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* fill adjacency */
1311b1a0a8a3SJed Brown           cnt = 0;
1312b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1313b1a0a8a3SJed Brown           row = ja + ia[i];
1314b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
13152fa5cd67SKarl Rupp             if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */
1316b1a0a8a3SJed Brown           }
1317b1a0a8a3SJed Brown           nnz += cnt;
1318b1a0a8a3SJed Brown           iia[i+1] = nnz;
1319b1a0a8a3SJed Brown         }
1320b1a0a8a3SJed Brown         /* Partitioning of the adjacency matrix */
13210298fd71SBarry Smith         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr);
1322b1a0a8a3SJed Brown         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
1323b1a0a8a3SJed Brown         ierr      = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
1324b1a0a8a3SJed Brown         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
1325b1a0a8a3SJed Brown         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
13266bf464f9SBarry Smith         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
1327b1a0a8a3SJed Brown         foundpart = PETSC_TRUE;
1328b1a0a8a3SJed Brown       }
1329b1a0a8a3SJed Brown       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1330b1a0a8a3SJed Brown     }
13316bf464f9SBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1332b1a0a8a3SJed Brown   }
1333b1a0a8a3SJed Brown   ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr);
1334b1a0a8a3SJed Brown   if (!foundpart) {
1335b1a0a8a3SJed Brown 
1336b1a0a8a3SJed Brown     /* Partitioning by contiguous chunks of rows */
1337b1a0a8a3SJed Brown 
1338b1a0a8a3SJed Brown     PetscInt mbs   = (rend-rstart)/bs;
1339b1a0a8a3SJed Brown     PetscInt start = rstart;
1340b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1341b1a0a8a3SJed Brown       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
1342b1a0a8a3SJed Brown       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
1343b1a0a8a3SJed Brown       start += count;
1344b1a0a8a3SJed Brown     }
1345b1a0a8a3SJed Brown 
1346b1a0a8a3SJed Brown   } else {
1347b1a0a8a3SJed Brown 
1348b1a0a8a3SJed Brown     /* Partitioning by adjacency of diagonal block  */
1349b1a0a8a3SJed Brown 
1350b1a0a8a3SJed Brown     const PetscInt *numbering;
1351b1a0a8a3SJed Brown     PetscInt       *count,nidx,*indices,*newidx,start=0;
1352b1a0a8a3SJed Brown     /* Get node count in each partition */
1353b1a0a8a3SJed Brown     ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr);
1354b1a0a8a3SJed Brown     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
1355b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1356b1a0a8a3SJed Brown       for (i=0; i<n; i++) count[i] *= bs;
1357b1a0a8a3SJed Brown     }
1358b1a0a8a3SJed Brown     /* Build indices from node numbering */
1359b1a0a8a3SJed Brown     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1360b1a0a8a3SJed Brown     ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr);
1361b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1362b1a0a8a3SJed Brown     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
1363b1a0a8a3SJed Brown     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
1364b1a0a8a3SJed Brown     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
1365b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1366b1a0a8a3SJed Brown       ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr);
13672fa5cd67SKarl Rupp       for (i=0; i<nidx; i++) {
13682fa5cd67SKarl Rupp         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
13692fa5cd67SKarl Rupp       }
1370b1a0a8a3SJed Brown       ierr    = PetscFree(indices);CHKERRQ(ierr);
1371b1a0a8a3SJed Brown       nidx   *= bs;
1372b1a0a8a3SJed Brown       indices = newidx;
1373b1a0a8a3SJed Brown     }
1374b1a0a8a3SJed Brown     /* Shift to get global indices */
1375b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] += rstart;
1376b1a0a8a3SJed Brown 
1377b1a0a8a3SJed Brown     /* Build the index sets for each block */
1378b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1379b1a0a8a3SJed Brown       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
1380b1a0a8a3SJed Brown       ierr   = ISSort(is[i]);CHKERRQ(ierr);
1381b1a0a8a3SJed Brown       start += count[i];
1382b1a0a8a3SJed Brown     }
1383b1a0a8a3SJed Brown 
1384b1a0a8a3SJed Brown     ierr = PetscFree(count);
1385b1a0a8a3SJed Brown     ierr = PetscFree(indices);
1386fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1387fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
1388b1a0a8a3SJed Brown   }
13896a4f0f73SDmitry Karpeev   *iis = is;
13906a4f0f73SDmitry Karpeev   if (!ois) PetscFunctionReturn(0);
13916a4f0f73SDmitry Karpeev   /*
13926a4f0f73SDmitry Karpeev    Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap
13936a4f0f73SDmitry Karpeev    has been requested, copy the inner subdomains over so they can be modified.
13946a4f0f73SDmitry Karpeev    */
13956a4f0f73SDmitry Karpeev   ierr = PetscMalloc(n*sizeof(IS),ois);CHKERRQ(ierr);
13966a4f0f73SDmitry Karpeev   for (i=0; i<n; ++i) {
13976a4f0f73SDmitry Karpeev     if (overlap > 0) { /* With positive overlap, (*iis)[i] will be modified */
13986a4f0f73SDmitry Karpeev       ierr = ISDuplicate((*iis)[i],(*ois)+i);CHKERRQ(ierr);
13996a4f0f73SDmitry Karpeev       ierr = ISCopy((*iis)[i],(*ois)[i]);CHKERRQ(ierr);
14006a4f0f73SDmitry Karpeev     } else {
14016a4f0f73SDmitry Karpeev       ierr      = PetscObjectReference((PetscObject)(*iis)[i]);CHKERRQ(ierr);
14026a4f0f73SDmitry Karpeev       (*ois)[i] = (*iis)[i];
14036a4f0f73SDmitry Karpeev     }
14046a4f0f73SDmitry Karpeev   }
14056a4f0f73SDmitry Karpeev   if (overlap > 0) {
14066a4f0f73SDmitry Karpeev     /* Extend the "overlapping" regions by a number of steps */
14076a4f0f73SDmitry Karpeev     ierr = MatIncreaseOverlap(A,n,*ois,overlap);CHKERRQ(ierr);
14086a4f0f73SDmitry Karpeev   }
1409b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1410b1a0a8a3SJed Brown }
1411b1a0a8a3SJed Brown 
1412b1a0a8a3SJed Brown #undef __FUNCT__
1413f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains"
1414b1a0a8a3SJed Brown /*@C
1415f746d493SDmitry Karpeev    PCGASMDestroySubdomains - Destroys the index sets created with
141606b43e7eSDmitry Karpeev    PCGASMCreateLocalSubdomains() or PCGASMCreateSubdomains2D. Should be
141706b43e7eSDmitry Karpeev    called after setting subdomains with PCGASMSetSubdomains().
1418b1a0a8a3SJed Brown 
1419b1a0a8a3SJed Brown    Collective
1420b1a0a8a3SJed Brown 
1421b1a0a8a3SJed Brown    Input Parameters:
1422b1a0a8a3SJed Brown +  n   - the number of index sets
14236a4f0f73SDmitry Karpeev .  iis - the array of inner subdomains,
14240298fd71SBarry Smith -  ois - the array of outer subdomains, can be NULL
1425b1a0a8a3SJed Brown 
14266a4f0f73SDmitry Karpeev    Level: intermediate
14276a4f0f73SDmitry Karpeev 
14286a4f0f73SDmitry Karpeev    Notes: this is merely a convenience subroutine that walks each list,
14296a4f0f73SDmitry Karpeev    destroys each IS on the list, and then frees the list.
1430b1a0a8a3SJed Brown 
1431f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1432b1a0a8a3SJed Brown 
143306b43e7eSDmitry Karpeev .seealso: PCGASMCreateLocalSubdomains(), PCGASMSetSubdomains()
1434b1a0a8a3SJed Brown @*/
14356a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMDestroySubdomains(PetscInt n, IS iis[], IS ois[])
1436b1a0a8a3SJed Brown {
1437b1a0a8a3SJed Brown   PetscInt       i;
1438b1a0a8a3SJed Brown   PetscErrorCode ierr;
14395fd66863SKarl Rupp 
1440b1a0a8a3SJed Brown   PetscFunctionBegin;
1441a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1442a35b7b57SDmitry Karpeev   if (iis) {
14436a4f0f73SDmitry Karpeev     PetscValidPointer(iis,2);
1444a35b7b57SDmitry Karpeev     for (i=0; i<n; i++) {
1445a35b7b57SDmitry Karpeev       ierr = ISDestroy(&iis[i]);CHKERRQ(ierr);
1446a35b7b57SDmitry Karpeev     }
14476a4f0f73SDmitry Karpeev     ierr = PetscFree(iis);CHKERRQ(ierr);
1448a35b7b57SDmitry Karpeev   }
14496a4f0f73SDmitry Karpeev   if (ois) {
1450a35b7b57SDmitry Karpeev     for (i=0; i<n; i++) {
1451a35b7b57SDmitry Karpeev       ierr = ISDestroy(&ois[i]);CHKERRQ(ierr);
1452a35b7b57SDmitry Karpeev     }
14536a4f0f73SDmitry Karpeev     ierr = PetscFree(ois);CHKERRQ(ierr);
1454b1a0a8a3SJed Brown   }
1455b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1456b1a0a8a3SJed Brown }
1457b1a0a8a3SJed Brown 
1458af538404SDmitry Karpeev 
1459af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \
1460af538404SDmitry Karpeev   {                                                                                                       \
1461af538404SDmitry Karpeev     PetscInt first_row = first/M, last_row = last/M+1;                                                     \
1462af538404SDmitry Karpeev     /*                                                                                                    \
1463af538404SDmitry Karpeev      Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners       \
1464af538404SDmitry Karpeev      of the bounding box of the intersection of the subdomain with the local ownership range (local       \
1465af538404SDmitry Karpeev      subdomain).                                                                                          \
1466af538404SDmitry Karpeev      Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows       \
1467af538404SDmitry Karpeev      of the intersection.                                                                                 \
1468af538404SDmitry Karpeev     */                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                             \
1469af538404SDmitry Karpeev     /* ylow_loc is the grid row containing the first element of the local sumbdomain */                   \
1470eec7e3faSDmitry Karpeev     *ylow_loc = PetscMax(first_row,ylow);                                                                    \
1471af538404SDmitry Karpeev     /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1472eec7e3faSDmitry Karpeev     *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft;                                                                            \
1473af538404SDmitry Karpeev     /* yhigh_loc is the grid row above the last local subdomain element */                                                                    \
1474eec7e3faSDmitry Karpeev     *yhigh_loc = PetscMin(last_row,yhigh);                                                                                                     \
1475af538404SDmitry Karpeev     /* xright is the offset of the end of the  local subdomain within its grid row (might actually be outside the local subdomain) */         \
1476eec7e3faSDmitry Karpeev     *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright;                                                                          \
1477af538404SDmitry Karpeev     /* Now compute the size of the local subdomain n. */ \
1478c3518bceSDmitry Karpeev     *n = 0;                                               \
1479eec7e3faSDmitry Karpeev     if (*ylow_loc < *yhigh_loc) {                           \
1480af538404SDmitry Karpeev       PetscInt width = xright-xleft;                     \
1481eec7e3faSDmitry Karpeev       *n += width*(*yhigh_loc-*ylow_loc-1);                 \
1482eec7e3faSDmitry Karpeev       *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \
1483eec7e3faSDmitry Karpeev       *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \
1484af538404SDmitry Karpeev     } \
1485af538404SDmitry Karpeev   }
1486af538404SDmitry Karpeev 
1487c3518bceSDmitry Karpeev 
1488eec7e3faSDmitry Karpeev 
1489eec7e3faSDmitry Karpeev #undef __FUNCT__
1490f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D"
1491b1a0a8a3SJed Brown /*@
1492f746d493SDmitry Karpeev    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1493b1a0a8a3SJed Brown    preconditioner for a two-dimensional problem on a regular grid.
1494b1a0a8a3SJed Brown 
1495af538404SDmitry Karpeev    Collective
1496b1a0a8a3SJed Brown 
1497b1a0a8a3SJed Brown    Input Parameters:
149806b43e7eSDmitry Karpeev +  M, N               - the global number of grid points in the x and y directions
1499af538404SDmitry Karpeev .  Mdomains, Ndomains - the global number of subdomains in the x and y directions
1500b1a0a8a3SJed Brown .  dof                - degrees of freedom per node
1501b1a0a8a3SJed Brown -  overlap            - overlap in mesh lines
1502b1a0a8a3SJed Brown 
1503b1a0a8a3SJed Brown    Output Parameters:
1504af538404SDmitry Karpeev +  Nsub - the number of local subdomains created
15056a4f0f73SDmitry Karpeev .  iis  - array of index sets defining inner (nonoverlapping) subdomains
15066a4f0f73SDmitry Karpeev -  ois  - array of index sets defining outer (overlapping, if overlap > 0) subdomains
1507b1a0a8a3SJed Brown 
1508b1a0a8a3SJed Brown 
1509b1a0a8a3SJed Brown    Level: advanced
1510b1a0a8a3SJed Brown 
1511f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid
1512b1a0a8a3SJed Brown 
15136a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1514f746d493SDmitry Karpeev           PCGASMSetOverlap()
1515b1a0a8a3SJed Brown @*/
15166a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMCreateSubdomains2D(PC pc, PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap, PetscInt *nsub,IS **iis,IS **ois)
1517b1a0a8a3SJed Brown {
1518b1a0a8a3SJed Brown   PetscErrorCode ierr;
15192bbb417fSDmitry Karpeev   PetscMPIInt    size, rank;
15202bbb417fSDmitry Karpeev   PetscInt       i, j;
15212bbb417fSDmitry Karpeev   PetscInt       maxheight, maxwidth;
15222bbb417fSDmitry Karpeev   PetscInt       xstart, xleft, xright, xleft_loc, xright_loc;
15232bbb417fSDmitry Karpeev   PetscInt       ystart, ylow,  yhigh,  ylow_loc,  yhigh_loc;
1524eec7e3faSDmitry Karpeev   PetscInt       x[2][2], y[2][2], n[2];
15252bbb417fSDmitry Karpeev   PetscInt       first, last;
15262bbb417fSDmitry Karpeev   PetscInt       nidx, *idx;
15272bbb417fSDmitry Karpeev   PetscInt       ii,jj,s,q,d;
1528af538404SDmitry Karpeev   PetscInt       k,kk;
15292bbb417fSDmitry Karpeev   PetscMPIInt    color;
15302bbb417fSDmitry Karpeev   MPI_Comm       comm, subcomm;
15316a4f0f73SDmitry Karpeev   IS             **xis = 0, **is = ois, **is_local = iis;
1532d34fcf5fSBarry Smith 
1533b1a0a8a3SJed Brown   PetscFunctionBegin;
15342bbb417fSDmitry Karpeev   ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr);
15352bbb417fSDmitry Karpeev   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
15362bbb417fSDmitry Karpeev   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
15372bbb417fSDmitry Karpeev   ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr);
1538d34fcf5fSBarry 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) "
15392bbb417fSDmitry Karpeev                                       "does not respect the number of degrees of freedom per grid point %D", first, last, dof);
1540d34fcf5fSBarry Smith 
1541af538404SDmitry Karpeev   /* Determine the number of domains with nonzero intersections with the local ownership range. */
15422bbb417fSDmitry Karpeev   s      = 0;
1543af538404SDmitry Karpeev   ystart = 0;
1544af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1545af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1546af538404SDmitry 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);
1547eec7e3faSDmitry Karpeev     /* Vertical domain limits with an overlap. */
1548eec7e3faSDmitry Karpeev     ylow   = PetscMax(ystart - overlap,0);
1549eec7e3faSDmitry Karpeev     yhigh  = PetscMin(ystart + maxheight + overlap,N);
1550b1a0a8a3SJed Brown     xstart = 0;
1551af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1552af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1553af538404SDmitry 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);
1554eec7e3faSDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1555eec7e3faSDmitry Karpeev       xleft  = PetscMax(xstart - overlap,0);
1556eec7e3faSDmitry Karpeev       xright = PetscMin(xstart + maxwidth + overlap,M);
1557af538404SDmitry Karpeev       /*
1558af538404SDmitry Karpeev          Determine whether this subdomain intersects this processor's ownership range of pc->pmat.
1559af538404SDmitry Karpeev       */
1560c3518bceSDmitry Karpeev       PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
15612fa5cd67SKarl Rupp       if (nidx) ++s;
1562af538404SDmitry Karpeev       xstart += maxwidth;
1563af538404SDmitry Karpeev     } /* for (i = 0; i < Mdomains; ++i) */
1564af538404SDmitry Karpeev     ystart += maxheight;
1565af538404SDmitry Karpeev   } /* for (j = 0; j < Ndomains; ++j) */
15662fa5cd67SKarl Rupp 
1567af538404SDmitry Karpeev   /* Now we can allocate the necessary number of ISs. */
1568af538404SDmitry Karpeev   *nsub  = s;
1569af538404SDmitry Karpeev   ierr   = PetscMalloc((*nsub)*sizeof(IS*),is);CHKERRQ(ierr);
1570af538404SDmitry Karpeev   ierr   = PetscMalloc((*nsub)*sizeof(IS*),is_local);CHKERRQ(ierr);
1571af538404SDmitry Karpeev   s      = 0;
1572af538404SDmitry Karpeev   ystart = 0;
1573af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1574af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1575af538404SDmitry 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);
1576af538404SDmitry Karpeev     /* Vertical domain limits with an overlap. */
1577eec7e3faSDmitry Karpeev     y[0][0] = PetscMax(ystart - overlap,0);
1578eec7e3faSDmitry Karpeev     y[0][1] = PetscMin(ystart + maxheight + overlap,N);
1579eec7e3faSDmitry Karpeev     /* Vertical domain limits without an overlap. */
1580eec7e3faSDmitry Karpeev     y[1][0] = ystart;
1581eec7e3faSDmitry Karpeev     y[1][1] = PetscMin(ystart + maxheight,N);
1582eec7e3faSDmitry Karpeev     xstart  = 0;
1583af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1584af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1585af538404SDmitry 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);
1586af538404SDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1587eec7e3faSDmitry Karpeev       x[0][0] = PetscMax(xstart - overlap,0);
1588eec7e3faSDmitry Karpeev       x[0][1] = PetscMin(xstart + maxwidth + overlap,M);
1589eec7e3faSDmitry Karpeev       /* Horizontal domain limits without an overlap. */
1590eec7e3faSDmitry Karpeev       x[1][0] = xstart;
1591eec7e3faSDmitry Karpeev       x[1][1] = PetscMin(xstart+maxwidth,M);
15922bbb417fSDmitry Karpeev       /*
15932bbb417fSDmitry Karpeev          Determine whether this domain intersects this processor's ownership range of pc->pmat.
15942bbb417fSDmitry Karpeev          Do this twice: first for the domains with overlaps, and once without.
15952bbb417fSDmitry Karpeev          During the first pass create the subcommunicators, and use them on the second pass as well.
15962bbb417fSDmitry Karpeev       */
15972bbb417fSDmitry Karpeev       for (q = 0; q < 2; ++q) {
15982bbb417fSDmitry Karpeev         /*
15992bbb417fSDmitry Karpeev           domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
16002bbb417fSDmitry Karpeev           according to whether the domain with an overlap or without is considered.
16012bbb417fSDmitry Karpeev         */
16022bbb417fSDmitry Karpeev         xleft = x[q][0]; xright = x[q][1];
16032bbb417fSDmitry Karpeev         ylow  = y[q][0]; yhigh  = y[q][1];
1604c3518bceSDmitry Karpeev         PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1605af538404SDmitry Karpeev         nidx *= dof;
1606eec7e3faSDmitry Karpeev         n[q]  = nidx;
16072bbb417fSDmitry Karpeev         /*
1608eec7e3faSDmitry Karpeev          Based on the counted number of indices in the local domain *with an overlap*,
1609af538404SDmitry Karpeev          construct a subcommunicator of all the processors supporting this domain.
1610eec7e3faSDmitry Karpeev          Observe that a domain with an overlap might have nontrivial local support,
1611eec7e3faSDmitry Karpeev          while the domain without an overlap might not.  Hence, the decision to participate
1612eec7e3faSDmitry Karpeev          in the subcommunicator must be based on the domain with an overlap.
16132bbb417fSDmitry Karpeev          */
16142bbb417fSDmitry Karpeev         if (q == 0) {
16152fa5cd67SKarl Rupp           if (nidx) color = 1;
16162fa5cd67SKarl Rupp           else color = MPI_UNDEFINED;
16172fa5cd67SKarl Rupp 
16182bbb417fSDmitry Karpeev           ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr);
16192bbb417fSDmitry Karpeev         }
1620af538404SDmitry Karpeev         /*
1621eec7e3faSDmitry Karpeev          Proceed only if the number of local indices *with an overlap* is nonzero.
1622af538404SDmitry Karpeev          */
1623eec7e3faSDmitry Karpeev         if (n[0]) {
16242fa5cd67SKarl Rupp           if (q == 0) xis = is;
1625af538404SDmitry Karpeev           if (q == 1) {
1626af538404SDmitry Karpeev             /*
1627eec7e3faSDmitry Karpeev              The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1628af538404SDmitry Karpeev              Moreover, if the overlap is zero, the two ISs are identical.
1629af538404SDmitry Karpeev              */
1630b1a0a8a3SJed Brown             if (overlap == 0) {
1631eec7e3faSDmitry Karpeev               (*is_local)[s] = (*is)[s];
16322bbb417fSDmitry Karpeev               ierr           = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr);
16332bbb417fSDmitry Karpeev               continue;
1634d34fcf5fSBarry Smith             } else {
16356a4f0f73SDmitry Karpeev               xis     = is_local;
1636eec7e3faSDmitry Karpeev               subcomm = ((PetscObject)(*is)[s])->comm;
16372bbb417fSDmitry Karpeev             }
1638af538404SDmitry Karpeev           } /* if (q == 1) */
16390298fd71SBarry Smith           idx  = NULL;
16402bbb417fSDmitry Karpeev           ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1641eec7e3faSDmitry Karpeev           if (nidx) {
16422bbb417fSDmitry Karpeev             k = 0;
16432bbb417fSDmitry Karpeev             for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1644af538404SDmitry Karpeev               PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft;
1645af538404SDmitry Karpeev               PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright;
1646af538404SDmitry Karpeev               kk = dof*(M*jj + x0);
16472bbb417fSDmitry Karpeev               for (ii=x0; ii<x1; ++ii) {
16482bbb417fSDmitry Karpeev                 for (d = 0; d < dof; ++d) {
16492bbb417fSDmitry Karpeev                   idx[k++] = kk++;
1650b1a0a8a3SJed Brown                 }
1651b1a0a8a3SJed Brown               }
1652b1a0a8a3SJed Brown             }
1653eec7e3faSDmitry Karpeev           }
16546a4f0f73SDmitry Karpeev           ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr);
1655eec7e3faSDmitry Karpeev         }/* if (n[0]) */
16562bbb417fSDmitry Karpeev       }/* for (q = 0; q < 2; ++q) */
16572fa5cd67SKarl Rupp       if (n[0]) ++s;
1658af538404SDmitry Karpeev       xstart += maxwidth;
1659af538404SDmitry Karpeev     } /* for (i = 0; i < Mdomains; ++i) */
16602bbb417fSDmitry Karpeev     ystart += maxheight;
1661af538404SDmitry Karpeev   } /* for (j = 0; j < Ndomains; ++j) */
1662b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1663b1a0a8a3SJed Brown }
1664b1a0a8a3SJed Brown 
1665b1a0a8a3SJed Brown #undef __FUNCT__
166606b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains"
1667b1a0a8a3SJed Brown /*@C
166806b43e7eSDmitry Karpeev     PCGASMGetSubdomains - Gets the subdomains supported on this processor
166906b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
1670b1a0a8a3SJed Brown 
1671b1a0a8a3SJed Brown     Not Collective
1672b1a0a8a3SJed Brown 
1673b1a0a8a3SJed Brown     Input Parameter:
1674b1a0a8a3SJed Brown .   pc - the preconditioner context
1675b1a0a8a3SJed Brown 
1676b1a0a8a3SJed Brown     Output Parameters:
1677b1a0a8a3SJed Brown +   n   - the number of subdomains for this processor (default value = 1)
16780298fd71SBarry Smith .   iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL)
16790298fd71SBarry Smith -   ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL)
1680b1a0a8a3SJed Brown 
1681b1a0a8a3SJed Brown 
1682b1a0a8a3SJed Brown     Notes:
16836a4f0f73SDmitry Karpeev     The user is responsible for destroying the ISs and freeing the returned arrays.
1684b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
1685b1a0a8a3SJed Brown 
1686b1a0a8a3SJed Brown     Level: advanced
1687b1a0a8a3SJed Brown 
168806b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz
1689b1a0a8a3SJed Brown 
16906a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
169106b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubmatrices()
1692b1a0a8a3SJed Brown @*/
16936a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[])
1694b1a0a8a3SJed Brown {
1695f746d493SDmitry Karpeev   PC_GASM        *osm;
1696b1a0a8a3SJed Brown   PetscErrorCode ierr;
1697b1a0a8a3SJed Brown   PetscBool      match;
16986a4f0f73SDmitry Karpeev   PetscInt       i;
16995fd66863SKarl Rupp 
1700b1a0a8a3SJed Brown   PetscFunctionBegin;
1701b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1702251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1703ce94432eSBarry Smith   if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1704f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1705ab3e923fSDmitry Karpeev   if (n) *n = osm->n;
17066a4f0f73SDmitry Karpeev   if (iis) {
17076a4f0f73SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(IS), iis);CHKERRQ(ierr);
17086a4f0f73SDmitry Karpeev   }
17096a4f0f73SDmitry Karpeev   if (ois) {
17106a4f0f73SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(IS), ois);CHKERRQ(ierr);
17116a4f0f73SDmitry Karpeev   }
17126a4f0f73SDmitry Karpeev   if (iis || ois) {
17136a4f0f73SDmitry Karpeev     for (i = 0; i < osm->n; ++i) {
17146a4f0f73SDmitry Karpeev       if (iis) (*iis)[i] = osm->iis[i];
17156a4f0f73SDmitry Karpeev       if (ois) (*ois)[i] = osm->ois[i];
17166a4f0f73SDmitry Karpeev     }
1717b1a0a8a3SJed Brown   }
1718b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1719b1a0a8a3SJed Brown }
1720b1a0a8a3SJed Brown 
1721b1a0a8a3SJed Brown #undef __FUNCT__
172206b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices"
1723b1a0a8a3SJed Brown /*@C
172406b43e7eSDmitry Karpeev     PCGASMGetSubmatrices - Gets the local submatrices (for this processor
1725b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1726b1a0a8a3SJed Brown 
1727b1a0a8a3SJed Brown     Not Collective
1728b1a0a8a3SJed Brown 
1729b1a0a8a3SJed Brown     Input Parameter:
1730b1a0a8a3SJed Brown .   pc - the preconditioner context
1731b1a0a8a3SJed Brown 
1732b1a0a8a3SJed Brown     Output Parameters:
1733b1a0a8a3SJed Brown +   n   - the number of matrices for this processor (default value = 1)
1734b1a0a8a3SJed Brown -   mat - the matrices
1735b1a0a8a3SJed Brown 
173606b43e7eSDmitry Karpeev     Notes: matrices returned by this routine have the same communicators as the index sets (IS)
173706b43e7eSDmitry Karpeev            used to define subdomains in PCGASMSetSubdomains(), or PETSC_COMM_SELF, if the
17386a4f0f73SDmitry Karpeev            subdomains were defined using PCGASMSetTotalSubdomains().
1739b1a0a8a3SJed Brown     Level: advanced
1740b1a0a8a3SJed Brown 
1741f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi
1742b1a0a8a3SJed Brown 
17436a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomain(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
174406b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains()
1745b1a0a8a3SJed Brown @*/
174606b43e7eSDmitry Karpeev PetscErrorCode  PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1747b1a0a8a3SJed Brown {
1748f746d493SDmitry Karpeev   PC_GASM        *osm;
1749b1a0a8a3SJed Brown   PetscErrorCode ierr;
1750b1a0a8a3SJed Brown   PetscBool      match;
1751b1a0a8a3SJed Brown 
1752b1a0a8a3SJed Brown   PetscFunctionBegin;
1753b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1754b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1755b1a0a8a3SJed Brown   if (mat) PetscValidPointer(mat,3);
1756ce94432eSBarry Smith   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1757251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1758ce94432eSBarry Smith   if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1759f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1760ab3e923fSDmitry Karpeev   if (n) *n = osm->n;
1761b1a0a8a3SJed Brown   if (mat) *mat = osm->pmat;
1762b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1763b1a0a8a3SJed Brown }
1764*d709ea83SDmitry Karpeev 
1765*d709ea83SDmitry Karpeev #undef __FUNCT__
1766*d709ea83SDmitry Karpeev #define __FUNCT__ "PCGASMSetDMSubdomains"
1767*d709ea83SDmitry Karpeev /*@
1768*d709ea83SDmitry Karpeev     PCGASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1769*d709ea83SDmitry Karpeev     Logically Collective
1770*d709ea83SDmitry Karpeev 
1771*d709ea83SDmitry Karpeev     Input Parameter:
1772*d709ea83SDmitry Karpeev +   pc  - the preconditioner
1773*d709ea83SDmitry Karpeev -   flg - boolean indicating whether to use subdomains defined by the DM
1774*d709ea83SDmitry Karpeev 
1775*d709ea83SDmitry Karpeev     Options Database Key:
1776*d709ea83SDmitry Karpeev .   -pc_gasm_dm_subdomains
1777*d709ea83SDmitry Karpeev 
1778*d709ea83SDmitry Karpeev     Level: intermediate
1779*d709ea83SDmitry Karpeev 
1780*d709ea83SDmitry Karpeev     Notes:
1781*d709ea83SDmitry Karpeev     PCGASMSetTotalSubdomains() and PCGASMSetOverlap() take precedence over PCGASMSetDMSubdomains(),
1782*d709ea83SDmitry Karpeev     so setting either of the first two effectively turns the latter off.
1783*d709ea83SDmitry Karpeev 
1784*d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1785*d709ea83SDmitry Karpeev 
1786*d709ea83SDmitry Karpeev .seealso: PCGASMGetDMSubdomains(), PCGASMSetTotalSubdomains(), PCGASMSetOverlap()
1787*d709ea83SDmitry Karpeev           PCGASMCreateSubdomains2D()
1788*d709ea83SDmitry Karpeev @*/
1789*d709ea83SDmitry Karpeev PetscErrorCode  PCGASMSetDMSubdomains(PC pc,PetscBool flg)
1790*d709ea83SDmitry Karpeev {
1791*d709ea83SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
1792*d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1793*d709ea83SDmitry Karpeev   PetscBool      match;
1794*d709ea83SDmitry Karpeev 
1795*d709ea83SDmitry Karpeev   PetscFunctionBegin;
1796*d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1797*d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc,flg,2);
1798*d709ea83SDmitry Karpeev   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1799*d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1800*d709ea83SDmitry Karpeev   if (match) {
1801*d709ea83SDmitry Karpeev     osm->dm_subdomains = flg;
1802*d709ea83SDmitry Karpeev   }
1803*d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1804*d709ea83SDmitry Karpeev }
1805*d709ea83SDmitry Karpeev 
1806*d709ea83SDmitry Karpeev #undef __FUNCT__
1807*d709ea83SDmitry Karpeev #define __FUNCT__ "PCGASMGetDMSubdomains"
1808*d709ea83SDmitry Karpeev /*@
1809*d709ea83SDmitry Karpeev     PCGASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1810*d709ea83SDmitry Karpeev     Not Collective
1811*d709ea83SDmitry Karpeev 
1812*d709ea83SDmitry Karpeev     Input Parameter:
1813*d709ea83SDmitry Karpeev .   pc  - the preconditioner
1814*d709ea83SDmitry Karpeev 
1815*d709ea83SDmitry Karpeev     Output Parameter:
1816*d709ea83SDmitry Karpeev .   flg - boolean indicating whether to use subdomains defined by the DM
1817*d709ea83SDmitry Karpeev 
1818*d709ea83SDmitry Karpeev     Level: intermediate
1819*d709ea83SDmitry Karpeev 
1820*d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz
1821*d709ea83SDmitry Karpeev 
1822*d709ea83SDmitry Karpeev .seealso: PCGASMSetDMSubdomains(), PCGASMSetTotalSubdomains(), PCGASMSetOverlap()
1823*d709ea83SDmitry Karpeev           PCGASMCreateSubdomains2D()
1824*d709ea83SDmitry Karpeev @*/
1825*d709ea83SDmitry Karpeev PetscErrorCode  PCGASMGetDMSubdomains(PC pc,PetscBool* flg)
1826*d709ea83SDmitry Karpeev {
1827*d709ea83SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
1828*d709ea83SDmitry Karpeev   PetscErrorCode ierr;
1829*d709ea83SDmitry Karpeev   PetscBool      match;
1830*d709ea83SDmitry Karpeev 
1831*d709ea83SDmitry Karpeev   PetscFunctionBegin;
1832*d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1833*d709ea83SDmitry Karpeev   PetscValidPointer(flg,2);
1834*d709ea83SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1835*d709ea83SDmitry Karpeev   if (match) {
1836*d709ea83SDmitry Karpeev     if (flg) *flg = osm->dm_subdomains;
1837*d709ea83SDmitry Karpeev   }
1838*d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1839*d709ea83SDmitry Karpeev }
1840