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