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