Actual source code: wb.c

  1: #define PETSCKSP_DLL


 4:  #include petscpc.h
 5:  #include petscmg.h
 6:  #include petscda.h
 7:  #include ../src/ksp/pc/impls/mg/mgimpl.h

  9: typedef struct {
 10:   DA           da;
 11:   PCExoticType type;
 12:   Mat          P;            /* the constructed interpolation matrix */
 13:   PetscTruth   directSolve;  /* use direct LU factorization to construct interpolation */
 14:   KSP          ksp;
 15: } PC_Exotic;

 17: const char *PCExoticTypes[] = {"face","wirebasket","PCExoticType","PC_Exotic",0};


 22: /*
 23:       DAGetWireBasketInterpolation - Gets the interpolation for a wirebasket based coarse space

 25: */
 26: PetscErrorCode DAGetWireBasketInterpolation(PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P)
 27: {
 28:   DA                     da = exotic->da;
 29:   PetscErrorCode         ierr;
 30:   PetscInt               dim,i,j,k,m,n,p,dof,Nint,Nface,Nwire,Nsurf,*Iint,*Isurf,cint = 0,csurf = 0,istart,jstart,kstart,*II,N,c = 0;
 31:   PetscInt               mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[26],*globals,Ng,*IIint,*IIsurf;
 32:   Mat                    Xint, Xsurf,Xint_tmp;
 33:   IS                     isint,issurf,is,row,col;
 34:   ISLocalToGlobalMapping ltg;
 35:   MPI_Comm               comm;
 36:   Mat                    A,Aii,Ais,Asi,*Aholder,iAii;
 37:   MatFactorInfo          info;
 38:   PetscScalar            *xsurf,*xint;
 39: #if defined(PETSC_USE_DEBUG_foo)
 40:   PetscScalar            tmp;
 41: #endif
 42:   PetscTable             ht;

 45:   DAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0);
 46:   if (dof != 1) SETERRQ(PETSC_ERR_SUP,"Only for single field problems");
 47:   if (dim != 3) SETERRQ(PETSC_ERR_SUP,"Only coded for 3d problems");
 48:   DAGetCorners(da,0,0,0,&m,&n,&p);
 49:   DAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);
 50:   istart = istart ? -1 : 0;
 51:   jstart = jstart ? -1 : 0;
 52:   kstart = kstart ? -1 : 0;

 54:   /* 
 55:     the columns of P are the interpolation of each coarse grid point (one for each vertex and edge) 
 56:     to all the local degrees of freedom (this includes the vertices, edges and faces).

 58:     Xint are the subset of the interpolation into the interior

 60:     Xface are the interpolation onto faces but not into the interior 

 62:     Xsurf are the interpolation onto the vertices and edges (the surfbasket) 
 63:                                         Xint
 64:     Symbolically one could write P = (  Xface  ) after interchanging the rows to match the natural ordering on the domain
 65:                                         Xsurf
 66:   */
 67:   N     = (m - istart)*(n - jstart)*(p - kstart);
 68:   Nint  = (m-2-istart)*(n-2-jstart)*(p-2-kstart);
 69:   Nface = 2*( (m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart) );
 70:   Nwire = 4*( (m-2-istart) + (n-2-jstart) + (p-2-kstart) ) + 8;
 71:   Nsurf = Nface + Nwire;
 72:   MatCreateSeqDense(MPI_COMM_SELF,Nint,26,PETSC_NULL,&Xint);
 73:   MatCreateSeqDense(MPI_COMM_SELF,Nsurf,26,PETSC_NULL,&Xsurf);
 74:   MatGetArray(Xsurf,&xsurf);

 76:   /*
 77:      Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of 
 78:      Xsurf will be all zero (thus making the coarse matrix singular). 
 79:   */
 80:   if (m-istart < 3) SETERRQ(PETSC_ERR_SUP,"Number of grid points per process in X direction must be at least 3");
 81:   if (n-jstart < 3) SETERRQ(PETSC_ERR_SUP,"Number of grid points per process in Y direction must be at least 3");
 82:   if (p-kstart < 3) SETERRQ(PETSC_ERR_SUP,"Number of grid points per process in Z direction must be at least 3");

 84:   cnt = 0;
 85:   xsurf[cnt++] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + Nsurf] = 1;} xsurf[cnt++ + 2*Nsurf] = 1;
 86:   for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 3*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 4*Nsurf] = 1;} xsurf[cnt++ + 5*Nsurf] = 1;}
 87:   xsurf[cnt++ + 6*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 7*Nsurf] = 1;} xsurf[cnt++ + 8*Nsurf] = 1;
 88:   for (k=1;k<p-1-kstart;k++) {
 89:     xsurf[cnt++ + 9*Nsurf] = 1;  for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 10*Nsurf] = 1;}  xsurf[cnt++ + 11*Nsurf] = 1;
 90:     for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 12*Nsurf] = 1; /* these are the interior nodes */ xsurf[cnt++ + 13*Nsurf] = 1;}
 91:     xsurf[cnt++ + 14*Nsurf] = 1;  for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 15*Nsurf] = 1;} xsurf[cnt++ + 16*Nsurf] = 1;
 92:   }
 93:   xsurf[cnt++ + 17*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 18*Nsurf] = 1;} xsurf[cnt++ + 19*Nsurf] = 1;
 94:   for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 20*Nsurf] = 1;  for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 21*Nsurf] = 1;} xsurf[cnt++ + 22*Nsurf] = 1;}
 95:   xsurf[cnt++ + 23*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 24*Nsurf] = 1;} xsurf[cnt++ + 25*Nsurf] = 1;

 97:   /* interpolations only sum to 1 when using direct solver */
 98: #if defined(PETSC_USE_DEBUG_foo)
 99:   for (i=0; i<Nsurf; i++) {
100:     tmp = 0.0;
101:     for (j=0; j<26; j++) {
102:       tmp += xsurf[i+j*Nsurf];
103:     }
104:     if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong Xsurf interpolation at i %D value %G",i,PetscAbsScalar(tmp));
105:   }
106: #endif
107:   MatRestoreArray(Xsurf,&xsurf);
108:   /* MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);*/


