Actual source code: daint.c

  2: #define PETSCDM_DLL
  3: 
 4:  #include ../src/dm/da/daimpl.h
 5:  #include petscmat.h


 10: /*
 11:       DAGetWireBasketInterpolation - Gets the interpolation for a wirebasket based coarse space

 13: */
 14: PetscErrorCode DAGetWireBasketInterpolation(DA da,Mat Aglobal,MatReuse reuse,Mat *P)
 15: {
 16:   PetscErrorCode         ierr;
 17:   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;
 18:   PetscInt               mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[26],*globals,Ng,*IIint,*IIsurf;
 19:   Mat                    Xint, Xsurf,Xint_tmp;
 20:   IS                     isint,issurf,is,row,col;
 21:   ISLocalToGlobalMapping ltg;
 22:   MPI_Comm               comm;
 23:   Mat                    A,Aii,Ais,Asi,*Aholder,iAii;
 24:   MatFactorInfo          info;
 25:   PetscScalar            *xsurf,*xint;
 26: #if defined(PETSC_USE_DEBUG)
 27:   PetscScalar            tmp;
 28: #endif
 29:   PetscTable             ht;

 32:   DAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0);
 33:   if (dof != 1) SETERRQ(PETSC_ERR_SUP,"Only for single field problems");
 34:   if (dim != 3) SETERRQ(PETSC_ERR_SUP,"Only coded for 3d problems");
 35:   DAGetCorners(da,0,0,0,&m,&n,&p);
 36:   DAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);
 37:   istart = istart ? -1 : 0;
 38:   jstart = jstart ? -1 : 0;
 39:   kstart = kstart ? -1 : 0;

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

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

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

 49:     Xsurf are the interpolation onto the vertices and edges (the surfbasket) 
 50:                                         Xint
 51:     Symbolically one could write P = (  Xface  ) after interchanging the rows to match the natural ordering on the domain
 52:                                         Xsurf
 53:   */
 54:   N     = (m - istart)*(n - jstart)*(p - kstart);
 55:   Nint  = (m-2-istart)*(n-2-jstart)*(p-2-kstart);
 56:   Nface = 2*( (m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart) );
 57:   Nwire = 4*( (m-2-istart) + (n-2-jstart) + (p-2-kstart) ) + 8;
 58:   Nsurf = Nface + Nwire;
 59:   MatCreateSeqDense(MPI_COMM_SELF,Nint,26,PETSC_NULL,&Xint);
 60:   MatCreateSeqDense(MPI_COMM_SELF,Nsurf,26,PETSC_NULL,&Xsurf);
 61:   MatGetArray(Xsurf,&xsurf);

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

 71:   cnt = 0;
 72:   xsurf[cnt++] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + Nsurf] = 1;} xsurf[cnt++ + 2*Nsurf] = 1;
 73:   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;}
 74:   xsurf[cnt++ + 6*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 7*Nsurf] = 1;} xsurf[cnt++ + 8*Nsurf] = 1;
 75:   for (k=1;k<p-1-kstart;k++) {
 76:     xsurf[cnt++ + 9*Nsurf] = 1;  for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 10*Nsurf] = 1;}  xsurf[cnt++ + 11*Nsurf] = 1;
 77:     for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 12*Nsurf] = 1; /* these are the interior nodes */ xsurf[cnt++ + 13*Nsurf] = 1;}
 78:     xsurf[cnt++ + 14*Nsurf] = 1;  for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 15*Nsurf] = 1;} xsurf[cnt++ + 16*Nsurf] = 1;
 79:   }
 80:   xsurf[cnt++ + 17*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 18*Nsurf] = 1;} xsurf[cnt++ + 19*Nsurf] = 1;
 81:   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;}
 82:   xsurf[cnt++ + 23*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 24*Nsurf] = 1;} xsurf[cnt++ + 25*Nsurf] = 1;

 84: #if defined(PETSC_USE_DEBUG)
 85:   for (i=0; i<Nsurf; i++) {
 86:     tmp = 0.0;
 87:     for (j=0; j<26; j++) {
 88:       tmp += xsurf[i+j*Nsurf];
 89:     }
 90:     if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong Xsurf interpolation at i %D value %G",i,PetscAbsScalar(tmp));
 91:   }
 92: #endif
 93:   MatRestoreArray(Xsurf,&xsurf);
 94:   /* MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);*/


 97:   /* 
 98:        I are the indices for all the needed vertices (in global numbering)
 99:        Iint are the indices for the interior values, I surf for the surface values
100:             (This is just for the part of the global matrix obtained with MatGetSubMatrix(), it 
101:              is NOT the local DA ordering.)
102:        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
103:   */
104: #define Endpoint(a,start,b) (a == 0 || a == (b-1-start))
105:   PetscMalloc3(N,PetscInt,&II,Nint,PetscInt,&Iint,Nsurf,PetscInt,&Isurf);
106:   PetscMalloc2(Nint,PetscInt,&IIint,Nsurf,PetscInt,&IIsurf);
107:   for (k=0; k<p-kstart; k++) {
108:     for (j=0; j<n-jstart; j++) {
109:       for (i=0; i<m-istart; i++) {
110:         II[c++] = i + j*mwidth + k*mwidth*nwidth;

112:         if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) {
113:           IIint[cint]  = i + j*mwidth + k*mwidth*nwidth;
114:           Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
115:         } else {
116:           IIsurf[csurf]  = i + j*mwidth + k*mwidth*nwidth;
117:           Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
118:         }
119:       }
120:     }
121:   }
122:   if (c != N) SETERRQ(PETSC_ERR_PLIB,"c != N");
123:   if (cint != Nint) SETERRQ(PETSC_ERR_PLIB,"cint != Nint");
124:   if (csurf != Nsurf) SETERRQ(PETSC_ERR_PLIB,"csurf != Nsurf");
125:   DAGetISLocalToGlobalMapping(da,&ltg);
126:   ISLocalToGlobalMappingApply(ltg,N,II,II);
127:   ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);
128:   ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);
129:   PetscObjectGetComm((PetscObject)da,&comm);
130:   ISCreateGeneral(comm,N,II,&is);
131:   ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,&isint);
132:   ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,&issurf);
133:   PetscFree3(II,Iint,Isurf);

