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