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