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