135:   MatGetSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);
136:   A    = *Aholder;
137:   PetscFree(Aholder);

139:   MatGetSubMatrix(A,isint,isint,PETSC_DECIDE,MAT_INITIAL_MATRIX,&Aii);
140:   MatGetSubMatrix(A,isint,issurf,PETSC_DECIDE,MAT_INITIAL_MATRIX,&Ais);
141:   MatGetSubMatrix(A,issurf,isint,PETSC_DECIDE,MAT_INITIAL_MATRIX,&Asi);

143:   /* 
144:      Solve for the interpolation onto the interior Xint
145:   */
146:   MatGetFactor(Aii,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&iAii);
147:   MatFactorInfoInitialize(&info);
148:   MatGetOrdering(Aii,MATORDERING_ND,&row,&col);
149:   MatLUFactorSymbolic(iAii,Aii,row,col,&info);
150:   ISDestroy(row);
151:   ISDestroy(col);
152:   MatLUFactorNumeric(iAii,Aii,&info);
153:   MatDuplicate(Xint,MAT_DO_NOT_COPY_VALUES,&Xint_tmp);
154:   MatMatMult(Ais,Xsurf,MAT_REUSE_MATRIX,PETSC_DETERMINE,&Xint_tmp);
155:   MatScale(Xint_tmp,-1.0);
156:   MatMatSolve(iAii,Xint_tmp,Xint);
157:   MatDestroy(Xint_tmp);
158:   MatDestroy(iAii);

