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