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