xref: /petsc/src/ksp/pc/impls/gasm/gasm.c (revision 6a4f0f735f30779341516db08481ffdda72a83fc)
1b1a0a8a3SJed Brown /*
2f746d493SDmitry Karpeev   This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation.
3*6a4f0f73SDmitry Karpeev   In this version each processor may intersect multiple subdomains and any subdomain may
4*6a4f0f73SDmitry Karpeev   intersect multiple processors.  Intersections of subdomains with processors are called *local
5*6a4f0f73SDmitry Karpeev   subdomains*.
6b1a0a8a3SJed Brown 
7*6a4f0f73SDmitry Karpeev        N    - total number of local subdomains on all processors  (set in PCGASMSetTotalSubdomains() or calculated in PCSetUp_GASM())
8*6a4f0f73SDmitry Karpeev        n    - actual number of local subdomains on this processor (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains())
9*6a4f0f73SDmitry 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. */
19*6a4f0f73SDmitry Karpeev   VecScatter gorestriction;       /* merged restriction to disjoint union of outer subdomains */
20*6a4f0f73SDmitry Karpeev   VecScatter girestriction;       /* merged restriction to disjoint union of inner subdomains */
21*6a4f0f73SDmitry Karpeev   IS         *ois;                /* index sets that define the outer (conceptually, overlapping) subdomains */
22*6a4f0f73SDmitry 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 */
25*6a4f0f73SDmitry 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) */
27*6a4f0f73SDmitry 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);
44*6a4f0f73SDmitry Karpeev   ierr = ISGetLocalSize(osm->ois[i], &nidx);CHKERRQ(ierr);
45eec7e3faSDmitry Karpeev   /*
46eec7e3faSDmitry Karpeev    No more than 15 characters per index plus a space.
47eec7e3faSDmitry Karpeev    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
48eec7e3faSDmitry Karpeev    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
49eec7e3faSDmitry Karpeev    For nidx == 0, the whole string 16 '\0'.
50eec7e3faSDmitry Karpeev    */
51eec7e3faSDmitry Karpeev   ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr);
5206b43e7eSDmitry Karpeev   ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer); CHKERRQ(ierr);
53*6a4f0f73SDmitry Karpeev   ierr = ISGetIndices(osm->ois[i], &idx);CHKERRQ(ierr);
54af538404SDmitry Karpeev   for(j = 0; j < nidx; ++j) {
55af538404SDmitry Karpeev     ierr = PetscViewerStringSPrintf(sviewer,"%D ", idx[j]); CHKERRQ(ierr);
56af538404SDmitry Karpeev   }
576bf464f9SBarry Smith   ierr = PetscViewerDestroy(&sviewer); CHKERRQ(ierr);
58*6a4f0f73SDmitry Karpeev   ierr = ISRestoreIndices(osm->ois[i],&idx);CHKERRQ(ierr);
59*6a4f0f73SDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n");CHKERRQ(ierr);
6006b43e7eSDmitry Karpeev   ierr = PetscViewerFlush(viewer);                                   CHKERRQ(ierr);
617b23a99aSBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
62af538404SDmitry Karpeev   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
63af538404SDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
647b23a99aSBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
65af538404SDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
6606b43e7eSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
67af538404SDmitry Karpeev   ierr = PetscFree(cidx);CHKERRQ(ierr);
68*6a4f0f73SDmitry Karpeev   if (osm->iis) {
69*6a4f0f73SDmitry Karpeev     ierr = ISGetLocalSize(osm->iis[i], &nidx);CHKERRQ(ierr);
70eec7e3faSDmitry Karpeev     /*
71eec7e3faSDmitry Karpeev      No more than 15 characters per index plus a space.
72eec7e3faSDmitry Karpeev      PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
73eec7e3faSDmitry Karpeev      in case nidx == 0. That will take care of the space for the trailing '\0' as well.
74eec7e3faSDmitry Karpeev      For nidx == 0, the whole string 16 '\0'.
75eec7e3faSDmitry Karpeev      */
76eec7e3faSDmitry Karpeev     ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr);
7706b43e7eSDmitry Karpeev     ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer); CHKERRQ(ierr);
78*6a4f0f73SDmitry Karpeev     ierr = ISGetIndices(osm->iis[i], &idx);CHKERRQ(ierr);
79af538404SDmitry Karpeev     for(j = 0; j < nidx; ++j) {
80af538404SDmitry Karpeev       ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr);
81af538404SDmitry Karpeev     }
82*6a4f0f73SDmitry Karpeev     ierr = ISRestoreIndices(osm->iis[i],&idx);CHKERRQ(ierr);
8306b43e7eSDmitry Karpeev     ierr = PetscViewerDestroy(&sviewer); CHKERRQ(ierr);
84*6a4f0f73SDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n");CHKERRQ(ierr);
85*6a4f0f73SDmitry Karpeev     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
86af538404SDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
87*6a4f0f73SDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
88*6a4f0f73SDmitry Karpeev     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
89af538404SDmitry Karpeev     ierr = PetscFree(cidx);CHKERRQ(ierr);
90af538404SDmitry Karpeev   }
9106b43e7eSDmitry Karpeev   PetscFunctionReturn(0);
9206b43e7eSDmitry Karpeev }
9306b43e7eSDmitry Karpeev 
9406b43e7eSDmitry Karpeev #undef __FUNCT__
9506b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains"
9606b43e7eSDmitry Karpeev static PetscErrorCode  PCGASMPrintSubdomains(PC pc)
9706b43e7eSDmitry Karpeev {
9806b43e7eSDmitry Karpeev   PC_GASM        *osm  = (PC_GASM*)pc->data;
9906b43e7eSDmitry Karpeev   const char     *prefix;
10006b43e7eSDmitry Karpeev   char           fname[PETSC_MAX_PATH_LEN+1];
10106b43e7eSDmitry Karpeev   PetscViewer    viewer;
10206b43e7eSDmitry Karpeev   PetscInt       i;
10306b43e7eSDmitry Karpeev   PetscBool      found;
10406b43e7eSDmitry Karpeev   PetscErrorCode ierr;
10506b43e7eSDmitry Karpeev 
10606b43e7eSDmitry Karpeev   PetscFunctionBegin;
10706b43e7eSDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
10806b43e7eSDmitry Karpeev   ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr);
10906b43e7eSDmitry Karpeev   if (!found) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
11006b43e7eSDmitry Karpeev   for (i=0;i<osm->n;++i) {
111*6a4f0f73SDmitry Karpeev     ierr = PetscViewerASCIIOpen(((PetscObject)((osm->ois)[i]))->comm,fname,&viewer);CHKERRQ(ierr);
11206b43e7eSDmitry Karpeev     ierr = PCGASMSubdomainView_Private(pc, i, viewer); CHKERRQ(ierr);
113af538404SDmitry Karpeev   }
114af538404SDmitry Karpeev   PetscFunctionReturn(0);
115af538404SDmitry Karpeev }
116af538404SDmitry Karpeev 
117af538404SDmitry Karpeev 
118af538404SDmitry Karpeev #undef __FUNCT__
119f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM"
120f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer)
121b1a0a8a3SJed Brown {
122f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
123af538404SDmitry Karpeev   const char     *prefix;
124b1a0a8a3SJed Brown   PetscErrorCode ierr;
125af538404SDmitry Karpeev   PetscMPIInt    rank, size;
126b1a0a8a3SJed Brown   PetscInt       i,bsz;
12706b43e7eSDmitry Karpeev   PetscBool      iascii,view_subdomains=PETSC_FALSE;
128b1a0a8a3SJed Brown   PetscViewer    sviewer;
12906b43e7eSDmitry Karpeev   PetscInt       count, l, gcount, *numbering, *permutation;
130*6a4f0f73SDmitry Karpeev   char overlap[256]     = "user-defined overlap";
131*6a4f0f73SDmitry Karpeev   char gsubdomains[256] = "unknown total number of subdomains";
13206b43e7eSDmitry Karpeev   char lsubdomains[256] = "unknown number of local  subdomains";
13306b43e7eSDmitry Karpeev   char msubdomains[256] = "unknown max number of local subdomains";
134b1a0a8a3SJed Brown   PetscFunctionBegin;
135af538404SDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm, &size);CHKERRQ(ierr);
136af538404SDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm, &rank);CHKERRQ(ierr);
13706b43e7eSDmitry Karpeev 
13806b43e7eSDmitry Karpeev 
13906b43e7eSDmitry Karpeev   ierr = PetscMalloc2(osm->n, PetscInt, &permutation, osm->n, PetscInt, &numbering); CHKERRQ(ierr);
14006b43e7eSDmitry Karpeev   for(i = 0; i < osm->n; ++i) permutation[i] = i;
141*6a4f0f73SDmitry Karpeev   ierr = PetscObjectsGetGlobalNumbering(((PetscObject)pc)->comm, osm->n, (PetscObject*)osm->ois, &gcount, numbering); CHKERRQ(ierr);
14206b43e7eSDmitry Karpeev   ierr = PetscSortIntWithPermutation(osm->n, numbering, permutation); CHKERRQ(ierr);
14306b43e7eSDmitry Karpeev 
14406b43e7eSDmitry Karpeev   if(osm->overlap >= 0) {
14506b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %D",osm->overlap);        CHKERRQ(ierr);
14606b43e7eSDmitry Karpeev   }
147*6a4f0f73SDmitry Karpeev   ierr = PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %D",gcount);      CHKERRQ(ierr);
14806b43e7eSDmitry Karpeev   if(osm->N > 0) {
14906b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(lsubdomains, sizeof(gsubdomains), "number of local subdomains = %D",osm->N);     CHKERRQ(ierr);
15006b43e7eSDmitry Karpeev   }
15106b43e7eSDmitry Karpeev   if(osm->nmax > 0){
15206b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr);
15306b43e7eSDmitry Karpeev   }
15406b43e7eSDmitry Karpeev 
155af538404SDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
15606b43e7eSDmitry Karpeev   ierr = PetscOptionsGetBool(prefix,"-pc_gasm_view_subdomains",&view_subdomains,PETSC_NULL);CHKERRQ(ierr);
15706b43e7eSDmitry Karpeev 
15806b43e7eSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii); CHKERRQ(ierr);
159b1a0a8a3SJed Brown   if (iascii) {
16006b43e7eSDmitry Karpeev     /*
16106b43e7eSDmitry Karpeev      Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer:
16206b43e7eSDmitry Karpeev      the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed
16306b43e7eSDmitry Karpeev      collectively on the comm.
16406b43e7eSDmitry Karpeev      */
16506b43e7eSDmitry Karpeev     ierr = PetscObjectName((PetscObject)viewer);                  CHKERRQ(ierr);
16606b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Generalized additive Schwarz:\n"); CHKERRQ(ierr);
16706b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Restriction/interpolation type: %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr);
16806b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",overlap);    CHKERRQ(ierr);
16906b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",gsubdomains);CHKERRQ(ierr);
17006b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",lsubdomains);CHKERRQ(ierr);
17106b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",msubdomains);CHKERRQ(ierr);
1727b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
173eec7e3faSDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d:%d] number of locally-supported subdomains = %D\n",(int)rank,(int)size,osm->n);CHKERRQ(ierr);
174af538404SDmitry Karpeev     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1757b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
176*6a4f0f73SDmitry Karpeev     /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */
17706b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Subdomain solver info is as follows:\n");CHKERRQ(ierr);
178b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
179b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
18006b43e7eSDmitry Karpeev     /* Make sure that everybody waits for the banner to be printed. */
18106b43e7eSDmitry Karpeev     ierr = MPI_Barrier(((PetscObject)viewer)->comm); CHKERRQ(ierr);
18206b43e7eSDmitry Karpeev     l = 0;
18306b43e7eSDmitry Karpeev     for(count = 0; count < gcount; ++count) {
18406b43e7eSDmitry Karpeev       /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
18506b43e7eSDmitry Karpeev       PetscMPIInt srank, ssize;
18606b43e7eSDmitry Karpeev       if(l<osm->n){
18706b43e7eSDmitry Karpeev         PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
18806b43e7eSDmitry Karpeev         if(numbering[d] == count) {
189*6a4f0f73SDmitry Karpeev           ierr = MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize); CHKERRQ(ierr);
190*6a4f0f73SDmitry Karpeev           ierr = MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank); CHKERRQ(ierr);
191*6a4f0f73SDmitry Karpeev           ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
192*6a4f0f73SDmitry Karpeev           ierr = ISGetLocalSize(osm->ois[d],&bsz);CHKERRQ(ierr);
1937b23a99aSBarry Smith           ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_TRUE);CHKERRQ(ierr);
19406b43e7eSDmitry 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);
195*6a4f0f73SDmitry Karpeev           ierr = PetscViewerFlush(sviewer); CHKERRQ(ierr);
196*6a4f0f73SDmitry Karpeev           ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_FALSE);CHKERRQ(ierr);
19706b43e7eSDmitry Karpeev           if(view_subdomains) {
19806b43e7eSDmitry Karpeev             ierr = PCGASMSubdomainView_Private(pc,d,sviewer); CHKERRQ(ierr);
19906b43e7eSDmitry Karpeev           }
200*6a4f0f73SDmitry Karpeev           if(!pc->setupcalled) {
201*6a4f0f73SDmitry Karpeev             PetscViewerASCIIPrintf(sviewer, "Solver not set up yet: PCSetUp() not yet called\n"); CHKERRQ(ierr);
202*6a4f0f73SDmitry Karpeev           }
203*6a4f0f73SDmitry Karpeev           else {
20406b43e7eSDmitry Karpeev             ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr);
205*6a4f0f73SDmitry Karpeev           }
20606b43e7eSDmitry Karpeev           ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
2077b23a99aSBarry Smith           ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
208*6a4f0f73SDmitry Karpeev           ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr);
20906b43e7eSDmitry Karpeev           ++l;
210b1a0a8a3SJed Brown         }
21106b43e7eSDmitry Karpeev       }
21206b43e7eSDmitry Karpeev       ierr = MPI_Barrier(((PetscObject)pc)->comm); CHKERRQ(ierr);
213b1a0a8a3SJed Brown     }
214b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
215b1a0a8a3SJed Brown   } else {
216f746d493SDmitry Karpeev     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCGASM",((PetscObject)viewer)->type_name);
217b1a0a8a3SJed Brown   }
21806b43e7eSDmitry Karpeev   ierr = PetscFree2(permutation,numbering); CHKERRQ(ierr);
219b1a0a8a3SJed Brown   PetscFunctionReturn(0);
220b1a0a8a3SJed Brown }
221b1a0a8a3SJed Brown 
222af538404SDmitry Karpeev 
223af538404SDmitry Karpeev 
224af538404SDmitry Karpeev 
225af538404SDmitry Karpeev 
226b1a0a8a3SJed Brown #undef __FUNCT__
227f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM"
228f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc)
229b1a0a8a3SJed Brown {
230f746d493SDmitry Karpeev   PC_GASM         *osm  = (PC_GASM*)pc->data;
231b1a0a8a3SJed Brown   PetscErrorCode ierr;
232b1a0a8a3SJed Brown   PetscBool      symset,flg;
23388389755SJed Brown   PetscInt       i;
234c10bc1a1SDmitry Karpeev   PetscMPIInt    rank, size;
235b1a0a8a3SJed Brown   MatReuse       scall = MAT_REUSE_MATRIX;
236b1a0a8a3SJed Brown   KSP            ksp;
237b1a0a8a3SJed Brown   PC             subpc;
238b1a0a8a3SJed Brown   const char     *prefix,*pprefix;
239f746d493SDmitry Karpeev   Vec            x,y;
240*6a4f0f73SDmitry Karpeev   PetscInt       oni;       /* Number of indices in the i-th local outer subdomain.               */
241*6a4f0f73SDmitry Karpeev   const PetscInt *oidxi;    /* Indices from the i-th subdomain local outer subdomain.             */
242*6a4f0f73SDmitry Karpeev   PetscInt       on;        /* Number of indices in the disjoint union of local outer subdomains. */
243*6a4f0f73SDmitry Karpeev   PetscInt       *oidx;     /* Indices in the disjoint union of local outer subdomains. */
244*6a4f0f73SDmitry Karpeev   IS             gois;      /* Disjoint union the global indices of outer subdomains.             */
245*6a4f0f73SDmitry Karpeev   IS             goid;      /* Identity IS of the size of the disjoint union of outer subdomains. */
246f746d493SDmitry Karpeev   PetscScalar    *gxarray, *gyarray;
247*6a4f0f73SDmitry Karpeev   PetscInt       gofirst;   /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy --
248*6a4f0f73SDmitry Karpeev                              over the disjoint union of outer subdomains. */
249fdc48646SDmitry Karpeev   DM             *domain_dm = PETSC_NULL;
250b1a0a8a3SJed Brown 
251b1a0a8a3SJed Brown   PetscFunctionBegin;
252c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
253c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
254b1a0a8a3SJed Brown   if (!pc->setupcalled) {
255b1a0a8a3SJed Brown 
256b1a0a8a3SJed Brown     if (!osm->type_set) {
257b1a0a8a3SJed Brown       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
258f746d493SDmitry Karpeev       if (symset && flg) { osm->type = PC_GASM_BASIC; }
259b1a0a8a3SJed Brown     }
260b1a0a8a3SJed Brown 
26106b43e7eSDmitry Karpeev     /*
26206b43e7eSDmitry Karpeev      If subdomains have been set, then the local number of subdomains, osm->n, is NOT PETSC_DECIDE and is at least 1.
26306b43e7eSDmitry 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.
26406b43e7eSDmitry Karpeev      */
26506b43e7eSDmitry Karpeev     if (osm->n == PETSC_DECIDE) {
26669ca1f37SDmitry Karpeev       /* no subdomains given */
26769ca1f37SDmitry Karpeev       /* try pc->dm first */
268fdc48646SDmitry Karpeev       if(pc->dm) {
269fdc48646SDmitry Karpeev         char      ddm_name[1024];
270fdc48646SDmitry Karpeev         DM        ddm;
271fdc48646SDmitry Karpeev         PetscBool flg;
272fdc48646SDmitry Karpeev         PetscInt     num_domains, d;
273fdc48646SDmitry Karpeev         char         **domain_names;
274fdc48646SDmitry Karpeev         IS           *domain_is;
275fdc48646SDmitry Karpeev         /* Allow the user to request a decomposition DM by name */
276fdc48646SDmitry Karpeev         ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr);
27706b43e7eSDmitry Karpeev         ierr = PetscOptionsString("-pc_gasm_decomposition","Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr);
278fdc48646SDmitry Karpeev         if(flg) {
279fdc48646SDmitry Karpeev           ierr = DMCreateDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr);
280fdc48646SDmitry Karpeev           if(!ddm) {
281fdc48646SDmitry Karpeev             SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown DM decomposition name %s", ddm_name);
282fdc48646SDmitry Karpeev           }
283fdc48646SDmitry Karpeev           ierr = PetscInfo(pc,"Using decomposition DM defined using options database\n");CHKERRQ(ierr);
284fdc48646SDmitry Karpeev           ierr = PCSetDM(pc,ddm); CHKERRQ(ierr);
285fdc48646SDmitry Karpeev         }
286fdc48646SDmitry Karpeev         ierr = DMCreateDecomposition(pc->dm, &num_domains, &domain_names, &domain_is, &domain_dm);    CHKERRQ(ierr);
28769ca1f37SDmitry Karpeev         if(num_domains) {
28806b43e7eSDmitry Karpeev           ierr = PCGASMSetSubdomains(pc, num_domains, domain_is, PETSC_NULL);CHKERRQ(ierr);
28969ca1f37SDmitry Karpeev         }
290fdc48646SDmitry Karpeev         for(d = 0; d < num_domains; ++d) {
291fdc48646SDmitry Karpeev           ierr = PetscFree(domain_names[d]); CHKERRQ(ierr);
292fdc48646SDmitry Karpeev           ierr = ISDestroy(&domain_is[d]);   CHKERRQ(ierr);
293fdc48646SDmitry Karpeev         }
294fdc48646SDmitry Karpeev         ierr = PetscFree(domain_names);CHKERRQ(ierr);
295fdc48646SDmitry Karpeev         ierr = PetscFree(domain_is);CHKERRQ(ierr);
296fdc48646SDmitry Karpeev       }
29706b43e7eSDmitry Karpeev       if (osm->n == PETSC_DECIDE) { /* still no subdomains; use one per processor */
29844fe18e5SDmitry Karpeev         osm->nmax = osm->n = 1;
299b1a0a8a3SJed Brown         ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
300f746d493SDmitry Karpeev         osm->N = size;
301fdc48646SDmitry Karpeev       }
30206b43e7eSDmitry Karpeev     }
30306b43e7eSDmitry Karpeev     if (osm->N == PETSC_DECIDE) {
304b1a0a8a3SJed Brown       PetscInt inwork[2], outwork[2];
305f746d493SDmitry Karpeev       /* determine global number of subdomains and the max number of local subdomains */
306f746d493SDmitry Karpeev       inwork[0] = inwork[1] = osm->n;
307b1a0a8a3SJed Brown       ierr = MPI_Allreduce(inwork,outwork,1,MPIU_2INT,PetscMaxSum_Op,((PetscObject)pc)->comm);CHKERRQ(ierr);
308f746d493SDmitry Karpeev       osm->nmax = outwork[0];
309f746d493SDmitry Karpeev       osm->N    = outwork[1];
310b1a0a8a3SJed Brown     }
311*6a4f0f73SDmitry Karpeev     if (!osm->iis){
31293c1ec32SDmitry Karpeev       /*
313*6a4f0f73SDmitry Karpeev        The local number of subdomains was set in PCGASMSetTotalSubdomains() or PCGASMSetSubdomains(),
314*6a4f0f73SDmitry Karpeev        but the actual subdomains have not been supplied (in PCGASMSetSubdomains()).
315*6a4f0f73SDmitry Karpeev        We create the requisite number of inner subdomains on PETSC_COMM_SELF (for now).
31693c1ec32SDmitry Karpeev        */
317*6a4f0f73SDmitry Karpeev       ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->overlap,osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr);
318b1a0a8a3SJed Brown     }
319*6a4f0f73SDmitry Karpeev 
320b1a0a8a3SJed Brown     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
321b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
322f746d493SDmitry Karpeev     ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&flg,PETSC_NULL);CHKERRQ(ierr);
323f746d493SDmitry Karpeev     if (flg) { ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); }
324b1a0a8a3SJed Brown 
325b1a0a8a3SJed Brown     if (osm->sort_indices) {
326f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
327*6a4f0f73SDmitry Karpeev         ierr = ISSort(osm->ois[i]);CHKERRQ(ierr);
328*6a4f0f73SDmitry Karpeev 	ierr = ISSort(osm->iis[i]);CHKERRQ(ierr);
329b1a0a8a3SJed Brown       }
330b1a0a8a3SJed Brown     }
331*6a4f0f73SDmitry Karpeev     /*
332*6a4f0f73SDmitry Karpeev      Merge the ISs, create merged vectors and restrictions.
333*6a4f0f73SDmitry Karpeev      */
334*6a4f0f73SDmitry Karpeev     /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */
335*6a4f0f73SDmitry Karpeev     on = 0;
336f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
337*6a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
338*6a4f0f73SDmitry Karpeev       on += oni;
339f746d493SDmitry Karpeev     }
340*6a4f0f73SDmitry Karpeev     ierr = PetscMalloc(on*sizeof(PetscInt), &oidx);CHKERRQ(ierr);
341*6a4f0f73SDmitry Karpeev     on = 0;
342f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
343*6a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
344*6a4f0f73SDmitry Karpeev       ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr);
345*6a4f0f73SDmitry Karpeev       ierr = PetscMemcpy(oidx+on, oidxi, sizeof(PetscInt)*oni);CHKERRQ(ierr);
346*6a4f0f73SDmitry Karpeev       ierr = ISRestoreIndices(osm->ois[i], &oidxi);CHKERRQ(ierr);
347*6a4f0f73SDmitry Karpeev       on += oni;
348f746d493SDmitry Karpeev     }
349*6a4f0f73SDmitry Karpeev     ierr = ISCreateGeneral(((PetscObject)(pc))->comm, on, oidx, PETSC_OWN_POINTER, &gois);CHKERRQ(ierr);
350f746d493SDmitry Karpeev     ierr = MatGetVecs(pc->pmat,&x,&y);CHKERRQ(ierr);
351*6a4f0f73SDmitry Karpeev     ierr = VecCreateMPI(((PetscObject)pc)->comm, on, PETSC_DECIDE, &osm->gx);CHKERRQ(ierr);
352f746d493SDmitry Karpeev     ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr);
353*6a4f0f73SDmitry Karpeev     ierr = VecGetOwnershipRange(osm->gx, &gofirst, PETSC_NULL);CHKERRQ(ierr);
354*6a4f0f73SDmitry Karpeev     ierr = ISCreateStride(((PetscObject)pc)->comm,on,gofirst,1, &goid);CHKERRQ(ierr);
355*6a4f0f73SDmitry Karpeev     ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr);
356*6a4f0f73SDmitry Karpeev     ierr = VecDestroy(&x); CHKERRQ(ierr);
357*6a4f0f73SDmitry Karpeev     /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */
358*6a4f0f73SDmitry Karpeev     { PetscInt       ini;     /* Number of indices the i-th a local inner subdomain. */
359*6a4f0f73SDmitry Karpeev       PetscInt       in;      /* Number of indices in the disjoint uniont of local inner subdomains. */
360*6a4f0f73SDmitry Karpeev       PetscInt       *iidx;   /* Global indices in the merged local inner subdomain. */
361*6a4f0f73SDmitry Karpeev       PetscInt       *ioidx;  /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
362*6a4f0f73SDmitry Karpeev       IS             giis;    /* IS for the disjoint union of inner subdomains. */
363*6a4f0f73SDmitry Karpeev       IS             giois;   /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */
364f746d493SDmitry Karpeev       /**/
365*6a4f0f73SDmitry Karpeev       in = 0;
366f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
367*6a4f0f73SDmitry Karpeev 	ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr);
368*6a4f0f73SDmitry Karpeev 	in += ini;
369f746d493SDmitry Karpeev       }
370*6a4f0f73SDmitry Karpeev       ierr = PetscMalloc(in*sizeof(PetscInt), &iidx); CHKERRQ(ierr);
371*6a4f0f73SDmitry Karpeev       ierr = PetscMalloc(in*sizeof(PetscInt), &ioidx);CHKERRQ(ierr);
372*6a4f0f73SDmitry Karpeev       ierr = VecGetOwnershipRange(osm->gx,&gofirst, PETSC_NULL); CHKERRQ(ierr);
373*6a4f0f73SDmitry Karpeev       in = 0;
374*6a4f0f73SDmitry Karpeev       on = 0;
375f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
376*6a4f0f73SDmitry Karpeev         const PetscInt *iidxi;        /* Global indices of the i-th local inner subdomain. */
377*6a4f0f73SDmitry Karpeev         ISLocalToGlobalMapping ltogi; /* Map from global to local indices of the i-th outer local subdomain. */
378*6a4f0f73SDmitry Karpeev         PetscInt       *ioidxi;       /* Local indices of the i-th local inner subdomain within the local outer subdomain. */
379*6a4f0f73SDmitry Karpeev         PetscInt       ioni;          /* Number of indices in ioidxi; if ioni != ini the inner subdomain is not a subdomain of the outer subdomain (error). */
380*6a4f0f73SDmitry Karpeev         PetscInt       k;
381*6a4f0f73SDmitry Karpeev 	ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr);
382*6a4f0f73SDmitry Karpeev 	ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
383*6a4f0f73SDmitry Karpeev 	ierr = ISGetIndices(osm->iis[i],&iidxi);CHKERRQ(ierr);
384*6a4f0f73SDmitry Karpeev 	ierr = PetscMemcpy(iidx+in, iidxi, sizeof(PetscInt)*ini);CHKERRQ(ierr);
385*6a4f0f73SDmitry Karpeev         ierr = ISLocalToGlobalMappingCreateIS(osm->ois[i],&ltogi);CHKERRQ(ierr);
386*6a4f0f73SDmitry Karpeev         ioidxi = ioidx+in;
387*6a4f0f73SDmitry Karpeev         ierr = ISGlobalToLocalMappingApply(ltogi,IS_GTOLM_DROP,ini,iidxi,&ioni,ioidxi);CHKERRQ(ierr);
388*6a4f0f73SDmitry Karpeev         ierr = ISLocalToGlobalMappingDestroy(&ltogi);CHKERRQ(ierr);
389*6a4f0f73SDmitry Karpeev 	ierr = ISRestoreIndices(osm->iis[i], &iidxi);CHKERRQ(ierr);
390*6a4f0f73SDmitry 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);
391*6a4f0f73SDmitry Karpeev         for(k = 0; k < ini; ++k) {
392*6a4f0f73SDmitry Karpeev           ioidxi[k] += gofirst+on;
393f746d493SDmitry Karpeev         }
394*6a4f0f73SDmitry Karpeev 	in += ini;
395*6a4f0f73SDmitry Karpeev         on += oni;
396f746d493SDmitry Karpeev       }
397*6a4f0f73SDmitry Karpeev       ierr = ISCreateGeneral(((PetscObject)pc)->comm, in, iidx,  PETSC_OWN_POINTER, &giis);CHKERRQ(ierr);
398*6a4f0f73SDmitry Karpeev       ierr = ISCreateGeneral(((PetscObject)pc)->comm, in, ioidx, PETSC_OWN_POINTER, &giois);CHKERRQ(ierr);
399*6a4f0f73SDmitry Karpeev       ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr);
400*6a4f0f73SDmitry Karpeev       ierr = VecDestroy(&y); CHKERRQ(ierr);
401*6a4f0f73SDmitry Karpeev       ierr = ISDestroy(&giis);  CHKERRQ(ierr);
402*6a4f0f73SDmitry Karpeev       ierr = ISDestroy(&giois); CHKERRQ(ierr);
403b1a0a8a3SJed Brown     }
404*6a4f0f73SDmitry Karpeev     ierr = ISDestroy(&goid);CHKERRQ(ierr);
405f746d493SDmitry Karpeev     /* Create the subdomain work vectors. */
406f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->x);CHKERRQ(ierr);
407f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->y);CHKERRQ(ierr);
408f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr);
409f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr);
410*6a4f0f73SDmitry Karpeev     for (i=0, on=0; i<osm->n; ++i, on += oni) {
411*6a4f0f73SDmitry Karpeev       PetscInt oNi;
412*6a4f0f73SDmitry Karpeev       ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr);
413*6a4f0f73SDmitry Karpeev       ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr);
414*6a4f0f73SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr);
415*6a4f0f73SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr);
416b1a0a8a3SJed Brown     }
417f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr);
418f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr);
419b1a0a8a3SJed Brown     /* Create the local solvers */
420f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(KSP *),&osm->ksp);CHKERRQ(ierr);
42144fe18e5SDmitry Karpeev     for (i=0; i<osm->n; i++) { /* KSPs are local */
422*6a4f0f73SDmitry Karpeev       ierr = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr);
423b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr);
424b1a0a8a3SJed Brown       ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
425b1a0a8a3SJed Brown       ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
426b1a0a8a3SJed Brown       ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
427b1a0a8a3SJed Brown       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
428b1a0a8a3SJed Brown       ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
429b1a0a8a3SJed Brown       ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
430b1a0a8a3SJed Brown       osm->ksp[i] = ksp;
431b1a0a8a3SJed Brown     }
432b1a0a8a3SJed Brown     scall = MAT_INITIAL_MATRIX;
433b1a0a8a3SJed Brown 
434*6a4f0f73SDmitry Karpeev   }/*if(!pc->setupcalled)*/
435*6a4f0f73SDmitry Karpeev   else {
436b1a0a8a3SJed Brown     /*
437b1a0a8a3SJed Brown        Destroy the blocks from the previous iteration
438b1a0a8a3SJed Brown     */
439b1a0a8a3SJed Brown     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
440f746d493SDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
441b1a0a8a3SJed Brown       scall = MAT_INITIAL_MATRIX;
442b1a0a8a3SJed Brown     }
443b1a0a8a3SJed Brown   }
444b1a0a8a3SJed Brown 
445b1a0a8a3SJed Brown   /*
446f746d493SDmitry Karpeev      Extract out the submatrices.
447b1a0a8a3SJed Brown   */
44881a5c4aaSDmitry Karpeev   if(size > 1) {
449*6a4f0f73SDmitry Karpeev     ierr = MatGetSubMatricesParallel(pc->pmat,osm->n,osm->ois, osm->ois,scall,&osm->pmat);CHKERRQ(ierr);
45081a5c4aaSDmitry Karpeev   }
45181a5c4aaSDmitry Karpeev   else {
452*6a4f0f73SDmitry Karpeev     ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->ois, osm->ois,scall,&osm->pmat);CHKERRQ(ierr);
45381a5c4aaSDmitry Karpeev   }
454b1a0a8a3SJed Brown   if (scall == MAT_INITIAL_MATRIX) {
455b1a0a8a3SJed Brown     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
456f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
457b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,osm->pmat[i]);CHKERRQ(ierr);
458b1a0a8a3SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
459b1a0a8a3SJed Brown     }
460b1a0a8a3SJed Brown   }
461b1a0a8a3SJed Brown 
462b1a0a8a3SJed Brown   /* Return control to the user so that the submatrices can be modified (e.g., to apply
463b1a0a8a3SJed Brown      different boundary conditions for the submatrices than for the global problem) */
464*6a4f0f73SDmitry Karpeev   ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
465b1a0a8a3SJed Brown 
466b1a0a8a3SJed Brown   /*
467*6a4f0f73SDmitry Karpeev      Loop over submatrices putting them into local ksps
468b1a0a8a3SJed Brown   */
469f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
470b1a0a8a3SJed Brown     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr);
471b1a0a8a3SJed Brown     if (!pc->setupcalled) {
472b1a0a8a3SJed Brown       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
473b1a0a8a3SJed Brown     }
474b1a0a8a3SJed Brown   }
475b1a0a8a3SJed Brown 
476b1a0a8a3SJed Brown   PetscFunctionReturn(0);
477b1a0a8a3SJed Brown }
478b1a0a8a3SJed Brown 
479b1a0a8a3SJed Brown #undef __FUNCT__
480f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM"
481f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
482b1a0a8a3SJed Brown {
483f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
484b1a0a8a3SJed Brown   PetscErrorCode ierr;
485b1a0a8a3SJed Brown   PetscInt       i;
486b1a0a8a3SJed Brown 
487b1a0a8a3SJed Brown   PetscFunctionBegin;
488f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
489b1a0a8a3SJed Brown     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
490b1a0a8a3SJed Brown   }
491b1a0a8a3SJed Brown   PetscFunctionReturn(0);
492b1a0a8a3SJed Brown }
493b1a0a8a3SJed Brown 
494b1a0a8a3SJed Brown #undef __FUNCT__
495f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM"
496f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y)
497b1a0a8a3SJed Brown {
498f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
499b1a0a8a3SJed Brown   PetscErrorCode ierr;
500f746d493SDmitry Karpeev   PetscInt       i;
501b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
502b1a0a8a3SJed Brown 
503b1a0a8a3SJed Brown   PetscFunctionBegin;
504b1a0a8a3SJed Brown   /*
505*6a4f0f73SDmitry Karpeev      Support for limiting the restriction or interpolation only to the inner
506b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
507b1a0a8a3SJed Brown   */
508f746d493SDmitry Karpeev   if(!(osm->type & PC_GASM_RESTRICT)) {
509b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
510f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
511*6a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
512b1a0a8a3SJed Brown   }
513*6a4f0f73SDmitry Karpeev   else {
514*6a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
515b1a0a8a3SJed Brown   }
516*6a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr);
517*6a4f0f73SDmitry Karpeev   if(!(osm->type & PC_GASM_RESTRICT)) {
518*6a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
519*6a4f0f73SDmitry Karpeev   }
520*6a4f0f73SDmitry Karpeev   else {
521*6a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
522*6a4f0f73SDmitry Karpeev   }
52386b96d36SDmitry Karpeev   /* do the subdomain solves */
52486b96d36SDmitry Karpeev   for (i=0; i<osm->n; ++i) {
525b1a0a8a3SJed Brown     ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
526b1a0a8a3SJed Brown   }
527*6a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(y); CHKERRQ(ierr);
528*6a4f0f73SDmitry Karpeev   if(!(osm->type & PC_GASM_INTERPOLATE)) {
529*6a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
530*6a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);  PetscFunctionReturn(0);
531*6a4f0f73SDmitry Karpeev   }
532*6a4f0f73SDmitry Karpeev   else {
533*6a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
534*6a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);  PetscFunctionReturn(0);
535*6a4f0f73SDmitry Karpeev   }
536b1a0a8a3SJed Brown }
537b1a0a8a3SJed Brown 
538b1a0a8a3SJed Brown #undef __FUNCT__
539f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM"
540f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y)
541b1a0a8a3SJed Brown {
542f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
543b1a0a8a3SJed Brown   PetscErrorCode ierr;
544f746d493SDmitry Karpeev   PetscInt       i;
545b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
546b1a0a8a3SJed Brown 
547b1a0a8a3SJed Brown   PetscFunctionBegin;
548b1a0a8a3SJed Brown   /*
549b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
550b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
551b1a0a8a3SJed Brown 
552f746d493SDmitry Karpeev      Note: these are reversed from the PCApply_GASM() because we are applying the
553b1a0a8a3SJed Brown      transpose of the three terms
554b1a0a8a3SJed Brown   */
555f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
556b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
557f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
558*6a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
559b1a0a8a3SJed Brown   }
560*6a4f0f73SDmitry Karpeev   else {
561*6a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
562b1a0a8a3SJed Brown   }
563*6a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr);
564*6a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
565*6a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
566*6a4f0f73SDmitry Karpeev   }
567*6a4f0f73SDmitry Karpeev   else {
568*6a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
569*6a4f0f73SDmitry Karpeev   }
570b1a0a8a3SJed Brown   /* do the local solves */
571ab3e923fSDmitry 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. */
572b1a0a8a3SJed Brown     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
573b1a0a8a3SJed Brown   }
574*6a4f0f73SDmitry Karpeev   ierr = VecZeroEntries(y); CHKERRQ(ierr);
575*6a4f0f73SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
576*6a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
577*6a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
578*6a4f0f73SDmitry Karpeev   }
579*6a4f0f73SDmitry Karpeev   else {
580*6a4f0f73SDmitry Karpeev     ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
581*6a4f0f73SDmitry Karpeev     ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
582*6a4f0f73SDmitry Karpeev   }
583*6a4f0f73SDmitry Karpeev 
584b1a0a8a3SJed Brown   PetscFunctionReturn(0);
585b1a0a8a3SJed Brown }
586b1a0a8a3SJed Brown 
587b1a0a8a3SJed Brown #undef __FUNCT__
588a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM"
589a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc)
590b1a0a8a3SJed Brown {
591f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
592b1a0a8a3SJed Brown   PetscErrorCode ierr;
593b1a0a8a3SJed Brown   PetscInt       i;
594b1a0a8a3SJed Brown 
595b1a0a8a3SJed Brown   PetscFunctionBegin;
596b1a0a8a3SJed Brown   if (osm->ksp) {
597f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
598a06653b4SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
599b1a0a8a3SJed Brown     }
600b1a0a8a3SJed Brown   }
601b1a0a8a3SJed Brown   if (osm->pmat) {
602f746d493SDmitry Karpeev     if (osm->n > 0) {
603ab3e923fSDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
604b1a0a8a3SJed Brown     }
605b1a0a8a3SJed Brown   }
606d34fcf5fSBarry Smith   if (osm->x) {
607f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
608fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
609fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
610b1a0a8a3SJed Brown     }
611d34fcf5fSBarry Smith   }
6126bf464f9SBarry Smith   ierr = VecDestroy(&osm->gx);CHKERRQ(ierr);
6136bf464f9SBarry Smith   ierr = VecDestroy(&osm->gy);CHKERRQ(ierr);
614ab3e923fSDmitry Karpeev 
615*6a4f0f73SDmitry Karpeev   ierr = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr);
616*6a4f0f73SDmitry Karpeev   ierr = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr);
617*6a4f0f73SDmitry Karpeev   if (osm->ois) {
618*6a4f0f73SDmitry Karpeev     ierr = PCGASMDestroySubdomains(osm->n,osm->ois,osm->iis);CHKERRQ(ierr);
619*6a4f0f73SDmitry Karpeev     osm->ois = 0;
620*6a4f0f73SDmitry Karpeev     osm->iis = 0;
621b6b0974bSDmitry Karpeev   }
622a06653b4SBarry Smith   PetscFunctionReturn(0);
623a06653b4SBarry Smith }
624a06653b4SBarry Smith 
625a06653b4SBarry Smith #undef __FUNCT__
626a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM"
627a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc)
628a06653b4SBarry Smith {
629a06653b4SBarry Smith   PC_GASM         *osm = (PC_GASM*)pc->data;
630a06653b4SBarry Smith   PetscErrorCode ierr;
631a06653b4SBarry Smith   PetscInt       i;
632a06653b4SBarry Smith 
633a06653b4SBarry Smith   PetscFunctionBegin;
634a06653b4SBarry Smith   ierr = PCReset_GASM(pc);CHKERRQ(ierr);
635a06653b4SBarry Smith   if (osm->ksp) {
636a06653b4SBarry Smith     for (i=0; i<osm->n; i++) {
6376bf464f9SBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
638a06653b4SBarry Smith     }
639a06653b4SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
640a06653b4SBarry Smith   }
641a06653b4SBarry Smith   ierr = PetscFree(osm->x);CHKERRQ(ierr);
642a06653b4SBarry Smith   ierr = PetscFree(osm->y);CHKERRQ(ierr);
643c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
644b1a0a8a3SJed Brown   PetscFunctionReturn(0);
645b1a0a8a3SJed Brown }
646b1a0a8a3SJed Brown 
647b1a0a8a3SJed Brown #undef __FUNCT__
648f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM"
649ab3e923fSDmitry Karpeev static PetscErrorCode PCSetFromOptions_GASM(PC pc) {
650f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
651b1a0a8a3SJed Brown   PetscErrorCode ierr;
652b1a0a8a3SJed Brown   PetscInt       blocks,ovl;
653b1a0a8a3SJed Brown   PetscBool      symset,flg;
654f746d493SDmitry Karpeev   PCGASMType      gasmtype;
655b1a0a8a3SJed Brown 
656b1a0a8a3SJed Brown   PetscFunctionBegin;
657b1a0a8a3SJed Brown   /* set the type to symmetric if matrix is symmetric */
658b1a0a8a3SJed Brown   if (!osm->type_set && pc->pmat) {
659b1a0a8a3SJed Brown     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
660f746d493SDmitry Karpeev     if (symset && flg) { osm->type = PC_GASM_BASIC; }
661b1a0a8a3SJed Brown   }
66244fe18e5SDmitry Karpeev   ierr = PetscOptionsHead("Generalized additive Schwarz options");CHKERRQ(ierr);
663*6a4f0f73SDmitry Karpeev     ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
664*6a4f0f73SDmitry Karpeev     osm->create_local = PETSC_TRUE;
665*6a4f0f73SDmitry 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);
666*6a4f0f73SDmitry Karpeev     if(!osm->create_local) SETERRQ(((PetscObject)pc)->comm, PETSC_ERR_SUP, "No support for autocreation of nonlocal subdomains yet.");
667*6a4f0f73SDmitry Karpeev 
668*6a4f0f73SDmitry Karpeev     if (flg) {ierr = PCGASMSetTotalSubdomains(pc,blocks,osm->create_local);CHKERRQ(ierr); }
66906b43e7eSDmitry Karpeev     ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
670f746d493SDmitry Karpeev     if (flg) {ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); }
671b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
672f746d493SDmitry Karpeev     ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr);
673f746d493SDmitry Karpeev     if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr); }
674b1a0a8a3SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
675b1a0a8a3SJed Brown   PetscFunctionReturn(0);
676b1a0a8a3SJed Brown }
677b1a0a8a3SJed Brown 
678b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/
679b1a0a8a3SJed Brown 
680b1a0a8a3SJed Brown EXTERN_C_BEGIN
681b1a0a8a3SJed Brown #undef __FUNCT__
68206b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM"
683*6a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS is_local[],IS is[])
684b1a0a8a3SJed Brown {
685f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
686b1a0a8a3SJed Brown   PetscErrorCode ierr;
687b1a0a8a3SJed Brown   PetscInt       i;
688b1a0a8a3SJed Brown 
689b1a0a8a3SJed Brown   PetscFunctionBegin;
69069ca1f37SDmitry Karpeev   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, n = %D",n);
691*6a4f0f73SDmitry Karpeev   if (pc->setupcalled && (n != osm->n || is_local)) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp().");
692b1a0a8a3SJed Brown 
693b1a0a8a3SJed Brown   if (!pc->setupcalled) {
69493c1ec32SDmitry Karpeev     osm->n            = n;
695*6a4f0f73SDmitry Karpeev     osm->ois           = 0;
696*6a4f0f73SDmitry Karpeev     osm->iis     = 0;
697b1a0a8a3SJed Brown     if (is) {
698b1a0a8a3SJed Brown       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
699b1a0a8a3SJed Brown     }
700b1a0a8a3SJed Brown     if (is_local) {
701b1a0a8a3SJed Brown       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
702b1a0a8a3SJed Brown     }
703*6a4f0f73SDmitry Karpeev     if (osm->ois) {
704*6a4f0f73SDmitry Karpeev       ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr);
705b1a0a8a3SJed Brown     }
706b1a0a8a3SJed Brown     if (is) {
707*6a4f0f73SDmitry Karpeev       ierr = PetscMalloc(n*sizeof(IS),&osm->ois);CHKERRQ(ierr);
708*6a4f0f73SDmitry Karpeev       for (i=0; i<n; i++) { osm->ois[i] = is[i]; }
709*6a4f0f73SDmitry Karpeev       /* Flag indicating that the user has set outer subdomains, so PCGASM should not increase their size. */
710b1a0a8a3SJed Brown       osm->overlap = -1;
711b1a0a8a3SJed Brown     }
712b1a0a8a3SJed Brown     if (is_local) {
713*6a4f0f73SDmitry Karpeev       ierr = PetscMalloc(n*sizeof(IS),&osm->iis);CHKERRQ(ierr);
714*6a4f0f73SDmitry Karpeev       for (i=0; i<n; i++) { osm->iis[i] = is_local[i]; }
715b1a0a8a3SJed Brown     }
716b1a0a8a3SJed Brown   }
717b1a0a8a3SJed Brown   PetscFunctionReturn(0);
718b1a0a8a3SJed Brown }
719b1a0a8a3SJed Brown EXTERN_C_END
720b1a0a8a3SJed Brown 
721b1a0a8a3SJed Brown EXTERN_C_BEGIN
722b1a0a8a3SJed Brown #undef __FUNCT__
723*6a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains_GASM"
724*6a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetTotalSubdomains_GASM(PC pc,PetscInt N, PetscBool create_local) {
725f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
726b1a0a8a3SJed Brown   PetscErrorCode ierr;
727b1a0a8a3SJed Brown   PetscMPIInt    rank,size;
728b1a0a8a3SJed Brown   PetscInt       n;
729*6a4f0f73SDmitry Karpeev   PetscInt       Nmin, Nmax;
730b1a0a8a3SJed Brown   PetscFunctionBegin;
731*6a4f0f73SDmitry Karpeev   if(!create_local) SETERRQ(((PetscObject)pc)->comm, PETSC_ERR_SUP, "No suppor for autocreation of nonlocal subdomains.");
732*6a4f0f73SDmitry Karpeev   if (N < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be > 0, N = %D",N);
733*6a4f0f73SDmitry Karpeev   ierr = MPI_Allreduce(&N,&Nmin,1,MPI_INT,MPI_MIN,((PetscObject)pc)->comm); CHKERRQ(ierr);
734*6a4f0f73SDmitry Karpeev   ierr = MPI_Allreduce(&N,&Nmax,1,MPI_INT,MPI_MAX,((PetscObject)pc)->comm); CHKERRQ(ierr);
735*6a4f0f73SDmitry Karpeev   if(Nmin != Nmax)
736*6a4f0f73SDmitry Karpeev     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);
737b1a0a8a3SJed Brown 
738*6a4f0f73SDmitry Karpeev   osm->create_local = create_local;
739b1a0a8a3SJed Brown   /*
740b1a0a8a3SJed Brown      Split the subdomains equally among all processors
741b1a0a8a3SJed Brown   */
742b1a0a8a3SJed Brown   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
743b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
744b1a0a8a3SJed Brown   n = N/size + ((N % size) > rank);
745*6a4f0f73SDmitry 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);
746f746d493SDmitry Karpeev   if (pc->setupcalled && n != osm->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before PCSetUp().");
747b1a0a8a3SJed Brown   if (!pc->setupcalled) {
748*6a4f0f73SDmitry Karpeev     if (osm->ois) {
749*6a4f0f73SDmitry Karpeev       ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr);
750b1a0a8a3SJed Brown     }
75193c1ec32SDmitry Karpeev     osm->N            = N;
752f746d493SDmitry Karpeev     osm->n            = n;
75306b43e7eSDmitry Karpeev     osm->nmax         = N/size + ((N%size)?1:0);
754*6a4f0f73SDmitry Karpeev     osm->ois           = 0;
755*6a4f0f73SDmitry Karpeev     osm->iis     = 0;
756b1a0a8a3SJed Brown   }
757b1a0a8a3SJed Brown   PetscFunctionReturn(0);
75806b43e7eSDmitry Karpeev }
759b1a0a8a3SJed Brown EXTERN_C_END
760b1a0a8a3SJed Brown 
761b1a0a8a3SJed Brown EXTERN_C_BEGIN
762b1a0a8a3SJed Brown #undef __FUNCT__
763f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM"
7647087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
765b1a0a8a3SJed Brown {
766f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
767b1a0a8a3SJed Brown 
768b1a0a8a3SJed Brown   PetscFunctionBegin;
769b1a0a8a3SJed Brown   if (ovl < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
770f746d493SDmitry Karpeev   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
771b1a0a8a3SJed Brown   if (!pc->setupcalled) {
772b1a0a8a3SJed Brown     osm->overlap = ovl;
773b1a0a8a3SJed Brown   }
774b1a0a8a3SJed Brown   PetscFunctionReturn(0);
775b1a0a8a3SJed Brown }
776b1a0a8a3SJed Brown EXTERN_C_END
777b1a0a8a3SJed Brown 
778b1a0a8a3SJed Brown EXTERN_C_BEGIN
779b1a0a8a3SJed Brown #undef __FUNCT__
780f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM"
7817087cfbeSBarry Smith PetscErrorCode  PCGASMSetType_GASM(PC pc,PCGASMType type)
782b1a0a8a3SJed Brown {
783f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
784b1a0a8a3SJed Brown 
785b1a0a8a3SJed Brown   PetscFunctionBegin;
786b1a0a8a3SJed Brown   osm->type     = type;
787b1a0a8a3SJed Brown   osm->type_set = PETSC_TRUE;
788b1a0a8a3SJed Brown   PetscFunctionReturn(0);
789b1a0a8a3SJed Brown }
790b1a0a8a3SJed Brown EXTERN_C_END
791b1a0a8a3SJed Brown 
792b1a0a8a3SJed Brown EXTERN_C_BEGIN
793b1a0a8a3SJed Brown #undef __FUNCT__
794f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM"
7957087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices_GASM(PC pc,PetscBool  doSort)
796b1a0a8a3SJed Brown {
797f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
798b1a0a8a3SJed Brown 
799b1a0a8a3SJed Brown   PetscFunctionBegin;
800b1a0a8a3SJed Brown   osm->sort_indices = doSort;
801b1a0a8a3SJed Brown   PetscFunctionReturn(0);
802b1a0a8a3SJed Brown }
803b1a0a8a3SJed Brown EXTERN_C_END
804b1a0a8a3SJed Brown 
805b1a0a8a3SJed Brown EXTERN_C_BEGIN
806b1a0a8a3SJed Brown #undef __FUNCT__
807f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM"
80844fe18e5SDmitry Karpeev /*
80944fe18e5SDmitry Karpeev    FIX: This routine might need to be modified once multiple ranks per subdomain are allowed.
81044fe18e5SDmitry Karpeev         In particular, it would upset the global subdomain number calculation.
81144fe18e5SDmitry Karpeev */
8127087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
813b1a0a8a3SJed Brown {
814f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
815b1a0a8a3SJed Brown   PetscErrorCode ierr;
816b1a0a8a3SJed Brown 
817b1a0a8a3SJed Brown   PetscFunctionBegin;
818ab3e923fSDmitry 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");
819b1a0a8a3SJed Brown 
820ab3e923fSDmitry Karpeev   if (n) {
821ab3e923fSDmitry Karpeev     *n = osm->n;
822b1a0a8a3SJed Brown   }
823ab3e923fSDmitry Karpeev   if (first) {
824ab3e923fSDmitry Karpeev     ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr);
825ab3e923fSDmitry Karpeev     *first -= osm->n;
826b1a0a8a3SJed Brown   }
827b1a0a8a3SJed Brown   if (ksp) {
828b1a0a8a3SJed Brown     /* Assume that local solves are now different; not necessarily
82906b43e7eSDmitry Karpeev        true, though!  This flag is used only for PCView_GASM() */
830b1a0a8a3SJed Brown     *ksp                   = osm->ksp;
831*6a4f0f73SDmitry Karpeev     osm->same_subdomain_solvers = PETSC_FALSE;
832b1a0a8a3SJed Brown   }
833b1a0a8a3SJed Brown   PetscFunctionReturn(0);
834ab3e923fSDmitry Karpeev }/* PCGASMGetSubKSP_GASM() */
835b1a0a8a3SJed Brown EXTERN_C_END
836b1a0a8a3SJed Brown 
837b1a0a8a3SJed Brown 
838b1a0a8a3SJed Brown #undef __FUNCT__
83906b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains"
840b1a0a8a3SJed Brown /*@C
84106b43e7eSDmitry Karpeev     PCGASMSetSubdomains - Sets the subdomains for this processor
84206b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
843b1a0a8a3SJed Brown 
844b1a0a8a3SJed Brown     Collective on PC
845b1a0a8a3SJed Brown 
846b1a0a8a3SJed Brown     Input Parameters:
847b1a0a8a3SJed Brown +   pc  - the preconditioner context
84806b43e7eSDmitry Karpeev .   n   - the number of subdomains for this processor
849*6a4f0f73SDmitry Karpeev .   iis - the index sets that define this processor's local inner subdomains
850b1a0a8a3SJed Brown          (or PETSC_NULL for PETSc to determine subdomains)
851*6a4f0f73SDmitry Karpeev -   ois- the index sets that define this processor's local outer subdomains
852*6a4f0f73SDmitry Karpeev          (or PETSC_NULL to use the same as iis)
853b1a0a8a3SJed Brown 
854b1a0a8a3SJed Brown     Notes:
85506b43e7eSDmitry Karpeev     The IS indices use the parallel, global numbering of the vector entries.
856*6a4f0f73SDmitry Karpeev     Inner subdomains are those where the correction is applied.
857*6a4f0f73SDmitry Karpeev     Outer subdomains are those where the residual necessary to obtain the
858*6a4f0f73SDmitry Karpeev     corrections is obtained (see PCGASMType for the use of inner/outer subdomains).
859*6a4f0f73SDmitry Karpeev     Both inner and outer subdomains can extend over several processors.
860*6a4f0f73SDmitry Karpeev     This processor's portion of a subdomain is known as a local subdomain.
861*6a4f0f73SDmitry Karpeev 
862*6a4f0f73SDmitry Karpeev     By default the GASM preconditioner uses 1 (local) subdomain per processor.
863*6a4f0f73SDmitry Karpeev     Use PCGASMSetTotalSubdomains() to set the total number of subdomains across
864*6a4f0f73SDmitry Karpeev     all processors that PCGASM will create automatically, and to specify whether
865*6a4f0f73SDmitry Karpeev     they should be local or not.
866*6a4f0f73SDmitry Karpeev 
867b1a0a8a3SJed Brown 
868b1a0a8a3SJed Brown     Level: advanced
869b1a0a8a3SJed Brown 
87006b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz
871b1a0a8a3SJed Brown 
872*6a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
87306b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
874b1a0a8a3SJed Brown @*/
875*6a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[])
876b1a0a8a3SJed Brown {
877b1a0a8a3SJed Brown   PetscErrorCode ierr;
878b1a0a8a3SJed Brown 
879b1a0a8a3SJed Brown   PetscFunctionBegin;
880b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
881*6a4f0f73SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr);
882b1a0a8a3SJed Brown   PetscFunctionReturn(0);
883b1a0a8a3SJed Brown }
884b1a0a8a3SJed Brown 
885b1a0a8a3SJed Brown #undef __FUNCT__
886*6a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains"
887b1a0a8a3SJed Brown /*@C
888*6a4f0f73SDmitry Karpeev     PCGASMSetTotalSubdomains - Sets the total number of subdomains to use in the generalized additive
889*6a4f0f73SDmitry Karpeev     Schwarz preconditioner.  The number of subdomains is cumulative across all processors in pc's
890*6a4f0f73SDmitry Karpeev     communicator. Either all or no processors in the PC communicator must call this routine with
891*6a4f0f73SDmitry Karpeev     the same N.  The subdomains will be created automatically during PCSetUp().
892b1a0a8a3SJed Brown 
893b1a0a8a3SJed Brown     Collective on PC
894b1a0a8a3SJed Brown 
895b1a0a8a3SJed Brown     Input Parameters:
896b1a0a8a3SJed Brown +   pc           - the preconditioner context
897*6a4f0f73SDmitry Karpeev .   N            - the total number of subdomains cumulative across all processors
898*6a4f0f73SDmitry Karpeev -   create_local - whether the subdomains to be created are to be local
899b1a0a8a3SJed Brown 
900b1a0a8a3SJed Brown     Options Database Key:
901*6a4f0f73SDmitry Karpeev     To set the total number of subdomains and let PCGASM autocreate them, rather than specify the index sets, use the following options:
902*6a4f0f73SDmitry Karpeev +    -pc_gasm_total_subdomains <n>                  - sets the total number of subdomains to be autocreated by PCGASM
903*6a4f0f73SDmitry Karpeev -    -pc_gasm_subdomains_create_local <true|false>  - whether autocreated subdomains should be local or not (default is true)
904b1a0a8a3SJed Brown 
90506b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.
906b1a0a8a3SJed Brown 
907*6a4f0f73SDmitry Karpeev 
908*6a4f0f73SDmitry Karpeev     Use PCGASMSetSubdomains() to set subdomains explicitly or to set different numbers
909*6a4f0f73SDmitry Karpeev     of subdomains per processor.
910b1a0a8a3SJed Brown 
911b1a0a8a3SJed Brown     Level: advanced
912b1a0a8a3SJed Brown 
913f746d493SDmitry Karpeev .keywords: PC, GASM, set, total, global, subdomains, additive Schwarz
914b1a0a8a3SJed Brown 
91506b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
916f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
917b1a0a8a3SJed Brown @*/
918*6a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMSetTotalSubdomains(PC pc,PetscInt N, PetscBool create_local)
919b1a0a8a3SJed Brown {
920b1a0a8a3SJed Brown   PetscErrorCode ierr;
921b1a0a8a3SJed Brown 
922b1a0a8a3SJed Brown   PetscFunctionBegin;
923b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
924*6a4f0f73SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetTotalSubdomains_C",(PC,PetscInt,PetscBool),(pc,N,create_local));CHKERRQ(ierr);
925b1a0a8a3SJed Brown   PetscFunctionReturn(0);
926b1a0a8a3SJed Brown }
927b1a0a8a3SJed Brown 
928b1a0a8a3SJed Brown #undef __FUNCT__
929f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap"
930b1a0a8a3SJed Brown /*@
931f746d493SDmitry Karpeev     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
932b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
933b1a0a8a3SJed Brown     PC communicator must call this routine.
934b1a0a8a3SJed Brown 
935b1a0a8a3SJed Brown     Logically Collective on PC
936b1a0a8a3SJed Brown 
937b1a0a8a3SJed Brown     Input Parameters:
938b1a0a8a3SJed Brown +   pc  - the preconditioner context
939b1a0a8a3SJed Brown -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
940b1a0a8a3SJed Brown 
941b1a0a8a3SJed Brown     Options Database Key:
94206b43e7eSDmitry Karpeev .   -pc_gasm_overlap <overlap> - Sets overlap
943b1a0a8a3SJed Brown 
944b1a0a8a3SJed Brown     Notes:
94506b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.  To use
946*6a4f0f73SDmitry Karpeev     multiple subdomain per perocessor, see PCGASMSetTotalSubdomains() or
947*6a4f0f73SDmitry Karpeev     PCGASMSetSubdomains() (and the option -pc_gasm_total_subdomains <n>).
948b1a0a8a3SJed Brown 
949b1a0a8a3SJed Brown     The overlap defaults to 1, so if one desires that no additional
950b1a0a8a3SJed Brown     overlap be computed beyond what may have been set with a call to
951*6a4f0f73SDmitry Karpeev     PCGASMSetTotalSubdomains() or PCGASMSetSubdomains(), then ovl
952b1a0a8a3SJed Brown     must be set to be 0.  In particular, if one does not explicitly set
95306b43e7eSDmitry Karpeev     the subdomains in application code, then all overlap would be computed
954f746d493SDmitry Karpeev     internally by PETSc, and using an overlap of 0 would result in an GASM
955b1a0a8a3SJed Brown     variant that is equivalent to the block Jacobi preconditioner.
956b1a0a8a3SJed Brown 
957b1a0a8a3SJed Brown     Note that one can define initial index sets with any overlap via
95806b43e7eSDmitry Karpeev     PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows
95906b43e7eSDmitry Karpeev     PETSc to extend that overlap further, if desired.
960b1a0a8a3SJed Brown 
961b1a0a8a3SJed Brown     Level: intermediate
962b1a0a8a3SJed Brown 
963f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap
964b1a0a8a3SJed Brown 
965*6a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
96606b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
967b1a0a8a3SJed Brown @*/
9687087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap(PC pc,PetscInt ovl)
969b1a0a8a3SJed Brown {
970b1a0a8a3SJed Brown   PetscErrorCode ierr;
971b1a0a8a3SJed Brown 
972b1a0a8a3SJed Brown   PetscFunctionBegin;
973b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
974b1a0a8a3SJed Brown   PetscValidLogicalCollectiveInt(pc,ovl,2);
975f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
976b1a0a8a3SJed Brown   PetscFunctionReturn(0);
977b1a0a8a3SJed Brown }
978b1a0a8a3SJed Brown 
979b1a0a8a3SJed Brown #undef __FUNCT__
980f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType"
981b1a0a8a3SJed Brown /*@
982f746d493SDmitry Karpeev     PCGASMSetType - Sets the type of restriction and interpolation used
983b1a0a8a3SJed Brown     for local problems in the additive Schwarz method.
984b1a0a8a3SJed Brown 
985b1a0a8a3SJed Brown     Logically Collective on PC
986b1a0a8a3SJed Brown 
987b1a0a8a3SJed Brown     Input Parameters:
988b1a0a8a3SJed Brown +   pc  - the preconditioner context
989f746d493SDmitry Karpeev -   type - variant of GASM, one of
990b1a0a8a3SJed Brown .vb
991f746d493SDmitry Karpeev       PC_GASM_BASIC       - full interpolation and restriction
992f746d493SDmitry Karpeev       PC_GASM_RESTRICT    - full restriction, local processor interpolation
993f746d493SDmitry Karpeev       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
994f746d493SDmitry Karpeev       PC_GASM_NONE        - local processor restriction and interpolation
995b1a0a8a3SJed Brown .ve
996b1a0a8a3SJed Brown 
997b1a0a8a3SJed Brown     Options Database Key:
998f746d493SDmitry Karpeev .   -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
999b1a0a8a3SJed Brown 
1000b1a0a8a3SJed Brown     Level: intermediate
1001b1a0a8a3SJed Brown 
1002f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1003b1a0a8a3SJed Brown 
1004*6a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1005f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1006b1a0a8a3SJed Brown @*/
10077087cfbeSBarry Smith PetscErrorCode  PCGASMSetType(PC pc,PCGASMType type)
1008b1a0a8a3SJed Brown {
1009b1a0a8a3SJed Brown   PetscErrorCode ierr;
1010b1a0a8a3SJed Brown 
1011b1a0a8a3SJed Brown   PetscFunctionBegin;
1012b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1013b1a0a8a3SJed Brown   PetscValidLogicalCollectiveEnum(pc,type,2);
1014f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr);
1015b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1016b1a0a8a3SJed Brown }
1017b1a0a8a3SJed Brown 
1018b1a0a8a3SJed Brown #undef __FUNCT__
1019f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices"
1020b1a0a8a3SJed Brown /*@
1021f746d493SDmitry Karpeev     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
1022b1a0a8a3SJed Brown 
1023b1a0a8a3SJed Brown     Logically Collective on PC
1024b1a0a8a3SJed Brown 
1025b1a0a8a3SJed Brown     Input Parameters:
1026b1a0a8a3SJed Brown +   pc  - the preconditioner context
1027b1a0a8a3SJed Brown -   doSort - sort the subdomain indices
1028b1a0a8a3SJed Brown 
1029b1a0a8a3SJed Brown     Level: intermediate
1030b1a0a8a3SJed Brown 
1031f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
1032b1a0a8a3SJed Brown 
1033*6a4f0f73SDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(),
1034f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1035b1a0a8a3SJed Brown @*/
10367087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices(PC pc,PetscBool  doSort)
1037b1a0a8a3SJed Brown {
1038b1a0a8a3SJed Brown   PetscErrorCode ierr;
1039b1a0a8a3SJed Brown 
1040b1a0a8a3SJed Brown   PetscFunctionBegin;
1041b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1042b1a0a8a3SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
1043f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
1044b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1045b1a0a8a3SJed Brown }
1046b1a0a8a3SJed Brown 
1047b1a0a8a3SJed Brown #undef __FUNCT__
1048f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP"
1049b1a0a8a3SJed Brown /*@C
1050f746d493SDmitry Karpeev    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
1051b1a0a8a3SJed Brown    this processor.
1052b1a0a8a3SJed Brown 
1053b1a0a8a3SJed Brown    Collective on PC iff first_local is requested
1054b1a0a8a3SJed Brown 
1055b1a0a8a3SJed Brown    Input Parameter:
1056b1a0a8a3SJed Brown .  pc - the preconditioner context
1057b1a0a8a3SJed Brown 
1058b1a0a8a3SJed Brown    Output Parameters:
1059b1a0a8a3SJed Brown +  n_local - the number of blocks on this processor or PETSC_NULL
1060b1a0a8a3SJed Brown .  first_local - the global number of the first block on this processor or PETSC_NULL,
1061b1a0a8a3SJed Brown                  all processors must request or all must pass PETSC_NULL
1062b1a0a8a3SJed Brown -  ksp - the array of KSP contexts
1063b1a0a8a3SJed Brown 
1064b1a0a8a3SJed Brown    Note:
1065f746d493SDmitry Karpeev    After PCGASMGetSubKSP() the array of KSPes is not to be freed
1066b1a0a8a3SJed Brown 
1067b1a0a8a3SJed Brown    Currently for some matrix implementations only 1 block per processor
1068b1a0a8a3SJed Brown    is supported.
1069b1a0a8a3SJed Brown 
1070f746d493SDmitry Karpeev    You must call KSPSetUp() before calling PCGASMGetSubKSP().
1071b1a0a8a3SJed Brown 
1072b1a0a8a3SJed Brown    Level: advanced
1073b1a0a8a3SJed Brown 
1074f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context
1075b1a0a8a3SJed Brown 
1076*6a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap(),
1077f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(),
1078b1a0a8a3SJed Brown @*/
10797087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1080b1a0a8a3SJed Brown {
1081b1a0a8a3SJed Brown   PetscErrorCode ierr;
1082b1a0a8a3SJed Brown 
1083b1a0a8a3SJed Brown   PetscFunctionBegin;
1084b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1085f746d493SDmitry Karpeev   ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
1086b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1087b1a0a8a3SJed Brown }
1088b1a0a8a3SJed Brown 
1089b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/
1090b1a0a8a3SJed Brown /*MC
1091f746d493SDmitry Karpeev    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1092b1a0a8a3SJed Brown            its own KSP object.
1093b1a0a8a3SJed Brown 
1094b1a0a8a3SJed Brown    Options Database Keys:
109506b43e7eSDmitry Karpeev +  -pc_gasm_total_block_count <n> - Sets total number of local subdomains (known as blocks) to be distributed among processors
109606b43e7eSDmitry Karpeev .  -pc_gasm_view_subdomains       - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view
109706b43e7eSDmitry Karpeev .  -pc_gasm_print_subdomains      - activates the printing of subdomain indices in PCSetUp()
109806b43e7eSDmitry Karpeev .  -pc_gasm_overlap <ovl>         - Sets overlap by which to (automatically) extend local subdomains
1099f746d493SDmitry Karpeev -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1100b1a0a8a3SJed Brown 
1101b1a0a8a3SJed Brown      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
1102f746d493SDmitry Karpeev       will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
1103f746d493SDmitry Karpeev       -pc_gasm_type basic to use the standard GASM.
1104b1a0a8a3SJed Brown 
1105b1a0a8a3SJed Brown    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
1106b1a0a8a3SJed Brown      than one processor. Defaults to one block per processor.
1107b1a0a8a3SJed Brown 
1108b1a0a8a3SJed Brown      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1109b1a0a8a3SJed Brown         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
1110b1a0a8a3SJed Brown 
1111f746d493SDmitry Karpeev      To set the options on the solvers separate for each block call PCGASMGetSubKSP()
1112b1a0a8a3SJed Brown          and set the options directly on the resulting KSP object (you can access its PC
1113b1a0a8a3SJed Brown          with KSPGetPC())
1114b1a0a8a3SJed Brown 
1115b1a0a8a3SJed Brown 
1116b1a0a8a3SJed Brown    Level: beginner
1117b1a0a8a3SJed Brown 
1118b1a0a8a3SJed Brown    Concepts: additive Schwarz method
1119b1a0a8a3SJed Brown 
1120b1a0a8a3SJed Brown     References:
1121b1a0a8a3SJed Brown     An additive variant of the Schwarz alternating method for the case of many subregions
1122b1a0a8a3SJed Brown     M Dryja, OB Widlund - Courant Institute, New York University Technical report
1123b1a0a8a3SJed Brown 
1124b1a0a8a3SJed Brown     Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1125b1a0a8a3SJed Brown     Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X.
1126b1a0a8a3SJed Brown 
1127b1a0a8a3SJed Brown .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
112806b43e7eSDmitry Karpeev            PCBJACOBI, PCGASMSetUseTrueLocal(), PCGASMGetSubKSP(), PCGASMSetSubdomains(),
1129*6a4f0f73SDmitry Karpeev            PCGASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType()
1130b1a0a8a3SJed Brown 
1131b1a0a8a3SJed Brown M*/
1132b1a0a8a3SJed Brown 
1133b1a0a8a3SJed Brown EXTERN_C_BEGIN
1134b1a0a8a3SJed Brown #undef __FUNCT__
1135f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM"
11367087cfbeSBarry Smith PetscErrorCode  PCCreate_GASM(PC pc)
1137b1a0a8a3SJed Brown {
1138b1a0a8a3SJed Brown   PetscErrorCode ierr;
1139f746d493SDmitry Karpeev   PC_GASM         *osm;
1140b1a0a8a3SJed Brown 
1141b1a0a8a3SJed Brown   PetscFunctionBegin;
1142f746d493SDmitry Karpeev   ierr = PetscNewLog(pc,PC_GASM,&osm);CHKERRQ(ierr);
1143ab3e923fSDmitry Karpeev   osm->N                 = PETSC_DECIDE;
114406b43e7eSDmitry Karpeev   osm->n                 = PETSC_DECIDE;
1145ab3e923fSDmitry Karpeev   osm->nmax              = 0;
1146b1a0a8a3SJed Brown   osm->overlap           = 1;
1147b1a0a8a3SJed Brown   osm->ksp               = 0;
1148*6a4f0f73SDmitry Karpeev   osm->gorestriction     = 0;
1149*6a4f0f73SDmitry Karpeev   osm->girestriction     = 0;
1150ab3e923fSDmitry Karpeev   osm->gx                = 0;
1151ab3e923fSDmitry Karpeev   osm->gy                = 0;
1152b1a0a8a3SJed Brown   osm->x                 = 0;
1153b1a0a8a3SJed Brown   osm->y                 = 0;
1154*6a4f0f73SDmitry Karpeev   osm->ois               = 0;
1155*6a4f0f73SDmitry Karpeev   osm->iis               = 0;
1156b1a0a8a3SJed Brown   osm->pmat              = 0;
1157f746d493SDmitry Karpeev   osm->type              = PC_GASM_RESTRICT;
1158*6a4f0f73SDmitry Karpeev   osm->same_subdomain_solvers = PETSC_TRUE;
1159b1a0a8a3SJed Brown   osm->sort_indices           = PETSC_TRUE;
1160b1a0a8a3SJed Brown 
1161b1a0a8a3SJed Brown   pc->data                   = (void*)osm;
1162f746d493SDmitry Karpeev   pc->ops->apply             = PCApply_GASM;
1163f746d493SDmitry Karpeev   pc->ops->applytranspose    = PCApplyTranspose_GASM;
1164f746d493SDmitry Karpeev   pc->ops->setup             = PCSetUp_GASM;
1165a06653b4SBarry Smith   pc->ops->reset             = PCReset_GASM;
1166f746d493SDmitry Karpeev   pc->ops->destroy           = PCDestroy_GASM;
1167f746d493SDmitry Karpeev   pc->ops->setfromoptions    = PCSetFromOptions_GASM;
1168f746d493SDmitry Karpeev   pc->ops->setuponblocks     = PCSetUpOnBlocks_GASM;
1169f746d493SDmitry Karpeev   pc->ops->view              = PCView_GASM;
1170b1a0a8a3SJed Brown   pc->ops->applyrichardson   = 0;
1171b1a0a8a3SJed Brown 
117206b43e7eSDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSubdomains_C","PCGASMSetSubdomains_GASM",
117306b43e7eSDmitry Karpeev                     PCGASMSetSubdomains_GASM);CHKERRQ(ierr);
1174*6a4f0f73SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetTotalSubdomains_C","PCGASMSetTotalSubdomains_GASM",
1175*6a4f0f73SDmitry Karpeev                     PCGASMSetTotalSubdomains_GASM);CHKERRQ(ierr);
1176f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetOverlap_C","PCGASMSetOverlap_GASM",
1177f746d493SDmitry Karpeev                     PCGASMSetOverlap_GASM);CHKERRQ(ierr);
1178f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetType_C","PCGASMSetType_GASM",
1179f746d493SDmitry Karpeev                     PCGASMSetType_GASM);CHKERRQ(ierr);
1180f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSortIndices_C","PCGASMSetSortIndices_GASM",
1181f746d493SDmitry Karpeev                     PCGASMSetSortIndices_GASM);CHKERRQ(ierr);
1182f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMGetSubKSP_C","PCGASMGetSubKSP_GASM",
1183f746d493SDmitry Karpeev                     PCGASMGetSubKSP_GASM);CHKERRQ(ierr);
1184b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1185b1a0a8a3SJed Brown }
1186b1a0a8a3SJed Brown EXTERN_C_END
1187b1a0a8a3SJed Brown 
1188b1a0a8a3SJed Brown 
1189b1a0a8a3SJed Brown #undef __FUNCT__
119006b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains"
1191b1a0a8a3SJed Brown /*@C
119206b43e7eSDmitry Karpeev    PCGASMCreateLocalSubdomains - Creates n local index sets for the overlapping
119306b43e7eSDmitry Karpeev    Schwarz preconditioner for a any problem based on its matrix.
1194b1a0a8a3SJed Brown 
1195b1a0a8a3SJed Brown    Collective
1196b1a0a8a3SJed Brown 
1197b1a0a8a3SJed Brown    Input Parameters:
1198b1a0a8a3SJed Brown +  A       - The global matrix operator
1199*6a4f0f73SDmitry Karpeev .  overlap - amount of overlap in outer subdomains
120006b43e7eSDmitry Karpeev -  n       - the number of local subdomains
1201b1a0a8a3SJed Brown 
1202b1a0a8a3SJed Brown    Output Parameters:
1203*6a4f0f73SDmitry Karpeev +  iis - the array of index sets defining the local inner subdomains (on which the correction is applied)
1204*6a4f0f73SDmitry Karpeev -  ois - the array of index sets defining the local outer subdomains (on which the residual is computed)
1205b1a0a8a3SJed Brown 
1206b1a0a8a3SJed Brown    Level: advanced
1207b1a0a8a3SJed Brown 
1208*6a4f0f73SDmitry Karpeev    Note: this generates n nonoverlapping local inner subdomains on PETSC_COMM_SELF;
1209*6a4f0f73SDmitry Karpeev          PCGASM will generate the overlap from these if you use them in PCGASMSetSubdomains() and set a
1210*6a4f0f73SDmitry Karpeev          nonzero overlap with PCGASMSetOverlap()
1211b1a0a8a3SJed Brown 
1212b1a0a8a3SJed Brown     In the Fortran version you must provide the array outis[] already allocated of length n.
1213b1a0a8a3SJed Brown 
1214f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1215b1a0a8a3SJed Brown 
121606b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains()
1217b1a0a8a3SJed Brown @*/
1218*6a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt overlap, PetscInt n, IS* iis[], IS* ois[])
1219b1a0a8a3SJed Brown {
1220b1a0a8a3SJed Brown   MatPartitioning           mpart;
1221b1a0a8a3SJed Brown   const char                *prefix;
122211bd1e4dSLisandro Dalcin   PetscErrorCode            (*f)(Mat,MatReuse,Mat*);
1223b1a0a8a3SJed Brown   PetscMPIInt               size;
1224b1a0a8a3SJed Brown   PetscInt                  i,j,rstart,rend,bs;
122511bd1e4dSLisandro Dalcin   PetscBool                 isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
1226b1a0a8a3SJed Brown   Mat                       Ad = PETSC_NULL, adj;
1227b1a0a8a3SJed Brown   IS                        ispart,isnumb,*is;
1228b1a0a8a3SJed Brown   PetscErrorCode            ierr;
1229b1a0a8a3SJed Brown 
1230b1a0a8a3SJed Brown   PetscFunctionBegin;
1231b1a0a8a3SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1232*6a4f0f73SDmitry Karpeev   PetscValidPointer(iis,4);
1233b1a0a8a3SJed Brown   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
1234b1a0a8a3SJed Brown 
1235b1a0a8a3SJed Brown   /* Get prefix, row distribution, and block size */
1236b1a0a8a3SJed Brown   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1237b1a0a8a3SJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1238b1a0a8a3SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1239b1a0a8a3SJed 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);
1240b1a0a8a3SJed Brown 
1241b1a0a8a3SJed Brown   /* Get diagonal block from matrix if possible */
1242b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1243b1a0a8a3SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr);
1244b1a0a8a3SJed Brown   if (f) {
124511bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
1246b1a0a8a3SJed Brown   } else if (size == 1) {
124711bd1e4dSLisandro Dalcin     Ad = A;
1248b1a0a8a3SJed Brown   }
1249b1a0a8a3SJed Brown   if (Ad) {
1250251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1251251f4c67SDmitry Karpeev     if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
1252b1a0a8a3SJed Brown   }
1253b1a0a8a3SJed Brown   if (Ad && n > 1) {
1254b1a0a8a3SJed Brown     PetscBool  match,done;
1255b1a0a8a3SJed Brown     /* Try to setup a good matrix partitioning if available */
1256b1a0a8a3SJed Brown     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
1257b1a0a8a3SJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
1258b1a0a8a3SJed Brown     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1259251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
1260b1a0a8a3SJed Brown     if (!match) {
1261251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
1262b1a0a8a3SJed Brown     }
1263b1a0a8a3SJed Brown     if (!match) { /* assume a "good" partitioner is available */
1264b1a0a8a3SJed Brown       PetscInt na,*ia,*ja;
1265b1a0a8a3SJed Brown       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1266b1a0a8a3SJed Brown       if (done) {
1267b1a0a8a3SJed Brown         /* Build adjacency matrix by hand. Unfortunately a call to
1268b1a0a8a3SJed Brown            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1269b1a0a8a3SJed Brown            remove the block-aij structure and we cannot expect
1270b1a0a8a3SJed Brown            MatPartitioning to split vertices as we need */
1271b1a0a8a3SJed Brown         PetscInt i,j,*row,len,nnz,cnt,*iia=0,*jja=0;
1272b1a0a8a3SJed Brown         nnz = 0;
1273b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* count number of nonzeros */
1274b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1275b1a0a8a3SJed Brown           row = ja + ia[i];
1276b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1277b1a0a8a3SJed Brown             if (row[j] == i) { /* don't count diagonal */
1278b1a0a8a3SJed Brown               len--; break;
1279b1a0a8a3SJed Brown             }
1280b1a0a8a3SJed Brown           }
1281b1a0a8a3SJed Brown           nnz += len;
1282b1a0a8a3SJed Brown         }
1283b1a0a8a3SJed Brown         ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr);
1284b1a0a8a3SJed Brown         ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr);
1285b1a0a8a3SJed Brown         nnz    = 0;
1286b1a0a8a3SJed Brown         iia[0] = 0;
1287b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* fill adjacency */
1288b1a0a8a3SJed Brown           cnt = 0;
1289b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1290b1a0a8a3SJed Brown           row = ja + ia[i];
1291b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1292b1a0a8a3SJed Brown             if (row[j] != i) { /* if not diagonal */
1293b1a0a8a3SJed Brown               jja[nnz+cnt++] = row[j];
1294b1a0a8a3SJed Brown             }
1295b1a0a8a3SJed Brown           }
1296b1a0a8a3SJed Brown           nnz += cnt;
1297b1a0a8a3SJed Brown           iia[i+1] = nnz;
1298b1a0a8a3SJed Brown         }
1299b1a0a8a3SJed Brown         /* Partitioning of the adjacency matrix */
1300b1a0a8a3SJed Brown         ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);CHKERRQ(ierr);
1301b1a0a8a3SJed Brown         ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
1302b1a0a8a3SJed Brown         ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
1303b1a0a8a3SJed Brown         ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
1304b1a0a8a3SJed Brown         ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
13056bf464f9SBarry Smith         ierr = MatDestroy(&adj);CHKERRQ(ierr);
1306b1a0a8a3SJed Brown         foundpart = PETSC_TRUE;
1307b1a0a8a3SJed Brown       }
1308b1a0a8a3SJed Brown       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1309b1a0a8a3SJed Brown     }
13106bf464f9SBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1311b1a0a8a3SJed Brown   }
1312b1a0a8a3SJed Brown   ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr);
1313b1a0a8a3SJed Brown   if (!foundpart) {
1314b1a0a8a3SJed Brown 
1315b1a0a8a3SJed Brown     /* Partitioning by contiguous chunks of rows */
1316b1a0a8a3SJed Brown 
1317b1a0a8a3SJed Brown     PetscInt mbs   = (rend-rstart)/bs;
1318b1a0a8a3SJed Brown     PetscInt start = rstart;
1319b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1320b1a0a8a3SJed Brown       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
1321b1a0a8a3SJed Brown       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
1322b1a0a8a3SJed Brown       start += count;
1323b1a0a8a3SJed Brown     }
1324b1a0a8a3SJed Brown 
1325b1a0a8a3SJed Brown   } else {
1326b1a0a8a3SJed Brown 
1327b1a0a8a3SJed Brown     /* Partitioning by adjacency of diagonal block  */
1328b1a0a8a3SJed Brown 
1329b1a0a8a3SJed Brown     const PetscInt *numbering;
1330b1a0a8a3SJed Brown     PetscInt       *count,nidx,*indices,*newidx,start=0;
1331b1a0a8a3SJed Brown     /* Get node count in each partition */
1332b1a0a8a3SJed Brown     ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr);
1333b1a0a8a3SJed Brown     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
1334b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1335b1a0a8a3SJed Brown       for (i=0; i<n; i++) count[i] *= bs;
1336b1a0a8a3SJed Brown     }
1337b1a0a8a3SJed Brown     /* Build indices from node numbering */
1338b1a0a8a3SJed Brown     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1339b1a0a8a3SJed Brown     ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr);
1340b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1341b1a0a8a3SJed Brown     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
1342b1a0a8a3SJed Brown     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
1343b1a0a8a3SJed Brown     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
1344b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1345b1a0a8a3SJed Brown       ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr);
1346b1a0a8a3SJed Brown       for (i=0; i<nidx; i++)
1347b1a0a8a3SJed Brown         for (j=0; j<bs; j++)
1348b1a0a8a3SJed Brown           newidx[i*bs+j] = indices[i]*bs + j;
1349b1a0a8a3SJed Brown       ierr = PetscFree(indices);CHKERRQ(ierr);
1350b1a0a8a3SJed Brown       nidx   *= bs;
1351b1a0a8a3SJed Brown       indices = newidx;
1352b1a0a8a3SJed Brown     }
1353b1a0a8a3SJed Brown     /* Shift to get global indices */
1354b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] += rstart;
1355b1a0a8a3SJed Brown 
1356b1a0a8a3SJed Brown     /* Build the index sets for each block */
1357b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1358b1a0a8a3SJed Brown       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
1359b1a0a8a3SJed Brown       ierr   = ISSort(is[i]);CHKERRQ(ierr);
1360b1a0a8a3SJed Brown       start += count[i];
1361b1a0a8a3SJed Brown     }
1362b1a0a8a3SJed Brown 
1363b1a0a8a3SJed Brown     ierr = PetscFree(count);
1364b1a0a8a3SJed Brown     ierr = PetscFree(indices);
1365fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1366fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
1367b1a0a8a3SJed Brown   }
1368*6a4f0f73SDmitry Karpeev   *iis = is;
1369*6a4f0f73SDmitry Karpeev   if(!ois) PetscFunctionReturn(0);
1370*6a4f0f73SDmitry Karpeev   /*
1371*6a4f0f73SDmitry Karpeev    Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap
1372*6a4f0f73SDmitry Karpeev    has been requested, copy the inner subdomains over so they can be modified.
1373*6a4f0f73SDmitry Karpeev    */
1374*6a4f0f73SDmitry Karpeev   ierr = PetscMalloc(n*sizeof(IS),ois);CHKERRQ(ierr);
1375*6a4f0f73SDmitry Karpeev   for (i=0; i<n; ++i) {
1376*6a4f0f73SDmitry Karpeev     if (overlap > 0) { /* With positive overlap, (*iis)[i] will be modified */
1377*6a4f0f73SDmitry Karpeev       ierr = ISDuplicate((*iis)[i],(*ois)+i);CHKERRQ(ierr);
1378*6a4f0f73SDmitry Karpeev       ierr = ISCopy((*iis)[i],(*ois)[i]);CHKERRQ(ierr);
1379*6a4f0f73SDmitry Karpeev     } else {
1380*6a4f0f73SDmitry Karpeev       ierr = PetscObjectReference((PetscObject)(*iis)[i]);CHKERRQ(ierr);
1381*6a4f0f73SDmitry Karpeev       (*ois)[i] = (*iis)[i];
1382*6a4f0f73SDmitry Karpeev     }
1383*6a4f0f73SDmitry Karpeev   }
1384*6a4f0f73SDmitry Karpeev   if (overlap > 0) {
1385*6a4f0f73SDmitry Karpeev     /* Extend the "overlapping" regions by a number of steps */
1386*6a4f0f73SDmitry Karpeev     ierr = MatIncreaseOverlap(A,n,*ois,overlap);CHKERRQ(ierr);
1387*6a4f0f73SDmitry Karpeev   }
1388b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1389b1a0a8a3SJed Brown }
1390b1a0a8a3SJed Brown 
1391b1a0a8a3SJed Brown #undef __FUNCT__
1392f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains"
1393b1a0a8a3SJed Brown /*@C
1394f746d493SDmitry Karpeev    PCGASMDestroySubdomains - Destroys the index sets created with
139506b43e7eSDmitry Karpeev    PCGASMCreateLocalSubdomains() or PCGASMCreateSubdomains2D. Should be
139606b43e7eSDmitry Karpeev    called after setting subdomains with PCGASMSetSubdomains().
1397b1a0a8a3SJed Brown 
1398b1a0a8a3SJed Brown    Collective
1399b1a0a8a3SJed Brown 
1400b1a0a8a3SJed Brown    Input Parameters:
1401b1a0a8a3SJed Brown +  n   - the number of index sets
1402*6a4f0f73SDmitry Karpeev .  iis - the array of inner subdomains,
1403*6a4f0f73SDmitry Karpeev -  ois - the array of outer subdomains, can be PETSC_NULL
1404b1a0a8a3SJed Brown 
1405*6a4f0f73SDmitry Karpeev    Level: intermediate
1406*6a4f0f73SDmitry Karpeev 
1407*6a4f0f73SDmitry Karpeev    Notes: this is merely a convenience subroutine that walks each list,
1408*6a4f0f73SDmitry Karpeev    destroys each IS on the list, and then frees the list.
1409b1a0a8a3SJed Brown 
1410f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1411b1a0a8a3SJed Brown 
141206b43e7eSDmitry Karpeev .seealso: PCGASMCreateLocalSubdomains(), PCGASMSetSubdomains()
1413b1a0a8a3SJed Brown @*/
1414*6a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMDestroySubdomains(PetscInt n, IS iis[], IS ois[])
1415b1a0a8a3SJed Brown {
1416b1a0a8a3SJed Brown   PetscInt       i;
1417b1a0a8a3SJed Brown   PetscErrorCode ierr;
1418b1a0a8a3SJed Brown   PetscFunctionBegin;
1419b1a0a8a3SJed Brown   if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"n must be > 0: n = %D",n);
1420*6a4f0f73SDmitry Karpeev   PetscValidPointer(iis,2);
1421*6a4f0f73SDmitry Karpeev   for (i=0; i<n; i++) { ierr = ISDestroy(&iis[i]);CHKERRQ(ierr); }
1422*6a4f0f73SDmitry Karpeev   ierr = PetscFree(iis);CHKERRQ(ierr);
1423*6a4f0f73SDmitry Karpeev   if (ois) {
1424*6a4f0f73SDmitry Karpeev     for (i=0; i<n; i++) { ierr = ISDestroy(&ois[i]);CHKERRQ(ierr); }
1425*6a4f0f73SDmitry Karpeev     ierr = PetscFree(ois);CHKERRQ(ierr);
1426b1a0a8a3SJed Brown   }
1427b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1428b1a0a8a3SJed Brown }
1429b1a0a8a3SJed Brown 
1430af538404SDmitry Karpeev 
1431af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \
1432af538404SDmitry Karpeev {                                                                                                       \
1433af538404SDmitry Karpeev  PetscInt first_row = first/M, last_row = last/M+1;                                                     \
1434af538404SDmitry Karpeev   /*                                                                                                    \
1435af538404SDmitry Karpeev    Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners       \
1436af538404SDmitry Karpeev    of the bounding box of the intersection of the subdomain with the local ownership range (local       \
1437af538404SDmitry Karpeev    subdomain).                                                                                          \
1438af538404SDmitry Karpeev    Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows       \
1439af538404SDmitry Karpeev    of the intersection.                                                                                 \
1440af538404SDmitry Karpeev   */                                                                                                    \
1441af538404SDmitry Karpeev   /* ylow_loc is the grid row containing the first element of the local sumbdomain */                   \
1442eec7e3faSDmitry Karpeev   *ylow_loc = PetscMax(first_row,ylow);                                                                    \
1443af538404SDmitry Karpeev   /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1444eec7e3faSDmitry Karpeev   *xleft_loc = *ylow_loc==first_row?PetscMax(first%M,xleft):xleft;                                                                            \
1445af538404SDmitry Karpeev   /* yhigh_loc is the grid row above the last local subdomain element */                                                                    \
1446eec7e3faSDmitry Karpeev   *yhigh_loc = PetscMin(last_row,yhigh);                                                                                                     \
1447af538404SDmitry Karpeev   /* xright is the offset of the end of the  local subdomain within its grid row (might actually be outside the local subdomain) */         \
1448eec7e3faSDmitry Karpeev   *xright_loc = *yhigh_loc==last_row?PetscMin(xright,last%M):xright;                                                                          \
1449af538404SDmitry Karpeev   /* Now compute the size of the local subdomain n. */ \
1450c3518bceSDmitry Karpeev   *n = 0;                                               \
1451eec7e3faSDmitry Karpeev   if(*ylow_loc < *yhigh_loc) {                           \
1452af538404SDmitry Karpeev     PetscInt width = xright-xleft;                     \
1453eec7e3faSDmitry Karpeev     *n += width*(*yhigh_loc-*ylow_loc-1);                 \
1454eec7e3faSDmitry Karpeev     *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \
1455eec7e3faSDmitry Karpeev     *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \
1456af538404SDmitry Karpeev   }\
1457af538404SDmitry Karpeev }
1458af538404SDmitry Karpeev 
1459c3518bceSDmitry Karpeev 
1460eec7e3faSDmitry Karpeev 
1461eec7e3faSDmitry Karpeev #undef __FUNCT__
1462f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D"
1463b1a0a8a3SJed Brown /*@
1464f746d493SDmitry Karpeev    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1465b1a0a8a3SJed Brown    preconditioner for a two-dimensional problem on a regular grid.
1466b1a0a8a3SJed Brown 
1467af538404SDmitry Karpeev    Collective
1468b1a0a8a3SJed Brown 
1469b1a0a8a3SJed Brown    Input Parameters:
147006b43e7eSDmitry Karpeev +  M, N               - the global number of grid points in the x and y directions
1471af538404SDmitry Karpeev .  Mdomains, Ndomains - the global number of subdomains in the x and y directions
1472b1a0a8a3SJed Brown .  dof                - degrees of freedom per node
1473b1a0a8a3SJed Brown -  overlap            - overlap in mesh lines
1474b1a0a8a3SJed Brown 
1475b1a0a8a3SJed Brown    Output Parameters:
1476af538404SDmitry Karpeev +  Nsub - the number of local subdomains created
1477*6a4f0f73SDmitry Karpeev .  iis  - array of index sets defining inner (nonoverlapping) subdomains
1478*6a4f0f73SDmitry Karpeev -  ois  - array of index sets defining outer (overlapping, if overlap > 0) subdomains
1479b1a0a8a3SJed Brown 
1480b1a0a8a3SJed Brown 
1481b1a0a8a3SJed Brown    Level: advanced
1482b1a0a8a3SJed Brown 
1483f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid
1484b1a0a8a3SJed Brown 
1485*6a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1486f746d493SDmitry Karpeev           PCGASMSetOverlap()
1487b1a0a8a3SJed Brown @*/
1488*6a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMCreateSubdomains2D(PC pc, PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap, PetscInt *nsub,IS **iis,IS **ois)
1489b1a0a8a3SJed Brown {
1490b1a0a8a3SJed Brown   PetscErrorCode ierr;
14912bbb417fSDmitry Karpeev   PetscMPIInt    size, rank;
14922bbb417fSDmitry Karpeev   PetscInt       i, j;
14932bbb417fSDmitry Karpeev   PetscInt       maxheight, maxwidth;
14942bbb417fSDmitry Karpeev   PetscInt       xstart, xleft, xright, xleft_loc, xright_loc;
14952bbb417fSDmitry Karpeev   PetscInt       ystart, ylow,  yhigh,  ylow_loc,  yhigh_loc;
1496eec7e3faSDmitry Karpeev   PetscInt       x[2][2], y[2][2], n[2];
14972bbb417fSDmitry Karpeev   PetscInt       first, last;
14982bbb417fSDmitry Karpeev   PetscInt       nidx, *idx;
14992bbb417fSDmitry Karpeev   PetscInt       ii,jj,s,q,d;
1500af538404SDmitry Karpeev   PetscInt       k,kk;
15012bbb417fSDmitry Karpeev   PetscMPIInt    color;
15022bbb417fSDmitry Karpeev   MPI_Comm       comm, subcomm;
1503*6a4f0f73SDmitry Karpeev   IS             **xis = 0, **is = ois, **is_local = iis;
1504d34fcf5fSBarry Smith 
1505b1a0a8a3SJed Brown   PetscFunctionBegin;
15062bbb417fSDmitry Karpeev   ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr);
15072bbb417fSDmitry Karpeev   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
15082bbb417fSDmitry Karpeev   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
15092bbb417fSDmitry Karpeev   ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr);
1510d34fcf5fSBarry 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) "
15112bbb417fSDmitry Karpeev 	     "does not respect the number of degrees of freedom per grid point %D", first, last, dof);
1512d34fcf5fSBarry Smith 
1513af538404SDmitry Karpeev   /* Determine the number of domains with nonzero intersections with the local ownership range. */
15142bbb417fSDmitry Karpeev   s = 0;
1515af538404SDmitry Karpeev   ystart = 0;
1516af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1517af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1518af538404SDmitry 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);
1519eec7e3faSDmitry Karpeev     /* Vertical domain limits with an overlap. */
1520eec7e3faSDmitry Karpeev     ylow = PetscMax(ystart - overlap,0);
1521eec7e3faSDmitry Karpeev     yhigh = PetscMin(ystart + maxheight + overlap,N);
1522b1a0a8a3SJed Brown     xstart = 0;
1523af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1524af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1525af538404SDmitry 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);
1526eec7e3faSDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1527eec7e3faSDmitry Karpeev       xleft   = PetscMax(xstart - overlap,0);
1528eec7e3faSDmitry Karpeev       xright  = PetscMin(xstart + maxwidth + overlap,M);
1529af538404SDmitry Karpeev       /*
1530af538404SDmitry Karpeev 	 Determine whether this subdomain intersects this processor's ownership range of pc->pmat.
1531af538404SDmitry Karpeev       */
1532c3518bceSDmitry Karpeev       PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1533af538404SDmitry Karpeev       if(nidx) {
1534af538404SDmitry Karpeev         ++s;
1535af538404SDmitry Karpeev       }
1536af538404SDmitry Karpeev       xstart += maxwidth;
1537af538404SDmitry Karpeev     }/* for(i = 0; i < Mdomains; ++i) */
1538af538404SDmitry Karpeev     ystart += maxheight;
1539af538404SDmitry Karpeev   }/* for(j = 0; j < Ndomains; ++j) */
1540af538404SDmitry Karpeev   /* Now we can allocate the necessary number of ISs. */
1541af538404SDmitry Karpeev   *nsub = s;
1542af538404SDmitry Karpeev   ierr = PetscMalloc((*nsub)*sizeof(IS*),is);CHKERRQ(ierr);
1543af538404SDmitry Karpeev   ierr = PetscMalloc((*nsub)*sizeof(IS*),is_local);CHKERRQ(ierr);
1544af538404SDmitry Karpeev   s = 0;
1545af538404SDmitry Karpeev   ystart = 0;
1546af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1547af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1548af538404SDmitry 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);
1549af538404SDmitry Karpeev     /* Vertical domain limits with an overlap. */
1550eec7e3faSDmitry Karpeev     y[0][0] = PetscMax(ystart - overlap,0);
1551eec7e3faSDmitry Karpeev     y[0][1] = PetscMin(ystart + maxheight + overlap,N);
1552eec7e3faSDmitry Karpeev     /* Vertical domain limits without an overlap. */
1553eec7e3faSDmitry Karpeev     y[1][0] = ystart;
1554eec7e3faSDmitry Karpeev     y[1][1] = PetscMin(ystart + maxheight,N);
1555eec7e3faSDmitry Karpeev     xstart = 0;
1556af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1557af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1558af538404SDmitry 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);
1559af538404SDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1560eec7e3faSDmitry Karpeev       x[0][0]  = PetscMax(xstart - overlap,0);
1561eec7e3faSDmitry Karpeev       x[0][1]  = PetscMin(xstart + maxwidth + overlap,M);
1562eec7e3faSDmitry Karpeev       /* Horizontal domain limits without an overlap. */
1563eec7e3faSDmitry Karpeev       x[1][0] = xstart;
1564eec7e3faSDmitry Karpeev       x[1][1] = PetscMin(xstart+maxwidth,M);
15652bbb417fSDmitry Karpeev       /*
15662bbb417fSDmitry Karpeev 	 Determine whether this domain intersects this processor's ownership range of pc->pmat.
15672bbb417fSDmitry Karpeev 	 Do this twice: first for the domains with overlaps, and once without.
15682bbb417fSDmitry Karpeev 	 During the first pass create the subcommunicators, and use them on the second pass as well.
15692bbb417fSDmitry Karpeev       */
15702bbb417fSDmitry Karpeev       for(q = 0; q < 2; ++q) {
15712bbb417fSDmitry Karpeev 	/*
15722bbb417fSDmitry Karpeev 	  domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
15732bbb417fSDmitry Karpeev 	  according to whether the domain with an overlap or without is considered.
15742bbb417fSDmitry Karpeev 	*/
15752bbb417fSDmitry Karpeev 	xleft = x[q][0]; xright = x[q][1];
15762bbb417fSDmitry Karpeev 	ylow  = y[q][0]; yhigh  = y[q][1];
1577c3518bceSDmitry Karpeev         PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1578af538404SDmitry Karpeev 	nidx *= dof;
1579eec7e3faSDmitry Karpeev         n[q] = nidx;
15802bbb417fSDmitry Karpeev         /*
1581eec7e3faSDmitry Karpeev          Based on the counted number of indices in the local domain *with an overlap*,
1582af538404SDmitry Karpeev          construct a subcommunicator of all the processors supporting this domain.
1583eec7e3faSDmitry Karpeev          Observe that a domain with an overlap might have nontrivial local support,
1584eec7e3faSDmitry Karpeev          while the domain without an overlap might not.  Hence, the decision to participate
1585eec7e3faSDmitry Karpeev          in the subcommunicator must be based on the domain with an overlap.
15862bbb417fSDmitry Karpeev          */
15872bbb417fSDmitry Karpeev 	if (q == 0) {
15882bbb417fSDmitry Karpeev 	  if(nidx) {
15892bbb417fSDmitry Karpeev 	    color = 1;
1590d34fcf5fSBarry Smith 	  } else {
15912bbb417fSDmitry Karpeev 	    color = MPI_UNDEFINED;
15922bbb417fSDmitry Karpeev 	  }
15932bbb417fSDmitry Karpeev 	  ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr);
15942bbb417fSDmitry Karpeev 	}
1595af538404SDmitry Karpeev         /*
1596eec7e3faSDmitry Karpeev          Proceed only if the number of local indices *with an overlap* is nonzero.
1597af538404SDmitry Karpeev          */
1598eec7e3faSDmitry Karpeev         if (n[0]) {
1599af538404SDmitry Karpeev           if(q == 0) {
1600*6a4f0f73SDmitry Karpeev             xis = is;
1601af538404SDmitry Karpeev           }
1602af538404SDmitry Karpeev           if (q == 1) {
1603af538404SDmitry Karpeev             /*
1604eec7e3faSDmitry Karpeev              The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1605af538404SDmitry Karpeev              Moreover, if the overlap is zero, the two ISs are identical.
1606af538404SDmitry Karpeev              */
1607b1a0a8a3SJed Brown             if (overlap == 0) {
1608eec7e3faSDmitry Karpeev               (*is_local)[s] = (*is)[s];
16092bbb417fSDmitry Karpeev               ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr);
16102bbb417fSDmitry Karpeev               continue;
1611d34fcf5fSBarry Smith             } else {
1612*6a4f0f73SDmitry Karpeev               xis = is_local;
1613eec7e3faSDmitry Karpeev               subcomm = ((PetscObject)(*is)[s])->comm;
16142bbb417fSDmitry Karpeev             }
1615af538404SDmitry Karpeev           }/* if(q == 1) */
1616eec7e3faSDmitry Karpeev           idx = PETSC_NULL;
16172bbb417fSDmitry Karpeev 	  ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1618eec7e3faSDmitry Karpeev           if(nidx) {
16192bbb417fSDmitry Karpeev             k    = 0;
16202bbb417fSDmitry Karpeev             for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1621af538404SDmitry Karpeev               PetscInt x0 = (jj==ylow_loc)?xleft_loc:xleft;
1622af538404SDmitry Karpeev               PetscInt x1 = (jj==yhigh_loc-1)?xright_loc:xright;
1623af538404SDmitry Karpeev               kk = dof*(M*jj + x0);
16242bbb417fSDmitry Karpeev               for (ii=x0; ii<x1; ++ii) {
16252bbb417fSDmitry Karpeev                 for(d = 0; d < dof; ++d) {
16262bbb417fSDmitry Karpeev                   idx[k++] = kk++;
1627b1a0a8a3SJed Brown                 }
1628b1a0a8a3SJed Brown               }
1629b1a0a8a3SJed Brown             }
1630eec7e3faSDmitry Karpeev           }
1631*6a4f0f73SDmitry Karpeev 	  ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr);
1632eec7e3faSDmitry Karpeev 	}/* if(n[0]) */
16332bbb417fSDmitry Karpeev       }/* for(q = 0; q < 2; ++q) */
1634eec7e3faSDmitry Karpeev       if(n[0]) {
16352bbb417fSDmitry Karpeev         ++s;
1636b1a0a8a3SJed Brown       }
1637af538404SDmitry Karpeev       xstart += maxwidth;
1638af538404SDmitry Karpeev     }/* for(i = 0; i < Mdomains; ++i) */
16392bbb417fSDmitry Karpeev     ystart += maxheight;
1640af538404SDmitry Karpeev   }/* for(j = 0; j < Ndomains; ++j) */
1641b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1642b1a0a8a3SJed Brown }
1643b1a0a8a3SJed Brown 
1644b1a0a8a3SJed Brown #undef __FUNCT__
164506b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains"
1646b1a0a8a3SJed Brown /*@C
164706b43e7eSDmitry Karpeev     PCGASMGetSubdomains - Gets the subdomains supported on this processor
164806b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
1649b1a0a8a3SJed Brown 
1650b1a0a8a3SJed Brown     Not Collective
1651b1a0a8a3SJed Brown 
1652b1a0a8a3SJed Brown     Input Parameter:
1653b1a0a8a3SJed Brown .   pc - the preconditioner context
1654b1a0a8a3SJed Brown 
1655b1a0a8a3SJed Brown     Output Parameters:
1656b1a0a8a3SJed Brown +   n   - the number of subdomains for this processor (default value = 1)
1657*6a4f0f73SDmitry Karpeev .   iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be PETSC_NULL)
1658*6a4f0f73SDmitry Karpeev -   ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be PETSC_NULL)
1659b1a0a8a3SJed Brown 
1660b1a0a8a3SJed Brown 
1661b1a0a8a3SJed Brown     Notes:
1662*6a4f0f73SDmitry Karpeev     The user is responsible for destroying the ISs and freeing the returned arrays.
1663b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
1664b1a0a8a3SJed Brown 
1665b1a0a8a3SJed Brown     Level: advanced
1666b1a0a8a3SJed Brown 
166706b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz
1668b1a0a8a3SJed Brown 
1669*6a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
167006b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubmatrices()
1671b1a0a8a3SJed Brown @*/
1672*6a4f0f73SDmitry Karpeev PetscErrorCode  PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[])
1673b1a0a8a3SJed Brown {
1674f746d493SDmitry Karpeev   PC_GASM         *osm;
1675b1a0a8a3SJed Brown   PetscErrorCode ierr;
1676b1a0a8a3SJed Brown   PetscBool      match;
1677*6a4f0f73SDmitry Karpeev   PetscInt       i;
1678b1a0a8a3SJed Brown   PetscFunctionBegin;
1679b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1680251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1681*6a4f0f73SDmitry Karpeev   if (!match)
1682*6a4f0f73SDmitry Karpeev     SETERRQ2(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1683f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1684ab3e923fSDmitry Karpeev   if (n)  *n  = osm->n;
1685*6a4f0f73SDmitry Karpeev   if(iis) {
1686*6a4f0f73SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(IS), iis); CHKERRQ(ierr);
1687*6a4f0f73SDmitry Karpeev   }
1688*6a4f0f73SDmitry Karpeev   if(ois) {
1689*6a4f0f73SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(IS), ois); CHKERRQ(ierr);
1690*6a4f0f73SDmitry Karpeev   }
1691*6a4f0f73SDmitry Karpeev   if(iis || ois) {
1692*6a4f0f73SDmitry Karpeev     for(i = 0; i < osm->n; ++i) {
1693*6a4f0f73SDmitry Karpeev       if(iis) (*iis)[i] = osm->iis[i];
1694*6a4f0f73SDmitry Karpeev       if(ois) (*ois)[i] = osm->ois[i];
1695*6a4f0f73SDmitry Karpeev     }
1696b1a0a8a3SJed Brown   }
1697b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1698b1a0a8a3SJed Brown }
1699b1a0a8a3SJed Brown 
1700b1a0a8a3SJed Brown #undef __FUNCT__
170106b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices"
1702b1a0a8a3SJed Brown /*@C
170306b43e7eSDmitry Karpeev     PCGASMGetSubmatrices - Gets the local submatrices (for this processor
1704b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1705b1a0a8a3SJed Brown 
1706b1a0a8a3SJed Brown     Not Collective
1707b1a0a8a3SJed Brown 
1708b1a0a8a3SJed Brown     Input Parameter:
1709b1a0a8a3SJed Brown .   pc - the preconditioner context
1710b1a0a8a3SJed Brown 
1711b1a0a8a3SJed Brown     Output Parameters:
1712b1a0a8a3SJed Brown +   n   - the number of matrices for this processor (default value = 1)
1713b1a0a8a3SJed Brown -   mat - the matrices
1714b1a0a8a3SJed Brown 
171506b43e7eSDmitry Karpeev     Notes: matrices returned by this routine have the same communicators as the index sets (IS)
171606b43e7eSDmitry Karpeev            used to define subdomains in PCGASMSetSubdomains(), or PETSC_COMM_SELF, if the
1717*6a4f0f73SDmitry Karpeev            subdomains were defined using PCGASMSetTotalSubdomains().
1718b1a0a8a3SJed Brown     Level: advanced
1719b1a0a8a3SJed Brown 
1720f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi
1721b1a0a8a3SJed Brown 
1722*6a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomain(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
172306b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains()
1724b1a0a8a3SJed Brown @*/
172506b43e7eSDmitry Karpeev PetscErrorCode  PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1726b1a0a8a3SJed Brown {
1727f746d493SDmitry Karpeev   PC_GASM         *osm;
1728b1a0a8a3SJed Brown   PetscErrorCode ierr;
1729b1a0a8a3SJed Brown   PetscBool      match;
1730b1a0a8a3SJed Brown 
1731b1a0a8a3SJed Brown   PetscFunctionBegin;
1732b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1733b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1734b1a0a8a3SJed Brown   if (mat) PetscValidPointer(mat,3);
1735b1a0a8a3SJed Brown   if (!pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1736251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
173706b43e7eSDmitry Karpeev   if (!match) SETERRQ2(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1738f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1739ab3e923fSDmitry Karpeev   if (n)   *n   = osm->n;
1740b1a0a8a3SJed Brown   if (mat) *mat = osm->pmat;
174106b43e7eSDmitry Karpeev 
1742b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1743b1a0a8a3SJed Brown }
1744