111:   /* 
112:        I are the indices for all the needed vertices (in global numbering)
113:        Iint are the indices for the interior values, I surf for the surface values
114:             (This is just for the part of the global matrix obtained with MatGetSubMatrix(), it 
115:              is NOT the local DA ordering.)
116:        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
117:   */
118: #define Endpoint(a,start,b) (a == 0 || a == (b-1-start))
119:   PetscMalloc3(N,PetscInt,&II,Nint,PetscInt,&Iint,Nsurf,PetscInt,&Isurf);
120:   PetscMalloc2(Nint,PetscInt,&IIint,Nsurf,PetscInt,&IIsurf);
121:   for (k=0; k<p-kstart; k++) {
122:     for (j=0; j<n-jstart; j++) {
123:       for (i=0; i<m-istart; i++) {
124:         II[c++] = i + j*mwidth + k*mwidth*nwidth;

126:         if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) {
127:           IIint[cint]  = i + j*mwidth + k*mwidth*nwidth;
128:           Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
129:         } else {
130:           IIsurf[csurf]  = i + j*mwidth + k*mwidth*nwidth;
131:           Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
132:         }
133:       }
134:     }
135:   }
136:   if (c != N) SETERRQ(PETSC_ERR_PLIB,"c != N");
137:   if (cint != Nint) SETERRQ(PETSC_ERR_PLIB,"cint != Nint");
138:   if (csurf != Nsurf) SETERRQ(PETSC_ERR_PLIB,"csurf != Nsurf");
139:   DAGetISLocalToGlobalMapping(da,&ltg);
140:   ISLocalToGlobalMappingApply(ltg,N,II,II);
141:   ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);
142:   ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);
143:   PetscObjectGetComm((PetscObject)da,&comm);
144:   ISCreateGeneral(comm,N,II,&is);
145:   ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,&isint);
146:   ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,&issurf);
147:   PetscFree3(II,Iint,Isurf);

149:   MatGetSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);
150:   A    = *Aholder;
151:   PetscFree(Aholder);

153:   MatGetSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);
154:   MatGetSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);
155:   MatGetSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);

157:   /* 
158:      Solve for the interpolation onto the interior Xint
159:   */
160:   MatDuplicate(Xint,MAT_DO_NOT_COPY_VALUES,&Xint_tmp);
161:   MatMatMult(Ais,Xsurf,MAT_REUSE_MATRIX,PETSC_DETERMINE,&Xint_tmp);
162:   MatScale(Xint_tmp,-1.0);
163:   if (exotic->directSolve) {
164:     MatGetFactor(Aii,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&iAii);
165:     MatFactorInfoInitialize(&info);
166:     MatGetOrdering(Aii,MATORDERING_ND,&row,&col);
167:     MatLUFactorSymbolic(iAii,Aii,row,col,&info);
168:     ISDestroy(row);
169:     ISDestroy(col);
170:     MatLUFactorNumeric(iAii,Aii,&info);
171:     MatMatSolve(iAii,Xint_tmp,Xint);
172:     MatDestroy(iAii);
173:   } else {
174:     Vec         b,x;
175:     PetscScalar *xint_tmp;

177:     MatGetArray(Xint,&xint);
178:     VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&x);
179:     MatGetArray(Xint_tmp,&xint_tmp);
180:     VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&b);
181:     KSPSetOperators(exotic->ksp,Aii,Aii,SAME_NONZERO_PATTERN);
182:     for (i=0; i<26; i++) {
183:       VecPlaceArray(x,xint+i*Nint);
184:       VecPlaceArray(b,xint_tmp+i*Nint);
185:       KSPSolve(exotic->ksp,b,x);
186:       VecResetArray(x);
187:       VecResetArray(b);
188:     }
189:     MatRestoreArray(Xint,&xint);
190:     MatRestoreArray(Xint_tmp,&xint_tmp);
191:     VecDestroy(x);
192:     VecDestroy(b);
193:   }
194:   MatDestroy(Xint_tmp);

