xref: /petsc/src/ksp/pc/impls/gasm/gasm.c (revision 167188bf39e6f1466bb0de974b5078c45f49fb54)
1b1a0a8a3SJed Brown #define PETSCKSP_DLL
2b1a0a8a3SJed Brown 
3b1a0a8a3SJed Brown /*
4f746d493SDmitry Karpeev   This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation.
544fe18e5SDmitry Karpeev   In this version each processor may have any number of subdomains and a subdomain may live on multiple
6f746d493SDmitry Karpeev   processors.
7b1a0a8a3SJed Brown 
8af538404SDmitry Karpeev        N    - total number of subdomains on all processors
9f746d493SDmitry Karpeev        n    - actual number of subdomains on this processor
10f746d493SDmitry Karpeev        nmax - maximum number of subdomains per processor
11b1a0a8a3SJed Brown */
12b1a0a8a3SJed Brown #include "private/pcimpl.h"     /*I "petscpc.h" I*/
13b1a0a8a3SJed Brown 
14b1a0a8a3SJed Brown typedef struct {
15f746d493SDmitry Karpeev   PetscInt   N,n,nmax;
16b1a0a8a3SJed Brown   PetscInt   overlap;             /* overlap requested by user */
17b1a0a8a3SJed Brown   KSP        *ksp;                /* linear solvers for each block */
18f746d493SDmitry Karpeev   Vec        gx,gy;               /* Merged work vectors */
19f746d493SDmitry Karpeev   Vec        *x,*y;               /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */
20f746d493SDmitry Karpeev   IS         gis, gis_local;      /* merged ISs */
21f746d493SDmitry Karpeev   VecScatter grestriction;        /* merged restriction */
22f746d493SDmitry Karpeev   VecScatter gprolongation;       /* merged prolongation */
23b1a0a8a3SJed Brown   IS         *is;                 /* index set that defines each overlapping subdomain */
24f746d493SDmitry Karpeev   IS         *is_local;           /* index set that defines each local subdomain (same as subdomain with the overlap removed); may be NULL */
25f746d493SDmitry Karpeev   Mat        *pmat;               /* subdomain block matrices */
26f746d493SDmitry Karpeev   PCGASMType  type;               /* use reduced interpolation, restriction or both */
27b1a0a8a3SJed Brown   PetscBool  type_set;            /* if user set this value (so won't change it for symmetric problems) */
28b1a0a8a3SJed Brown   PetscBool  same_local_solves;   /* flag indicating whether all local solvers are same */
29b1a0a8a3SJed Brown   PetscBool  sort_indices;        /* flag to sort subdomain indices */
30f746d493SDmitry Karpeev } PC_GASM;
31b1a0a8a3SJed Brown 
32b1a0a8a3SJed Brown #undef __FUNCT__
33af538404SDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains"
34af538404SDmitry Karpeev static PetscErrorCode  PCGASMPrintSubdomains(PC pc)
35af538404SDmitry Karpeev {
36af538404SDmitry Karpeev   PC_GASM         *osm  = (PC_GASM*)pc->data;
37af538404SDmitry Karpeev   const char     *prefix;
38af538404SDmitry Karpeev   char           fname[PETSC_MAX_PATH_LEN+1];
39af538404SDmitry Karpeev   PetscViewer    viewer, sviewer;
40af538404SDmitry Karpeev   PetscInt       i,j,nidx;
41af538404SDmitry Karpeev   const PetscInt *idx;
42af538404SDmitry Karpeev   char           *cidx;
43af538404SDmitry Karpeev   PetscBool      found;
44eec7e3faSDmitry Karpeev   PetscMPIInt    size, rank;
45af538404SDmitry Karpeev   PetscErrorCode ierr;
46af538404SDmitry Karpeev 
47af538404SDmitry Karpeev   PetscFunctionBegin;
48eec7e3faSDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm, &size); CHKERRQ(ierr);
49eec7e3faSDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm, &rank); CHKERRQ(ierr);
50af538404SDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
51af538404SDmitry Karpeev   ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr);
52af538404SDmitry Karpeev   if (!found) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
53af538404SDmitry Karpeev   for (i=0;i<osm->n;++i) {
54af538404SDmitry Karpeev     ierr = PetscViewerASCIIOpen(((PetscObject)((osm->is)[i]))->comm,fname,&viewer);CHKERRQ(ierr);
55af538404SDmitry Karpeev     ierr = ISGetLocalSize(osm->is[i], &nidx); CHKERRQ(ierr);
56eec7e3faSDmitry Karpeev     /*
57eec7e3faSDmitry Karpeev      No more than 15 characters per index plus a space.
58eec7e3faSDmitry Karpeev      PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
59eec7e3faSDmitry Karpeev      in case nidx == 0. That will take care of the space for the trailing '\0' as well.
60eec7e3faSDmitry Karpeev      For nidx == 0, the whole string 16 '\0'.
61eec7e3faSDmitry Karpeev      */
62eec7e3faSDmitry Karpeev     ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx); CHKERRQ(ierr);
63af538404SDmitry Karpeev     ierr = ISGetIndices(osm->is[i], &idx);    CHKERRQ(ierr);
64eec7e3faSDmitry Karpeev     ierr = PetscViewerStringOpen(((PetscObject)(osm->is[i]))->comm, cidx, 16*(nidx+1), &sviewer); CHKERRQ(ierr);
65af538404SDmitry Karpeev     for(j = 0; j < nidx; ++j) {
66af538404SDmitry Karpeev       ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]); CHKERRQ(ierr);
67af538404SDmitry Karpeev     }
68af538404SDmitry Karpeev     ierr = PetscViewerDestroy(sviewer); CHKERRQ(ierr);
69af538404SDmitry Karpeev     ierr = ISRestoreIndices(osm->is[i],&idx); CHKERRQ(ierr);
70eec7e3faSDmitry Karpeev 
71af538404SDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer, "Subdomain with overlap\n"); CHKERRQ(ierr);
72af538404SDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx); CHKERRQ(ierr);
73af538404SDmitry Karpeev     ierr = PetscViewerFlush(viewer); CHKERRQ(ierr);
74af538404SDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer, "\n"); CHKERRQ(ierr);
75af538404SDmitry Karpeev     ierr = PetscFree(cidx); CHKERRQ(ierr);
76af538404SDmitry Karpeev     ierr = PetscViewerDestroy(viewer); CHKERRQ(ierr);
77af538404SDmitry Karpeev     if (osm->is_local) {
78af538404SDmitry Karpeev       ierr = PetscViewerASCIIOpen(((PetscObject)((osm->is)[i]))->comm,fname,&viewer);CHKERRQ(ierr);
79af538404SDmitry Karpeev       ierr = ISGetLocalSize(osm->is_local[i], &nidx); CHKERRQ(ierr);
80eec7e3faSDmitry Karpeev       /*
81eec7e3faSDmitry Karpeev        No more than 15 characters per index plus a space.
82eec7e3faSDmitry Karpeev        PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx,
83eec7e3faSDmitry Karpeev        in case nidx == 0. That will take care of the space for the trailing '\0' as well.
84eec7e3faSDmitry Karpeev        For nidx == 0, the whole string 16 '\0'.
85eec7e3faSDmitry Karpeev        */
86eec7e3faSDmitry Karpeev       ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx); CHKERRQ(ierr);
87af538404SDmitry Karpeev       ierr = ISGetIndices(osm->is_local[i], &idx);    CHKERRQ(ierr);
88eec7e3faSDmitry Karpeev       ierr = PetscViewerStringOpen(((PetscObject)(osm->is_local[i]))->comm, cidx, 16*(nidx+1), &sviewer); CHKERRQ(ierr);
89af538404SDmitry Karpeev       for(j = 0; j < nidx; ++j) {
90af538404SDmitry Karpeev         ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]); CHKERRQ(ierr);
91af538404SDmitry Karpeev       }
92af538404SDmitry Karpeev       ierr = PetscViewerDestroy(sviewer); CHKERRQ(ierr);
93af538404SDmitry Karpeev       ierr = ISRestoreIndices(osm->is_local[i],&idx); CHKERRQ(ierr);
94eec7e3faSDmitry Karpeev 
95af538404SDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Subdomain without overlap\n"); CHKERRQ(ierr);
96af538404SDmitry Karpeev       ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx); CHKERRQ(ierr);
97af538404SDmitry Karpeev       ierr = PetscViewerFlush(viewer); CHKERRQ(ierr);
98af538404SDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "\n"); CHKERRQ(ierr);
99af538404SDmitry Karpeev       ierr = PetscFree(cidx); CHKERRQ(ierr);
100af538404SDmitry Karpeev       ierr = PetscViewerDestroy(viewer); CHKERRQ(ierr);
101af538404SDmitry Karpeev     }
102af538404SDmitry Karpeev   }
103af538404SDmitry Karpeev   PetscFunctionReturn(0);
104af538404SDmitry Karpeev }
105af538404SDmitry Karpeev 
106af538404SDmitry Karpeev 
107af538404SDmitry Karpeev #undef __FUNCT__
108f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM"
109f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer)
110b1a0a8a3SJed Brown {
111f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
112af538404SDmitry Karpeev   const char     *prefix;
113b1a0a8a3SJed Brown   PetscErrorCode ierr;
114af538404SDmitry Karpeev   PetscMPIInt    rank, size;
115b1a0a8a3SJed Brown   PetscInt       i,bsz;
116af538404SDmitry Karpeev   PetscBool      iascii,isstring, print_subdomains;
117b1a0a8a3SJed Brown   PetscViewer    sviewer;
118b1a0a8a3SJed Brown 
119b1a0a8a3SJed Brown 
120b1a0a8a3SJed Brown   PetscFunctionBegin;
121af538404SDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm, &size); CHKERRQ(ierr);
122af538404SDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm, &rank); CHKERRQ(ierr);
123b1a0a8a3SJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
124b1a0a8a3SJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
125af538404SDmitry Karpeev   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
126af538404SDmitry Karpeev   ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&print_subdomains,PETSC_NULL);CHKERRQ(ierr);
127b1a0a8a3SJed Brown   if (iascii) {
128af538404SDmitry Karpeev     char overlaps[256] = "user-defined overlap",subdomains[256] = "total subdomains set";
129b1a0a8a3SJed Brown     if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof overlaps,"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);}
130af538404SDmitry Karpeev     if (osm->nmax > 0)     {ierr = PetscSNPrintf(subdomains,sizeof subdomains,"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr);}
131eec7e3faSDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d:%d] number of locally-supported subdomains = %D\n",(int)rank,(int)size,osm->n);CHKERRQ(ierr);
132af538404SDmitry Karpeev     ierr = PetscViewerFlush(viewer); CHKERRQ(ierr);
133af538404SDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"  Generalized additive Schwarz: %s, %s\n",subdomains,overlaps);CHKERRQ(ierr);
134f746d493SDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"  Generalized additive Schwarz: restriction/interpolation type - %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr);
135b1a0a8a3SJed Brown     ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
136b1a0a8a3SJed Brown     if (osm->same_local_solves) {
137b1a0a8a3SJed Brown       if (osm->ksp) {
138af538404SDmitry Karpeev         ierr = PetscViewerASCIIPrintf(viewer,"  Local solve is same for all subdomains, in the following KSP and PC objects:\n");CHKERRQ(ierr);
139b1a0a8a3SJed Brown         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
140b1a0a8a3SJed Brown         if (!rank) {
141b1a0a8a3SJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
142b1a0a8a3SJed Brown           ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);
143b1a0a8a3SJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
144b1a0a8a3SJed Brown         }
145b1a0a8a3SJed Brown         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
146b1a0a8a3SJed Brown       }
147b1a0a8a3SJed Brown     } else {
148b1a0a8a3SJed Brown       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
149af538404SDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each subdomain is in the following KSP and PC objects:\n");CHKERRQ(ierr);
150b1a0a8a3SJed Brown       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
151b1a0a8a3SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
152f746d493SDmitry Karpeev       for (i=0; i<osm->nmax; i++) {
153b1a0a8a3SJed Brown         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
154f746d493SDmitry Karpeev         if (i < osm->n) {
155b1a0a8a3SJed Brown           ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr);
156af538404SDmitry Karpeev           ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d:%d] local subdomain number %D, size = %D\n",(int)rank,(int)size,i,bsz);CHKERRQ(ierr);
157b1a0a8a3SJed Brown           ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr);
158b1a0a8a3SJed Brown           ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
159b1a0a8a3SJed Brown         }
160b1a0a8a3SJed Brown         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
161b1a0a8a3SJed Brown       }
162b1a0a8a3SJed Brown       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
163b1a0a8a3SJed Brown       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
164b1a0a8a3SJed Brown     }
165b1a0a8a3SJed Brown   } else if (isstring) {
166af538404SDmitry Karpeev     ierr = PetscViewerStringSPrintf(viewer," subdomains=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCGASMTypes[osm->type]);CHKERRQ(ierr);
167b1a0a8a3SJed Brown     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
168b1a0a8a3SJed Brown     if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);}
169b1a0a8a3SJed Brown     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
170b1a0a8a3SJed Brown   } else {
171f746d493SDmitry Karpeev     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCGASM",((PetscObject)viewer)->type_name);
172b1a0a8a3SJed Brown   }
173b1a0a8a3SJed Brown   ierr = PetscViewerFlush(viewer); CHKERRQ(ierr);
174af538404SDmitry Karpeev   if(print_subdomains) {
175af538404SDmitry Karpeev     ierr = PCGASMPrintSubdomains(pc); CHKERRQ(ierr);
176af538404SDmitry Karpeev   }
177b1a0a8a3SJed Brown   PetscFunctionReturn(0);
178b1a0a8a3SJed Brown }
179b1a0a8a3SJed Brown 
180af538404SDmitry Karpeev 
181af538404SDmitry Karpeev 
182af538404SDmitry Karpeev 
183af538404SDmitry Karpeev 
184b1a0a8a3SJed Brown #undef __FUNCT__
185f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM"
186f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc)
187b1a0a8a3SJed Brown {
188f746d493SDmitry Karpeev   PC_GASM         *osm  = (PC_GASM*)pc->data;
189b1a0a8a3SJed Brown   PetscErrorCode ierr;
190b1a0a8a3SJed Brown   PetscBool      symset,flg;
19188389755SJed Brown   PetscInt       i;
192c10bc1a1SDmitry Karpeev   PetscMPIInt    rank, size;
193b1a0a8a3SJed Brown   MatReuse       scall = MAT_REUSE_MATRIX;
194b1a0a8a3SJed Brown   KSP            ksp;
195b1a0a8a3SJed Brown   PC             subpc;
196b1a0a8a3SJed Brown   const char     *prefix,*pprefix;
197f746d493SDmitry Karpeev   PetscInt       dn;       /* Number of indices in a single subdomain assigned to this processor. */
198f746d493SDmitry Karpeev   const PetscInt *didx;    /* Indices from a single subdomain assigned to this processor. */
199f746d493SDmitry Karpeev   PetscInt       ddn;      /* Number of indices in all subdomains assigned to this processor. */
200f746d493SDmitry Karpeev   PetscInt       *ddidx;   /* Indices of all subdomains assigned to this processor. */
201f746d493SDmitry Karpeev   IS             gid;      /* Identity IS of the size of all subdomains assigned to this processor. */
202f746d493SDmitry Karpeev   Vec            x,y;
203f746d493SDmitry Karpeev   PetscScalar    *gxarray, *gyarray;
2046d98aee9SDmitry Karpeev   PetscInt       gfirst, glast;
205b1a0a8a3SJed Brown 
206b1a0a8a3SJed Brown   PetscFunctionBegin;
207c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
208c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
209b1a0a8a3SJed Brown   if (!pc->setupcalled) {
210b1a0a8a3SJed Brown 
211b1a0a8a3SJed Brown     if (!osm->type_set) {
212b1a0a8a3SJed Brown       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
213f746d493SDmitry Karpeev       if (symset && flg) { osm->type = PC_GASM_BASIC; }
214b1a0a8a3SJed Brown     }
215b1a0a8a3SJed Brown 
21644fe18e5SDmitry Karpeev     if (osm->N == PETSC_DECIDE && osm->n < 1) {
217b1a0a8a3SJed Brown       /* no subdomains given, use one per processor */
21844fe18e5SDmitry Karpeev       osm->nmax = osm->n = 1;
219b1a0a8a3SJed Brown       ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
220f746d493SDmitry Karpeev       osm->N = size;
221f746d493SDmitry Karpeev     } else if (osm->N == PETSC_DECIDE) {
222b1a0a8a3SJed Brown       PetscInt inwork[2], outwork[2];
223f746d493SDmitry Karpeev       /* determine global number of subdomains and the max number of local subdomains */
224f746d493SDmitry Karpeev       inwork[0] = inwork[1] = osm->n;
225b1a0a8a3SJed Brown       ierr = MPI_Allreduce(inwork,outwork,1,MPIU_2INT,PetscMaxSum_Op,((PetscObject)pc)->comm);CHKERRQ(ierr);
226f746d493SDmitry Karpeev       osm->nmax = outwork[0];
227f746d493SDmitry Karpeev       osm->N    = outwork[1];
228b1a0a8a3SJed Brown     }
229b1a0a8a3SJed Brown     if (!osm->is){ /* create the index sets */
230f746d493SDmitry Karpeev       ierr = PCGASMCreateSubdomains(pc->pmat,osm->n,&osm->is);CHKERRQ(ierr);
231b1a0a8a3SJed Brown     }
23293c1ec32SDmitry Karpeev     if (!osm->is_local) {
23393c1ec32SDmitry Karpeev       /*
23493c1ec32SDmitry Karpeev 	 This indicates that osm->is should define a nonoverlapping decomposition
23593c1ec32SDmitry Karpeev 	 (there is no way to really guarantee that if subdomains are set by the user through PCGASMSetLocalSubdomains,
23693c1ec32SDmitry Karpeev 	  but the assumption is that either the user does the right thing, or subdomains in ossm->is have been created
23793c1ec32SDmitry Karpeev 	  via PCGASMCreateSubdomains, which guarantees a nonoverlapping decomposition).
23893c1ec32SDmitry Karpeev 	 Therefore, osm->is will be used to define osm->is_local.
23993c1ec32SDmitry Karpeev 	 If a nonzero overlap has been requested by the user, then osm->is will be expanded and will overlap,
24093c1ec32SDmitry Karpeev 	 so osm->is_local should obtain a copy of osm->is while they are still (presumably) nonoverlapping.
24193c1ec32SDmitry Karpeev 	 Otherwise (no overlap has been requested), osm->is_local are simply aliases for osm->is.
24293c1ec32SDmitry Karpeev       */
243f746d493SDmitry Karpeev       ierr = PetscMalloc(osm->n*sizeof(IS),&osm->is_local);CHKERRQ(ierr);
244f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
245b1a0a8a3SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
246b1a0a8a3SJed Brown           ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr);
247b1a0a8a3SJed Brown           ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr);
248b1a0a8a3SJed Brown         } else {
249b1a0a8a3SJed Brown           ierr = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr);
250b1a0a8a3SJed Brown           osm->is_local[i] = osm->is[i];
251b1a0a8a3SJed Brown         }
252b1a0a8a3SJed Brown       }
253b1a0a8a3SJed Brown     }
25493c1ec32SDmitry Karpeev     /* Beyond this point osm->is_local is not null. */
25593c1ec32SDmitry Karpeev     if (osm->overlap > 0) {
25693c1ec32SDmitry Karpeev       /* Extend the "overlapping" regions by a number of steps */
25793c1ec32SDmitry Karpeev       ierr = MatIncreaseOverlap(pc->pmat,osm->n,osm->is,osm->overlap);CHKERRQ(ierr);
25893c1ec32SDmitry Karpeev     }
259b1a0a8a3SJed Brown     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
260b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
261f746d493SDmitry Karpeev     ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&flg,PETSC_NULL);CHKERRQ(ierr);
262f746d493SDmitry Karpeev     if (flg) { ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); }
263b1a0a8a3SJed Brown 
264b1a0a8a3SJed Brown     if (osm->sort_indices) {
265f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
266b1a0a8a3SJed Brown         ierr = ISSort(osm->is[i]);CHKERRQ(ierr);
267b1a0a8a3SJed Brown 	ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr);
268b1a0a8a3SJed Brown       }
269b1a0a8a3SJed Brown     }
270f746d493SDmitry Karpeev     /* Merge the ISs, create merged vectors and scatter contexts. */
27193c1ec32SDmitry Karpeev     /* Restriction ISs. */
272f746d493SDmitry Karpeev     ddn = 0;
273f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
274f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn);          CHKERRQ(ierr);
275f746d493SDmitry Karpeev       ddn += dn;
276f746d493SDmitry Karpeev     }
277f746d493SDmitry Karpeev     ierr = PetscMalloc(ddn*sizeof(PetscInt), &ddidx); CHKERRQ(ierr);
278f746d493SDmitry Karpeev     ddn = 0;
279f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
280f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn); CHKERRQ(ierr);
281f746d493SDmitry Karpeev       ierr = ISGetIndices(osm->is[i],&didx); CHKERRQ(ierr);
282f746d493SDmitry Karpeev       ierr = PetscMemcpy(ddidx+ddn, didx, sizeof(PetscInt)*dn); CHKERRQ(ierr);
283f746d493SDmitry Karpeev       ierr = ISRestoreIndices(osm->is[i], &didx); CHKERRQ(ierr);
28493c1ec32SDmitry Karpeev       ddn += dn;
285f746d493SDmitry Karpeev     }
286f746d493SDmitry Karpeev     ierr = ISCreateGeneral(((PetscObject)(pc))->comm, ddn, ddidx, PETSC_OWN_POINTER, &osm->gis); CHKERRQ(ierr);
287f746d493SDmitry Karpeev     ierr = MatGetVecs(pc->pmat,&x,&y);CHKERRQ(ierr);
288f746d493SDmitry Karpeev     ierr = VecCreateMPI(((PetscObject)pc)->comm, ddn, PETSC_DECIDE, &osm->gx); CHKERRQ(ierr);
289f746d493SDmitry Karpeev     ierr = VecDuplicate(osm->gx,&osm->gy);                                     CHKERRQ(ierr);
2906d98aee9SDmitry Karpeev     ierr = VecGetOwnershipRange(osm->gx, &gfirst, &glast);                     CHKERRQ(ierr);
2916d98aee9SDmitry Karpeev     ierr = ISCreateStride(((PetscObject)pc)->comm,ddn,gfirst,1, &gid);              CHKERRQ(ierr);
292f746d493SDmitry Karpeev     ierr = VecScatterCreate(x,osm->gis,osm->gx,gid, &(osm->grestriction));     CHKERRQ(ierr);
293f746d493SDmitry Karpeev     ierr = ISDestroy(gid);                                                     CHKERRQ(ierr);
29493c1ec32SDmitry Karpeev     /* Prolongation ISs */
29581a5c4aaSDmitry Karpeev     { PetscInt       dn_local;       /* Number of indices in the local part of a single domain assigned to this processor. */
296f746d493SDmitry Karpeev       const PetscInt *didx_local;    /* Global indices from the local part of a single domain assigned to this processor. */
297f746d493SDmitry Karpeev       PetscInt       ddn_local;      /* Number of indices in the local part of the disjoint union all domains assigned to this processor. */
298f746d493SDmitry Karpeev       PetscInt       *ddidx_local;   /* Global indices of the local part of the disjoint union of all domains assigned to this processor. */
299f746d493SDmitry Karpeev       /**/
30093c1ec32SDmitry Karpeev       ISLocalToGlobalMapping ltog;          /* Map from global to local indices on the disjoint union of subdomains: "local" ind's run from 0 to ddn-1. */
30193c1ec32SDmitry Karpeev       PetscInt              *ddidx_llocal;  /* Mapped local indices of the disjoint union of local parts of subdomains. */
302f746d493SDmitry 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. */
303f746d493SDmitry Karpeev       IS                     gis_llocal;    /* IS with ddidx_llocal indices. */
304f746d493SDmitry Karpeev       PetscInt               j;
305f746d493SDmitry Karpeev       ddn_local = 0;
306f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
307f746d493SDmitry Karpeev 	ierr = ISGetLocalSize(osm->is_local[i],&dn_local); CHKERRQ(ierr);
308f746d493SDmitry Karpeev 	ddn_local += dn_local;
309f746d493SDmitry Karpeev       }
310f746d493SDmitry Karpeev       ierr = PetscMalloc(ddn_local*sizeof(PetscInt), &ddidx_local); CHKERRQ(ierr);
311f746d493SDmitry Karpeev       ddn_local = 0;
312f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
313f746d493SDmitry Karpeev 	ierr = ISGetLocalSize(osm->is_local[i],&dn_local); CHKERRQ(ierr);
314f746d493SDmitry Karpeev 	ierr = ISGetIndices(osm->is_local[i],&didx_local); CHKERRQ(ierr);
315f746d493SDmitry Karpeev 	ierr = PetscMemcpy(ddidx_local+ddn_local, didx_local, sizeof(PetscInt)*dn_local); CHKERRQ(ierr);
316f746d493SDmitry Karpeev 	ierr = ISRestoreIndices(osm->is_local[i], &didx_local); CHKERRQ(ierr);
317f746d493SDmitry Karpeev 	ddn_local += dn_local;
318f746d493SDmitry Karpeev       }
319ab3e923fSDmitry Karpeev       ierr = PetscMalloc(sizeof(PetscInt)*ddn_local, &ddidx_llocal); CHKERRQ(ierr);
320f746d493SDmitry Karpeev       ierr = ISCreateGeneral(((PetscObject)pc)->comm, ddn_local, ddidx_local, PETSC_OWN_POINTER, &(osm->gis_local)); CHKERRQ(ierr);
321f746d493SDmitry Karpeev       ierr = ISLocalToGlobalMappingCreateIS(osm->gis,&ltog);CHKERRQ(ierr);
32281a5c4aaSDmitry Karpeev       ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,ddn_local,ddidx_local,&ddn_llocal,ddidx_llocal);CHKERRQ(ierr);
323b1a0a8a3SJed Brown       ierr = ISLocalToGlobalMappingDestroy(ltog);CHKERRQ(ierr);
32493c1ec32SDmitry 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);
325f746d493SDmitry Karpeev       /* Now convert these localized indices into the global indices into the merged output vector. */
3266d98aee9SDmitry Karpeev       ierr = VecGetOwnershipRange(osm->gy, &gfirst, &glast); CHKERRQ(ierr);
327f746d493SDmitry Karpeev       for(j=0; j < ddn_llocal; ++j) {
3286d98aee9SDmitry Karpeev 	ddidx_llocal[j] += gfirst;
329f746d493SDmitry Karpeev       }
330f746d493SDmitry Karpeev       ierr = ISCreateGeneral(PETSC_COMM_SELF,ddn_llocal,ddidx_llocal,PETSC_OWN_POINTER,&gis_llocal); CHKERRQ(ierr);
331f746d493SDmitry Karpeev       ierr = VecScatterCreate(y,osm->gis_local,osm->gy,gis_llocal,&osm->gprolongation);              CHKERRQ(ierr);
332f746d493SDmitry Karpeev       ierr = ISDestroy(gis_llocal);CHKERRQ(ierr);
333b1a0a8a3SJed Brown     }
334f746d493SDmitry Karpeev     /* Create the subdomain work vectors. */
335f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->x);CHKERRQ(ierr);
336f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->y);CHKERRQ(ierr);
3376d98aee9SDmitry Karpeev     ierr = VecGetOwnershipRange(osm->gx, &gfirst, &glast);CHKERRQ(ierr);
338f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gx, &gxarray);                                     CHKERRQ(ierr);
339f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gy, &gyarray);                                     CHKERRQ(ierr);
340f746d493SDmitry Karpeev     for (i=0, ddn=0; i<osm->n; ++i, ddn += dn) {
34186b96d36SDmitry Karpeev       PetscInt dN;
342f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn);CHKERRQ(ierr);
34386b96d36SDmitry Karpeev       ierr = ISGetSize(osm->is[i],&dN);CHKERRQ(ierr);
34486b96d36SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->is[i]))->comm,dn,dN,gxarray+ddn,&osm->x[i]);CHKERRQ(ierr);
34586b96d36SDmitry Karpeev       ierr = VecCreateMPIWithArray(((PetscObject)(osm->is[i]))->comm,dn,dN,gyarray+ddn,&osm->y[i]);CHKERRQ(ierr);
346b1a0a8a3SJed Brown     }
347f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gx, &gxarray); CHKERRQ(ierr);
348f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gy, &gyarray); CHKERRQ(ierr);
349f746d493SDmitry Karpeev     ierr = VecDestroy(x); CHKERRQ(ierr);
350f746d493SDmitry Karpeev     ierr = VecDestroy(y); CHKERRQ(ierr);
351b1a0a8a3SJed Brown     /* Create the local solvers */
352f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(KSP *),&osm->ksp);CHKERRQ(ierr);
35344fe18e5SDmitry Karpeev     for (i=0; i<osm->n; i++) { /* KSPs are local */
35486b96d36SDmitry Karpeev       ierr = KSPCreate(((PetscObject)(osm->is[i]))->comm,&ksp);CHKERRQ(ierr);
355b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr);
356b1a0a8a3SJed Brown       ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
357b1a0a8a3SJed Brown       ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
358b1a0a8a3SJed Brown       ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
359b1a0a8a3SJed Brown       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
360b1a0a8a3SJed Brown       ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
361b1a0a8a3SJed Brown       ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
362b1a0a8a3SJed Brown       osm->ksp[i] = ksp;
363b1a0a8a3SJed Brown     }
364b1a0a8a3SJed Brown     scall = MAT_INITIAL_MATRIX;
365b1a0a8a3SJed Brown 
366b1a0a8a3SJed Brown   } else {
367b1a0a8a3SJed Brown     /*
368b1a0a8a3SJed Brown        Destroy the blocks from the previous iteration
369b1a0a8a3SJed Brown     */
370b1a0a8a3SJed Brown     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
371f746d493SDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
372b1a0a8a3SJed Brown       scall = MAT_INITIAL_MATRIX;
373b1a0a8a3SJed Brown     }
374b1a0a8a3SJed Brown   }
375b1a0a8a3SJed Brown 
376b1a0a8a3SJed Brown   /*
377f746d493SDmitry Karpeev      Extract out the submatrices.
378b1a0a8a3SJed Brown   */
37981a5c4aaSDmitry Karpeev   if(size > 1) {
3806d42056aSDmitry Karpeev     ierr = MatGetSubMatricesParallel(pc->pmat,osm->n,osm->is, osm->is,scall,&osm->pmat);CHKERRQ(ierr);
38181a5c4aaSDmitry Karpeev   }
38281a5c4aaSDmitry Karpeev   else {
38381a5c4aaSDmitry Karpeev     ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->is, osm->is,scall,&osm->pmat);CHKERRQ(ierr);
38481a5c4aaSDmitry Karpeev   }
385b1a0a8a3SJed Brown   if (scall == MAT_INITIAL_MATRIX) {
386b1a0a8a3SJed Brown     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
387f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
388b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,osm->pmat[i]);CHKERRQ(ierr);
389b1a0a8a3SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
390b1a0a8a3SJed Brown     }
391b1a0a8a3SJed Brown   }
392b1a0a8a3SJed Brown 
393b1a0a8a3SJed Brown   /* Return control to the user so that the submatrices can be modified (e.g., to apply
394b1a0a8a3SJed Brown      different boundary conditions for the submatrices than for the global problem) */
395c10bc1a1SDmitry Karpeev   ierr = PCModifySubMatrices(pc,osm->n,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
396b1a0a8a3SJed Brown 
397b1a0a8a3SJed Brown   /*
398f746d493SDmitry Karpeev      Loop over submatrices putting them into local ksp
399b1a0a8a3SJed Brown   */
400f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
401b1a0a8a3SJed Brown     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr);
402b1a0a8a3SJed Brown     if (!pc->setupcalled) {
403b1a0a8a3SJed Brown       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
404b1a0a8a3SJed Brown     }
405b1a0a8a3SJed Brown   }
406b1a0a8a3SJed Brown 
407b1a0a8a3SJed Brown   PetscFunctionReturn(0);
408b1a0a8a3SJed Brown }
409b1a0a8a3SJed Brown 
410b1a0a8a3SJed Brown #undef __FUNCT__
411f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM"
412f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
413b1a0a8a3SJed Brown {
414f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
415b1a0a8a3SJed Brown   PetscErrorCode ierr;
416b1a0a8a3SJed Brown   PetscInt       i;
417b1a0a8a3SJed Brown 
418b1a0a8a3SJed Brown   PetscFunctionBegin;
419f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
420b1a0a8a3SJed Brown     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
421b1a0a8a3SJed Brown   }
422b1a0a8a3SJed Brown   PetscFunctionReturn(0);
423b1a0a8a3SJed Brown }
424b1a0a8a3SJed Brown 
425b1a0a8a3SJed Brown #undef __FUNCT__
426f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM"
427f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y)
428b1a0a8a3SJed Brown {
429f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
430b1a0a8a3SJed Brown   PetscErrorCode ierr;
431f746d493SDmitry Karpeev   PetscInt       i;
432b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
433b1a0a8a3SJed Brown 
434b1a0a8a3SJed Brown   PetscFunctionBegin;
435b1a0a8a3SJed Brown   /*
436b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
437b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
438b1a0a8a3SJed Brown   */
439f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
440b1a0a8a3SJed Brown     forward = SCATTER_FORWARD_LOCAL;
441b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
442f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx); CHKERRQ(ierr);
443b1a0a8a3SJed Brown   }
444f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
445b1a0a8a3SJed Brown     reverse = SCATTER_REVERSE_LOCAL;
446b1a0a8a3SJed Brown   }
447b1a0a8a3SJed Brown 
448f746d493SDmitry Karpeev   ierr = VecScatterBegin(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
449b1a0a8a3SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
450f746d493SDmitry Karpeev   ierr = VecScatterEnd(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
45186b96d36SDmitry Karpeev   /* do the subdomain solves */
45286b96d36SDmitry Karpeev   for (i=0; i<osm->n; ++i) {
453b1a0a8a3SJed Brown     ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
454b1a0a8a3SJed Brown   }
455f746d493SDmitry Karpeev   ierr = VecScatterBegin(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
456f746d493SDmitry Karpeev   ierr = VecScatterEnd(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
457b1a0a8a3SJed Brown   PetscFunctionReturn(0);
458b1a0a8a3SJed Brown }
459b1a0a8a3SJed Brown 
460b1a0a8a3SJed Brown #undef __FUNCT__
461f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM"
462f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y)
463b1a0a8a3SJed Brown {
464f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
465b1a0a8a3SJed Brown   PetscErrorCode ierr;
466f746d493SDmitry Karpeev   PetscInt       i;
467b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
468b1a0a8a3SJed Brown 
469b1a0a8a3SJed Brown   PetscFunctionBegin;
470b1a0a8a3SJed Brown   /*
471b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
472b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
473b1a0a8a3SJed Brown 
474f746d493SDmitry Karpeev      Note: these are reversed from the PCApply_GASM() because we are applying the
475b1a0a8a3SJed Brown      transpose of the three terms
476b1a0a8a3SJed Brown   */
477f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
478b1a0a8a3SJed Brown     forward = SCATTER_FORWARD_LOCAL;
479b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
480f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
481b1a0a8a3SJed Brown   }
482f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
483b1a0a8a3SJed Brown     reverse = SCATTER_REVERSE_LOCAL;
484b1a0a8a3SJed Brown   }
485b1a0a8a3SJed Brown 
486ab3e923fSDmitry Karpeev   ierr = VecScatterBegin(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
487b1a0a8a3SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
488ab3e923fSDmitry Karpeev   ierr = VecScatterEnd(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
489b1a0a8a3SJed Brown   /* do the local solves */
490ab3e923fSDmitry 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. */
491b1a0a8a3SJed Brown     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
492b1a0a8a3SJed Brown   }
493ab3e923fSDmitry Karpeev   ierr = VecScatterBegin(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
494ab3e923fSDmitry Karpeev   ierr = VecScatterEnd(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
495b1a0a8a3SJed Brown   PetscFunctionReturn(0);
496b1a0a8a3SJed Brown }
497b1a0a8a3SJed Brown 
498b1a0a8a3SJed Brown #undef __FUNCT__
499f746d493SDmitry Karpeev #define __FUNCT__ "PCDestroy_GASM"
500f746d493SDmitry Karpeev static PetscErrorCode PCDestroy_GASM(PC pc)
501b1a0a8a3SJed Brown {
502f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
503b1a0a8a3SJed Brown   PetscErrorCode ierr;
504b1a0a8a3SJed Brown   PetscInt       i;
505b1a0a8a3SJed Brown 
506b1a0a8a3SJed Brown   PetscFunctionBegin;
507b1a0a8a3SJed Brown   if (osm->ksp) {
508f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
509b1a0a8a3SJed Brown       ierr = KSPDestroy(osm->ksp[i]);CHKERRQ(ierr);
510b1a0a8a3SJed Brown     }
511b1a0a8a3SJed Brown     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
512b1a0a8a3SJed Brown   }
513b1a0a8a3SJed Brown   if (osm->pmat) {
514f746d493SDmitry Karpeev     if (osm->n > 0) {
515ab3e923fSDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
516b1a0a8a3SJed Brown     }
517b1a0a8a3SJed Brown   }
518f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
519b1a0a8a3SJed Brown     ierr = VecDestroy(osm->x[i]);CHKERRQ(ierr);
520b1a0a8a3SJed Brown     ierr = VecDestroy(osm->y[i]);CHKERRQ(ierr);
521b1a0a8a3SJed Brown   }
522b1a0a8a3SJed Brown   ierr = PetscFree(osm->x);CHKERRQ(ierr);
523b1a0a8a3SJed Brown   ierr = PetscFree(osm->y);CHKERRQ(ierr);
524f746d493SDmitry Karpeev   ierr = VecDestroy(osm->gx); CHKERRQ(ierr);
525ab3e923fSDmitry Karpeev   ierr = VecDestroy(osm->gy); CHKERRQ(ierr);
526ab3e923fSDmitry Karpeev 
527ab3e923fSDmitry Karpeev   ierr = VecScatterDestroy(osm->grestriction); CHKERRQ(ierr);
528ab3e923fSDmitry Karpeev   ierr = VecScatterDestroy(osm->gprolongation); CHKERRQ(ierr);
529f746d493SDmitry Karpeev 
530b1a0a8a3SJed Brown   if (osm->is) {
531ab3e923fSDmitry Karpeev     ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
532b1a0a8a3SJed Brown   }
533*167188bfSDmitry Karpeev 
534*167188bfSDmitry Karpeev   if(osm->gis) {
535*167188bfSDmitry Karpeev     ierr = ISDestroy(osm->gis); CHKERRQ(ierr);
536*167188bfSDmitry Karpeev   }
537*167188bfSDmitry Karpeev   if(osm->gis_local) {
538*167188bfSDmitry Karpeev     ierr = ISDestroy(osm->gis_local); CHKERRQ(ierr);
539*167188bfSDmitry Karpeev   }
540c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
541b1a0a8a3SJed Brown   PetscFunctionReturn(0);
542b1a0a8a3SJed Brown }
543b1a0a8a3SJed Brown 
544b1a0a8a3SJed Brown #undef __FUNCT__
545f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM"
546ab3e923fSDmitry Karpeev static PetscErrorCode PCSetFromOptions_GASM(PC pc) {
547f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
548b1a0a8a3SJed Brown   PetscErrorCode ierr;
549b1a0a8a3SJed Brown   PetscInt       blocks,ovl;
550b1a0a8a3SJed Brown   PetscBool      symset,flg;
551f746d493SDmitry Karpeev   PCGASMType      gasmtype;
552b1a0a8a3SJed Brown 
553b1a0a8a3SJed Brown   PetscFunctionBegin;
554b1a0a8a3SJed Brown   /* set the type to symmetric if matrix is symmetric */
555b1a0a8a3SJed Brown   if (!osm->type_set && pc->pmat) {
556b1a0a8a3SJed Brown     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
557f746d493SDmitry Karpeev     if (symset && flg) { osm->type = PC_GASM_BASIC; }
558b1a0a8a3SJed Brown   }
55944fe18e5SDmitry Karpeev   ierr = PetscOptionsHead("Generalized additive Schwarz options");CHKERRQ(ierr);
560f746d493SDmitry Karpeev     ierr = PetscOptionsInt("-pc_gasm_blocks","Number of subdomains","PCGASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
561ab3e923fSDmitry Karpeev     if (flg) {ierr = PCGASMSetTotalSubdomains(pc,blocks);CHKERRQ(ierr); }
562f746d493SDmitry Karpeev     ierr = PetscOptionsInt("-pc_gasm_overlap","Number of grid points overlap","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
563f746d493SDmitry Karpeev     if (flg) {ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); }
564b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
565f746d493SDmitry Karpeev     ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr);
566f746d493SDmitry Karpeev     if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr); }
567b1a0a8a3SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
568b1a0a8a3SJed Brown   PetscFunctionReturn(0);
569b1a0a8a3SJed Brown }
570b1a0a8a3SJed Brown 
571b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/
572b1a0a8a3SJed Brown 
573b1a0a8a3SJed Brown EXTERN_C_BEGIN
574b1a0a8a3SJed Brown #undef __FUNCT__
575f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetLocalSubdomains_GASM"
5767087cfbeSBarry Smith PetscErrorCode  PCGASMSetLocalSubdomains_GASM(PC pc,PetscInt n,IS is[],IS is_local[])
577b1a0a8a3SJed Brown {
578f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
579b1a0a8a3SJed Brown   PetscErrorCode ierr;
580b1a0a8a3SJed Brown   PetscInt       i;
581b1a0a8a3SJed Brown 
582b1a0a8a3SJed Brown   PetscFunctionBegin;
583b1a0a8a3SJed Brown   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
584f746d493SDmitry Karpeev   if (pc->setupcalled && (n != osm->n || is)) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetLocalSubdomains() should be called before calling PCSetUp().");
585b1a0a8a3SJed Brown 
586b1a0a8a3SJed Brown   if (!pc->setupcalled) {
58793c1ec32SDmitry Karpeev     osm->n            = n;
58893c1ec32SDmitry Karpeev     osm->is           = 0;
58993c1ec32SDmitry Karpeev     osm->is_local     = 0;
590b1a0a8a3SJed Brown     if (is) {
591b1a0a8a3SJed Brown       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
592b1a0a8a3SJed Brown     }
593b1a0a8a3SJed Brown     if (is_local) {
594b1a0a8a3SJed Brown       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
595b1a0a8a3SJed Brown     }
596b1a0a8a3SJed Brown     if (osm->is) {
597ab3e923fSDmitry Karpeev       ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
598b1a0a8a3SJed Brown     }
599b1a0a8a3SJed Brown     if (is) {
600b1a0a8a3SJed Brown       ierr = PetscMalloc(n*sizeof(IS),&osm->is);CHKERRQ(ierr);
601b1a0a8a3SJed Brown       for (i=0; i<n; i++) { osm->is[i] = is[i]; }
602f746d493SDmitry Karpeev       /* Flag indicating that the user has set overlapping subdomains so PCGASM should not increase their size. */
603b1a0a8a3SJed Brown       osm->overlap = -1;
604b1a0a8a3SJed Brown     }
605b1a0a8a3SJed Brown     if (is_local) {
606b1a0a8a3SJed Brown       ierr = PetscMalloc(n*sizeof(IS),&osm->is_local);CHKERRQ(ierr);
607b1a0a8a3SJed Brown       for (i=0; i<n; i++) { osm->is_local[i] = is_local[i]; }
608b1a0a8a3SJed Brown     }
609b1a0a8a3SJed Brown   }
610b1a0a8a3SJed Brown   PetscFunctionReturn(0);
611b1a0a8a3SJed Brown }
612b1a0a8a3SJed Brown EXTERN_C_END
613b1a0a8a3SJed Brown 
614b1a0a8a3SJed Brown EXTERN_C_BEGIN
615b1a0a8a3SJed Brown #undef __FUNCT__
616f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains_GASM"
6177087cfbeSBarry Smith PetscErrorCode  PCGASMSetTotalSubdomains_GASM(PC pc,PetscInt N) {
618f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
619b1a0a8a3SJed Brown   PetscErrorCode ierr;
620b1a0a8a3SJed Brown   PetscMPIInt    rank,size;
621b1a0a8a3SJed Brown   PetscInt       n;
622b1a0a8a3SJed Brown 
623b1a0a8a3SJed Brown   PetscFunctionBegin;
624b1a0a8a3SJed Brown   if (N < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
625b1a0a8a3SJed Brown 
626b1a0a8a3SJed Brown   /*
627b1a0a8a3SJed Brown      Split the subdomains equally among all processors
628b1a0a8a3SJed Brown   */
629b1a0a8a3SJed Brown   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
630b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
631b1a0a8a3SJed Brown   n = N/size + ((N % size) > rank);
632b1a0a8a3SJed 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);
633f746d493SDmitry Karpeev   if (pc->setupcalled && n != osm->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before PCSetUp().");
634b1a0a8a3SJed Brown   if (!pc->setupcalled) {
635b1a0a8a3SJed Brown     if (osm->is) {
636ab3e923fSDmitry Karpeev       ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
637b1a0a8a3SJed Brown     }
63893c1ec32SDmitry Karpeev     osm->N            = N;
639f746d493SDmitry Karpeev     osm->n            = n;
640b1a0a8a3SJed Brown     osm->is           = 0;
641b1a0a8a3SJed Brown     osm->is_local     = 0;
642b1a0a8a3SJed Brown   }
643b1a0a8a3SJed Brown   PetscFunctionReturn(0);
644ab3e923fSDmitry Karpeev }/* PCGASMSetTotalSubdomains_GASM() */
645b1a0a8a3SJed Brown EXTERN_C_END
646b1a0a8a3SJed Brown 
647b1a0a8a3SJed Brown EXTERN_C_BEGIN
648b1a0a8a3SJed Brown #undef __FUNCT__
649f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM"
6507087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
651b1a0a8a3SJed Brown {
652f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
653b1a0a8a3SJed Brown 
654b1a0a8a3SJed Brown   PetscFunctionBegin;
655b1a0a8a3SJed Brown   if (ovl < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
656f746d493SDmitry Karpeev   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
657b1a0a8a3SJed Brown   if (!pc->setupcalled) {
658b1a0a8a3SJed Brown     osm->overlap = ovl;
659b1a0a8a3SJed Brown   }
660b1a0a8a3SJed Brown   PetscFunctionReturn(0);
661b1a0a8a3SJed Brown }
662b1a0a8a3SJed Brown EXTERN_C_END
663b1a0a8a3SJed Brown 
664b1a0a8a3SJed Brown EXTERN_C_BEGIN
665b1a0a8a3SJed Brown #undef __FUNCT__
666f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM"
6677087cfbeSBarry Smith PetscErrorCode  PCGASMSetType_GASM(PC pc,PCGASMType type)
668b1a0a8a3SJed Brown {
669f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
670b1a0a8a3SJed Brown 
671b1a0a8a3SJed Brown   PetscFunctionBegin;
672b1a0a8a3SJed Brown   osm->type     = type;
673b1a0a8a3SJed Brown   osm->type_set = PETSC_TRUE;
674b1a0a8a3SJed Brown   PetscFunctionReturn(0);
675b1a0a8a3SJed Brown }
676b1a0a8a3SJed Brown EXTERN_C_END
677b1a0a8a3SJed Brown 
678b1a0a8a3SJed Brown EXTERN_C_BEGIN
679b1a0a8a3SJed Brown #undef __FUNCT__
680f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM"
6817087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices_GASM(PC pc,PetscBool  doSort)
682b1a0a8a3SJed Brown {
683f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
684b1a0a8a3SJed Brown 
685b1a0a8a3SJed Brown   PetscFunctionBegin;
686b1a0a8a3SJed Brown   osm->sort_indices = doSort;
687b1a0a8a3SJed Brown   PetscFunctionReturn(0);
688b1a0a8a3SJed Brown }
689b1a0a8a3SJed Brown EXTERN_C_END
690b1a0a8a3SJed Brown 
691b1a0a8a3SJed Brown EXTERN_C_BEGIN
692b1a0a8a3SJed Brown #undef __FUNCT__
693f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM"
69444fe18e5SDmitry Karpeev /*
69544fe18e5SDmitry Karpeev    FIX: This routine might need to be modified once multiple ranks per subdomain are allowed.
69644fe18e5SDmitry Karpeev         In particular, it would upset the global subdomain number calculation.
69744fe18e5SDmitry Karpeev */
6987087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
699b1a0a8a3SJed Brown {
700f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
701b1a0a8a3SJed Brown   PetscErrorCode ierr;
702b1a0a8a3SJed Brown 
703b1a0a8a3SJed Brown   PetscFunctionBegin;
704ab3e923fSDmitry 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");
705b1a0a8a3SJed Brown 
706ab3e923fSDmitry Karpeev   if (n) {
707ab3e923fSDmitry Karpeev     *n = osm->n;
708b1a0a8a3SJed Brown   }
709ab3e923fSDmitry Karpeev   if (first) {
710ab3e923fSDmitry Karpeev     ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr);
711ab3e923fSDmitry Karpeev     *first -= osm->n;
712b1a0a8a3SJed Brown   }
713b1a0a8a3SJed Brown   if (ksp) {
714b1a0a8a3SJed Brown     /* Assume that local solves are now different; not necessarily
715f746d493SDmitry Karpeev        true though!  This flag is used only for PCView_GASM() */
716b1a0a8a3SJed Brown     *ksp                   = osm->ksp;
717b1a0a8a3SJed Brown     osm->same_local_solves = PETSC_FALSE;
718b1a0a8a3SJed Brown   }
719b1a0a8a3SJed Brown   PetscFunctionReturn(0);
720ab3e923fSDmitry Karpeev }/* PCGASMGetSubKSP_GASM() */
721b1a0a8a3SJed Brown EXTERN_C_END
722b1a0a8a3SJed Brown 
723b1a0a8a3SJed Brown 
724b1a0a8a3SJed Brown #undef __FUNCT__
725f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetLocalSubdomains"
726b1a0a8a3SJed Brown /*@C
727f746d493SDmitry Karpeev     PCGASMSetLocalSubdomains - Sets the local subdomains (for this processor
728b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
729b1a0a8a3SJed Brown 
730b1a0a8a3SJed Brown     Collective on PC
731b1a0a8a3SJed Brown 
732b1a0a8a3SJed Brown     Input Parameters:
733b1a0a8a3SJed Brown +   pc - the preconditioner context
734b1a0a8a3SJed Brown .   n - the number of subdomains for this processor (default value = 1)
735b1a0a8a3SJed Brown .   is - the index set that defines the subdomains for this processor
736b1a0a8a3SJed Brown          (or PETSC_NULL for PETSc to determine subdomains)
737b1a0a8a3SJed Brown -   is_local - the index sets that define the local part of the subdomains for this processor
738b1a0a8a3SJed Brown          (or PETSC_NULL to use the default of 1 subdomain per process)
739b1a0a8a3SJed Brown 
740b1a0a8a3SJed Brown     Notes:
741b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
742b1a0a8a3SJed Brown 
743f746d493SDmitry Karpeev     By default the GASM preconditioner uses 1 block per processor.
744b1a0a8a3SJed Brown 
745f746d493SDmitry Karpeev     Use PCGASMSetTotalSubdomains() to set the subdomains for all processors.
746b1a0a8a3SJed Brown 
747b1a0a8a3SJed Brown     Level: advanced
748b1a0a8a3SJed Brown 
749f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz
750b1a0a8a3SJed Brown 
751f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
752f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetLocalSubdomains()
753b1a0a8a3SJed Brown @*/
7547087cfbeSBarry Smith PetscErrorCode  PCGASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
755b1a0a8a3SJed Brown {
756b1a0a8a3SJed Brown   PetscErrorCode ierr;
757b1a0a8a3SJed Brown 
758b1a0a8a3SJed Brown   PetscFunctionBegin;
759b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
760f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr);
761b1a0a8a3SJed Brown   PetscFunctionReturn(0);
762b1a0a8a3SJed Brown }
763b1a0a8a3SJed Brown 
764b1a0a8a3SJed Brown #undef __FUNCT__
765f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains"
766b1a0a8a3SJed Brown /*@C
767f746d493SDmitry Karpeev     PCGASMSetTotalSubdomains - Sets the subdomains for all processor for the
768b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
769b1a0a8a3SJed Brown     PC communicator must call this routine, with the same index sets.
770b1a0a8a3SJed Brown 
771b1a0a8a3SJed Brown     Collective on PC
772b1a0a8a3SJed Brown 
773b1a0a8a3SJed Brown     Input Parameters:
774b1a0a8a3SJed Brown +   pc - the preconditioner context
775b1a0a8a3SJed Brown .   n - the number of subdomains for all processors
776b1a0a8a3SJed Brown .   is - the index sets that define the subdomains for all processor
777b1a0a8a3SJed Brown          (or PETSC_NULL for PETSc to determine subdomains)
778b1a0a8a3SJed Brown -   is_local - the index sets that define the local part of the subdomains for this processor
779b1a0a8a3SJed Brown          (or PETSC_NULL to use the default of 1 subdomain per process)
780b1a0a8a3SJed Brown 
781b1a0a8a3SJed Brown     Options Database Key:
782b1a0a8a3SJed Brown     To set the total number of subdomain blocks rather than specify the
783b1a0a8a3SJed Brown     index sets, use the option
784f746d493SDmitry Karpeev .    -pc_gasm_blocks <blks> - Sets total blocks
785b1a0a8a3SJed Brown 
786b1a0a8a3SJed Brown     Notes:
787b1a0a8a3SJed Brown     Currently you cannot use this to set the actual subdomains with the argument is.
788b1a0a8a3SJed Brown 
789f746d493SDmitry Karpeev     By default the GASM preconditioner uses 1 block per processor.
790b1a0a8a3SJed Brown 
791b1a0a8a3SJed Brown     These index sets cannot be destroyed until after completion of the
792f746d493SDmitry Karpeev     linear solves for which the GASM preconditioner is being used.
793b1a0a8a3SJed Brown 
794f746d493SDmitry Karpeev     Use PCGASMSetLocalSubdomains() to set local subdomains.
795b1a0a8a3SJed Brown 
796b1a0a8a3SJed Brown     Level: advanced
797b1a0a8a3SJed Brown 
798f746d493SDmitry Karpeev .keywords: PC, GASM, set, total, global, subdomains, additive Schwarz
799b1a0a8a3SJed Brown 
800f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
801f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
802b1a0a8a3SJed Brown @*/
8037087cfbeSBarry Smith PetscErrorCode  PCGASMSetTotalSubdomains(PC pc,PetscInt N)
804b1a0a8a3SJed Brown {
805b1a0a8a3SJed Brown   PetscErrorCode ierr;
806b1a0a8a3SJed Brown 
807b1a0a8a3SJed Brown   PetscFunctionBegin;
808b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
809f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetTotalSubdomains_C",(PC,PetscInt),(pc,N));CHKERRQ(ierr);
810b1a0a8a3SJed Brown   PetscFunctionReturn(0);
811b1a0a8a3SJed Brown }
812b1a0a8a3SJed Brown 
813b1a0a8a3SJed Brown #undef __FUNCT__
814f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap"
815b1a0a8a3SJed Brown /*@
816f746d493SDmitry Karpeev     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
817b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
818b1a0a8a3SJed Brown     PC communicator must call this routine.
819b1a0a8a3SJed Brown 
820b1a0a8a3SJed Brown     Logically Collective on PC
821b1a0a8a3SJed Brown 
822b1a0a8a3SJed Brown     Input Parameters:
823b1a0a8a3SJed Brown +   pc  - the preconditioner context
824b1a0a8a3SJed Brown -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
825b1a0a8a3SJed Brown 
826b1a0a8a3SJed Brown     Options Database Key:
827f746d493SDmitry Karpeev .   -pc_gasm_overlap <ovl> - Sets overlap
828b1a0a8a3SJed Brown 
829b1a0a8a3SJed Brown     Notes:
830f746d493SDmitry Karpeev     By default the GASM preconditioner uses 1 block per processor.  To use
831f746d493SDmitry Karpeev     multiple blocks per perocessor, see PCGASMSetTotalSubdomains() and
832f746d493SDmitry Karpeev     PCGASMSetLocalSubdomains() (and the option -pc_gasm_blocks <blks>).
833b1a0a8a3SJed Brown 
834b1a0a8a3SJed Brown     The overlap defaults to 1, so if one desires that no additional
835b1a0a8a3SJed Brown     overlap be computed beyond what may have been set with a call to
836f746d493SDmitry Karpeev     PCGASMSetTotalSubdomains() or PCGASMSetLocalSubdomains(), then ovl
837b1a0a8a3SJed Brown     must be set to be 0.  In particular, if one does not explicitly set
838b1a0a8a3SJed Brown     the subdomains an application code, then all overlap would be computed
839f746d493SDmitry Karpeev     internally by PETSc, and using an overlap of 0 would result in an GASM
840b1a0a8a3SJed Brown     variant that is equivalent to the block Jacobi preconditioner.
841b1a0a8a3SJed Brown 
842b1a0a8a3SJed Brown     Note that one can define initial index sets with any overlap via
843f746d493SDmitry Karpeev     PCGASMSetTotalSubdomains() or PCGASMSetLocalSubdomains(); the routine
844f746d493SDmitry Karpeev     PCGASMSetOverlap() merely allows PETSc to extend that overlap further
845b1a0a8a3SJed Brown     if desired.
846b1a0a8a3SJed Brown 
847b1a0a8a3SJed Brown     Level: intermediate
848b1a0a8a3SJed Brown 
849f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap
850b1a0a8a3SJed Brown 
851f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetLocalSubdomains(), PCGASMGetSubKSP(),
852f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetLocalSubdomains()
853b1a0a8a3SJed Brown @*/
8547087cfbeSBarry Smith PetscErrorCode  PCGASMSetOverlap(PC pc,PetscInt ovl)
855b1a0a8a3SJed Brown {
856b1a0a8a3SJed Brown   PetscErrorCode ierr;
857b1a0a8a3SJed Brown 
858b1a0a8a3SJed Brown   PetscFunctionBegin;
859b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
860b1a0a8a3SJed Brown   PetscValidLogicalCollectiveInt(pc,ovl,2);
861f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
862b1a0a8a3SJed Brown   PetscFunctionReturn(0);
863b1a0a8a3SJed Brown }
864b1a0a8a3SJed Brown 
865b1a0a8a3SJed Brown #undef __FUNCT__
866f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType"
867b1a0a8a3SJed Brown /*@
868f746d493SDmitry Karpeev     PCGASMSetType - Sets the type of restriction and interpolation used
869b1a0a8a3SJed Brown     for local problems in the additive Schwarz method.
870b1a0a8a3SJed Brown 
871b1a0a8a3SJed Brown     Logically Collective on PC
872b1a0a8a3SJed Brown 
873b1a0a8a3SJed Brown     Input Parameters:
874b1a0a8a3SJed Brown +   pc  - the preconditioner context
875f746d493SDmitry Karpeev -   type - variant of GASM, one of
876b1a0a8a3SJed Brown .vb
877f746d493SDmitry Karpeev       PC_GASM_BASIC       - full interpolation and restriction
878f746d493SDmitry Karpeev       PC_GASM_RESTRICT    - full restriction, local processor interpolation
879f746d493SDmitry Karpeev       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
880f746d493SDmitry Karpeev       PC_GASM_NONE        - local processor restriction and interpolation
881b1a0a8a3SJed Brown .ve
882b1a0a8a3SJed Brown 
883b1a0a8a3SJed Brown     Options Database Key:
884f746d493SDmitry Karpeev .   -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
885b1a0a8a3SJed Brown 
886b1a0a8a3SJed Brown     Level: intermediate
887b1a0a8a3SJed Brown 
888f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
889b1a0a8a3SJed Brown 
890f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(),
891f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
892b1a0a8a3SJed Brown @*/
8937087cfbeSBarry Smith PetscErrorCode  PCGASMSetType(PC pc,PCGASMType type)
894b1a0a8a3SJed Brown {
895b1a0a8a3SJed Brown   PetscErrorCode ierr;
896b1a0a8a3SJed Brown 
897b1a0a8a3SJed Brown   PetscFunctionBegin;
898b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
899b1a0a8a3SJed Brown   PetscValidLogicalCollectiveEnum(pc,type,2);
900f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr);
901b1a0a8a3SJed Brown   PetscFunctionReturn(0);
902b1a0a8a3SJed Brown }
903b1a0a8a3SJed Brown 
904b1a0a8a3SJed Brown #undef __FUNCT__
905f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices"
906b1a0a8a3SJed Brown /*@
907f746d493SDmitry Karpeev     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
908b1a0a8a3SJed Brown 
909b1a0a8a3SJed Brown     Logically Collective on PC
910b1a0a8a3SJed Brown 
911b1a0a8a3SJed Brown     Input Parameters:
912b1a0a8a3SJed Brown +   pc  - the preconditioner context
913b1a0a8a3SJed Brown -   doSort - sort the subdomain indices
914b1a0a8a3SJed Brown 
915b1a0a8a3SJed Brown     Level: intermediate
916b1a0a8a3SJed Brown 
917f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
918b1a0a8a3SJed Brown 
919f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(),
920f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
921b1a0a8a3SJed Brown @*/
9227087cfbeSBarry Smith PetscErrorCode  PCGASMSetSortIndices(PC pc,PetscBool  doSort)
923b1a0a8a3SJed Brown {
924b1a0a8a3SJed Brown   PetscErrorCode ierr;
925b1a0a8a3SJed Brown 
926b1a0a8a3SJed Brown   PetscFunctionBegin;
927b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
928b1a0a8a3SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
929f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
930b1a0a8a3SJed Brown   PetscFunctionReturn(0);
931b1a0a8a3SJed Brown }
932b1a0a8a3SJed Brown 
933b1a0a8a3SJed Brown #undef __FUNCT__
934f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP"
935b1a0a8a3SJed Brown /*@C
936f746d493SDmitry Karpeev    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
937b1a0a8a3SJed Brown    this processor.
938b1a0a8a3SJed Brown 
939b1a0a8a3SJed Brown    Collective on PC iff first_local is requested
940b1a0a8a3SJed Brown 
941b1a0a8a3SJed Brown    Input Parameter:
942b1a0a8a3SJed Brown .  pc - the preconditioner context
943b1a0a8a3SJed Brown 
944b1a0a8a3SJed Brown    Output Parameters:
945b1a0a8a3SJed Brown +  n_local - the number of blocks on this processor or PETSC_NULL
946b1a0a8a3SJed Brown .  first_local - the global number of the first block on this processor or PETSC_NULL,
947b1a0a8a3SJed Brown                  all processors must request or all must pass PETSC_NULL
948b1a0a8a3SJed Brown -  ksp - the array of KSP contexts
949b1a0a8a3SJed Brown 
950b1a0a8a3SJed Brown    Note:
951f746d493SDmitry Karpeev    After PCGASMGetSubKSP() the array of KSPes is not to be freed
952b1a0a8a3SJed Brown 
953b1a0a8a3SJed Brown    Currently for some matrix implementations only 1 block per processor
954b1a0a8a3SJed Brown    is supported.
955b1a0a8a3SJed Brown 
956f746d493SDmitry Karpeev    You must call KSPSetUp() before calling PCGASMGetSubKSP().
957b1a0a8a3SJed Brown 
958b1a0a8a3SJed Brown    Level: advanced
959b1a0a8a3SJed Brown 
960f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context
961b1a0a8a3SJed Brown 
962f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMSetOverlap(),
963f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(),
964b1a0a8a3SJed Brown @*/
9657087cfbeSBarry Smith PetscErrorCode  PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
966b1a0a8a3SJed Brown {
967b1a0a8a3SJed Brown   PetscErrorCode ierr;
968b1a0a8a3SJed Brown 
969b1a0a8a3SJed Brown   PetscFunctionBegin;
970b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
971f746d493SDmitry Karpeev   ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
972b1a0a8a3SJed Brown   PetscFunctionReturn(0);
973b1a0a8a3SJed Brown }
974b1a0a8a3SJed Brown 
975b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/
976b1a0a8a3SJed Brown /*MC
977f746d493SDmitry Karpeev    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
978b1a0a8a3SJed Brown            its own KSP object.
979b1a0a8a3SJed Brown 
980b1a0a8a3SJed Brown    Options Database Keys:
981f746d493SDmitry Karpeev +  -pc_gasm_truelocal - Activates PCGGASMSetUseTrueLocal()
982f746d493SDmitry Karpeev .  -pc_gasm_blocks <blks> - Sets total blocks
983f746d493SDmitry Karpeev .  -pc_gasm_overlap <ovl> - Sets overlap
984f746d493SDmitry Karpeev -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
985b1a0a8a3SJed Brown 
986b1a0a8a3SJed Brown      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
987f746d493SDmitry Karpeev       will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
988f746d493SDmitry Karpeev       -pc_gasm_type basic to use the standard GASM.
989b1a0a8a3SJed Brown 
990b1a0a8a3SJed Brown    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
991b1a0a8a3SJed Brown      than one processor. Defaults to one block per processor.
992b1a0a8a3SJed Brown 
993b1a0a8a3SJed Brown      To set options on the solvers for each block append -sub_ to all the KSP, and PC
994b1a0a8a3SJed Brown         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
995b1a0a8a3SJed Brown 
996f746d493SDmitry Karpeev      To set the options on the solvers separate for each block call PCGASMGetSubKSP()
997b1a0a8a3SJed Brown          and set the options directly on the resulting KSP object (you can access its PC
998b1a0a8a3SJed Brown          with KSPGetPC())
999b1a0a8a3SJed Brown 
1000b1a0a8a3SJed Brown 
1001b1a0a8a3SJed Brown    Level: beginner
1002b1a0a8a3SJed Brown 
1003b1a0a8a3SJed Brown    Concepts: additive Schwarz method
1004b1a0a8a3SJed Brown 
1005b1a0a8a3SJed Brown     References:
1006b1a0a8a3SJed Brown     An additive variant of the Schwarz alternating method for the case of many subregions
1007b1a0a8a3SJed Brown     M Dryja, OB Widlund - Courant Institute, New York University Technical report
1008b1a0a8a3SJed Brown 
1009b1a0a8a3SJed Brown     Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1010b1a0a8a3SJed Brown     Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X.
1011b1a0a8a3SJed Brown 
1012b1a0a8a3SJed Brown .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1013f746d493SDmitry Karpeev            PCBJACOBI, PCGASMSetUseTrueLocal(), PCGASMGetSubKSP(), PCGASMSetLocalSubdomains(),
1014f746d493SDmitry Karpeev            PCGASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType()
1015b1a0a8a3SJed Brown 
1016b1a0a8a3SJed Brown M*/
1017b1a0a8a3SJed Brown 
1018b1a0a8a3SJed Brown EXTERN_C_BEGIN
1019b1a0a8a3SJed Brown #undef __FUNCT__
1020f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM"
10217087cfbeSBarry Smith PetscErrorCode  PCCreate_GASM(PC pc)
1022b1a0a8a3SJed Brown {
1023b1a0a8a3SJed Brown   PetscErrorCode ierr;
1024f746d493SDmitry Karpeev   PC_GASM         *osm;
1025b1a0a8a3SJed Brown 
1026b1a0a8a3SJed Brown   PetscFunctionBegin;
1027f746d493SDmitry Karpeev   ierr = PetscNewLog(pc,PC_GASM,&osm);CHKERRQ(ierr);
1028ab3e923fSDmitry Karpeev   osm->N                 = PETSC_DECIDE;
1029ab3e923fSDmitry Karpeev   osm->n                 = 0;
1030ab3e923fSDmitry Karpeev   osm->nmax              = 0;
1031b1a0a8a3SJed Brown   osm->overlap           = 1;
1032b1a0a8a3SJed Brown   osm->ksp               = 0;
1033ab3e923fSDmitry Karpeev   osm->grestriction      = 0;
1034ab3e923fSDmitry Karpeev   osm->gprolongation     = 0;
1035ab3e923fSDmitry Karpeev   osm->gx                = 0;
1036ab3e923fSDmitry Karpeev   osm->gy                = 0;
1037b1a0a8a3SJed Brown   osm->x                 = 0;
1038b1a0a8a3SJed Brown   osm->y                 = 0;
1039b1a0a8a3SJed Brown   osm->is                = 0;
1040b1a0a8a3SJed Brown   osm->is_local          = 0;
1041b1a0a8a3SJed Brown   osm->pmat              = 0;
1042f746d493SDmitry Karpeev   osm->type              = PC_GASM_RESTRICT;
1043b1a0a8a3SJed Brown   osm->same_local_solves = PETSC_TRUE;
1044b1a0a8a3SJed Brown   osm->sort_indices      = PETSC_TRUE;
1045b1a0a8a3SJed Brown 
1046b1a0a8a3SJed Brown   pc->data                   = (void*)osm;
1047f746d493SDmitry Karpeev   pc->ops->apply             = PCApply_GASM;
1048f746d493SDmitry Karpeev   pc->ops->applytranspose    = PCApplyTranspose_GASM;
1049f746d493SDmitry Karpeev   pc->ops->setup             = PCSetUp_GASM;
1050f746d493SDmitry Karpeev   pc->ops->destroy           = PCDestroy_GASM;
1051f746d493SDmitry Karpeev   pc->ops->setfromoptions    = PCSetFromOptions_GASM;
1052f746d493SDmitry Karpeev   pc->ops->setuponblocks     = PCSetUpOnBlocks_GASM;
1053f746d493SDmitry Karpeev   pc->ops->view              = PCView_GASM;
1054b1a0a8a3SJed Brown   pc->ops->applyrichardson   = 0;
1055b1a0a8a3SJed Brown 
1056f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetLocalSubdomains_C","PCGASMSetLocalSubdomains_GASM",
1057f746d493SDmitry Karpeev                     PCGASMSetLocalSubdomains_GASM);CHKERRQ(ierr);
1058f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetTotalSubdomains_C","PCGASMSetTotalSubdomains_GASM",
1059f746d493SDmitry Karpeev                     PCGASMSetTotalSubdomains_GASM);CHKERRQ(ierr);
1060f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetOverlap_C","PCGASMSetOverlap_GASM",
1061f746d493SDmitry Karpeev                     PCGASMSetOverlap_GASM);CHKERRQ(ierr);
1062f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetType_C","PCGASMSetType_GASM",
1063f746d493SDmitry Karpeev                     PCGASMSetType_GASM);CHKERRQ(ierr);
1064f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSortIndices_C","PCGASMSetSortIndices_GASM",
1065f746d493SDmitry Karpeev                     PCGASMSetSortIndices_GASM);CHKERRQ(ierr);
1066f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMGetSubKSP_C","PCGASMGetSubKSP_GASM",
1067f746d493SDmitry Karpeev                     PCGASMGetSubKSP_GASM);CHKERRQ(ierr);
1068b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1069b1a0a8a3SJed Brown }
1070b1a0a8a3SJed Brown EXTERN_C_END
1071b1a0a8a3SJed Brown 
1072b1a0a8a3SJed Brown 
1073b1a0a8a3SJed Brown #undef __FUNCT__
1074f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains"
1075b1a0a8a3SJed Brown /*@C
1076f746d493SDmitry Karpeev    PCGASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
1077b1a0a8a3SJed Brown    preconditioner for a any problem on a general grid.
1078b1a0a8a3SJed Brown 
1079b1a0a8a3SJed Brown    Collective
1080b1a0a8a3SJed Brown 
1081b1a0a8a3SJed Brown    Input Parameters:
1082b1a0a8a3SJed Brown +  A - The global matrix operator
1083b1a0a8a3SJed Brown -  n - the number of local blocks
1084b1a0a8a3SJed Brown 
1085b1a0a8a3SJed Brown    Output Parameters:
1086b1a0a8a3SJed Brown .  outis - the array of index sets defining the subdomains
1087b1a0a8a3SJed Brown 
1088b1a0a8a3SJed Brown    Level: advanced
1089b1a0a8a3SJed Brown 
1090f746d493SDmitry Karpeev    Note: this generates nonoverlapping subdomains; the PCGASM will generate the overlap
1091f746d493SDmitry Karpeev     from these if you use PCGASMSetLocalSubdomains()
1092b1a0a8a3SJed Brown 
1093b1a0a8a3SJed Brown     In the Fortran version you must provide the array outis[] already allocated of length n.
1094b1a0a8a3SJed Brown 
1095f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1096b1a0a8a3SJed Brown 
1097f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMDestroySubdomains()
1098b1a0a8a3SJed Brown @*/
10997087cfbeSBarry Smith PetscErrorCode  PCGASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
1100b1a0a8a3SJed Brown {
1101b1a0a8a3SJed Brown   MatPartitioning           mpart;
1102b1a0a8a3SJed Brown   const char                *prefix;
1103b1a0a8a3SJed Brown   PetscErrorCode            (*f)(Mat,PetscBool *,MatReuse,Mat*);
1104b1a0a8a3SJed Brown   PetscMPIInt               size;
1105b1a0a8a3SJed Brown   PetscInt                  i,j,rstart,rend,bs;
1106b1a0a8a3SJed Brown   PetscBool                 iscopy = PETSC_FALSE,isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
1107b1a0a8a3SJed Brown   Mat                       Ad = PETSC_NULL, adj;
1108b1a0a8a3SJed Brown   IS                        ispart,isnumb,*is;
1109b1a0a8a3SJed Brown   PetscErrorCode            ierr;
1110b1a0a8a3SJed Brown 
1111b1a0a8a3SJed Brown   PetscFunctionBegin;
1112b1a0a8a3SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1113b1a0a8a3SJed Brown   PetscValidPointer(outis,3);
1114b1a0a8a3SJed Brown   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
1115b1a0a8a3SJed Brown 
1116b1a0a8a3SJed Brown   /* Get prefix, row distribution, and block size */
1117b1a0a8a3SJed Brown   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1118b1a0a8a3SJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1119b1a0a8a3SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1120b1a0a8a3SJed 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);
1121b1a0a8a3SJed Brown 
1122b1a0a8a3SJed Brown   /* Get diagonal block from matrix if possible */
1123b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1124b1a0a8a3SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr);
1125b1a0a8a3SJed Brown   if (f) {
1126b1a0a8a3SJed Brown     ierr = (*f)(A,&iscopy,MAT_INITIAL_MATRIX,&Ad);CHKERRQ(ierr);
1127b1a0a8a3SJed Brown   } else if (size == 1) {
1128b1a0a8a3SJed Brown     iscopy = PETSC_FALSE; Ad = A;
1129b1a0a8a3SJed Brown   } else {
1130b1a0a8a3SJed Brown     iscopy = PETSC_FALSE; Ad = PETSC_NULL;
1131b1a0a8a3SJed Brown   }
1132b1a0a8a3SJed Brown   if (Ad) {
1133b1a0a8a3SJed Brown     ierr = PetscTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1134b1a0a8a3SJed Brown     if (!isbaij) {ierr = PetscTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
1135b1a0a8a3SJed Brown   }
1136b1a0a8a3SJed Brown   if (Ad && n > 1) {
1137b1a0a8a3SJed Brown     PetscBool  match,done;
1138b1a0a8a3SJed Brown     /* Try to setup a good matrix partitioning if available */
1139b1a0a8a3SJed Brown     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
1140b1a0a8a3SJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
1141b1a0a8a3SJed Brown     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1142b1a0a8a3SJed Brown     ierr = PetscTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
1143b1a0a8a3SJed Brown     if (!match) {
1144b1a0a8a3SJed Brown       ierr = PetscTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
1145b1a0a8a3SJed Brown     }
1146b1a0a8a3SJed Brown     if (!match) { /* assume a "good" partitioner is available */
1147b1a0a8a3SJed Brown       PetscInt na,*ia,*ja;
1148b1a0a8a3SJed Brown       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1149b1a0a8a3SJed Brown       if (done) {
1150b1a0a8a3SJed Brown 	/* Build adjacency matrix by hand. Unfortunately a call to
1151b1a0a8a3SJed Brown 	   MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1152b1a0a8a3SJed Brown 	   remove the block-aij structure and we cannot expect
1153b1a0a8a3SJed Brown 	   MatPartitioning to split vertices as we need */
1154b1a0a8a3SJed Brown 	PetscInt i,j,*row,len,nnz,cnt,*iia=0,*jja=0;
1155b1a0a8a3SJed Brown 	nnz = 0;
1156b1a0a8a3SJed Brown 	for (i=0; i<na; i++) { /* count number of nonzeros */
1157b1a0a8a3SJed Brown 	  len = ia[i+1] - ia[i];
1158b1a0a8a3SJed Brown 	  row = ja + ia[i];
1159b1a0a8a3SJed Brown 	  for (j=0; j<len; j++) {
1160b1a0a8a3SJed Brown 	    if (row[j] == i) { /* don't count diagonal */
1161b1a0a8a3SJed Brown 	      len--; break;
1162b1a0a8a3SJed Brown 	    }
1163b1a0a8a3SJed Brown 	  }
1164b1a0a8a3SJed Brown 	  nnz += len;
1165b1a0a8a3SJed Brown 	}
1166b1a0a8a3SJed Brown 	ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr);
1167b1a0a8a3SJed Brown 	ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr);
1168b1a0a8a3SJed Brown 	nnz    = 0;
1169b1a0a8a3SJed Brown 	iia[0] = 0;
1170b1a0a8a3SJed Brown 	for (i=0; i<na; i++) { /* fill adjacency */
1171b1a0a8a3SJed Brown 	  cnt = 0;
1172b1a0a8a3SJed Brown 	  len = ia[i+1] - ia[i];
1173b1a0a8a3SJed Brown 	  row = ja + ia[i];
1174b1a0a8a3SJed Brown 	  for (j=0; j<len; j++) {
1175b1a0a8a3SJed Brown 	    if (row[j] != i) { /* if not diagonal */
1176b1a0a8a3SJed Brown 	      jja[nnz+cnt++] = row[j];
1177b1a0a8a3SJed Brown 	    }
1178b1a0a8a3SJed Brown 	  }
1179b1a0a8a3SJed Brown 	  nnz += cnt;
1180b1a0a8a3SJed Brown 	  iia[i+1] = nnz;
1181b1a0a8a3SJed Brown 	}
1182b1a0a8a3SJed Brown 	/* Partitioning of the adjacency matrix */
1183b1a0a8a3SJed Brown 	ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);CHKERRQ(ierr);
1184b1a0a8a3SJed Brown 	ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
1185b1a0a8a3SJed Brown 	ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
1186b1a0a8a3SJed Brown 	ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
1187b1a0a8a3SJed Brown 	ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
1188b1a0a8a3SJed Brown 	ierr = MatDestroy(adj);CHKERRQ(ierr);
1189b1a0a8a3SJed Brown 	foundpart = PETSC_TRUE;
1190b1a0a8a3SJed Brown       }
1191b1a0a8a3SJed Brown       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1192b1a0a8a3SJed Brown     }
1193b1a0a8a3SJed Brown     ierr = MatPartitioningDestroy(mpart);CHKERRQ(ierr);
1194b1a0a8a3SJed Brown   }
1195b1a0a8a3SJed Brown   if (iscopy) {ierr = MatDestroy(Ad);CHKERRQ(ierr);}
1196b1a0a8a3SJed Brown 
1197b1a0a8a3SJed Brown   ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr);
1198b1a0a8a3SJed Brown   *outis = is;
1199b1a0a8a3SJed Brown 
1200b1a0a8a3SJed Brown   if (!foundpart) {
1201b1a0a8a3SJed Brown 
1202b1a0a8a3SJed Brown     /* Partitioning by contiguous chunks of rows */
1203b1a0a8a3SJed Brown 
1204b1a0a8a3SJed Brown     PetscInt mbs   = (rend-rstart)/bs;
1205b1a0a8a3SJed Brown     PetscInt start = rstart;
1206b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1207b1a0a8a3SJed Brown       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
1208b1a0a8a3SJed Brown       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
1209b1a0a8a3SJed Brown       start += count;
1210b1a0a8a3SJed Brown     }
1211b1a0a8a3SJed Brown 
1212b1a0a8a3SJed Brown   } else {
1213b1a0a8a3SJed Brown 
1214b1a0a8a3SJed Brown     /* Partitioning by adjacency of diagonal block  */
1215b1a0a8a3SJed Brown 
1216b1a0a8a3SJed Brown     const PetscInt *numbering;
1217b1a0a8a3SJed Brown     PetscInt       *count,nidx,*indices,*newidx,start=0;
1218b1a0a8a3SJed Brown     /* Get node count in each partition */
1219b1a0a8a3SJed Brown     ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr);
1220b1a0a8a3SJed Brown     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
1221b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1222b1a0a8a3SJed Brown       for (i=0; i<n; i++) count[i] *= bs;
1223b1a0a8a3SJed Brown     }
1224b1a0a8a3SJed Brown     /* Build indices from node numbering */
1225b1a0a8a3SJed Brown     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1226b1a0a8a3SJed Brown     ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr);
1227b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1228b1a0a8a3SJed Brown     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
1229b1a0a8a3SJed Brown     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
1230b1a0a8a3SJed Brown     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
1231b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1232b1a0a8a3SJed Brown       ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr);
1233b1a0a8a3SJed Brown       for (i=0; i<nidx; i++)
1234b1a0a8a3SJed Brown 	for (j=0; j<bs; j++)
1235b1a0a8a3SJed Brown 	  newidx[i*bs+j] = indices[i]*bs + j;
1236b1a0a8a3SJed Brown       ierr = PetscFree(indices);CHKERRQ(ierr);
1237b1a0a8a3SJed Brown       nidx   *= bs;
1238b1a0a8a3SJed Brown       indices = newidx;
1239b1a0a8a3SJed Brown     }
1240b1a0a8a3SJed Brown     /* Shift to get global indices */
1241b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] += rstart;
1242b1a0a8a3SJed Brown 
1243b1a0a8a3SJed Brown     /* Build the index sets for each block */
1244b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1245b1a0a8a3SJed Brown       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
1246b1a0a8a3SJed Brown       ierr   = ISSort(is[i]);CHKERRQ(ierr);
1247b1a0a8a3SJed Brown       start += count[i];
1248b1a0a8a3SJed Brown     }
1249b1a0a8a3SJed Brown 
1250b1a0a8a3SJed Brown     ierr = PetscFree(count);
1251b1a0a8a3SJed Brown     ierr = PetscFree(indices);
1252b1a0a8a3SJed Brown     ierr = ISDestroy(isnumb);CHKERRQ(ierr);
1253b1a0a8a3SJed Brown     ierr = ISDestroy(ispart);CHKERRQ(ierr);
1254b1a0a8a3SJed Brown   }
1255b1a0a8a3SJed Brown 
1256b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1257b1a0a8a3SJed Brown }
1258b1a0a8a3SJed Brown 
1259b1a0a8a3SJed Brown #undef __FUNCT__
1260f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains"
1261b1a0a8a3SJed Brown /*@C
1262f746d493SDmitry Karpeev    PCGASMDestroySubdomains - Destroys the index sets created with
1263f746d493SDmitry Karpeev    PCGASMCreateSubdomains(). Should be called after setting subdomains
1264f746d493SDmitry Karpeev    with PCGASMSetLocalSubdomains().
1265b1a0a8a3SJed Brown 
1266b1a0a8a3SJed Brown    Collective
1267b1a0a8a3SJed Brown 
1268b1a0a8a3SJed Brown    Input Parameters:
1269b1a0a8a3SJed Brown +  n - the number of index sets
1270b1a0a8a3SJed Brown .  is - the array of index sets
1271b1a0a8a3SJed Brown -  is_local - the array of local index sets, can be PETSC_NULL
1272b1a0a8a3SJed Brown 
1273b1a0a8a3SJed Brown    Level: advanced
1274b1a0a8a3SJed Brown 
1275f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1276b1a0a8a3SJed Brown 
1277f746d493SDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetLocalSubdomains()
1278b1a0a8a3SJed Brown @*/
12797087cfbeSBarry Smith PetscErrorCode  PCGASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
1280b1a0a8a3SJed Brown {
1281b1a0a8a3SJed Brown   PetscInt       i;
1282b1a0a8a3SJed Brown   PetscErrorCode ierr;
1283b1a0a8a3SJed Brown   PetscFunctionBegin;
1284b1a0a8a3SJed Brown   if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"n must be > 0: n = %D",n);
1285b1a0a8a3SJed Brown   PetscValidPointer(is,2);
1286b1a0a8a3SJed Brown   for (i=0; i<n; i++) { ierr = ISDestroy(is[i]);CHKERRQ(ierr); }
1287b1a0a8a3SJed Brown   ierr = PetscFree(is);CHKERRQ(ierr);
1288b1a0a8a3SJed Brown   if (is_local) {
1289b1a0a8a3SJed Brown     PetscValidPointer(is_local,3);
1290b1a0a8a3SJed Brown     for (i=0; i<n; i++) { ierr = ISDestroy(is_local[i]);CHKERRQ(ierr); }
1291b1a0a8a3SJed Brown     ierr = PetscFree(is_local);CHKERRQ(ierr);
1292b1a0a8a3SJed Brown   }
1293b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1294b1a0a8a3SJed Brown }
1295b1a0a8a3SJed Brown 
1296af538404SDmitry Karpeev 
1297af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \
1298af538404SDmitry Karpeev {                                                                                                       \
1299af538404SDmitry Karpeev  PetscInt first_row = first/M, last_row = last/M+1;                                                     \
1300af538404SDmitry Karpeev   /*                                                                                                    \
1301af538404SDmitry Karpeev    Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners       \
1302af538404SDmitry Karpeev    of the bounding box of the intersection of the subdomain with the local ownership range (local       \
1303af538404SDmitry Karpeev    subdomain).                                                                                          \
1304af538404SDmitry Karpeev    Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows       \
1305af538404SDmitry Karpeev    of the intersection.                                                                                 \
1306af538404SDmitry Karpeev   */                                                                                                    \
1307af538404SDmitry Karpeev   /* ylow_loc is the grid row containing the first element of the local sumbdomain */                   \
1308eec7e3faSDmitry Karpeev   *ylow_loc = PetscMax(first_row,ylow);                                                                    \
1309af538404SDmitry Karpeev   /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \
1310eec7e3faSDmitry Karpeev   *xleft_loc = *ylow_loc==first_row?PetscMax(first%M,xleft):xleft;                                                                            \
1311af538404SDmitry Karpeev   /* yhigh_loc is the grid row above the last local subdomain element */                                                                    \
1312eec7e3faSDmitry Karpeev   *yhigh_loc = PetscMin(last_row,yhigh);                                                                                                     \
1313af538404SDmitry Karpeev   /* xright is the offset of the end of the  local subdomain within its grid row (might actually be outside the local subdomain) */         \
1314eec7e3faSDmitry Karpeev   *xright_loc = *yhigh_loc==last_row?PetscMin(xright,last%M):xright;                                                                          \
1315af538404SDmitry Karpeev   /* Now compute the size of the local subdomain n. */ \
1316c3518bceSDmitry Karpeev   *n = 0;                                               \
1317eec7e3faSDmitry Karpeev   if(*ylow_loc < *yhigh_loc) {                           \
1318af538404SDmitry Karpeev     PetscInt width = xright-xleft;                     \
1319eec7e3faSDmitry Karpeev     *n += width*(*yhigh_loc-*ylow_loc-1);                 \
1320eec7e3faSDmitry Karpeev     *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \
1321eec7e3faSDmitry Karpeev     *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \
1322af538404SDmitry Karpeev   }\
1323af538404SDmitry Karpeev }
1324af538404SDmitry Karpeev 
1325c3518bceSDmitry Karpeev 
1326eec7e3faSDmitry Karpeev 
1327eec7e3faSDmitry Karpeev #undef __FUNCT__
1328f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D"
1329b1a0a8a3SJed Brown /*@
1330f746d493SDmitry Karpeev    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1331b1a0a8a3SJed Brown    preconditioner for a two-dimensional problem on a regular grid.
1332b1a0a8a3SJed Brown 
1333af538404SDmitry Karpeev    Collective
1334b1a0a8a3SJed Brown 
1335b1a0a8a3SJed Brown    Input Parameters:
1336af538404SDmitry Karpeev +  M, N - the global number of mesh points in the x and y directions
1337af538404SDmitry Karpeev .  Mdomains, Ndomains - the global number of subdomains in the x and y directions
1338b1a0a8a3SJed Brown .  dof - degrees of freedom per node
1339b1a0a8a3SJed Brown -  overlap - overlap in mesh lines
1340b1a0a8a3SJed Brown 
1341b1a0a8a3SJed Brown    Output Parameters:
1342af538404SDmitry Karpeev +  Nsub - the number of local subdomains created
1343b1a0a8a3SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
1344b1a0a8a3SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
1345b1a0a8a3SJed Brown 
1346b1a0a8a3SJed Brown    Note:
1347b1a0a8a3SJed Brown    Presently PCAMSCreateSubdomains2d() is valid only for sequential
1348b1a0a8a3SJed Brown    preconditioners.  More general related routines are
1349f746d493SDmitry Karpeev    PCGASMSetTotalSubdomains() and PCGASMSetLocalSubdomains().
1350b1a0a8a3SJed Brown 
1351b1a0a8a3SJed Brown    Level: advanced
1352b1a0a8a3SJed Brown 
1353f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid
1354b1a0a8a3SJed Brown 
1355f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetLocalSubdomains(), PCGASMGetSubKSP(),
1356f746d493SDmitry Karpeev           PCGASMSetOverlap()
1357b1a0a8a3SJed Brown @*/
13587087cfbeSBarry Smith PetscErrorCode  PCGASMCreateSubdomains2D(PC pc, PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap, PetscInt *nsub,IS **is,IS **is_local)
1359b1a0a8a3SJed Brown {
1360b1a0a8a3SJed Brown   PetscErrorCode ierr;
13612bbb417fSDmitry Karpeev   PetscMPIInt size, rank;
13622bbb417fSDmitry Karpeev   PetscInt i, j;
13632bbb417fSDmitry Karpeev   PetscInt maxheight, maxwidth;
13642bbb417fSDmitry Karpeev   PetscInt xstart, xleft, xright, xleft_loc, xright_loc;
13652bbb417fSDmitry Karpeev   PetscInt ystart, ylow,  yhigh,  ylow_loc,  yhigh_loc;
1366eec7e3faSDmitry Karpeev   PetscInt x[2][2], y[2][2], n[2];
13672bbb417fSDmitry Karpeev   PetscInt first, last;
13682bbb417fSDmitry Karpeev   PetscInt nidx, *idx;
13692bbb417fSDmitry Karpeev   PetscInt ii,jj,s,q,d;
1370af538404SDmitry Karpeev   PetscInt k,kk;
13712bbb417fSDmitry Karpeev   PetscMPIInt color;
13722bbb417fSDmitry Karpeev   MPI_Comm comm, subcomm;
1373af538404SDmitry Karpeev 
13742bbb417fSDmitry Karpeev   IS **iis;
1375b1a0a8a3SJed Brown   PetscFunctionBegin;
13762bbb417fSDmitry Karpeev   ierr = PetscObjectGetComm((PetscObject)pc, &comm);     CHKERRQ(ierr);
13772bbb417fSDmitry Karpeev   ierr = MPI_Comm_size(comm, &size);                     CHKERRQ(ierr);
13782bbb417fSDmitry Karpeev   ierr = MPI_Comm_rank(comm, &rank);                     CHKERRQ(ierr);
13792bbb417fSDmitry Karpeev   ierr = MatGetOwnershipRange(pc->pmat, &first, &last);  CHKERRQ(ierr);
13802bbb417fSDmitry Karpeev   if(first%dof || last%dof) {
13812bbb417fSDmitry Karpeev     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Matrix row partitioning unsuitable for domain decomposition: local row range (%D,%D) "
13822bbb417fSDmitry Karpeev 	     "does not respect the number of degrees of freedom per grid point %D", first, last, dof);
13832bbb417fSDmitry Karpeev   }
1384af538404SDmitry Karpeev   /* Determine the number of domains with nonzero intersections with the local ownership range. */
13852bbb417fSDmitry Karpeev   s = 0;
1386af538404SDmitry Karpeev   ystart = 0;
1387af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1388af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1389af538404SDmitry 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);
1390eec7e3faSDmitry Karpeev     /* Vertical domain limits with an overlap. */
1391eec7e3faSDmitry Karpeev     ylow = PetscMax(ystart - overlap,0);
1392eec7e3faSDmitry Karpeev     yhigh = PetscMin(ystart + maxheight + overlap,N);
1393b1a0a8a3SJed Brown     xstart = 0;
1394af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1395af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1396af538404SDmitry 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);
1397eec7e3faSDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1398eec7e3faSDmitry Karpeev       xleft   = PetscMax(xstart - overlap,0);
1399eec7e3faSDmitry Karpeev       xright  = PetscMin(xstart + maxwidth + overlap,M);
1400af538404SDmitry Karpeev       /*
1401af538404SDmitry Karpeev 	 Determine whether this subdomain intersects this processor's ownership range of pc->pmat.
1402af538404SDmitry Karpeev       */
1403c3518bceSDmitry Karpeev       PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1404af538404SDmitry Karpeev       if(nidx) {
1405af538404SDmitry Karpeev         ++s;
1406af538404SDmitry Karpeev       }
1407af538404SDmitry Karpeev       xstart += maxwidth;
1408af538404SDmitry Karpeev     }/* for(i = 0; i < Mdomains; ++i) */
1409af538404SDmitry Karpeev     ystart += maxheight;
1410af538404SDmitry Karpeev   }/* for(j = 0; j < Ndomains; ++j) */
1411af538404SDmitry Karpeev   /* Now we can allocate the necessary number of ISs. */
1412af538404SDmitry Karpeev   *nsub = s;
1413af538404SDmitry Karpeev   ierr = PetscMalloc((*nsub)*sizeof(IS*),is);         CHKERRQ(ierr);
1414af538404SDmitry Karpeev   ierr = PetscMalloc((*nsub)*sizeof(IS*),is_local);   CHKERRQ(ierr);
1415af538404SDmitry Karpeev   s = 0;
1416af538404SDmitry Karpeev   ystart = 0;
1417af538404SDmitry Karpeev   for (j=0; j<Ndomains; ++j) {
1418af538404SDmitry Karpeev     maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */
1419af538404SDmitry 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);
1420af538404SDmitry Karpeev     /* Vertical domain limits with an overlap. */
1421eec7e3faSDmitry Karpeev     y[0][0] = PetscMax(ystart - overlap,0);
1422eec7e3faSDmitry Karpeev     y[0][1] = PetscMin(ystart + maxheight + overlap,N);
1423eec7e3faSDmitry Karpeev     /* Vertical domain limits without an overlap. */
1424eec7e3faSDmitry Karpeev     y[1][0] = ystart;
1425eec7e3faSDmitry Karpeev     y[1][1] = PetscMin(ystart + maxheight,N);
1426eec7e3faSDmitry Karpeev     xstart = 0;
1427af538404SDmitry Karpeev     for (i=0; i<Mdomains; ++i) {
1428af538404SDmitry Karpeev       maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */
1429af538404SDmitry 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);
1430af538404SDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1431eec7e3faSDmitry Karpeev       x[0][0]  = PetscMax(xstart - overlap,0);
1432eec7e3faSDmitry Karpeev       x[0][1]  = PetscMin(xstart + maxwidth + overlap,M);
1433eec7e3faSDmitry Karpeev       /* Horizontal domain limits without an overlap. */
1434eec7e3faSDmitry Karpeev       x[1][0] = xstart;
1435eec7e3faSDmitry Karpeev       x[1][1] = PetscMin(xstart+maxwidth,M);
14362bbb417fSDmitry Karpeev       /*
14372bbb417fSDmitry Karpeev 	 Determine whether this domain intersects this processor's ownership range of pc->pmat.
14382bbb417fSDmitry Karpeev 	 Do this twice: first for the domains with overlaps, and once without.
14392bbb417fSDmitry Karpeev 	 During the first pass create the subcommunicators, and use them on the second pass as well.
14402bbb417fSDmitry Karpeev       */
14412bbb417fSDmitry Karpeev       for(q = 0; q < 2; ++q) {
14422bbb417fSDmitry Karpeev 	/*
14432bbb417fSDmitry Karpeev 	  domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
14442bbb417fSDmitry Karpeev 	  according to whether the domain with an overlap or without is considered.
14452bbb417fSDmitry Karpeev 	*/
14462bbb417fSDmitry Karpeev 	xleft = x[q][0]; xright = x[q][1];
14472bbb417fSDmitry Karpeev 	ylow  = y[q][0]; yhigh  = y[q][1];
1448c3518bceSDmitry Karpeev         PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx));
1449af538404SDmitry Karpeev 	nidx *= dof;
1450eec7e3faSDmitry Karpeev         n[q] = nidx;
14512bbb417fSDmitry Karpeev         /*
1452eec7e3faSDmitry Karpeev          Based on the counted number of indices in the local domain *with an overlap*,
1453af538404SDmitry Karpeev          construct a subcommunicator of all the processors supporting this domain.
1454eec7e3faSDmitry Karpeev          Observe that a domain with an overlap might have nontrivial local support,
1455eec7e3faSDmitry Karpeev          while the domain without an overlap might not.  Hence, the decision to participate
1456eec7e3faSDmitry Karpeev          in the subcommunicator must be based on the domain with an overlap.
14572bbb417fSDmitry Karpeev          */
14582bbb417fSDmitry Karpeev 	if(q == 0) {
14592bbb417fSDmitry Karpeev 	  if(nidx) {
14602bbb417fSDmitry Karpeev 	    color = 1;
1461b1a0a8a3SJed Brown 	  }
14622bbb417fSDmitry Karpeev 	  else {
14632bbb417fSDmitry Karpeev 	    color = MPI_UNDEFINED;
14642bbb417fSDmitry Karpeev 	  }
14652bbb417fSDmitry Karpeev 	  ierr = MPI_Comm_split(comm, color, rank, &subcomm); CHKERRQ(ierr);
14662bbb417fSDmitry Karpeev 	}
1467af538404SDmitry Karpeev         /*
1468eec7e3faSDmitry Karpeev          Proceed only if the number of local indices *with an overlap* is nonzero.
1469af538404SDmitry Karpeev          */
1470eec7e3faSDmitry Karpeev         if(n[0]) {
1471af538404SDmitry Karpeev           if(q == 0) {
1472eec7e3faSDmitry Karpeev             iis = is;
1473af538404SDmitry Karpeev           }
1474af538404SDmitry Karpeev           if(q == 1) {
1475af538404SDmitry Karpeev             /*
1476eec7e3faSDmitry Karpeev              The IS for the no-overlap subdomain shares a communicator with the overlapping domain.
1477af538404SDmitry Karpeev              Moreover, if the overlap is zero, the two ISs are identical.
1478af538404SDmitry Karpeev              */
1479b1a0a8a3SJed Brown             if(overlap == 0) {
1480eec7e3faSDmitry Karpeev               (*is_local)[s] = (*is)[s];
14812bbb417fSDmitry Karpeev               ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr);
14822bbb417fSDmitry Karpeev               continue;
14832bbb417fSDmitry Karpeev             }
1484af538404SDmitry Karpeev             else {
1485eec7e3faSDmitry Karpeev               iis = is_local;
1486eec7e3faSDmitry Karpeev               subcomm = ((PetscObject)(*is)[s])->comm;
14872bbb417fSDmitry Karpeev             }
1488af538404SDmitry Karpeev           }/* if(q == 1) */
1489eec7e3faSDmitry Karpeev           idx = PETSC_NULL;
14902bbb417fSDmitry Karpeev 	  ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1491eec7e3faSDmitry Karpeev           if(nidx) {
14922bbb417fSDmitry Karpeev             k    = 0;
14932bbb417fSDmitry Karpeev             for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1494af538404SDmitry Karpeev               PetscInt x0 = (jj==ylow_loc)?xleft_loc:xleft;
1495af538404SDmitry Karpeev               PetscInt x1 = (jj==yhigh_loc-1)?xright_loc:xright;
1496af538404SDmitry Karpeev               kk = dof*(M*jj + x0);
14972bbb417fSDmitry Karpeev               for (ii=x0; ii<x1; ++ii) {
14982bbb417fSDmitry Karpeev                 for(d = 0; d < dof; ++d) {
14992bbb417fSDmitry Karpeev                   idx[k++] = kk++;
1500b1a0a8a3SJed Brown                 }
1501b1a0a8a3SJed Brown               }
1502b1a0a8a3SJed Brown             }
1503eec7e3faSDmitry Karpeev           }
15042bbb417fSDmitry Karpeev 	  ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*iis)+s);CHKERRQ(ierr);
1505eec7e3faSDmitry Karpeev 	}/* if(n[0]) */
15062bbb417fSDmitry Karpeev       }/* for(q = 0; q < 2; ++q) */
1507eec7e3faSDmitry Karpeev       if(n[0]) {
15082bbb417fSDmitry Karpeev         ++s;
1509b1a0a8a3SJed Brown       }
1510af538404SDmitry Karpeev       xstart += maxwidth;
1511af538404SDmitry Karpeev     }/* for(i = 0; i < Mdomains; ++i) */
15122bbb417fSDmitry Karpeev     ystart += maxheight;
1513af538404SDmitry Karpeev   }/* for(j = 0; j < Ndomains; ++j) */
1514b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1515b1a0a8a3SJed Brown }
1516b1a0a8a3SJed Brown 
1517b1a0a8a3SJed Brown #undef __FUNCT__
1518f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetLocalSubdomains"
1519b1a0a8a3SJed Brown /*@C
1520f746d493SDmitry Karpeev     PCGASMGetLocalSubdomains - Gets the local subdomains (for this processor
1521b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1522b1a0a8a3SJed Brown 
1523b1a0a8a3SJed Brown     Not Collective
1524b1a0a8a3SJed Brown 
1525b1a0a8a3SJed Brown     Input Parameter:
1526b1a0a8a3SJed Brown .   pc - the preconditioner context
1527b1a0a8a3SJed Brown 
1528b1a0a8a3SJed Brown     Output Parameters:
1529b1a0a8a3SJed Brown +   n - the number of subdomains for this processor (default value = 1)
1530b1a0a8a3SJed Brown .   is - the index sets that define the subdomains for this processor
1531b1a0a8a3SJed Brown -   is_local - the index sets that define the local part of the subdomains for this processor (can be PETSC_NULL)
1532b1a0a8a3SJed Brown 
1533b1a0a8a3SJed Brown 
1534b1a0a8a3SJed Brown     Notes:
1535b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
1536b1a0a8a3SJed Brown 
1537b1a0a8a3SJed Brown     Level: advanced
1538b1a0a8a3SJed Brown 
1539f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz
1540b1a0a8a3SJed Brown 
1541f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
1542f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetLocalSubdomains(), PCGASMGetLocalSubmatrices()
1543b1a0a8a3SJed Brown @*/
15447087cfbeSBarry Smith PetscErrorCode  PCGASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
1545b1a0a8a3SJed Brown {
1546f746d493SDmitry Karpeev   PC_GASM         *osm;
1547b1a0a8a3SJed Brown   PetscErrorCode ierr;
1548b1a0a8a3SJed Brown   PetscBool      match;
1549b1a0a8a3SJed Brown 
1550b1a0a8a3SJed Brown   PetscFunctionBegin;
1551b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1552b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1553b1a0a8a3SJed Brown   if (is) PetscValidPointer(is,3);
1554f746d493SDmitry Karpeev   ierr = PetscTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1555b1a0a8a3SJed Brown   if (!match) {
1556b1a0a8a3SJed Brown     if (n)  *n  = 0;
1557b1a0a8a3SJed Brown     if (is) *is = PETSC_NULL;
1558b1a0a8a3SJed Brown   } else {
1559f746d493SDmitry Karpeev     osm = (PC_GASM*)pc->data;
1560ab3e923fSDmitry Karpeev     if (n)  *n  = osm->n;
1561b1a0a8a3SJed Brown     if (is) *is = osm->is;
1562b1a0a8a3SJed Brown     if (is_local) *is_local = osm->is_local;
1563b1a0a8a3SJed Brown   }
1564b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1565b1a0a8a3SJed Brown }
1566b1a0a8a3SJed Brown 
1567b1a0a8a3SJed Brown #undef __FUNCT__
1568f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetLocalSubmatrices"
1569b1a0a8a3SJed Brown /*@C
1570f746d493SDmitry Karpeev     PCGASMGetLocalSubmatrices - Gets the local submatrices (for this processor
1571b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1572b1a0a8a3SJed Brown 
1573b1a0a8a3SJed Brown     Not Collective
1574b1a0a8a3SJed Brown 
1575b1a0a8a3SJed Brown     Input Parameter:
1576b1a0a8a3SJed Brown .   pc - the preconditioner context
1577b1a0a8a3SJed Brown 
1578b1a0a8a3SJed Brown     Output Parameters:
1579b1a0a8a3SJed Brown +   n - the number of matrices for this processor (default value = 1)
1580b1a0a8a3SJed Brown -   mat - the matrices
1581b1a0a8a3SJed Brown 
1582b1a0a8a3SJed Brown 
1583b1a0a8a3SJed Brown     Level: advanced
1584b1a0a8a3SJed Brown 
1585f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi
1586b1a0a8a3SJed Brown 
1587f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
1588f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetLocalSubdomains(), PCGASMGetLocalSubdomains()
1589b1a0a8a3SJed Brown @*/
15907087cfbeSBarry Smith PetscErrorCode  PCGASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1591b1a0a8a3SJed Brown {
1592f746d493SDmitry Karpeev   PC_GASM         *osm;
1593b1a0a8a3SJed Brown   PetscErrorCode ierr;
1594b1a0a8a3SJed Brown   PetscBool      match;
1595b1a0a8a3SJed Brown 
1596b1a0a8a3SJed Brown   PetscFunctionBegin;
1597b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1598b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1599b1a0a8a3SJed Brown   if (mat) PetscValidPointer(mat,3);
1600b1a0a8a3SJed Brown   if (!pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1601f746d493SDmitry Karpeev   ierr = PetscTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1602b1a0a8a3SJed Brown   if (!match) {
1603b1a0a8a3SJed Brown     if (n)   *n   = 0;
1604b1a0a8a3SJed Brown     if (mat) *mat = PETSC_NULL;
1605b1a0a8a3SJed Brown   } else {
1606f746d493SDmitry Karpeev     osm = (PC_GASM*)pc->data;
1607ab3e923fSDmitry Karpeev     if (n)   *n   = osm->n;
1608b1a0a8a3SJed Brown     if (mat) *mat = osm->pmat;
1609b1a0a8a3SJed Brown   }
1610b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1611b1a0a8a3SJed Brown }
1612