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