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