196: #if defined(PETSC_USE_DEBUG_foo)
197:   MatGetArray(Xint,&xint);
198:   for (i=0; i<Nint; i++) {
199:     tmp = 0.0;
200:     for (j=0; j<26; j++) {
201:       tmp += xint[i+j*Nint];
202:     }
203:     if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong Xint interpolation at i %D value %G",i,PetscAbsScalar(tmp));
204:   }
205:   MatRestoreArray(Xint,&xint);
206:   /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD); */
207: #endif


210:   /*         total vertices             total faces                                  total edges */
211:   Ntotal = (mp + 1)*(np + 1)*(pp + 1) + mp*np*(pp+1) + mp*pp*(np+1) + np*pp*(mp+1) + mp*(np+1)*(pp+1) + np*(mp+1)*(pp+1) +  pp*(mp+1)*(np+1);

213:   /*
214:       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points 
215:   */
216:   cnt = 0;
217:   gl[cnt++] = 0;  { gl[cnt++] = 1;} gl[cnt++] = m-istart-1;
218:   { gl[cnt++] = mwidth;  { gl[cnt++] = mwidth+1;} gl[cnt++] = mwidth + m-istart-1;}
219:   gl[cnt++] = mwidth*(n-jstart-1);  { gl[cnt++] = mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*(n-jstart-1) + m-istart-1;
220:   {
221:     gl[cnt++] = mwidth*nwidth;  { gl[cnt++] = mwidth*nwidth+1;}  gl[cnt++] = mwidth*nwidth+ m-istart-1;
222:     { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;}
223:     gl[cnt++] = mwidth*nwidth+ mwidth*(n-jstart-1);   { gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1) + m-istart-1;
224:   }
225:   gl[cnt++] = mwidth*nwidth*(p-kstart-1); { gl[cnt++] = mwidth*nwidth*(p-kstart-1)+1;} gl[cnt++] = mwidth*nwidth*(p-kstart-1) +  m-istart-1;
226:   { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth;   { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth+1;} gl[cnt++] = mwidth*nwidth*(p-kstart-1)+mwidth+m-istart-1;}
227:   gl[cnt++] = mwidth*nwidth*(p-kstart-1) +  mwidth*(n-jstart-1);  { gl[cnt++] = mwidth*nwidth*(p-kstart-1)+ mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth*(n-jstart-1) + m-istart-1;

229:   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
230:   /* convert that to global numbering and get them on all processes */
231:   ISLocalToGlobalMappingApply(ltg,26,gl,gl);
232:   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
233:   PetscMalloc(26*mp*np*pp*sizeof(PetscInt),&globals);
234:   MPI_Allgather(gl,26,MPIU_INT,globals,26,MPIU_INT,((PetscObject)da)->comm);

236:   /* Number the coarse grid points from 0 to Ntotal */
237:   PetscTableCreate(Ntotal/3,&ht);
238:   for (i=0; i<26*mp*np*pp; i++){
239:     PetscTableAddCount(ht,globals[i]+1);
240:   }
241:   PetscTableGetCount(ht,&cnt);
242:   if (cnt != Ntotal) SETERRQ2(PETSC_ERR_PLIB,"Hash table size %D not equal to total number coarse grid points %D",cnt,Ntotal);
243:   PetscFree(globals);
244:   for (i=0; i<26; i++) {
245:     PetscTableFind(ht,gl[i]+1,&gl[i]);
246:     gl[i]--;
247:   }
248:   PetscTableDestroy(ht);
249:   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */

251:   /* construct global interpolation matrix */
252:   MatGetLocalSize(Aglobal,&Ng,PETSC_NULL);
253:   if (reuse == MAT_INITIAL_MATRIX) {
254:     MatCreateMPIAIJ(((PetscObject)da)->comm,Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,PETSC_NULL,Nint+Nsurf,PETSC_NULL,P);
255:   } else {
256:     MatZeroEntries(*P);
257:   }
258:   MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);
259:   MatGetArray(Xint,&xint);
260:   MatSetValues(*P,Nint,IIint,26,gl,xint,INSERT_VALUES);
261:   MatRestoreArray(Xint,&xint);
262:   MatGetArray(Xsurf,&xsurf);
263:   MatSetValues(*P,Nsurf,IIsurf,26,gl,xsurf,INSERT_VALUES);
264:   MatRestoreArray(Xsurf,&xsurf);
265:   MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);
266:   MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);
267:   PetscFree2(IIint,IIsurf);

