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