160: #if defined(PETSC_USE_DEBUG)
161:   MatGetArray(Xint,&xint);
162:   for (i=0; i<Nint; i++) {
163:     tmp = 0.0;
164:     for (j=0; j<26; j++) {
165:       tmp += xint[i+j*Nint];
166:     }
167:     if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong Xint interpolation at i %D value %G",i,PetscAbsScalar(tmp));
168:   }
169:   MatRestoreArray(Xint,&xint);
170:   /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD); */
171: #endif


174:   /*         total vertices             total faces                                  total edges */
175:   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);

177:   /*
178:       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points 
179:   */
180:   cnt = 0;
181:   gl[cnt++] = 0;  { gl[cnt++] = 1;} gl[cnt++] = m-istart-1;
182:   { gl[cnt++] = mwidth;  { gl[cnt++] = mwidth+1;} gl[cnt++] = mwidth + m-istart-1;}
183:   gl[cnt++] = mwidth*(n-jstart-1);  { gl[cnt++] = mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*(n-jstart-1) + m-istart-1;
184:   {
185:     gl[cnt++] = mwidth*nwidth;  { gl[cnt++] = mwidth*nwidth+1;}  gl[cnt++] = mwidth*nwidth+ m-istart-1;
186:     { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;}
187:     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;
188:   }
189:   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;
190:   { 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;}
191:   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;

193:   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
194:   /* convert that to global numbering and get them on all processes */
195:   ISLocalToGlobalMappingApply(ltg,26,gl,gl);
196:   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */
197:   PetscMalloc(26*mp*np*pp*sizeof(PetscInt),&globals);
198:   MPI_Allgather(gl,26,MPIU_INT,globals,26,MPIU_INT,((PetscObject)da)->comm);

200:   /* Number the coarse grid points from 0 to Ntotal */
201:   PetscTableCreate(Ntotal/3,&ht);
202:   for (i=0; i<26*mp*np*pp; i++){
203:     PetscTableAddCount(ht,globals[i]+1);
204:   }
205:   PetscTableGetCount(ht,&cnt);
206:   if (cnt != Ntotal) SETERRQ2(PETSC_ERR_PLIB,"Hash table size %D not equal to total number coarse grid points %D",cnt,Ntotal);
207:   PetscFree(globals);
208:   for (i=0; i<26; i++) {
209:     PetscTableFind(ht,gl[i]+1,&gl[i]);
210:     gl[i]--;
211:   }
212:   PetscTableDestroy(ht);
213:   /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */

215:   /* construct global interpolation matrix */
216:   MatGetLocalSize(Aglobal,&Ng,PETSC_NULL);
217:   if (reuse == MAT_INITIAL_MATRIX) {
218:     MatCreateMPIAIJ(((PetscObject)da)->comm,Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,PETSC_NULL,Nint+Nsurf,PETSC_NULL,P);
219:   } else {
220:     MatZeroEntries(*P);
221:   }
222:   MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);
223:   MatGetArray(Xint,&xint);
224:   MatSetValues(*P,Nint,IIint,26,gl,xint,INSERT_VALUES);
225:   MatRestoreArray(Xint,&xint);
226:   MatGetArray(Xsurf,&xsurf);
227:   MatSetValues(*P,Nsurf,IIsurf,26,gl,xsurf,INSERT_VALUES);
228:   MatRestoreArray(Xsurf,&xsurf);
229:   MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);
230:   MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);
231:   PetscFree2(IIint,IIsurf);