269: #if defined(PETSC_USE_DEBUG_foo)
270:   {
271:     Vec         x,y;
272:     PetscScalar *yy;
273:     VecCreateMPI(((PetscObject)da)->comm,Ng,PETSC_DETERMINE,&y);
274:     VecCreateMPI(((PetscObject)da)->comm,PETSC_DETERMINE,Ntotal,&x);
275:     VecSet(x,1.0);
276:     MatMult(*P,x,y);
277:     VecGetArray(y,&yy);
278:     for (i=0; i<Ng; i++) {
279:       if (PetscAbsScalar(yy[i]-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong p interpolation at i %D value %G",i,PetscAbsScalar(yy[i]));
280:     }
281:     VecRestoreArray(y,&yy);
282:     VecDestroy(x);
283:     VecDestroy(y);
284:   }
285: #endif
286: 
287:   MatDestroy(Aii);
288:   MatDestroy(Ais);
289:   MatDestroy(Asi);
290:   MatDestroy(A);
291:   ISDestroy(is);
292:   ISDestroy(isint);
293:   ISDestroy(issurf);
294:   MatDestroy(Xint);
295:   MatDestroy(Xsurf);
296:   return(0);
297: }

301: /*
302:       DAGetFaceInterpolation - Gets the interpolation for a face based coarse space

304: */
305: PetscErrorCode DAGetFaceInterpolation(PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P)
306: {
307:   DA                     da = exotic->da;
308:   PetscErrorCode         ierr;
309:   PetscInt               dim,i,j,k,m,n,p,dof,Nint,Nface,Nwire,Nsurf,*Iint,*Isurf,cint = 0,csurf = 0,istart,jstart,kstart,*II,N,c = 0;
310:   PetscInt               mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[6],*globals,Ng,*IIint,*IIsurf;
311:   Mat                    Xint, Xsurf,Xint_tmp;
312:   IS                     isint,issurf,is,row,col;
313:   ISLocalToGlobalMapping ltg;
314:   MPI_Comm               comm;
315:   Mat                    A,Aii,Ais,Asi,*Aholder,iAii;
316:   MatFactorInfo          info;
317:   PetscScalar            *xsurf,*xint;
318: #if defined(PETSC_USE_DEBUG_foo)
319:   PetscScalar            tmp;
320: #endif
321:   PetscTable             ht;

324:   DAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0);
325:   if (dof != 1) SETERRQ(PETSC_ERR_SUP,"Only for single field problems");
326:   if (dim != 3) SETERRQ(PETSC_ERR_SUP,"Only coded for 3d problems");
327:   DAGetCorners(da,0,0,0,&m,&n,&p);
328:   DAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);
329:   istart = istart ? -1 : 0;
330:   jstart = jstart ? -1 : 0;
331:   kstart = kstart ? -1 : 0;

333:   /* 
334:     the columns of P are the interpolation of each coarse grid point (one for each vertex and edge) 
335:     to all the local degrees of freedom (this includes the vertices, edges and faces).

337:     Xint are the subset of the interpolation into the interior

339:     Xface are the interpolation onto faces but not into the interior 

341:     Xsurf are the interpolation onto the vertices and edges (the surfbasket) 
342:                                         Xint
343:     Symbolically one could write P = (  Xface  ) after interchanging the rows to match the natural ordering on the domain
344:                                         Xsurf
345:   */
346:   N     = (m - istart)*(n - jstart)*(p - kstart);
347:   Nint  = (m-2-istart)*(n-2-jstart)*(p-2-kstart);
348:   Nface = 2*( (m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart) );
349:   Nwire = 4*( (m-2-istart) + (n-2-jstart) + (p-2-kstart) ) + 8;
350:   Nsurf = Nface + Nwire;
351:   MatCreateSeqDense(MPI_COMM_SELF,Nint,6,PETSC_NULL,&Xint);
352:   MatCreateSeqDense(MPI_COMM_SELF,Nsurf,6,PETSC_NULL,&Xsurf);
353:   MatGetArray(Xsurf,&xsurf);

355:   /*
356:      Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of 
357:      Xsurf will be all zero (thus making the coarse matrix singular). 
358:   */
359:   if (m-istart < 3) SETERRQ(PETSC_ERR_SUP,"Number of grid points per process in X direction must be at least 3");
360:   if (n-jstart < 3) SETERRQ(PETSC_ERR_SUP,"Number of grid points per process in Y direction must be at least 3");
361:   if (p-kstart < 3) SETERRQ(PETSC_ERR_SUP,"Number of grid points per process in Z direction must be at least 3");

363:   cnt = 0;
364:   for (j=1;j<n-1-jstart;j++) {  for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 0*Nsurf] = 1;} }
365:    for (k=1;k<p-1-kstart;k++) {
366:     for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 1*Nsurf] = 1;}
367:     for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 2*Nsurf] = 1; /* these are the interior nodes */ xsurf[cnt++ + 3*Nsurf] = 1;}
368:     for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 4*Nsurf] = 1;}
369:   }
370:   for (j=1;j<n-1-jstart;j++) {for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 5*Nsurf] = 1;} }

372: #if defined(PETSC_USE_DEBUG_foo)
373:   for (i=0; i<Nsurf; i++) {
374:     tmp = 0.0;
375:     for (j=0; j<6; j++) {
376:       tmp += xsurf[i+j*Nsurf];
377:     }
378:     if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong Xsurf interpolation at i %D value %G",i,PetscAbsScalar(tmp));
379:   }
380: #endif
381:   MatRestoreArray(Xsurf,&xsurf);
382:   /* MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);*/


