xref: /petsc/src/ksp/pc/impls/gasm/gasm.c (revision 06b43e7e4d6b1e790b8e5a0a4e79ed97fab2d6a9)
1b1a0a8a3SJed Brown /*
2f746d493SDmitry Karpeev   This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation.
344fe18e5SDmitry Karpeev   In this version each processor may have any number of subdomains and a subdomain may live on multiple
4f746d493SDmitry Karpeev   processors.
5b1a0a8a3SJed Brown 
6*06b43e7eSDmitry Karpeev        N    - total number of local subdomains on all processors  (set in PCGASMSetTotalBlockCount() or calculated in PCSetUp_GASM())
7*06b43e7eSDmitry Karpeev        n    - actual number of local subdomains on this processor (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalBlockCount())
8*06b43e7eSDmitry Karpeev        nmax - maximum number of local subdomains per processor    (calculated in PCGASMSetTotalBlockCount() or in PCSetUp_GASM())
9b1a0a8a3SJed Brown */
10b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>     /*I "petscpc.h" I*/
11b1a0a8a3SJed Brown 
12b1a0a8a3SJed Brown typedef struct {
13f746d493SDmitry Karpeev   PetscInt   N,n,nmax;
14b1a0a8a3SJed Brown   PetscInt   overlap;             /* overlap requested by user */
15b1a0a8a3SJed Brown   KSP        *ksp;                /* linear solvers for each block */
16f746d493SDmitry Karpeev   Vec        gx,gy;               /* Merged work vectors */
17f746d493SDmitry Karpeev   Vec        *x,*y;               /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */
18f746d493SDmitry Karpeev   IS         gis, gis_local;      /* merged ISs */
19f746d493SDmitry Karpeev   VecScatter grestriction;        /* merged restriction */
20f746d493SDmitry Karpeev   VecScatter gprolongation;       /* merged prolongation */
21b1a0a8a3SJed Brown   IS         *is;                 /* index set that defines each overlapping subdomain */
22f746d493SDmitry Karpeev   IS         *is_local;           /* index set that defines each local subdomain (same as subdomain with the overlap removed); may be NULL */
23f746d493SDmitry Karpeev   Mat        *pmat;               /* subdomain block matrices */
24f746d493SDmitry Karpeev   PCGASMType  type;               /* use reduced interpolation, restriction or both */
25b1a0a8a3SJed Brown   PetscBool  type_set;            /* if user set this value (so won't change it for symmetric problems) */
26b1a0a8a3SJed Brown   PetscBool  same_local_solves;   /* flag indicating whether all local solvers are same */
27b1a0a8a3SJed Brown   PetscBool  sort_indices;        /* flag to sort subdomain indices */
28f746d493SDmitry Karpeev } PC_GASM;
29b1a0a8a3SJed Brown 
30b1a0a8a3SJed Brown #undef __FUNCT__
31*06b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSubdomainView_Private"
32*06b43e7eSDmitry Karpeev static PetscErrorCode  PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer)
33af538404SDmitry Karpeev {
34af538404SDmitry Karpeev   PC_GASM        *osm  = (PC_GASM*)pc->data;
35*06b43e7eSDmitry Karpeev   PetscInt       j,nidx;
36af538404SDmitry Karpeev   const PetscInt *idx;
37*06b43e7eSDmitry Karpeev   PetscViewer    sviewer;
38af538404SDmitry Karpeev   char           *cidx;
39af538404SDmitry Karpeev   PetscErrorCode ierr;
40af538404SDmitry Karpeev 
41af538404SDmitry Karpeev   PetscFunctionBegin;
42*06b43e7eSDmitry 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);
43af538404SDmitry Karpeev   ierr = ISGetLocalSize(osm->is[i], &nidx);CHKERRQ(ierr);
44eec7e3faSDmitry Karpeev   /*
45eec7e3faSDmitry Karpeev    No more than 15 characters per index plus a space.
46eec7e3faSDmitry Karpeev    PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
47eec7e3faSDmitry Karpeev    in case nidx == 0. That will take care of the space for the trailing '\0' as well.
48eec7e3faSDmitry Karpeev    For nidx == 0, the whole string 16 '\0'.
49eec7e3faSDmitry Karpeev    */
50eec7e3faSDmitry Karpeev   ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr);
51*06b43e7eSDmitry Karpeev   ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer); CHKERRQ(ierr);
52af538404SDmitry Karpeev   ierr = ISGetIndices(osm->is[i], &idx);CHKERRQ(ierr);
53af538404SDmitry Karpeev   for(j = 0; j < nidx; ++j) {
54af538404SDmitry Karpeev     ierr = PetscViewerStringSPrintf(sviewer,"%D ", idx[j]); CHKERRQ(ierr);
55af538404SDmitry Karpeev   }
566bf464f9SBarry Smith   ierr = PetscViewerDestroy(&sviewer); CHKERRQ(ierr);
57af538404SDmitry Karpeev   ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr);
58*06b43e7eSDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "Subdomain with overlap:\n");CHKERRQ(ierr);
59*06b43e7eSDmitry Karpeev   ierr = PetscViewerFlush(viewer);                                   CHKERRQ(ierr);
607b23a99aSBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
61af538404SDmitry Karpeev   ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
62af538404SDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
637b23a99aSBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
64af538404SDmitry Karpeev   ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
65*06b43e7eSDmitry Karpeev   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
66af538404SDmitry Karpeev   ierr = PetscFree(cidx);CHKERRQ(ierr);
67af538404SDmitry Karpeev   if (osm->is_local) {
68af538404SDmitry Karpeev     ierr = ISGetLocalSize(osm->is_local[i], &nidx);CHKERRQ(ierr);
69eec7e3faSDmitry Karpeev     /*
70eec7e3faSDmitry Karpeev      No more than 15 characters per index plus a space.
71eec7e3faSDmitry Karpeev      PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
72eec7e3faSDmitry Karpeev      in case nidx == 0. That will take care of the space for the trailing '\0' as well.
73eec7e3faSDmitry Karpeev      For nidx == 0, the whole string 16 '\0'.
74eec7e3faSDmitry Karpeev      */
75eec7e3faSDmitry Karpeev     ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr);
76*06b43e7eSDmitry Karpeev     ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer); CHKERRQ(ierr);
77af538404SDmitry Karpeev     ierr = ISGetIndices(osm->is_local[i], &idx);CHKERRQ(ierr);
78af538404SDmitry Karpeev     for(j = 0; j < nidx; ++j) {
79af538404SDmitry Karpeev       ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr);
80af538404SDmitry Karpeev     }
81af538404SDmitry Karpeev     ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
82*06b43e7eSDmitry Karpeev     ierr = PetscViewerDestroy(&sviewer); CHKERRQ(ierr);
83*06b43e7eSDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer, "Subdomain without overlap:\n");CHKERRQ(ierr);
84af538404SDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr);
85*06b43e7eSDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer, "\n");CHKERRQ(ierr);
86af538404SDmitry Karpeev     ierr = PetscFree(cidx);CHKERRQ(ierr);
87af538404SDmitry Karpeev   }
88*06b43e7eSDmitry Karpeev   PetscFunctionReturn(0);
89*06b43e7eSDmitry Karpeev }
90*06b43e7eSDmitry Karpeev 
91*06b43e7eSDmitry Karpeev #undef __FUNCT__
92*06b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains"
93*06b43e7eSDmitry Karpeev static PetscErrorCode  PCGASMPrintSubdomains(PC pc)
94*06b43e7eSDmitry Karpeev {
95*06b43e7eSDmitry Karpeev   PC_GASM        *osm  = (PC_GASM*)pc->data;
96*06b43e7eSDmitry Karpeev   const char     *prefix;
97*06b43e7eSDmitry Karpeev   char           fname[PETSC_MAX_PATH_LEN+1];
98*06b43e7eSDmitry Karpeev   PetscViewer    viewer;
99*06b43e7eSDmitry Karpeev   PetscInt       i;
100*06b43e7eSDmitry Karpeev   PetscBool      found;
101*06b43e7eSDmitry Karpeev   PetscErrorCode ierr;
102*06b43e7eSDmitry Karpeev 
103*06b43e7eSDmitry Karpeev   PetscFunctionBegin;
104*06b43e7eSDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
105*06b43e7eSDmitry Karpeev   ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr);
106*06b43e7eSDmitry Karpeev   if (!found) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
107*06b43e7eSDmitry Karpeev   for (i=0;i<osm->n;++i) {
108*06b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIOpen(((PetscObject)((osm->is)[i]))->comm,fname,&viewer);CHKERRQ(ierr);
109*06b43e7eSDmitry Karpeev     ierr = PCGASMSubdomainView_Private(pc, i, viewer); CHKERRQ(ierr);
110af538404SDmitry Karpeev   }
111af538404SDmitry Karpeev   PetscFunctionReturn(0);
112af538404SDmitry Karpeev }
113af538404SDmitry Karpeev 
114af538404SDmitry Karpeev 
115af538404SDmitry Karpeev #undef __FUNCT__
116f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM"
117f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer)
118b1a0a8a3SJed Brown {
119f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
120af538404SDmitry Karpeev   const char     *prefix;
121b1a0a8a3SJed Brown   PetscErrorCode ierr;
122af538404SDmitry Karpeev   PetscMPIInt    rank, size;
123b1a0a8a3SJed Brown   PetscInt       i,bsz;
124*06b43e7eSDmitry Karpeev   PetscBool      iascii,view_subdomains=PETSC_FALSE;
125b1a0a8a3SJed Brown   PetscViewer    sviewer;
126*06b43e7eSDmitry Karpeev   PetscInt       count, l, gcount, *numbering, *permutation;
127*06b43e7eSDmitry Karpeev   char overlap[256]     = "no additional overlap requested";
128*06b43e7eSDmitry Karpeev   char gsubdomains[256] = "unknown number of global subdomains";
129*06b43e7eSDmitry Karpeev   char lsubdomains[256] = "unknown number of local  subdomains";
130*06b43e7eSDmitry Karpeev   char msubdomains[256] = "unknown max number of local subdomains";
131b1a0a8a3SJed Brown   PetscFunctionBegin;
132af538404SDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm, &size);CHKERRQ(ierr);
133af538404SDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm, &rank);CHKERRQ(ierr);
134*06b43e7eSDmitry Karpeev 
135*06b43e7eSDmitry Karpeev 
136*06b43e7eSDmitry Karpeev   ierr = PetscMalloc2(osm->n, PetscInt, &permutation, osm->n, PetscInt, &numbering); CHKERRQ(ierr);
137*06b43e7eSDmitry Karpeev   for(i = 0; i < osm->n; ++i) permutation[i] = i;
138*06b43e7eSDmitry Karpeev   ierr = PetscObjectListGetGlobalNumbering(((PetscObject)pc)->comm, osm->n, (PetscObject*)osm->is, &gcount, numbering); CHKERRQ(ierr);
139*06b43e7eSDmitry Karpeev   ierr = PetscSortIntWithPermutation(osm->n, numbering, permutation); CHKERRQ(ierr);
140*06b43e7eSDmitry Karpeev 
141*06b43e7eSDmitry Karpeev   if(osm->overlap >= 0) {
142*06b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %D",osm->overlap);        CHKERRQ(ierr);
143*06b43e7eSDmitry Karpeev   }
144*06b43e7eSDmitry Karpeev   ierr = PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "number of global subdomains = %D",gcount);      CHKERRQ(ierr);
145*06b43e7eSDmitry Karpeev   if(osm->N > 0) {
146*06b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(lsubdomains, sizeof(gsubdomains), "number of local subdomains = %D",osm->N);     CHKERRQ(ierr);
147*06b43e7eSDmitry Karpeev   }
148*06b43e7eSDmitry Karpeev   if(osm->nmax > 0){
149*06b43e7eSDmitry Karpeev     ierr = PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr);
150*06b43e7eSDmitry Karpeev   }
151*06b43e7eSDmitry Karpeev 
152af538404SDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
153*06b43e7eSDmitry Karpeev   ierr = PetscOptionsGetBool(prefix,"-pc_gasm_view_subdomains",&view_subdomains,PETSC_NULL);CHKERRQ(ierr);
154*06b43e7eSDmitry Karpeev 
155*06b43e7eSDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii); CHKERRQ(ierr);
156b1a0a8a3SJed Brown   if (iascii) {
157*06b43e7eSDmitry Karpeev     /*
158*06b43e7eSDmitry Karpeev      Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer:
159*06b43e7eSDmitry Karpeev      the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed
160*06b43e7eSDmitry Karpeev      collectively on the comm.
161*06b43e7eSDmitry Karpeev      */
162*06b43e7eSDmitry Karpeev     ierr = PetscObjectName((PetscObject)viewer);                  CHKERRQ(ierr);
163*06b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Generalized additive Schwarz:\n"); CHKERRQ(ierr);
164*06b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Restriction/interpolation type: %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr);
165*06b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",overlap);    CHKERRQ(ierr);
166*06b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",gsubdomains);CHKERRQ(ierr);
167*06b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",lsubdomains);CHKERRQ(ierr);
168*06b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"%s\n",msubdomains);CHKERRQ(ierr);
1697b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
170eec7e3faSDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d:%d] number of locally-supported subdomains = %D\n",(int)rank,(int)size,osm->n);CHKERRQ(ierr);
171af538404SDmitry Karpeev     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1727b23a99aSBarry Smith     ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
173*06b43e7eSDmitry Karpeev     /* Cannot take advantage of osm->same_local_solves without a global numbering of subdomains. */
174*06b43e7eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Subdomain solver info is as follows:\n");CHKERRQ(ierr);
175b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
176b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
177*06b43e7eSDmitry Karpeev     /* Make sure that everybody waits for the banner to be printed. */
178*06b43e7eSDmitry Karpeev     ierr = MPI_Barrier(((PetscObject)viewer)->comm); CHKERRQ(ierr);
179*06b43e7eSDmitry Karpeev     l = 0;
180*06b43e7eSDmitry Karpeev     for(count = 0; count < gcount; ++count) {
181*06b43e7eSDmitry Karpeev       /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */
182*06b43e7eSDmitry Karpeev       PetscMPIInt srank, ssize;
183*06b43e7eSDmitry Karpeev       if(l<osm->n){
184*06b43e7eSDmitry Karpeev         PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */
185*06b43e7eSDmitry Karpeev         if(numbering[d] == count) {
186*06b43e7eSDmitry Karpeev           ierr = MPI_Comm_size(((PetscObject)osm->is[d])->comm, &ssize); CHKERRQ(ierr);
187*06b43e7eSDmitry Karpeev           ierr = MPI_Comm_rank(((PetscObject)osm->is[d])->comm, &srank); CHKERRQ(ierr);
188*06b43e7eSDmitry Karpeev           ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->is[d])->comm, &sviewer);CHKERRQ(ierr);
189*06b43e7eSDmitry Karpeev           ierr = ISGetLocalSize(osm->is[d],&bsz);CHKERRQ(ierr);
1907b23a99aSBarry Smith           ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_TRUE);CHKERRQ(ierr);
191*06b43e7eSDmitry 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);
192*06b43e7eSDmitry Karpeev           ierr = PetscViewerFlush(viewer); CHKERRQ(ierr);
193*06b43e7eSDmitry Karpeev           if(view_subdomains) {
194*06b43e7eSDmitry Karpeev             ierr = PCGASMSubdomainView_Private(pc,d,sviewer); CHKERRQ(ierr);
195*06b43e7eSDmitry Karpeev           }
196*06b43e7eSDmitry Karpeev           ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr);
197*06b43e7eSDmitry Karpeev           ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
1987b23a99aSBarry Smith           ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr);
1997b23a99aSBarry Smith           ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_FALSE);CHKERRQ(ierr);
200*06b43e7eSDmitry Karpeev           ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->is[d])->comm, &sviewer);CHKERRQ(ierr);
201*06b43e7eSDmitry Karpeev           ++l;
202b1a0a8a3SJed Brown         }
203*06b43e7eSDmitry Karpeev       }
204*06b43e7eSDmitry Karpeev       ierr = MPI_Barrier(((PetscObject)pc)->comm); CHKERRQ(ierr);
205b1a0a8a3SJed Brown     }
206b1a0a8a3SJed Brown     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
207b1a0a8a3SJed Brown   } else {
208f746d493SDmitry Karpeev     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCGASM",((PetscObject)viewer)->type_name);
209b1a0a8a3SJed Brown   }
210*06b43e7eSDmitry Karpeev   ierr = PetscFree2(permutation,numbering); CHKERRQ(ierr);
211b1a0a8a3SJed Brown   PetscFunctionReturn(0);
212b1a0a8a3SJed Brown }
213b1a0a8a3SJed Brown 
214af538404SDmitry Karpeev 
215af538404SDmitry Karpeev 
216af538404SDmitry Karpeev 
217af538404SDmitry Karpeev 
218b1a0a8a3SJed Brown #undef __FUNCT__
219f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM"
220f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc)
221b1a0a8a3SJed Brown {
222f746d493SDmitry Karpeev   PC_GASM         *osm  = (PC_GASM*)pc->data;
223b1a0a8a3SJed Brown   PetscErrorCode ierr;
224b1a0a8a3SJed Brown   PetscBool      symset,flg;
22588389755SJed Brown   PetscInt       i;
226c10bc1a1SDmitry Karpeev   PetscMPIInt    rank, size;
227b1a0a8a3SJed Brown   MatReuse       scall = MAT_REUSE_MATRIX;
228b1a0a8a3SJed Brown   KSP            ksp;
229b1a0a8a3SJed Brown   PC             subpc;
230b1a0a8a3SJed Brown   const char     *prefix,*pprefix;
231f746d493SDmitry Karpeev   PetscInt       dn;       /* Number of indices in a single subdomain assigned to this processor. */
232f746d493SDmitry Karpeev   const PetscInt *didx;    /* Indices from a single subdomain assigned to this processor. */
233f746d493SDmitry Karpeev   PetscInt       ddn;      /* Number of indices in all subdomains assigned to this processor. */
234f746d493SDmitry Karpeev   PetscInt       *ddidx;   /* Indices of all subdomains assigned to this processor. */
235f746d493SDmitry Karpeev   IS             gid;      /* Identity IS of the size of all subdomains assigned to this processor. */
236f746d493SDmitry Karpeev   Vec            x,y;
237f746d493SDmitry Karpeev   PetscScalar    *gxarray, *gyarray;
2386d98aee9SDmitry Karpeev   PetscInt       gfirst, glast;
239fdc48646SDmitry Karpeev   DM             *domain_dm = PETSC_NULL;
240b1a0a8a3SJed Brown 
241b1a0a8a3SJed Brown   PetscFunctionBegin;
242c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
243c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
244b1a0a8a3SJed Brown   if (!pc->setupcalled) {
245b1a0a8a3SJed Brown 
246b1a0a8a3SJed Brown     if (!osm->type_set) {
247b1a0a8a3SJed Brown       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
248f746d493SDmitry Karpeev       if (symset && flg) { osm->type = PC_GASM_BASIC; }
249b1a0a8a3SJed Brown     }
250b1a0a8a3SJed Brown 
251*06b43e7eSDmitry Karpeev     /*
252*06b43e7eSDmitry Karpeev      If subdomains have been set, then the local number of subdomains, osm->n, is NOT PETSC_DECIDE and is at least 1.
253*06b43e7eSDmitry 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.
254*06b43e7eSDmitry Karpeev      */
255*06b43e7eSDmitry Karpeev     if (osm->n == PETSC_DECIDE) {
25669ca1f37SDmitry Karpeev       /* no subdomains given */
25769ca1f37SDmitry Karpeev       /* try pc->dm first */
258fdc48646SDmitry Karpeev       if(pc->dm) {
259fdc48646SDmitry Karpeev         char      ddm_name[1024];
260fdc48646SDmitry Karpeev         DM        ddm;
261fdc48646SDmitry Karpeev         PetscBool flg;
262fdc48646SDmitry Karpeev         PetscInt     num_domains, d;
263fdc48646SDmitry Karpeev         char         **domain_names;
264fdc48646SDmitry Karpeev         IS           *domain_is;
265fdc48646SDmitry Karpeev         /* Allow the user to request a decomposition DM by name */
266fdc48646SDmitry Karpeev         ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr);
267*06b43e7eSDmitry Karpeev         ierr = PetscOptionsString("-pc_gasm_decomposition","Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr);
268fdc48646SDmitry Karpeev         if(flg) {
269fdc48646SDmitry Karpeev           ierr = DMCreateDecompositionDM(pc->dm, ddm_name, &ddm); CHKERRQ(ierr);
270fdc48646SDmitry Karpeev           if(!ddm) {
271fdc48646SDmitry Karpeev             SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown DM decomposition name %s", ddm_name);
272fdc48646SDmitry Karpeev           }
273fdc48646SDmitry Karpeev           ierr = PetscInfo(pc,"Using decomposition DM defined using options database\n");CHKERRQ(ierr);
274fdc48646SDmitry Karpeev           ierr = PCSetDM(pc,ddm); CHKERRQ(ierr);
275fdc48646SDmitry Karpeev         }
276fdc48646SDmitry Karpeev         ierr = DMCreateDecomposition(pc->dm, &num_domains, &domain_names, &domain_is, &domain_dm);    CHKERRQ(ierr);
27769ca1f37SDmitry Karpeev         if(num_domains) {
278*06b43e7eSDmitry Karpeev           ierr = PCGASMSetSubdomains(pc, num_domains, domain_is, PETSC_NULL);CHKERRQ(ierr);
27969ca1f37SDmitry Karpeev         }
280fdc48646SDmitry Karpeev         for(d = 0; d < num_domains; ++d) {
281fdc48646SDmitry Karpeev           ierr = PetscFree(domain_names[d]); CHKERRQ(ierr);
282fdc48646SDmitry Karpeev           ierr = ISDestroy(&domain_is[d]);   CHKERRQ(ierr);
283fdc48646SDmitry Karpeev         }
284fdc48646SDmitry Karpeev         ierr = PetscFree(domain_names);CHKERRQ(ierr);
285fdc48646SDmitry Karpeev         ierr = PetscFree(domain_is);CHKERRQ(ierr);
286fdc48646SDmitry Karpeev       }
287*06b43e7eSDmitry Karpeev       if (osm->n == PETSC_DECIDE) { /* still no subdomains; use one per processor */
28844fe18e5SDmitry Karpeev         osm->nmax = osm->n = 1;
289b1a0a8a3SJed Brown         ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
290f746d493SDmitry Karpeev         osm->N = size;
291fdc48646SDmitry Karpeev       }
292*06b43e7eSDmitry Karpeev     }
293*06b43e7eSDmitry Karpeev     if (osm->N == PETSC_DECIDE) {
294b1a0a8a3SJed Brown       PetscInt inwork[2], outwork[2];
295f746d493SDmitry Karpeev       /* determine global number of subdomains and the max number of local subdomains */
296f746d493SDmitry Karpeev       inwork[0] = inwork[1] = osm->n;
297b1a0a8a3SJed Brown       ierr = MPI_Allreduce(inwork,outwork,1,MPIU_2INT,PetscMaxSum_Op,((PetscObject)pc)->comm);CHKERRQ(ierr);
298f746d493SDmitry Karpeev       osm->nmax = outwork[0];
299f746d493SDmitry Karpeev       osm->N    = outwork[1];
300b1a0a8a3SJed Brown     }
301*06b43e7eSDmitry Karpeev     if (!osm->is){
302*06b43e7eSDmitry Karpeev       /* The local number of subdomains was set in PCGASMSetTotalBlockCount() or PCGASMSetSubdomains(),
303*06b43e7eSDmitry Karpeev        but the actual subdomains have not been supplied (in PCGASMSetSubdomains()).
304*06b43e7eSDmitry Karpeev        We create the requisite number of subdomains on PETSC_COMM_SELF. */
305*06b43e7eSDmitry Karpeev       ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->is,&osm->is_local);CHKERRQ(ierr);
306b1a0a8a3SJed Brown     }
30793c1ec32SDmitry Karpeev     if (!osm->is_local) {
30893c1ec32SDmitry Karpeev       /*
30993c1ec32SDmitry Karpeev 	 This indicates that osm->is should define a nonoverlapping decomposition
310*06b43e7eSDmitry Karpeev 	 (there is no way to really guarantee that if subdomains are set by the user through PCGASMSetSubdomains,
311*06b43e7eSDmitry Karpeev 	  but the assumption is that either the user does the right thing, or subdomains in osm->is have been created
312*06b43e7eSDmitry Karpeev 	  via PCGASMCreateLocalSubdomains(), which guarantees a nonoverlapping decomposition).
31393c1ec32SDmitry Karpeev 	 Therefore, osm->is will be used to define osm->is_local.
31493c1ec32SDmitry Karpeev 	 If a nonzero overlap has been requested by the user, then osm->is will be expanded and will overlap,
31593c1ec32SDmitry Karpeev 	 so osm->is_local should obtain a copy of osm->is while they are still (presumably) nonoverlapping.
31693c1ec32SDmitry Karpeev 	 Otherwise (no overlap has been requested), osm->is_local are simply aliases for osm->is.
31793c1ec32SDmitry Karpeev       */
318f746d493SDmitry Karpeev       ierr = PetscMalloc(osm->n*sizeof(IS),&osm->is_local);CHKERRQ(ierr);
319f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
320b1a0a8a3SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
321b1a0a8a3SJed Brown           ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr);
322b1a0a8a3SJed Brown           ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr);
323b1a0a8a3SJed Brown         } else {
324b1a0a8a3SJed Brown           ierr = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr);
325b1a0a8a3SJed Brown           osm->is_local[i] = osm->is[i];
326b1a0a8a3SJed Brown         }
327b1a0a8a3SJed Brown       }
328b1a0a8a3SJed Brown     }
32993c1ec32SDmitry Karpeev     /* Beyond this point osm->is_local is not null. */
33093c1ec32SDmitry Karpeev     if (osm->overlap > 0) {
33193c1ec32SDmitry Karpeev       /* Extend the "overlapping" regions by a number of steps */
33293c1ec32SDmitry Karpeev       ierr = MatIncreaseOverlap(pc->pmat,osm->n,osm->is,osm->overlap);CHKERRQ(ierr);
33393c1ec32SDmitry Karpeev     }
334b1a0a8a3SJed Brown     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
335b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
336f746d493SDmitry Karpeev     ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&flg,PETSC_NULL);CHKERRQ(ierr);
337f746d493SDmitry Karpeev     if (flg) { ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); }
338b1a0a8a3SJed Brown 
339b1a0a8a3SJed Brown     if (osm->sort_indices) {
340f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
341b1a0a8a3SJed Brown         ierr = ISSort(osm->is[i]);CHKERRQ(ierr);
342b1a0a8a3SJed Brown 	ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr);
343b1a0a8a3SJed Brown       }
344b1a0a8a3SJed Brown     }
345f746d493SDmitry Karpeev     /* Merge the ISs, create merged vectors and scatter contexts. */
34693c1ec32SDmitry Karpeev     /* Restriction ISs. */
347f746d493SDmitry Karpeev     ddn = 0;
348f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
349f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn);CHKERRQ(ierr);
350f746d493SDmitry Karpeev       ddn += dn;
351f746d493SDmitry Karpeev     }
352f746d493SDmitry Karpeev     ierr = PetscMalloc(ddn*sizeof(PetscInt), &ddidx);CHKERRQ(ierr);
353f746d493SDmitry Karpeev     ddn = 0;
354f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
355f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn);CHKERRQ(ierr);
356f746d493SDmitry Karpeev       ierr = ISGetIndices(osm->is[i],&didx);CHKERRQ(ierr);
357f746d493SDmitry Karpeev       ierr = PetscMemcpy(ddidx+ddn, didx, sizeof(PetscInt)*dn);CHKERRQ(ierr);
358f746d493SDmitry Karpeev       ierr = ISRestoreIndices(osm->is[i], &didx);CHKERRQ(ierr);
35993c1ec32SDmitry Karpeev       ddn += dn;
360f746d493SDmitry Karpeev     }
361f746d493SDmitry Karpeev     ierr = ISCreateGeneral(((PetscObject)(pc))->comm, ddn, ddidx, PETSC_OWN_POINTER, &osm->gis);CHKERRQ(ierr);
362f746d493SDmitry Karpeev     ierr = MatGetVecs(pc->pmat,&x,&y);CHKERRQ(ierr);
363f746d493SDmitry Karpeev     ierr = VecCreateMPI(((PetscObject)pc)->comm, ddn, PETSC_DECIDE, &osm->gx);CHKERRQ(ierr);
364f746d493SDmitry Karpeev     ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr);
3656d98aee9SDmitry Karpeev     ierr = VecGetOwnershipRange(osm->gx, &gfirst, &glast);CHKERRQ(ierr);
3666d98aee9SDmitry Karpeev     ierr = ISCreateStride(((PetscObject)pc)->comm,ddn,gfirst,1, &gid);CHKERRQ(ierr);
367f746d493SDmitry Karpeev     ierr = VecScatterCreate(x,osm->gis,osm->gx,gid, &(osm->grestriction));CHKERRQ(ierr);
368fcfd50ebSBarry Smith     ierr = ISDestroy(&gid);CHKERRQ(ierr);
36993c1ec32SDmitry Karpeev     /* Prolongation ISs */
37081a5c4aaSDmitry Karpeev     { PetscInt       dn_local;       /* Number of indices in the local part of a single domain assigned to this processor. */
371f746d493SDmitry Karpeev       const PetscInt *didx_local;    /* Global indices from the local part of a single domain assigned to this processor. */
372f746d493SDmitry Karpeev       PetscInt       ddn_local;      /* Number of indices in the local part of the disjoint union all domains assigned to this processor. */
373f746d493SDmitry Karpeev       PetscInt       *ddidx_local;   /* Global indices of the local part of the disjoint union of all domains assigned to this processor. */
374f746d493SDmitry Karpeev       /**/
37593c1ec32SDmitry Karpeev       ISLocalToGlobalMapping ltog;          /* Map from global to local indices on the disjoint union of subdomains: "local" ind's run from 0 to ddn-1. */
37693c1ec32SDmitry Karpeev       PetscInt              *ddidx_llocal;  /* Mapped local indices of the disjoint union of local parts of subdomains. */
377f746d493SDmitry Karpeev       PetscInt               ddn_llocal;    /* Number of indices in ddidx_llocal; must equal ddn_local, or else gis_local is not a sub-IS of gis. */
378f746d493SDmitry Karpeev       IS                     gis_llocal;    /* IS with ddidx_llocal indices. */
379f746d493SDmitry Karpeev       PetscInt               j;
380f746d493SDmitry Karpeev       ddn_local = 0;
381f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
382f746d493SDmitry Karpeev 	ierr = ISGetLocalSize(osm->is_local[i],&dn_local);CHKERRQ(ierr);
383f746d493SDmitry Karpeev 	ddn_local += dn_local;
384f746d493SDmitry Karpeev       }
385f746d493SDmitry Karpeev       ierr = PetscMalloc(ddn_local*sizeof(PetscInt), &ddidx_local);CHKERRQ(ierr);
386f746d493SDmitry Karpeev       ddn_local = 0;
387f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
388f746d493SDmitry Karpeev 	ierr = ISGetLocalSize(osm->is_local[i],&dn_local);CHKERRQ(ierr);
389f746d493SDmitry Karpeev 	ierr = ISGetIndices(osm->is_local[i],&didx_local);CHKERRQ(ierr);
390f746d493SDmitry Karpeev 	ierr = PetscMemcpy(ddidx_local+ddn_local, didx_local, sizeof(PetscInt)*dn_local);CHKERRQ(ierr);
391f746d493SDmitry Karpeev 	ierr = ISRestoreIndices(osm->is_local[i], &didx_local);CHKERRQ(ierr);
392f746d493SDmitry Karpeev 	ddn_local += dn_local;
393f746d493SDmitry Karpeev       }
394ab3e923fSDmitry Karpeev       ierr = PetscMalloc(sizeof(PetscInt)*ddn_local, &ddidx_llocal);CHKERRQ(ierr);
395f746d493SDmitry Karpeev       ierr = ISCreateGeneral(((PetscObject)pc)->comm, ddn_local, ddidx_local, PETSC_OWN_POINTER, &(osm->gis_local));CHKERRQ(ierr);
396f746d493SDmitry Karpeev       ierr = ISLocalToGlobalMappingCreateIS(osm->gis,&ltog);CHKERRQ(ierr);
39781a5c4aaSDmitry Karpeev       ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,ddn_local,ddidx_local,&ddn_llocal,ddidx_llocal);CHKERRQ(ierr);
3986bf464f9SBarry Smith       ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
39993c1ec32SDmitry Karpeev       if (ddn_llocal != ddn_local) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gis_local contains %D indices outside of gis", ddn_llocal - ddn_local);
400f746d493SDmitry Karpeev       /* Now convert these localized indices into the global indices into the merged output vector. */
4016d98aee9SDmitry Karpeev       ierr = VecGetOwnershipRange(osm->gy, &gfirst, &glast);CHKERRQ(ierr);
402f746d493SDmitry Karpeev       for(j=0; j < ddn_llocal; ++j) {
4036d98aee9SDmitry Karpeev 	ddidx_llocal[j] += gfirst;
404f746d493SDmitry Karpeev       }
405f746d493SDmitry Karpeev       ierr = ISCreateGeneral(PETSC_COMM_SELF,ddn_llocal,ddidx_llocal,PETSC_OWN_POINTER,&gis_llocal);CHKERRQ(ierr);
406f746d493SDmitry Karpeev       ierr = VecScatterCreate(y,osm->gis_local,osm->gy,gis_llocal,&osm->gprolongation);CHKERRQ(ierr);
407fcfd50ebSBarry Smith       ierr = ISDestroy(&gis_llocal);CHKERRQ(ierr);
408b1a0a8a3SJed Brown     }
409f746d493SDmitry Karpeev     /* Create the subdomain work vectors. */
410f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->x);CHKERRQ(ierr);
411f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->y);CHKERRQ(ierr);
4126d98aee9SDmitry Karpeev     ierr = VecGetOwnershipRange(osm->gx, &gfirst, &glast);CHKERRQ(ierr);
413f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr);
414f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr);
415f746d493SDmitry Karpeev     for (i=0, ddn=0; i<osm->n; ++i, ddn += dn) {
41686b96d36SDmitry Karpeev       PetscInt dN;
417f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn);CHKERRQ(ierr);
41886b96d36SDmitry Karpeev       ierr = ISGetSize(osm->is[i],&dN);CHKERRQ(ierr);
419778a2246SBarry Smith       ierr = VecCreateMPIWithArray(((PetscObject)(osm->is[i]))->comm,1,dn,dN,gxarray+ddn,&osm->x[i]);CHKERRQ(ierr);
420778a2246SBarry Smith       ierr = VecCreateMPIWithArray(((PetscObject)(osm->is[i]))->comm,1,dn,dN,gyarray+ddn,&osm->y[i]);CHKERRQ(ierr);
421b1a0a8a3SJed Brown     }
422f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr);
423f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr);
4246bf464f9SBarry Smith     ierr = VecDestroy(&x);CHKERRQ(ierr);
4256bf464f9SBarry Smith     ierr = VecDestroy(&y);CHKERRQ(ierr);
426b1a0a8a3SJed Brown     /* Create the local solvers */
427f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(KSP *),&osm->ksp);CHKERRQ(ierr);
42844fe18e5SDmitry Karpeev     for (i=0; i<osm->n; i++) { /* KSPs are local */
42986b96d36SDmitry Karpeev       ierr = KSPCreate(((PetscObject)(osm->is[i]))->comm,&ksp);CHKERRQ(ierr);
430b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr);
431b1a0a8a3SJed Brown       ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
432b1a0a8a3SJed Brown       ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
433b1a0a8a3SJed Brown       ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
434b1a0a8a3SJed Brown       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
435b1a0a8a3SJed Brown       ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
436b1a0a8a3SJed Brown       ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
437b1a0a8a3SJed Brown       osm->ksp[i] = ksp;
438b1a0a8a3SJed Brown     }
439b1a0a8a3SJed Brown     scall = MAT_INITIAL_MATRIX;
440b1a0a8a3SJed Brown 
441b1a0a8a3SJed Brown   } else {
442b1a0a8a3SJed Brown     /*
443b1a0a8a3SJed Brown        Destroy the blocks from the previous iteration
444b1a0a8a3SJed Brown     */
445b1a0a8a3SJed Brown     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
446f746d493SDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
447b1a0a8a3SJed Brown       scall = MAT_INITIAL_MATRIX;
448b1a0a8a3SJed Brown     }
449b1a0a8a3SJed Brown   }
450b1a0a8a3SJed Brown 
451b1a0a8a3SJed Brown   /*
452f746d493SDmitry Karpeev      Extract out the submatrices.
453b1a0a8a3SJed Brown   */
45481a5c4aaSDmitry Karpeev   if(size > 1) {
4556d42056aSDmitry Karpeev     ierr = MatGetSubMatricesParallel(pc->pmat,osm->n,osm->is, osm->is,scall,&osm->pmat);CHKERRQ(ierr);
45681a5c4aaSDmitry Karpeev   }
45781a5c4aaSDmitry Karpeev   else {
45881a5c4aaSDmitry Karpeev     ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->is, osm->is,scall,&osm->pmat);CHKERRQ(ierr);
45981a5c4aaSDmitry Karpeev   }
460b1a0a8a3SJed Brown   if (scall == MAT_INITIAL_MATRIX) {
461b1a0a8a3SJed Brown     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
462f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
463b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,osm->pmat[i]);CHKERRQ(ierr);
464b1a0a8a3SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
465b1a0a8a3SJed Brown     }
466b1a0a8a3SJed Brown   }
467b1a0a8a3SJed Brown 
468b1a0a8a3SJed Brown   /* Return control to the user so that the submatrices can be modified (e.g., to apply
469b1a0a8a3SJed Brown      different boundary conditions for the submatrices than for the global problem) */
470c10bc1a1SDmitry Karpeev   ierr = PCModifySubMatrices(pc,osm->n,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
471b1a0a8a3SJed Brown 
472b1a0a8a3SJed Brown   /*
473f746d493SDmitry Karpeev      Loop over submatrices putting them into local ksp
474b1a0a8a3SJed Brown   */
475f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
476b1a0a8a3SJed Brown     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr);
477b1a0a8a3SJed Brown     if (!pc->setupcalled) {
478b1a0a8a3SJed Brown       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
479b1a0a8a3SJed Brown     }
480b1a0a8a3SJed Brown   }
481b1a0a8a3SJed Brown 
482b1a0a8a3SJed Brown   PetscFunctionReturn(0);
483b1a0a8a3SJed Brown }
484b1a0a8a3SJed Brown 
485b1a0a8a3SJed Brown #undef __FUNCT__
486f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM"
487f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
488b1a0a8a3SJed Brown {
489f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
490b1a0a8a3SJed Brown   PetscErrorCode ierr;
491b1a0a8a3SJed Brown   PetscInt       i;
492b1a0a8a3SJed Brown 
493b1a0a8a3SJed Brown   PetscFunctionBegin;
494f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
495b1a0a8a3SJed Brown     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
496b1a0a8a3SJed Brown   }
497b1a0a8a3SJed Brown   PetscFunctionReturn(0);
498b1a0a8a3SJed Brown }
499b1a0a8a3SJed Brown 
500b1a0a8a3SJed Brown #undef __FUNCT__
501f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM"
502f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y)
503b1a0a8a3SJed Brown {
504f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
505b1a0a8a3SJed Brown   PetscErrorCode ierr;
506f746d493SDmitry Karpeev   PetscInt       i;
507b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
508b1a0a8a3SJed Brown 
509b1a0a8a3SJed Brown   PetscFunctionBegin;
510b1a0a8a3SJed Brown   /*
511b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
512b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
513b1a0a8a3SJed Brown   */
514f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
515b1a0a8a3SJed Brown     forward = SCATTER_FORWARD_LOCAL;
516b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
517f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
518b1a0a8a3SJed Brown   }
519f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
520b1a0a8a3SJed Brown     reverse = SCATTER_REVERSE_LOCAL;
521b1a0a8a3SJed Brown   }
522b1a0a8a3SJed Brown 
523f746d493SDmitry Karpeev   ierr = VecScatterBegin(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
524b1a0a8a3SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
525f746d493SDmitry Karpeev   ierr = VecScatterEnd(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
52686b96d36SDmitry Karpeev   /* do the subdomain solves */
52786b96d36SDmitry Karpeev   for (i=0; i<osm->n; ++i) {
528b1a0a8a3SJed Brown     ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
529b1a0a8a3SJed Brown   }
530f746d493SDmitry Karpeev   ierr = VecScatterBegin(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
531f746d493SDmitry Karpeev   ierr = VecScatterEnd(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
532b1a0a8a3SJed Brown   PetscFunctionReturn(0);
533b1a0a8a3SJed Brown }
534b1a0a8a3SJed Brown 
535b1a0a8a3SJed Brown #undef __FUNCT__
536f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM"
537f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y)
538b1a0a8a3SJed Brown {
539f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
540b1a0a8a3SJed Brown   PetscErrorCode ierr;
541f746d493SDmitry Karpeev   PetscInt       i;
542b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
543b1a0a8a3SJed Brown 
544b1a0a8a3SJed Brown   PetscFunctionBegin;
545b1a0a8a3SJed Brown   /*
546b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
547b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
548b1a0a8a3SJed Brown 
549f746d493SDmitry Karpeev      Note: these are reversed from the PCApply_GASM() because we are applying the
550b1a0a8a3SJed Brown      transpose of the three terms
551b1a0a8a3SJed Brown   */
552f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
553b1a0a8a3SJed Brown     forward = SCATTER_FORWARD_LOCAL;
554b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
555f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
556b1a0a8a3SJed Brown   }
557f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
558b1a0a8a3SJed Brown     reverse = SCATTER_REVERSE_LOCAL;
559b1a0a8a3SJed Brown   }
560b1a0a8a3SJed Brown 
561ab3e923fSDmitry Karpeev   ierr = VecScatterBegin(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
562b1a0a8a3SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
563ab3e923fSDmitry Karpeev   ierr = VecScatterEnd(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
564b1a0a8a3SJed Brown   /* do the local solves */
565ab3e923fSDmitry 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. */
566b1a0a8a3SJed Brown     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
567b1a0a8a3SJed Brown   }
568ab3e923fSDmitry Karpeev   ierr = VecScatterBegin(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
569ab3e923fSDmitry Karpeev   ierr = VecScatterEnd(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
570b1a0a8a3SJed Brown   PetscFunctionReturn(0);
571b1a0a8a3SJed Brown }
572b1a0a8a3SJed Brown 
573b1a0a8a3SJed Brown #undef __FUNCT__
574a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM"
575a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc)
576b1a0a8a3SJed Brown {
577f746d493SDmitry Karpeev   PC_GASM        *osm = (PC_GASM*)pc->data;
578b1a0a8a3SJed Brown   PetscErrorCode ierr;
579b1a0a8a3SJed Brown   PetscInt       i;
580b1a0a8a3SJed Brown 
581b1a0a8a3SJed Brown   PetscFunctionBegin;
582b1a0a8a3SJed Brown   if (osm->ksp) {
583f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
584a06653b4SBarry Smith       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
585b1a0a8a3SJed Brown     }
586b1a0a8a3SJed Brown   }
587b1a0a8a3SJed Brown   if (osm->pmat) {
588f746d493SDmitry Karpeev     if (osm->n > 0) {
589ab3e923fSDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
590b1a0a8a3SJed Brown     }
591b1a0a8a3SJed Brown   }
592d34fcf5fSBarry Smith   if (osm->x) {
593f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
594fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
595fcfd50ebSBarry Smith       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
596b1a0a8a3SJed Brown     }
597d34fcf5fSBarry Smith   }
5986bf464f9SBarry Smith   ierr = VecDestroy(&osm->gx);CHKERRQ(ierr);
5996bf464f9SBarry Smith   ierr = VecDestroy(&osm->gy);CHKERRQ(ierr);
600ab3e923fSDmitry Karpeev 
6016bf464f9SBarry Smith   ierr = VecScatterDestroy(&osm->grestriction);CHKERRQ(ierr);
6026bf464f9SBarry Smith   ierr = VecScatterDestroy(&osm->gprolongation);CHKERRQ(ierr);
603b6b0974bSDmitry Karpeev   if (osm->is) {
604b6b0974bSDmitry Karpeev     ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
605b6b0974bSDmitry Karpeev     osm->is = 0;
606b6b0974bSDmitry Karpeev     osm->is_local = 0;
607b6b0974bSDmitry Karpeev   }
608fcfd50ebSBarry Smith   ierr = ISDestroy(&osm->gis);CHKERRQ(ierr);
609fcfd50ebSBarry Smith   ierr = ISDestroy(&osm->gis_local);CHKERRQ(ierr);
610a06653b4SBarry Smith   PetscFunctionReturn(0);
611a06653b4SBarry Smith }
612a06653b4SBarry Smith 
613a06653b4SBarry Smith #undef __FUNCT__
614a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM"
615a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc)
616a06653b4SBarry Smith {
617a06653b4SBarry Smith   PC_GASM         *osm = (PC_GASM*)pc->data;
618a06653b4SBarry Smith   PetscErrorCode ierr;
619a06653b4SBarry Smith   PetscInt       i;
620a06653b4SBarry Smith 
621a06653b4SBarry Smith   PetscFunctionBegin;
622a06653b4SBarry Smith   ierr = PCReset_GASM(pc);CHKERRQ(ierr);
623a06653b4SBarry Smith   if (osm->ksp) {
624a06653b4SBarry Smith     for (i=0; i<osm->n; i++) {
6256bf464f9SBarry Smith       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
626a06653b4SBarry Smith     }
627a06653b4SBarry Smith     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
628a06653b4SBarry Smith   }
629a06653b4SBarry Smith   ierr = PetscFree(osm->x);CHKERRQ(ierr);
630a06653b4SBarry Smith   ierr = PetscFree(osm->y);CHKERRQ(ierr);
631c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
632b1a0a8a3SJed Brown   PetscFunctionReturn(0);
633b1a0a8a3SJed Brown }
634b1a0a8a3SJed Brown 
635b1a0a8a3SJed Brown #undef __FUNCT__
636f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM"
637ab3e923fSDmitry Karpeev static PetscErrorCode PCSetFromOptions_GASM(PC pc) {
638f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
639b1a0a8a3SJed Brown   PetscErrorCode ierr;
640b1a0a8a3SJed Brown   PetscInt       blocks,ovl;
641b1a0a8a3SJed Brown   PetscBool      symset,flg;
642f746d493SDmitry Karpeev   PCGASMType      gasmtype;
643b1a0a8a3SJed Brown 
644b1a0a8a3SJed Brown   PetscFunctionBegin;
645b1a0a8a3SJed Brown   /* set the type to symmetric if matrix is symmetric */
646b1a0a8a3SJed Brown   if (!osm->type_set && pc->pmat) {
647b1a0a8a3SJed Brown     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
648f746d493SDmitry Karpeev     if (symset && flg) { osm->type = PC_GASM_BASIC; }
649b1a0a8a3SJed Brown   }
65044fe18e5SDmitry Karpeev   ierr = PetscOptionsHead("Generalized additive Schwarz options");CHKERRQ(ierr);
651*06b43e7eSDmitry Karpeev     ierr = PetscOptionsInt("-pc_gasm_total_block_count","Total number of local subdomains across communicator","PCGASMSetTotalBlockCount",osm->n,&blocks,&flg);CHKERRQ(ierr);
652*06b43e7eSDmitry Karpeev     if (flg) {ierr = PCGASMSetTotalBlockCount(pc,blocks);CHKERRQ(ierr); }
653*06b43e7eSDmitry Karpeev     ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
654f746d493SDmitry Karpeev     if (flg) {ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); }
655b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
656f746d493SDmitry Karpeev     ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr);
657f746d493SDmitry Karpeev     if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr); }
658b1a0a8a3SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
659b1a0a8a3SJed Brown   PetscFunctionReturn(0);
660b1a0a8a3SJed Brown }
661b1a0a8a3SJed Brown 
662b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/
663b1a0a8a3SJed Brown 
664b1a0a8a3SJed Brown EXTERN_C_BEGIN
665b1a0a8a3SJed Brown #undef __FUNCT__
666*06b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM"
667*06b43e7eSDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS is[],IS is_local[])
668b1a0a8a3SJed Brown {
669f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
670b1a0a8a3SJed Brown   PetscErrorCode ierr;
671b1a0a8a3SJed Brown   PetscInt       i;
672b1a0a8a3SJed Brown 
673b1a0a8a3SJed Brown   PetscFunctionBegin;
67469ca1f37SDmitry Karpeev   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, n = %D",n);
675*06b43e7eSDmitry Karpeev   if (pc->setupcalled && (n != osm->n || is)) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp().");
676b1a0a8a3SJed Brown 
677b1a0a8a3SJed Brown   if (!pc->setupcalled) {
67893c1ec32SDmitry Karpeev     osm->n            = n;
67993c1ec32SDmitry Karpeev     osm->is           = 0;
68093c1ec32SDmitry Karpeev     osm->is_local     = 0;
681b1a0a8a3SJed Brown     if (is) {
682b1a0a8a3SJed Brown       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
683b1a0a8a3SJed Brown     }
684b1a0a8a3SJed Brown     if (is_local) {
685b1a0a8a3SJed Brown       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
686b1a0a8a3SJed Brown     }
687b1a0a8a3SJed Brown     if (osm->is) {
688ab3e923fSDmitry Karpeev       ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
689b1a0a8a3SJed Brown     }
690b1a0a8a3SJed Brown     if (is) {
691b1a0a8a3SJed Brown       ierr = PetscMalloc(n*sizeof(IS),&osm->is);CHKERRQ(ierr);
692b1a0a8a3SJed Brown       for (i=0; i<n; i++) { osm->is[i] = is[i]; }
693f746d493SDmitry Karpeev       /* Flag indicating that the user has set overlapping subdomains so PCGASM should not increase their size. */
694b1a0a8a3SJed Brown       osm->overlap = -1;
695b1a0a8a3SJed Brown     }
696b1a0a8a3SJed Brown     if (is_local) {
697b1a0a8a3SJed Brown       ierr = PetscMalloc(n*sizeof(IS),&osm->is_local);CHKERRQ(ierr);
698b1a0a8a3SJed Brown       for (i=0; i<n; i++) { osm->is_local[i] = is_local[i]; }
699b1a0a8a3SJed Brown     }
700b1a0a8a3SJed Brown   }
701b1a0a8a3SJed Brown   PetscFunctionReturn(0);
702b1a0a8a3SJed Brown }
703b1a0a8a3SJed Brown EXTERN_C_END
704b1a0a8a3SJed Brown 
705b1a0a8a3SJed Brown EXTERN_C_BEGIN
706b1a0a8a3SJed Brown #undef __FUNCT__
707*06b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalBlockCount_GASM"
708*06b43e7eSDmitry Karpeev PetscErrorCode  PCGASMSetTotalBlockCount_GASM(PC pc,PetscInt N) {
709f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
710b1a0a8a3SJed Brown   PetscErrorCode ierr;
711b1a0a8a3SJed Brown   PetscMPIInt    rank,size;
712b1a0a8a3SJed Brown   PetscInt       n;
713b1a0a8a3SJed Brown 
714b1a0a8a3SJed Brown   PetscFunctionBegin;
715*06b43e7eSDmitry Karpeev   if (N < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of total number of subdomains must be > 0, N = %D",N);
716b1a0a8a3SJed Brown 
717b1a0a8a3SJed Brown   /*
718b1a0a8a3SJed Brown      Split the subdomains equally among all processors
719b1a0a8a3SJed Brown   */
720b1a0a8a3SJed Brown   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
721b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
722b1a0a8a3SJed Brown   n = N/size + ((N % size) > rank);
723b1a0a8a3SJed Brown   if (!n) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Process %d must have at least one block: total processors %d total blocks %D",(int)rank,(int)size,N);
724f746d493SDmitry Karpeev   if (pc->setupcalled && n != osm->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before PCSetUp().");
725b1a0a8a3SJed Brown   if (!pc->setupcalled) {
726b1a0a8a3SJed Brown     if (osm->is) {
727ab3e923fSDmitry Karpeev       ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
728b1a0a8a3SJed Brown     }
72993c1ec32SDmitry Karpeev     osm->N            = N;
730f746d493SDmitry Karpeev     osm->n            = n;
731*06b43e7eSDmitry Karpeev     osm->nmax         = N/size + ((N%size)?1:0);
732b1a0a8a3SJed Brown     osm->is           = 0;
733b1a0a8a3SJed Brown     osm->is_local     = 0;
734b1a0a8a3SJed Brown   }
735b1a0a8a3SJed Brown   PetscFunctionReturn(0);
736*06b43e7eSDmitry Karpeev }
737b1a0a8a3SJed Brown EXTERN_C_END
738b1a0a8a3SJed Brown 
739b1a0a8a3SJed Brown EXTERN_C_BEGIN
740b1a0a8a3SJed Brown #undef __FUNCT__
741f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM"
7427087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
743b1a0a8a3SJed Brown {
744f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
745b1a0a8a3SJed Brown 
746b1a0a8a3SJed Brown   PetscFunctionBegin;
747b1a0a8a3SJed Brown   if (ovl < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
748f746d493SDmitry Karpeev   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
749b1a0a8a3SJed Brown   if (!pc->setupcalled) {
750b1a0a8a3SJed Brown     osm->overlap = ovl;
751b1a0a8a3SJed Brown   }
752b1a0a8a3SJed Brown   PetscFunctionReturn(0);
753b1a0a8a3SJed Brown }
754b1a0a8a3SJed Brown EXTERN_C_END
755b1a0a8a3SJed Brown 
756b1a0a8a3SJed Brown EXTERN_C_BEGIN
757b1a0a8a3SJed Brown #undef __FUNCT__
758f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM"
7597087cfbeSBarry Smith PetscErrorCode  PCGASMSetType_GASM(PC pc,PCGASMType type)
760b1a0a8a3SJed Brown {
761f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
762b1a0a8a3SJed Brown 
763b1a0a8a3SJed Brown   PetscFunctionBegin;
764b1a0a8a3SJed Brown   osm->type     = type;
765b1a0a8a3SJed Brown   osm->type_set = PETSC_TRUE;
766b1a0a8a3SJed Brown   PetscFunctionReturn(0);
767b1a0a8a3SJed Brown }
768b1a0a8a3SJed Brown EXTERN_C_END
769b1a0a8a3SJed Brown 
770b1a0a8a3SJed Brown EXTERN_C_BEGIN
771b1a0a8a3SJed Brown #undef __FUNCT__
772f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM"
7737087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices_GASM(PC pc,PetscBool  doSort)
774b1a0a8a3SJed Brown {
775f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
776b1a0a8a3SJed Brown 
777b1a0a8a3SJed Brown   PetscFunctionBegin;
778b1a0a8a3SJed Brown   osm->sort_indices = doSort;
779b1a0a8a3SJed Brown   PetscFunctionReturn(0);
780b1a0a8a3SJed Brown }
781b1a0a8a3SJed Brown EXTERN_C_END
782b1a0a8a3SJed Brown 
783b1a0a8a3SJed Brown EXTERN_C_BEGIN
784b1a0a8a3SJed Brown #undef __FUNCT__
785f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM"
78644fe18e5SDmitry Karpeev /*
78744fe18e5SDmitry Karpeev    FIX: This routine might need to be modified once multiple ranks per subdomain are allowed.
78844fe18e5SDmitry Karpeev         In particular, it would upset the global subdomain number calculation.
78944fe18e5SDmitry Karpeev */
7907087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
791b1a0a8a3SJed Brown {
792f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
793b1a0a8a3SJed Brown   PetscErrorCode ierr;
794b1a0a8a3SJed Brown 
795b1a0a8a3SJed Brown   PetscFunctionBegin;
796ab3e923fSDmitry 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");
797b1a0a8a3SJed Brown 
798ab3e923fSDmitry Karpeev   if (n) {
799ab3e923fSDmitry Karpeev     *n = osm->n;
800b1a0a8a3SJed Brown   }
801ab3e923fSDmitry Karpeev   if (first) {
802ab3e923fSDmitry Karpeev     ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr);
803ab3e923fSDmitry Karpeev     *first -= osm->n;
804b1a0a8a3SJed Brown   }
805b1a0a8a3SJed Brown   if (ksp) {
806b1a0a8a3SJed Brown     /* Assume that local solves are now different; not necessarily
807*06b43e7eSDmitry Karpeev        true, though!  This flag is used only for PCView_GASM() */
808b1a0a8a3SJed Brown     *ksp                   = osm->ksp;
809b1a0a8a3SJed Brown     osm->same_local_solves = PETSC_FALSE;
810b1a0a8a3SJed Brown   }
811b1a0a8a3SJed Brown   PetscFunctionReturn(0);
812ab3e923fSDmitry Karpeev }/* PCGASMGetSubKSP_GASM() */
813b1a0a8a3SJed Brown EXTERN_C_END
814b1a0a8a3SJed Brown 
815b1a0a8a3SJed Brown 
816b1a0a8a3SJed Brown #undef __FUNCT__
817*06b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains"
818b1a0a8a3SJed Brown /*@C
819*06b43e7eSDmitry Karpeev     PCGASMSetSubdomains - Sets the subdomains for this processor
820*06b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
821b1a0a8a3SJed Brown 
822b1a0a8a3SJed Brown     Collective on PC
823b1a0a8a3SJed Brown 
824b1a0a8a3SJed Brown     Input Parameters:
825b1a0a8a3SJed Brown +   pc - the preconditioner context
826*06b43e7eSDmitry Karpeev .   n  - the number of subdomains for this processor
827*06b43e7eSDmitry Karpeev .   is - the index sets that define this processor's part of the subdomains (on which the correction is obtained)
828b1a0a8a3SJed Brown          (or PETSC_NULL for PETSc to determine subdomains)
829*06b43e7eSDmitry Karpeev -   is_local - the index sets that define the local part of this processor's part of the subdomains (on which the correction is applied)
830*06b43e7eSDmitry Karpeev          (or PETSC_NULL to use the same as is)
831b1a0a8a3SJed Brown 
832b1a0a8a3SJed Brown     Notes:
833*06b43e7eSDmitry Karpeev     The IS indices use the parallel, global numbering of the vector entries.
834*06b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.
835*06b43e7eSDmitry Karpeev     Use PCGASMSetTotalBlockCount() to set the total number of subdomains across all processors
836*06b43e7eSDmitry Karpeev     and to allow PETSc to assign them to processor and to compute the subdomains from that.
837b1a0a8a3SJed Brown 
838b1a0a8a3SJed Brown     Level: advanced
839b1a0a8a3SJed Brown 
840*06b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz
841b1a0a8a3SJed Brown 
842*06b43e7eSDmitry Karpeev .seealso: PCGASMSetTotalBlockCount(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
843*06b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
844b1a0a8a3SJed Brown @*/
845*06b43e7eSDmitry Karpeev PetscErrorCode  PCGASMSetSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
846b1a0a8a3SJed Brown {
847b1a0a8a3SJed Brown   PetscErrorCode ierr;
848b1a0a8a3SJed Brown 
849b1a0a8a3SJed Brown   PetscFunctionBegin;
850b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
851*06b43e7eSDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr);
852b1a0a8a3SJed Brown   PetscFunctionReturn(0);
853b1a0a8a3SJed Brown }
854b1a0a8a3SJed Brown 
855b1a0a8a3SJed Brown #undef __FUNCT__
856*06b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalBlockCount"
857b1a0a8a3SJed Brown /*@C
858*06b43e7eSDmitry Karpeev     PCGASMSetTotalBlockCount - Sets the total number of subdomains
859*06b43e7eSDmitry Karpeev     cumulative across all processors across the communicator of a generalized
860b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
861*06b43e7eSDmitry Karpeev     PC communicator must call this routine, with the same number.
862b1a0a8a3SJed Brown 
863b1a0a8a3SJed Brown     Collective on PC
864b1a0a8a3SJed Brown 
865b1a0a8a3SJed Brown     Input Parameters:
866b1a0a8a3SJed Brown +   pc - the preconditioner context
867*06b43e7eSDmitry Karpeev -   n  - the total number of subdomains cumulative across all processors
868b1a0a8a3SJed Brown 
869b1a0a8a3SJed Brown     Options Database Key:
870b1a0a8a3SJed Brown     To set the total number of subdomain blocks rather than specify the
871b1a0a8a3SJed Brown     index sets, use the option
872*06b43e7eSDmitry Karpeev .    -pc_gasm_total_block_count <n> - Sets the total number of local subdomains (known as blocks) to be distributed among the processors
873b1a0a8a3SJed Brown 
874*06b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.
875b1a0a8a3SJed Brown 
876*06b43e7eSDmitry Karpeev     Use PCGASMSetSubdomains() to set local subdomains.
877b1a0a8a3SJed Brown 
878b1a0a8a3SJed Brown     Level: advanced
879b1a0a8a3SJed Brown 
880f746d493SDmitry Karpeev .keywords: PC, GASM, set, total, global, subdomains, additive Schwarz
881b1a0a8a3SJed Brown 
882*06b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
883f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
884b1a0a8a3SJed Brown @*/
885*06b43e7eSDmitry Karpeev PetscErrorCode  PCGASMSetTotalBlockCount(PC pc,PetscInt N)
886b1a0a8a3SJed Brown {
887b1a0a8a3SJed Brown   PetscErrorCode ierr;
888b1a0a8a3SJed Brown 
889b1a0a8a3SJed Brown   PetscFunctionBegin;
890b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
891*06b43e7eSDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetTotalBlockCount_C",(PC,PetscInt),(pc,N));CHKERRQ(ierr);
892b1a0a8a3SJed Brown   PetscFunctionReturn(0);
893b1a0a8a3SJed Brown }
894b1a0a8a3SJed Brown 
895b1a0a8a3SJed Brown #undef __FUNCT__
896f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap"
897b1a0a8a3SJed Brown /*@
898f746d493SDmitry Karpeev     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
899b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
900b1a0a8a3SJed Brown     PC communicator must call this routine.
901b1a0a8a3SJed Brown 
902b1a0a8a3SJed Brown     Logically Collective on PC
903b1a0a8a3SJed Brown 
904b1a0a8a3SJed Brown     Input Parameters:
905b1a0a8a3SJed Brown +   pc  - the preconditioner context
906b1a0a8a3SJed Brown -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
907b1a0a8a3SJed Brown 
908b1a0a8a3SJed Brown     Options Database Key:
909*06b43e7eSDmitry Karpeev .   -pc_gasm_overlap <overlap> - Sets overlap
910b1a0a8a3SJed Brown 
911b1a0a8a3SJed Brown     Notes:
912*06b43e7eSDmitry Karpeev     By default the GASM preconditioner uses 1 subdomain per processor.  To use
913*06b43e7eSDmitry Karpeev     multiple subdomain per perocessor, see PCGASMSetTotalBlockCount() or
914*06b43e7eSDmitry Karpeev     PCGASMSetSubdomains() (and the option -pc_gasm_subdomain_count <n>).
915b1a0a8a3SJed Brown 
916b1a0a8a3SJed Brown     The overlap defaults to 1, so if one desires that no additional
917b1a0a8a3SJed Brown     overlap be computed beyond what may have been set with a call to
918*06b43e7eSDmitry Karpeev     PCGASMSetTotalBlockCount() or PCGASMSetSubdomains(), then ovl
919b1a0a8a3SJed Brown     must be set to be 0.  In particular, if one does not explicitly set
920*06b43e7eSDmitry Karpeev     the subdomains in application code, then all overlap would be computed
921f746d493SDmitry Karpeev     internally by PETSc, and using an overlap of 0 would result in an GASM
922b1a0a8a3SJed Brown     variant that is equivalent to the block Jacobi preconditioner.
923b1a0a8a3SJed Brown 
924b1a0a8a3SJed Brown     Note that one can define initial index sets with any overlap via
925*06b43e7eSDmitry Karpeev     PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows
926*06b43e7eSDmitry Karpeev     PETSc to extend that overlap further, if desired.
927b1a0a8a3SJed Brown 
928b1a0a8a3SJed Brown     Level: intermediate
929b1a0a8a3SJed Brown 
930f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap
931b1a0a8a3SJed Brown 
932*06b43e7eSDmitry Karpeev .seealso: PCGASMSetTotalBlockCount(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
933*06b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetSubdomains()
934b1a0a8a3SJed Brown @*/
9357087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap(PC pc,PetscInt ovl)
936b1a0a8a3SJed Brown {
937b1a0a8a3SJed Brown   PetscErrorCode ierr;
938b1a0a8a3SJed Brown 
939b1a0a8a3SJed Brown   PetscFunctionBegin;
940b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
941b1a0a8a3SJed Brown   PetscValidLogicalCollectiveInt(pc,ovl,2);
942f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
943b1a0a8a3SJed Brown   PetscFunctionReturn(0);
944b1a0a8a3SJed Brown }
945b1a0a8a3SJed Brown 
946b1a0a8a3SJed Brown #undef __FUNCT__
947f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType"
948b1a0a8a3SJed Brown /*@
949f746d493SDmitry Karpeev     PCGASMSetType - Sets the type of restriction and interpolation used
950b1a0a8a3SJed Brown     for local problems in the additive Schwarz method.
951b1a0a8a3SJed Brown 
952b1a0a8a3SJed Brown     Logically Collective on PC
953b1a0a8a3SJed Brown 
954b1a0a8a3SJed Brown     Input Parameters:
955b1a0a8a3SJed Brown +   pc  - the preconditioner context
956f746d493SDmitry Karpeev -   type - variant of GASM, one of
957b1a0a8a3SJed Brown .vb
958f746d493SDmitry Karpeev       PC_GASM_BASIC       - full interpolation and restriction
959f746d493SDmitry Karpeev       PC_GASM_RESTRICT    - full restriction, local processor interpolation
960f746d493SDmitry Karpeev       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
961f746d493SDmitry Karpeev       PC_GASM_NONE        - local processor restriction and interpolation
962b1a0a8a3SJed Brown .ve
963b1a0a8a3SJed Brown 
964b1a0a8a3SJed Brown     Options Database Key:
965f746d493SDmitry Karpeev .   -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
966b1a0a8a3SJed Brown 
967b1a0a8a3SJed Brown     Level: intermediate
968b1a0a8a3SJed Brown 
969f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
970b1a0a8a3SJed Brown 
971*06b43e7eSDmitry Karpeev .seealso: PCGASMSetTotalBlockCount(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
972f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
973b1a0a8a3SJed Brown @*/
9747087cfbeSBarry Smith PetscErrorCode  PCGASMSetType(PC pc,PCGASMType type)
975b1a0a8a3SJed Brown {
976b1a0a8a3SJed Brown   PetscErrorCode ierr;
977b1a0a8a3SJed Brown 
978b1a0a8a3SJed Brown   PetscFunctionBegin;
979b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
980b1a0a8a3SJed Brown   PetscValidLogicalCollectiveEnum(pc,type,2);
981f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr);
982b1a0a8a3SJed Brown   PetscFunctionReturn(0);
983b1a0a8a3SJed Brown }
984b1a0a8a3SJed Brown 
985b1a0a8a3SJed Brown #undef __FUNCT__
986f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices"
987b1a0a8a3SJed Brown /*@
988f746d493SDmitry Karpeev     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
989b1a0a8a3SJed Brown 
990b1a0a8a3SJed Brown     Logically Collective on PC
991b1a0a8a3SJed Brown 
992b1a0a8a3SJed Brown     Input Parameters:
993b1a0a8a3SJed Brown +   pc  - the preconditioner context
994b1a0a8a3SJed Brown -   doSort - sort the subdomain indices
995b1a0a8a3SJed Brown 
996b1a0a8a3SJed Brown     Level: intermediate
997b1a0a8a3SJed Brown 
998f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
999b1a0a8a3SJed Brown 
1000*06b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetTotalBlockCount(), PCGASMGetSubKSP(),
1001f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
1002b1a0a8a3SJed Brown @*/
10037087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices(PC pc,PetscBool  doSort)
1004b1a0a8a3SJed Brown {
1005b1a0a8a3SJed Brown   PetscErrorCode ierr;
1006b1a0a8a3SJed Brown 
1007b1a0a8a3SJed Brown   PetscFunctionBegin;
1008b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1009b1a0a8a3SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
1010f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
1011b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1012b1a0a8a3SJed Brown }
1013b1a0a8a3SJed Brown 
1014b1a0a8a3SJed Brown #undef __FUNCT__
1015f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP"
1016b1a0a8a3SJed Brown /*@C
1017f746d493SDmitry Karpeev    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
1018b1a0a8a3SJed Brown    this processor.
1019b1a0a8a3SJed Brown 
1020b1a0a8a3SJed Brown    Collective on PC iff first_local is requested
1021b1a0a8a3SJed Brown 
1022b1a0a8a3SJed Brown    Input Parameter:
1023b1a0a8a3SJed Brown .  pc - the preconditioner context
1024b1a0a8a3SJed Brown 
1025b1a0a8a3SJed Brown    Output Parameters:
1026b1a0a8a3SJed Brown +  n_local - the number of blocks on this processor or PETSC_NULL
1027b1a0a8a3SJed Brown .  first_local - the global number of the first block on this processor or PETSC_NULL,
1028b1a0a8a3SJed Brown                  all processors must request or all must pass PETSC_NULL
1029b1a0a8a3SJed Brown -  ksp - the array of KSP contexts
1030b1a0a8a3SJed Brown 
1031b1a0a8a3SJed Brown    Note:
1032f746d493SDmitry Karpeev    After PCGASMGetSubKSP() the array of KSPes is not to be freed
1033b1a0a8a3SJed Brown 
1034b1a0a8a3SJed Brown    Currently for some matrix implementations only 1 block per processor
1035b1a0a8a3SJed Brown    is supported.
1036b1a0a8a3SJed Brown 
1037f746d493SDmitry Karpeev    You must call KSPSetUp() before calling PCGASMGetSubKSP().
1038b1a0a8a3SJed Brown 
1039b1a0a8a3SJed Brown    Level: advanced
1040b1a0a8a3SJed Brown 
1041f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context
1042b1a0a8a3SJed Brown 
1043*06b43e7eSDmitry Karpeev .seealso: PCGASMSetTotalBlockCount(), PCGASMSetSubdomains(), PCGASMSetOverlap(),
1044f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(),
1045b1a0a8a3SJed Brown @*/
10467087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1047b1a0a8a3SJed Brown {
1048b1a0a8a3SJed Brown   PetscErrorCode ierr;
1049b1a0a8a3SJed Brown 
1050b1a0a8a3SJed Brown   PetscFunctionBegin;
1051b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1052f746d493SDmitry Karpeev   ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
1053b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1054b1a0a8a3SJed Brown }
1055b1a0a8a3SJed Brown 
1056b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/
1057b1a0a8a3SJed Brown /*MC
1058f746d493SDmitry Karpeev    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1059b1a0a8a3SJed Brown            its own KSP object.
1060b1a0a8a3SJed Brown 
1061b1a0a8a3SJed Brown    Options Database Keys:
1062*06b43e7eSDmitry Karpeev +  -pc_gasm_total_block_count <n> - Sets total number of local subdomains (known as blocks) to be distributed among processors
1063*06b43e7eSDmitry Karpeev .  -pc_gasm_view_subdomains       - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view
1064*06b43e7eSDmitry Karpeev .  -pc_gasm_print_subdomains      - activates the printing of subdomain indices in PCSetUp()
1065*06b43e7eSDmitry Karpeev .  -pc_gasm_overlap <ovl>         - Sets overlap by which to (automatically) extend local subdomains
1066f746d493SDmitry Karpeev -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
1067b1a0a8a3SJed Brown 
1068b1a0a8a3SJed Brown      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
1069f746d493SDmitry Karpeev       will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
1070f746d493SDmitry Karpeev       -pc_gasm_type basic to use the standard GASM.
1071b1a0a8a3SJed Brown 
1072b1a0a8a3SJed Brown    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
1073b1a0a8a3SJed Brown      than one processor. Defaults to one block per processor.
1074b1a0a8a3SJed Brown 
1075b1a0a8a3SJed Brown      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1076b1a0a8a3SJed Brown         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
1077b1a0a8a3SJed Brown 
1078f746d493SDmitry Karpeev      To set the options on the solvers separate for each block call PCGASMGetSubKSP()
1079b1a0a8a3SJed Brown          and set the options directly on the resulting KSP object (you can access its PC
1080b1a0a8a3SJed Brown          with KSPGetPC())
1081b1a0a8a3SJed Brown 
1082b1a0a8a3SJed Brown 
1083b1a0a8a3SJed Brown    Level: beginner
1084b1a0a8a3SJed Brown 
1085b1a0a8a3SJed Brown    Concepts: additive Schwarz method
1086b1a0a8a3SJed Brown 
1087b1a0a8a3SJed Brown     References:
1088b1a0a8a3SJed Brown     An additive variant of the Schwarz alternating method for the case of many subregions
1089b1a0a8a3SJed Brown     M Dryja, OB Widlund - Courant Institute, New York University Technical report
1090b1a0a8a3SJed Brown 
1091b1a0a8a3SJed Brown     Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1092b1a0a8a3SJed Brown     Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X.
1093b1a0a8a3SJed Brown 
1094b1a0a8a3SJed Brown .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1095*06b43e7eSDmitry Karpeev            PCBJACOBI, PCGASMSetUseTrueLocal(), PCGASMGetSubKSP(), PCGASMSetSubdomains(),
1096*06b43e7eSDmitry Karpeev            PCGASMSetTotalBlockCount(), PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType()
1097b1a0a8a3SJed Brown 
1098b1a0a8a3SJed Brown M*/
1099b1a0a8a3SJed Brown 
1100b1a0a8a3SJed Brown EXTERN_C_BEGIN
1101b1a0a8a3SJed Brown #undef __FUNCT__
1102f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM"
11037087cfbeSBarry Smith PetscErrorCode  PCCreate_GASM(PC pc)
1104b1a0a8a3SJed Brown {
1105b1a0a8a3SJed Brown   PetscErrorCode ierr;
1106f746d493SDmitry Karpeev   PC_GASM         *osm;
1107b1a0a8a3SJed Brown 
1108b1a0a8a3SJed Brown   PetscFunctionBegin;
1109f746d493SDmitry Karpeev   ierr = PetscNewLog(pc,PC_GASM,&osm);CHKERRQ(ierr);
1110ab3e923fSDmitry Karpeev   osm->N                 = PETSC_DECIDE;
1111*06b43e7eSDmitry Karpeev   osm->n                 = PETSC_DECIDE;
1112ab3e923fSDmitry Karpeev   osm->nmax              = 0;
1113b1a0a8a3SJed Brown   osm->overlap           = 1;
1114b1a0a8a3SJed Brown   osm->ksp               = 0;
1115ab3e923fSDmitry Karpeev   osm->grestriction      = 0;
1116ab3e923fSDmitry Karpeev   osm->gprolongation     = 0;
1117ab3e923fSDmitry Karpeev   osm->gx                = 0;
1118ab3e923fSDmitry Karpeev   osm->gy                = 0;
1119b1a0a8a3SJed Brown   osm->x                 = 0;
1120b1a0a8a3SJed Brown   osm->y                 = 0;
1121b1a0a8a3SJed Brown   osm->is                = 0;
1122b1a0a8a3SJed Brown   osm->is_local          = 0;
1123b1a0a8a3SJed Brown   osm->pmat              = 0;
1124f746d493SDmitry Karpeev   osm->type              = PC_GASM_RESTRICT;
1125b1a0a8a3SJed Brown   osm->same_local_solves = PETSC_TRUE;
1126b1a0a8a3SJed Brown   osm->sort_indices      = PETSC_TRUE;
1127b1a0a8a3SJed Brown 
1128b1a0a8a3SJed Brown   pc->data                   = (void*)osm;
1129f746d493SDmitry Karpeev   pc->ops->apply             = PCApply_GASM;
1130f746d493SDmitry Karpeev   pc->ops->applytranspose    = PCApplyTranspose_GASM;
1131f746d493SDmitry Karpeev   pc->ops->setup             = PCSetUp_GASM;
1132a06653b4SBarry Smith   pc->ops->reset             = PCReset_GASM;
1133f746d493SDmitry Karpeev   pc->ops->destroy           = PCDestroy_GASM;
1134f746d493SDmitry Karpeev   pc->ops->setfromoptions    = PCSetFromOptions_GASM;
1135f746d493SDmitry Karpeev   pc->ops->setuponblocks     = PCSetUpOnBlocks_GASM;
1136f746d493SDmitry Karpeev   pc->ops->view              = PCView_GASM;
1137b1a0a8a3SJed Brown   pc->ops->applyrichardson   = 0;
1138b1a0a8a3SJed Brown 
1139*06b43e7eSDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSubdomains_C","PCGASMSetSubdomains_GASM",
1140*06b43e7eSDmitry Karpeev                     PCGASMSetSubdomains_GASM);CHKERRQ(ierr);
1141*06b43e7eSDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetTotalBlockCount_C","PCGASMSetTotalBlockCount_GASM",
1142*06b43e7eSDmitry Karpeev                     PCGASMSetTotalBlockCount_GASM);CHKERRQ(ierr);
1143f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetOverlap_C","PCGASMSetOverlap_GASM",
1144f746d493SDmitry Karpeev                     PCGASMSetOverlap_GASM);CHKERRQ(ierr);
1145f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetType_C","PCGASMSetType_GASM",
1146f746d493SDmitry Karpeev                     PCGASMSetType_GASM);CHKERRQ(ierr);
1147f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSortIndices_C","PCGASMSetSortIndices_GASM",
1148f746d493SDmitry Karpeev                     PCGASMSetSortIndices_GASM);CHKERRQ(ierr);
1149f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMGetSubKSP_C","PCGASMGetSubKSP_GASM",
1150f746d493SDmitry Karpeev                     PCGASMGetSubKSP_GASM);CHKERRQ(ierr);
1151b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1152b1a0a8a3SJed Brown }
1153b1a0a8a3SJed Brown EXTERN_C_END
1154b1a0a8a3SJed Brown 
1155b1a0a8a3SJed Brown 
1156b1a0a8a3SJed Brown #undef __FUNCT__
1157*06b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains"
1158b1a0a8a3SJed Brown /*@C
1159*06b43e7eSDmitry Karpeev    PCGASMCreateLocalSubdomains - Creates n local index sets for the overlapping
1160*06b43e7eSDmitry Karpeev    Schwarz preconditioner for a any problem based on its matrix.
1161b1a0a8a3SJed Brown 
1162b1a0a8a3SJed Brown    Collective
1163b1a0a8a3SJed Brown 
1164b1a0a8a3SJed Brown    Input Parameters:
1165b1a0a8a3SJed Brown +  A - The global matrix operator
1166*06b43e7eSDmitry Karpeev -  n - the number of local subdomains
1167b1a0a8a3SJed Brown 
1168b1a0a8a3SJed Brown    Output Parameters:
1169*06b43e7eSDmitry Karpeev +  outis       - the array of index sets defining the local subdomains (on which the correction is computed)
1170*06b43e7eSDmitry Karpeev -  outis_local - the array of index sets defining the local portions of the subdomains (on which the correction is applied)
1171b1a0a8a3SJed Brown 
1172b1a0a8a3SJed Brown    Level: advanced
1173b1a0a8a3SJed Brown 
1174*06b43e7eSDmitry Karpeev    Note: this generates n nonoverlapping subdomains on PETSC_COMM_SELF;
1175*06b43e7eSDmitry Karpeev          PCGASM will generate the overlap from these if you use them
1176*06b43e7eSDmitry Karpeev          in PCGASMSetSubdomains() and set a nonzero overlap with
1177*06b43e7eSDmitry Karpeev          PCGASMSetOverlap()
1178b1a0a8a3SJed Brown 
1179b1a0a8a3SJed Brown     In the Fortran version you must provide the array outis[] already allocated of length n.
1180b1a0a8a3SJed Brown 
1181f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1182b1a0a8a3SJed Brown 
1183*06b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains()
1184b1a0a8a3SJed Brown @*/
1185*06b43e7eSDmitry Karpeev PetscErrorCode  PCGASMCreateLocalSubdomains(Mat A, PetscInt n, IS* outis[], IS*outis_local[])
1186b1a0a8a3SJed Brown {
1187b1a0a8a3SJed Brown   MatPartitioning           mpart;
1188b1a0a8a3SJed Brown   const char                *prefix;
118911bd1e4dSLisandro Dalcin   PetscErrorCode            (*f)(Mat,MatReuse,Mat*);
1190b1a0a8a3SJed Brown   PetscMPIInt               size;
1191b1a0a8a3SJed Brown   PetscInt                  i,j,rstart,rend,bs;
119211bd1e4dSLisandro Dalcin   PetscBool                 isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
1193b1a0a8a3SJed Brown   Mat                       Ad = PETSC_NULL, adj;
1194b1a0a8a3SJed Brown   IS                        ispart,isnumb,*is;
1195b1a0a8a3SJed Brown   PetscErrorCode            ierr;
1196b1a0a8a3SJed Brown 
1197b1a0a8a3SJed Brown   PetscFunctionBegin;
1198b1a0a8a3SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1199b1a0a8a3SJed Brown   PetscValidPointer(outis,3);
1200b1a0a8a3SJed Brown   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
1201b1a0a8a3SJed Brown 
1202b1a0a8a3SJed Brown   /* Get prefix, row distribution, and block size */
1203b1a0a8a3SJed Brown   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1204b1a0a8a3SJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1205b1a0a8a3SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1206b1a0a8a3SJed 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);
1207b1a0a8a3SJed Brown 
1208b1a0a8a3SJed Brown   /* Get diagonal block from matrix if possible */
1209b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1210b1a0a8a3SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr);
1211b1a0a8a3SJed Brown   if (f) {
121211bd1e4dSLisandro Dalcin     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
1213b1a0a8a3SJed Brown   } else if (size == 1) {
121411bd1e4dSLisandro Dalcin     Ad = A;
1215b1a0a8a3SJed Brown   }
1216b1a0a8a3SJed Brown   if (Ad) {
1217251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1218251f4c67SDmitry Karpeev     if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
1219b1a0a8a3SJed Brown   }
1220b1a0a8a3SJed Brown   if (Ad && n > 1) {
1221b1a0a8a3SJed Brown     PetscBool  match,done;
1222b1a0a8a3SJed Brown     /* Try to setup a good matrix partitioning if available */
1223b1a0a8a3SJed Brown     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
1224b1a0a8a3SJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
1225b1a0a8a3SJed Brown     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1226251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
1227b1a0a8a3SJed Brown     if (!match) {
1228251f4c67SDmitry Karpeev       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
1229b1a0a8a3SJed Brown     }
1230b1a0a8a3SJed Brown     if (!match) { /* assume a "good" partitioner is available */
1231b1a0a8a3SJed Brown       PetscInt na,*ia,*ja;
1232b1a0a8a3SJed Brown       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1233b1a0a8a3SJed Brown       if (done) {
1234b1a0a8a3SJed Brown         /* Build adjacency matrix by hand. Unfortunately a call to
1235b1a0a8a3SJed Brown            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1236b1a0a8a3SJed Brown            remove the block-aij structure and we cannot expect
1237b1a0a8a3SJed Brown            MatPartitioning to split vertices as we need */
1238b1a0a8a3SJed Brown         PetscInt i,j,*row,len,nnz,cnt,*iia=0,*jja=0;
1239b1a0a8a3SJed Brown         nnz = 0;
1240b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* count number of nonzeros */
1241b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1242b1a0a8a3SJed Brown           row = ja + ia[i];
1243b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1244b1a0a8a3SJed Brown             if (row[j] == i) { /* don't count diagonal */
1245b1a0a8a3SJed Brown               len--; break;
1246b1a0a8a3SJed Brown             }
1247b1a0a8a3SJed Brown           }
1248b1a0a8a3SJed Brown           nnz += len;
1249b1a0a8a3SJed Brown         }
1250b1a0a8a3SJed Brown         ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr);
1251b1a0a8a3SJed Brown         ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr);
1252b1a0a8a3SJed Brown         nnz    = 0;
1253b1a0a8a3SJed Brown         iia[0] = 0;
1254b1a0a8a3SJed Brown         for (i=0; i<na; i++) { /* fill adjacency */
1255b1a0a8a3SJed Brown           cnt = 0;
1256b1a0a8a3SJed Brown           len = ia[i+1] - ia[i];
1257b1a0a8a3SJed Brown           row = ja + ia[i];
1258b1a0a8a3SJed Brown           for (j=0; j<len; j++) {
1259b1a0a8a3SJed Brown             if (row[j] != i) { /* if not diagonal */
1260b1a0a8a3SJed Brown               jja[nnz+cnt++] = row[j];
1261b1a0a8a3SJed Brown             }
1262b1a0a8a3SJed Brown           }
1263b1a0a8a3SJed Brown           nnz += cnt;
1264b1a0a8a3SJed Brown           iia[i+1] = nnz;
1265b1a0a8a3SJed Brown         }
1266b1a0a8a3SJed Brown         /* Partitioning of the adjacency matrix */
1267b1a0a8a3SJed Brown         ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);CHKERRQ(ierr);
1268b1a0a8a3SJed Brown         ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
1269b1a0a8a3SJed Brown         ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
1270b1a0a8a3SJed Brown         ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
1271b1a0a8a3SJed Brown         ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
12726bf464f9SBarry Smith         ierr = MatDestroy(&adj);CHKERRQ(ierr);
1273b1a0a8a3SJed Brown         foundpart = PETSC_TRUE;
1274b1a0a8a3SJed Brown       }
1275b1a0a8a3SJed Brown       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1276b1a0a8a3SJed Brown     }
12776bf464f9SBarry Smith     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1278b1a0a8a3SJed Brown   }
1279b1a0a8a3SJed Brown 
1280b1a0a8a3SJed Brown   ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr);
1281b1a0a8a3SJed Brown   *outis = is;
1282b1a0a8a3SJed Brown 
1283b1a0a8a3SJed Brown   if (!foundpart) {
1284b1a0a8a3SJed Brown 
1285b1a0a8a3SJed Brown     /* Partitioning by contiguous chunks of rows */
1286b1a0a8a3SJed Brown 
1287b1a0a8a3SJed Brown     PetscInt mbs   = (rend-rstart)/bs;
1288b1a0a8a3SJed Brown     PetscInt start = rstart;
1289b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1290b1a0a8a3SJed Brown       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
1291b1a0a8a3SJed Brown       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
1292b1a0a8a3SJed Brown       start += count;
1293b1a0a8a3SJed Brown     }
1294b1a0a8a3SJed Brown 
1295b1a0a8a3SJed Brown   } else {
1296b1a0a8a3SJed Brown 
1297b1a0a8a3SJed Brown     /* Partitioning by adjacency of diagonal block  */
1298b1a0a8a3SJed Brown 
1299b1a0a8a3SJed Brown     const PetscInt *numbering;
1300b1a0a8a3SJed Brown     PetscInt       *count,nidx,*indices,*newidx,start=0;
1301b1a0a8a3SJed Brown     /* Get node count in each partition */
1302b1a0a8a3SJed Brown     ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr);
1303b1a0a8a3SJed Brown     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
1304b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1305b1a0a8a3SJed Brown       for (i=0; i<n; i++) count[i] *= bs;
1306b1a0a8a3SJed Brown     }
1307b1a0a8a3SJed Brown     /* Build indices from node numbering */
1308b1a0a8a3SJed Brown     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1309b1a0a8a3SJed Brown     ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr);
1310b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1311b1a0a8a3SJed Brown     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
1312b1a0a8a3SJed Brown     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
1313b1a0a8a3SJed Brown     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
1314b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1315b1a0a8a3SJed Brown       ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr);
1316b1a0a8a3SJed Brown       for (i=0; i<nidx; i++)
1317b1a0a8a3SJed Brown         for (j=0; j<bs; j++)
1318b1a0a8a3SJed Brown           newidx[i*bs+j] = indices[i]*bs + j;
1319b1a0a8a3SJed Brown       ierr = PetscFree(indices);CHKERRQ(ierr);
1320b1a0a8a3SJed Brown       nidx   *= bs;
1321b1a0a8a3SJed Brown       indices = newidx;
1322b1a0a8a3SJed Brown     }
1323b1a0a8a3SJed Brown     /* Shift to get global indices */
1324b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] += rstart;
1325b1a0a8a3SJed Brown 
1326b1a0a8a3SJed Brown     /* Build the index sets for each block */
1327b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1328b1a0a8a3SJed Brown       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
1329b1a0a8a3SJed Brown       ierr   = ISSort(is[i]);CHKERRQ(ierr);
1330b1a0a8a3SJed Brown       start += count[i];
1331b1a0a8a3SJed Brown     }
1332b1a0a8a3SJed Brown 
1333b1a0a8a3SJed Brown     ierr = PetscFree(count);
1334b1a0a8a3SJed Brown     ierr = PetscFree(indices);
1335fcfd50ebSBarry Smith     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1336fcfd50ebSBarry Smith     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
1337b1a0a8a3SJed Brown   }
1338b1a0a8a3SJed Brown 
1339b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1340b1a0a8a3SJed Brown }
1341b1a0a8a3SJed Brown 
1342b1a0a8a3SJed Brown #undef __FUNCT__
1343f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains"
1344b1a0a8a3SJed Brown /*@C
1345f746d493SDmitry Karpeev    PCGASMDestroySubdomains - Destroys the index sets created with
1346*06b43e7eSDmitry Karpeev    PCGASMCreateLocalSubdomains() or PCGASMCreateSubdomains2D. Should be
1347*06b43e7eSDmitry Karpeev    called after setting subdomains with PCGASMSetSubdomains().
1348b1a0a8a3SJed Brown 
1349b1a0a8a3SJed Brown    Collective
1350b1a0a8a3SJed Brown 
1351b1a0a8a3SJed Brown    Input Parameters:
1352b1a0a8a3SJed Brown +  n - the number of index sets
1353b1a0a8a3SJed Brown .  is - the array of index sets
1354b1a0a8a3SJed Brown -  is_local - the array of local index sets, can be PETSC_NULL
1355b1a0a8a3SJed Brown 
1356b1a0a8a3SJed Brown    Level: advanced
1357b1a0a8a3SJed Brown 
1358f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1359b1a0a8a3SJed Brown 
1360*06b43e7eSDmitry Karpeev .seealso: PCGASMCreateLocalSubdomains(), PCGASMSetSubdomains()
1361b1a0a8a3SJed Brown @*/
13627087cfbeSBarry Smith PetscErrorCode  PCGASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
1363b1a0a8a3SJed Brown {
1364b1a0a8a3SJed Brown   PetscInt       i;
1365b1a0a8a3SJed Brown   PetscErrorCode ierr;
1366b1a0a8a3SJed Brown   PetscFunctionBegin;
1367b1a0a8a3SJed Brown   if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"n must be > 0: n = %D",n);
1368b1a0a8a3SJed Brown   PetscValidPointer(is,2);
1369fcfd50ebSBarry Smith   for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); }
1370b1a0a8a3SJed Brown   ierr = PetscFree(is);CHKERRQ(ierr);
1371b1a0a8a3SJed Brown   if (is_local) {
1372b1a0a8a3SJed Brown     PetscValidPointer(is_local,3);
1373fcfd50ebSBarry Smith     for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); }
1374b1a0a8a3SJed Brown     ierr = PetscFree(is_local);CHKERRQ(ierr);
1375b1a0a8a3SJed Brown   }
1376b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1377b1a0a8a3SJed Brown }
1378b1a0a8a3SJed Brown 
1379af538404SDmitry Karpeev 
1380af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \
1381af538404SDmitry Karpeev {                                                                                                       \
1382af538404SDmitry Karpeev  PetscInt first_row = first/M, last_row = last/M+1;                                                     \
1383af538404SDmitry Karpeev   /*                                                                                                    \
1384af538404SDmitry Karpeev    Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners       \
1385af538404SDmitry Karpeev    of the bounding box of the intersection of the subdomain with the local ownership range (local       \
1386af538404SDmitry Karpeev    subdomain).                                                                                          \
1387af538404SDmitry Karpeev    Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows       \
1388af538404SDmitry Karpeev    of the intersection.                                                                                 \
1389af538404SDmitry Karpeev   */                                                                                                    \
1390af538404SDmitry Karpeev   /* ylow_loc is the grid row containing the first element of the local sumbdomain */                   \
1391eec7e3faSDmitry Karpeev   *ylow_loc = PetscMax(first_row,ylow);                                                                    \
1392af538404SDmitry Karpeev   /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1393eec7e3faSDmitry Karpeev   *xleft_loc = *ylow_loc==first_row?PetscMax(first%M,xleft):xleft;                                                                            \
1394af538404SDmitry Karpeev   /* yhigh_loc is the grid row above the last local subdomain element */                                                                    \
1395eec7e3faSDmitry Karpeev   *yhigh_loc = PetscMin(last_row,yhigh);                                                                                                     \
1396af538404SDmitry Karpeev   /* xright is the offset of the end of the  local subdomain within its grid row (might actually be outside the local subdomain) */         \
1397eec7e3faSDmitry Karpeev   *xright_loc = *yhigh_loc==last_row?PetscMin(xright,last%M):xright;                                                                          \
1398af538404SDmitry Karpeev   /* Now compute the size of the local subdomain n. */ \
1399c3518bceSDmitry Karpeev   *n = 0;                                               \
1400eec7e3faSDmitry Karpeev   if(*ylow_loc < *yhigh_loc) {                           \
1401af538404SDmitry Karpeev     PetscInt width = xright-xleft;                     \
1402eec7e3faSDmitry Karpeev     *n += width*(*yhigh_loc-*ylow_loc-1);                 \
1403eec7e3faSDmitry Karpeev     *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \
1404eec7e3faSDmitry Karpeev     *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \
1405af538404SDmitry Karpeev   }\
1406af538404SDmitry Karpeev }
1407af538404SDmitry Karpeev 
1408c3518bceSDmitry Karpeev 
1409eec7e3faSDmitry Karpeev 
1410eec7e3faSDmitry Karpeev #undef __FUNCT__
1411f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D"
1412b1a0a8a3SJed Brown /*@
1413f746d493SDmitry Karpeev    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1414b1a0a8a3SJed Brown    preconditioner for a two-dimensional problem on a regular grid.
1415b1a0a8a3SJed Brown 
1416af538404SDmitry Karpeev    Collective
1417b1a0a8a3SJed Brown 
1418b1a0a8a3SJed Brown    Input Parameters:
1419*06b43e7eSDmitry Karpeev +  M, N               - the global number of grid points in the x and y directions
1420af538404SDmitry Karpeev .  Mdomains, Ndomains - the global number of subdomains in the x and y directions
1421b1a0a8a3SJed Brown .  dof                - degrees of freedom per node
1422b1a0a8a3SJed Brown -  overlap            - overlap in mesh lines
1423b1a0a8a3SJed Brown 
1424b1a0a8a3SJed Brown    Output Parameters:
1425af538404SDmitry Karpeev +  Nsub - the number of local subdomains created
1426b1a0a8a3SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
1427b1a0a8a3SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
1428b1a0a8a3SJed Brown 
1429b1a0a8a3SJed Brown 
1430b1a0a8a3SJed Brown    Level: advanced
1431b1a0a8a3SJed Brown 
1432f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid
1433b1a0a8a3SJed Brown 
1434*06b43e7eSDmitry Karpeev .seealso: PCGASMSetTotalBlockCount(), PCGASMSetSubdomains(), PCGASMGetSubKSP(),
1435f746d493SDmitry Karpeev           PCGASMSetOverlap()
1436b1a0a8a3SJed Brown @*/
14377087cfbeSBarry Smith PetscErrorCode  PCGASMCreateSubdomains2D(PC pc, PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap, PetscInt *nsub,IS **is,IS **is_local)
1438b1a0a8a3SJed Brown {
1439b1a0a8a3SJed Brown   PetscErrorCode ierr;
14402bbb417fSDmitry Karpeev   PetscMPIInt    size, rank;
14412bbb417fSDmitry Karpeev   PetscInt       i, j;
14422bbb417fSDmitry Karpeev   PetscInt       maxheight, maxwidth;
14432bbb417fSDmitry Karpeev   PetscInt       xstart, xleft, xright, xleft_loc, xright_loc;
14442bbb417fSDmitry Karpeev   PetscInt       ystart, ylow,  yhigh,  ylow_loc,  yhigh_loc;
1445eec7e3faSDmitry Karpeev   PetscInt       x[2][2], y[2][2], n[2];
14462bbb417fSDmitry Karpeev   PetscInt       first, last;
14472bbb417fSDmitry Karpeev   PetscInt       nidx, *idx;
14482bbb417fSDmitry Karpeev   PetscInt       ii,jj,s,q,d;
1449af538404SDmitry Karpeev   PetscInt       k,kk;
14502bbb417fSDmitry Karpeev   PetscMPIInt    color;
14512bbb417fSDmitry Karpeev   MPI_Comm       comm, subcomm;
14525d38c402SBarry Smith   IS             **iis = 0;
1453d34fcf5fSBarry Smith 
1454b1a0a8a3SJed Brown   PetscFunctionBegin;
14552bbb417fSDmitry Karpeev   ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr);
14562bbb417fSDmitry Karpeev   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
14572bbb417fSDmitry Karpeev   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
14582bbb417fSDmitry Karpeev   ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr);
1459d34fcf5fSBarry 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) "
14602bbb417fSDmitry Karpeev 	     "does not respect the number of degrees of freedom per grid point %D", first, last, dof);
1461d34fcf5fSBarry Smith 
1462af538404SDmitry Karpeev   /* Determine the number of domains with nonzero intersections with the local ownership range. */
14632bbb417fSDmitry Karpeev   s = 0;
1464af538404SDmitry Karpeev   ystart = 0;
1465af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1466af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1467af538404SDmitry 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);
1468eec7e3faSDmitry Karpeev     /* Vertical domain limits with an overlap. */
1469eec7e3faSDmitry Karpeev     ylow = PetscMax(ystart - overlap,0);
1470eec7e3faSDmitry Karpeev     yhigh = PetscMin(ystart + maxheight + overlap,N);
1471b1a0a8a3SJed Brown     xstart = 0;
1472af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1473af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1474af538404SDmitry 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);
1475eec7e3faSDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1476eec7e3faSDmitry Karpeev       xleft   = PetscMax(xstart - overlap,0);
1477eec7e3faSDmitry Karpeev       xright  = PetscMin(xstart + maxwidth + overlap,M);
1478af538404SDmitry Karpeev       /*
1479af538404SDmitry Karpeev 	 Determine whether this subdomain intersects this processor's ownership range of pc->pmat.
1480af538404SDmitry Karpeev       */
1481c3518bceSDmitry Karpeev       PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1482af538404SDmitry Karpeev       if(nidx) {
1483af538404SDmitry Karpeev         ++s;
1484af538404SDmitry Karpeev       }
1485af538404SDmitry Karpeev       xstart += maxwidth;
1486af538404SDmitry Karpeev     }/* for(i = 0; i < Mdomains; ++i) */
1487af538404SDmitry Karpeev     ystart += maxheight;
1488af538404SDmitry Karpeev   }/* for(j = 0; j < Ndomains; ++j) */
1489af538404SDmitry Karpeev   /* Now we can allocate the necessary number of ISs. */
1490af538404SDmitry Karpeev   *nsub = s;
1491af538404SDmitry Karpeev   ierr = PetscMalloc((*nsub)*sizeof(IS*),is);CHKERRQ(ierr);
1492af538404SDmitry Karpeev   ierr = PetscMalloc((*nsub)*sizeof(IS*),is_local);CHKERRQ(ierr);
1493af538404SDmitry Karpeev   s = 0;
1494af538404SDmitry Karpeev   ystart = 0;
1495af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1496af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1497af538404SDmitry 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);
1498af538404SDmitry Karpeev     /* Vertical domain limits with an overlap. */
1499eec7e3faSDmitry Karpeev     y[0][0] = PetscMax(ystart - overlap,0);
1500eec7e3faSDmitry Karpeev     y[0][1] = PetscMin(ystart + maxheight + overlap,N);
1501eec7e3faSDmitry Karpeev     /* Vertical domain limits without an overlap. */
1502eec7e3faSDmitry Karpeev     y[1][0] = ystart;
1503eec7e3faSDmitry Karpeev     y[1][1] = PetscMin(ystart + maxheight,N);
1504eec7e3faSDmitry Karpeev     xstart = 0;
1505af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1506af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1507af538404SDmitry 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);
1508af538404SDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1509eec7e3faSDmitry Karpeev       x[0][0]  = PetscMax(xstart - overlap,0);
1510eec7e3faSDmitry Karpeev       x[0][1]  = PetscMin(xstart + maxwidth + overlap,M);
1511eec7e3faSDmitry Karpeev       /* Horizontal domain limits without an overlap. */
1512eec7e3faSDmitry Karpeev       x[1][0] = xstart;
1513eec7e3faSDmitry Karpeev       x[1][1] = PetscMin(xstart+maxwidth,M);
15142bbb417fSDmitry Karpeev       /*
15152bbb417fSDmitry Karpeev 	 Determine whether this domain intersects this processor's ownership range of pc->pmat.
15162bbb417fSDmitry Karpeev 	 Do this twice: first for the domains with overlaps, and once without.
15172bbb417fSDmitry Karpeev 	 During the first pass create the subcommunicators, and use them on the second pass as well.
15182bbb417fSDmitry Karpeev       */
15192bbb417fSDmitry Karpeev       for(q = 0; q < 2; ++q) {
15202bbb417fSDmitry Karpeev 	/*
15212bbb417fSDmitry Karpeev 	  domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
15222bbb417fSDmitry Karpeev 	  according to whether the domain with an overlap or without is considered.
15232bbb417fSDmitry Karpeev 	*/
15242bbb417fSDmitry Karpeev 	xleft = x[q][0]; xright = x[q][1];
15252bbb417fSDmitry Karpeev 	ylow  = y[q][0]; yhigh  = y[q][1];
1526c3518bceSDmitry Karpeev         PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1527af538404SDmitry Karpeev 	nidx *= dof;
1528eec7e3faSDmitry Karpeev         n[q] = nidx;
15292bbb417fSDmitry Karpeev         /*
1530eec7e3faSDmitry Karpeev          Based on the counted number of indices in the local domain *with an overlap*,
1531af538404SDmitry Karpeev          construct a subcommunicator of all the processors supporting this domain.
1532eec7e3faSDmitry Karpeev          Observe that a domain with an overlap might have nontrivial local support,
1533eec7e3faSDmitry Karpeev          while the domain without an overlap might not.  Hence, the decision to participate
1534eec7e3faSDmitry Karpeev          in the subcommunicator must be based on the domain with an overlap.
15352bbb417fSDmitry Karpeev          */
15362bbb417fSDmitry Karpeev 	if (q == 0) {
15372bbb417fSDmitry Karpeev 	  if(nidx) {
15382bbb417fSDmitry Karpeev 	    color = 1;
1539d34fcf5fSBarry Smith 	  } else {
15402bbb417fSDmitry Karpeev 	    color = MPI_UNDEFINED;
15412bbb417fSDmitry Karpeev 	  }
15422bbb417fSDmitry Karpeev 	  ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr);
15432bbb417fSDmitry Karpeev 	}
1544af538404SDmitry Karpeev         /*
1545eec7e3faSDmitry Karpeev          Proceed only if the number of local indices *with an overlap* is nonzero.
1546af538404SDmitry Karpeev          */
1547eec7e3faSDmitry Karpeev         if (n[0]) {
1548af538404SDmitry Karpeev           if(q == 0) {
1549eec7e3faSDmitry Karpeev             iis = is;
1550af538404SDmitry Karpeev           }
1551af538404SDmitry Karpeev           if (q == 1) {
1552af538404SDmitry Karpeev             /*
1553eec7e3faSDmitry Karpeev              The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1554af538404SDmitry Karpeev              Moreover, if the overlap is zero, the two ISs are identical.
1555af538404SDmitry Karpeev              */
1556b1a0a8a3SJed Brown             if (overlap == 0) {
1557eec7e3faSDmitry Karpeev               (*is_local)[s] = (*is)[s];
15582bbb417fSDmitry Karpeev               ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr);
15592bbb417fSDmitry Karpeev               continue;
1560d34fcf5fSBarry Smith             } else {
1561eec7e3faSDmitry Karpeev               iis = is_local;
1562eec7e3faSDmitry Karpeev               subcomm = ((PetscObject)(*is)[s])->comm;
15632bbb417fSDmitry Karpeev             }
1564af538404SDmitry Karpeev           }/* if(q == 1) */
1565eec7e3faSDmitry Karpeev           idx = PETSC_NULL;
15662bbb417fSDmitry Karpeev 	  ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1567eec7e3faSDmitry Karpeev           if(nidx) {
15682bbb417fSDmitry Karpeev             k    = 0;
15692bbb417fSDmitry Karpeev             for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1570af538404SDmitry Karpeev               PetscInt x0 = (jj==ylow_loc)?xleft_loc:xleft;
1571af538404SDmitry Karpeev               PetscInt x1 = (jj==yhigh_loc-1)?xright_loc:xright;
1572af538404SDmitry Karpeev               kk = dof*(M*jj + x0);
15732bbb417fSDmitry Karpeev               for (ii=x0; ii<x1; ++ii) {
15742bbb417fSDmitry Karpeev                 for(d = 0; d < dof; ++d) {
15752bbb417fSDmitry Karpeev                   idx[k++] = kk++;
1576b1a0a8a3SJed Brown                 }
1577b1a0a8a3SJed Brown               }
1578b1a0a8a3SJed Brown             }
1579eec7e3faSDmitry Karpeev           }
15802bbb417fSDmitry Karpeev 	  ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*iis)+s);CHKERRQ(ierr);
1581eec7e3faSDmitry Karpeev 	}/* if(n[0]) */
15822bbb417fSDmitry Karpeev       }/* for(q = 0; q < 2; ++q) */
1583eec7e3faSDmitry Karpeev       if(n[0]) {
15842bbb417fSDmitry Karpeev         ++s;
1585b1a0a8a3SJed Brown       }
1586af538404SDmitry Karpeev       xstart += maxwidth;
1587af538404SDmitry Karpeev     }/* for(i = 0; i < Mdomains; ++i) */
15882bbb417fSDmitry Karpeev     ystart += maxheight;
1589af538404SDmitry Karpeev   }/* for(j = 0; j < Ndomains; ++j) */
1590b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1591b1a0a8a3SJed Brown }
1592b1a0a8a3SJed Brown 
1593b1a0a8a3SJed Brown #undef __FUNCT__
1594*06b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains"
1595b1a0a8a3SJed Brown /*@C
1596*06b43e7eSDmitry Karpeev     PCGASMGetSubdomains - Gets the subdomains supported on this processor
1597*06b43e7eSDmitry Karpeev     for the additive Schwarz preconditioner.
1598b1a0a8a3SJed Brown 
1599b1a0a8a3SJed Brown     Not Collective
1600b1a0a8a3SJed Brown 
1601b1a0a8a3SJed Brown     Input Parameter:
1602b1a0a8a3SJed Brown .   pc - the preconditioner context
1603b1a0a8a3SJed Brown 
1604b1a0a8a3SJed Brown     Output Parameters:
1605b1a0a8a3SJed Brown +   n        - the number of subdomains for this processor (default value = 1)
1606*06b43e7eSDmitry Karpeev .   is       - the index sets that define the subdomains supported on this processor
1607b1a0a8a3SJed Brown -   is_local - the index sets that define the local part of the subdomains for this processor (can be PETSC_NULL)
1608b1a0a8a3SJed Brown 
1609b1a0a8a3SJed Brown 
1610b1a0a8a3SJed Brown     Notes:
1611b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
1612b1a0a8a3SJed Brown 
1613b1a0a8a3SJed Brown     Level: advanced
1614b1a0a8a3SJed Brown 
1615*06b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz
1616b1a0a8a3SJed Brown 
1617*06b43e7eSDmitry Karpeev .seealso: PCGASMSetTotalBlockCount(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
1618*06b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubmatrices()
1619b1a0a8a3SJed Brown @*/
1620*06b43e7eSDmitry Karpeev PetscErrorCode  PCGASMGetSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
1621b1a0a8a3SJed Brown {
1622f746d493SDmitry Karpeev   PC_GASM         *osm;
1623b1a0a8a3SJed Brown   PetscErrorCode ierr;
1624b1a0a8a3SJed Brown   PetscBool      match;
1625b1a0a8a3SJed Brown 
1626b1a0a8a3SJed Brown   PetscFunctionBegin;
1627b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1628b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1629b1a0a8a3SJed Brown   if (is) PetscValidPointer(is,3);
1630251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1631b1a0a8a3SJed Brown   if (!match) {
1632b1a0a8a3SJed Brown     if (n)  *n  = 0;
1633b1a0a8a3SJed Brown     if (is) *is = PETSC_NULL;
1634b1a0a8a3SJed Brown   } else {
1635f746d493SDmitry Karpeev     osm = (PC_GASM*)pc->data;
1636ab3e923fSDmitry Karpeev     if (n)  *n  = osm->n;
1637b1a0a8a3SJed Brown     if (is) *is = osm->is;
1638b1a0a8a3SJed Brown     if (is_local) *is_local = osm->is_local;
1639b1a0a8a3SJed Brown   }
1640b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1641b1a0a8a3SJed Brown }
1642b1a0a8a3SJed Brown 
1643b1a0a8a3SJed Brown #undef __FUNCT__
1644*06b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices"
1645b1a0a8a3SJed Brown /*@C
1646*06b43e7eSDmitry Karpeev     PCGASMGetSubmatrices - Gets the local submatrices (for this processor
1647b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1648b1a0a8a3SJed Brown 
1649b1a0a8a3SJed Brown     Not Collective
1650b1a0a8a3SJed Brown 
1651b1a0a8a3SJed Brown     Input Parameter:
1652b1a0a8a3SJed Brown .   pc - the preconditioner context
1653b1a0a8a3SJed Brown 
1654b1a0a8a3SJed Brown     Output Parameters:
1655b1a0a8a3SJed Brown +   n   - the number of matrices for this processor (default value = 1)
1656b1a0a8a3SJed Brown -   mat - the matrices
1657b1a0a8a3SJed Brown 
1658*06b43e7eSDmitry Karpeev     Notes: matrices returned by this routine have the same communicators as the index sets (IS)
1659*06b43e7eSDmitry Karpeev            used to define subdomains in PCGASMSetSubdomains(), or PETSC_COMM_SELF, if the
1660*06b43e7eSDmitry Karpeev            subdomains were defined using PCGASMSetTotalBlockCount().
1661b1a0a8a3SJed Brown     Level: advanced
1662b1a0a8a3SJed Brown 
1663f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi
1664b1a0a8a3SJed Brown 
1665*06b43e7eSDmitry Karpeev .seealso: PCGASMSetTotalBlockCount(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
1666*06b43e7eSDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains()
1667b1a0a8a3SJed Brown @*/
1668*06b43e7eSDmitry Karpeev PetscErrorCode  PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1669b1a0a8a3SJed Brown {
1670f746d493SDmitry Karpeev   PC_GASM         *osm;
1671b1a0a8a3SJed Brown   PetscErrorCode ierr;
1672b1a0a8a3SJed Brown   PetscBool      match;
1673b1a0a8a3SJed Brown 
1674b1a0a8a3SJed Brown   PetscFunctionBegin;
1675b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1676b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1677b1a0a8a3SJed Brown   if (mat) PetscValidPointer(mat,3);
1678b1a0a8a3SJed Brown   if (!pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1679251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1680*06b43e7eSDmitry Karpeev   if (!match) SETERRQ2(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name);
1681f746d493SDmitry Karpeev   osm = (PC_GASM*)pc->data;
1682ab3e923fSDmitry Karpeev   if (n)   *n   = osm->n;
1683b1a0a8a3SJed Brown   if (mat) *mat = osm->pmat;
1684*06b43e7eSDmitry Karpeev 
1685b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1686b1a0a8a3SJed Brown }
1687