233: #if defined(PETSC_USE_DEBUG)
234:   {
235:     Vec         x,y;
236:     PetscScalar *yy;
237:     VecCreateMPI(((PetscObject)da)->comm,Ng,PETSC_DETERMINE,&y);
238:     VecCreateMPI(((PetscObject)da)->comm,PETSC_DETERMINE,Ntotal,&x);
239:     VecSet(x,1.0);
240:     MatMult(*P,x,y);
241:     VecGetArray(y,&yy);
242:     for (i=0; i<Ng; i++) {
243:       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]));
244:     }
245:     VecRestoreArray(y,&yy);
246:     VecDestroy(x);
247:     VecDestroy(y);
248:   }
249: #endif
250: 
251:   MatDestroy(Aii);
252:   MatDestroy(Ais);
253:   MatDestroy(Asi);
254:   MatDestroy(A);
255:   ISDestroy(is);
256:   ISDestroy(isint);
257:   ISDestroy(issurf);
258:   MatDestroy(Xint);
259:   MatDestroy(Xsurf);
260:   return(0);
261: }

265: /*
266:       DAGetFaceInterpolation - Gets the interpolation for a face based coarse space

268: */
269: PetscErrorCode DAGetFaceInterpolation(DA da,Mat Aglobal,MatReuse reuse,Mat *P)
270: {
271:   PetscErrorCode         ierr;
272:   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;
273:   PetscInt               mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[6],*globals,Ng,*IIint,*IIsurf;
274:   Mat                    Xint, Xsurf,Xint_tmp;
275:   IS                     isint,issurf,is,row,col;
276:   ISLocalToGlobalMapping ltg;
277:   MPI_Comm               comm;
278:   Mat                    A,Aii,Ais,Asi,*Aholder,iAii;
279:   MatFactorInfo          info;
280:   PetscScalar            *xsurf,*xint;
281: #if defined(PETSC_USE_DEBUG_foo)
282:   PetscScalar            tmp;
283: #endif
284:   PetscTable             ht;

287:   DAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0);
288:   if (dof != 1) SETERRQ(PETSC_ERR_SUP,"Only for single field problems");
289:   if (dim != 3) SETERRQ(PETSC_ERR_SUP,"Only coded for 3d problems");
290:   DAGetCorners(da,0,0,0,&m,&n,&p);
291:   DAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);
292:   istart = istart ? -1 : 0;
293:   jstart = jstart ? -1 : 0;
294:   kstart = kstart ? -1 : 0;

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

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

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

304:     Xsurf are the interpolation onto the vertices and edges (the surfbasket) 
305:                                         Xint
306:     Symbolically one could write P = (  Xface  ) after interchanging the rows to match the natural ordering on the domain
307:                                         Xsurf
308:   */
309:   N     = (m - istart)*(n - jstart)*(p - kstart);
310:   Nint  = (m-2-istart)*(n-2-jstart)*(p-2-kstart);
311:   Nface = 2*( (m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart) );
312:   Nwire = 4*( (m-2-istart) + (n-2-jstart) + (p-2-kstart) ) + 8;
313:   Nsurf = Nface + Nwire;
314:   MatCreateSeqDense(MPI_COMM_SELF,Nint,6,PETSC_NULL,&Xint);
315:   MatCreateSeqDense(MPI_COMM_SELF,Nsurf,6,PETSC_NULL,&Xsurf);
316:   MatGetArray(Xsurf,&xsurf);

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

326:   cnt = 0;
327:   for (j=1;j<n-1-jstart;j++) {  for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 0*Nsurf] = 1;} }
328:    for (k=1;k<p-1-kstart;k++) {
329:     for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 1*Nsurf] = 1;}
330:     for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 2*Nsurf] = 1; /* these are the interior nodes */ xsurf[cnt++ + 3*Nsurf] = 1;}
331:     for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 4*Nsurf] = 1;}
332:   }
333:   for (j=1;j<n-1-jstart;j++) {for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 5*Nsurf] = 1;} }

335: #if defined(PETSC_USE_DEBUG_foo)
336:   for (i=0; i<Nsurf; i++) {
337:     tmp = 0.0;
338:     for (j=0; j<6; j++) {
339:       tmp += xsurf[i+j*Nsurf];
340:     }
341:     if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong Xsurf interpolation at i %D value %G",i,PetscAbsScalar(tmp));
342:   }
343: #endif
344:   MatRestoreArray(Xsurf,&xsurf);
345:   /* MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);*/


