xref: /petsc/src/ksp/pc/impls/asm/asm.c (revision 2fa40bb9206b96114faa7cb222621ec184d31cd2)
1 /*
2   This file defines an additive Schwarz preconditioner for any Mat implementation.
3 
4   Note that each processor may have any number of subdomains. But in order to
5   deal easily with the VecScatter(), we treat each processor as if it has the
6   same number of subdomains.
7 
8        n - total number of true subdomains on all processors
9        n_local_true - actual number of subdomains on this processor
10        n_local = maximum over all processors of n_local_true
11 */
12 
13 #include <petsc/private/pcasmimpl.h> /*I "petscpc.h" I*/
14 
15 static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer)
16 {
17   PC_ASM            *osm = (PC_ASM*)pc->data;
18   PetscErrorCode    ierr;
19   PetscMPIInt       rank;
20   PetscInt          i,bsz;
21   PetscBool         iascii,isstring;
22   PetscViewer       sviewer;
23   PetscViewerFormat format;
24   const char        *prefix;
25 
26   PetscFunctionBegin;
27   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
28   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
29   if (iascii) {
30     char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set";
31     if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof(overlaps),"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);}
32     if (osm->n > 0) {ierr = PetscSNPrintf(blocks,sizeof(blocks),"total subdomain blocks = %D",osm->n);CHKERRQ(ierr);}
33     ierr = PetscViewerASCIIPrintf(viewer,"  %s, %s\n",blocks,overlaps);CHKERRQ(ierr);
34     ierr = PetscViewerASCIIPrintf(viewer,"  restriction/interpolation type - %s\n",PCASMTypes[osm->type]);CHKERRQ(ierr);
35     if (osm->dm_subdomains) {ierr = PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: using DM to define subdomains\n");CHKERRQ(ierr);}
36     if (osm->loctype != PC_COMPOSITE_ADDITIVE) {ierr = PetscViewerASCIIPrintf(viewer,"  Additive Schwarz: local solve composition type - %s\n",PCCompositeTypes[osm->loctype]);CHKERRQ(ierr);}
37     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRMPI(ierr);
38     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
39     if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
40       if (osm->ksp) {
41         ierr = PetscViewerASCIIPrintf(viewer,"  Local solver information for first block is in the following KSP and PC objects on rank 0:\n");CHKERRQ(ierr);
42         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
43         ierr = PetscViewerASCIIPrintf(viewer,"  Use -%sksp_view ::ascii_info_detail to display information for all blocks\n",prefix?prefix:"");CHKERRQ(ierr);
44         ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
45         if (rank == 0) {
46           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
47           ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);
48           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
49         }
50         ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
51       }
52     } else {
53       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
54       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,osm->n_local_true);CHKERRQ(ierr);
55       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
56       ierr = PetscViewerASCIIPrintf(viewer,"  Local solver information for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr);
57       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
58       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
59       ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
60       for (i=0; i<osm->n_local_true; i++) {
61         ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr);
62         ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr);
63         ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr);
64         ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
65       }
66       ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
67       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
68       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
69       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
70     }
71   } else if (isstring) {
72     ierr = PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCASMTypes[osm->type]);CHKERRQ(ierr);
73     ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
74     if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);}
75     ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
76   }
77   PetscFunctionReturn(0);
78 }
79 
80 static PetscErrorCode PCASMPrintSubdomains(PC pc)
81 {
82   PC_ASM         *osm = (PC_ASM*)pc->data;
83   const char     *prefix;
84   char           fname[PETSC_MAX_PATH_LEN+1];
85   PetscViewer    viewer, sviewer;
86   char           *s;
87   PetscInt       i,j,nidx;
88   const PetscInt *idx;
89   PetscMPIInt    rank, size;
90   PetscErrorCode ierr;
91 
92   PetscFunctionBegin;
93   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRMPI(ierr);
94   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRMPI(ierr);
95   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
96   ierr = PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,sizeof(fname),NULL);CHKERRQ(ierr);
97   if (fname[0] == 0) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); };
98   ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr);
99   for (i=0; i<osm->n_local; i++) {
100     if (i < osm->n_local_true) {
101       ierr = ISGetLocalSize(osm->is[i],&nidx);CHKERRQ(ierr);
102       ierr = ISGetIndices(osm->is[i],&idx);CHKERRQ(ierr);
103       /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
104 #define len  16*(nidx+1)+512
105       ierr = PetscMalloc1(len,&s);CHKERRQ(ierr);
106       ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer);CHKERRQ(ierr);
107 #undef len
108       ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D with overlap:\n", rank, size, i);CHKERRQ(ierr);
109       for (j=0; j<nidx; j++) {
110         ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr);
111       }
112       ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr);
113       ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr);
114       ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
115       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
116       ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr);
117       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
118       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
119       ierr = PetscFree(s);CHKERRQ(ierr);
120       if (osm->is_local) {
121         /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
122 #define len  16*(nidx+1)+512
123         ierr = PetscMalloc1(len, &s);CHKERRQ(ierr);
124         ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer);CHKERRQ(ierr);
125 #undef len
126         ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D without overlap:\n", rank, size, i);CHKERRQ(ierr);
127         ierr = ISGetLocalSize(osm->is_local[i],&nidx);CHKERRQ(ierr);
128         ierr = ISGetIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
129         for (j=0; j<nidx; j++) {
130           ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr);
131         }
132         ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr);
133         ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr);
134         ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr);
135         ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
136         ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr);
137         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
138         ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
139         ierr = PetscFree(s);CHKERRQ(ierr);
140       }
141     } else {
142       /* Participate in collective viewer calls. */
143       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
144       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
145       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
146       /* Assume either all ranks have is_local or none do. */
147       if (osm->is_local) {
148         ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
149         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
150         ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
151       }
152     }
153   }
154   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
155   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
156   PetscFunctionReturn(0);
157 }
158 
159 static PetscErrorCode PCSetUp_ASM(PC pc)
160 {
161   PC_ASM         *osm = (PC_ASM*)pc->data;
162   PetscErrorCode ierr;
163   PetscBool      flg;
164   PetscInt       i,m,m_local;
165   MatReuse       scall = MAT_REUSE_MATRIX;
166   IS             isl;
167   KSP            ksp;
168   PC             subpc;
169   const char     *prefix,*pprefix;
170   Vec            vec;
171   DM             *domain_dm = NULL;
172 
173   PetscFunctionBegin;
174   if (!pc->setupcalled) {
175     PetscInt m;
176 
177     /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */
178     if (osm->n_local_true == PETSC_DECIDE) {
179       /* no subdomains given */
180       /* try pc->dm first, if allowed */
181       if (osm->dm_subdomains && pc->dm) {
182         PetscInt  num_domains, d;
183         char      **domain_names;
184         IS        *inner_domain_is, *outer_domain_is;
185         ierr = DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm);CHKERRQ(ierr);
186         osm->overlap = -1; /* We do not want to increase the overlap of the IS.
187                               A future improvement of this code might allow one to use
188                               DM-defined subdomains and also increase the overlap,
189                               but that is not currently supported */
190         if (num_domains) {
191           ierr = PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is);CHKERRQ(ierr);
192         }
193         for (d = 0; d < num_domains; ++d) {
194           if (domain_names)    {ierr = PetscFree(domain_names[d]);CHKERRQ(ierr);}
195           if (inner_domain_is) {ierr = ISDestroy(&inner_domain_is[d]);CHKERRQ(ierr);}
196           if (outer_domain_is) {ierr = ISDestroy(&outer_domain_is[d]);CHKERRQ(ierr);}
197         }
198         ierr = PetscFree(domain_names);CHKERRQ(ierr);
199         ierr = PetscFree(inner_domain_is);CHKERRQ(ierr);
200         ierr = PetscFree(outer_domain_is);CHKERRQ(ierr);
201       }
202       if (osm->n_local_true == PETSC_DECIDE) {
203         /* still no subdomains; use one subdomain per processor */
204         osm->n_local_true = 1;
205       }
206     }
207     { /* determine the global and max number of subdomains */
208       struct {PetscInt max,sum;} inwork,outwork;
209       PetscMPIInt size;
210 
211       inwork.max   = osm->n_local_true;
212       inwork.sum   = osm->n_local_true;
213       ierr         = MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc));CHKERRMPI(ierr);
214       osm->n_local = outwork.max;
215       osm->n       = outwork.sum;
216 
217       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRMPI(ierr);
218       if (outwork.max == 1 && outwork.sum == size) {
219         /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */
220         ierr = MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE);CHKERRQ(ierr);
221       }
222     }
223     if (!osm->is) { /* create the index sets */
224       ierr = PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);CHKERRQ(ierr);
225     }
226     if (osm->n_local_true > 1 && !osm->is_local) {
227       ierr = PetscMalloc1(osm->n_local_true,&osm->is_local);CHKERRQ(ierr);
228       for (i=0; i<osm->n_local_true; i++) {
229         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
230           ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr);
231           ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr);
232         } else {
233           ierr             = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr);
234           osm->is_local[i] = osm->is[i];
235         }
236       }
237     }
238     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
239     flg  = PETSC_FALSE;
240     ierr = PetscOptionsGetBool(NULL,prefix,"-pc_asm_print_subdomains",&flg,NULL);CHKERRQ(ierr);
241     if (flg) { ierr = PCASMPrintSubdomains(pc);CHKERRQ(ierr); }
242 
243     if (osm->overlap > 0) {
244       /* Extend the "overlapping" regions by a number of steps */
245       ierr = MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);CHKERRQ(ierr);
246     }
247     if (osm->sort_indices) {
248       for (i=0; i<osm->n_local_true; i++) {
249         ierr = ISSort(osm->is[i]);CHKERRQ(ierr);
250         if (osm->is_local) {
251           ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr);
252         }
253       }
254     }
255 
256     if (!osm->ksp) {
257       /* Create the local solvers */
258       ierr = PetscMalloc1(osm->n_local_true,&osm->ksp);CHKERRQ(ierr);
259       if (domain_dm) {
260         ierr = PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n");CHKERRQ(ierr);
261       }
262       for (i=0; i<osm->n_local_true; i++) {
263         ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr);
264         ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr);
265         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr);
266         ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr);
267         ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
268         ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr);
269         ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
270         ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr);
271         ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr);
272         if (domain_dm) {
273           ierr = KSPSetDM(ksp, domain_dm[i]);CHKERRQ(ierr);
274           ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr);
275           ierr = DMDestroy(&domain_dm[i]);CHKERRQ(ierr);
276         }
277         osm->ksp[i] = ksp;
278       }
279       if (domain_dm) {
280         ierr = PetscFree(domain_dm);CHKERRQ(ierr);
281       }
282     }
283 
284     ierr = ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis);CHKERRQ(ierr);
285     ierr = ISSortRemoveDups(osm->lis);CHKERRQ(ierr);
286     ierr = ISGetLocalSize(osm->lis, &m);CHKERRQ(ierr);
287 
288     scall = MAT_INITIAL_MATRIX;
289   } else {
290     /*
291        Destroy the blocks from the previous iteration
292     */
293     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
294       ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
295       scall = MAT_INITIAL_MATRIX;
296     }
297   }
298 
299   /* Destroy previous submatrices of a different type than pc->pmat since MAT_REUSE_MATRIX won't work in that case */
300   if ((scall == MAT_REUSE_MATRIX) && osm->sub_mat_type) {
301     if (osm->n_local_true > 0) {
302       ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
303     }
304     scall = MAT_INITIAL_MATRIX;
305   }
306 
307   /*
308      Extract out the submatrices
309   */
310   ierr = MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr);
311   if (scall == MAT_INITIAL_MATRIX) {
312     ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr);
313     for (i=0; i<osm->n_local_true; i++) {
314       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr);
315       ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr);
316     }
317   }
318 
319   /* Convert the types of the submatrices (if needbe) */
320   if (osm->sub_mat_type) {
321     for (i=0; i<osm->n_local_true; i++) {
322       ierr = MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]));CHKERRQ(ierr);
323     }
324   }
325 
326   if (!pc->setupcalled) {
327     VecType vtype;
328 
329     /* Create the local work vectors (from the local matrices) and scatter contexts */
330     ierr = MatCreateVecs(pc->pmat,&vec,NULL);CHKERRQ(ierr);
331 
332     if (osm->is_local && (osm->type == PC_ASM_INTERPOLATE || osm->type == PC_ASM_NONE)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot use interpolate or none PCASMType if is_local was provided to PCASMSetLocalSubdomains()");
333     if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) {
334       ierr = PetscMalloc1(osm->n_local_true,&osm->lprolongation);CHKERRQ(ierr);
335     }
336     ierr = PetscMalloc1(osm->n_local_true,&osm->lrestriction);CHKERRQ(ierr);
337     ierr = PetscMalloc1(osm->n_local_true,&osm->x);CHKERRQ(ierr);
338     ierr = PetscMalloc1(osm->n_local_true,&osm->y);CHKERRQ(ierr);
339 
340     ierr = ISGetLocalSize(osm->lis,&m);CHKERRQ(ierr);
341     ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
342     ierr = MatGetVecType(osm->pmat[0],&vtype);CHKERRQ(ierr);
343     ierr = VecCreate(PETSC_COMM_SELF,&osm->lx);CHKERRQ(ierr);
344     ierr = VecSetSizes(osm->lx,m,m);CHKERRQ(ierr);
345     ierr = VecSetType(osm->lx,vtype);CHKERRQ(ierr);
346     ierr = VecDuplicate(osm->lx, &osm->ly);CHKERRQ(ierr);
347     ierr = VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->restriction);CHKERRQ(ierr);
348     ierr = ISDestroy(&isl);CHKERRQ(ierr);
349 
350     for (i=0; i<osm->n_local_true; ++i) {
351       ISLocalToGlobalMapping ltog;
352       IS                     isll;
353       const PetscInt         *idx_is;
354       PetscInt               *idx_lis,nout;
355 
356       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
357       ierr = MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);CHKERRQ(ierr);
358       ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr);
359 
360       /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */
361       ierr = ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);CHKERRQ(ierr);
362       ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr);
363       ierr = ISGetIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
364       ierr = PetscMalloc1(m,&idx_lis);CHKERRQ(ierr);
365       ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis);CHKERRQ(ierr);
366       if (nout != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis");
367       ierr = ISRestoreIndices(osm->is[i], &idx_is);CHKERRQ(ierr);
368       ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
369       ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
370       ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr);
371       ierr = VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lrestriction[i]);CHKERRQ(ierr);
372       ierr = ISDestroy(&isll);CHKERRQ(ierr);
373       ierr = ISDestroy(&isl);CHKERRQ(ierr);
374       if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overlapping is_local[i] entries */
375         ISLocalToGlobalMapping ltog;
376         IS                     isll,isll_local;
377         const PetscInt         *idx_local;
378         PetscInt               *idx1, *idx2, nout;
379 
380         ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr);
381         ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
382 
383         ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],&ltog);CHKERRQ(ierr);
384         ierr = PetscMalloc1(m_local,&idx1);CHKERRQ(ierr);
385         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1);CHKERRQ(ierr);
386         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
387         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is");
388         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr);
389 
390         ierr = ISLocalToGlobalMappingCreateIS(osm->lis,&ltog);CHKERRQ(ierr);
391         ierr = PetscMalloc1(m_local,&idx2);CHKERRQ(ierr);
392         ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2);CHKERRQ(ierr);
393         ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
394         if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis");
395         ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local);CHKERRQ(ierr);
396 
397         ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr);
398         ierr = VecScatterCreate(osm->y[i],isll,osm->ly,isll_local,&osm->lprolongation[i]);CHKERRQ(ierr);
399 
400         ierr = ISDestroy(&isll);CHKERRQ(ierr);
401         ierr = ISDestroy(&isll_local);CHKERRQ(ierr);
402       }
403     }
404     ierr = VecDestroy(&vec);CHKERRQ(ierr);
405   }
406 
407   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
408     IS      *cis;
409     PetscInt c;
410 
411     ierr = PetscMalloc1(osm->n_local_true, &cis);CHKERRQ(ierr);
412     for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis;
413     ierr = MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);CHKERRQ(ierr);
414     ierr = PetscFree(cis);CHKERRQ(ierr);
415   }
416 
417   /* Return control to the user so that the submatrices can be modified (e.g., to apply
418      different boundary conditions for the submatrices than for the global problem) */
419   ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr);
420 
421   /*
422      Loop over subdomains putting them into local ksp
423   */
424   ierr = KSPGetOptionsPrefix(osm->ksp[0],&prefix);CHKERRQ(ierr);
425   for (i=0; i<osm->n_local_true; i++) {
426     ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr);
427     ierr = MatSetOptionsPrefix(osm->pmat[i],prefix);CHKERRQ(ierr);
428     if (!pc->setupcalled) {
429       ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr);
430     }
431   }
432   PetscFunctionReturn(0);
433 }
434 
435 static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc)
436 {
437   PC_ASM             *osm = (PC_ASM*)pc->data;
438   PetscErrorCode     ierr;
439   PetscInt           i;
440   KSPConvergedReason reason;
441 
442   PetscFunctionBegin;
443   for (i=0; i<osm->n_local_true; i++) {
444     ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr);
445     ierr = KSPGetConvergedReason(osm->ksp[i],&reason);CHKERRQ(ierr);
446     if (reason == KSP_DIVERGED_PC_FAILED) {
447       pc->failedreason = PC_SUBPC_ERROR;
448     }
449   }
450   PetscFunctionReturn(0);
451 }
452 
453 static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y)
454 {
455   PC_ASM         *osm = (PC_ASM*)pc->data;
456   PetscErrorCode ierr;
457   PetscInt       i,n_local_true = osm->n_local_true;
458   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
459 
460   PetscFunctionBegin;
461   /*
462      support for limiting the restriction or interpolation to only local
463      subdomain values (leaving the other values 0).
464   */
465   if (!(osm->type & PC_ASM_RESTRICT)) {
466     forward = SCATTER_FORWARD_LOCAL;
467     /* have to zero the work RHS since scatter may leave some slots empty */
468     ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr);
469   }
470   if (!(osm->type & PC_ASM_INTERPOLATE)) {
471     reverse = SCATTER_REVERSE_LOCAL;
472   }
473 
474   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) {
475     /* zero the global and the local solutions */
476     ierr = VecSet(y, 0.0);CHKERRQ(ierr);
477     ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
478 
479     /* copy the global RHS to local RHS including the ghost nodes */
480     ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
481     ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
482 
483     /* restrict local RHS to the overlapping 0-block RHS */
484     ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
485     ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
486 
487     /* do the local solves */
488     for (i = 0; i < n_local_true; ++i) {
489 
490       /* solve the overlapping i-block */
491       ierr = PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i],0);CHKERRQ(ierr);
492       ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
493       ierr = KSPCheckSolve(osm->ksp[i], pc, osm->y[i]);CHKERRQ(ierr);
494       ierr = PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0);CHKERRQ(ierr);
495 
496       if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution (only for restrictive additive) */
497         ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
498         ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
499       } else { /* interpolate the overlapping i-block solution to the local solution */
500         ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
501         ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
502       }
503 
504       if (i < n_local_true-1) {
505         /* restrict local RHS to the overlapping (i+1)-block RHS */
506         ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
507         ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
508 
509         if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
510           /* update the overlapping (i+1)-block RHS using the current local solution */
511           ierr = MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]);CHKERRQ(ierr);
512           ierr = VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1]);CHKERRQ(ierr);
513         }
514       }
515     }
516     /* add the local solution to the global solution including the ghost nodes */
517     ierr = VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
518     ierr = VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
519   } else SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
520   PetscFunctionReturn(0);
521 }
522 
523 static PetscErrorCode PCMatApply_ASM(PC pc,Mat X,Mat Y)
524 {
525   PC_ASM         *osm = (PC_ASM*)pc->data;
526   Mat            Z,W;
527   Vec            x;
528   PetscInt       i,m,N;
529   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
530   PetscErrorCode ierr;
531 
532   PetscFunctionBegin;
533   if (osm->n_local_true > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not yet implemented");
534   /*
535      support for limiting the restriction or interpolation to only local
536      subdomain values (leaving the other values 0).
537   */
538   if (!(osm->type & PC_ASM_RESTRICT)) {
539     forward = SCATTER_FORWARD_LOCAL;
540     /* have to zero the work RHS since scatter may leave some slots empty */
541     ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr);
542   }
543   if (!(osm->type & PC_ASM_INTERPOLATE)) {
544     reverse = SCATTER_REVERSE_LOCAL;
545   }
546   ierr = VecGetLocalSize(osm->x[0], &m);CHKERRQ(ierr);
547   ierr = MatGetSize(X, NULL, &N);CHKERRQ(ierr);
548   ierr = MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &Z);CHKERRQ(ierr);
549   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) {
550     /* zero the global and the local solutions */
551     ierr = MatZeroEntries(Y);CHKERRQ(ierr);
552     ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
553 
554     for (i = 0; i < N; ++i) {
555       ierr = MatDenseGetColumnVecRead(X, i, &x);CHKERRQ(ierr);
556       /* copy the global RHS to local RHS including the ghost nodes */
557       ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
558       ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
559       ierr = MatDenseRestoreColumnVecRead(X, i, &x);CHKERRQ(ierr);
560 
561       ierr = MatDenseGetColumnVecWrite(Z, i, &x);CHKERRQ(ierr);
562       /* restrict local RHS to the overlapping 0-block RHS */
563       ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward);CHKERRQ(ierr);
564       ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward);CHKERRQ(ierr);
565       ierr = MatDenseRestoreColumnVecWrite(Z, i, &x);CHKERRQ(ierr);
566     }
567     ierr = MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &W);CHKERRQ(ierr);
568     /* solve the overlapping 0-block */
569     ierr = PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0);CHKERRQ(ierr);
570     ierr = KSPMatSolve(osm->ksp[0], Z, W);CHKERRQ(ierr);
571     ierr = KSPCheckSolve(osm->ksp[0], pc, NULL);CHKERRQ(ierr);
572     ierr = PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[0], Z, W,0);CHKERRQ(ierr);
573     ierr = MatDestroy(&Z);CHKERRQ(ierr);
574 
575     for (i = 0; i < N; ++i) {
576       ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
577       ierr = MatDenseGetColumnVecRead(W, i, &x);CHKERRQ(ierr);
578       if (osm->lprolongation) { /* interpolate the non-overlapping 0-block solution to the local solution (only for restrictive additive) */
579         ierr = VecScatterBegin(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
580         ierr = VecScatterEnd(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
581       } else { /* interpolate the overlapping 0-block solution to the local solution */
582         ierr = VecScatterBegin(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
583         ierr = VecScatterEnd(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
584       }
585       ierr = MatDenseRestoreColumnVecRead(W, i, &x);CHKERRQ(ierr);
586 
587       ierr = MatDenseGetColumnVecWrite(Y, i, &x);CHKERRQ(ierr);
588       /* add the local solution to the global solution including the ghost nodes */
589       ierr = VecScatterBegin(osm->restriction, osm->ly, x, ADD_VALUES, reverse);CHKERRQ(ierr);
590       ierr = VecScatterEnd(osm->restriction, osm->ly, x, ADD_VALUES, reverse);CHKERRQ(ierr);
591       ierr = MatDenseRestoreColumnVecWrite(Y, i, &x);CHKERRQ(ierr);
592     }
593     ierr = MatDestroy(&W);CHKERRQ(ierr);
594   } else SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
595   PetscFunctionReturn(0);
596 }
597 
598 static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y)
599 {
600   PC_ASM         *osm = (PC_ASM*)pc->data;
601   PetscErrorCode ierr;
602   PetscInt       i,n_local_true = osm->n_local_true;
603   ScatterMode    forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE;
604 
605   PetscFunctionBegin;
606   /*
607      Support for limiting the restriction or interpolation to only local
608      subdomain values (leaving the other values 0).
609 
610      Note: these are reversed from the PCApply_ASM() because we are applying the
611      transpose of the three terms
612   */
613 
614   if (!(osm->type & PC_ASM_INTERPOLATE)) {
615     forward = SCATTER_FORWARD_LOCAL;
616     /* have to zero the work RHS since scatter may leave some slots empty */
617     ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr);
618   }
619   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;
620 
621   /* zero the global and the local solutions */
622   ierr = VecSet(y, 0.0);CHKERRQ(ierr);
623   ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr);
624 
625   /* Copy the global RHS to local RHS including the ghost nodes */
626   ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
627   ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr);
628 
629   /* Restrict local RHS to the overlapping 0-block RHS */
630   ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
631   ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr);
632 
633   /* do the local solves */
634   for (i = 0; i < n_local_true; ++i) {
635 
636     /* solve the overlapping i-block */
637     ierr = PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
638     ierr = KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr);
639     ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr);
640     ierr = PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr);
641 
642     if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution */
643       ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
644       ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr);
645     } else { /* interpolate the overlapping i-block solution to the local solution */
646       ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
647       ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr);
648     }
649 
650     if (i < n_local_true-1) {
651       /* Restrict local RHS to the overlapping (i+1)-block RHS */
652       ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
653       ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr);
654     }
655   }
656   /* Add the local solution to the global solution including the ghost nodes */
657   ierr = VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
658   ierr = VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr);
659 
660   PetscFunctionReturn(0);
661 
662 }
663 
664 static PetscErrorCode PCReset_ASM(PC pc)
665 {
666   PC_ASM         *osm = (PC_ASM*)pc->data;
667   PetscErrorCode ierr;
668   PetscInt       i;
669 
670   PetscFunctionBegin;
671   if (osm->ksp) {
672     for (i=0; i<osm->n_local_true; i++) {
673       ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr);
674     }
675   }
676   if (osm->pmat) {
677     if (osm->n_local_true > 0) {
678       ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr);
679     }
680   }
681   if (osm->lrestriction) {
682     ierr = VecScatterDestroy(&osm->restriction);CHKERRQ(ierr);
683     for (i=0; i<osm->n_local_true; i++) {
684       ierr = VecScatterDestroy(&osm->lrestriction[i]);CHKERRQ(ierr);
685       if (osm->lprolongation) {ierr = VecScatterDestroy(&osm->lprolongation[i]);CHKERRQ(ierr);}
686       ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr);
687       ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr);
688     }
689     ierr = PetscFree(osm->lrestriction);CHKERRQ(ierr);
690     if (osm->lprolongation) {ierr = PetscFree(osm->lprolongation);CHKERRQ(ierr);}
691     ierr = PetscFree(osm->x);CHKERRQ(ierr);
692     ierr = PetscFree(osm->y);CHKERRQ(ierr);
693 
694   }
695   ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
696   ierr = ISDestroy(&osm->lis);CHKERRQ(ierr);
697   ierr = VecDestroy(&osm->lx);CHKERRQ(ierr);
698   ierr = VecDestroy(&osm->ly);CHKERRQ(ierr);
699   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
700     ierr = MatDestroyMatrices(osm->n_local_true, &osm->lmats);CHKERRQ(ierr);
701   }
702 
703   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
704 
705   osm->is       = NULL;
706   osm->is_local = NULL;
707   PetscFunctionReturn(0);
708 }
709 
710 static PetscErrorCode PCDestroy_ASM(PC pc)
711 {
712   PC_ASM         *osm = (PC_ASM*)pc->data;
713   PetscErrorCode ierr;
714   PetscInt       i;
715 
716   PetscFunctionBegin;
717   ierr = PCReset_ASM(pc);CHKERRQ(ierr);
718   if (osm->ksp) {
719     for (i=0; i<osm->n_local_true; i++) {
720       ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr);
721     }
722     ierr = PetscFree(osm->ksp);CHKERRQ(ierr);
723   }
724   ierr = PetscFree(pc->data);CHKERRQ(ierr);
725 
726   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",NULL);CHKERRQ(ierr);
727   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",NULL);CHKERRQ(ierr);
728   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",NULL);CHKERRQ(ierr);
729   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",NULL);CHKERRQ(ierr);
730   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",NULL);CHKERRQ(ierr);
731   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",NULL);CHKERRQ(ierr);
732   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",NULL);CHKERRQ(ierr);
733   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",NULL);CHKERRQ(ierr);
734   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",NULL);CHKERRQ(ierr);
735   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",NULL);CHKERRQ(ierr);
736   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",NULL);CHKERRQ(ierr);
737   PetscFunctionReturn(0);
738 }
739 
740 static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc)
741 {
742   PC_ASM         *osm = (PC_ASM*)pc->data;
743   PetscErrorCode ierr;
744   PetscInt       blocks,ovl;
745   PetscBool      flg;
746   PCASMType      asmtype;
747   PCCompositeType loctype;
748   char           sub_mat_type[256];
749 
750   PetscFunctionBegin;
751   ierr = PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");CHKERRQ(ierr);
752   ierr = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr);
753   ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr);
754   if (flg) {
755     ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr);
756     osm->dm_subdomains = PETSC_FALSE;
757   }
758   ierr = PetscOptionsInt("-pc_asm_local_blocks","Number of local subdomains","PCASMSetLocalSubdomains",osm->n_local_true,&blocks,&flg);CHKERRQ(ierr);
759   if (flg) {
760     ierr = PCASMSetLocalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr);
761     osm->dm_subdomains = PETSC_FALSE;
762   }
763   ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr);
764   if (flg) {
765     ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr);
766     osm->dm_subdomains = PETSC_FALSE;
767   }
768   flg  = PETSC_FALSE;
769   ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
770   if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); }
771   flg  = PETSC_FALSE;
772   ierr = PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr);
773   if (flg) {ierr = PCASMSetLocalType(pc,loctype);CHKERRQ(ierr); }
774   ierr = PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);CHKERRQ(ierr);
775   if (flg) {
776     ierr = PCASMSetSubMatType(pc,sub_mat_type);CHKERRQ(ierr);
777   }
778   ierr = PetscOptionsTail();CHKERRQ(ierr);
779   PetscFunctionReturn(0);
780 }
781 
782 /*------------------------------------------------------------------------------------*/
783 
784 static PetscErrorCode  PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[])
785 {
786   PC_ASM         *osm = (PC_ASM*)pc->data;
787   PetscErrorCode ierr;
788   PetscInt       i;
789 
790   PetscFunctionBegin;
791   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n);
792   if (pc->setupcalled && (n != osm->n_local_true || is)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetLocalSubdomains() should be called before calling PCSetUp().");
793 
794   if (!pc->setupcalled) {
795     if (is) {
796       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);}
797     }
798     if (is_local) {
799       for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);}
800     }
801     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
802 
803     osm->n_local_true = n;
804     osm->is           = NULL;
805     osm->is_local     = NULL;
806     if (is) {
807       ierr = PetscMalloc1(n,&osm->is);CHKERRQ(ierr);
808       for (i=0; i<n; i++) osm->is[i] = is[i];
809       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
810       osm->overlap = -1;
811     }
812     if (is_local) {
813       ierr = PetscMalloc1(n,&osm->is_local);CHKERRQ(ierr);
814       for (i=0; i<n; i++) osm->is_local[i] = is_local[i];
815       if (!is) {
816         ierr = PetscMalloc1(osm->n_local_true,&osm->is);CHKERRQ(ierr);
817         for (i=0; i<osm->n_local_true; i++) {
818           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
819             ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr);
820             ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr);
821           } else {
822             ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr);
823             osm->is[i] = osm->is_local[i];
824           }
825         }
826       }
827     }
828   }
829   PetscFunctionReturn(0);
830 }
831 
832 static PetscErrorCode  PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local)
833 {
834   PC_ASM         *osm = (PC_ASM*)pc->data;
835   PetscErrorCode ierr;
836   PetscMPIInt    rank,size;
837   PetscInt       n;
838 
839   PetscFunctionBegin;
840   if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N);
841   if (is || is_local) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Use PCASMSetLocalSubdomains() to set specific index sets\n\they cannot be set globally yet.");
842 
843   /*
844      Split the subdomains equally among all processors
845   */
846   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRMPI(ierr);
847   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRMPI(ierr);
848   n    = N/size + ((N % size) > rank);
849   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);
850   if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp().");
851   if (!pc->setupcalled) {
852     ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr);
853 
854     osm->n_local_true = n;
855     osm->is           = NULL;
856     osm->is_local     = NULL;
857   }
858   PetscFunctionReturn(0);
859 }
860 
861 static PetscErrorCode  PCASMSetOverlap_ASM(PC pc,PetscInt ovl)
862 {
863   PC_ASM *osm = (PC_ASM*)pc->data;
864 
865   PetscFunctionBegin;
866   if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested");
867   if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp().");
868   if (!pc->setupcalled) osm->overlap = ovl;
869   PetscFunctionReturn(0);
870 }
871 
872 static PetscErrorCode  PCASMSetType_ASM(PC pc,PCASMType type)
873 {
874   PC_ASM *osm = (PC_ASM*)pc->data;
875 
876   PetscFunctionBegin;
877   osm->type     = type;
878   osm->type_set = PETSC_TRUE;
879   PetscFunctionReturn(0);
880 }
881 
882 static PetscErrorCode  PCASMGetType_ASM(PC pc,PCASMType *type)
883 {
884   PC_ASM *osm = (PC_ASM*)pc->data;
885 
886   PetscFunctionBegin;
887   *type = osm->type;
888   PetscFunctionReturn(0);
889 }
890 
891 static PetscErrorCode  PCASMSetLocalType_ASM(PC pc, PCCompositeType type)
892 {
893   PC_ASM *osm = (PC_ASM *) pc->data;
894 
895   PetscFunctionBegin;
896   if (type != PC_COMPOSITE_ADDITIVE && type != PC_COMPOSITE_MULTIPLICATIVE) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Only supports additive or multiplicative as the local type");
897   osm->loctype = type;
898   PetscFunctionReturn(0);
899 }
900 
901 static PetscErrorCode  PCASMGetLocalType_ASM(PC pc, PCCompositeType *type)
902 {
903   PC_ASM *osm = (PC_ASM *) pc->data;
904 
905   PetscFunctionBegin;
906   *type = osm->loctype;
907   PetscFunctionReturn(0);
908 }
909 
910 static PetscErrorCode  PCASMSetSortIndices_ASM(PC pc,PetscBool  doSort)
911 {
912   PC_ASM *osm = (PC_ASM*)pc->data;
913 
914   PetscFunctionBegin;
915   osm->sort_indices = doSort;
916   PetscFunctionReturn(0);
917 }
918 
919 static PetscErrorCode  PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp)
920 {
921   PC_ASM         *osm = (PC_ASM*)pc->data;
922   PetscErrorCode ierr;
923 
924   PetscFunctionBegin;
925   if (osm->n_local_true < 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Need to call PCSetUp() on PC (or KSPSetUp() on the outer KSP object) before calling here");
926 
927   if (n_local) *n_local = osm->n_local_true;
928   if (first_local) {
929     ierr          = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRMPI(ierr);
930     *first_local -= osm->n_local_true;
931   }
932   if (ksp) *ksp   = osm->ksp;
933   PetscFunctionReturn(0);
934 }
935 
936 static PetscErrorCode  PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type)
937 {
938   PC_ASM         *osm = (PC_ASM*)pc->data;
939 
940   PetscFunctionBegin;
941   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
942   PetscValidPointer(sub_mat_type,2);
943   *sub_mat_type = osm->sub_mat_type;
944   PetscFunctionReturn(0);
945 }
946 
947 static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type)
948 {
949   PetscErrorCode    ierr;
950   PC_ASM            *osm = (PC_ASM*)pc->data;
951 
952   PetscFunctionBegin;
953   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
954   ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr);
955   ierr = PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);CHKERRQ(ierr);
956   PetscFunctionReturn(0);
957 }
958 
959 /*@C
960     PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner.
961 
962     Collective on pc
963 
964     Input Parameters:
965 +   pc - the preconditioner context
966 .   n - the number of subdomains for this processor (default value = 1)
967 .   is - the index set that defines the subdomains for this processor
968          (or NULL for PETSc to determine subdomains)
969 -   is_local - the index sets that define the local part of the subdomains for this processor, not used unless PCASMType is PC_ASM_RESTRICT
970          (or NULL to not provide these)
971 
972     Options Database Key:
973     To set the total number of subdomain blocks rather than specify the
974     index sets, use the option
975 .    -pc_asm_local_blocks <blks> - Sets local blocks
976 
977     Notes:
978     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
979 
980     By default the ASM preconditioner uses 1 block per processor.
981 
982     Use PCASMSetTotalSubdomains() to set the subdomains for all processors.
983 
984     If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated
985     back to form the global solution (this is the standard restricted additive Schwarz method)
986 
987     If the is_local is provided and PCASMType is PC_ASM_INTERPOLATE or PC_ASM_NONE then an error is generated since there is
988     no code to handle that case.
989 
990     Level: advanced
991 
992 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
993           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType()
994 @*/
995 PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
996 {
997   PetscErrorCode ierr;
998 
999   PetscFunctionBegin;
1000   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1001   ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr);
1002   PetscFunctionReturn(0);
1003 }
1004 
1005 /*@C
1006     PCASMSetTotalSubdomains - Sets the subdomains for all processors for the
1007     additive Schwarz preconditioner.  Either all or no processors in the
1008     PC communicator must call this routine, with the same index sets.
1009 
1010     Collective on pc
1011 
1012     Input Parameters:
1013 +   pc - the preconditioner context
1014 .   N  - the number of subdomains for all processors
1015 .   is - the index sets that define the subdomains for all processors
1016          (or NULL to ask PETSc to determine the subdomains)
1017 -   is_local - the index sets that define the local part of the subdomains for this processor
1018          (or NULL to not provide this information)
1019 
1020     Options Database Key:
1021     To set the total number of subdomain blocks rather than specify the
1022     index sets, use the option
1023 .    -pc_asm_blocks <blks> - Sets total blocks
1024 
1025     Notes:
1026     Currently you cannot use this to set the actual subdomains with the argument is or is_local.
1027 
1028     By default the ASM preconditioner uses 1 block per processor.
1029 
1030     These index sets cannot be destroyed until after completion of the
1031     linear solves for which the ASM preconditioner is being used.
1032 
1033     Use PCASMSetLocalSubdomains() to set local subdomains.
1034 
1035     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
1036 
1037     Level: advanced
1038 
1039 .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1040           PCASMCreateSubdomains2D()
1041 @*/
1042 PetscErrorCode  PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[])
1043 {
1044   PetscErrorCode ierr;
1045 
1046   PetscFunctionBegin;
1047   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1048   ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr);
1049   PetscFunctionReturn(0);
1050 }
1051 
1052 /*@
1053     PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
1054     additive Schwarz preconditioner.  Either all or no processors in the
1055     PC communicator must call this routine.
1056 
1057     Logically Collective on pc
1058 
1059     Input Parameters:
1060 +   pc  - the preconditioner context
1061 -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
1062 
1063     Options Database Key:
1064 .   -pc_asm_overlap <ovl> - Sets overlap
1065 
1066     Notes:
1067     By default the ASM preconditioner uses 1 block per processor.  To use
1068     multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
1069     PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).
1070 
1071     The overlap defaults to 1, so if one desires that no additional
1072     overlap be computed beyond what may have been set with a call to
1073     PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
1074     must be set to be 0.  In particular, if one does not explicitly set
1075     the subdomains an application code, then all overlap would be computed
1076     internally by PETSc, and using an overlap of 0 would result in an ASM
1077     variant that is equivalent to the block Jacobi preconditioner.
1078 
1079     The default algorithm used by PETSc to increase overlap is fast, but not scalable,
1080     use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
1081 
1082     Note that one can define initial index sets with any overlap via
1083     PCASMSetLocalSubdomains(); the routine
1084     PCASMSetOverlap() merely allows PETSc to extend that overlap further
1085     if desired.
1086 
1087     Level: intermediate
1088 
1089 .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
1090           PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap()
1091 @*/
1092 PetscErrorCode  PCASMSetOverlap(PC pc,PetscInt ovl)
1093 {
1094   PetscErrorCode ierr;
1095 
1096   PetscFunctionBegin;
1097   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1098   PetscValidLogicalCollectiveInt(pc,ovl,2);
1099   ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr);
1100   PetscFunctionReturn(0);
1101 }
1102 
1103 /*@
1104     PCASMSetType - Sets the type of restriction and interpolation used
1105     for local problems in the additive Schwarz method.
1106 
1107     Logically Collective on pc
1108 
1109     Input Parameters:
1110 +   pc  - the preconditioner context
1111 -   type - variant of ASM, one of
1112 .vb
1113       PC_ASM_BASIC       - full interpolation and restriction
1114       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1115       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1116       PC_ASM_NONE        - local processor restriction and interpolation
1117 .ve
1118 
1119     Options Database Key:
1120 .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
1121 
1122     Notes:
1123     if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected
1124     to limit the local processor interpolation
1125 
1126     Level: intermediate
1127 
1128 .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1129           PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType()
1130 @*/
1131 PetscErrorCode  PCASMSetType(PC pc,PCASMType type)
1132 {
1133   PetscErrorCode ierr;
1134 
1135   PetscFunctionBegin;
1136   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1137   PetscValidLogicalCollectiveEnum(pc,type,2);
1138   ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr);
1139   PetscFunctionReturn(0);
1140 }
1141 
1142 /*@
1143     PCASMGetType - Gets the type of restriction and interpolation used
1144     for local problems in the additive Schwarz method.
1145 
1146     Logically Collective on pc
1147 
1148     Input Parameter:
1149 .   pc  - the preconditioner context
1150 
1151     Output Parameter:
1152 .   type - variant of ASM, one of
1153 
1154 .vb
1155       PC_ASM_BASIC       - full interpolation and restriction
1156       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1157       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1158       PC_ASM_NONE        - local processor restriction and interpolation
1159 .ve
1160 
1161     Options Database Key:
1162 .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
1163 
1164     Level: intermediate
1165 
1166 .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1167           PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType()
1168 @*/
1169 PetscErrorCode  PCASMGetType(PC pc,PCASMType *type)
1170 {
1171   PetscErrorCode ierr;
1172 
1173   PetscFunctionBegin;
1174   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1175   ierr = PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));CHKERRQ(ierr);
1176   PetscFunctionReturn(0);
1177 }
1178 
1179 /*@
1180   PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method.
1181 
1182   Logically Collective on pc
1183 
1184   Input Parameters:
1185 + pc  - the preconditioner context
1186 - type - type of composition, one of
1187 .vb
1188   PC_COMPOSITE_ADDITIVE       - local additive combination
1189   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
1190 .ve
1191 
1192   Options Database Key:
1193 . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
1194 
1195   Level: intermediate
1196 
1197 .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
1198 @*/
1199 PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type)
1200 {
1201   PetscErrorCode ierr;
1202 
1203   PetscFunctionBegin;
1204   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1205   PetscValidLogicalCollectiveEnum(pc, type, 2);
1206   ierr = PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));CHKERRQ(ierr);
1207   PetscFunctionReturn(0);
1208 }
1209 
1210 /*@
1211   PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method.
1212 
1213   Logically Collective on pc
1214 
1215   Input Parameter:
1216 . pc  - the preconditioner context
1217 
1218   Output Parameter:
1219 . type - type of composition, one of
1220 .vb
1221   PC_COMPOSITE_ADDITIVE       - local additive combination
1222   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
1223 .ve
1224 
1225   Options Database Key:
1226 . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
1227 
1228   Level: intermediate
1229 
1230 .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType
1231 @*/
1232 PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type)
1233 {
1234   PetscErrorCode ierr;
1235 
1236   PetscFunctionBegin;
1237   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1238   PetscValidPointer(type, 2);
1239   ierr = PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));CHKERRQ(ierr);
1240   PetscFunctionReturn(0);
1241 }
1242 
1243 /*@
1244     PCASMSetSortIndices - Determines whether subdomain indices are sorted.
1245 
1246     Logically Collective on pc
1247 
1248     Input Parameters:
1249 +   pc  - the preconditioner context
1250 -   doSort - sort the subdomain indices
1251 
1252     Level: intermediate
1253 
1254 .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(),
1255           PCASMCreateSubdomains2D()
1256 @*/
1257 PetscErrorCode  PCASMSetSortIndices(PC pc,PetscBool doSort)
1258 {
1259   PetscErrorCode ierr;
1260 
1261   PetscFunctionBegin;
1262   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1263   PetscValidLogicalCollectiveBool(pc,doSort,2);
1264   ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr);
1265   PetscFunctionReturn(0);
1266 }
1267 
1268 /*@C
1269    PCASMGetSubKSP - Gets the local KSP contexts for all blocks on
1270    this processor.
1271 
1272    Collective on pc iff first_local is requested
1273 
1274    Input Parameter:
1275 .  pc - the preconditioner context
1276 
1277    Output Parameters:
1278 +  n_local - the number of blocks on this processor or NULL
1279 .  first_local - the global number of the first block on this processor or NULL,
1280                  all processors must request or all must pass NULL
1281 -  ksp - the array of KSP contexts
1282 
1283    Note:
1284    After PCASMGetSubKSP() the array of KSPes is not to be freed.
1285 
1286    You must call KSPSetUp() before calling PCASMGetSubKSP().
1287 
1288    Fortran note:
1289    The output argument 'ksp' must be an array of sufficient length or PETSC_NULL_KSP. The latter can be used to learn the necessary length.
1290 
1291    Level: advanced
1292 
1293 .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(),
1294           PCASMCreateSubdomains2D(),
1295 @*/
1296 PetscErrorCode  PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[])
1297 {
1298   PetscErrorCode ierr;
1299 
1300   PetscFunctionBegin;
1301   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1302   ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr);
1303   PetscFunctionReturn(0);
1304 }
1305 
1306 /* -------------------------------------------------------------------------------------*/
1307 /*MC
1308    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
1309            its own KSP object.
1310 
1311    Options Database Keys:
1312 +  -pc_asm_blocks <blks> - Sets total blocks
1313 .  -pc_asm_overlap <ovl> - Sets overlap
1314 .  -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict
1315 -  -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive
1316 
1317      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
1318       will get a different convergence rate due to the default option of -pc_asm_type restrict. Use
1319       -pc_asm_type basic to use the standard ASM.
1320 
1321    Notes:
1322     Each processor can have one or more blocks, but a block cannot be shared by more
1323      than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor.
1324 
1325      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1326         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
1327 
1328      To set the options on the solvers separate for each block call PCASMGetSubKSP()
1329          and set the options directly on the resulting KSP object (you can access its PC
1330          with KSPGetPC())
1331 
1332    Level: beginner
1333 
1334     References:
1335 +   1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
1336      Courant Institute, New York University Technical report
1337 -   2. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
1338     Cambridge University Press.
1339 
1340 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1341            PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType()
1342            PCASMSetTotalSubdomains(), PCSetModifySubMatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType
1343 
1344 M*/
1345 
1346 PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc)
1347 {
1348   PetscErrorCode ierr;
1349   PC_ASM         *osm;
1350 
1351   PetscFunctionBegin;
1352   ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr);
1353 
1354   osm->n                 = PETSC_DECIDE;
1355   osm->n_local           = 0;
1356   osm->n_local_true      = PETSC_DECIDE;
1357   osm->overlap           = 1;
1358   osm->ksp               = NULL;
1359   osm->restriction       = NULL;
1360   osm->lprolongation     = NULL;
1361   osm->lrestriction      = NULL;
1362   osm->x                 = NULL;
1363   osm->y                 = NULL;
1364   osm->is                = NULL;
1365   osm->is_local          = NULL;
1366   osm->mat               = NULL;
1367   osm->pmat              = NULL;
1368   osm->type              = PC_ASM_RESTRICT;
1369   osm->loctype           = PC_COMPOSITE_ADDITIVE;
1370   osm->sort_indices      = PETSC_TRUE;
1371   osm->dm_subdomains     = PETSC_FALSE;
1372   osm->sub_mat_type      = NULL;
1373 
1374   pc->data                 = (void*)osm;
1375   pc->ops->apply           = PCApply_ASM;
1376   pc->ops->matapply        = PCMatApply_ASM;
1377   pc->ops->applytranspose  = PCApplyTranspose_ASM;
1378   pc->ops->setup           = PCSetUp_ASM;
1379   pc->ops->reset           = PCReset_ASM;
1380   pc->ops->destroy         = PCDestroy_ASM;
1381   pc->ops->setfromoptions  = PCSetFromOptions_ASM;
1382   pc->ops->setuponblocks   = PCSetUpOnBlocks_ASM;
1383   pc->ops->view            = PCView_ASM;
1384   pc->ops->applyrichardson = NULL;
1385 
1386   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr);
1387   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr);
1388   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);CHKERRQ(ierr);
1389   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);CHKERRQ(ierr);
1390   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);CHKERRQ(ierr);
1391   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);CHKERRQ(ierr);
1392   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);CHKERRQ(ierr);
1393   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);CHKERRQ(ierr);
1394   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);CHKERRQ(ierr);
1395   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);CHKERRQ(ierr);
1396   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);CHKERRQ(ierr);
1397   PetscFunctionReturn(0);
1398 }
1399 
1400 /*@C
1401    PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
1402    preconditioner for any problem on a general grid.
1403 
1404    Collective
1405 
1406    Input Parameters:
1407 +  A - The global matrix operator
1408 -  n - the number of local blocks
1409 
1410    Output Parameters:
1411 .  outis - the array of index sets defining the subdomains
1412 
1413    Level: advanced
1414 
1415    Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap
1416     from these if you use PCASMSetLocalSubdomains()
1417 
1418     In the Fortran version you must provide the array outis[] already allocated of length n.
1419 
1420 .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains()
1421 @*/
1422 PetscErrorCode  PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[])
1423 {
1424   MatPartitioning mpart;
1425   const char      *prefix;
1426   PetscInt        i,j,rstart,rend,bs;
1427   PetscBool       hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE;
1428   Mat             Ad     = NULL, adj;
1429   IS              ispart,isnumb,*is;
1430   PetscErrorCode  ierr;
1431 
1432   PetscFunctionBegin;
1433   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1434   PetscValidPointer(outis,3);
1435   if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);
1436 
1437   /* Get prefix, row distribution, and block size */
1438   ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr);
1439   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
1440   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1441   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);
1442 
1443   /* Get diagonal block from matrix if possible */
1444   ierr = MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr);
1445   if (hasop) {
1446     ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr);
1447   }
1448   if (Ad) {
1449     ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr);
1450     if (!isbaij) {ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);}
1451   }
1452   if (Ad && n > 1) {
1453     PetscBool match,done;
1454     /* Try to setup a good matrix partitioning if available */
1455     ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr);
1456     ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr);
1457     ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr);
1458     ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr);
1459     if (!match) {
1460       ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr);
1461     }
1462     if (!match) { /* assume a "good" partitioner is available */
1463       PetscInt       na;
1464       const PetscInt *ia,*ja;
1465       ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1466       if (done) {
1467         /* Build adjacency matrix by hand. Unfortunately a call to
1468            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
1469            remove the block-aij structure and we cannot expect
1470            MatPartitioning to split vertices as we need */
1471         PetscInt       i,j,len,nnz,cnt,*iia=NULL,*jja=NULL;
1472         const PetscInt *row;
1473         nnz = 0;
1474         for (i=0; i<na; i++) { /* count number of nonzeros */
1475           len = ia[i+1] - ia[i];
1476           row = ja + ia[i];
1477           for (j=0; j<len; j++) {
1478             if (row[j] == i) { /* don't count diagonal */
1479               len--; break;
1480             }
1481           }
1482           nnz += len;
1483         }
1484         ierr   = PetscMalloc1(na+1,&iia);CHKERRQ(ierr);
1485         ierr   = PetscMalloc1(nnz,&jja);CHKERRQ(ierr);
1486         nnz    = 0;
1487         iia[0] = 0;
1488         for (i=0; i<na; i++) { /* fill adjacency */
1489           cnt = 0;
1490           len = ia[i+1] - ia[i];
1491           row = ja + ia[i];
1492           for (j=0; j<len; j++) {
1493             if (row[j] != i) { /* if not diagonal */
1494               jja[nnz+cnt++] = row[j];
1495             }
1496           }
1497           nnz     += cnt;
1498           iia[i+1] = nnz;
1499         }
1500         /* Partitioning of the adjacency matrix */
1501         ierr      = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr);
1502         ierr      = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr);
1503         ierr      = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr);
1504         ierr      = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr);
1505         ierr      = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr);
1506         ierr      = MatDestroy(&adj);CHKERRQ(ierr);
1507         foundpart = PETSC_TRUE;
1508       }
1509       ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr);
1510     }
1511     ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr);
1512   }
1513 
1514   ierr   = PetscMalloc1(n,&is);CHKERRQ(ierr);
1515   *outis = is;
1516 
1517   if (!foundpart) {
1518 
1519     /* Partitioning by contiguous chunks of rows */
1520 
1521     PetscInt mbs   = (rend-rstart)/bs;
1522     PetscInt start = rstart;
1523     for (i=0; i<n; i++) {
1524       PetscInt count = (mbs/n + ((mbs % n) > i)) * bs;
1525       ierr   = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr);
1526       start += count;
1527     }
1528 
1529   } else {
1530 
1531     /* Partitioning by adjacency of diagonal block  */
1532 
1533     const PetscInt *numbering;
1534     PetscInt       *count,nidx,*indices,*newidx,start=0;
1535     /* Get node count in each partition */
1536     ierr = PetscMalloc1(n,&count);CHKERRQ(ierr);
1537     ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr);
1538     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1539       for (i=0; i<n; i++) count[i] *= bs;
1540     }
1541     /* Build indices from node numbering */
1542     ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr);
1543     ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr);
1544     for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */
1545     ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr);
1546     ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr);
1547     ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr);
1548     if (isbaij && bs > 1) { /* adjust for the block-aij case */
1549       ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr);
1550       for (i=0; i<nidx; i++) {
1551         for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j;
1552       }
1553       ierr    = PetscFree(indices);CHKERRQ(ierr);
1554       nidx   *= bs;
1555       indices = newidx;
1556     }
1557     /* Shift to get global indices */
1558     for (i=0; i<nidx; i++) indices[i] += rstart;
1559 
1560     /* Build the index sets for each block */
1561     for (i=0; i<n; i++) {
1562       ierr   = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr);
1563       ierr   = ISSort(is[i]);CHKERRQ(ierr);
1564       start += count[i];
1565     }
1566 
1567     ierr = PetscFree(count);CHKERRQ(ierr);
1568     ierr = PetscFree(indices);CHKERRQ(ierr);
1569     ierr = ISDestroy(&isnumb);CHKERRQ(ierr);
1570     ierr = ISDestroy(&ispart);CHKERRQ(ierr);
1571 
1572   }
1573   PetscFunctionReturn(0);
1574 }
1575 
1576 /*@C
1577    PCASMDestroySubdomains - Destroys the index sets created with
1578    PCASMCreateSubdomains(). Should be called after setting subdomains
1579    with PCASMSetLocalSubdomains().
1580 
1581    Collective
1582 
1583    Input Parameters:
1584 +  n - the number of index sets
1585 .  is - the array of index sets
1586 -  is_local - the array of local index sets, can be NULL
1587 
1588    Level: advanced
1589 
1590 .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains()
1591 @*/
1592 PetscErrorCode  PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[])
1593 {
1594   PetscInt       i;
1595   PetscErrorCode ierr;
1596 
1597   PetscFunctionBegin;
1598   if (n <= 0) PetscFunctionReturn(0);
1599   if (is) {
1600     PetscValidPointer(is,2);
1601     for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); }
1602     ierr = PetscFree(is);CHKERRQ(ierr);
1603   }
1604   if (is_local) {
1605     PetscValidPointer(is_local,3);
1606     for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); }
1607     ierr = PetscFree(is_local);CHKERRQ(ierr);
1608   }
1609   PetscFunctionReturn(0);
1610 }
1611 
1612 /*@
1613    PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1614    preconditioner for a two-dimensional problem on a regular grid.
1615 
1616    Not Collective
1617 
1618    Input Parameters:
1619 +  m   - the number of mesh points in the x direction
1620 .  n   - the number of mesh points in the y direction
1621 .  M   - the number of subdomains in the x direction
1622 .  N   - the number of subdomains in the y direction
1623 .  dof - degrees of freedom per node
1624 -  overlap - overlap in mesh lines
1625 
1626    Output Parameters:
1627 +  Nsub - the number of subdomains created
1628 .  is - array of index sets defining overlapping (if overlap > 0) subdomains
1629 -  is_local - array of index sets defining non-overlapping subdomains
1630 
1631    Note:
1632    Presently PCAMSCreateSubdomains2d() is valid only for sequential
1633    preconditioners.  More general related routines are
1634    PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().
1635 
1636    Level: advanced
1637 
1638 .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(),
1639           PCASMSetOverlap()
1640 @*/
1641 PetscErrorCode  PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local)
1642 {
1643   PetscInt       i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer;
1644   PetscErrorCode ierr;
1645   PetscInt       nidx,*idx,loc,ii,jj,count;
1646 
1647   PetscFunctionBegin;
1648   if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," ");
1649 
1650   *Nsub     = N*M;
1651   ierr      = PetscMalloc1(*Nsub,is);CHKERRQ(ierr);
1652   ierr      = PetscMalloc1(*Nsub,is_local);CHKERRQ(ierr);
1653   ystart    = 0;
1654   loc_outer = 0;
1655   for (i=0; i<N; i++) {
1656     height = n/N + ((n % N) > i); /* height of subdomain */
1657     if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n");
1658     yleft  = ystart - overlap; if (yleft < 0) yleft = 0;
1659     yright = ystart + height + overlap; if (yright > n) yright = n;
1660     xstart = 0;
1661     for (j=0; j<M; j++) {
1662       width = m/M + ((m % M) > j); /* width of subdomain */
1663       if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m");
1664       xleft  = xstart - overlap; if (xleft < 0) xleft = 0;
1665       xright = xstart + width + overlap; if (xright > m) xright = m;
1666       nidx   = (xright - xleft)*(yright - yleft);
1667       ierr   = PetscMalloc1(nidx,&idx);CHKERRQ(ierr);
1668       loc    = 0;
1669       for (ii=yleft; ii<yright; ii++) {
1670         count = m*ii + xleft;
1671         for (jj=xleft; jj<xright; jj++) idx[loc++] = count++;
1672       }
1673       ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr);
1674       if (overlap == 0) {
1675         ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr);
1676 
1677         (*is_local)[loc_outer] = (*is)[loc_outer];
1678       } else {
1679         for (loc=0,ii=ystart; ii<ystart+height; ii++) {
1680           for (jj=xstart; jj<xstart+width; jj++) {
1681             idx[loc++] = m*ii + jj;
1682           }
1683         }
1684         ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr);
1685       }
1686       ierr    = PetscFree(idx);CHKERRQ(ierr);
1687       xstart += width;
1688       loc_outer++;
1689     }
1690     ystart += height;
1691   }
1692   for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); }
1693   PetscFunctionReturn(0);
1694 }
1695 
1696 /*@C
1697     PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
1698     only) for the additive Schwarz preconditioner.
1699 
1700     Not Collective
1701 
1702     Input Parameter:
1703 .   pc - the preconditioner context
1704 
1705     Output Parameters:
1706 +   n - if requested, the number of subdomains for this processor (default value = 1)
1707 .   is - if requested, the index sets that define the subdomains for this processor
1708 -   is_local - if requested, the index sets that define the local part of the subdomains for this processor (can be NULL)
1709 
1710     Notes:
1711     The IS numbering is in the parallel, global numbering of the vector.
1712 
1713     Level: advanced
1714 
1715 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1716           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices()
1717 @*/
1718 PetscErrorCode  PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[])
1719 {
1720   PC_ASM         *osm = (PC_ASM*)pc->data;
1721   PetscErrorCode ierr;
1722   PetscBool      match;
1723 
1724   PetscFunctionBegin;
1725   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1726   if (n) PetscValidIntPointer(n,2);
1727   if (is) PetscValidPointer(is,3);
1728   if (is_local) PetscValidPointer(is_local,4);
1729   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1730   if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM");
1731   if (n) *n = osm->n_local_true;
1732   if (is) *is = osm->is;
1733   if (is_local) *is_local = osm->is_local;
1734   PetscFunctionReturn(0);
1735 }
1736 
1737 /*@C
1738     PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
1739     only) for the additive Schwarz preconditioner.
1740 
1741     Not Collective
1742 
1743     Input Parameter:
1744 .   pc - the preconditioner context
1745 
1746     Output Parameters:
1747 +   n - if requested, the number of matrices for this processor (default value = 1)
1748 -   mat - if requested, the matrices
1749 
1750     Level: advanced
1751 
1752     Notes:
1753     Call after PCSetUp() (or KSPSetUp()) but before PCApply() and before PCSetUpOnBlocks())
1754 
1755            Usually one would use PCSetModifySubMatrices() to change the submatrices in building the preconditioner.
1756 
1757 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(),
1758           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubMatrices()
1759 @*/
1760 PetscErrorCode  PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[])
1761 {
1762   PC_ASM         *osm;
1763   PetscErrorCode ierr;
1764   PetscBool      match;
1765 
1766   PetscFunctionBegin;
1767   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1768   if (n) PetscValidIntPointer(n,2);
1769   if (mat) PetscValidPointer(mat,3);
1770   if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp().");
1771   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1772   if (!match) {
1773     if (n) *n = 0;
1774     if (mat) *mat = NULL;
1775   } else {
1776     osm = (PC_ASM*)pc->data;
1777     if (n) *n = osm->n_local_true;
1778     if (mat) *mat = osm->pmat;
1779   }
1780   PetscFunctionReturn(0);
1781 }
1782 
1783 /*@
1784     PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1785 
1786     Logically Collective
1787 
1788     Input Parameters:
1789 +   pc  - the preconditioner
1790 -   flg - boolean indicating whether to use subdomains defined by the DM
1791 
1792     Options Database Key:
1793 .   -pc_asm_dm_subdomains
1794 
1795     Level: intermediate
1796 
1797     Notes:
1798     PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(),
1799     so setting either of the first two effectively turns the latter off.
1800 
1801 .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1802           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1803 @*/
1804 PetscErrorCode  PCASMSetDMSubdomains(PC pc,PetscBool flg)
1805 {
1806   PC_ASM         *osm = (PC_ASM*)pc->data;
1807   PetscErrorCode ierr;
1808   PetscBool      match;
1809 
1810   PetscFunctionBegin;
1811   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1812   PetscValidLogicalCollectiveBool(pc,flg,2);
1813   if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC.");
1814   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1815   if (match) {
1816     osm->dm_subdomains = flg;
1817   }
1818   PetscFunctionReturn(0);
1819 }
1820 
1821 /*@
1822     PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1823     Not Collective
1824 
1825     Input Parameter:
1826 .   pc  - the preconditioner
1827 
1828     Output Parameter:
1829 .   flg - boolean indicating whether to use subdomains defined by the DM
1830 
1831     Level: intermediate
1832 
1833 .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap()
1834           PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains()
1835 @*/
1836 PetscErrorCode  PCASMGetDMSubdomains(PC pc,PetscBool* flg)
1837 {
1838   PC_ASM         *osm = (PC_ASM*)pc->data;
1839   PetscErrorCode ierr;
1840   PetscBool      match;
1841 
1842   PetscFunctionBegin;
1843   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1844   PetscValidBoolPointer(flg,2);
1845   ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr);
1846   if (match) *flg = osm->dm_subdomains;
1847   else *flg = PETSC_FALSE;
1848   PetscFunctionReturn(0);
1849 }
1850 
1851 /*@
1852      PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string.
1853 
1854    Not Collective
1855 
1856    Input Parameter:
1857 .  pc - the PC
1858 
1859    Output Parameter:
1860 .  -pc_asm_sub_mat_type - name of matrix type
1861 
1862    Level: advanced
1863 
1864 .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
1865 @*/
1866 PetscErrorCode  PCASMGetSubMatType(PC pc,MatType *sub_mat_type)
1867 {
1868   PetscErrorCode ierr;
1869 
1870   PetscFunctionBegin;
1871   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1872   ierr = PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));CHKERRQ(ierr);
1873   PetscFunctionReturn(0);
1874 }
1875 
1876 /*@
1877      PCASMSetSubMatType - Set the type of matrix used for ASM subsolves
1878 
1879    Collective on Mat
1880 
1881    Input Parameters:
1882 +  pc             - the PC object
1883 -  sub_mat_type   - matrix type
1884 
1885    Options Database Key:
1886 .  -pc_asm_sub_mat_type  <sub_mat_type> - Sets the matrix type used for subsolves, for example, seqaijviennacl. If you specify a base name like aijviennacl, the corresponding sequential type is assumed.
1887 
1888    Notes:
1889    See "${PETSC_DIR}/include/petscmat.h" for available types
1890 
1891   Level: advanced
1892 
1893 .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat
1894 @*/
1895 PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type)
1896 {
1897   PetscErrorCode ierr;
1898 
1899   PetscFunctionBegin;
1900   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1901   ierr = PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));CHKERRQ(ierr);
1902   PetscFunctionReturn(0);
1903 }
1904