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