348:   /* 
349:        I are the indices for all the needed vertices (in global numbering)
350:        Iint are the indices for the interior values, I surf for the surface values
351:             (This is just for the part of the global matrix obtained with MatGetSubMatrix(), it 
352:              is NOT the local DA ordering.)
353:        IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering
354:   */
355: #define Endpoint(a,start,b) (a == 0 || a == (b-1-start))
356:   PetscMalloc3(N,PetscInt,&II,Nint,PetscInt,&Iint,Nsurf,PetscInt,&Isurf);
357:   PetscMalloc2(Nint,PetscInt,&IIint,Nsurf,PetscInt,&IIsurf);
358:   for (k=0; k<p-kstart; k++) {
359:     for (j=0; j<n-jstart; j++) {
360:       for (i=0; i<m-istart; i++) {
361:         II[c++] = i + j*mwidth + k*mwidth*nwidth;

363:         if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) {
364:           IIint[cint]  = i + j*mwidth + k*mwidth*nwidth;
365:           Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
366:         } else {
367:           IIsurf[csurf]  = i + j*mwidth + k*mwidth*nwidth;
368:           Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart);
369:         }
370:       }
371:     }
372:   }
373:   if (c != N) SETERRQ(PETSC_ERR_PLIB,"c != N");
374:   if (cint != Nint) SETERRQ(PETSC_ERR_PLIB,"cint != Nint");
375:   if (csurf != Nsurf) SETERRQ(PETSC_ERR_PLIB,"csurf != Nsurf");
376:   DAGetISLocalToGlobalMapping(da,&ltg);
377:   ISLocalToGlobalMappingApply(ltg,N,II,II);
378:   ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);
379:   ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);
380:   PetscObjectGetComm((PetscObject)da,&comm);
381:   ISCreateGeneral(comm,N,II,&is);
382:   ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,&isint);
383:   ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,&issurf);
384:   PetscFree3(II,Iint,Isurf);

386:   MatGetSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);
387:   A    = *Aholder;
388:   PetscFree(Aholder);

390:   MatGetSubMatrix(A,isint,isint,PETSC_DECIDE,MAT_INITIAL_MATRIX,&Aii);
391:   MatGetSubMatrix(A,isint,issurf,PETSC_DECIDE,MAT_INITIAL_MATRIX,&Ais);
392:   MatGetSubMatrix(A,issurf,isint,PETSC_DECIDE,MAT_INITIAL_MATRIX,&Asi);

394:   /* 
395:      Solve for the interpolation onto the interior Xint
396:   */
397:   MatGetFactor(Aii,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&iAii);
398:   MatFactorInfoInitialize(&info);
399:   MatGetOrdering(Aii,MATORDERING_ND,&row,&col);
400:   MatLUFactorSymbolic(iAii,Aii,row,col,&info);
401:   ISDestroy(row);
402:   ISDestroy(col);
403:   MatLUFactorNumeric(iAii,Aii,&info);
404:   MatDuplicate(Xint,MAT_DO_NOT_COPY_VALUES,&Xint_tmp);
405:   MatMatMult(Ais,Xsurf,MAT_REUSE_MATRIX,PETSC_DETERMINE,&Xint_tmp);
406:   MatScale(Xint_tmp,-1.0);
407:   MatMatSolve(iAii,Xint_tmp,Xint);
408:   MatDestroy(Xint_tmp);
409:   MatDestroy(iAii);

411: #if defined(PETSC_USE_DEBUG_foo)
412:   MatGetArray(Xint,&xint);
413:   for (i=0; i<Nint; i++) {
414:     tmp = 0.0;
415:     for (j=0; j<6; j++) {
416:       tmp += xint[i+j*Nint];
417:     }
418:     if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong Xint interpolation at i %D value %G",i,PetscAbsScalar(tmp));
419:   }
420:   MatRestoreArray(Xint,&xint);
421:   /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD); */
422: #endif