385:   /* 
386:        I are the indices for all the needed vertices (in global numbering)
387:        Iint are the indices for the interior values, I surf for the surface values
388:             (This is just for the part of the global matrix obtained with MatGetSubMatrix(), it 
389:              is NOT the local DA ordering.)
390:        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
391:   */
392: #define Endpoint(a,start,b) (a == 0 || a == (b-1-start))
393:   PetscMalloc3(N,PetscInt,&II,Nint,PetscInt,&Iint,Nsurf,PetscInt,&Isurf);
394:   PetscMalloc2(Nint,PetscInt,&IIint,Nsurf,PetscInt,&IIsurf);
395:   for (k=0; k<p-kstart; k++) {
396:     for (j=0; j<n-jstart; j++) {
397:       for (i=0; i<m-istart; i++) {
398:         II[c++] = i + j*mwidth + k*mwidth*nwidth;

400:         if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) {
401:           IIint[cint]  = i + j*mwidth + k*mwidth*nwidth;
402:           Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
403:         } else {
404:           IIsurf[csurf]  = i + j*mwidth + k*mwidth*nwidth;
405:           Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
406:         }
407:       }
408:     }
409:   }
410:   if (c != N) SETERRQ(PETSC_ERR_PLIB,"c != N");
411:   if (cint != Nint) SETERRQ(PETSC_ERR_PLIB,"cint != Nint");
412:   if (csurf != Nsurf) SETERRQ(PETSC_ERR_PLIB,"csurf != Nsurf");
413:   DAGetISLocalToGlobalMapping(da,&ltg);
414:   ISLocalToGlobalMappingApply(ltg,N,II,II);
415:   ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);
416:   ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);
417:   PetscObjectGetComm((PetscObject)da,&comm);
418:   ISCreateGeneral(comm,N,II,&is);
419:   ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,&isint);
420:   ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,&issurf);
421:   PetscFree3(II,Iint,Isurf);

423:   MatGetSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);
424:   A    = *Aholder;
425:   PetscFree(Aholder);

427:   MatGetSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);
428:   MatGetSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);
429:   MatGetSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);

431:   /* 
432:      Solve for the interpolation onto the interior Xint
433:   */
434:   MatDuplicate(Xint,MAT_DO_NOT_COPY_VALUES,&Xint_tmp);
435:   MatMatMult(Ais,Xsurf,MAT_REUSE_MATRIX,PETSC_DETERMINE,&Xint_tmp);
436:   MatScale(Xint_tmp,-1.0);

438:   if (exotic->directSolve) {
439:     MatGetFactor(Aii,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&iAii);
440:     MatFactorInfoInitialize(&info);
441:     MatGetOrdering(Aii,MATORDERING_ND,&row,&col);
442:     MatLUFactorSymbolic(iAii,Aii,row,col,&info);
443:     ISDestroy(row);
444:     ISDestroy(col);
445:     MatLUFactorNumeric(iAii,Aii,&info);
446:     MatMatSolve(iAii,Xint_tmp,Xint);
447:     MatDestroy(iAii);
448:   } else {
449:     Vec         b,x;
450:     PetscScalar *xint_tmp;

452:     MatGetArray(Xint,&xint);
453:     VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&x);
454:     MatGetArray(Xint_tmp,&xint_tmp);
455:     VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&b);
456:     KSPSetOperators(exotic->ksp,Aii,Aii,SAME_NONZERO_PATTERN);
457:     for (i=0; i<6; i++) {
458:       VecPlaceArray(x,xint+i*Nint);
459:       VecPlaceArray(b,xint_tmp+i*Nint);
460:       KSPSolve(exotic->ksp,b,x);
461:       VecResetArray(x);
462:       VecResetArray(b);
463:     }
464:     MatRestoreArray(Xint,&xint);
465:     MatRestoreArray(Xint_tmp,&xint_tmp);
466:     VecDestroy(x);
467:     VecDestroy(b);
468:   }
469:   MatDestroy(Xint_tmp);

471: #if defined(PETSC_USE_DEBUG_foo)
472:   MatGetArray(Xint,&xint);
473:   for (i=0; i<Nint; i++) {
474:     tmp = 0.0;
475:     for (j=0; j<6; j++) {
476:       tmp += xint[i+j*Nint];
477:     }
478:     if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong Xint interpolation at i %D value %G",i,PetscAbsScalar(tmp));
479:   }
480:   MatRestoreArray(Xint,&xint);
481:   /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD); */
482: #endif


485:   /*         total faces    */
486:   Ntotal =  mp*np*(pp+1) + mp*pp*(np+1) + np*pp*(mp+1);

488:   /*
489:       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points 
490:   */
491:   cnt = 0;
492:   { gl[cnt++] = mwidth+1;}
493:   {
494:     { gl[cnt++] = mwidth*nwidth+1;}
495:     { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;}
496:     { gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1)+1;}
497:   }
498:   { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth+1;}

