xref: /petsc/src/ksp/pc/impls/gasm/gasm.c (revision 2bbb417fbd5e323b52bae6f285fe756e6fb6203a)
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 
8f746d493SDmitry Karpeev        N    - total number of true 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__
33f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM"
34f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer)
35b1a0a8a3SJed Brown {
36f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
37b1a0a8a3SJed Brown   PetscErrorCode ierr;
38b1a0a8a3SJed Brown   PetscMPIInt    rank;
39b1a0a8a3SJed Brown   PetscInt       i,bsz;
40b1a0a8a3SJed Brown   PetscBool      iascii,isstring;
41b1a0a8a3SJed Brown   PetscViewer    sviewer;
42b1a0a8a3SJed Brown 
43b1a0a8a3SJed Brown 
44b1a0a8a3SJed Brown   PetscFunctionBegin;
45b1a0a8a3SJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
46b1a0a8a3SJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
47b1a0a8a3SJed Brown   if (iascii) {
48b1a0a8a3SJed Brown     char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set";
49b1a0a8a3SJed Brown     if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof overlaps,"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);}
5044fe18e5SDmitry Karpeev     if (osm->nmax > 0)     {ierr = PetscSNPrintf(blocks,sizeof blocks,"max number of  subdomain blocks = %D",osm->nmax);CHKERRQ(ierr);}
51f746d493SDmitry Karpeev     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,osm->n);CHKERRQ(ierr);
52f746d493SDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"  Generalized additive Schwarz: %s, %s\n",blocks,overlaps);CHKERRQ(ierr);
53f746d493SDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"  Generalized additive Schwarz: restriction/interpolation type - %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr);
54b1a0a8a3SJed Brown     ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
55b1a0a8a3SJed Brown     if (osm->same_local_solves) {
56b1a0a8a3SJed Brown       if (osm->ksp) {
57b1a0a8a3SJed Brown         ierr = PetscViewerASCIIPrintf(viewer,"  Local solve is same for all blocks, in the following KSP and PC objects:\n");CHKERRQ(ierr);
58b1a0a8a3SJed Brown         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
59b1a0a8a3SJed Brown         if (!rank) {
60b1a0a8a3SJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
61b1a0a8a3SJed Brown           ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);
62b1a0a8a3SJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
63b1a0a8a3SJed Brown         }
64b1a0a8a3SJed Brown         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
65b1a0a8a3SJed Brown       }
66b1a0a8a3SJed Brown     } else {
67b1a0a8a3SJed Brown       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
68b1a0a8a3SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr);
69b1a0a8a3SJed Brown       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
70b1a0a8a3SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
71f746d493SDmitry Karpeev       for (i=0; i<osm->nmax; i++) {
72b1a0a8a3SJed Brown         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
73f746d493SDmitry Karpeev         if (i < osm->n) {
74b1a0a8a3SJed Brown           ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr);
75b1a0a8a3SJed Brown           ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr);
76b1a0a8a3SJed Brown           ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr);
77b1a0a8a3SJed Brown           ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
78b1a0a8a3SJed Brown         }
79b1a0a8a3SJed Brown         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
80b1a0a8a3SJed Brown       }
81b1a0a8a3SJed Brown       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
82b1a0a8a3SJed Brown       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
83b1a0a8a3SJed Brown     }
84b1a0a8a3SJed Brown   } else if (isstring) {
85f746d493SDmitry Karpeev     ierr = PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCGASMTypes[osm->type]);CHKERRQ(ierr);
86b1a0a8a3SJed Brown     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
87b1a0a8a3SJed Brown     if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);}
88b1a0a8a3SJed Brown     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
89b1a0a8a3SJed Brown   } else {
90f746d493SDmitry Karpeev     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCGASM",((PetscObject)viewer)->type_name);
91b1a0a8a3SJed Brown   }
92b1a0a8a3SJed Brown   PetscFunctionReturn(0);
93b1a0a8a3SJed Brown }
94b1a0a8a3SJed Brown 
95b1a0a8a3SJed Brown #undef __FUNCT__
96f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains"
97f746d493SDmitry Karpeev static PetscErrorCode PCGASMPrintSubdomains(PC pc)
98b1a0a8a3SJed Brown {
99f746d493SDmitry Karpeev   PC_GASM         *osm  = (PC_GASM*)pc->data;
100b1a0a8a3SJed Brown   const char     *prefix;
101b1a0a8a3SJed Brown   char           fname[PETSC_MAX_PATH_LEN+1];
102b1a0a8a3SJed Brown   PetscViewer    viewer;
103b1a0a8a3SJed Brown   PetscInt       i,j,nidx;
104b1a0a8a3SJed Brown   const PetscInt *idx;
105b1a0a8a3SJed Brown   PetscErrorCode ierr;
106b1a0a8a3SJed Brown 
107b1a0a8a3SJed Brown   PetscFunctionBegin;
108b1a0a8a3SJed Brown   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
109f746d493SDmitry Karpeev   ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,PETSC_NULL);CHKERRQ(ierr);
110b1a0a8a3SJed Brown   if (fname[0] == 0) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
111b1a0a8a3SJed Brown   ierr = PetscViewerASCIIOpen(((PetscObject)pc)->comm,fname,&viewer);CHKERRQ(ierr);
112f746d493SDmitry Karpeev   for (i=0;i<osm->n;++i) {
113b1a0a8a3SJed Brown     ierr = ISGetLocalSize(osm->is[i],&nidx);CHKERRQ(ierr);
114b1a0a8a3SJed Brown     ierr = ISGetIndices(osm->is[i],&idx);CHKERRQ(ierr);
115b1a0a8a3SJed Brown     for (j=0; j<nidx; j++) {
116b1a0a8a3SJed Brown       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%D ",idx[j]);CHKERRQ(ierr);
117b1a0a8a3SJed Brown     }
118b1a0a8a3SJed Brown     ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr);
119b1a0a8a3SJed Brown     ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
120b1a0a8a3SJed Brown     if (osm->is_local) {
121b1a0a8a3SJed Brown       ierr = ISGetLocalSize(osm->is_local[i],&nidx);CHKERRQ(ierr);
122b1a0a8a3SJed Brown       ierr = ISGetIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
123b1a0a8a3SJed Brown       for (j=0; j<nidx; j++) {
124b1a0a8a3SJed Brown         ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%D ",idx[j]);CHKERRQ(ierr);
125b1a0a8a3SJed Brown       }
126b1a0a8a3SJed Brown       ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
127b1a0a8a3SJed Brown       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr);
128b1a0a8a3SJed Brown     }
129b1a0a8a3SJed Brown   }
130b1a0a8a3SJed Brown   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
131b1a0a8a3SJed Brown   ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
132b1a0a8a3SJed Brown   PetscFunctionReturn(0);
133b1a0a8a3SJed Brown }
134b1a0a8a3SJed Brown 
135b1a0a8a3SJed Brown #undef __FUNCT__
136f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM"
137f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc)
138b1a0a8a3SJed Brown {
139f746d493SDmitry Karpeev   PC_GASM         *osm  = (PC_GASM*)pc->data;
140b1a0a8a3SJed Brown   PetscErrorCode ierr;
141b1a0a8a3SJed Brown   PetscBool      symset,flg;
14288389755SJed Brown   PetscInt       i;
143c10bc1a1SDmitry Karpeev   PetscMPIInt    rank, size;
144b1a0a8a3SJed Brown   MatReuse       scall = MAT_REUSE_MATRIX;
145b1a0a8a3SJed Brown   KSP            ksp;
146b1a0a8a3SJed Brown   PC             subpc;
147b1a0a8a3SJed Brown   const char     *prefix,*pprefix;
148f746d493SDmitry Karpeev   PetscInt       dn;       /* Number of indices in a single subdomain assigned to this processor. */
149f746d493SDmitry Karpeev   const PetscInt *didx;    /* Indices from a single subdomain assigned to this processor. */
150f746d493SDmitry Karpeev   PetscInt       ddn;      /* Number of indices in all subdomains assigned to this processor. */
151f746d493SDmitry Karpeev   PetscInt       *ddidx;   /* Indices of all subdomains assigned to this processor. */
152f746d493SDmitry Karpeev   IS             gid;      /* Identity IS of the size of all subdomains assigned to this processor. */
153f746d493SDmitry Karpeev   Vec            x,y;
154f746d493SDmitry Karpeev   PetscScalar    *gxarray, *gyarray;
1556d98aee9SDmitry Karpeev   PetscInt       gfirst, glast;
156b1a0a8a3SJed Brown 
157b1a0a8a3SJed Brown   PetscFunctionBegin;
158c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
159c10bc1a1SDmitry Karpeev   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
160b1a0a8a3SJed Brown   if (!pc->setupcalled) {
161b1a0a8a3SJed Brown 
162b1a0a8a3SJed Brown     if (!osm->type_set) {
163b1a0a8a3SJed Brown       ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
164f746d493SDmitry Karpeev       if (symset && flg) { osm->type = PC_GASM_BASIC; }
165b1a0a8a3SJed Brown     }
166b1a0a8a3SJed Brown 
16744fe18e5SDmitry Karpeev     if (osm->N == PETSC_DECIDE && osm->n < 1) {
168b1a0a8a3SJed Brown       /* no subdomains given, use one per processor */
16944fe18e5SDmitry Karpeev       osm->nmax = osm->n = 1;
170b1a0a8a3SJed Brown       ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
171f746d493SDmitry Karpeev       osm->N = size;
172f746d493SDmitry Karpeev     } else if (osm->N == PETSC_DECIDE) {
173b1a0a8a3SJed Brown       PetscInt inwork[2], outwork[2];
174f746d493SDmitry Karpeev       /* determine global number of subdomains and the max number of local subdomains */
175f746d493SDmitry Karpeev       inwork[0] = inwork[1] = osm->n;
176b1a0a8a3SJed Brown       ierr = MPI_Allreduce(inwork,outwork,1,MPIU_2INT,PetscMaxSum_Op,((PetscObject)pc)->comm);CHKERRQ(ierr);
177f746d493SDmitry Karpeev       osm->nmax = outwork[0];
178f746d493SDmitry Karpeev       osm->N    = outwork[1];
179b1a0a8a3SJed Brown     }
180b1a0a8a3SJed Brown     if (!osm->is){ /* create the index sets */
181f746d493SDmitry Karpeev       ierr = PCGASMCreateSubdomains(pc->pmat,osm->n,&osm->is);CHKERRQ(ierr);
182b1a0a8a3SJed Brown     }
18393c1ec32SDmitry Karpeev     if (!osm->is_local) {
18493c1ec32SDmitry Karpeev       /*
18593c1ec32SDmitry Karpeev 	 This indicates that osm->is should define a nonoverlapping decomposition
18693c1ec32SDmitry Karpeev 	 (there is no way to really guarantee that if subdomains are set by the user through PCGASMSetLocalSubdomains,
18793c1ec32SDmitry Karpeev 	  but the assumption is that either the user does the right thing, or subdomains in ossm->is have been created
18893c1ec32SDmitry Karpeev 	  via PCGASMCreateSubdomains, which guarantees a nonoverlapping decomposition).
18993c1ec32SDmitry Karpeev 	 Therefore, osm->is will be used to define osm->is_local.
19093c1ec32SDmitry Karpeev 	 If a nonzero overlap has been requested by the user, then osm->is will be expanded and will overlap,
19193c1ec32SDmitry Karpeev 	 so osm->is_local should obtain a copy of osm->is while they are still (presumably) nonoverlapping.
19293c1ec32SDmitry Karpeev 	 Otherwise (no overlap has been requested), osm->is_local are simply aliases for osm->is.
19393c1ec32SDmitry Karpeev       */
194f746d493SDmitry Karpeev       ierr = PetscMalloc(osm->n*sizeof(IS),&osm->is_local);CHKERRQ(ierr);
195f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
196b1a0a8a3SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
197b1a0a8a3SJed Brown           ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr);
198b1a0a8a3SJed Brown           ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr);
199b1a0a8a3SJed Brown         } else {
200b1a0a8a3SJed Brown           ierr = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr);
201b1a0a8a3SJed Brown           osm->is_local[i] = osm->is[i];
202b1a0a8a3SJed Brown         }
203b1a0a8a3SJed Brown       }
204b1a0a8a3SJed Brown     }
20593c1ec32SDmitry Karpeev     /* Beyond this point osm->is_local is not null. */
20693c1ec32SDmitry Karpeev     if (osm->overlap > 0) {
20793c1ec32SDmitry Karpeev       /* Extend the "overlapping" regions by a number of steps */
20893c1ec32SDmitry Karpeev       ierr = MatIncreaseOverlap(pc->pmat,osm->n,osm->is,osm->overlap);CHKERRQ(ierr);
20993c1ec32SDmitry Karpeev     }
210b1a0a8a3SJed Brown     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
211b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
212f746d493SDmitry Karpeev     ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&flg,PETSC_NULL);CHKERRQ(ierr);
213f746d493SDmitry Karpeev     if (flg) { ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); }
214b1a0a8a3SJed Brown 
215b1a0a8a3SJed Brown     if (osm->sort_indices) {
216f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
217b1a0a8a3SJed Brown         ierr = ISSort(osm->is[i]);CHKERRQ(ierr);
218b1a0a8a3SJed Brown 	ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr);
219b1a0a8a3SJed Brown       }
220b1a0a8a3SJed Brown     }
221f746d493SDmitry Karpeev     /* Merge the ISs, create merged vectors and scatter contexts. */
22293c1ec32SDmitry Karpeev     /* Restriction ISs. */
223f746d493SDmitry Karpeev     ddn = 0;
224f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
225f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn);          CHKERRQ(ierr);
226f746d493SDmitry Karpeev       ddn += dn;
227f746d493SDmitry Karpeev     }
228f746d493SDmitry Karpeev     ierr = PetscMalloc(ddn*sizeof(PetscInt), &ddidx); CHKERRQ(ierr);
229f746d493SDmitry Karpeev     ddn = 0;
230f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
231f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn); CHKERRQ(ierr);
232f746d493SDmitry Karpeev       ierr = ISGetIndices(osm->is[i],&didx); CHKERRQ(ierr);
233f746d493SDmitry Karpeev       ierr = PetscMemcpy(ddidx+ddn, didx, sizeof(PetscInt)*dn); CHKERRQ(ierr);
234f746d493SDmitry Karpeev       ierr = ISRestoreIndices(osm->is[i], &didx); CHKERRQ(ierr);
23593c1ec32SDmitry Karpeev       ddn += dn;
236f746d493SDmitry Karpeev     }
237f746d493SDmitry Karpeev     ierr = ISCreateGeneral(((PetscObject)(pc))->comm, ddn, ddidx, PETSC_OWN_POINTER, &osm->gis); CHKERRQ(ierr);
238f746d493SDmitry Karpeev     ierr = MatGetVecs(pc->pmat,&x,&y);CHKERRQ(ierr);
239f746d493SDmitry Karpeev     ierr = VecCreateMPI(((PetscObject)pc)->comm, ddn, PETSC_DECIDE, &osm->gx); CHKERRQ(ierr);
240f746d493SDmitry Karpeev     ierr = VecDuplicate(osm->gx,&osm->gy);                                     CHKERRQ(ierr);
2416d98aee9SDmitry Karpeev     ierr = VecGetOwnershipRange(osm->gx, &gfirst, &glast);                     CHKERRQ(ierr);
2426d98aee9SDmitry Karpeev     ierr = ISCreateStride(((PetscObject)pc)->comm,ddn,gfirst,1, &gid);              CHKERRQ(ierr);
243f746d493SDmitry Karpeev     ierr = VecScatterCreate(x,osm->gis,osm->gx,gid, &(osm->grestriction));     CHKERRQ(ierr);
244f746d493SDmitry Karpeev     ierr = ISDestroy(gid);                                                     CHKERRQ(ierr);
24593c1ec32SDmitry Karpeev     /* Prolongation ISs */
24681a5c4aaSDmitry Karpeev     { PetscInt       dn_local;       /* Number of indices in the local part of a single domain assigned to this processor. */
247f746d493SDmitry Karpeev       const PetscInt *didx_local;    /* Global indices from the local part of a single domain assigned to this processor. */
248f746d493SDmitry Karpeev       PetscInt       ddn_local;      /* Number of indices in the local part of the disjoint union all domains assigned to this processor. */
249f746d493SDmitry Karpeev       PetscInt       *ddidx_local;   /* Global indices of the local part of the disjoint union of all domains assigned to this processor. */
250f746d493SDmitry Karpeev       /**/
25193c1ec32SDmitry Karpeev       ISLocalToGlobalMapping ltog;          /* Map from global to local indices on the disjoint union of subdomains: "local" ind's run from 0 to ddn-1. */
25293c1ec32SDmitry Karpeev       PetscInt              *ddidx_llocal;  /* Mapped local indices of the disjoint union of local parts of subdomains. */
253f746d493SDmitry 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. */
254f746d493SDmitry Karpeev       IS                     gis_llocal;    /* IS with ddidx_llocal indices. */
255f746d493SDmitry Karpeev       PetscInt               j;
256f746d493SDmitry Karpeev       ddn_local = 0;
257f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
258f746d493SDmitry Karpeev 	ierr = ISGetLocalSize(osm->is_local[i],&dn_local); CHKERRQ(ierr);
259f746d493SDmitry Karpeev 	ddn_local += dn_local;
260f746d493SDmitry Karpeev       }
261f746d493SDmitry Karpeev       ierr = PetscMalloc(ddn_local*sizeof(PetscInt), &ddidx_local); CHKERRQ(ierr);
262f746d493SDmitry Karpeev       ddn_local = 0;
263f746d493SDmitry Karpeev       for (i=0; i<osm->n; i++) {
264f746d493SDmitry Karpeev 	ierr = ISGetLocalSize(osm->is_local[i],&dn_local); CHKERRQ(ierr);
265f746d493SDmitry Karpeev 	ierr = ISGetIndices(osm->is_local[i],&didx_local); CHKERRQ(ierr);
266f746d493SDmitry Karpeev 	ierr = PetscMemcpy(ddidx_local+ddn_local, didx_local, sizeof(PetscInt)*dn_local); CHKERRQ(ierr);
267f746d493SDmitry Karpeev 	ierr = ISRestoreIndices(osm->is_local[i], &didx_local); CHKERRQ(ierr);
268f746d493SDmitry Karpeev 	ddn_local += dn_local;
269f746d493SDmitry Karpeev       }
270ab3e923fSDmitry Karpeev       ierr = PetscMalloc(sizeof(PetscInt)*ddn_local, &ddidx_llocal); CHKERRQ(ierr);
271f746d493SDmitry Karpeev       ierr = ISCreateGeneral(((PetscObject)pc)->comm, ddn_local, ddidx_local, PETSC_OWN_POINTER, &(osm->gis_local)); CHKERRQ(ierr);
272f746d493SDmitry Karpeev       ierr = ISLocalToGlobalMappingCreateIS(osm->gis,&ltog);CHKERRQ(ierr);
27381a5c4aaSDmitry Karpeev       ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,ddn_local,ddidx_local,&ddn_llocal,ddidx_llocal);CHKERRQ(ierr);
274b1a0a8a3SJed Brown       ierr = ISLocalToGlobalMappingDestroy(ltog);CHKERRQ(ierr);
27593c1ec32SDmitry 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);
276f746d493SDmitry Karpeev       /* Now convert these localized indices into the global indices into the merged output vector. */
2776d98aee9SDmitry Karpeev       ierr = VecGetOwnershipRange(osm->gy, &gfirst, &glast); CHKERRQ(ierr);
278f746d493SDmitry Karpeev       for(j=0; j < ddn_llocal; ++j) {
2796d98aee9SDmitry Karpeev 	ddidx_llocal[j] += gfirst;
280f746d493SDmitry Karpeev       }
281f746d493SDmitry Karpeev       ierr = ISCreateGeneral(PETSC_COMM_SELF,ddn_llocal,ddidx_llocal,PETSC_OWN_POINTER,&gis_llocal); CHKERRQ(ierr);
282f746d493SDmitry Karpeev       ierr = VecScatterCreate(y,osm->gis_local,osm->gy,gis_llocal,&osm->gprolongation);              CHKERRQ(ierr);
283f746d493SDmitry Karpeev       ierr = ISDestroy(gis_llocal);CHKERRQ(ierr);
284b1a0a8a3SJed Brown     }
285f746d493SDmitry Karpeev     /* Create the subdomain work vectors. */
286f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->x);CHKERRQ(ierr);
287f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->y);CHKERRQ(ierr);
2886d98aee9SDmitry Karpeev     ierr = VecGetOwnershipRange(osm->gx, &gfirst, &glast);CHKERRQ(ierr);
289f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gx, &gxarray);                                     CHKERRQ(ierr);
290f746d493SDmitry Karpeev     ierr = VecGetArray(osm->gy, &gyarray);                                     CHKERRQ(ierr);
291f746d493SDmitry Karpeev     for (i=0, ddn=0; i<osm->n; ++i, ddn += dn) {
292f746d493SDmitry Karpeev       ierr = ISGetLocalSize(osm->is[i],&dn);CHKERRQ(ierr);
293f746d493SDmitry Karpeev       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dn,gxarray+ddn,&osm->x[i]);CHKERRQ(ierr);
294f746d493SDmitry Karpeev       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dn,gyarray+ddn,&osm->y[i]);CHKERRQ(ierr);
295b1a0a8a3SJed Brown     }
296f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gx, &gxarray); CHKERRQ(ierr);
297f746d493SDmitry Karpeev     ierr = VecRestoreArray(osm->gy, &gyarray); CHKERRQ(ierr);
298f746d493SDmitry Karpeev     ierr = VecDestroy(x); CHKERRQ(ierr);
299f746d493SDmitry Karpeev     ierr = VecDestroy(y); CHKERRQ(ierr);
300b1a0a8a3SJed Brown     /* Create the local solvers */
301f746d493SDmitry Karpeev     ierr = PetscMalloc(osm->n*sizeof(KSP *),&osm->ksp);CHKERRQ(ierr);
30244fe18e5SDmitry Karpeev     for (i=0; i<osm->n; i++) { /* KSPs are local */
303b1a0a8a3SJed Brown       ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
304b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr);
305b1a0a8a3SJed Brown       ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
306b1a0a8a3SJed Brown       ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
307b1a0a8a3SJed Brown       ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
308b1a0a8a3SJed Brown       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
309b1a0a8a3SJed Brown       ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
310b1a0a8a3SJed Brown       ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
311b1a0a8a3SJed Brown       osm->ksp[i] = ksp;
312b1a0a8a3SJed Brown     }
313b1a0a8a3SJed Brown     scall = MAT_INITIAL_MATRIX;
314b1a0a8a3SJed Brown 
315b1a0a8a3SJed Brown   } else {
316b1a0a8a3SJed Brown     /*
317b1a0a8a3SJed Brown        Destroy the blocks from the previous iteration
318b1a0a8a3SJed Brown     */
319b1a0a8a3SJed Brown     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
320f746d493SDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
321b1a0a8a3SJed Brown       scall = MAT_INITIAL_MATRIX;
322b1a0a8a3SJed Brown     }
323b1a0a8a3SJed Brown   }
324b1a0a8a3SJed Brown 
325b1a0a8a3SJed Brown   /*
326f746d493SDmitry Karpeev      Extract out the submatrices.
327b1a0a8a3SJed Brown   */
32881a5c4aaSDmitry Karpeev   if(size > 1) {
3296d42056aSDmitry Karpeev     ierr = MatGetSubMatricesParallel(pc->pmat,osm->n,osm->is, osm->is,scall,&osm->pmat);CHKERRQ(ierr);
33081a5c4aaSDmitry Karpeev   }
33181a5c4aaSDmitry Karpeev   else {
33281a5c4aaSDmitry Karpeev     ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->is, osm->is,scall,&osm->pmat);CHKERRQ(ierr);
33381a5c4aaSDmitry Karpeev   }
334b1a0a8a3SJed Brown   if (scall == MAT_INITIAL_MATRIX) {
335b1a0a8a3SJed Brown     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
336f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
337b1a0a8a3SJed Brown       ierr = PetscLogObjectParent(pc,osm->pmat[i]);CHKERRQ(ierr);
338b1a0a8a3SJed Brown       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
339b1a0a8a3SJed Brown     }
340b1a0a8a3SJed Brown   }
341b1a0a8a3SJed Brown 
342b1a0a8a3SJed Brown   /* Return control to the user so that the submatrices can be modified (e.g., to apply
343b1a0a8a3SJed Brown      different boundary conditions for the submatrices than for the global problem) */
344c10bc1a1SDmitry Karpeev   ierr = PCModifySubMatrices(pc,osm->n,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
345b1a0a8a3SJed Brown 
346b1a0a8a3SJed Brown   /*
347f746d493SDmitry Karpeev      Loop over submatrices putting them into local ksp
348b1a0a8a3SJed Brown   */
349f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
350b1a0a8a3SJed Brown     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr);
351b1a0a8a3SJed Brown     if (!pc->setupcalled) {
352b1a0a8a3SJed Brown       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
353b1a0a8a3SJed Brown     }
354b1a0a8a3SJed Brown   }
355b1a0a8a3SJed Brown 
356b1a0a8a3SJed Brown   PetscFunctionReturn(0);
357b1a0a8a3SJed Brown }
358b1a0a8a3SJed Brown 
359b1a0a8a3SJed Brown #undef __FUNCT__
360f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM"
361f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc)
362b1a0a8a3SJed Brown {
363f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
364b1a0a8a3SJed Brown   PetscErrorCode ierr;
365b1a0a8a3SJed Brown   PetscInt       i;
366b1a0a8a3SJed Brown 
367b1a0a8a3SJed Brown   PetscFunctionBegin;
368f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
369b1a0a8a3SJed Brown     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
370b1a0a8a3SJed Brown   }
371b1a0a8a3SJed Brown   PetscFunctionReturn(0);
372b1a0a8a3SJed Brown }
373b1a0a8a3SJed Brown 
374b1a0a8a3SJed Brown #undef __FUNCT__
375f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM"
376f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y)
377b1a0a8a3SJed Brown {
378f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
379b1a0a8a3SJed Brown   PetscErrorCode ierr;
380f746d493SDmitry Karpeev   PetscInt       i;
381b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
382b1a0a8a3SJed Brown 
383b1a0a8a3SJed Brown   PetscFunctionBegin;
384b1a0a8a3SJed Brown   /*
385b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
386b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
387b1a0a8a3SJed Brown   */
388f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
389b1a0a8a3SJed Brown     forward = SCATTER_FORWARD_LOCAL;
390b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
391f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx); CHKERRQ(ierr);
392b1a0a8a3SJed Brown   }
393f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
394b1a0a8a3SJed Brown     reverse = SCATTER_REVERSE_LOCAL;
395b1a0a8a3SJed Brown   }
396b1a0a8a3SJed Brown 
397f746d493SDmitry Karpeev   ierr = VecScatterBegin(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
398b1a0a8a3SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
399f746d493SDmitry Karpeev   ierr = VecScatterEnd(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
400b1a0a8a3SJed Brown   /* do the local solves */
401f746d493SDmitry 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. */
402b1a0a8a3SJed Brown     ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
403b1a0a8a3SJed Brown   }
404f746d493SDmitry Karpeev   ierr = VecScatterBegin(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
405f746d493SDmitry Karpeev   ierr = VecScatterEnd(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
406b1a0a8a3SJed Brown   PetscFunctionReturn(0);
407b1a0a8a3SJed Brown }
408b1a0a8a3SJed Brown 
409b1a0a8a3SJed Brown #undef __FUNCT__
410f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM"
411f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y)
412b1a0a8a3SJed Brown {
413f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
414b1a0a8a3SJed Brown   PetscErrorCode ierr;
415f746d493SDmitry Karpeev   PetscInt       i;
416b1a0a8a3SJed Brown   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
417b1a0a8a3SJed Brown 
418b1a0a8a3SJed Brown   PetscFunctionBegin;
419b1a0a8a3SJed Brown   /*
420b1a0a8a3SJed Brown      Support for limiting the restriction or interpolation to only local
421b1a0a8a3SJed Brown      subdomain values (leaving the other values 0).
422b1a0a8a3SJed Brown 
423f746d493SDmitry Karpeev      Note: these are reversed from the PCApply_GASM() because we are applying the
424b1a0a8a3SJed Brown      transpose of the three terms
425b1a0a8a3SJed Brown   */
426f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_INTERPOLATE)) {
427b1a0a8a3SJed Brown     forward = SCATTER_FORWARD_LOCAL;
428b1a0a8a3SJed Brown     /* have to zero the work RHS since scatter may leave some slots empty */
429f746d493SDmitry Karpeev     ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr);
430b1a0a8a3SJed Brown   }
431f746d493SDmitry Karpeev   if (!(osm->type & PC_GASM_RESTRICT)) {
432b1a0a8a3SJed Brown     reverse = SCATTER_REVERSE_LOCAL;
433b1a0a8a3SJed Brown   }
434b1a0a8a3SJed Brown 
435ab3e923fSDmitry Karpeev   ierr = VecScatterBegin(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
436b1a0a8a3SJed Brown   ierr = VecZeroEntries(y);CHKERRQ(ierr);
437ab3e923fSDmitry Karpeev   ierr = VecScatterEnd(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr);
438b1a0a8a3SJed Brown   /* do the local solves */
439ab3e923fSDmitry 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. */
440b1a0a8a3SJed Brown     ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr);
441b1a0a8a3SJed Brown   }
442ab3e923fSDmitry Karpeev   ierr = VecScatterBegin(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
443ab3e923fSDmitry Karpeev   ierr = VecScatterEnd(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr);
444b1a0a8a3SJed Brown   PetscFunctionReturn(0);
445b1a0a8a3SJed Brown }
446b1a0a8a3SJed Brown 
447b1a0a8a3SJed Brown #undef __FUNCT__
448f746d493SDmitry Karpeev #define __FUNCT__ "PCDestroy_GASM"
449f746d493SDmitry Karpeev static PetscErrorCode PCDestroy_GASM(PC pc)
450b1a0a8a3SJed Brown {
451f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
452b1a0a8a3SJed Brown   PetscErrorCode ierr;
453b1a0a8a3SJed Brown   PetscInt       i;
454b1a0a8a3SJed Brown 
455b1a0a8a3SJed Brown   PetscFunctionBegin;
456b1a0a8a3SJed Brown   if (osm->ksp) {
457f746d493SDmitry Karpeev     for (i=0; i<osm->n; i++) {
458b1a0a8a3SJed Brown       ierr = KSPDestroy(osm->ksp[i]);CHKERRQ(ierr);
459b1a0a8a3SJed Brown     }
460b1a0a8a3SJed Brown     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
461b1a0a8a3SJed Brown   }
462b1a0a8a3SJed Brown   if (osm->pmat) {
463f746d493SDmitry Karpeev     if (osm->n > 0) {
464ab3e923fSDmitry Karpeev       ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr);
465b1a0a8a3SJed Brown     }
466b1a0a8a3SJed Brown   }
467f746d493SDmitry Karpeev   for (i=0; i<osm->n; i++) {
468b1a0a8a3SJed Brown     ierr = VecDestroy(osm->x[i]);CHKERRQ(ierr);
469b1a0a8a3SJed Brown     ierr = VecDestroy(osm->y[i]);CHKERRQ(ierr);
470b1a0a8a3SJed Brown   }
471b1a0a8a3SJed Brown   ierr = PetscFree(osm->x);CHKERRQ(ierr);
472b1a0a8a3SJed Brown   ierr = PetscFree(osm->y);CHKERRQ(ierr);
473f746d493SDmitry Karpeev   ierr = VecDestroy(osm->gx); CHKERRQ(ierr);
474ab3e923fSDmitry Karpeev   ierr = VecDestroy(osm->gy); CHKERRQ(ierr);
475ab3e923fSDmitry Karpeev 
476ab3e923fSDmitry Karpeev   ierr = VecScatterDestroy(osm->grestriction); CHKERRQ(ierr);
477ab3e923fSDmitry Karpeev   ierr = VecScatterDestroy(osm->gprolongation); CHKERRQ(ierr);
478f746d493SDmitry Karpeev 
479b1a0a8a3SJed Brown   if (osm->is) {
480ab3e923fSDmitry Karpeev     ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
481b1a0a8a3SJed Brown   }
482b1a0a8a3SJed Brown   ierr = PetscFree(osm);CHKERRQ(ierr);
483b1a0a8a3SJed Brown   PetscFunctionReturn(0);
484b1a0a8a3SJed Brown }
485b1a0a8a3SJed Brown 
486b1a0a8a3SJed Brown #undef __FUNCT__
487f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM"
488ab3e923fSDmitry Karpeev static PetscErrorCode PCSetFromOptions_GASM(PC pc) {
489f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
490b1a0a8a3SJed Brown   PetscErrorCode ierr;
491b1a0a8a3SJed Brown   PetscInt       blocks,ovl;
492b1a0a8a3SJed Brown   PetscBool      symset,flg;
493f746d493SDmitry Karpeev   PCGASMType      gasmtype;
494b1a0a8a3SJed Brown 
495b1a0a8a3SJed Brown   PetscFunctionBegin;
496b1a0a8a3SJed Brown   /* set the type to symmetric if matrix is symmetric */
497b1a0a8a3SJed Brown   if (!osm->type_set && pc->pmat) {
498b1a0a8a3SJed Brown     ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr);
499f746d493SDmitry Karpeev     if (symset && flg) { osm->type = PC_GASM_BASIC; }
500b1a0a8a3SJed Brown   }
50144fe18e5SDmitry Karpeev   ierr = PetscOptionsHead("Generalized additive Schwarz options");CHKERRQ(ierr);
502f746d493SDmitry Karpeev     ierr = PetscOptionsInt("-pc_gasm_blocks","Number of subdomains","PCGASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
503ab3e923fSDmitry Karpeev     if (flg) {ierr = PCGASMSetTotalSubdomains(pc,blocks);CHKERRQ(ierr); }
504f746d493SDmitry Karpeev     ierr = PetscOptionsInt("-pc_gasm_overlap","Number of grid points overlap","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
505f746d493SDmitry Karpeev     if (flg) {ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); }
506b1a0a8a3SJed Brown     flg  = PETSC_FALSE;
507f746d493SDmitry Karpeev     ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr);
508f746d493SDmitry Karpeev     if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr); }
509b1a0a8a3SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
510b1a0a8a3SJed Brown   PetscFunctionReturn(0);
511b1a0a8a3SJed Brown }
512b1a0a8a3SJed Brown 
513b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/
514b1a0a8a3SJed Brown 
515b1a0a8a3SJed Brown EXTERN_C_BEGIN
516b1a0a8a3SJed Brown #undef __FUNCT__
517f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetLocalSubdomains_GASM"
518f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetLocalSubdomains_GASM(PC pc,PetscInt n,IS is[],IS is_local[])
519b1a0a8a3SJed Brown {
520f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
521b1a0a8a3SJed Brown   PetscErrorCode ierr;
522b1a0a8a3SJed Brown   PetscInt       i;
523b1a0a8a3SJed Brown 
524b1a0a8a3SJed Brown   PetscFunctionBegin;
525b1a0a8a3SJed Brown   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
526f746d493SDmitry Karpeev   if (pc->setupcalled && (n != osm->n || is)) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetLocalSubdomains() should be called before calling PCSetUp().");
527b1a0a8a3SJed Brown 
528b1a0a8a3SJed Brown   if (!pc->setupcalled) {
52993c1ec32SDmitry Karpeev     osm->n            = n;
53093c1ec32SDmitry Karpeev     osm->is           = 0;
53193c1ec32SDmitry Karpeev     osm->is_local     = 0;
532b1a0a8a3SJed Brown     if (is) {
533b1a0a8a3SJed Brown       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
534b1a0a8a3SJed Brown     }
535b1a0a8a3SJed Brown     if (is_local) {
536b1a0a8a3SJed Brown       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
537b1a0a8a3SJed Brown     }
538b1a0a8a3SJed Brown     if (osm->is) {
539ab3e923fSDmitry Karpeev       ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
540b1a0a8a3SJed Brown     }
541b1a0a8a3SJed Brown     if (is) {
542b1a0a8a3SJed Brown       ierr = PetscMalloc(n*sizeof(IS),&osm->is);CHKERRQ(ierr);
543b1a0a8a3SJed Brown       for (i=0; i<n; i++) { osm->is[i] = is[i]; }
544f746d493SDmitry Karpeev       /* Flag indicating that the user has set overlapping subdomains so PCGASM should not increase their size. */
545b1a0a8a3SJed Brown       osm->overlap = -1;
546b1a0a8a3SJed Brown     }
547b1a0a8a3SJed Brown     if (is_local) {
548b1a0a8a3SJed Brown       ierr = PetscMalloc(n*sizeof(IS),&osm->is_local);CHKERRQ(ierr);
549b1a0a8a3SJed Brown       for (i=0; i<n; i++) { osm->is_local[i] = is_local[i]; }
550b1a0a8a3SJed Brown     }
551b1a0a8a3SJed Brown   }
552b1a0a8a3SJed Brown   PetscFunctionReturn(0);
553b1a0a8a3SJed Brown }
554b1a0a8a3SJed Brown EXTERN_C_END
555b1a0a8a3SJed Brown 
556b1a0a8a3SJed Brown EXTERN_C_BEGIN
557b1a0a8a3SJed Brown #undef __FUNCT__
558f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains_GASM"
559ab3e923fSDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetTotalSubdomains_GASM(PC pc,PetscInt N) {
560f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
561b1a0a8a3SJed Brown   PetscErrorCode ierr;
562b1a0a8a3SJed Brown   PetscMPIInt    rank,size;
563b1a0a8a3SJed Brown   PetscInt       n;
564b1a0a8a3SJed Brown 
565b1a0a8a3SJed Brown   PetscFunctionBegin;
566b1a0a8a3SJed Brown   if (N < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
567b1a0a8a3SJed Brown 
568b1a0a8a3SJed Brown   /*
569b1a0a8a3SJed Brown      Split the subdomains equally among all processors
570b1a0a8a3SJed Brown   */
571b1a0a8a3SJed Brown   ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr);
572b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
573b1a0a8a3SJed Brown   n = N/size + ((N % size) > rank);
574b1a0a8a3SJed 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);
575f746d493SDmitry Karpeev   if (pc->setupcalled && n != osm->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before PCSetUp().");
576b1a0a8a3SJed Brown   if (!pc->setupcalled) {
577b1a0a8a3SJed Brown     if (osm->is) {
578ab3e923fSDmitry Karpeev       ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr);
579b1a0a8a3SJed Brown     }
58093c1ec32SDmitry Karpeev     osm->N            = N;
581f746d493SDmitry Karpeev     osm->n            = n;
582b1a0a8a3SJed Brown     osm->is           = 0;
583b1a0a8a3SJed Brown     osm->is_local     = 0;
584b1a0a8a3SJed Brown   }
585b1a0a8a3SJed Brown   PetscFunctionReturn(0);
586ab3e923fSDmitry Karpeev }/* PCGASMSetTotalSubdomains_GASM() */
587b1a0a8a3SJed Brown EXTERN_C_END
588b1a0a8a3SJed Brown 
589b1a0a8a3SJed Brown EXTERN_C_BEGIN
590b1a0a8a3SJed Brown #undef __FUNCT__
591f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM"
592f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetOverlap_GASM(PC pc,PetscInt ovl)
593b1a0a8a3SJed Brown {
594f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
595b1a0a8a3SJed Brown 
596b1a0a8a3SJed Brown   PetscFunctionBegin;
597b1a0a8a3SJed Brown   if (ovl < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
598f746d493SDmitry Karpeev   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp().");
599b1a0a8a3SJed Brown   if (!pc->setupcalled) {
600b1a0a8a3SJed Brown     osm->overlap = ovl;
601b1a0a8a3SJed Brown   }
602b1a0a8a3SJed Brown   PetscFunctionReturn(0);
603b1a0a8a3SJed Brown }
604b1a0a8a3SJed Brown EXTERN_C_END
605b1a0a8a3SJed Brown 
606b1a0a8a3SJed Brown EXTERN_C_BEGIN
607b1a0a8a3SJed Brown #undef __FUNCT__
608f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM"
609f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetType_GASM(PC pc,PCGASMType type)
610b1a0a8a3SJed Brown {
611f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
612b1a0a8a3SJed Brown 
613b1a0a8a3SJed Brown   PetscFunctionBegin;
614b1a0a8a3SJed Brown   osm->type     = type;
615b1a0a8a3SJed Brown   osm->type_set = PETSC_TRUE;
616b1a0a8a3SJed Brown   PetscFunctionReturn(0);
617b1a0a8a3SJed Brown }
618b1a0a8a3SJed Brown EXTERN_C_END
619b1a0a8a3SJed Brown 
620b1a0a8a3SJed Brown EXTERN_C_BEGIN
621b1a0a8a3SJed Brown #undef __FUNCT__
622f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM"
623f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetSortIndices_GASM(PC pc,PetscBool  doSort)
624b1a0a8a3SJed Brown {
625f746d493SDmitry Karpeev   PC_GASM *osm = (PC_GASM*)pc->data;
626b1a0a8a3SJed Brown 
627b1a0a8a3SJed Brown   PetscFunctionBegin;
628b1a0a8a3SJed Brown   osm->sort_indices = doSort;
629b1a0a8a3SJed Brown   PetscFunctionReturn(0);
630b1a0a8a3SJed Brown }
631b1a0a8a3SJed Brown EXTERN_C_END
632b1a0a8a3SJed Brown 
633b1a0a8a3SJed Brown EXTERN_C_BEGIN
634b1a0a8a3SJed Brown #undef __FUNCT__
635f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM"
63644fe18e5SDmitry Karpeev /*
63744fe18e5SDmitry Karpeev    FIX: This routine might need to be modified once multiple ranks per subdomain are allowed.
63844fe18e5SDmitry Karpeev         In particular, it would upset the global subdomain number calculation.
63944fe18e5SDmitry Karpeev */
640ab3e923fSDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp)
641b1a0a8a3SJed Brown {
642f746d493SDmitry Karpeev   PC_GASM         *osm = (PC_GASM*)pc->data;
643b1a0a8a3SJed Brown   PetscErrorCode ierr;
644b1a0a8a3SJed Brown 
645b1a0a8a3SJed Brown   PetscFunctionBegin;
646ab3e923fSDmitry 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");
647b1a0a8a3SJed Brown 
648ab3e923fSDmitry Karpeev   if (n) {
649ab3e923fSDmitry Karpeev     *n = osm->n;
650b1a0a8a3SJed Brown   }
651ab3e923fSDmitry Karpeev   if (first) {
652ab3e923fSDmitry Karpeev     ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr);
653ab3e923fSDmitry Karpeev     *first -= osm->n;
654b1a0a8a3SJed Brown   }
655b1a0a8a3SJed Brown   if (ksp) {
656b1a0a8a3SJed Brown     /* Assume that local solves are now different; not necessarily
657f746d493SDmitry Karpeev        true though!  This flag is used only for PCView_GASM() */
658b1a0a8a3SJed Brown     *ksp                   = osm->ksp;
659b1a0a8a3SJed Brown     osm->same_local_solves = PETSC_FALSE;
660b1a0a8a3SJed Brown   }
661b1a0a8a3SJed Brown   PetscFunctionReturn(0);
662ab3e923fSDmitry Karpeev }/* PCGASMGetSubKSP_GASM() */
663b1a0a8a3SJed Brown EXTERN_C_END
664b1a0a8a3SJed Brown 
665b1a0a8a3SJed Brown 
666b1a0a8a3SJed Brown #undef __FUNCT__
667f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetLocalSubdomains"
668b1a0a8a3SJed Brown /*@C
669f746d493SDmitry Karpeev     PCGASMSetLocalSubdomains - Sets the local subdomains (for this processor
670b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
671b1a0a8a3SJed Brown 
672b1a0a8a3SJed Brown     Collective on PC
673b1a0a8a3SJed Brown 
674b1a0a8a3SJed Brown     Input Parameters:
675b1a0a8a3SJed Brown +   pc - the preconditioner context
676b1a0a8a3SJed Brown .   n - the number of subdomains for this processor (default value = 1)
677b1a0a8a3SJed Brown .   is - the index set that defines the subdomains for this processor
678b1a0a8a3SJed Brown          (or PETSC_NULL for PETSc to determine subdomains)
679b1a0a8a3SJed Brown -   is_local - the index sets that define the local part of the subdomains for this processor
680b1a0a8a3SJed Brown          (or PETSC_NULL to use the default of 1 subdomain per process)
681b1a0a8a3SJed Brown 
682b1a0a8a3SJed Brown     Notes:
683b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
684b1a0a8a3SJed Brown 
685f746d493SDmitry Karpeev     By default the GASM preconditioner uses 1 block per processor.
686b1a0a8a3SJed Brown 
687f746d493SDmitry Karpeev     Use PCGASMSetTotalSubdomains() to set the subdomains for all processors.
688b1a0a8a3SJed Brown 
689b1a0a8a3SJed Brown     Level: advanced
690b1a0a8a3SJed Brown 
691f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz
692b1a0a8a3SJed Brown 
693f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
694f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetLocalSubdomains()
695b1a0a8a3SJed Brown @*/
696f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
697b1a0a8a3SJed Brown {
698b1a0a8a3SJed Brown   PetscErrorCode ierr;
699b1a0a8a3SJed Brown 
700b1a0a8a3SJed Brown   PetscFunctionBegin;
701b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
702f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr);
703b1a0a8a3SJed Brown   PetscFunctionReturn(0);
704b1a0a8a3SJed Brown }
705b1a0a8a3SJed Brown 
706b1a0a8a3SJed Brown #undef __FUNCT__
707f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains"
708b1a0a8a3SJed Brown /*@C
709f746d493SDmitry Karpeev     PCGASMSetTotalSubdomains - Sets the subdomains for all processor for the
710b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
711b1a0a8a3SJed Brown     PC communicator must call this routine, with the same index sets.
712b1a0a8a3SJed Brown 
713b1a0a8a3SJed Brown     Collective on PC
714b1a0a8a3SJed Brown 
715b1a0a8a3SJed Brown     Input Parameters:
716b1a0a8a3SJed Brown +   pc - the preconditioner context
717b1a0a8a3SJed Brown .   n - the number of subdomains for all processors
718b1a0a8a3SJed Brown .   is - the index sets that define the subdomains for all processor
719b1a0a8a3SJed Brown          (or PETSC_NULL for PETSc to determine subdomains)
720b1a0a8a3SJed Brown -   is_local - the index sets that define the local part of the subdomains for this processor
721b1a0a8a3SJed Brown          (or PETSC_NULL to use the default of 1 subdomain per process)
722b1a0a8a3SJed Brown 
723b1a0a8a3SJed Brown     Options Database Key:
724b1a0a8a3SJed Brown     To set the total number of subdomain blocks rather than specify the
725b1a0a8a3SJed Brown     index sets, use the option
726f746d493SDmitry Karpeev .    -pc_gasm_blocks <blks> - Sets total blocks
727b1a0a8a3SJed Brown 
728b1a0a8a3SJed Brown     Notes:
729b1a0a8a3SJed Brown     Currently you cannot use this to set the actual subdomains with the argument is.
730b1a0a8a3SJed Brown 
731f746d493SDmitry Karpeev     By default the GASM preconditioner uses 1 block per processor.
732b1a0a8a3SJed Brown 
733b1a0a8a3SJed Brown     These index sets cannot be destroyed until after completion of the
734f746d493SDmitry Karpeev     linear solves for which the GASM preconditioner is being used.
735b1a0a8a3SJed Brown 
736f746d493SDmitry Karpeev     Use PCGASMSetLocalSubdomains() to set local subdomains.
737b1a0a8a3SJed Brown 
738b1a0a8a3SJed Brown     Level: advanced
739b1a0a8a3SJed Brown 
740f746d493SDmitry Karpeev .keywords: PC, GASM, set, total, global, subdomains, additive Schwarz
741b1a0a8a3SJed Brown 
742f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
743f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
744b1a0a8a3SJed Brown @*/
745f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetTotalSubdomains(PC pc,PetscInt N)
746b1a0a8a3SJed Brown {
747b1a0a8a3SJed Brown   PetscErrorCode ierr;
748b1a0a8a3SJed Brown 
749b1a0a8a3SJed Brown   PetscFunctionBegin;
750b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
751f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetTotalSubdomains_C",(PC,PetscInt),(pc,N));CHKERRQ(ierr);
752b1a0a8a3SJed Brown   PetscFunctionReturn(0);
753b1a0a8a3SJed Brown }
754b1a0a8a3SJed Brown 
755b1a0a8a3SJed Brown #undef __FUNCT__
756f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap"
757b1a0a8a3SJed Brown /*@
758f746d493SDmitry Karpeev     PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the
759b1a0a8a3SJed Brown     additive Schwarz preconditioner.  Either all or no processors in the
760b1a0a8a3SJed Brown     PC communicator must call this routine.
761b1a0a8a3SJed Brown 
762b1a0a8a3SJed Brown     Logically Collective on PC
763b1a0a8a3SJed Brown 
764b1a0a8a3SJed Brown     Input Parameters:
765b1a0a8a3SJed Brown +   pc  - the preconditioner context
766b1a0a8a3SJed Brown -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
767b1a0a8a3SJed Brown 
768b1a0a8a3SJed Brown     Options Database Key:
769f746d493SDmitry Karpeev .   -pc_gasm_overlap <ovl> - Sets overlap
770b1a0a8a3SJed Brown 
771b1a0a8a3SJed Brown     Notes:
772f746d493SDmitry Karpeev     By default the GASM preconditioner uses 1 block per processor.  To use
773f746d493SDmitry Karpeev     multiple blocks per perocessor, see PCGASMSetTotalSubdomains() and
774f746d493SDmitry Karpeev     PCGASMSetLocalSubdomains() (and the option -pc_gasm_blocks <blks>).
775b1a0a8a3SJed Brown 
776b1a0a8a3SJed Brown     The overlap defaults to 1, so if one desires that no additional
777b1a0a8a3SJed Brown     overlap be computed beyond what may have been set with a call to
778f746d493SDmitry Karpeev     PCGASMSetTotalSubdomains() or PCGASMSetLocalSubdomains(), then ovl
779b1a0a8a3SJed Brown     must be set to be 0.  In particular, if one does not explicitly set
780b1a0a8a3SJed Brown     the subdomains an application code, then all overlap would be computed
781f746d493SDmitry Karpeev     internally by PETSc, and using an overlap of 0 would result in an GASM
782b1a0a8a3SJed Brown     variant that is equivalent to the block Jacobi preconditioner.
783b1a0a8a3SJed Brown 
784b1a0a8a3SJed Brown     Note that one can define initial index sets with any overlap via
785f746d493SDmitry Karpeev     PCGASMSetTotalSubdomains() or PCGASMSetLocalSubdomains(); the routine
786f746d493SDmitry Karpeev     PCGASMSetOverlap() merely allows PETSc to extend that overlap further
787b1a0a8a3SJed Brown     if desired.
788b1a0a8a3SJed Brown 
789b1a0a8a3SJed Brown     Level: intermediate
790b1a0a8a3SJed Brown 
791f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap
792b1a0a8a3SJed Brown 
793f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetLocalSubdomains(), PCGASMGetSubKSP(),
794f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMGetLocalSubdomains()
795b1a0a8a3SJed Brown @*/
796f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetOverlap(PC pc,PetscInt ovl)
797b1a0a8a3SJed Brown {
798b1a0a8a3SJed Brown   PetscErrorCode ierr;
799b1a0a8a3SJed Brown 
800b1a0a8a3SJed Brown   PetscFunctionBegin;
801b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
802b1a0a8a3SJed Brown   PetscValidLogicalCollectiveInt(pc,ovl,2);
803f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
804b1a0a8a3SJed Brown   PetscFunctionReturn(0);
805b1a0a8a3SJed Brown }
806b1a0a8a3SJed Brown 
807b1a0a8a3SJed Brown #undef __FUNCT__
808f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType"
809b1a0a8a3SJed Brown /*@
810f746d493SDmitry Karpeev     PCGASMSetType - Sets the type of restriction and interpolation used
811b1a0a8a3SJed Brown     for local problems in the additive Schwarz method.
812b1a0a8a3SJed Brown 
813b1a0a8a3SJed Brown     Logically Collective on PC
814b1a0a8a3SJed Brown 
815b1a0a8a3SJed Brown     Input Parameters:
816b1a0a8a3SJed Brown +   pc  - the preconditioner context
817f746d493SDmitry Karpeev -   type - variant of GASM, one of
818b1a0a8a3SJed Brown .vb
819f746d493SDmitry Karpeev       PC_GASM_BASIC       - full interpolation and restriction
820f746d493SDmitry Karpeev       PC_GASM_RESTRICT    - full restriction, local processor interpolation
821f746d493SDmitry Karpeev       PC_GASM_INTERPOLATE - full interpolation, local processor restriction
822f746d493SDmitry Karpeev       PC_GASM_NONE        - local processor restriction and interpolation
823b1a0a8a3SJed Brown .ve
824b1a0a8a3SJed Brown 
825b1a0a8a3SJed Brown     Options Database Key:
826f746d493SDmitry Karpeev .   -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
827b1a0a8a3SJed Brown 
828b1a0a8a3SJed Brown     Level: intermediate
829b1a0a8a3SJed Brown 
830f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
831b1a0a8a3SJed Brown 
832f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(),
833f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
834b1a0a8a3SJed Brown @*/
835f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetType(PC pc,PCGASMType type)
836b1a0a8a3SJed Brown {
837b1a0a8a3SJed Brown   PetscErrorCode ierr;
838b1a0a8a3SJed Brown 
839b1a0a8a3SJed Brown   PetscFunctionBegin;
840b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
841b1a0a8a3SJed Brown   PetscValidLogicalCollectiveEnum(pc,type,2);
842f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr);
843b1a0a8a3SJed Brown   PetscFunctionReturn(0);
844b1a0a8a3SJed Brown }
845b1a0a8a3SJed Brown 
846b1a0a8a3SJed Brown #undef __FUNCT__
847f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices"
848b1a0a8a3SJed Brown /*@
849f746d493SDmitry Karpeev     PCGASMSetSortIndices - Determines whether subdomain indices are sorted.
850b1a0a8a3SJed Brown 
851b1a0a8a3SJed Brown     Logically Collective on PC
852b1a0a8a3SJed Brown 
853b1a0a8a3SJed Brown     Input Parameters:
854b1a0a8a3SJed Brown +   pc  - the preconditioner context
855b1a0a8a3SJed Brown -   doSort - sort the subdomain indices
856b1a0a8a3SJed Brown 
857b1a0a8a3SJed Brown     Level: intermediate
858b1a0a8a3SJed Brown 
859f746d493SDmitry Karpeev .keywords: PC, GASM, set, type
860b1a0a8a3SJed Brown 
861f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(),
862f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D()
863b1a0a8a3SJed Brown @*/
864f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetSortIndices(PC pc,PetscBool  doSort)
865b1a0a8a3SJed Brown {
866b1a0a8a3SJed Brown   PetscErrorCode ierr;
867b1a0a8a3SJed Brown 
868b1a0a8a3SJed Brown   PetscFunctionBegin;
869b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
870b1a0a8a3SJed Brown   PetscValidLogicalCollectiveBool(pc,doSort,2);
871f746d493SDmitry Karpeev   ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
872b1a0a8a3SJed Brown   PetscFunctionReturn(0);
873b1a0a8a3SJed Brown }
874b1a0a8a3SJed Brown 
875b1a0a8a3SJed Brown #undef __FUNCT__
876f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP"
877b1a0a8a3SJed Brown /*@C
878f746d493SDmitry Karpeev    PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on
879b1a0a8a3SJed Brown    this processor.
880b1a0a8a3SJed Brown 
881b1a0a8a3SJed Brown    Collective on PC iff first_local is requested
882b1a0a8a3SJed Brown 
883b1a0a8a3SJed Brown    Input Parameter:
884b1a0a8a3SJed Brown .  pc - the preconditioner context
885b1a0a8a3SJed Brown 
886b1a0a8a3SJed Brown    Output Parameters:
887b1a0a8a3SJed Brown +  n_local - the number of blocks on this processor or PETSC_NULL
888b1a0a8a3SJed Brown .  first_local - the global number of the first block on this processor or PETSC_NULL,
889b1a0a8a3SJed Brown                  all processors must request or all must pass PETSC_NULL
890b1a0a8a3SJed Brown -  ksp - the array of KSP contexts
891b1a0a8a3SJed Brown 
892b1a0a8a3SJed Brown    Note:
893f746d493SDmitry Karpeev    After PCGASMGetSubKSP() the array of KSPes is not to be freed
894b1a0a8a3SJed Brown 
895b1a0a8a3SJed Brown    Currently for some matrix implementations only 1 block per processor
896b1a0a8a3SJed Brown    is supported.
897b1a0a8a3SJed Brown 
898f746d493SDmitry Karpeev    You must call KSPSetUp() before calling PCGASMGetSubKSP().
899b1a0a8a3SJed Brown 
900b1a0a8a3SJed Brown    Level: advanced
901b1a0a8a3SJed Brown 
902f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context
903b1a0a8a3SJed Brown 
904f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMSetOverlap(),
905f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(),
906b1a0a8a3SJed Brown @*/
907f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
908b1a0a8a3SJed Brown {
909b1a0a8a3SJed Brown   PetscErrorCode ierr;
910b1a0a8a3SJed Brown 
911b1a0a8a3SJed Brown   PetscFunctionBegin;
912b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
913f746d493SDmitry Karpeev   ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
914b1a0a8a3SJed Brown   PetscFunctionReturn(0);
915b1a0a8a3SJed Brown }
916b1a0a8a3SJed Brown 
917b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/
918b1a0a8a3SJed Brown /*MC
919f746d493SDmitry Karpeev    PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
920b1a0a8a3SJed Brown            its own KSP object.
921b1a0a8a3SJed Brown 
922b1a0a8a3SJed Brown    Options Database Keys:
923f746d493SDmitry Karpeev +  -pc_gasm_truelocal - Activates PCGGASMSetUseTrueLocal()
924f746d493SDmitry Karpeev .  -pc_gasm_blocks <blks> - Sets total blocks
925f746d493SDmitry Karpeev .  -pc_gasm_overlap <ovl> - Sets overlap
926f746d493SDmitry Karpeev -  -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type
927b1a0a8a3SJed Brown 
928b1a0a8a3SJed Brown      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
929f746d493SDmitry Karpeev       will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use
930f746d493SDmitry Karpeev       -pc_gasm_type basic to use the standard GASM.
931b1a0a8a3SJed Brown 
932b1a0a8a3SJed Brown    Notes: Each processor can have one or more blocks, but a block cannot be shared by more
933b1a0a8a3SJed Brown      than one processor. Defaults to one block per processor.
934b1a0a8a3SJed Brown 
935b1a0a8a3SJed Brown      To set options on the solvers for each block append -sub_ to all the KSP, and PC
936b1a0a8a3SJed Brown         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
937b1a0a8a3SJed Brown 
938f746d493SDmitry Karpeev      To set the options on the solvers separate for each block call PCGASMGetSubKSP()
939b1a0a8a3SJed Brown          and set the options directly on the resulting KSP object (you can access its PC
940b1a0a8a3SJed Brown          with KSPGetPC())
941b1a0a8a3SJed Brown 
942b1a0a8a3SJed Brown 
943b1a0a8a3SJed Brown    Level: beginner
944b1a0a8a3SJed Brown 
945b1a0a8a3SJed Brown    Concepts: additive Schwarz method
946b1a0a8a3SJed Brown 
947b1a0a8a3SJed Brown     References:
948b1a0a8a3SJed Brown     An additive variant of the Schwarz alternating method for the case of many subregions
949b1a0a8a3SJed Brown     M Dryja, OB Widlund - Courant Institute, New York University Technical report
950b1a0a8a3SJed Brown 
951b1a0a8a3SJed Brown     Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
952b1a0a8a3SJed Brown     Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X.
953b1a0a8a3SJed Brown 
954b1a0a8a3SJed Brown .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
955f746d493SDmitry Karpeev            PCBJACOBI, PCGASMSetUseTrueLocal(), PCGASMGetSubKSP(), PCGASMSetLocalSubdomains(),
956f746d493SDmitry Karpeev            PCGASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType()
957b1a0a8a3SJed Brown 
958b1a0a8a3SJed Brown M*/
959b1a0a8a3SJed Brown 
960b1a0a8a3SJed Brown EXTERN_C_BEGIN
961b1a0a8a3SJed Brown #undef __FUNCT__
962f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM"
963f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_GASM(PC pc)
964b1a0a8a3SJed Brown {
965b1a0a8a3SJed Brown   PetscErrorCode ierr;
966f746d493SDmitry Karpeev   PC_GASM         *osm;
967b1a0a8a3SJed Brown 
968b1a0a8a3SJed Brown   PetscFunctionBegin;
969f746d493SDmitry Karpeev   ierr = PetscNewLog(pc,PC_GASM,&osm);CHKERRQ(ierr);
970ab3e923fSDmitry Karpeev   osm->N                 = PETSC_DECIDE;
971ab3e923fSDmitry Karpeev   osm->n                 = 0;
972ab3e923fSDmitry Karpeev   osm->nmax              = 0;
973b1a0a8a3SJed Brown   osm->overlap           = 1;
974b1a0a8a3SJed Brown   osm->ksp               = 0;
975ab3e923fSDmitry Karpeev   osm->grestriction      = 0;
976ab3e923fSDmitry Karpeev   osm->gprolongation     = 0;
977ab3e923fSDmitry Karpeev   osm->gx                = 0;
978ab3e923fSDmitry Karpeev   osm->gy                = 0;
979b1a0a8a3SJed Brown   osm->x                 = 0;
980b1a0a8a3SJed Brown   osm->y                 = 0;
981b1a0a8a3SJed Brown   osm->is                = 0;
982b1a0a8a3SJed Brown   osm->is_local          = 0;
983b1a0a8a3SJed Brown   osm->pmat              = 0;
984f746d493SDmitry Karpeev   osm->type              = PC_GASM_RESTRICT;
985b1a0a8a3SJed Brown   osm->same_local_solves = PETSC_TRUE;
986b1a0a8a3SJed Brown   osm->sort_indices      = PETSC_TRUE;
987b1a0a8a3SJed Brown 
988b1a0a8a3SJed Brown   pc->data                   = (void*)osm;
989f746d493SDmitry Karpeev   pc->ops->apply             = PCApply_GASM;
990f746d493SDmitry Karpeev   pc->ops->applytranspose    = PCApplyTranspose_GASM;
991f746d493SDmitry Karpeev   pc->ops->setup             = PCSetUp_GASM;
992f746d493SDmitry Karpeev   pc->ops->destroy           = PCDestroy_GASM;
993f746d493SDmitry Karpeev   pc->ops->setfromoptions    = PCSetFromOptions_GASM;
994f746d493SDmitry Karpeev   pc->ops->setuponblocks     = PCSetUpOnBlocks_GASM;
995f746d493SDmitry Karpeev   pc->ops->view              = PCView_GASM;
996b1a0a8a3SJed Brown   pc->ops->applyrichardson   = 0;
997b1a0a8a3SJed Brown 
998f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetLocalSubdomains_C","PCGASMSetLocalSubdomains_GASM",
999f746d493SDmitry Karpeev                     PCGASMSetLocalSubdomains_GASM);CHKERRQ(ierr);
1000f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetTotalSubdomains_C","PCGASMSetTotalSubdomains_GASM",
1001f746d493SDmitry Karpeev                     PCGASMSetTotalSubdomains_GASM);CHKERRQ(ierr);
1002f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetOverlap_C","PCGASMSetOverlap_GASM",
1003f746d493SDmitry Karpeev                     PCGASMSetOverlap_GASM);CHKERRQ(ierr);
1004f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetType_C","PCGASMSetType_GASM",
1005f746d493SDmitry Karpeev                     PCGASMSetType_GASM);CHKERRQ(ierr);
1006f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSortIndices_C","PCGASMSetSortIndices_GASM",
1007f746d493SDmitry Karpeev                     PCGASMSetSortIndices_GASM);CHKERRQ(ierr);
1008f746d493SDmitry Karpeev   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMGetSubKSP_C","PCGASMGetSubKSP_GASM",
1009f746d493SDmitry Karpeev                     PCGASMGetSubKSP_GASM);CHKERRQ(ierr);
1010b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1011b1a0a8a3SJed Brown }
1012b1a0a8a3SJed Brown EXTERN_C_END
1013b1a0a8a3SJed Brown 
1014b1a0a8a3SJed Brown 
1015b1a0a8a3SJed Brown #undef __FUNCT__
1016f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains"
1017b1a0a8a3SJed Brown /*@C
1018f746d493SDmitry Karpeev    PCGASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
1019b1a0a8a3SJed Brown    preconditioner for a any problem on a general grid.
1020b1a0a8a3SJed Brown 
1021b1a0a8a3SJed Brown    Collective
1022b1a0a8a3SJed Brown 
1023b1a0a8a3SJed Brown    Input Parameters:
1024b1a0a8a3SJed Brown +  A - The global matrix operator
1025b1a0a8a3SJed Brown -  n - the number of local blocks
1026b1a0a8a3SJed Brown 
1027b1a0a8a3SJed Brown    Output Parameters:
1028b1a0a8a3SJed Brown .  outis - the array of index sets defining the subdomains
1029b1a0a8a3SJed Brown 
1030b1a0a8a3SJed Brown    Level: advanced
1031b1a0a8a3SJed Brown 
1032f746d493SDmitry Karpeev    Note: this generates nonoverlapping subdomains; the PCGASM will generate the overlap
1033f746d493SDmitry Karpeev     from these if you use PCGASMSetLocalSubdomains()
1034b1a0a8a3SJed Brown 
1035b1a0a8a3SJed Brown     In the Fortran version you must provide the array outis[] already allocated of length n.
1036b1a0a8a3SJed Brown 
1037f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1038b1a0a8a3SJed Brown 
1039f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMDestroySubdomains()
1040b1a0a8a3SJed Brown @*/
1041f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
1042b1a0a8a3SJed Brown {
1043b1a0a8a3SJed Brown   MatPartitioning           mpart;
1044b1a0a8a3SJed Brown   const char                *prefix;
1045b1a0a8a3SJed Brown   PetscErrorCode            (*f)(Mat,PetscBool *,MatReuse,Mat*);
1046b1a0a8a3SJed Brown   PetscMPIInt               size;
1047b1a0a8a3SJed Brown   PetscInt                  i,j,rstart,rend,bs;
1048b1a0a8a3SJed Brown   PetscBool                 iscopy = PETSC_FALSE,isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
1049b1a0a8a3SJed Brown   Mat                       Ad = PETSC_NULL, adj;
1050b1a0a8a3SJed Brown   IS                        ispart,isnumb,*is;
1051b1a0a8a3SJed Brown   PetscErrorCode            ierr;
1052b1a0a8a3SJed Brown 
1053b1a0a8a3SJed Brown   PetscFunctionBegin;
1054b1a0a8a3SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1055b1a0a8a3SJed Brown   PetscValidPointer(outis,3);
1056b1a0a8a3SJed Brown   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
1057b1a0a8a3SJed Brown 
1058b1a0a8a3SJed Brown   /* Get prefix, row distribution, and block size */
1059b1a0a8a3SJed Brown   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1060b1a0a8a3SJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1061b1a0a8a3SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1062b1a0a8a3SJed 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);
1063b1a0a8a3SJed Brown 
1064b1a0a8a3SJed Brown   /* Get diagonal block from matrix if possible */
1065b1a0a8a3SJed Brown   ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);
1066b1a0a8a3SJed Brown   ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr);
1067b1a0a8a3SJed Brown   if (f) {
1068b1a0a8a3SJed Brown     ierr = (*f)(A,&iscopy,MAT_INITIAL_MATRIX,&Ad);CHKERRQ(ierr);
1069b1a0a8a3SJed Brown   } else if (size == 1) {
1070b1a0a8a3SJed Brown     iscopy = PETSC_FALSE; Ad = A;
1071b1a0a8a3SJed Brown   } else {
1072b1a0a8a3SJed Brown     iscopy = PETSC_FALSE; Ad = PETSC_NULL;
1073b1a0a8a3SJed Brown   }
1074b1a0a8a3SJed Brown   if (Ad) {
1075b1a0a8a3SJed Brown     ierr = PetscTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1076b1a0a8a3SJed Brown     if (!isbaij) {ierr = PetscTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
1077b1a0a8a3SJed Brown   }
1078b1a0a8a3SJed Brown   if (Ad && n > 1) {
1079b1a0a8a3SJed Brown     PetscBool  match,done;
1080b1a0a8a3SJed Brown     /* Try to setup a good matrix partitioning if available */
1081b1a0a8a3SJed Brown     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
1082b1a0a8a3SJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
1083b1a0a8a3SJed Brown     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1084b1a0a8a3SJed Brown     ierr = PetscTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
1085b1a0a8a3SJed Brown     if (!match) {
1086b1a0a8a3SJed Brown       ierr = PetscTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
1087b1a0a8a3SJed Brown     }
1088b1a0a8a3SJed Brown     if (!match) { /* assume a "good" partitioner is available */
1089b1a0a8a3SJed Brown       PetscInt na,*ia,*ja;
1090b1a0a8a3SJed Brown       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1091b1a0a8a3SJed Brown       if (done) {
1092b1a0a8a3SJed Brown 	/* Build adjacency matrix by hand. Unfortunately a call to
1093b1a0a8a3SJed Brown 	   MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1094b1a0a8a3SJed Brown 	   remove the block-aij structure and we cannot expect
1095b1a0a8a3SJed Brown 	   MatPartitioning to split vertices as we need */
1096b1a0a8a3SJed Brown 	PetscInt i,j,*row,len,nnz,cnt,*iia=0,*jja=0;
1097b1a0a8a3SJed Brown 	nnz = 0;
1098b1a0a8a3SJed Brown 	for (i=0; i<na; i++) { /* count number of nonzeros */
1099b1a0a8a3SJed Brown 	  len = ia[i+1] - ia[i];
1100b1a0a8a3SJed Brown 	  row = ja + ia[i];
1101b1a0a8a3SJed Brown 	  for (j=0; j<len; j++) {
1102b1a0a8a3SJed Brown 	    if (row[j] == i) { /* don't count diagonal */
1103b1a0a8a3SJed Brown 	      len--; break;
1104b1a0a8a3SJed Brown 	    }
1105b1a0a8a3SJed Brown 	  }
1106b1a0a8a3SJed Brown 	  nnz += len;
1107b1a0a8a3SJed Brown 	}
1108b1a0a8a3SJed Brown 	ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr);
1109b1a0a8a3SJed Brown 	ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr);
1110b1a0a8a3SJed Brown 	nnz    = 0;
1111b1a0a8a3SJed Brown 	iia[0] = 0;
1112b1a0a8a3SJed Brown 	for (i=0; i<na; i++) { /* fill adjacency */
1113b1a0a8a3SJed Brown 	  cnt = 0;
1114b1a0a8a3SJed Brown 	  len = ia[i+1] - ia[i];
1115b1a0a8a3SJed Brown 	  row = ja + ia[i];
1116b1a0a8a3SJed Brown 	  for (j=0; j<len; j++) {
1117b1a0a8a3SJed Brown 	    if (row[j] != i) { /* if not diagonal */
1118b1a0a8a3SJed Brown 	      jja[nnz+cnt++] = row[j];
1119b1a0a8a3SJed Brown 	    }
1120b1a0a8a3SJed Brown 	  }
1121b1a0a8a3SJed Brown 	  nnz += cnt;
1122b1a0a8a3SJed Brown 	  iia[i+1] = nnz;
1123b1a0a8a3SJed Brown 	}
1124b1a0a8a3SJed Brown 	/* Partitioning of the adjacency matrix */
1125b1a0a8a3SJed Brown 	ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);CHKERRQ(ierr);
1126b1a0a8a3SJed Brown 	ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
1127b1a0a8a3SJed Brown 	ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
1128b1a0a8a3SJed Brown 	ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
1129b1a0a8a3SJed Brown 	ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
1130b1a0a8a3SJed Brown 	ierr = MatDestroy(adj);CHKERRQ(ierr);
1131b1a0a8a3SJed Brown 	foundpart = PETSC_TRUE;
1132b1a0a8a3SJed Brown       }
1133b1a0a8a3SJed Brown       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1134b1a0a8a3SJed Brown     }
1135b1a0a8a3SJed Brown     ierr = MatPartitioningDestroy(mpart);CHKERRQ(ierr);
1136b1a0a8a3SJed Brown   }
1137b1a0a8a3SJed Brown   if (iscopy) {ierr = MatDestroy(Ad);CHKERRQ(ierr);}
1138b1a0a8a3SJed Brown 
1139b1a0a8a3SJed Brown   ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr);
1140b1a0a8a3SJed Brown   *outis = is;
1141b1a0a8a3SJed Brown 
1142b1a0a8a3SJed Brown   if (!foundpart) {
1143b1a0a8a3SJed Brown 
1144b1a0a8a3SJed Brown     /* Partitioning by contiguous chunks of rows */
1145b1a0a8a3SJed Brown 
1146b1a0a8a3SJed Brown     PetscInt mbs   = (rend-rstart)/bs;
1147b1a0a8a3SJed Brown     PetscInt start = rstart;
1148b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1149b1a0a8a3SJed Brown       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
1150b1a0a8a3SJed Brown       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
1151b1a0a8a3SJed Brown       start += count;
1152b1a0a8a3SJed Brown     }
1153b1a0a8a3SJed Brown 
1154b1a0a8a3SJed Brown   } else {
1155b1a0a8a3SJed Brown 
1156b1a0a8a3SJed Brown     /* Partitioning by adjacency of diagonal block  */
1157b1a0a8a3SJed Brown 
1158b1a0a8a3SJed Brown     const PetscInt *numbering;
1159b1a0a8a3SJed Brown     PetscInt       *count,nidx,*indices,*newidx,start=0;
1160b1a0a8a3SJed Brown     /* Get node count in each partition */
1161b1a0a8a3SJed Brown     ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr);
1162b1a0a8a3SJed Brown     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
1163b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1164b1a0a8a3SJed Brown       for (i=0; i<n; i++) count[i] *= bs;
1165b1a0a8a3SJed Brown     }
1166b1a0a8a3SJed Brown     /* Build indices from node numbering */
1167b1a0a8a3SJed Brown     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1168b1a0a8a3SJed Brown     ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr);
1169b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1170b1a0a8a3SJed Brown     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
1171b1a0a8a3SJed Brown     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
1172b1a0a8a3SJed Brown     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
1173b1a0a8a3SJed Brown     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1174b1a0a8a3SJed Brown       ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr);
1175b1a0a8a3SJed Brown       for (i=0; i<nidx; i++)
1176b1a0a8a3SJed Brown 	for (j=0; j<bs; j++)
1177b1a0a8a3SJed Brown 	  newidx[i*bs+j] = indices[i]*bs + j;
1178b1a0a8a3SJed Brown       ierr = PetscFree(indices);CHKERRQ(ierr);
1179b1a0a8a3SJed Brown       nidx   *= bs;
1180b1a0a8a3SJed Brown       indices = newidx;
1181b1a0a8a3SJed Brown     }
1182b1a0a8a3SJed Brown     /* Shift to get global indices */
1183b1a0a8a3SJed Brown     for (i=0; i<nidx; i++) indices[i] += rstart;
1184b1a0a8a3SJed Brown 
1185b1a0a8a3SJed Brown     /* Build the index sets for each block */
1186b1a0a8a3SJed Brown     for (i=0; i<n; i++) {
1187b1a0a8a3SJed Brown       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
1188b1a0a8a3SJed Brown       ierr   = ISSort(is[i]);CHKERRQ(ierr);
1189b1a0a8a3SJed Brown       start += count[i];
1190b1a0a8a3SJed Brown     }
1191b1a0a8a3SJed Brown 
1192b1a0a8a3SJed Brown     ierr = PetscFree(count);
1193b1a0a8a3SJed Brown     ierr = PetscFree(indices);
1194b1a0a8a3SJed Brown     ierr = ISDestroy(isnumb);CHKERRQ(ierr);
1195b1a0a8a3SJed Brown     ierr = ISDestroy(ispart);CHKERRQ(ierr);
1196b1a0a8a3SJed Brown 
1197b1a0a8a3SJed Brown   }
1198b1a0a8a3SJed Brown 
1199b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1200b1a0a8a3SJed Brown }
1201b1a0a8a3SJed Brown 
1202b1a0a8a3SJed Brown #undef __FUNCT__
1203f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains"
1204b1a0a8a3SJed Brown /*@C
1205f746d493SDmitry Karpeev    PCGASMDestroySubdomains - Destroys the index sets created with
1206f746d493SDmitry Karpeev    PCGASMCreateSubdomains(). Should be called after setting subdomains
1207f746d493SDmitry Karpeev    with PCGASMSetLocalSubdomains().
1208b1a0a8a3SJed Brown 
1209b1a0a8a3SJed Brown    Collective
1210b1a0a8a3SJed Brown 
1211b1a0a8a3SJed Brown    Input Parameters:
1212b1a0a8a3SJed Brown +  n - the number of index sets
1213b1a0a8a3SJed Brown .  is - the array of index sets
1214b1a0a8a3SJed Brown -  is_local - the array of local index sets, can be PETSC_NULL
1215b1a0a8a3SJed Brown 
1216b1a0a8a3SJed Brown    Level: advanced
1217b1a0a8a3SJed Brown 
1218f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid
1219b1a0a8a3SJed Brown 
1220f746d493SDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetLocalSubdomains()
1221b1a0a8a3SJed Brown @*/
1222f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
1223b1a0a8a3SJed Brown {
1224b1a0a8a3SJed Brown   PetscInt       i;
1225b1a0a8a3SJed Brown   PetscErrorCode ierr;
1226b1a0a8a3SJed Brown   PetscFunctionBegin;
1227b1a0a8a3SJed Brown   if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"n must be > 0: n = %D",n);
1228b1a0a8a3SJed Brown   PetscValidPointer(is,2);
1229b1a0a8a3SJed Brown   for (i=0; i<n; i++) { ierr = ISDestroy(is[i]);CHKERRQ(ierr); }
1230b1a0a8a3SJed Brown   ierr = PetscFree(is);CHKERRQ(ierr);
1231b1a0a8a3SJed Brown   if (is_local) {
1232b1a0a8a3SJed Brown     PetscValidPointer(is_local,3);
1233b1a0a8a3SJed Brown     for (i=0; i<n; i++) { ierr = ISDestroy(is_local[i]);CHKERRQ(ierr); }
1234b1a0a8a3SJed Brown     ierr = PetscFree(is_local);CHKERRQ(ierr);
1235b1a0a8a3SJed Brown   }
1236b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1237b1a0a8a3SJed Brown }
1238b1a0a8a3SJed Brown 
1239b1a0a8a3SJed Brown #undef __FUNCT__
1240f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D"
1241b1a0a8a3SJed Brown /*@
1242f746d493SDmitry Karpeev    PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1243b1a0a8a3SJed Brown    preconditioner for a two-dimensional problem on a regular grid.
1244b1a0a8a3SJed Brown 
1245b1a0a8a3SJed Brown    Not Collective
1246b1a0a8a3SJed Brown 
1247b1a0a8a3SJed Brown    Input Parameters:
1248b1a0a8a3SJed Brown +  m, n - the number of mesh points in the x and y directions
1249b1a0a8a3SJed Brown .  M, N - the number of subdomains in the x and y directions
1250b1a0a8a3SJed Brown .  dof - degrees of freedom per node
1251b1a0a8a3SJed Brown -  overlap - overlap in mesh lines
1252b1a0a8a3SJed Brown 
1253b1a0a8a3SJed Brown    Output Parameters:
1254b1a0a8a3SJed Brown +  Nsub - the number of subdomains created
1255b1a0a8a3SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
1256b1a0a8a3SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
1257b1a0a8a3SJed Brown 
1258b1a0a8a3SJed Brown    Note:
1259b1a0a8a3SJed Brown    Presently PCAMSCreateSubdomains2d() is valid only for sequential
1260b1a0a8a3SJed Brown    preconditioners.  More general related routines are
1261f746d493SDmitry Karpeev    PCGASMSetTotalSubdomains() and PCGASMSetLocalSubdomains().
1262b1a0a8a3SJed Brown 
1263b1a0a8a3SJed Brown    Level: advanced
1264b1a0a8a3SJed Brown 
1265f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid
1266b1a0a8a3SJed Brown 
1267f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetLocalSubdomains(), PCGASMGetSubKSP(),
1268f746d493SDmitry Karpeev           PCGASMSetOverlap()
1269b1a0a8a3SJed Brown @*/
1270*2bbb417fSDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMCreateSubdomains2D(PC pc, PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local)
1271b1a0a8a3SJed Brown {
1272b1a0a8a3SJed Brown   PetscErrorCode ierr;
1273*2bbb417fSDmitry Karpeev   PetscMPIInt size, rank;
1274*2bbb417fSDmitry Karpeev   PetscInt i, j;
1275*2bbb417fSDmitry Karpeev   PetscInt maxheight, maxwidth;
1276*2bbb417fSDmitry Karpeev   PetscInt xstart, xleft, xright, xleft_loc, xright_loc;
1277*2bbb417fSDmitry Karpeev   PetscInt ystart, ylow,  yhigh,  ylow_loc,  yhigh_loc;
1278*2bbb417fSDmitry Karpeev   PetscInt x[2][2], y[2][2];
1279*2bbb417fSDmitry Karpeev   PetscInt first, last;
1280*2bbb417fSDmitry Karpeev   PetscInt nidx, *idx;
1281*2bbb417fSDmitry Karpeev   PetscInt ii,jj,s,q,d;
1282*2bbb417fSDmitry Karpeev   PetscMPIInt color;
1283*2bbb417fSDmitry Karpeev   MPI_Comm comm, subcomm;
1284*2bbb417fSDmitry Karpeev   IS **iis;
1285b1a0a8a3SJed Brown   PetscFunctionBegin;
1286*2bbb417fSDmitry Karpeev   ierr = PetscObjectGetComm((PetscObject)pc, &comm);     CHKERRQ(ierr);
1287*2bbb417fSDmitry Karpeev   ierr = MPI_Comm_size(comm, &size);                     CHKERRQ(ierr);
1288*2bbb417fSDmitry Karpeev   ierr = MPI_Comm_rank(comm, &rank);                     CHKERRQ(ierr);
1289*2bbb417fSDmitry Karpeev   ierr = MatGetOwnershipRange(pc->pmat, &first, &last);  CHKERRQ(ierr);
1290*2bbb417fSDmitry Karpeev   if(first%dof || last%dof) {
1291*2bbb417fSDmitry Karpeev     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Matrix row partitioning unsuitable for domain decomposition: local row range (%D,%D) "
1292*2bbb417fSDmitry Karpeev 	     "does not respect the number of degrees of freedom per grid point %D", first, last, dof);
1293*2bbb417fSDmitry Karpeev   }
1294*2bbb417fSDmitry Karpeev   *Nsub = M*N;
1295b1a0a8a3SJed Brown   ierr = PetscMalloc((*Nsub)*sizeof(IS*),is);         CHKERRQ(ierr);
1296b1a0a8a3SJed Brown   ierr = PetscMalloc((*Nsub)*sizeof(IS*),is_local);   CHKERRQ(ierr);
1297b1a0a8a3SJed Brown   ystart = 0;
1298*2bbb417fSDmitry Karpeev   s = 0;
1299*2bbb417fSDmitry Karpeev   for (j=0; j<N; j++) {
1300*2bbb417fSDmitry Karpeev     maxheight = n/N + ((n % N) > j); /* Maximal height of subdomain */
1301*2bbb417fSDmitry Karpeev     if (maxheight < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains for mesh height %D", N, n);
1302*2bbb417fSDmitry Karpeev     /* Vertical domain limits with an overlap. */
1303*2bbb417fSDmitry Karpeev     y[0][1] = PetscMin(ystart + maxheight + overlap,n);
1304*2bbb417fSDmitry Karpeev     y[0][0] = PetscMax(ystart - overlap,0);
1305*2bbb417fSDmitry Karpeev     /* Vertical domain limits without an overlap. */
1306*2bbb417fSDmitry Karpeev     y[1][0] = ystart;
1307*2bbb417fSDmitry Karpeev     y[1][1] = PetscMin(ystart + maxheight,n);
1308b1a0a8a3SJed Brown     xstart = 0;
1309*2bbb417fSDmitry Karpeev     for (i=0; i<M; i++) {
1310*2bbb417fSDmitry Karpeev       maxwidth = m/M + ((m % M) > i); /* Maximal width of subdomain */
1311*2bbb417fSDmitry Karpeev       if (maxwidth < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains for mesh width %D", M, m);
1312*2bbb417fSDmitry Karpeev       /* Horizontal domain limits with an overlap. */
1313*2bbb417fSDmitry Karpeev       x[0][0]  = PetscMax(xstart - overlap,0);
1314*2bbb417fSDmitry Karpeev       x[0][1]  = PetscMin(xstart + maxwidth + overlap,m);
1315*2bbb417fSDmitry Karpeev       /* Horizontal domain limits without an overlap. */
1316*2bbb417fSDmitry Karpeev       x[1][0] = xstart;
1317*2bbb417fSDmitry Karpeev       x[1][1] = PetscMin(xstart+maxwidth,m);
1318*2bbb417fSDmitry Karpeev       /*
1319*2bbb417fSDmitry Karpeev 	 Determine whether this domain intersects this processor's ownership range of pc->pmat.
1320*2bbb417fSDmitry Karpeev 	 Do this twice: first for the domains with overlaps, and once without.
1321*2bbb417fSDmitry Karpeev 	 During the first pass create the subcommunicators, and use them on the second pass as well.
1322*2bbb417fSDmitry Karpeev       */
1323*2bbb417fSDmitry Karpeev       /* ylow_loc is the grid row containing first */
1324*2bbb417fSDmitry Karpeev       ylow_loc = first/m;
1325*2bbb417fSDmitry Karpeev       /* xleft_loc is the offset of first within the yylow-th row */
1326*2bbb417fSDmitry Karpeev       xleft_loc = first%m;
1327*2bbb417fSDmitry Karpeev       /* yhigh_loc is the grid row above the one containing last */
1328*2bbb417fSDmitry Karpeev       yhigh_loc = last/m + 1;
1329*2bbb417fSDmitry Karpeev       /* xright is one above the offset of last within the yhigh_loc-th row */
1330*2bbb417fSDmitry Karpeev       xright_loc = last%m + 1;
1331*2bbb417fSDmitry Karpeev       for(q = 0; q < 2; ++q) {
1332*2bbb417fSDmitry Karpeev 	/*
1333*2bbb417fSDmitry Karpeev 	  domain limits, (xleft, xright) and (ylow, yheigh) are adjusted
1334*2bbb417fSDmitry Karpeev 	  according to whether the domain with an overlap or without is considered.
1335*2bbb417fSDmitry Karpeev 	*/
1336*2bbb417fSDmitry Karpeev 	xleft = x[q][0]; xright = x[q][1];
1337*2bbb417fSDmitry Karpeev 	ylow  = y[q][0]; yhigh  = y[q][1];
1338*2bbb417fSDmitry Karpeev 	/*
1339*2bbb417fSDmitry Karpeev 	  (ylow_loc,xleft) and (yhigh_loc,xright) are the corners of the bounding box of
1340*2bbb417fSDmitry Karpeev 	  the intersection of the domain with the local ownership range;
1341*2bbb417fSDmitry Karpeev 	  xleft_loc and xright_loc are the lower and upper bounds in the first and last rows
1342*2bbb417fSDmitry Karpeev 	  of the intersection;
1343*2bbb417fSDmitry Karpeev 	  nidx is the size of the intersection.
1344*2bbb417fSDmitry Karpeev 	*/
1345*2bbb417fSDmitry Karpeev 	nidx = 0;
1346*2bbb417fSDmitry Karpeev 	if(ylow < yhigh_loc && yhigh > ylow_loc) {
1347*2bbb417fSDmitry Karpeev 	  PetscInt width = xright-xleft;
1348*2bbb417fSDmitry Karpeev 	  nidx += width*PetscMax(yhigh_loc-ylow_loc-1,0);
1349*2bbb417fSDmitry Karpeev 	  nidx += PetscMin(xright_loc-xright,width);
1350*2bbb417fSDmitry Karpeev 	  nidx -= PetscMax(xleft_loc-xleft, width);
1351b1a0a8a3SJed Brown 	}
1352*2bbb417fSDmitry Karpeev 	if(q == 0) {
1353*2bbb417fSDmitry Karpeev 	  iis = is;
1354*2bbb417fSDmitry Karpeev 	  if(nidx) {
1355*2bbb417fSDmitry Karpeev 	    color = 1;
1356b1a0a8a3SJed Brown 	  }
1357*2bbb417fSDmitry Karpeev 	  else {
1358*2bbb417fSDmitry Karpeev 	    color = MPI_UNDEFINED;
1359*2bbb417fSDmitry Karpeev 	  }
1360*2bbb417fSDmitry Karpeev 	  ierr = MPI_Comm_split(comm, color, rank, &subcomm); CHKERRQ(ierr);
1361*2bbb417fSDmitry Karpeev 	}
1362*2bbb417fSDmitry Karpeev 	else {
1363b1a0a8a3SJed Brown 	  if(overlap == 0) {
1364*2bbb417fSDmitry Karpeev 	    ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr);
1365*2bbb417fSDmitry Karpeev 	    (*is_local)[s] = (*is)[s];
1366*2bbb417fSDmitry Karpeev 	    continue;
1367*2bbb417fSDmitry Karpeev 	  }
1368*2bbb417fSDmitry Karpeev 	  iis = is_local;
1369*2bbb417fSDmitry Karpeev 	  subcomm = ((PetscObject)(*is)[s])->comm;
1370*2bbb417fSDmitry Karpeev 	}
1371*2bbb417fSDmitry Karpeev 	if(nidx) {
1372*2bbb417fSDmitry Karpeev 	  PetscInt k,kk;
1373*2bbb417fSDmitry Karpeev 	  ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr);
1374*2bbb417fSDmitry Karpeev 	  k    = 0;
1375*2bbb417fSDmitry Karpeev 	  for (jj=ylow_loc; jj<yhigh_loc; ++jj) {
1376*2bbb417fSDmitry Karpeev 	    PetscInt x0 = (jj==ylow_loc)?PetscMax(xleft_loc,xleft):xleft;
1377*2bbb417fSDmitry Karpeev 	    PetscInt x1 = (jj==yhigh_loc)?PetscMin(xright_loc,xright):xright;
1378*2bbb417fSDmitry Karpeev 	    kk = dof*(m*jj + x0);
1379*2bbb417fSDmitry Karpeev 	    for (ii=x0; ii<x1; ++ii) {
1380*2bbb417fSDmitry Karpeev 	      for(d = 0; d < dof; ++d) {
1381*2bbb417fSDmitry Karpeev 		idx[k++] = kk++;
1382b1a0a8a3SJed Brown 	      }
1383b1a0a8a3SJed Brown 	    }
1384b1a0a8a3SJed Brown 	  }
1385*2bbb417fSDmitry Karpeev 	  ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*iis)+s);CHKERRQ(ierr);
1386b1a0a8a3SJed Brown 	}
1387*2bbb417fSDmitry Karpeev       }/* for(q = 0; q < 2; ++q) */
1388*2bbb417fSDmitry Karpeev       xstart += maxwidth;
1389*2bbb417fSDmitry Karpeev       ++s;
1390b1a0a8a3SJed Brown     }
1391*2bbb417fSDmitry Karpeev     ystart += maxheight;
1392*2bbb417fSDmitry Karpeev   }
1393b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1394b1a0a8a3SJed Brown }
1395b1a0a8a3SJed Brown 
1396b1a0a8a3SJed Brown #undef __FUNCT__
1397f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetLocalSubdomains"
1398b1a0a8a3SJed Brown /*@C
1399f746d493SDmitry Karpeev     PCGASMGetLocalSubdomains - Gets the local subdomains (for this processor
1400b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1401b1a0a8a3SJed Brown 
1402b1a0a8a3SJed Brown     Not Collective
1403b1a0a8a3SJed Brown 
1404b1a0a8a3SJed Brown     Input Parameter:
1405b1a0a8a3SJed Brown .   pc - the preconditioner context
1406b1a0a8a3SJed Brown 
1407b1a0a8a3SJed Brown     Output Parameters:
1408b1a0a8a3SJed Brown +   n - the number of subdomains for this processor (default value = 1)
1409b1a0a8a3SJed Brown .   is - the index sets that define the subdomains for this processor
1410b1a0a8a3SJed Brown -   is_local - the index sets that define the local part of the subdomains for this processor (can be PETSC_NULL)
1411b1a0a8a3SJed Brown 
1412b1a0a8a3SJed Brown 
1413b1a0a8a3SJed Brown     Notes:
1414b1a0a8a3SJed Brown     The IS numbering is in the parallel, global numbering of the vector.
1415b1a0a8a3SJed Brown 
1416b1a0a8a3SJed Brown     Level: advanced
1417b1a0a8a3SJed Brown 
1418f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz
1419b1a0a8a3SJed Brown 
1420f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
1421f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetLocalSubdomains(), PCGASMGetLocalSubmatrices()
1422b1a0a8a3SJed Brown @*/
1423f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
1424b1a0a8a3SJed Brown {
1425f746d493SDmitry Karpeev   PC_GASM         *osm;
1426b1a0a8a3SJed Brown   PetscErrorCode ierr;
1427b1a0a8a3SJed Brown   PetscBool      match;
1428b1a0a8a3SJed Brown 
1429b1a0a8a3SJed Brown   PetscFunctionBegin;
1430b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1431b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1432b1a0a8a3SJed Brown   if (is) PetscValidPointer(is,3);
1433f746d493SDmitry Karpeev   ierr = PetscTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1434b1a0a8a3SJed Brown   if (!match) {
1435b1a0a8a3SJed Brown     if (n)  *n  = 0;
1436b1a0a8a3SJed Brown     if (is) *is = PETSC_NULL;
1437b1a0a8a3SJed Brown   } else {
1438f746d493SDmitry Karpeev     osm = (PC_GASM*)pc->data;
1439ab3e923fSDmitry Karpeev     if (n)  *n  = osm->n;
1440b1a0a8a3SJed Brown     if (is) *is = osm->is;
1441b1a0a8a3SJed Brown     if (is_local) *is_local = osm->is_local;
1442b1a0a8a3SJed Brown   }
1443b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1444b1a0a8a3SJed Brown }
1445b1a0a8a3SJed Brown 
1446b1a0a8a3SJed Brown #undef __FUNCT__
1447f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetLocalSubmatrices"
1448b1a0a8a3SJed Brown /*@C
1449f746d493SDmitry Karpeev     PCGASMGetLocalSubmatrices - Gets the local submatrices (for this processor
1450b1a0a8a3SJed Brown     only) for the additive Schwarz preconditioner.
1451b1a0a8a3SJed Brown 
1452b1a0a8a3SJed Brown     Not Collective
1453b1a0a8a3SJed Brown 
1454b1a0a8a3SJed Brown     Input Parameter:
1455b1a0a8a3SJed Brown .   pc - the preconditioner context
1456b1a0a8a3SJed Brown 
1457b1a0a8a3SJed Brown     Output Parameters:
1458b1a0a8a3SJed Brown +   n - the number of matrices for this processor (default value = 1)
1459b1a0a8a3SJed Brown -   mat - the matrices
1460b1a0a8a3SJed Brown 
1461b1a0a8a3SJed Brown 
1462b1a0a8a3SJed Brown     Level: advanced
1463b1a0a8a3SJed Brown 
1464f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi
1465b1a0a8a3SJed Brown 
1466f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(),
1467f746d493SDmitry Karpeev           PCGASMCreateSubdomains2D(), PCGASMSetLocalSubdomains(), PCGASMGetLocalSubdomains()
1468b1a0a8a3SJed Brown @*/
1469f746d493SDmitry Karpeev PetscErrorCode PETSCKSP_DLLEXPORT PCGASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1470b1a0a8a3SJed Brown {
1471f746d493SDmitry Karpeev   PC_GASM         *osm;
1472b1a0a8a3SJed Brown   PetscErrorCode ierr;
1473b1a0a8a3SJed Brown   PetscBool      match;
1474b1a0a8a3SJed Brown 
1475b1a0a8a3SJed Brown   PetscFunctionBegin;
1476b1a0a8a3SJed Brown   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1477b1a0a8a3SJed Brown   PetscValidIntPointer(n,2);
1478b1a0a8a3SJed Brown   if (mat) PetscValidPointer(mat,3);
1479b1a0a8a3SJed Brown   if (!pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp().");
1480f746d493SDmitry Karpeev   ierr = PetscTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr);
1481b1a0a8a3SJed Brown   if (!match) {
1482b1a0a8a3SJed Brown     if (n)   *n   = 0;
1483b1a0a8a3SJed Brown     if (mat) *mat = PETSC_NULL;
1484b1a0a8a3SJed Brown   } else {
1485f746d493SDmitry Karpeev     osm = (PC_GASM*)pc->data;
1486ab3e923fSDmitry Karpeev     if (n)   *n   = osm->n;
1487b1a0a8a3SJed Brown     if (mat) *mat = osm->pmat;
1488b1a0a8a3SJed Brown   }
1489b1a0a8a3SJed Brown   PetscFunctionReturn(0);
1490b1a0a8a3SJed Brown }
1491