425:   /*         total faces    */
426:   Ntotal =  mp*np*(pp+1) + mp*pp*(np+1) + np*pp*(mp+1);

428:   /*
429:       For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points 
430:   */
431:   cnt = 0;
432:   { gl[cnt++] = mwidth+1;}
433:   {
434:     { gl[cnt++] = mwidth*nwidth+1;}
435:     { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;}
436:     { gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1)+1;}
437:   }
438:   { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth+1;}

440:   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
441:   /* convert that to global numbering and get them on all processes */
442:   ISLocalToGlobalMappingApply(ltg,6,gl,gl);
443:   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */
444:   PetscMalloc(6*mp*np*pp*sizeof(PetscInt),&globals);
445:   MPI_Allgather(gl,6,MPIU_INT,globals,6,MPIU_INT,((PetscObject)da)->comm);

447:   /* Number the coarse grid points from 0 to Ntotal */
448:   PetscTableCreate(Ntotal/3,&ht);
449:   for (i=0; i<6*mp*np*pp; i++){
450:     PetscTableAddCount(ht,globals[i]+1);
451:   }
452:   PetscTableGetCount(ht,&cnt);
453:   if (cnt != Ntotal) SETERRQ2(PETSC_ERR_PLIB,"Hash table size %D not equal to total number coarse grid points %D",cnt,Ntotal);
454:   PetscFree(globals);
455:   for (i=0; i<6; i++) {
456:     PetscTableFind(ht,gl[i]+1,&gl[i]);
457:     gl[i]--;
458:   }
459:   PetscTableDestroy(ht);
460:   /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */

462:   /* construct global interpolation matrix */
463:   MatGetLocalSize(Aglobal,&Ng,PETSC_NULL);
464:   if (reuse == MAT_INITIAL_MATRIX) {
465:     MatCreateMPIAIJ(((PetscObject)da)->comm,Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,PETSC_NULL,Nint,PETSC_NULL,P);
466:   } else {
467:     MatZeroEntries(*P);
468:   }
469:   MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);
470:   MatGetArray(Xint,&xint);
471:   MatSetValues(*P,Nint,IIint,6,gl,xint,INSERT_VALUES);
472:   MatRestoreArray(Xint,&xint);
473:   MatGetArray(Xsurf,&xsurf);
474:   MatSetValues(*P,Nsurf,IIsurf,6,gl,xsurf,INSERT_VALUES);
475:   MatRestoreArray(Xsurf,&xsurf);
476:   MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);
477:   MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);
478:   PetscFree2(IIint,IIsurf);


481: #if defined(PETSC_USE_DEBUG_foo)
482:   {
483:     Vec         x,y;
484:     PetscScalar *yy;
485:     VecCreateMPI(((PetscObject)da)->comm,Ng,PETSC_DETERMINE,&y);
486:     VecCreateMPI(((PetscObject)da)->comm,PETSC_DETERMINE,Ntotal,&x);
487:     VecSet(x,1.0);
488:     MatMult(*P,x,y);
489:     VecGetArray(y,&yy);
490:     for (i=0; i<Ng; i++) {
491:       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]));
492:     }
493:     VecRestoreArray(y,&yy);
494:     VecDestroy(x);
495:     VecDestroy(y);
496:   }
497: #endif
498: 
499:   MatDestroy(Aii);
500:   MatDestroy(Ais);
501:   MatDestroy(Asi);
502:   MatDestroy(A);
503:   ISDestroy(is);
504:   ISDestroy(isint);
505:   ISDestroy(issurf);
506:   MatDestroy(Xint);
507:   MatDestroy(Xsurf);
508:   return(0);
509: }