500:   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
501:   /* convert that to global numbering and get them on all processes */
502:   ISLocalToGlobalMappingApply(ltg,6,gl,gl);
503:   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
504:   PetscMalloc(6*mp*np*pp*sizeof(PetscInt),&globals);
505:   MPI_Allgather(gl,6,MPIU_INT,globals,6,MPIU_INT,((PetscObject)da)->comm);

507:   /* Number the coarse grid points from 0 to Ntotal */
508:   PetscTableCreate(Ntotal/3,&ht);
509:   for (i=0; i<6*mp*np*pp; i++){
510:     PetscTableAddCount(ht,globals[i]+1);
511:   }
512:   PetscTableGetCount(ht,&cnt);
513:   if (cnt != Ntotal) SETERRQ2(PETSC_ERR_PLIB,"Hash table size %D not equal to total number coarse grid points %D",cnt,Ntotal);
514:   PetscFree(globals);
515:   for (i=0; i<6; i++) {
516:     PetscTableFind(ht,gl[i]+1,&gl[i]);
517:     gl[i]--;
518:   }
519:   PetscTableDestroy(ht);
520:   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */

522:   /* construct global interpolation matrix */
523:   MatGetLocalSize(Aglobal,&Ng,PETSC_NULL);
524:   if (reuse == MAT_INITIAL_MATRIX) {
525:     MatCreateMPIAIJ(((PetscObject)da)->comm,Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,PETSC_NULL,Nint,PETSC_NULL,P);
526:   } else {
527:     MatZeroEntries(*P);
528:   }
529:   MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);
530:   MatGetArray(Xint,&xint);
531:   MatSetValues(*P,Nint,IIint,6,gl,xint,INSERT_VALUES);
532:   MatRestoreArray(Xint,&xint);
533:   MatGetArray(Xsurf,&xsurf);
534:   MatSetValues(*P,Nsurf,IIsurf,6,gl,xsurf,INSERT_VALUES);
535:   MatRestoreArray(Xsurf,&xsurf);
536:   MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);
537:   MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);
538:   PetscFree2(IIint,IIsurf);


541: #if defined(PETSC_USE_DEBUG_foo)
542:   {
543:     Vec         x,y;
544:     PetscScalar *yy;
545:     VecCreateMPI(((PetscObject)da)->comm,Ng,PETSC_DETERMINE,&y);
546:     VecCreateMPI(((PetscObject)da)->comm,PETSC_DETERMINE,Ntotal,&x);
547:     VecSet(x,1.0);
548:     MatMult(*P,x,y);
549:     VecGetArray(y,&yy);
550:     for (i=0; i<Ng; i++) {
551:       if (PetscAbsScalar(yy[i]-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong p interpolation at i %D value %G",i,PetscAbsScalar(yy[i]));
552:     }
553:     VecRestoreArray(y,&yy);
554:     VecDestroy(x);
555:     VecDestroy(y);
556:   }
557: #endif
558: 
559:   MatDestroy(Aii);
560:   MatDestroy(Ais);
561:   MatDestroy(Asi);
562:   MatDestroy(A);
563:   ISDestroy(is);
564:   ISDestroy(isint);
565:   ISDestroy(issurf);
566:   MatDestroy(Xint);
567:   MatDestroy(Xsurf);
568:   return(0);
569: }


574: /*@
575:    PCExoticSetType - Sets the type of coarse grid interpolation to use

577:    Collective on PC

579:    Input Parameters:
580: +  pc - the preconditioner context
581: -  type - either PC_EXOTIC_FACE or PC_EXOTIC_WIREBASKET (defaults to face)

583:    Notes: The face based interpolation has 1 degree of freedom per face and ignores the 
584:      edge and vertex values completely in the coarse problem. For any seven point
585:      stencil the interpolation of a constant on all faces into the interior is that constant.

587:      The wirebasket interpolation has 1 degree of freedom per vertex, per edge and 
588:      per face. A constant on the subdomain boundary is interpolated as that constant
589:      in the interior of the domain. 

591:      The coarse grid matrix is obtained via the Galerkin computation A_c = R A R^T, hence 
592:      if A is nonsingular A_c is also nonsingular.

594:      Both interpolations are suitable for only scalar problems.

596:    Level: intermediate


599: .seealso: PCEXOTIC, PCExoticType()
600: @*/
601: PetscErrorCode  PCExoticSetType(PC pc,PCExoticType type)
602: {
603:   PetscErrorCode ierr,(*f)(PC,PCExoticType);

607:   PetscObjectQueryFunction((PetscObject)pc,"PCExoticSetType_C",(void (**)(void))&f);
608:   if (f) {
609:     (*f)(pc,type);
610:   }
611:   return(0);
612: }

616: PetscErrorCode  PCExoticSetType_Exotic(PC pc,PCExoticType type)
617: {
618:   PC_MG     *mg = (PC_MG*)pc->data;
619:   PC_Exotic *ctx = (PC_Exotic*) mg->innerctx;

622:   ctx->type = type;
623:   return(0);
624: }

628: PetscErrorCode PCSetUp_Exotic(PC pc)
629: {
631:   Mat            A;
632:   PC_MG          *mg = (PC_MG*)pc->data;
633:   PC_Exotic      *ex = (PC_Exotic*) mg->innerctx;
634:   MatReuse       reuse = (ex->P) ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX;

637:   PCGetOperators(pc,PETSC_NULL,&A,PETSC_NULL);
638:   if (ex->type == PC_EXOTIC_FACE) {
639:     DAGetFaceInterpolation(ex,A,reuse,&ex->P);
640:   } else if (ex->type == PC_EXOTIC_WIREBASKET) {
641:     DAGetWireBasketInterpolation(ex,A,reuse,&ex->P);
642:   } else SETERRQ1(PETSC_ERR_PLIB,"Unknown exotic coarse space %d",ex->type);
643:   PCMGSetInterpolation(pc,1,ex->P);
644:   PCSetUp_MG(pc);
645:   return(0);
646: }

650: PetscErrorCode PCDestroy_Exotic(PC pc)
651: {
653:   PC_MG          *mg = (PC_MG*)pc->data;
654:   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;

657:   if (ctx->da) {DADestroy(ctx->da);}
658:   if (ctx->P) {MatDestroy(ctx->P);}
659:   if (ctx->ksp) {KSPDestroy(ctx->ksp);}
660:   PetscFree(ctx);
661:   PCDestroy_MG(pc);
662:   return(0);
663: }

667: PetscErrorCode PCSetUp_Exotic_Error(PC pc)
668: {
670:   SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"You are using the Exotic preconditioner but never called PCSetDA()");
671:   return(0);
672: }

