xref: /petsc/src/ksp/pc/impls/gasm/gasm.c (revision 11aeaf0a669438c21ebbbd10be78fc394d59507d)
1b1a0a8a3SJed Brown /*
2f746d493SDmitry Karpeev   This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation.
36a4f0f73SDmitry Karpeev   In this version each processor may intersect multiple subdomains and any subdomain may
46a4f0f73SDmitry Karpeev   intersect multiple processors.  Intersections of subdomains with processors are called *local
56a4f0f73SDmitry Karpeev   subdomains*.
6b1a0a8a3SJed Brown 
76a4f0f73SDmitry Karpeev        N    - total number of local subdomains on all processors  (set in PCGASMSetTotalSubdomains() or calculated in PCSetUp_GASM())
86a4f0f73SDmitry Karpeev        n    - actual number of local subdomains on this processor (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains())
96a4f0f73SDmitry Karpeev        nmax - maximum number of local subdomains per processor    (calculated in PCGASMSetTotalSubdomains() or in PCSetUp_GASM())
10b1a0a8a3SJed Brown */
11b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>     /*I "petscpc.h" I*/
12b1a0a8a3SJed Brown 
13b1a0a8a3SJed Brown typedef struct {
14f746d493SDmitry Karpeev   PetscInt   N,n,nmax;
15b1a0a8a3SJed Brown   PetscInt   overlap;             /* overlap requested by user */
16b1a0a8a3SJed Brown   KSP        *ksp;                /* linear solvers for each block */
17f746d493SDmitry Karpeev   Vec        gx,gy;               /* Merged work vectors */
18f746d493SDmitry Karpeev   Vec        *x,*y;               /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */
196a4f0f73SDmitry Karpeev   VecScatter gorestriction;       /* merged restriction to disjoint union of outer subdomains */
206a4f0f73SDmitry Karpeev   VecScatter girestriction;       /* merged restriction to disjoint union of inner subdomains */
216a4f0f73SDmitry Karpeev   IS         *ois;                /* index sets that define the outer (conceptually, overlapping) subdomains */
226a4f0f73SDmitry Karpeev   IS         *iis;                /* index sets that define the inner (conceptually, nonoverlapping) subdomains */
23f746d493SDmitry Karpeev   Mat        *pmat;               /* subdomain block matrices */
24f746d493SDmitry Karpeev   PCGASMType type;                /* use reduced interpolation, restriction or both */
256a4f0f73SDmitry Karpeev   PetscBool  create_local;           /* whether the autocreated subdomains are local or not. */
26b1a0a8a3SJed Brown   PetscBool  type_set;               /* if user set this value (so won't change it for symmetric problems) */
276a4f0f73SDmitry Karpeev   PetscBool  same_subdomain_solvers; /* flag indicating whether all local solvers are same */
28b1a0a8a3SJed Brown   PetscBool  sort_indices;           /* flag to sort subdomain indices */
29f746d493SDmitry Karpeev } PC_GASM;
30b1a0a8a3SJed Brown 
31b1a0a8a3SJed Brown #undef __FUNCT__
3206b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSubdomainView_Private"
3306b43e7eSDmitry Karpeev static PetscErrorCode  PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer)
34af538404SDmitry Karpeev {
35af538404SDmitry Karpeev   PC_GASM        *osm  = (PC_GASM*)pc->data;
3606b43e7eSDmitry Karpeev   PetscInt       j,nidx;
37af538404SDmitry Karpeev   const PetscInt *idx;
3806b43e7eSDmitry Karpeev   PetscViewer    sviewer;
39af538404SDmitry Karpeev   char           *cidx;
40af538404SDmitry Karpeev   PetscErrorCode ierr;
41af538404SDmitry Karpeev 
42af538404SDmitry Karpeev   PetscFunctionBegin;
4306b43e7eSDmitry Karpeev   if (i < 0 || i > osm->n) SETERRQ2(((PetscObject)viewer)->comm, PETSC_ERR_ARG_WRONG, "Invalid subdomain %D: must nonnegative and less than %D", i, osm->n);
444bde246dSDmitry Karpeev   /* Inner subdomains. */
454bde246dSDmitry Karpeev   ierr = ISGetLocalSize(osm->iis[i], &nidx);CHKERRQ(ierr);
464bde246dSDmitry Karpeev   /*
474bde246dSDmitry Karpeev    No more than 15 characters per index plus a space.
484bde246dSDmitry Karpeev    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
494bde246dSDmitry Karpeev    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
504bde246dSDmitry Karpeev    For nidx == 0, the whole string 16 '\0'.
514bde246dSDmitry Karpeev    */
524bde246dSDmitry Karpeev   ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr);
534bde246dSDmitry Karpeev   ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr);
544bde246dSDmitry Karpeev   ierr = ISGetIndices(osm->iis[i], &idx);CHKERRQ(ierr);
554bde246dSDmitry Karpeev   for (j = 0; j < nidx; ++j) {
564bde246dSDmitry Karpeev     ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr);
574bde246dSDmitry Karpeev   }
584bde246dSDmitry Karpeev   ierr = ISRestoreIndices(osm->iis[i],&idx);CHKERRQ(ierr);
594bde246dSDmitry Karpeev   ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
604bde246dSDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n");CHKERRQ(ierr);
614bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
624bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
634bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
644bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
654bde246dSDmitry Karpeev   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
664bde246dSDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
674bde246dSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
684bde246dSDmitry Karpeev   ierr = PetscFree(cidx);CHKERRQ(ierr);
694bde246dSDmitry Karpeev   /* Outer subdomains. */
706a4f0f73SDmitry Karpeev   ierr = ISGetLocalSize(osm->ois[i], &nidx);CHKERRQ(ierr);
71eec7e3faSDmitry Karpeev   /*
72eec7e3faSDmitry Karpeev    No more than 15 characters per index plus a space.
73eec7e3faSDmitry Karpeev    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
74eec7e3faSDmitry Karpeev    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
75eec7e3faSDmitry Karpeev    For nidx == 0, the whole string 16 '\0'.
76eec7e3faSDmitry Karpeev    */
77eec7e3faSDmitry Karpeev   ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr);
7806b43e7eSDmitry Karpeev   ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr);
796a4f0f73SDmitry Karpeev   ierr = ISGetIndices(osm->ois[i], &idx);CHKERRQ(ierr);
80af538404SDmitry Karpeev   for (j = 0; j < nidx; ++j) {
81af538404SDmitry Karpeev     ierr = PetscViewerStringSPrintf(sviewer,"%D ", idx[j]);CHKERRQ(ierr);
82af538404SDmitry Karpeev   }
836bf464f9SBarry Smith   ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
846a4f0f73SDmitry Karpeev   ierr = ISRestoreIndices(osm->ois[i],&idx);CHKERRQ(ierr);
856a4f0f73SDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n");CHKERRQ(ierr);
8606b43e7eSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
877b23a99aSBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
88af538404SDmitry Karpeev   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
89af538404SDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
907b23a99aSBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
91af538404SDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
9206b43e7eSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
93af538404SDmitry Karpeev   ierr = PetscFree(cidx);CHKERRQ(ierr);
944bde246dSDmitry Karpeev 
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;
1074bde246dSDmitry Karpeev   PetscViewer    viewer, sviewer = PETSC_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;
1134bde246dSDmitry Karpeev   ierr = PetscObjectsGetGlobalNumbering(((PetscObject)pc)->comm, 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); };
1184bde246dSDmitry Karpeev   ierr = PetscViewerASCIIOpen(((PetscObject)pc)->comm,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     }
1364bde246dSDmitry Karpeev     ierr = MPI_Barrier(((PetscObject)pc)->comm);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";
160*11aeaf0aSBarry Smith 
161b1a0a8a3SJed Brown   PetscFunctionBegin;
162af538404SDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm, &size);CHKERRQ(ierr);
163af538404SDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm, &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;
1666a4f0f73SDmitry Karpeev   ierr = PetscObjectsGetGlobalNumbering(((PetscObject)pc)->comm, 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);
18106b43e7eSDmitry Karpeev   ierr = PetscOptionsGetBool(prefix,"-pc_gasm_view_subdomains",&view_subdomains,PETSC_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. */
20606b43e7eSDmitry Karpeev     ierr = MPI_Barrier(((PetscObject)viewer)->comm);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);
2276a4f0f73SDmitry Karpeev           }
2286a4f0f73SDmitry Karpeev           else {
22906b43e7eSDmitry Karpeev             ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr);
2306a4f0f73SDmitry Karpeev           }
23106b43e7eSDmitry Karpeev           ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
2327b23a99aSBarry Smith           ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
2336a4f0f73SDmitry Karpeev           ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
23406b43e7eSDmitry Karpeev           ++l;
235b1a0a8a3SJed Brown         }
23606b43e7eSDmitry Karpeev       }
23706b43e7eSDmitry Karpeev       ierr = MPI_Barrier(((PetscObject)pc)->comm);CHKERRQ(ierr);
238b1a0a8a3SJed Brown     }
239b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
240b1a0a8a3SJed Brown   }
24106b43e7eSDmitry Karpeev   ierr = PetscFree2(permutation,numbering);CHKERRQ(ierr);
242b1a0a8a3SJed Brown   PetscFunctionReturn(0);
243b1a0a8a3SJed Brown }
244b1a0a8a3SJed Brown 
245af538404SDmitry Karpeev 
246af538404SDmitry Karpeev 
247af538404SDmitry Karpeev 
248af538404SDmitry Karpeev 
249b1a0a8a3SJed Brown #undef __FUNCT__
250f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM"
251f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc)
252b1a0a8a3SJed Brown {
253f746d493SDmitry Karpeev   PC_GASM         *osm  = (PC_GASM*)pc->data;
254b1a0a8a3SJed Brown   PetscErrorCode ierr;
255b1a0a8a3SJed Brown   PetscBool      symset,flg;
25688389755SJed Brown   PetscInt       i;
257c10bc1a1SDmitry Karpeev   PetscMPIInt    rank, size;
258b1a0a8a3SJed Brown   MatReuse       scall = MAT_REUSE_MATRIX;
259b1a0a8a3SJed Brown   KSP            ksp;
260b1a0a8a3SJed Brown   PC             subpc;
261b1a0a8a3SJed Brown   const char     *prefix,*pprefix;
262f746d493SDmitry Karpeev   Vec            x,y;
2636a4f0f73SDmitry Karpeev   PetscInt       oni;       /* Number of indices in the i-th local outer subdomain.               */
2646a4f0f73SDmitry Karpeev   const PetscInt *oidxi;    /* Indices from the i-th subdomain local outer subdomain.             */
2656a4f0f73SDmitry Karpeev   PetscInt       on;        /* Number of indices in the disjoint union of local outer subdomains. */
2666a4f0f73SDmitry Karpeev   PetscInt       *oidx;     /* Indices in the disjoint union of local outer subdomains. */
2676a4f0f73SDmitry Karpeev   IS             gois;      /* Disjoint union the global indices of outer subdomains.             */
2686a4f0f73SDmitry Karpeev   IS             goid;      /* Identity IS of the size of the disjoint union of outer subdomains. */
269f746d493SDmitry Karpeev   PetscScalar    *gxarray, *gyarray;
2706a4f0f73SDmitry Karpeev   PetscInt       gofirst;   /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy --
2716a4f0f73SDmitry Karpeev                              over the disjoint union of outer subdomains. */
272a35b7b57SDmitry Karpeev   DM             *subdomain_dm = PETSC_NULL;
273b1a0a8a3SJed Brown 
274b1a0a8a3SJed Brown   PetscFunctionBegin;
275c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
276c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
277b1a0a8a3SJed Brown   if (!pc->setupcalled) {
278b1a0a8a3SJed Brown 
279b1a0a8a3SJed Brown     if (!osm->type_set) {
280b1a0a8a3SJed Brown       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
281f746d493SDmitry Karpeev       if (symset && flg) { osm->type = PC_GASM_BASIC; }
282b1a0a8a3SJed Brown     }
283b1a0a8a3SJed Brown 
28406b43e7eSDmitry Karpeev     /*
28506b43e7eSDmitry Karpeev      If subdomains have been set, then the local number of subdomains, osm->n, is NOT PETSC_DECIDE and is at least 1.
28606b43e7eSDmitry 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.
28706b43e7eSDmitry Karpeev      */
28806b43e7eSDmitry Karpeev     if (osm->n == PETSC_DECIDE) {
28969ca1f37SDmitry Karpeev       /* no subdomains given */
29069ca1f37SDmitry Karpeev       /* try pc->dm first */
291fdc48646SDmitry Karpeev       if (pc->dm) {
292fdc48646SDmitry Karpeev         char      ddm_name[1024];
293fdc48646SDmitry Karpeev         DM        ddm;
294fdc48646SDmitry Karpeev         PetscBool flg;
295a35b7b57SDmitry Karpeev         PetscInt     num_subdomains, d;
296a35b7b57SDmitry Karpeev         char         **subdomain_names;
297a35b7b57SDmitry Karpeev         IS           *inner_subdomain_is, *outer_subdomain_is;
298fdc48646SDmitry Karpeev         /* Allow the user to request a decomposition DM by name */
299fdc48646SDmitry Karpeev         ierr = PetscStrncpy(ddm_name, "", 1024);CHKERRQ(ierr);
300a35b7b57SDmitry Karpeev         ierr = PetscOptionsString("-pc_gasm_decomposition","Name of the DM defining the decomposition", "PCSetDM",ddm_name,ddm_name,1024,&flg);CHKERRQ(ierr);
301fdc48646SDmitry Karpeev         if (flg) {
30216621825SDmitry Karpeev           ierr = DMCreateDomainDecompositionDM(pc->dm, ddm_name, &ddm);CHKERRQ(ierr);
303f23aa3ddSBarry Smith           if (!ddm) SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown DM decomposition name %s", ddm_name);
304fdc48646SDmitry Karpeev           ierr = PetscInfo(pc,"Using decomposition DM defined using options database\n");CHKERRQ(ierr);
305fdc48646SDmitry Karpeev           ierr = PCSetDM(pc,ddm);CHKERRQ(ierr);
306fdc48646SDmitry Karpeev         }
307a35b7b57SDmitry Karpeev         ierr = DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);CHKERRQ(ierr);
308a35b7b57SDmitry Karpeev         if (num_subdomains) {
309a35b7b57SDmitry Karpeev           ierr = PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);CHKERRQ(ierr);
31069ca1f37SDmitry Karpeev         }
311a35b7b57SDmitry Karpeev         for (d = 0; d < num_subdomains; ++d) {
312a35b7b57SDmitry Karpeev           if (subdomain_names)    {ierr = PetscFree(subdomain_names[d]);CHKERRQ(ierr);}
313a35b7b57SDmitry Karpeev           if (inner_subdomain_is) {ierr = ISDestroy(&inner_subdomain_is[d]);CHKERRQ(ierr);}
314a35b7b57SDmitry Karpeev           if (outer_subdomain_is) {ierr = ISDestroy(&outer_subdomain_is[d]);CHKERRQ(ierr);}
315fdc48646SDmitry Karpeev         }
316a35b7b57SDmitry Karpeev         ierr = PetscFree(subdomain_names);CHKERRQ(ierr);
317a35b7b57SDmitry Karpeev         ierr = PetscFree(inner_subdomain_is);CHKERRQ(ierr);
318a35b7b57SDmitry Karpeev         ierr = PetscFree(outer_subdomain_is);CHKERRQ(ierr);
319fdc48646SDmitry Karpeev       }
32006b43e7eSDmitry Karpeev       if (osm->n == PETSC_DECIDE) { /* still no subdomains; use one per processor */
32144fe18e5SDmitry Karpeev         osm->nmax = osm->n = 1;
322b1a0a8a3SJed Brown         ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
323f746d493SDmitry Karpeev         osm->N = size;
324fdc48646SDmitry Karpeev       }
32506b43e7eSDmitry Karpeev     }
326a35b7b57SDmitry Karpeev     if (!osm->iis) {
327a35b7b57SDmitry Karpeev       /*
328a35b7b57SDmitry Karpeev        The local number of subdomains was set just above, or in PCGASMSetTotalSubdomains(), or in PCGASMSetSubdomains(),
329a35b7b57SDmitry Karpeev        but the actual subdomains have not been supplied (in PCGASMSetSubdomains()).
330a35b7b57SDmitry Karpeev        We create the requisite number of inner subdomains on PETSC_COMM_SELF (for now).
331a35b7b57SDmitry Karpeev        */
332a35b7b57SDmitry Karpeev       ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->overlap,osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr);
333a35b7b57SDmitry Karpeev     }
33406b43e7eSDmitry Karpeev     if (osm->N == PETSC_DECIDE) {
3356ac3741eSJed Brown       struct {PetscInt max,sum;} inwork,outwork;
336f746d493SDmitry Karpeev       /* determine global number of subdomains and the max number of local subdomains */
3376ac3741eSJed Brown       inwork.max = osm->n;
3386ac3741eSJed Brown       inwork.sum = osm->n;
3396ac3741eSJed Brown       ierr = MPI_Allreduce(&inwork,&outwork,1,MPIU_2INT,PetscMaxSum_Op,((PetscObject)pc)->comm);CHKERRQ(ierr);
3406ac3741eSJed Brown       osm->nmax = outwork.max;
3416ac3741eSJed Brown       osm->N    = outwork.sum;
342b1a0a8a3SJed Brown     }
3436a4f0f73SDmitry Karpeev 
344b1a0a8a3SJed Brown     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
345b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
346f746d493SDmitry Karpeev     ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&flg,PETSC_NULL);CHKERRQ(ierr);
347f746d493SDmitry Karpeev     if (flg) { ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); }
348b1a0a8a3SJed Brown 
349b1a0a8a3SJed Brown     if (osm->sort_indices) {
350f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
3516a4f0f73SDmitry Karpeev         ierr = ISSort(osm->ois[i]);CHKERRQ(ierr);
3526a4f0f73SDmitry Karpeev         ierr = ISSort(osm->iis[i]);CHKERRQ(ierr);
353b1a0a8a3SJed Brown       }
354b1a0a8a3SJed Brown     }
3556a4f0f73SDmitry Karpeev     /*
3566a4f0f73SDmitry Karpeev      Merge the ISs, create merged vectors and restrictions.
3576a4f0f73SDmitry Karpeev      */
3586a4f0f73SDmitry Karpeev     /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */
3596a4f0f73SDmitry Karpeev     on = 0;
360f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
3616a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
3626a4f0f73SDmitry Karpeev       on += oni;
363f746d493SDmitry Karpeev     }
3646a4f0f73SDmitry Karpeev     ierr = PetscMalloc(on*sizeof(PetscInt), &oidx);CHKERRQ(ierr);
3656a4f0f73SDmitry Karpeev     on = 0;
366f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
3676a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
3686a4f0f73SDmitry Karpeev       ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr);
3696a4f0f73SDmitry Karpeev       ierr = PetscMemcpy(oidx+on, oidxi, sizeof(PetscInt)*oni);CHKERRQ(ierr);
3706a4f0f73SDmitry Karpeev       ierr = ISRestoreIndices(osm->ois[i], &oidxi);CHKERRQ(ierr);
3716a4f0f73SDmitry Karpeev       on += oni;
372f746d493SDmitry Karpeev     }
3736a4f0f73SDmitry Karpeev     ierr = ISCreateGeneral(((PetscObject)(pc))->comm, on, oidx, PETSC_OWN_POINTER, &gois);CHKERRQ(ierr);
374f746d493SDmitry Karpeev     ierr = MatGetVecs(pc->pmat,&x,&y);CHKERRQ(ierr);
3756a4f0f73SDmitry Karpeev     ierr = VecCreateMPI(((PetscObject)pc)->comm, on, PETSC_DECIDE, &osm->gx);CHKERRQ(ierr);
376f746d493SDmitry Karpeev     ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr);
3776a4f0f73SDmitry Karpeev     ierr = VecGetOwnershipRange(osm->gx, &gofirst, PETSC_NULL);CHKERRQ(ierr);
3786a4f0f73SDmitry Karpeev     ierr = ISCreateStride(((PetscObject)pc)->comm,on,gofirst,1, &goid);CHKERRQ(ierr);
3796a4f0f73SDmitry Karpeev     ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr);
3806a4f0f73SDmitry Karpeev     ierr = VecDestroy(&x);CHKERRQ(ierr);
3817105d80bSDmitry Karpeev     ierr = ISDestroy(&gois);CHKERRQ(ierr);
3826a4f0f73SDmitry Karpeev     /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */
3836a4f0f73SDmitry Karpeev     { PetscInt       ini;     /* Number of indices the i-th a local inner subdomain. */
3846a4f0f73SDmitry Karpeev       PetscInt       in;      /* Number of indices in the disjoint uniont of local inner subdomains. */
3856a4f0f73SDmitry Karpeev       PetscInt       *iidx;   /* Global indices in the merged local inner subdomain. */
3866a4f0f73SDmitry Karpeev       PetscInt       *ioidx;  /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
3876a4f0f73SDmitry Karpeev       IS             giis;    /* IS for the disjoint union of inner subdomains. */
3886a4f0f73SDmitry Karpeev       IS             giois;   /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
389f746d493SDmitry Karpeev       /**/
3906a4f0f73SDmitry Karpeev       in = 0;
391f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
3926a4f0f73SDmitry Karpeev         ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr);
3936a4f0f73SDmitry Karpeev         in += ini;
394f746d493SDmitry Karpeev       }
3956a4f0f73SDmitry Karpeev       ierr = PetscMalloc(in*sizeof(PetscInt), &iidx);CHKERRQ(ierr);
3966a4f0f73SDmitry Karpeev       ierr = PetscMalloc(in*sizeof(PetscInt), &ioidx);CHKERRQ(ierr);
3976a4f0f73SDmitry Karpeev       ierr = VecGetOwnershipRange(osm->gx,&gofirst, PETSC_NULL);CHKERRQ(ierr);
3986a4f0f73SDmitry Karpeev       in = 0;
3996a4f0f73SDmitry Karpeev       on = 0;
400f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
4016a4f0f73SDmitry Karpeev         const PetscInt *iidxi;        /* Global indices of the i-th local inner subdomain. */
4026a4f0f73SDmitry Karpeev         ISLocalToGlobalMapping ltogi; /* Map from global to local indices of the i-th outer local subdomain. */
4036a4f0f73SDmitry Karpeev         PetscInt       *ioidxi;       /* Local indices of the i-th local inner subdomain within the local outer subdomain. */
4046a4f0f73SDmitry Karpeev         PetscInt       ioni;          /* Number of indices in ioidxi; if ioni != ini the inner subdomain is not a subdomain of the outer subdomain (error). */
4056a4f0f73SDmitry Karpeev         PetscInt       k;
4066a4f0f73SDmitry Karpeev         ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr);
4076a4f0f73SDmitry Karpeev         ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
4086a4f0f73SDmitry Karpeev         ierr = ISGetIndices(osm->iis[i],&iidxi);CHKERRQ(ierr);
4096a4f0f73SDmitry Karpeev         ierr = PetscMemcpy(iidx+in, iidxi, sizeof(PetscInt)*ini);CHKERRQ(ierr);
4106a4f0f73SDmitry Karpeev         ierr = ISLocalToGlobalMappingCreateIS(osm->ois[i],&ltogi);CHKERRQ(ierr);
4116a4f0f73SDmitry Karpeev         ioidxi = ioidx+in;
4126a4f0f73SDmitry Karpeev         ierr = ISGlobalToLocalMappingApply(ltogi,IS_GTOLM_DROP,ini,iidxi,&ioni,ioidxi);CHKERRQ(ierr);
4136a4f0f73SDmitry Karpeev         ierr = ISLocalToGlobalMappingDestroy(&ltogi);CHKERRQ(ierr);
4146a4f0f73SDmitry Karpeev         ierr = ISRestoreIndices(osm->iis[i], &iidxi);CHKERRQ(ierr);
4156a4f0f73SDmitry 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);
4166a4f0f73SDmitry Karpeev         for (k = 0; k < ini; ++k) {
4176a4f0f73SDmitry Karpeev           ioidxi[k] += gofirst+on;
418f746d493SDmitry Karpeev         }
4196a4f0f73SDmitry Karpeev         in += ini;
4206a4f0f73SDmitry Karpeev         on += oni;
421f746d493SDmitry Karpeev       }
4226a4f0f73SDmitry Karpeev       ierr = ISCreateGeneral(((PetscObject)pc)->comm, in, iidx,  PETSC_OWN_POINTER, &giis);CHKERRQ(ierr);
4236a4f0f73SDmitry Karpeev       ierr = ISCreateGeneral(((PetscObject)pc)->comm, in, ioidx, PETSC_OWN_POINTER, &giois);CHKERRQ(ierr);
4246a4f0f73SDmitry Karpeev       ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr);
4256a4f0f73SDmitry Karpeev       ierr = VecDestroy(&y);CHKERRQ(ierr);
4266a4f0f73SDmitry Karpeev       ierr = ISDestroy(&giis);CHKERRQ(ierr);
4276a4f0f73SDmitry Karpeev       ierr = ISDestroy(&giois);CHKERRQ(ierr);
428b1a0a8a3SJed Brown     }
4296a4f0f73SDmitry Karpeev     ierr = ISDestroy(&goid);CHKERRQ(ierr);
430f746d493SDmitry Karpeev     /* Create the subdomain work vectors. */
431f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->x);CHKERRQ(ierr);
432f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->y);CHKERRQ(ierr);
433f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr);
434f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr);
4356a4f0f73SDmitry Karpeev     for (i=0, on=0; i<osm->n; ++i, on += oni) {
4366a4f0f73SDmitry Karpeev       PetscInt oNi;
4376a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
4386a4f0f73SDmitry Karpeev       ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr);
4396a4f0f73SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr);
4406a4f0f73SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr);
441b1a0a8a3SJed Brown     }
442f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr);
443f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr);
444b1a0a8a3SJed Brown     /* Create the local solvers */
445f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(KSP *),&osm->ksp);CHKERRQ(ierr);
44644fe18e5SDmitry Karpeev     for (i=0; i<osm->n; i++) { /* KSPs are local */
4476a4f0f73SDmitry Karpeev       ierr = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr);
448b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr);
449b1a0a8a3SJed Brown       ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
450b1a0a8a3SJed Brown       ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
451b1a0a8a3SJed Brown       ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
452b1a0a8a3SJed Brown       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
453b1a0a8a3SJed Brown       ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
454b1a0a8a3SJed Brown       ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
455b1a0a8a3SJed Brown       osm->ksp[i] = ksp;
456b1a0a8a3SJed Brown     }
457b1a0a8a3SJed Brown     scall = MAT_INITIAL_MATRIX;
458b1a0a8a3SJed Brown 
4596a4f0f73SDmitry Karpeev   }/*if (!pc->setupcalled)*/
4606a4f0f73SDmitry Karpeev   else {
461b1a0a8a3SJed Brown     /*
462b1a0a8a3SJed Brown        Destroy the blocks from the previous iteration
463b1a0a8a3SJed Brown     */
464b1a0a8a3SJed Brown     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
465f746d493SDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
466b1a0a8a3SJed Brown       scall = MAT_INITIAL_MATRIX;
467b1a0a8a3SJed Brown     }
468b1a0a8a3SJed Brown   }
469b1a0a8a3SJed Brown 
470b1a0a8a3SJed Brown   /*
471f746d493SDmitry Karpeev      Extract out the submatrices.
472b1a0a8a3SJed Brown   */
47381a5c4aaSDmitry Karpeev   if (size > 1) {
4746a4f0f73SDmitry Karpeev     ierr = MatGetSubMatricesParallel(pc->pmat,osm->n,osm->ois, osm->ois,scall,&osm->pmat);CHKERRQ(ierr);
47581a5c4aaSDmitry Karpeev   }
47681a5c4aaSDmitry Karpeev   else {
4776a4f0f73SDmitry Karpeev     ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->ois, osm->ois,scall,&osm->pmat);CHKERRQ(ierr);
47881a5c4aaSDmitry Karpeev   }
479b1a0a8a3SJed Brown   if (scall == MAT_INITIAL_MATRIX) {
480b1a0a8a3SJed Brown     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
481f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
482b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,osm->pmat[i]);CHKERRQ(ierr);
483b1a0a8a3SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
484b1a0a8a3SJed Brown     }
485b1a0a8a3SJed Brown   }
486b1a0a8a3SJed Brown 
487b1a0a8a3SJed Brown   /* Return control to the user so that the submatrices can be modified (e.g., to apply
488b1a0a8a3SJed Brown      different boundary conditions for the submatrices than for the global problem) */
4896a4f0f73SDmitry Karpeev   ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
490b1a0a8a3SJed Brown 
491b1a0a8a3SJed Brown   /*
4926a4f0f73SDmitry Karpeev      Loop over submatrices putting them into local ksps
493b1a0a8a3SJed Brown   */
494f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
495b1a0a8a3SJed Brown     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr);
496b1a0a8a3SJed Brown     if (!pc->setupcalled) {
497b1a0a8a3SJed Brown       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
498b1a0a8a3SJed Brown     }
499b1a0a8a3SJed Brown   }
500b1a0a8a3SJed Brown 
501b1a0a8a3SJed Brown   PetscFunctionReturn(0);
502b1a0a8a3SJed Brown }
503b1a0a8a3SJed Brown 
504b1a0a8a3SJed Brown #undef __FUNCT__
505f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM"
506f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
507b1a0a8a3SJed Brown {
508f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
509b1a0a8a3SJed Brown   PetscErrorCode ierr;
510b1a0a8a3SJed Brown   PetscInt       i;
511b1a0a8a3SJed Brown 
512b1a0a8a3SJed Brown   PetscFunctionBegin;
513f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
514b1a0a8a3SJed Brown     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
515b1a0a8a3SJed Brown   }
516b1a0a8a3SJed Brown   PetscFunctionReturn(0);
517b1a0a8a3SJed Brown }
518b1a0a8a3SJed Brown 
519b1a0a8a3SJed Brown #undef __FUNCT__
520f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM"
521f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y)
522b1a0a8a3SJed Brown {
523f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
524b1a0a8a3SJed Brown   PetscErrorCode ierr;
525f746d493SDmitry Karpeev   PetscInt       i;
526b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
527b1a0a8a3SJed Brown 
528b1a0a8a3SJed Brown   PetscFunctionBegin;
529b1a0a8a3SJed Brown   /*
5306a4f0f73SDmitry Karpeev      Support for limiting the restriction or interpolation only to the inner
531b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
532b1a0a8a3SJed Brown   */
533f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
534b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
535f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
5366a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
537b1a0a8a3SJed Brown   }
5386a4f0f73SDmitry Karpeev   else {
5396a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
540b1a0a8a3SJed Brown   }
5416a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr);
5426a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
5436a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5446a4f0f73SDmitry Karpeev   }
5456a4f0f73SDmitry Karpeev   else {
5466a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5476a4f0f73SDmitry Karpeev   }
54886b96d36SDmitry Karpeev   /* do the subdomain solves */
54986b96d36SDmitry Karpeev   for (i=0; i<osm->n; ++i) {
550b1a0a8a3SJed Brown     ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
551b1a0a8a3SJed Brown   }
5526a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(y);CHKERRQ(ierr);
5536a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
5546a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5556a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);  PetscFunctionReturn(0);
5566a4f0f73SDmitry Karpeev   }
5576a4f0f73SDmitry Karpeev   else {
5586a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
5596a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);  PetscFunctionReturn(0);
5606a4f0f73SDmitry Karpeev   }
561b1a0a8a3SJed Brown }
562b1a0a8a3SJed Brown 
563b1a0a8a3SJed Brown #undef __FUNCT__
564f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM"
565f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y)
566b1a0a8a3SJed Brown {
567f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
568b1a0a8a3SJed Brown   PetscErrorCode ierr;
569f746d493SDmitry Karpeev   PetscInt       i;
570b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
571b1a0a8a3SJed Brown 
572b1a0a8a3SJed Brown   PetscFunctionBegin;
573b1a0a8a3SJed Brown   /*
574b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
575b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
576b1a0a8a3SJed Brown 
577f746d493SDmitry Karpeev      Note: these are reversed from the PCApply_GASM() because we are applying the
578b1a0a8a3SJed Brown      transpose of the three terms
579b1a0a8a3SJed Brown   */
580f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
581b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
582f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
5836a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
584b1a0a8a3SJed Brown   }
5856a4f0f73SDmitry Karpeev   else {
5866a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
587b1a0a8a3SJed Brown   }
5886a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr);
5896a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
5906a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5916a4f0f73SDmitry Karpeev   }
5926a4f0f73SDmitry Karpeev   else {
5936a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
5946a4f0f73SDmitry Karpeev   }
595b1a0a8a3SJed Brown   /* do the local solves */
596ab3e923fSDmitry 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. */
597b1a0a8a3SJed Brown     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
598b1a0a8a3SJed Brown   }
5996a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(y);CHKERRQ(ierr);
6006a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
6016a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6026a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6036a4f0f73SDmitry Karpeev   }
6046a4f0f73SDmitry Karpeev   else {
6056a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6066a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
6076a4f0f73SDmitry Karpeev   }
6086a4f0f73SDmitry Karpeev 
609b1a0a8a3SJed Brown   PetscFunctionReturn(0);
610b1a0a8a3SJed Brown }
611b1a0a8a3SJed Brown 
612b1a0a8a3SJed Brown #undef __FUNCT__
613a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM"
614a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc)
615b1a0a8a3SJed Brown {
616f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
617b1a0a8a3SJed Brown   PetscErrorCode ierr;
618b1a0a8a3SJed Brown   PetscInt       i;
619b1a0a8a3SJed Brown 
620b1a0a8a3SJed Brown   PetscFunctionBegin;
621b1a0a8a3SJed Brown   if (osm->ksp) {
622f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
623a06653b4SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
624b1a0a8a3SJed Brown     }
625b1a0a8a3SJed Brown   }
626b1a0a8a3SJed Brown   if (osm->pmat) {
627f746d493SDmitry Karpeev     if (osm->n > 0) {
628ab3e923fSDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
629b1a0a8a3SJed Brown     }
630b1a0a8a3SJed Brown   }
631d34fcf5fSBarry Smith   if (osm->x) {
632f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
633fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
634fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
635b1a0a8a3SJed Brown     }
636d34fcf5fSBarry Smith   }
6376bf464f9SBarry Smith   ierr = VecDestroy(&osm->gx);CHKERRQ(ierr);
6386bf464f9SBarry Smith   ierr = VecDestroy(&osm->gy);CHKERRQ(ierr);
639ab3e923fSDmitry Karpeev 
6406a4f0f73SDmitry Karpeev   ierr = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr);
6416a4f0f73SDmitry Karpeev   ierr = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr);
6426a4f0f73SDmitry Karpeev   ierr = PCGASMDestroySubdomains(osm->n,osm->ois,osm->iis);CHKERRQ(ierr);
6436a4f0f73SDmitry Karpeev   osm->ois = 0;
6446a4f0f73SDmitry Karpeev   osm->iis = 0;
645a06653b4SBarry Smith   PetscFunctionReturn(0);
646a06653b4SBarry Smith }
647a06653b4SBarry Smith 
648a06653b4SBarry Smith #undef __FUNCT__
649a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM"
650a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc)
651a06653b4SBarry Smith {
652a06653b4SBarry Smith   PC_GASM         *osm = (PC_GASM*)pc->data;
653a06653b4SBarry Smith   PetscErrorCode ierr;
654a06653b4SBarry Smith   PetscInt       i;
655a06653b4SBarry Smith 
656a06653b4SBarry Smith   PetscFunctionBegin;
657a06653b4SBarry Smith   ierr = PCReset_GASM(pc);CHKERRQ(ierr);
658a06653b4SBarry Smith   if (osm->ksp) {
659a06653b4SBarry Smith     for (i=0; i<osm->n; i++) {
6606bf464f9SBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
661a06653b4SBarry Smith     }
662a06653b4SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
663a06653b4SBarry Smith   }
664a06653b4SBarry Smith   ierr = PetscFree(osm->x);CHKERRQ(ierr);
665a06653b4SBarry Smith   ierr = PetscFree(osm->y);CHKERRQ(ierr);
666c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
667b1a0a8a3SJed Brown   PetscFunctionReturn(0);
668b1a0a8a3SJed Brown }
669b1a0a8a3SJed Brown 
670b1a0a8a3SJed Brown #undef __FUNCT__
671f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM"
672a6dfd86eSKarl Rupp static PetscErrorCode PCSetFromOptions_GASM(PC pc)
673a6dfd86eSKarl Rupp {
674f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
675b1a0a8a3SJed Brown   PetscErrorCode ierr;
676b1a0a8a3SJed Brown   PetscInt       blocks,ovl;
677b1a0a8a3SJed Brown   PetscBool      symset,flg;
678f746d493SDmitry Karpeev   PCGASMType      gasmtype;
679b1a0a8a3SJed Brown 
680b1a0a8a3SJed Brown   PetscFunctionBegin;
681b1a0a8a3SJed Brown   /* set the type to symmetric if matrix is symmetric */
682b1a0a8a3SJed Brown   if (!osm->type_set && pc->pmat) {
683b1a0a8a3SJed Brown     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
684f746d493SDmitry Karpeev     if (symset && flg) { osm->type = PC_GASM_BASIC; }
685b1a0a8a3SJed Brown   }
68644fe18e5SDmitry Karpeev   ierr = PetscOptionsHead("Generalized additive Schwarz options");CHKERRQ(ierr);
6876a4f0f73SDmitry Karpeev     ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
6886a4f0f73SDmitry Karpeev     osm->create_local = PETSC_TRUE;
6896a4f0f73SDmitry Karpeev     ierr = PetscOptionsBool("-pc_gasm_subdomains_create_local","Whether to make autocreated subdomains local (true by default)","PCGASMSetTotalSubdomains",osm->create_local,&osm->create_local,&flg);CHKERRQ(ierr);
6906a4f0f73SDmitry Karpeev     if (!osm->create_local) SETERRQ(((PetscObject)pc)->comm, PETSC_ERR_SUP, "No support for autocreation of nonlocal subdomains yet.");
6916a4f0f73SDmitry Karpeev 
6926a4f0f73SDmitry Karpeev     if (flg) {ierr = PCGASMSetTotalSubdomains(pc,blocks,osm->create_local);CHKERRQ(ierr); }
69306b43e7eSDmitry Karpeev     ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
694f746d493SDmitry Karpeev     if (flg) {ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); }
695b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
696f746d493SDmitry Karpeev     ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr);
697f746d493SDmitry Karpeev     if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr); }
698b1a0a8a3SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
699b1a0a8a3SJed Brown   PetscFunctionReturn(0);
700b1a0a8a3SJed Brown }
701b1a0a8a3SJed Brown 
702b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/
703b1a0a8a3SJed Brown 
704b1a0a8a3SJed Brown EXTERN_C_BEGIN
705b1a0a8a3SJed Brown #undef __FUNCT__
70606b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM"
707a35b7b57SDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[])
708b1a0a8a3SJed Brown {
709f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
710b1a0a8a3SJed Brown   PetscErrorCode ierr;
711b1a0a8a3SJed Brown   PetscInt       i;
712b1a0a8a3SJed Brown 
713b1a0a8a3SJed Brown   PetscFunctionBegin;
71469ca1f37SDmitry Karpeev   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, n = %D",n);
715a35b7b57SDmitry Karpeev   if (pc->setupcalled && (n != osm->n || iis || ois)) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp().");
716b1a0a8a3SJed Brown 
717b1a0a8a3SJed Brown   if (!pc->setupcalled) {
71893c1ec32SDmitry Karpeev     osm->n       = n;
7196a4f0f73SDmitry Karpeev     osm->ois     = 0;
7206a4f0f73SDmitry Karpeev     osm->iis     = 0;
721a35b7b57SDmitry Karpeev     if (ois) {
722a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);}
723b1a0a8a3SJed Brown     }
724a35b7b57SDmitry Karpeev     if (iis) {
725a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);}
726b1a0a8a3SJed Brown     }
7276a4f0f73SDmitry Karpeev     ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr);
728a35b7b57SDmitry Karpeev     if (ois) {
7296a4f0f73SDmitry Karpeev       ierr = PetscMalloc(n*sizeof(IS),&osm->ois);CHKERRQ(ierr);
730a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) { osm->ois[i] = ois[i]; }
7316a4f0f73SDmitry Karpeev       /* Flag indicating that the user has set outer subdomains, so PCGASM should not increase their size. */
732b1a0a8a3SJed Brown       osm->overlap = -1;
733a35b7b57SDmitry Karpeev       if (!iis) {
7346a4f0f73SDmitry Karpeev         ierr = PetscMalloc(n*sizeof(IS),&osm->iis);CHKERRQ(ierr);
735a35b7b57SDmitry Karpeev         for (i=0; i<n; i++) {
736a35b7b57SDmitry Karpeev           for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);}
737a35b7b57SDmitry Karpeev           osm->iis[i] = ois[i];
738a35b7b57SDmitry Karpeev         }
739a35b7b57SDmitry Karpeev       }
740a35b7b57SDmitry Karpeev     }
741a35b7b57SDmitry Karpeev     if (iis) {
742a35b7b57SDmitry Karpeev       ierr = PetscMalloc(n*sizeof(IS),&osm->iis);CHKERRQ(ierr);
743a35b7b57SDmitry Karpeev       for (i=0; i<n; i++) { osm->iis[i] = iis[i]; }
744a35b7b57SDmitry Karpeev       if (!ois) {
745a35b7b57SDmitry Karpeev         ierr = PetscMalloc(n*sizeof(IS),&osm->ois);CHKERRQ(ierr);
746a35b7b57SDmitry Karpeev         for (i=0; i<n; i++) {
747a35b7b57SDmitry Karpeev           for (i=0; i<n; i++) {
748a35b7b57SDmitry Karpeev             ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);
749a35b7b57SDmitry Karpeev             osm->ois[i] = iis[i];
750a35b7b57SDmitry Karpeev           }
751a35b7b57SDmitry Karpeev         }
752a35b7b57SDmitry Karpeev         if (osm->overlap > 0) {
753a35b7b57SDmitry Karpeev           /* Extend the "overlapping" regions by a number of steps */
754a35b7b57SDmitry Karpeev           ierr = MatIncreaseOverlap(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr);
755a35b7b57SDmitry Karpeev         }
756a35b7b57SDmitry Karpeev       }
757b1a0a8a3SJed Brown     }
758b1a0a8a3SJed Brown   }
759b1a0a8a3SJed Brown   PetscFunctionReturn(0);
760b1a0a8a3SJed Brown }
761b1a0a8a3SJed Brown EXTERN_C_END
762b1a0a8a3SJed Brown 
763b1a0a8a3SJed Brown EXTERN_C_BEGIN
764b1a0a8a3SJed Brown #undef __FUNCT__
7656a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains_GASM"
7660adebc6cSBarry Smith PetscErrorCode  PCGASMSetTotalSubdomains_GASM(PC pc,PetscInt N, PetscBool create_local)
7670adebc6cSBarry Smith {
768f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
769b1a0a8a3SJed Brown   PetscErrorCode ierr;
770b1a0a8a3SJed Brown   PetscMPIInt    rank,size;
771b1a0a8a3SJed Brown   PetscInt       n;
7726a4f0f73SDmitry Karpeev   PetscInt       Nmin, Nmax;
773b1a0a8a3SJed Brown   PetscFunctionBegin;
7746a4f0f73SDmitry Karpeev   if (!create_local) SETERRQ(((PetscObject)pc)->comm, PETSC_ERR_SUP, "No suppor for autocreation of nonlocal subdomains.");
7756a4f0f73SDmitry Karpeev   if (N < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be > 0, N = %D",N);
776adcec1e5SJed Brown   ierr = MPI_Allreduce(&N,&Nmin,1,MPIU_INT,MPIU_MIN,((PetscObject)pc)->comm);CHKERRQ(ierr);
777adcec1e5SJed Brown   ierr = MPI_Allreduce(&N,&Nmax,1,MPIU_INT,MPIU_MAX,((PetscObject)pc)->comm);CHKERRQ(ierr);
778f23aa3ddSBarry Smith   if (Nmin != Nmax) SETERRQ2(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONG, "All processors must use the same number of subdomains.  min(N) = %D != %D = max(N)", Nmin, Nmax);
779b1a0a8a3SJed Brown 
7806a4f0f73SDmitry Karpeev   osm->create_local = create_local;
781b1a0a8a3SJed Brown   /*
782b1a0a8a3SJed Brown      Split the subdomains equally among all processors
783b1a0a8a3SJed Brown   */
784b1a0a8a3SJed Brown   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
785b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
786b1a0a8a3SJed Brown   n = N/size + ((N % size) > rank);
7876a4f0f73SDmitry 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);
788f746d493SDmitry Karpeev   if (pc->setupcalled && n != osm->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before PCSetUp().");
789b1a0a8a3SJed Brown   if (!pc->setupcalled) {
7906a4f0f73SDmitry Karpeev     ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr);
79193c1ec32SDmitry Karpeev     osm->N            = N;
792f746d493SDmitry Karpeev     osm->n            = n;
79306b43e7eSDmitry Karpeev     osm->nmax         = N/size + ((N%size)?1:0);
7946a4f0f73SDmitry Karpeev     osm->ois           = 0;
7956a4f0f73SDmitry Karpeev     osm->iis     = 0;
796b1a0a8a3SJed Brown   }
797b1a0a8a3SJed Brown   PetscFunctionReturn(0);
79806b43e7eSDmitry Karpeev }
799b1a0a8a3SJed Brown EXTERN_C_END
800b1a0a8a3SJed Brown 
801b1a0a8a3SJed Brown EXTERN_C_BEGIN
802b1a0a8a3SJed Brown #undef __FUNCT__
803f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM"
8047087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
805b1a0a8a3SJed Brown {
806f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
807b1a0a8a3SJed Brown 
808b1a0a8a3SJed Brown   PetscFunctionBegin;
809b1a0a8a3SJed Brown   if (ovl < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
810f746d493SDmitry Karpeev   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
811b1a0a8a3SJed Brown   if (!pc->setupcalled) {
812b1a0a8a3SJed Brown     osm->overlap = ovl;
813b1a0a8a3SJed Brown   }
814b1a0a8a3SJed Brown   PetscFunctionReturn(0);
815b1a0a8a3SJed Brown }
816b1a0a8a3SJed Brown EXTERN_C_END
817b1a0a8a3SJed Brown 
818b1a0a8a3SJed Brown EXTERN_C_BEGIN
819b1a0a8a3SJed Brown #undef __FUNCT__
820f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM"
8217087cfbeSBarry Smith PetscErrorCode  PCGASMSetType_GASM(PC pc,PCGASMType type)
822b1a0a8a3SJed Brown {
823f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
824b1a0a8a3SJed Brown 
825b1a0a8a3SJed Brown   PetscFunctionBegin;
826b1a0a8a3SJed Brown   osm->type     = type;
827b1a0a8a3SJed Brown   osm->type_set = PETSC_TRUE;
828b1a0a8a3SJed Brown   PetscFunctionReturn(0);
829b1a0a8a3SJed Brown }
830b1a0a8a3SJed Brown EXTERN_C_END
831b1a0a8a3SJed Brown 
832b1a0a8a3SJed Brown EXTERN_C_BEGIN
833b1a0a8a3SJed Brown #undef __FUNCT__
834f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM"
8357087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices_GASM(PC pc,PetscBool  doSort)
836b1a0a8a3SJed Brown {
837f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
838b1a0a8a3SJed Brown 
839b1a0a8a3SJed Brown   PetscFunctionBegin;
840b1a0a8a3SJed Brown   osm->sort_indices = doSort;
841b1a0a8a3SJed Brown   PetscFunctionReturn(0);
842b1a0a8a3SJed Brown }
843b1a0a8a3SJed Brown EXTERN_C_END
844b1a0a8a3SJed Brown 
845b1a0a8a3SJed Brown EXTERN_C_BEGIN
846b1a0a8a3SJed Brown #undef __FUNCT__
847f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM"
84844fe18e5SDmitry Karpeev /*
84944fe18e5SDmitry Karpeev    FIX: This routine might need to be modified once multiple ranks per subdomain are allowed.
85044fe18e5SDmitry Karpeev         In particular, it would upset the global subdomain number calculation.
85144fe18e5SDmitry Karpeev */
8527087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
853b1a0a8a3SJed Brown {
854f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
855b1a0a8a3SJed Brown   PetscErrorCode ierr;
856b1a0a8a3SJed Brown 
857b1a0a8a3SJed Brown   PetscFunctionBegin;
858ab3e923fSDmitry Karpeev   if (osm->n < 1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here");
859b1a0a8a3SJed Brown 
860ab3e923fSDmitry Karpeev   if (n) {
861ab3e923fSDmitry Karpeev     *n = osm->n;
862b1a0a8a3SJed Brown   }
863ab3e923fSDmitry Karpeev   if (first) {
864ab3e923fSDmitry Karpeev     ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr);
865ab3e923fSDmitry Karpeev     *first -= osm->n;
866b1a0a8a3SJed Brown   }
867b1a0a8a3SJed Brown   if (ksp) {
868b1a0a8a3SJed Brown     /* Assume that local solves are now different; not necessarily
86906b43e7eSDmitry Karpeev        true, though!  This flag is used only for PCView_GASM() */
870b1a0a8a3SJed Brown     *ksp                   = osm->ksp;
8716a4f0f73SDmitry Karpeev     osm->same_subdomain_solvers = PETSC_FALSE;
872b1a0a8a3SJed Brown   }
873b1a0a8a3SJed Brown   PetscFunctionReturn(0);
874ab3e923fSDmitry Karpeev }/* PCGASMGetSubKSP_GASM() */
875b1a0a8a3SJed Brown EXTERN_C_END
876b1a0a8a3SJed Brown 
877b1a0a8a3SJed Brown 
878b1a0a8a3SJed Brown #undef __FUNCT__
87906b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains"
880b1a0a8a3SJed Brown /*@C
88106b43e7eSDmitry Karpeev     PCGASMSetSubdomains - Sets the subdomains for this processor
88206b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
883b1a0a8a3SJed Brown 
884b1a0a8a3SJed Brown     Collective on PC
885b1a0a8a3SJed Brown 
886b1a0a8a3SJed Brown     Input Parameters:
887b1a0a8a3SJed Brown +   pc  - the preconditioner context
88806b43e7eSDmitry Karpeev .   n   - the number of subdomains for this processor
8896a4f0f73SDmitry Karpeev .   iis - the index sets that define this processor's local inner subdomains
890b1a0a8a3SJed Brown          (or PETSC_NULL for PETSc to determine subdomains)
8916a4f0f73SDmitry Karpeev -   ois- the index sets that define this processor's local outer subdomains
8926a4f0f73SDmitry Karpeev          (or PETSC_NULL to use the same as iis)
893b1a0a8a3SJed Brown 
894b1a0a8a3SJed Brown     Notes:
89506b43e7eSDmitry Karpeev     The IS indices use the parallel, global numbering of the vector entries.
8966a4f0f73SDmitry Karpeev     Inner subdomains are those where the correction is applied.
8976a4f0f73SDmitry Karpeev     Outer subdomains are those where the residual necessary to obtain the
8986a4f0f73SDmitry Karpeev     corrections is obtained (see PCGASMType for the use of inner/outer subdomains).
8996a4f0f73SDmitry Karpeev     Both inner and outer subdomains can extend over several processors.
9006a4f0f73SDmitry Karpeev     This processor's portion of a subdomain is known as a local subdomain.
9016a4f0f73SDmitry Karpeev 
9026a4f0f73SDmitry Karpeev     By default the GASM preconditioner uses 1 (local) subdomain per processor.
9036a4f0f73SDmitry Karpeev     Use PCGASMSetTotalSubdomains() to set the total number of subdomains across
9046a4f0f73SDmitry Karpeev     all processors that PCGASM will create automatically, and to specify whether
9056a4f0f73SDmitry Karpeev     they should be local or not.
9066a4f0f73SDmitry Karpeev 
907b1a0a8a3SJed Brown 
908b1a0a8a3SJed Brown     Level: advanced
909b1a0a8a3SJed Brown 
91006b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz
911b1a0a8a3SJed Brown 
9126a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
91306b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
914b1a0a8a3SJed Brown @*/
9156a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[])
916b1a0a8a3SJed Brown {
917b1a0a8a3SJed Brown   PetscErrorCode ierr;
918b1a0a8a3SJed Brown 
919b1a0a8a3SJed Brown   PetscFunctionBegin;
920b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9216a4f0f73SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr);
922b1a0a8a3SJed Brown   PetscFunctionReturn(0);
923b1a0a8a3SJed Brown }
924b1a0a8a3SJed Brown 
925b1a0a8a3SJed Brown #undef __FUNCT__
9266a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains"
927b1a0a8a3SJed Brown /*@C
9286a4f0f73SDmitry Karpeev     PCGASMSetTotalSubdomains - Sets the total number of subdomains to use in the generalized additive
9296a4f0f73SDmitry Karpeev     Schwarz preconditioner.  The number of subdomains is cumulative across all processors in pc's
9306a4f0f73SDmitry Karpeev     communicator. Either all or no processors in the PC communicator must call this routine with
9316a4f0f73SDmitry Karpeev     the same N.  The subdomains will be created automatically during PCSetUp().
932b1a0a8a3SJed Brown 
933b1a0a8a3SJed Brown     Collective on PC
934b1a0a8a3SJed Brown 
935b1a0a8a3SJed Brown     Input Parameters:
936b1a0a8a3SJed Brown +   pc           - the preconditioner context
9376a4f0f73SDmitry Karpeev .   N            - the total number of subdomains cumulative across all processors
9386a4f0f73SDmitry Karpeev -   create_local - whether the subdomains to be created are to be local
939b1a0a8a3SJed Brown 
940b1a0a8a3SJed Brown     Options Database Key:
9416a4f0f73SDmitry Karpeev     To set the total number of subdomains and let PCGASM autocreate them, rather than specify the index sets, use the following options:
9426a4f0f73SDmitry Karpeev +    -pc_gasm_total_subdomains <n>                  - sets the total number of subdomains to be autocreated by PCGASM
9436a4f0f73SDmitry Karpeev -    -pc_gasm_subdomains_create_local <true|false>  - whether autocreated subdomains should be local or not (default is true)
944b1a0a8a3SJed Brown 
94506b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.
946b1a0a8a3SJed Brown 
9476a4f0f73SDmitry Karpeev 
9486a4f0f73SDmitry Karpeev     Use PCGASMSetSubdomains() to set subdomains explicitly or to set different numbers
9496a4f0f73SDmitry Karpeev     of subdomains per processor.
950b1a0a8a3SJed Brown 
951b1a0a8a3SJed Brown     Level: advanced
952b1a0a8a3SJed Brown 
953f746d493SDmitry Karpeev .keywords: PC, GASM, set, total, global, subdomains, additive Schwarz
954b1a0a8a3SJed Brown 
95506b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
956f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
957b1a0a8a3SJed Brown @*/
9586a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetTotalSubdomains(PC pc,PetscInt N, PetscBool create_local)
959b1a0a8a3SJed Brown {
960b1a0a8a3SJed Brown   PetscErrorCode ierr;
961b1a0a8a3SJed Brown 
962b1a0a8a3SJed Brown   PetscFunctionBegin;
963b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9646a4f0f73SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetTotalSubdomains_C",(PC,PetscInt,PetscBool),(pc,N,create_local));CHKERRQ(ierr);
965b1a0a8a3SJed Brown   PetscFunctionReturn(0);
966b1a0a8a3SJed Brown }
967b1a0a8a3SJed Brown 
968b1a0a8a3SJed Brown #undef __FUNCT__
969f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap"
970b1a0a8a3SJed Brown /*@
971f746d493SDmitry Karpeev     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
972b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
973b1a0a8a3SJed Brown     PC communicator must call this routine.
974b1a0a8a3SJed Brown 
975b1a0a8a3SJed Brown     Logically Collective on PC
976b1a0a8a3SJed Brown 
977b1a0a8a3SJed Brown     Input Parameters:
978b1a0a8a3SJed Brown +   pc  - the preconditioner context
979b1a0a8a3SJed Brown -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
980b1a0a8a3SJed Brown 
981b1a0a8a3SJed Brown     Options Database Key:
98206b43e7eSDmitry Karpeev .   -pc_gasm_overlap <overlap> - Sets overlap
983b1a0a8a3SJed Brown 
984b1a0a8a3SJed Brown     Notes:
98506b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.  To use
9866a4f0f73SDmitry Karpeev     multiple subdomain per perocessor, see PCGASMSetTotalSubdomains() or
9876a4f0f73SDmitry Karpeev     PCGASMSetSubdomains() (and the option -pc_gasm_total_subdomains <n>).
988b1a0a8a3SJed Brown 
989b1a0a8a3SJed Brown     The overlap defaults to 1, so if one desires that no additional
990b1a0a8a3SJed Brown     overlap be computed beyond what may have been set with a call to
9916a4f0f73SDmitry Karpeev     PCGASMSetTotalSubdomains() or PCGASMSetSubdomains(), then ovl
992b1a0a8a3SJed Brown     must be set to be 0.  In particular, if one does not explicitly set
99306b43e7eSDmitry Karpeev     the subdomains in application code, then all overlap would be computed
994f746d493SDmitry Karpeev     internally by PETSc, and using an overlap of 0 would result in an GASM
995b1a0a8a3SJed Brown     variant that is equivalent to the block Jacobi preconditioner.
996b1a0a8a3SJed Brown 
997b1a0a8a3SJed Brown     Note that one can define initial index sets with any overlap via
99806b43e7eSDmitry Karpeev     PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows
99906b43e7eSDmitry Karpeev     PETSc to extend that overlap further, if desired.
1000b1a0a8a3SJed Brown 
1001b1a0a8a3SJed Brown     Level: intermediate
1002b1a0a8a3SJed Brown 
1003f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap
1004b1a0a8a3SJed Brown 
10056a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
100606b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
1007b1a0a8a3SJed Brown @*/
10087087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap(PC pc,PetscInt ovl)
1009b1a0a8a3SJed Brown {
1010b1a0a8a3SJed Brown   PetscErrorCode ierr;
1011b1a0a8a3SJed Brown 
1012b1a0a8a3SJed Brown   PetscFunctionBegin;
1013b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1014b1a0a8a3SJed Brown   PetscValidLogicalCollectiveInt(pc,ovl,2);
1015f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
1016b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1017b1a0a8a3SJed Brown }
1018b1a0a8a3SJed Brown 
1019b1a0a8a3SJed Brown #undef __FUNCT__
1020f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType"
1021b1a0a8a3SJed Brown /*@
1022f746d493SDmitry Karpeev     PCGASMSetType - Sets the type of restriction and interpolation used
1023b1a0a8a3SJed Brown     for local problems in the additive Schwarz method.
1024b1a0a8a3SJed Brown 
1025b1a0a8a3SJed Brown     Logically Collective on PC
1026b1a0a8a3SJed Brown 
1027b1a0a8a3SJed Brown     Input Parameters:
1028b1a0a8a3SJed Brown +   pc  - the preconditioner context
1029f746d493SDmitry Karpeev -   type - variant of GASM, one of
1030b1a0a8a3SJed Brown .vb
1031f746d493SDmitry Karpeev       PC_GASM_BASIC       - full interpolation and restriction
1032f746d493SDmitry Karpeev       PC_GASM_RESTRICT    - full restriction, local processor interpolation
1033f746d493SDmitry Karpeev       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
1034f746d493SDmitry Karpeev       PC_GASM_NONE        - local processor restriction and interpolation
1035b1a0a8a3SJed Brown .ve
1036b1a0a8a3SJed Brown 
1037b1a0a8a3SJed Brown     Options Database Key:
1038f746d493SDmitry Karpeev .   -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1039b1a0a8a3SJed Brown 
1040b1a0a8a3SJed Brown     Level: intermediate
1041b1a0a8a3SJed Brown 
1042f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1043b1a0a8a3SJed Brown 
10446a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1045f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1046b1a0a8a3SJed Brown @*/
10477087cfbeSBarry Smith PetscErrorCode  PCGASMSetType(PC pc,PCGASMType type)
1048b1a0a8a3SJed Brown {
1049b1a0a8a3SJed Brown   PetscErrorCode ierr;
1050b1a0a8a3SJed Brown 
1051b1a0a8a3SJed Brown   PetscFunctionBegin;
1052b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1053b1a0a8a3SJed Brown   PetscValidLogicalCollectiveEnum(pc,type,2);
1054f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr);
1055b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1056b1a0a8a3SJed Brown }
1057b1a0a8a3SJed Brown 
1058b1a0a8a3SJed Brown #undef __FUNCT__
1059f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices"
1060b1a0a8a3SJed Brown /*@
1061f746d493SDmitry Karpeev     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
1062b1a0a8a3SJed Brown 
1063b1a0a8a3SJed Brown     Logically Collective on PC
1064b1a0a8a3SJed Brown 
1065b1a0a8a3SJed Brown     Input Parameters:
1066b1a0a8a3SJed Brown +   pc  - the preconditioner context
1067b1a0a8a3SJed Brown -   doSort - sort the subdomain indices
1068b1a0a8a3SJed Brown 
1069b1a0a8a3SJed Brown     Level: intermediate
1070b1a0a8a3SJed Brown 
1071f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1072b1a0a8a3SJed Brown 
10736a4f0f73SDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(),
1074f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1075b1a0a8a3SJed Brown @*/
10767087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices(PC pc,PetscBool  doSort)
1077b1a0a8a3SJed Brown {
1078b1a0a8a3SJed Brown   PetscErrorCode ierr;
1079b1a0a8a3SJed Brown 
1080b1a0a8a3SJed Brown   PetscFunctionBegin;
1081b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1082b1a0a8a3SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
1083f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
1084b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1085b1a0a8a3SJed Brown }
1086b1a0a8a3SJed Brown 
1087b1a0a8a3SJed Brown #undef __FUNCT__
1088f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP"
1089b1a0a8a3SJed Brown /*@C
1090f746d493SDmitry Karpeev    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
1091b1a0a8a3SJed Brown    this processor.
1092b1a0a8a3SJed Brown 
1093b1a0a8a3SJed Brown    Collective on PC iff first_local is requested
1094b1a0a8a3SJed Brown 
1095b1a0a8a3SJed Brown    Input Parameter:
1096b1a0a8a3SJed Brown .  pc - the preconditioner context
1097b1a0a8a3SJed Brown 
1098b1a0a8a3SJed Brown    Output Parameters:
1099b1a0a8a3SJed Brown +  n_local - the number of blocks on this processor or PETSC_NULL
1100b1a0a8a3SJed Brown .  first_local - the global number of the first block on this processor or PETSC_NULL,
1101b1a0a8a3SJed Brown                  all processors must request or all must pass PETSC_NULL
1102b1a0a8a3SJed Brown -  ksp - the array of KSP contexts
1103b1a0a8a3SJed Brown 
1104b1a0a8a3SJed Brown    Note:
1105f746d493SDmitry Karpeev    After PCGASMGetSubKSP() the array of KSPes is not to be freed
1106b1a0a8a3SJed Brown 
1107b1a0a8a3SJed Brown    Currently for some matrix implementations only 1 block per processor
1108b1a0a8a3SJed Brown    is supported.
1109b1a0a8a3SJed Brown 
1110f746d493SDmitry Karpeev    You must call KSPSetUp() before calling PCGASMGetSubKSP().
1111b1a0a8a3SJed Brown 
1112b1a0a8a3SJed Brown    Level: advanced
1113b1a0a8a3SJed Brown 
1114f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context
1115b1a0a8a3SJed Brown 
11166a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap(),
1117f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(),
1118b1a0a8a3SJed Brown @*/
11197087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1120b1a0a8a3SJed Brown {
1121b1a0a8a3SJed Brown   PetscErrorCode ierr;
1122b1a0a8a3SJed Brown 
1123b1a0a8a3SJed Brown   PetscFunctionBegin;
1124b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1125f746d493SDmitry Karpeev   ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
1126b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1127b1a0a8a3SJed Brown }
1128b1a0a8a3SJed Brown 
1129b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/
1130b1a0a8a3SJed Brown /*MC
1131f746d493SDmitry Karpeev    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1132b1a0a8a3SJed Brown            its own KSP object.
1133b1a0a8a3SJed Brown 
1134b1a0a8a3SJed Brown    Options Database Keys:
113506b43e7eSDmitry Karpeev +  -pc_gasm_total_block_count <n> - Sets total number of local subdomains (known as blocks) to be distributed among processors
113606b43e7eSDmitry Karpeev .  -pc_gasm_view_subdomains       - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view
113706b43e7eSDmitry Karpeev .  -pc_gasm_print_subdomains      - activates the printing of subdomain indices in PCSetUp()
113806b43e7eSDmitry Karpeev .  -pc_gasm_overlap <ovl>         - Sets overlap by which to (automatically) extend local subdomains
1139f746d493SDmitry Karpeev -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1140b1a0a8a3SJed Brown 
1141b1a0a8a3SJed Brown      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
1142f746d493SDmitry Karpeev       will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
1143f746d493SDmitry Karpeev       -pc_gasm_type basic to use the standard GASM.
1144b1a0a8a3SJed Brown 
1145b1a0a8a3SJed Brown    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
1146b1a0a8a3SJed Brown      than one processor. Defaults to one block per processor.
1147b1a0a8a3SJed Brown 
1148b1a0a8a3SJed Brown      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1149b1a0a8a3SJed Brown         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
1150b1a0a8a3SJed Brown 
1151f746d493SDmitry Karpeev      To set the options on the solvers separate for each block call PCGASMGetSubKSP()
1152b1a0a8a3SJed Brown          and set the options directly on the resulting KSP object (you can access its PC
1153b1a0a8a3SJed Brown          with KSPGetPC())
1154b1a0a8a3SJed Brown 
1155b1a0a8a3SJed Brown 
1156b1a0a8a3SJed Brown    Level: beginner
1157b1a0a8a3SJed Brown 
1158b1a0a8a3SJed Brown    Concepts: additive Schwarz method
1159b1a0a8a3SJed Brown 
1160b1a0a8a3SJed Brown     References:
1161b1a0a8a3SJed Brown     An additive variant of the Schwarz alternating method for the case of many subregions
1162b1a0a8a3SJed Brown     M Dryja, OB Widlund - Courant Institute, New York University Technical report
1163b1a0a8a3SJed Brown 
1164b1a0a8a3SJed Brown     Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1165b1a0a8a3SJed Brown     Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X.
1166b1a0a8a3SJed Brown 
1167b1a0a8a3SJed Brown .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
116806b43e7eSDmitry Karpeev            PCBJACOBI, PCGASMSetUseTrueLocal(), PCGASMGetSubKSP(), PCGASMSetSubdomains(),
11696a4f0f73SDmitry Karpeev            PCGASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType()
1170b1a0a8a3SJed Brown 
1171b1a0a8a3SJed Brown M*/
1172b1a0a8a3SJed Brown 
1173b1a0a8a3SJed Brown EXTERN_C_BEGIN
1174b1a0a8a3SJed Brown #undef __FUNCT__
1175f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM"
11767087cfbeSBarry Smith PetscErrorCode  PCCreate_GASM(PC pc)
1177b1a0a8a3SJed Brown {
1178b1a0a8a3SJed Brown   PetscErrorCode ierr;
1179f746d493SDmitry Karpeev   PC_GASM         *osm;
1180b1a0a8a3SJed Brown 
1181b1a0a8a3SJed Brown   PetscFunctionBegin;
1182f746d493SDmitry Karpeev   ierr = PetscNewLog(pc,PC_GASM,&osm);CHKERRQ(ierr);
1183ab3e923fSDmitry Karpeev   osm->N                 = PETSC_DECIDE;
118406b43e7eSDmitry Karpeev   osm->n                 = PETSC_DECIDE;
1185ab3e923fSDmitry Karpeev   osm->nmax              = 0;
1186b1a0a8a3SJed Brown   osm->overlap           = 1;
1187b1a0a8a3SJed Brown   osm->ksp               = 0;
11886a4f0f73SDmitry Karpeev   osm->gorestriction     = 0;
11896a4f0f73SDmitry Karpeev   osm->girestriction     = 0;
1190ab3e923fSDmitry Karpeev   osm->gx                = 0;
1191ab3e923fSDmitry Karpeev   osm->gy                = 0;
1192b1a0a8a3SJed Brown   osm->x                 = 0;
1193b1a0a8a3SJed Brown   osm->y                 = 0;
11946a4f0f73SDmitry Karpeev   osm->ois               = 0;
11956a4f0f73SDmitry Karpeev   osm->iis               = 0;
1196b1a0a8a3SJed Brown   osm->pmat              = 0;
1197f746d493SDmitry Karpeev   osm->type              = PC_GASM_RESTRICT;
11986a4f0f73SDmitry Karpeev   osm->same_subdomain_solvers = PETSC_TRUE;
1199b1a0a8a3SJed Brown   osm->sort_indices           = PETSC_TRUE;
1200b1a0a8a3SJed Brown 
1201b1a0a8a3SJed Brown   pc->data                   = (void*)osm;
1202f746d493SDmitry Karpeev   pc->ops->apply             = PCApply_GASM;
1203f746d493SDmitry Karpeev   pc->ops->applytranspose    = PCApplyTranspose_GASM;
1204f746d493SDmitry Karpeev   pc->ops->setup             = PCSetUp_GASM;
1205a06653b4SBarry Smith   pc->ops->reset             = PCReset_GASM;
1206f746d493SDmitry Karpeev   pc->ops->destroy           = PCDestroy_GASM;
1207f746d493SDmitry Karpeev   pc->ops->setfromoptions    = PCSetFromOptions_GASM;
1208f746d493SDmitry Karpeev   pc->ops->setuponblocks     = PCSetUpOnBlocks_GASM;
1209f746d493SDmitry Karpeev   pc->ops->view              = PCView_GASM;
1210b1a0a8a3SJed Brown   pc->ops->applyrichardson   = 0;
1211b1a0a8a3SJed Brown 
121206b43e7eSDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSubdomains_C","PCGASMSetSubdomains_GASM",
121306b43e7eSDmitry Karpeev                     PCGASMSetSubdomains_GASM);CHKERRQ(ierr);
12146a4f0f73SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetTotalSubdomains_C","PCGASMSetTotalSubdomains_GASM",
12156a4f0f73SDmitry Karpeev                     PCGASMSetTotalSubdomains_GASM);CHKERRQ(ierr);
1216f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetOverlap_C","PCGASMSetOverlap_GASM",
1217f746d493SDmitry Karpeev                     PCGASMSetOverlap_GASM);CHKERRQ(ierr);
1218f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetType_C","PCGASMSetType_GASM",
1219f746d493SDmitry Karpeev                     PCGASMSetType_GASM);CHKERRQ(ierr);
1220f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSortIndices_C","PCGASMSetSortIndices_GASM",
1221f746d493SDmitry Karpeev                     PCGASMSetSortIndices_GASM);CHKERRQ(ierr);
1222f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMGetSubKSP_C","PCGASMGetSubKSP_GASM",
1223f746d493SDmitry Karpeev                     PCGASMGetSubKSP_GASM);CHKERRQ(ierr);
1224b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1225b1a0a8a3SJed Brown }
1226b1a0a8a3SJed Brown EXTERN_C_END
1227b1a0a8a3SJed Brown 
1228b1a0a8a3SJed Brown 
1229b1a0a8a3SJed Brown #undef __FUNCT__
123006b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains"
1231b1a0a8a3SJed Brown /*@C
123206b43e7eSDmitry Karpeev    PCGASMCreateLocalSubdomains - Creates n local index sets for the overlapping
123306b43e7eSDmitry Karpeev    Schwarz preconditioner for a any problem based on its matrix.
1234b1a0a8a3SJed Brown 
1235b1a0a8a3SJed Brown    Collective
1236b1a0a8a3SJed Brown 
1237b1a0a8a3SJed Brown    Input Parameters:
1238b1a0a8a3SJed Brown +  A       - The global matrix operator
12396a4f0f73SDmitry Karpeev .  overlap - amount of overlap in outer subdomains
124006b43e7eSDmitry Karpeev -  n       - the number of local subdomains
1241b1a0a8a3SJed Brown 
1242b1a0a8a3SJed Brown    Output Parameters:
12436a4f0f73SDmitry Karpeev +  iis - the array of index sets defining the local inner subdomains (on which the correction is applied)
12446a4f0f73SDmitry Karpeev -  ois - the array of index sets defining the local outer subdomains (on which the residual is computed)
1245b1a0a8a3SJed Brown 
1246b1a0a8a3SJed Brown    Level: advanced
1247b1a0a8a3SJed Brown 
12486a4f0f73SDmitry Karpeev    Note: this generates n nonoverlapping local inner subdomains on PETSC_COMM_SELF;
12496a4f0f73SDmitry Karpeev          PCGASM will generate the overlap from these if you use them in PCGASMSetSubdomains() and set a
12506a4f0f73SDmitry Karpeev          nonzero overlap with PCGASMSetOverlap()
1251b1a0a8a3SJed Brown 
1252b1a0a8a3SJed Brown     In the Fortran version you must provide the array outis[] already allocated of length n.
1253b1a0a8a3SJed Brown 
1254f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1255b1a0a8a3SJed Brown 
125606b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains()
1257b1a0a8a3SJed Brown @*/
12586a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt overlap, PetscInt n, IS* iis[], IS* ois[])
1259b1a0a8a3SJed Brown {
1260b1a0a8a3SJed Brown   MatPartitioning           mpart;
1261b1a0a8a3SJed Brown   const char                *prefix;
126211bd1e4dSLisandro Dalcin   PetscErrorCode            (*f)(Mat,MatReuse,Mat*);
1263b1a0a8a3SJed Brown   PetscMPIInt               size;
1264b1a0a8a3SJed Brown   PetscInt                  i,j,rstart,rend,bs;
126511bd1e4dSLisandro Dalcin   PetscBool                 isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
1266b1a0a8a3SJed Brown   Mat                       Ad = PETSC_NULL, adj;
1267b1a0a8a3SJed Brown   IS                        ispart,isnumb,*is;
1268b1a0a8a3SJed Brown   PetscErrorCode            ierr;
1269b1a0a8a3SJed Brown 
1270b1a0a8a3SJed Brown   PetscFunctionBegin;
1271b1a0a8a3SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
12726a4f0f73SDmitry Karpeev   PetscValidPointer(iis,4);
1273b1a0a8a3SJed Brown   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
1274b1a0a8a3SJed Brown 
1275b1a0a8a3SJed Brown   /* Get prefix, row distribution, and block size */
1276b1a0a8a3SJed Brown   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1277b1a0a8a3SJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1278b1a0a8a3SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1279b1a0a8a3SJed 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);
1280b1a0a8a3SJed Brown 
1281b1a0a8a3SJed Brown   /* Get diagonal block from matrix if possible */
1282b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1283b1a0a8a3SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr);
1284b1a0a8a3SJed Brown   if (f) {
128511bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
1286b1a0a8a3SJed Brown   } else if (size == 1) {
128711bd1e4dSLisandro Dalcin     Ad = A;
1288b1a0a8a3SJed Brown   }
1289b1a0a8a3SJed Brown   if (Ad) {
1290251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1291251f4c67SDmitry Karpeev     if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
1292b1a0a8a3SJed Brown   }
1293b1a0a8a3SJed Brown   if (Ad && n > 1) {
1294b1a0a8a3SJed Brown     PetscBool  match,done;
1295b1a0a8a3SJed Brown     /* Try to setup a good matrix partitioning if available */
1296b1a0a8a3SJed Brown     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
1297b1a0a8a3SJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
1298b1a0a8a3SJed Brown     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1299251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
1300b1a0a8a3SJed Brown     if (!match) {
1301251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
1302b1a0a8a3SJed Brown     }
1303b1a0a8a3SJed Brown     if (!match) { /* assume a "good" partitioner is available */
13041a83f524SJed Brown       PetscInt na;
13051a83f524SJed Brown       const PetscInt *ia,*ja;
1306b1a0a8a3SJed Brown       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1307b1a0a8a3SJed Brown       if (done) {
1308b1a0a8a3SJed Brown         /* Build adjacency matrix by hand. Unfortunately a call to
1309b1a0a8a3SJed Brown            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1310b1a0a8a3SJed Brown            remove the block-aij structure and we cannot expect
1311b1a0a8a3SJed Brown            MatPartitioning to split vertices as we need */
13121a83f524SJed Brown         PetscInt i,j,len,nnz,cnt,*iia=0,*jja=0;
13131a83f524SJed Brown         const PetscInt *row;
1314b1a0a8a3SJed Brown         nnz = 0;
1315b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* count number of nonzeros */
1316b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1317b1a0a8a3SJed Brown           row = ja + ia[i];
1318b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1319b1a0a8a3SJed Brown             if (row[j] == i) { /* don't count diagonal */
1320b1a0a8a3SJed Brown               len--; break;
1321b1a0a8a3SJed Brown             }
1322b1a0a8a3SJed Brown           }
1323b1a0a8a3SJed Brown           nnz += len;
1324b1a0a8a3SJed Brown         }
1325b1a0a8a3SJed Brown         ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr);
1326b1a0a8a3SJed Brown         ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr);
1327b1a0a8a3SJed Brown         nnz    = 0;
1328b1a0a8a3SJed Brown         iia[0] = 0;
1329b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* fill adjacency */
1330b1a0a8a3SJed Brown           cnt = 0;
1331b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1332b1a0a8a3SJed Brown           row = ja + ia[i];
1333b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1334b1a0a8a3SJed Brown             if (row[j] != i) { /* if not diagonal */
1335b1a0a8a3SJed Brown               jja[nnz+cnt++] = row[j];
1336b1a0a8a3SJed Brown             }
1337b1a0a8a3SJed Brown           }
1338b1a0a8a3SJed Brown           nnz += cnt;
1339b1a0a8a3SJed Brown           iia[i+1] = nnz;
1340b1a0a8a3SJed Brown         }
1341b1a0a8a3SJed Brown         /* Partitioning of the adjacency matrix */
1342b1a0a8a3SJed Brown         ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);CHKERRQ(ierr);
1343b1a0a8a3SJed Brown         ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
1344b1a0a8a3SJed Brown         ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
1345b1a0a8a3SJed Brown         ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
1346b1a0a8a3SJed Brown         ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
13476bf464f9SBarry Smith         ierr = MatDestroy(&adj);CHKERRQ(ierr);
1348b1a0a8a3SJed Brown         foundpart = PETSC_TRUE;
1349b1a0a8a3SJed Brown       }
1350b1a0a8a3SJed Brown       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1351b1a0a8a3SJed Brown     }
13526bf464f9SBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1353b1a0a8a3SJed Brown   }
1354b1a0a8a3SJed Brown   ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr);
1355b1a0a8a3SJed Brown   if (!foundpart) {
1356b1a0a8a3SJed Brown 
1357b1a0a8a3SJed Brown     /* Partitioning by contiguous chunks of rows */
1358b1a0a8a3SJed Brown 
1359b1a0a8a3SJed Brown     PetscInt mbs   = (rend-rstart)/bs;
1360b1a0a8a3SJed Brown     PetscInt start = rstart;
1361b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1362b1a0a8a3SJed Brown       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
1363b1a0a8a3SJed Brown       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
1364b1a0a8a3SJed Brown       start += count;
1365b1a0a8a3SJed Brown     }
1366b1a0a8a3SJed Brown 
1367b1a0a8a3SJed Brown   } else {
1368b1a0a8a3SJed Brown 
1369b1a0a8a3SJed Brown     /* Partitioning by adjacency of diagonal block  */
1370b1a0a8a3SJed Brown 
1371b1a0a8a3SJed Brown     const PetscInt *numbering;
1372b1a0a8a3SJed Brown     PetscInt       *count,nidx,*indices,*newidx,start=0;
1373b1a0a8a3SJed Brown     /* Get node count in each partition */
1374b1a0a8a3SJed Brown     ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr);
1375b1a0a8a3SJed Brown     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
1376b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1377b1a0a8a3SJed Brown       for (i=0; i<n; i++) count[i] *= bs;
1378b1a0a8a3SJed Brown     }
1379b1a0a8a3SJed Brown     /* Build indices from node numbering */
1380b1a0a8a3SJed Brown     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1381b1a0a8a3SJed Brown     ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr);
1382b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1383b1a0a8a3SJed Brown     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
1384b1a0a8a3SJed Brown     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
1385b1a0a8a3SJed Brown     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
1386b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1387b1a0a8a3SJed Brown       ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr);
1388b1a0a8a3SJed Brown       for (i=0; i<nidx; i++)
1389b1a0a8a3SJed Brown         for (j=0; j<bs; j++)
1390b1a0a8a3SJed Brown           newidx[i*bs+j] = indices[i]*bs + j;
1391b1a0a8a3SJed Brown       ierr = PetscFree(indices);CHKERRQ(ierr);
1392b1a0a8a3SJed Brown       nidx   *= bs;
1393b1a0a8a3SJed Brown       indices = newidx;
1394b1a0a8a3SJed Brown     }
1395b1a0a8a3SJed Brown     /* Shift to get global indices */
1396b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] += rstart;
1397b1a0a8a3SJed Brown 
1398b1a0a8a3SJed Brown     /* Build the index sets for each block */
1399b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1400b1a0a8a3SJed Brown       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
1401b1a0a8a3SJed Brown       ierr   = ISSort(is[i]);CHKERRQ(ierr);
1402b1a0a8a3SJed Brown       start += count[i];
1403b1a0a8a3SJed Brown     }
1404b1a0a8a3SJed Brown 
1405b1a0a8a3SJed Brown     ierr = PetscFree(count);
1406b1a0a8a3SJed Brown     ierr = PetscFree(indices);
1407fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1408fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
1409b1a0a8a3SJed Brown   }
14106a4f0f73SDmitry Karpeev   *iis = is;
14116a4f0f73SDmitry Karpeev   if (!ois) PetscFunctionReturn(0);
14126a4f0f73SDmitry Karpeev   /*
14136a4f0f73SDmitry Karpeev    Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap
14146a4f0f73SDmitry Karpeev    has been requested, copy the inner subdomains over so they can be modified.
14156a4f0f73SDmitry Karpeev    */
14166a4f0f73SDmitry Karpeev   ierr = PetscMalloc(n*sizeof(IS),ois);CHKERRQ(ierr);
14176a4f0f73SDmitry Karpeev   for (i=0; i<n; ++i) {
14186a4f0f73SDmitry Karpeev     if (overlap > 0) { /* With positive overlap, (*iis)[i] will be modified */
14196a4f0f73SDmitry Karpeev       ierr = ISDuplicate((*iis)[i],(*ois)+i);CHKERRQ(ierr);
14206a4f0f73SDmitry Karpeev       ierr = ISCopy((*iis)[i],(*ois)[i]);CHKERRQ(ierr);
14216a4f0f73SDmitry Karpeev     } else {
14226a4f0f73SDmitry Karpeev       ierr = PetscObjectReference((PetscObject)(*iis)[i]);CHKERRQ(ierr);
14236a4f0f73SDmitry Karpeev       (*ois)[i] = (*iis)[i];
14246a4f0f73SDmitry Karpeev     }
14256a4f0f73SDmitry Karpeev   }
14266a4f0f73SDmitry Karpeev   if (overlap > 0) {
14276a4f0f73SDmitry Karpeev     /* Extend the "overlapping" regions by a number of steps */
14286a4f0f73SDmitry Karpeev     ierr = MatIncreaseOverlap(A,n,*ois,overlap);CHKERRQ(ierr);
14296a4f0f73SDmitry Karpeev   }
1430b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1431b1a0a8a3SJed Brown }
1432b1a0a8a3SJed Brown 
1433b1a0a8a3SJed Brown #undef __FUNCT__
1434f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains"
1435b1a0a8a3SJed Brown /*@C
1436f746d493SDmitry Karpeev    PCGASMDestroySubdomains - Destroys the index sets created with
143706b43e7eSDmitry Karpeev    PCGASMCreateLocalSubdomains() or PCGASMCreateSubdomains2D. Should be
143806b43e7eSDmitry Karpeev    called after setting subdomains with PCGASMSetSubdomains().
1439b1a0a8a3SJed Brown 
1440b1a0a8a3SJed Brown    Collective
1441b1a0a8a3SJed Brown 
1442b1a0a8a3SJed Brown    Input Parameters:
1443b1a0a8a3SJed Brown +  n   - the number of index sets
14446a4f0f73SDmitry Karpeev .  iis - the array of inner subdomains,
14456a4f0f73SDmitry Karpeev -  ois - the array of outer subdomains, can be PETSC_NULL
1446b1a0a8a3SJed Brown 
14476a4f0f73SDmitry Karpeev    Level: intermediate
14486a4f0f73SDmitry Karpeev 
14496a4f0f73SDmitry Karpeev    Notes: this is merely a convenience subroutine that walks each list,
14506a4f0f73SDmitry Karpeev    destroys each IS on the list, and then frees the list.
1451b1a0a8a3SJed Brown 
1452f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1453b1a0a8a3SJed Brown 
145406b43e7eSDmitry Karpeev .seealso: PCGASMCreateLocalSubdomains(), PCGASMSetSubdomains()
1455b1a0a8a3SJed Brown @*/
14566a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMDestroySubdomains(PetscInt n, IS iis[], IS ois[])
1457b1a0a8a3SJed Brown {
1458b1a0a8a3SJed Brown   PetscInt       i;
1459b1a0a8a3SJed Brown   PetscErrorCode ierr;
1460b1a0a8a3SJed Brown   PetscFunctionBegin;
1461a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1462a35b7b57SDmitry Karpeev   if (iis) {
14636a4f0f73SDmitry Karpeev     PetscValidPointer(iis,2);
1464a35b7b57SDmitry Karpeev     for (i=0; i<n; i++) {
1465a35b7b57SDmitry Karpeev       ierr = ISDestroy(&iis[i]);CHKERRQ(ierr);
1466a35b7b57SDmitry Karpeev     }
14676a4f0f73SDmitry Karpeev     ierr = PetscFree(iis);CHKERRQ(ierr);
1468a35b7b57SDmitry Karpeev   }
14696a4f0f73SDmitry Karpeev   if (ois) {
1470a35b7b57SDmitry Karpeev     for (i=0; i<n; i++) {
1471a35b7b57SDmitry Karpeev       ierr = ISDestroy(&ois[i]);CHKERRQ(ierr);
1472a35b7b57SDmitry Karpeev     }
14736a4f0f73SDmitry Karpeev     ierr = PetscFree(ois);CHKERRQ(ierr);
1474b1a0a8a3SJed Brown   }
1475b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1476b1a0a8a3SJed Brown }
1477b1a0a8a3SJed Brown 
1478af538404SDmitry Karpeev 
1479af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \
1480af538404SDmitry Karpeev {                                                                                                       \
1481af538404SDmitry Karpeev  PetscInt first_row = first/M, last_row = last/M+1;                                                     \
1482af538404SDmitry Karpeev   /*                                                                                                    \
1483af538404SDmitry Karpeev    Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners       \
1484af538404SDmitry Karpeev    of the bounding box of the intersection of the subdomain with the local ownership range (local       \
1485af538404SDmitry Karpeev    subdomain).                                                                                          \
1486af538404SDmitry Karpeev    Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows       \
1487af538404SDmitry Karpeev    of the intersection.                                                                                 \
1488af538404SDmitry Karpeev   */                                                                                                    \
1489af538404SDmitry Karpeev   /* ylow_loc is the grid row containing the first element of the local sumbdomain */                   \
1490eec7e3faSDmitry Karpeev   *ylow_loc = PetscMax(first_row,ylow);                                                                    \
1491af538404SDmitry Karpeev   /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1492eec7e3faSDmitry Karpeev   *xleft_loc = *ylow_loc==first_row?PetscMax(first%M,xleft):xleft;                                                                            \
1493af538404SDmitry Karpeev   /* yhigh_loc is the grid row above the last local subdomain element */                                                                    \
1494eec7e3faSDmitry Karpeev   *yhigh_loc = PetscMin(last_row,yhigh);                                                                                                     \
1495af538404SDmitry Karpeev   /* xright is the offset of the end of the  local subdomain within its grid row (might actually be outside the local subdomain) */         \
1496eec7e3faSDmitry Karpeev   *xright_loc = *yhigh_loc==last_row?PetscMin(xright,last%M):xright;                                                                          \
1497af538404SDmitry Karpeev   /* Now compute the size of the local subdomain n. */ \
1498c3518bceSDmitry Karpeev   *n = 0;                                               \
1499eec7e3faSDmitry Karpeev   if (*ylow_loc < *yhigh_loc) {                           \
1500af538404SDmitry Karpeev     PetscInt width = xright-xleft;                     \
1501eec7e3faSDmitry Karpeev     *n += width*(*yhigh_loc-*ylow_loc-1);                 \
1502eec7e3faSDmitry Karpeev     *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \
1503eec7e3faSDmitry Karpeev     *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \
1504af538404SDmitry Karpeev   }\
1505af538404SDmitry Karpeev }
1506af538404SDmitry Karpeev 
1507c3518bceSDmitry Karpeev 
1508eec7e3faSDmitry Karpeev 
1509eec7e3faSDmitry Karpeev #undef __FUNCT__
1510f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D"
1511b1a0a8a3SJed Brown /*@
1512f746d493SDmitry Karpeev    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1513b1a0a8a3SJed Brown    preconditioner for a two-dimensional problem on a regular grid.
1514b1a0a8a3SJed Brown 
1515af538404SDmitry Karpeev    Collective
1516b1a0a8a3SJed Brown 
1517b1a0a8a3SJed Brown    Input Parameters:
151806b43e7eSDmitry Karpeev +  M, N               - the global number of grid points in the x and y directions
1519af538404SDmitry Karpeev .  Mdomains, Ndomains - the global number of subdomains in the x and y directions
1520b1a0a8a3SJed Brown .  dof                - degrees of freedom per node
1521b1a0a8a3SJed Brown -  overlap            - overlap in mesh lines
1522b1a0a8a3SJed Brown 
1523b1a0a8a3SJed Brown    Output Parameters:
1524af538404SDmitry Karpeev +  Nsub - the number of local subdomains created
15256a4f0f73SDmitry Karpeev .  iis  - array of index sets defining inner (nonoverlapping) subdomains
15266a4f0f73SDmitry Karpeev -  ois  - array of index sets defining outer (overlapping, if overlap > 0) subdomains
1527b1a0a8a3SJed Brown 
1528b1a0a8a3SJed Brown 
1529b1a0a8a3SJed Brown    Level: advanced
1530b1a0a8a3SJed Brown 
1531f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid
1532b1a0a8a3SJed Brown 
15336a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1534f746d493SDmitry Karpeev           PCGASMSetOverlap()
1535b1a0a8a3SJed Brown @*/
15366a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMCreateSubdomains2D(PC pc, PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap, PetscInt *nsub,IS **iis,IS **ois)
1537b1a0a8a3SJed Brown {
1538b1a0a8a3SJed Brown   PetscErrorCode ierr;
15392bbb417fSDmitry Karpeev   PetscMPIInt    size, rank;
15402bbb417fSDmitry Karpeev   PetscInt       i, j;
15412bbb417fSDmitry Karpeev   PetscInt       maxheight, maxwidth;
15422bbb417fSDmitry Karpeev   PetscInt       xstart, xleft, xright, xleft_loc, xright_loc;
15432bbb417fSDmitry Karpeev   PetscInt       ystart, ylow,  yhigh,  ylow_loc,  yhigh_loc;
1544eec7e3faSDmitry Karpeev   PetscInt       x[2][2], y[2][2], n[2];
15452bbb417fSDmitry Karpeev   PetscInt       first, last;
15462bbb417fSDmitry Karpeev   PetscInt       nidx, *idx;
15472bbb417fSDmitry Karpeev   PetscInt       ii,jj,s,q,d;
1548af538404SDmitry Karpeev   PetscInt       k,kk;
15492bbb417fSDmitry Karpeev   PetscMPIInt    color;
15502bbb417fSDmitry Karpeev   MPI_Comm       comm, subcomm;
15516a4f0f73SDmitry Karpeev   IS             **xis = 0, **is = ois, **is_local = iis;
1552d34fcf5fSBarry Smith 
1553b1a0a8a3SJed Brown   PetscFunctionBegin;
15542bbb417fSDmitry Karpeev   ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr);
15552bbb417fSDmitry Karpeev   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
15562bbb417fSDmitry Karpeev   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
15572bbb417fSDmitry Karpeev   ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr);
1558d34fcf5fSBarry 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) "
15592bbb417fSDmitry Karpeev              "does not respect the number of degrees of freedom per grid point %D", first, last, dof);
1560d34fcf5fSBarry Smith 
1561af538404SDmitry Karpeev   /* Determine the number of domains with nonzero intersections with the local ownership range. */
15622bbb417fSDmitry Karpeev   s = 0;
1563af538404SDmitry Karpeev   ystart = 0;
1564af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1565af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1566af538404SDmitry 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);
1567eec7e3faSDmitry Karpeev     /* Vertical domain limits with an overlap. */
1568eec7e3faSDmitry Karpeev     ylow = PetscMax(ystart - overlap,0);
1569eec7e3faSDmitry Karpeev     yhigh = PetscMin(ystart + maxheight + overlap,N);
1570b1a0a8a3SJed Brown     xstart = 0;
1571af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1572af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1573af538404SDmitry 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);
1574eec7e3faSDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1575eec7e3faSDmitry Karpeev       xleft   = PetscMax(xstart - overlap,0);
1576eec7e3faSDmitry Karpeev       xright  = PetscMin(xstart + maxwidth + overlap,M);
1577af538404SDmitry Karpeev       /*
1578af538404SDmitry Karpeev          Determine whether this subdomain intersects this processor's ownership range of pc->pmat.
1579af538404SDmitry Karpeev       */
1580c3518bceSDmitry Karpeev       PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1581af538404SDmitry Karpeev       if (nidx) {
1582af538404SDmitry Karpeev         ++s;
1583af538404SDmitry Karpeev       }
1584af538404SDmitry Karpeev       xstart += maxwidth;
1585af538404SDmitry Karpeev     }/* for (i = 0; i < Mdomains; ++i) */
1586af538404SDmitry Karpeev     ystart += maxheight;
1587af538404SDmitry Karpeev   }/* for (j = 0; j < Ndomains; ++j) */
1588af538404SDmitry Karpeev   /* Now we can allocate the necessary number of ISs. */
1589af538404SDmitry Karpeev   *nsub = s;
1590af538404SDmitry Karpeev   ierr = PetscMalloc((*nsub)*sizeof(IS*),is);CHKERRQ(ierr);
1591af538404SDmitry Karpeev   ierr = PetscMalloc((*nsub)*sizeof(IS*),is_local);CHKERRQ(ierr);
1592af538404SDmitry Karpeev   s = 0;
1593af538404SDmitry Karpeev   ystart = 0;
1594af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1595af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1596af538404SDmitry 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);
1597af538404SDmitry Karpeev     /* Vertical domain limits with an overlap. */
1598eec7e3faSDmitry Karpeev     y[0][0] = PetscMax(ystart - overlap,0);
1599eec7e3faSDmitry Karpeev     y[0][1] = PetscMin(ystart + maxheight + overlap,N);
1600eec7e3faSDmitry Karpeev     /* Vertical domain limits without an overlap. */
1601eec7e3faSDmitry Karpeev     y[1][0] = ystart;
1602eec7e3faSDmitry Karpeev     y[1][1] = PetscMin(ystart + maxheight,N);
1603eec7e3faSDmitry Karpeev     xstart = 0;
1604af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1605af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1606af538404SDmitry 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);
1607af538404SDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1608eec7e3faSDmitry Karpeev       x[0][0]  = PetscMax(xstart - overlap,0);
1609eec7e3faSDmitry Karpeev       x[0][1]  = PetscMin(xstart + maxwidth + overlap,M);
1610eec7e3faSDmitry Karpeev       /* Horizontal domain limits without an overlap. */
1611eec7e3faSDmitry Karpeev       x[1][0] = xstart;
1612eec7e3faSDmitry Karpeev       x[1][1] = PetscMin(xstart+maxwidth,M);
16132bbb417fSDmitry Karpeev       /*
16142bbb417fSDmitry Karpeev          Determine whether this domain intersects this processor's ownership range of pc->pmat.
16152bbb417fSDmitry Karpeev          Do this twice: first for the domains with overlaps, and once without.
16162bbb417fSDmitry Karpeev          During the first pass create the subcommunicators, and use them on the second pass as well.
16172bbb417fSDmitry Karpeev       */
16182bbb417fSDmitry Karpeev       for (q = 0; q < 2; ++q) {
16192bbb417fSDmitry Karpeev         /*
16202bbb417fSDmitry Karpeev           domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
16212bbb417fSDmitry Karpeev           according to whether the domain with an overlap or without is considered.
16222bbb417fSDmitry Karpeev         */
16232bbb417fSDmitry Karpeev         xleft = x[q][0]; xright = x[q][1];
16242bbb417fSDmitry Karpeev         ylow  = y[q][0]; yhigh  = y[q][1];
1625c3518bceSDmitry Karpeev         PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1626af538404SDmitry Karpeev         nidx *= dof;
1627eec7e3faSDmitry Karpeev         n[q] = nidx;
16282bbb417fSDmitry Karpeev         /*
1629eec7e3faSDmitry Karpeev          Based on the counted number of indices in the local domain *with an overlap*,
1630af538404SDmitry Karpeev          construct a subcommunicator of all the processors supporting this domain.
1631eec7e3faSDmitry Karpeev          Observe that a domain with an overlap might have nontrivial local support,
1632eec7e3faSDmitry Karpeev          while the domain without an overlap might not.  Hence, the decision to participate
1633eec7e3faSDmitry Karpeev          in the subcommunicator must be based on the domain with an overlap.
16342bbb417fSDmitry Karpeev          */
16352bbb417fSDmitry Karpeev         if (q == 0) {
16362bbb417fSDmitry Karpeev           if (nidx) {
16372bbb417fSDmitry Karpeev             color = 1;
1638d34fcf5fSBarry Smith           } else {
16392bbb417fSDmitry Karpeev             color = MPI_UNDEFINED;
16402bbb417fSDmitry Karpeev           }
16412bbb417fSDmitry Karpeev           ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr);
16422bbb417fSDmitry Karpeev         }
1643af538404SDmitry Karpeev         /*
1644eec7e3faSDmitry Karpeev          Proceed only if the number of local indices *with an overlap* is nonzero.
1645af538404SDmitry Karpeev          */
1646eec7e3faSDmitry Karpeev         if (n[0]) {
1647af538404SDmitry Karpeev           if (q == 0) {
16486a4f0f73SDmitry Karpeev             xis = is;
1649af538404SDmitry Karpeev           }
1650af538404SDmitry Karpeev           if (q == 1) {
1651af538404SDmitry Karpeev             /*
1652eec7e3faSDmitry Karpeev              The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1653af538404SDmitry Karpeev              Moreover, if the overlap is zero, the two ISs are identical.
1654af538404SDmitry Karpeev              */
1655b1a0a8a3SJed Brown             if (overlap == 0) {
1656eec7e3faSDmitry Karpeev               (*is_local)[s] = (*is)[s];
16572bbb417fSDmitry Karpeev               ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr);
16582bbb417fSDmitry Karpeev               continue;
1659d34fcf5fSBarry Smith             } else {
16606a4f0f73SDmitry Karpeev               xis = is_local;
1661eec7e3faSDmitry Karpeev               subcomm = ((PetscObject)(*is)[s])->comm;
16622bbb417fSDmitry Karpeev             }
1663af538404SDmitry Karpeev           }/* if (q == 1) */
1664eec7e3faSDmitry Karpeev           idx = PETSC_NULL;
16652bbb417fSDmitry Karpeev           ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1666eec7e3faSDmitry Karpeev           if (nidx) {
16672bbb417fSDmitry Karpeev             k    = 0;
16682bbb417fSDmitry Karpeev             for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1669af538404SDmitry Karpeev               PetscInt x0 = (jj==ylow_loc)?xleft_loc:xleft;
1670af538404SDmitry Karpeev               PetscInt x1 = (jj==yhigh_loc-1)?xright_loc:xright;
1671af538404SDmitry Karpeev               kk = dof*(M*jj + x0);
16722bbb417fSDmitry Karpeev               for (ii=x0; ii<x1; ++ii) {
16732bbb417fSDmitry Karpeev                 for (d = 0; d < dof; ++d) {
16742bbb417fSDmitry Karpeev                   idx[k++] = kk++;
1675b1a0a8a3SJed Brown                 }
1676b1a0a8a3SJed Brown               }
1677b1a0a8a3SJed Brown             }
1678eec7e3faSDmitry Karpeev           }
16796a4f0f73SDmitry Karpeev           ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr);
1680eec7e3faSDmitry Karpeev         }/* if (n[0]) */
16812bbb417fSDmitry Karpeev       }/* for (q = 0; q < 2; ++q) */
1682eec7e3faSDmitry Karpeev       if (n[0]) {
16832bbb417fSDmitry Karpeev         ++s;
1684b1a0a8a3SJed Brown       }
1685af538404SDmitry Karpeev       xstart += maxwidth;
1686af538404SDmitry Karpeev     }/* for (i = 0; i < Mdomains; ++i) */
16872bbb417fSDmitry Karpeev     ystart += maxheight;
1688af538404SDmitry Karpeev   }/* for (j = 0; j < Ndomains; ++j) */
1689b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1690b1a0a8a3SJed Brown }
1691b1a0a8a3SJed Brown 
1692b1a0a8a3SJed Brown #undef __FUNCT__
169306b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains"
1694b1a0a8a3SJed Brown /*@C
169506b43e7eSDmitry Karpeev     PCGASMGetSubdomains - Gets the subdomains supported on this processor
169606b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
1697b1a0a8a3SJed Brown 
1698b1a0a8a3SJed Brown     Not Collective
1699b1a0a8a3SJed Brown 
1700b1a0a8a3SJed Brown     Input Parameter:
1701b1a0a8a3SJed Brown .   pc - the preconditioner context
1702b1a0a8a3SJed Brown 
1703b1a0a8a3SJed Brown     Output Parameters:
1704b1a0a8a3SJed Brown +   n   - the number of subdomains for this processor (default value = 1)
17056a4f0f73SDmitry Karpeev .   iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be PETSC_NULL)
17066a4f0f73SDmitry Karpeev -   ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be PETSC_NULL)
1707b1a0a8a3SJed Brown 
1708b1a0a8a3SJed Brown 
1709b1a0a8a3SJed Brown     Notes:
17106a4f0f73SDmitry Karpeev     The user is responsible for destroying the ISs and freeing the returned arrays.
1711b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
1712b1a0a8a3SJed Brown 
1713b1a0a8a3SJed Brown     Level: advanced
1714b1a0a8a3SJed Brown 
171506b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz
1716b1a0a8a3SJed Brown 
17176a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
171806b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubmatrices()
1719b1a0a8a3SJed Brown @*/
17206a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[])
1721b1a0a8a3SJed Brown {
1722f746d493SDmitry Karpeev   PC_GASM         *osm;
1723b1a0a8a3SJed Brown   PetscErrorCode ierr;
1724b1a0a8a3SJed Brown   PetscBool      match;
17256a4f0f73SDmitry Karpeev   PetscInt       i;
1726b1a0a8a3SJed Brown   PetscFunctionBegin;
1727b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1728251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1729f23aa3ddSBarry Smith   if (!match) SETERRQ2(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1730f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1731ab3e923fSDmitry Karpeev   if (n)  *n  = osm->n;
17326a4f0f73SDmitry Karpeev   if (iis) {
17336a4f0f73SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(IS), iis);CHKERRQ(ierr);
17346a4f0f73SDmitry Karpeev   }
17356a4f0f73SDmitry Karpeev   if (ois) {
17366a4f0f73SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(IS), ois);CHKERRQ(ierr);
17376a4f0f73SDmitry Karpeev   }
17386a4f0f73SDmitry Karpeev   if (iis || ois) {
17396a4f0f73SDmitry Karpeev     for (i = 0; i < osm->n; ++i) {
17406a4f0f73SDmitry Karpeev       if (iis) (*iis)[i] = osm->iis[i];
17416a4f0f73SDmitry Karpeev       if (ois) (*ois)[i] = osm->ois[i];
17426a4f0f73SDmitry Karpeev     }
1743b1a0a8a3SJed Brown   }
1744b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1745b1a0a8a3SJed Brown }
1746b1a0a8a3SJed Brown 
1747b1a0a8a3SJed Brown #undef __FUNCT__
174806b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices"
1749b1a0a8a3SJed Brown /*@C
175006b43e7eSDmitry Karpeev     PCGASMGetSubmatrices - Gets the local submatrices (for this processor
1751b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1752b1a0a8a3SJed Brown 
1753b1a0a8a3SJed Brown     Not Collective
1754b1a0a8a3SJed Brown 
1755b1a0a8a3SJed Brown     Input Parameter:
1756b1a0a8a3SJed Brown .   pc - the preconditioner context
1757b1a0a8a3SJed Brown 
1758b1a0a8a3SJed Brown     Output Parameters:
1759b1a0a8a3SJed Brown +   n   - the number of matrices for this processor (default value = 1)
1760b1a0a8a3SJed Brown -   mat - the matrices
1761b1a0a8a3SJed Brown 
176206b43e7eSDmitry Karpeev     Notes: matrices returned by this routine have the same communicators as the index sets (IS)
176306b43e7eSDmitry Karpeev            used to define subdomains in PCGASMSetSubdomains(), or PETSC_COMM_SELF, if the
17646a4f0f73SDmitry Karpeev            subdomains were defined using PCGASMSetTotalSubdomains().
1765b1a0a8a3SJed Brown     Level: advanced
1766b1a0a8a3SJed Brown 
1767f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi
1768b1a0a8a3SJed Brown 
17696a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomain(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
177006b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains()
1771b1a0a8a3SJed Brown @*/
177206b43e7eSDmitry Karpeev PetscErrorCode  PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1773b1a0a8a3SJed Brown {
1774f746d493SDmitry Karpeev   PC_GASM         *osm;
1775b1a0a8a3SJed Brown   PetscErrorCode ierr;
1776b1a0a8a3SJed Brown   PetscBool      match;
1777b1a0a8a3SJed Brown 
1778b1a0a8a3SJed Brown   PetscFunctionBegin;
1779b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1780b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1781b1a0a8a3SJed Brown   if (mat) PetscValidPointer(mat,3);
1782b1a0a8a3SJed Brown   if (!pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1783251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
178406b43e7eSDmitry Karpeev   if (!match) SETERRQ2(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1785f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1786ab3e923fSDmitry Karpeev   if (n)   *n   = osm->n;
1787b1a0a8a3SJed Brown   if (mat) *mat = osm->pmat;
178806b43e7eSDmitry Karpeev 
1789b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1790b1a0a8a3SJed Brown }
1791