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