676: /*@
677:    PCSetDA - Sets the DA that is to be used by the PCEXOTIC or certain other preconditioners

679:    Collective on PC

681:    Input Parameters:
682: +  pc - the preconditioner context
683: -  da - the da

685:    Level: intermediate


688: .seealso: PCEXOTIC, PCExoticType()
689: @*/
690: PetscErrorCode  PCSetDA(PC pc,DA da)
691: {
692:   PetscErrorCode ierr,(*f)(PC,DA);

697:   PetscObjectQueryFunction((PetscObject)pc,"PCSetDA_C",(void (**)(void))&f);
698:   if (f) {
699:     (*f)(pc,da);
700:   }
701:   return(0);
702: }


707: PetscErrorCode  PCSetDA_Exotic(PC pc,DA da)
708: {
710:   PC_MG          *mg = (PC_MG*)pc->data;
711:   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;

714:   ctx->da = da;
715:   pc->ops->setup = PCSetUp_Exotic;
716:   PetscObjectReference((PetscObject)da);
717:   return(0);
718: }

722: PetscErrorCode PCView_Exotic(PC pc,PetscViewer viewer)
723: {
724:   PC_MG          *mg = (PC_MG*)pc->data;
726:   PetscTruth     iascii;
727:   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;

730:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
731:   if (iascii) {
732:     PetscViewerASCIIPrintf(viewer,"    Exotic type = %s\n",PCExoticTypes[ctx->type]);
733:     if (ctx->directSolve) {
734:       PetscViewerASCIIPrintf(viewer,"      Using direct solver to construct interpolation\n");
735:     } else {
736:       PetscViewer sviewer;
737:       PetscMPIInt rank;

739:       PetscViewerASCIIPrintf(viewer,"      Using iterative solver to construct interpolation\n");
740:       PetscViewerASCIIPushTab(viewer);
741:       PetscViewerASCIIPushTab(viewer);  /* should not need to push this twice? */
742:       PetscViewerGetSingleton(viewer,&sviewer);
743:       MPI_Comm_rank(((PetscObject)pc)->comm,&rank);
744:       if (!rank) {
745:         KSPView(ctx->ksp,sviewer);
746:       }
747:       PetscViewerRestoreSingleton(viewer,&sviewer);
748:       PetscViewerASCIIPopTab(viewer);
749:       PetscViewerASCIIPopTab(viewer);
750:     }
751:   }
752:   PCView_MG(pc,viewer);
753:   return(0);
754: }

758: PetscErrorCode PCSetFromOptions_Exotic(PC pc)
759: {
761:   PetscTruth     flg;
762:   PC_MG          *mg = (PC_MG*)pc->data;
763:   PCExoticType   mgctype;
764:   PC_Exotic      *ctx = (PC_Exotic*) mg->innerctx;

767:   PetscOptionsHead("Exotic coarse space options");
768:     PetscOptionsEnum("-pc_exotic_type","face or wirebasket","PCExoticSetType",PCExoticTypes,(PetscEnum)ctx->type,(PetscEnum*)&mgctype,&flg);
769:     if (flg) {
770:       PCExoticSetType(pc,mgctype);
771:     }
772:     PetscOptionsTruth("-pc_exotic_direct_solver","use direct solver to construct interpolation","None",ctx->directSolve,&ctx->directSolve,PETSC_NULL);
773:     if (!ctx->directSolve) {
774:       if (!ctx->ksp) {
775:         const char *prefix;
776:         KSPCreate(PETSC_COMM_SELF,&ctx->ksp);
777:         PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)pc,1);
778:         PetscLogObjectParent(pc,ctx->ksp);
779:         PCGetOptionsPrefix(pc,&prefix);
780:         KSPSetOptionsPrefix(ctx->ksp,prefix);
781:         KSPAppendOptionsPrefix(ctx->ksp,"exotic_");
782:       }
783:       KSPSetFromOptions(ctx->ksp);
784:     }
785:   PetscOptionsTail();
786:   return(0);
787: }


