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