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