790: /*MC
791:      PCEXOTIC - Two level overlapping Schwarz preconditioner with exotic (non-standard) coarse grid spaces

793:      This uses the PCMG infrastructure restricted to two levels and the face and wirebasket based coarse
794:    grid spaces. These coarse grid spaces originate in the work of Bramble, Pasciak  and Schatz, "The Construction
795:    of Preconditioners for Elliptic Problems by Substructing IV", Mathematics of Computation, volume 53 pages 1--24, 1989.
796:    They were generalized slightly in "Domain Decomposition Method for Linear Elasticity", Ph. D. thesis, Barry Smith,
797:    New York University, 1990. They were then explored in great detail in Dryja, Smith, Widlund, "Schwarz Analysis
798:    of Iterative Substructuring Methods for Elliptic Problems in Three Dimensions, SIAM Journal on Numerical 
799:    Analysis, volume 31. pages 1662-1694, 1994. These were developed in the context of iterative substructuring preconditioners.
800:    They were then ingeniously applied as coarse grid spaces for overlapping Schwarz methods by Dohrmann and Widlund.
801:    They refer to them as GDSW (generalized Dryja, Smith, Widlund preconditioners). See, for example, 
802:    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Extending theory for domain decomposition algorithms to irregular subdomains. In Ulrich Langer, Marco
803:    Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings
804:    of the 17th International Conference on Domain Decomposition Methods in
805:    Science and Engineering, held in Strobl, Austria, July 3-7, 2006, number 60 in
806:    Springer-Verlag, Lecture Notes in Computational Science and Engineering, pages 255-261, 2007.
807:    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. A family of energy min-
808:    imizing coarse spaces for overlapping Schwarz preconditioners. In Ulrich Langer,
809:    Marco Discacciati, David Keyes, OlofWidlund, andWalter Zulehner, editors, Proceedings
810:    of the 17th International Conference on Domain Decomposition Methods
811:    in Science and Engineering, held in Strobl, Austria, July 3-7, 2006, number 60 in
812:    Springer-Verlag, Lecture Notes in Computational Science and Engineering, pages 247-254, 2007
813:    Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Domain decomposition
814:    for less regular subdomains: Overlapping Schwarz in two dimensions. SIAM J.
815:    Numer. Anal., 46(4):2153-2168, 2008.
816:    Clark R. Dohrmann and Olof B. Widlund. An overlapping Schwarz
817:    algorithm for almost incompressible elasticity. Technical Report
818:    TR2008-912, Department of Computer Science, Courant Institute
819:    of Mathematical Sciences, New York University, May 2008. URL:
820:    http://cs.nyu.edu/csweb/Research/TechReports/TR2008-912/TR2008-912.pdf

822:    Options Database: The usual PCMG options are supported, such as -mg_levels_pc_type <type> -mg_coarse_pc_type <type>
823:       -pc_mg_type <type>

825:    Level: advanced

827: .seealso:  PCMG, PCSetDA(), PCExoticType, PCExoticSetType()
828: M*/

833: PetscErrorCode  PCCreate_Exotic(PC pc)
834: {
836:   PC_Exotic      *ex;
837:   PC_MG          *mg;

840:   /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
841:   if (pc->ops->destroy) {  (*pc->ops->destroy)(pc); pc->data = 0;}
842:   PetscStrfree(((PetscObject)pc)->type_name);
843:   ((PetscObject)pc)->type_name = 0;

845:   PCSetType(pc,PCMG);
846:   PCMGSetLevels(pc,2,PETSC_NULL);
847:   PCMGSetGalerkin(pc);
848:   PetscNew(PC_Exotic,&ex);\
849:   ex->type = PC_EXOTIC_FACE;
850:   mg = (PC_MG*) pc->data;
851:   mg->innerctx = ex;


854:   pc->ops->setfromoptions = PCSetFromOptions_Exotic;
855:   pc->ops->view           = PCView_Exotic;
856:   pc->ops->destroy        = PCDestroy_Exotic;
857:   pc->ops->setup          = PCSetUp_Exotic_Error;
858:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCExoticSetType_C","PCExoticSetType_Exotic",PCExoticSetType_Exotic);
859:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSetDA_C","PCSetDA_Exotic",PCSetDA_Exotic);
860:   return(0);
861: }