Actual source code: dainterp.c

  1: /*$Id: dainterp.c,v 1.25 2001/08/07 03:04:39 balay Exp $*/
  2: 
  3: /*
  4:   Code for interpolating between grids represented by DAs
  5: */

 7:  #include src/dm/da/daimpl.h
 8:  #include petscmg.h

 10: int DMGetInterpolationScale(DM dac,DM daf,Mat mat,Vec *scale)
 11: {
 12:   int    ierr;
 13:   Vec    fine;
 14:   PetscScalar one = 1.0;

 17:   DMCreateGlobalVector(daf,&fine);
 18:   DMCreateGlobalVector(dac,scale);
 19:   VecSet(&one,fine);
 20:   MatRestrict(mat,fine,*scale);
 21:   VecDestroy(fine);
 22:   VecReciprocal(*scale);
 23:   return(0);
 24: }

 26: int DAGetInterpolation_1D_Q1(DA dac,DA daf,Mat *A)
 27: {
 28:   int            ierr,i,i_start,m_f,Mx,*idx_f;
 29:   int            m_ghost,*idx_c,m_ghost_c;
 30:   int            row,col,i_start_ghost,mx,m_c,nc,ratio;
 31:   int            i_c,i_start_c,i_start_ghost_c,cols[2],dof;
 32:   PetscScalar    v[2],x;
 33:   Mat            mat;
 34:   DAPeriodicType pt;

 37:   DAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&pt,0);
 38:   DAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0);
 39:   if (pt == DA_XPERIODIC) {
 40:     ratio = mx/Mx;
 41:     if (ratio*Mx != mx) SETERRQ2(1,"Ratio between levels: mx/Mx  must be integer: mx %d Mx %d",mx,Mx);
 42:   } else {
 43:     ratio = (mx-1)/(Mx-1);
 44:     if (ratio*(Mx-1) != mx-1) SETERRQ2(1,"Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %d Mx %d",mx,Mx);
 45:   }

 47:   DAGetCorners(daf,&i_start,0,0,&m_f,0,0);
 48:   DAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);
 49:   DAGetGlobalIndices(daf,PETSC_NULL,&idx_f);

 51:   DAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);
 52:   DAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);
 53:   DAGetGlobalIndices(dac,PETSC_NULL,&idx_c);

 55:   /* create interpolation matrix */
 56:   MatCreateMPIAIJ(dac->comm,m_f,m_c,mx,Mx,2,0,0,0,&mat);
 57:   if (!DAXPeriodic(pt)){MatSetOption(mat,MAT_COLUMNS_SORTED);}

 59:   /* loop over local fine grid nodes setting interpolation for those*/
 60:   for (i=i_start; i<i_start+m_f; i++) {
 61:     /* convert to local "natural" numbering and then to PETSc global numbering */
 62:     row    = idx_f[dof*(i-i_start_ghost)]/dof;

 64:     i_c = (i/ratio);    /* coarse grid node to left of fine grid node */
 65:     if (i_c < i_start_ghost_c) SETERRQ3(1,"Processor's coarse DA must lie over fine DAn
 66:     i_start %d i_c %d i_start_ghost_c %d",i_start,i_c,i_start_ghost_c);

 68:     /* 
 69:          Only include those interpolation points that are truly 
 70:          nonzero. Note this is very important for final grid lines
 71:          in x direction; since they have no right neighbor
 72:     */
 73:     x  = ((double)(i - i_c*ratio))/((double)ratio);
 74:     /* printf("i j %d %d %g %gn",i,j,x,y); */
 75:     nc = 0;
 76:       /* one left and below; or we are right on it */
 77:     col      = dof*(i_c-i_start_ghost_c);
 78:     cols[nc] = idx_c[col]/dof;
 79:     v[nc++]  = - x + 1.0;
 80:     /* one right? */
 81:     if (i_c*ratio != i) {
 82:       cols[nc] = idx_c[col+dof]/dof;
 83:       v[nc++]  = x;
 84:     }
 85:     MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);
 86:   }
 87:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
 88:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
 89:   MatCreateMAIJ(mat,dof,A);
 90:   MatDestroy(mat);
 91:   PetscLogFlops(5*m_f);
 92:   return(0);
 93: }

 95: int DAGetInterpolation_1D_Q0(DA dac,DA daf,Mat *A)
 96: {
 97:   int            ierr,i,i_start,m_f,Mx,*idx_f;
 98:   int            m_ghost,*idx_c,m_ghost_c;
 99:   int            row,col,i_start_ghost,mx,m_c,nc,ratio;
100:   int            i_c,i_start_c,i_start_ghost_c,cols[2],dof;
101:   PetscScalar    v[2],x;
102:   Mat            mat;
103:   DAPeriodicType pt;

106:   DAGetInfo(dac,0,&Mx,0,0,0,0,0,0,0,&pt,0);
107:   DAGetInfo(daf,0,&mx,0,0,0,0,0,&dof,0,0,0);
108:   if (pt == DA_XPERIODIC) {
109:     ratio = mx/Mx;
110:     if (ratio*Mx != mx) SETERRQ2(1,"Ratio between levels: mx/Mx  must be integer: mx %d Mx %d",mx,Mx);
111:   } else {
112:     ratio = (mx-1)/(Mx-1);
113:     if (ratio*(Mx-1) != mx-1) SETERRQ2(1,"Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %d Mx %d",mx,Mx);
114:   }

116:   DAGetCorners(daf,&i_start,0,0,&m_f,0,0);
117:   DAGetGhostCorners(daf,&i_start_ghost,0,0,&m_ghost,0,0);
118:   DAGetGlobalIndices(daf,PETSC_NULL,&idx_f);

120:   DAGetCorners(dac,&i_start_c,0,0,&m_c,0,0);
121:   DAGetGhostCorners(dac,&i_start_ghost_c,0,0,&m_ghost_c,0,0);
122:   DAGetGlobalIndices(dac,PETSC_NULL,&idx_c);

124:   /* create interpolation matrix */
125:   MatCreateMPIAIJ(dac->comm,m_f,m_c,mx,Mx,2,0,0,0,&mat);
126:   if (!DAXPeriodic(pt)) {MatSetOption(mat,MAT_COLUMNS_SORTED);}

128:   /* loop over local fine grid nodes setting interpolation for those*/
129:   for (i=i_start; i<i_start+m_f; i++) {
130:     /* convert to local "natural" numbering and then to PETSc global numbering */
131:     row    = idx_f[dof*(i-i_start_ghost)]/dof;

133:     i_c = (i/ratio);    /* coarse grid node to left of fine grid node */

135:     /* 
136:          Only include those interpolation points that are truly 
137:          nonzero. Note this is very important for final grid lines
138:          in x direction; since they have no right neighbor
139:     */
140:     x  = ((double)(i - i_c*ratio))/((double)ratio);
141:     /* printf("i j %d %d %g %gn",i,j,x,y); */
142:     nc = 0;
143:       /* one left and below; or we are right on it */
144:     col      = dof*(i_c-i_start_ghost_c);
145:     cols[nc] = idx_c[col]/dof;
146:     v[nc++]  = - x + 1.0;
147:     /* one right? */
148:     if (i_c*ratio != i) {
149:       cols[nc] = idx_c[col+dof]/dof;
150:       v[nc++]  = x;
151:     }
152:     MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);
153:   }
154:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
155:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
156:   MatCreateMAIJ(mat,dof,A);
157:   MatDestroy(mat);
158:   PetscLogFlops(5*m_f);
159:   return(0);
160: }


163: /*   dof degree of freedom per node, nonperiodic */
164: int DAGetInterpolation_2D_Q1(DA dac,DA daf,Mat *A)
165: {
166:   int            ierr,i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof;
167:   int            m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,*dnz,*onz;
168:   int            row,col,i_start_ghost,j_start_ghost,cols[4],mx,m_c,my,nc,ratioi,ratioj;
169:   int            i_c,j_c,i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c;
170:   int            size_c,size_f,rank_f,col_shift,col_scale;
171:   PetscScalar    v[4],x,y;
172:   Mat            mat;
173:   DAPeriodicType pt;


177:   DAGetInfo(dac,0,&Mx,&My,0,0,0,0,0,0,&pt,0);
178:   DAGetInfo(daf,0,&mx,&my,0,0,0,0,&dof,0,0,0);
179:   if (DAXPeriodic(pt)){
180:     ratioi = mx/Mx;
181:     if (ratioi*Mx != mx) SETERRQ2(1,"Ratio between levels: mx/Mx  must be integer: mx %d Mx %d",mx,Mx);
182:   } else {
183:     ratioi = (mx-1)/(Mx-1);
184:     if (ratioi*(Mx-1) != mx-1) SETERRQ2(1,"Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %d Mx %d",mx,Mx);
185:   }
186:   if (DAYPeriodic(pt)){
187:     ratioj = my/My;
188:     if (ratioj*My != my) SETERRQ2(1,"Ratio between levels: my/My  must be integer: my %d My %d",my,My);
189:   } else {
190:     ratioj = (my-1)/(My-1);
191:     if (ratioj*(My-1) != my-1) SETERRQ2(1,"Ratio between levels: (my - 1)/(My - 1) must be integer: my %d My %d",my,My);
192:   }


195:   DAGetCorners(daf,&i_start,&j_start,0,&m_f,&n_f,0);
196:   DAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,0,&m_ghost,&n_ghost,0);
197:   DAGetGlobalIndices(daf,PETSC_NULL,&idx_f);

199:   DAGetCorners(dac,&i_start_c,&j_start_c,0,&m_c,&n_c,0);
200:   DAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,0,&m_ghost_c,&n_ghost_c,0);
201:   DAGetGlobalIndices(dac,PETSC_NULL,&idx_c);

203:   /*
204:      Used for handling a coarse DA that lives on 1/4 the processors of the fine DA.
205:      The coarse vector is then duplicated 4 times (each time it lives on 1/4 of the 
206:      processors). It's effective length is hence 4 times its normal length, this is
207:      why the col_scale is multiplied by the interpolation matrix column sizes.
208:      sol_shift allows each set of 1/4 processors do its own interpolation using ITS
209:      copy of the coarse vector. A bit of a hack but you do better.

211:      In the standard case when size_f == size_c col_scale == 1 and col_shift == 0
212:   */
213:   MPI_Comm_size(dac->comm,&size_c);
214:   MPI_Comm_size(daf->comm,&size_f);
215:   MPI_Comm_rank(daf->comm,&rank_f);
216:   col_scale = size_f/size_c;
217:   col_shift = Mx*My*(rank_f/size_c);

219:   MatPreallocateInitialize(daf->comm,m_f*n_f,col_scale*m_c*n_c,dnz,onz);
220:   for (j=j_start; j<j_start+n_f; j++) {
221:     for (i=i_start; i<i_start+m_f; i++) {
222:       /* convert to local "natural" numbering and then to PETSc global numbering */
223:       row    = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;

225:       i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
226:       j_c = (j/ratioj);    /* coarse grid node below fine grid node */

228:       if (j_c < j_start_ghost_c) SETERRQ3(1,"Processor's coarse DA must lie over fine DAn
229:     j_start %d j_c %d j_start_ghost_c %d",j_start,j_c,j_start_ghost_c);
230:       if (i_c < i_start_ghost_c) SETERRQ3(1,"Processor's coarse DA must lie over fine DAn
231:     i_start %d i_c %d i_start_ghost_c %d",i_start,i_c,i_start_ghost_c);

233:       /* 
234:          Only include those interpolation points that are truly 
235:          nonzero. Note this is very important for final grid lines
236:          in x and y directions; since they have no right/top neighbors
237:       */
238:       nc = 0;
239:       /* one left and below; or we are right on it */
240:       col        = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
241:       cols[nc++] = col_shift + idx_c[col]/dof;
242:       /* one right and below */
243:       if (i_c*ratioi != i) {
244:         cols[nc++] = col_shift + idx_c[col+dof]/dof;
245:       }
246:       /* one left and above */
247:       if (j_c*ratioj != j) {
248:         cols[nc++] = col_shift + idx_c[col+m_ghost_c*dof]/dof;
249:       }
250:       /* one right and above */
251:       if (j_c*ratioi != j && i_c*ratioj != i) {
252:         cols[nc++] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof;
253:       }
254:       MatPreallocateSet(row,nc,cols,dnz,onz);
255:     }
256:   }
257:   MatCreateMPIAIJ(daf->comm,m_f*n_f,col_scale*m_c*n_c,mx*my,col_scale*Mx*My,0,dnz,0,onz,&mat);
258:   MatPreallocateFinalize(dnz,onz);
259:   if (!DAXPeriodic(pt) && !DAYPeriodic(pt)) {MatSetOption(mat,MAT_COLUMNS_SORTED);}

261:   /* loop over local fine grid nodes setting interpolation for those*/
262:   for (j=j_start; j<j_start+n_f; j++) {
263:     for (i=i_start; i<i_start+m_f; i++) {
264:       /* convert to local "natural" numbering and then to PETSc global numbering */
265:       row    = idx_f[dof*(m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;

267:       i_c = (i/ratioi);    /* coarse grid node to left of fine grid node */
268:       j_c = (j/ratioj);    /* coarse grid node below fine grid node */

270:       /* 
271:          Only include those interpolation points that are truly 
272:          nonzero. Note this is very important for final grid lines
273:          in x and y directions; since they have no right/top neighbors
274:       */
275:       x  = ((double)(i - i_c*ratioi))/((double)ratioi);
276:       y  = ((double)(j - j_c*ratioj))/((double)ratioj);
277:       /* printf("i j %d %d %g %gn",i,j,x,y); */
278:       nc = 0;
279:       /* one left and below; or we are right on it */
280:       col      = dof*(m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
281:       cols[nc] = col_shift + idx_c[col]/dof;
282:       v[nc++]  = x*y - x - y + 1.0;
283:       /* one right and below */
284:       if (i_c*ratioi != i) {
285:         cols[nc] = col_shift + idx_c[col+dof]/dof;
286:         v[nc++]  = -x*y + x;
287:       }
288:       /* one left and above */
289:       if (j_c*ratioj != j) {
290:         cols[nc] = col_shift + idx_c[col+m_ghost_c*dof]/dof;
291:         v[nc++]  = -x*y + y;
292:       }
293:       /* one right and above */
294:       if (j_c*ratioj != j && i_c*ratioi != i) {
295:         cols[nc] = col_shift + idx_c[col+(m_ghost_c+1)*dof]/dof;
296:         v[nc++]  = x*y;
297:       }
298:       MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);
299:     }
300:   }
301:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
302:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
303:   MatCreateMAIJ(mat,dof,A);
304:   MatDestroy(mat);
305:   PetscLogFlops(13*m_f*n_f);
306:   return(0);
307: }


310: /*   dof degree of freedom per node, nonperiodic */
311: int DAGetInterpolation_3D_Q1(DA dac,DA daf,Mat *A)
312: {
313:   int            ierr,i,j,i_start,j_start,m_f,n_f,Mx,My,*idx_f,dof,l;
314:   int            m_ghost,n_ghost,*idx_c,m_ghost_c,n_ghost_c,Mz,mz;
315:   int            row,col,i_start_ghost,j_start_ghost,cols[8],mx,m_c,my,nc,ratioi,ratioj,ratiok;
316:   int            i_c,j_c,i_start_c,j_start_c,n_c,i_start_ghost_c,j_start_ghost_c;
317:   int            l_start,p_f,l_start_ghost,p_ghost,l_start_c,p_c;
318:   int            l_start_ghost_c,p_ghost_c,l_c,*dnz,*onz;
319:   PetscScalar    v[8],x,y,z;
320:   Mat            mat;
321:   DAPeriodicType pt;

324:   DAGetInfo(dac,0,&Mx,&My,&Mz,0,0,0,0,0,&pt,0);
325:   DAGetInfo(daf,0,&mx,&my,&mz,0,0,0,&dof,0,0,0);
326:   if (DAXPeriodic(pt)){
327:     ratioi = mx/Mx;
328:     if (ratioi*Mx != mx) SETERRQ2(1,"Ratio between levels: mx/Mx  must be integer: mx %d Mx %d",mx,Mx);
329:   } else {
330:     ratioi = (mx-1)/(Mx-1);
331:     if (ratioi*(Mx-1) != mx-1) SETERRQ2(1,"Ratio between levels: (mx - 1)/(Mx - 1) must be integer: mx %d Mx %d",mx,Mx);
332:   }
333:   if (DAYPeriodic(pt)){
334:     ratioj = my/My;
335:     if (ratioj*My != my) SETERRQ2(1,"Ratio between levels: my/My  must be integer: my %d My %d",my,My);
336:   } else {
337:     ratioj = (my-1)/(My-1);
338:     if (ratioj*(My-1) != my-1) SETERRQ2(1,"Ratio between levels: (my - 1)/(My - 1) must be integer: my %d My %d",my,My);
339:   }
340:   if (DAZPeriodic(pt)){
341:     ratiok = mz/Mz;
342:     if (ratiok*Mz != mz) SETERRQ2(1,"Ratio between levels: mz/Mz  must be integer: mz %d Mz %d",mz,Mz);
343:   } else {
344:     ratiok = (mz-1)/(Mz-1);
345:     if (ratiok*(Mz-1) != mz-1) SETERRQ2(1,"Ratio between levels: (mz - 1)/(Mz - 1) must be integer: mz %d Mz %d",mz,Mz);
346:   }

348:   DAGetCorners(daf,&i_start,&j_start,&l_start,&m_f,&n_f,&p_f);
349:   DAGetGhostCorners(daf,&i_start_ghost,&j_start_ghost,&l_start_ghost,&m_ghost,&n_ghost,&p_ghost);
350:   DAGetGlobalIndices(daf,PETSC_NULL,&idx_f);

352:   DAGetCorners(dac,&i_start_c,&j_start_c,&l_start_c,&m_c,&n_c,&p_c);
353:   DAGetGhostCorners(dac,&i_start_ghost_c,&j_start_ghost_c,&l_start_ghost_c,&m_ghost_c,&n_ghost_c,&p_ghost_c);
354:   DAGetGlobalIndices(dac,PETSC_NULL,&idx_c);

356:   /* create interpolation matrix, determining exact preallocation */
357:   MatPreallocateInitialize(dac->comm,m_f*n_f*p_f,m_c*n_c*p_c,dnz,onz);
358:   /* loop over local fine grid nodes counting interpolating points */
359:   for (l=l_start; l<l_start+p_f; l++) {
360:     for (j=j_start; j<j_start+n_f; j++) {
361:       for (i=i_start; i<i_start+m_f; i++) {
362:         /* convert to local "natural" numbering and then to PETSc global numbering */
363:         row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;
364:         i_c = (i/ratioi);
365:         j_c = (j/ratioj);
366:         l_c = (l/ratiok);
367:         if (l_c < l_start_ghost_c) SETERRQ3(1,"Processor's coarse DA must lie over fine DAn
368:           l_start %d l_c %d l_start_ghost_c %d",l_start,l_c,l_start_ghost_c);
369:         if (j_c < j_start_ghost_c) SETERRQ3(1,"Processor's coarse DA must lie over fine DAn
370:           j_start %d j_c %d j_start_ghost_c %d",j_start,j_c,j_start_ghost_c);
371:         if (i_c < i_start_ghost_c) SETERRQ3(1,"Processor's coarse DA must lie over fine DAn
372:           i_start %d i_c %d i_start_ghost_c %d",i_start,i_c,i_start_ghost_c);

374:         /* 
375:          Only include those interpolation points that are truly 
376:          nonzero. Note this is very important for final grid lines
377:          in x and y directions; since they have no right/top neighbors
378:         */
379:         nc       = 0;
380:         col      = dof*(m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c) + m_ghost_c*(j_c-j_start_ghost_c) + (i_c-i_start_ghost_c));
381:         cols[nc++] = idx_c[col]/dof;
382:         if (i_c*ratioi != i) {
383:           cols[nc++] = idx_c[col+dof]/dof;
384:         }
385:         if (j_c*ratioj != j) {
386:           cols[nc++] = idx_c[col+m_ghost_c*dof]/dof;
387:         }
388:         if (l_c*ratiok != l) {
389:           cols[nc++] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof;
390:         }
391:         if (j_c*ratioj != j && i_c*ratioi != i) {
392:           cols[nc++] = idx_c[col+(m_ghost_c+1)*dof]/dof;
393:         }
394:         if (j_c*ratioj != j && l_c*ratiok != l) {
395:           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof;
396:         }
397:         if (i_c*ratioi != i && l_c*ratiok != l) {
398:           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof;
399:         }
400:         if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) {
401:           cols[nc++] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof;
402:         }
403:         MatPreallocateSet(row,nc,cols,dnz,onz);
404:       }
405:     }
406:   }
407:   MatCreateMPIAIJ(dac->comm,m_f*n_f*p_f,m_c*n_c*p_c,mx*my*mz,Mx*My*Mz,0,dnz,0,onz,&mat);
408:   MatPreallocateFinalize(dnz,onz);
409:   if (!DAXPeriodic(pt) && !DAYPeriodic(pt) && !DAZPeriodic(pt)) {MatSetOption(mat,MAT_COLUMNS_SORTED);}

411:   /* loop over local fine grid nodes setting interpolation for those*/
412:   for (l=l_start; l<l_start+p_f; l++) {
413:     for (j=j_start; j<j_start+n_f; j++) {
414:       for (i=i_start; i<i_start+m_f; i++) {
415:         /* convert to local "natural" numbering and then to PETSc global numbering */
416:         row = idx_f[dof*(m_ghost*n_ghost*(l-l_start_ghost) + m_ghost*(j-j_start_ghost) + (i-i_start_ghost))]/dof;

418:         i_c = (i/ratioi);
419:         j_c = (j/ratioj);
420:         l_c = (l/ratiok);

422:         /* 
423:            Only include those interpolation points that are truly 
424:            nonzero. Note this is very important for final grid lines
425:            in x and y directions; since they have no right/top neighbors
426:         */
427:         x  = ((double)(i - i_c*ratioi))/((double)ratioi);
428:         y  = ((double)(j - j_c*ratioj))/((double)ratioj);
429:         z  = ((double)(l - l_c*ratiok))/((double)ratiok);
430:         /* printf("i j l %d %d %d %g %g %gn",i,j,l,x,y,z); */
431:         nc = 0;
432:         /* one left and below; or we are right on it */
433:         col      = dof*(m_ghost_c*n_ghost_c*(l_c-l_start_ghost_c)+m_ghost_c*(j_c-j_start_ghost_c)+(i_c-i_start_ghost_c));

435:         cols[nc] = idx_c[col]/dof;
436:         v[nc++]  = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.));

438:         if (i_c*ratioi != i) {
439:           cols[nc] = idx_c[col+dof]/dof;
440:           v[nc++]  = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. - (2.0*z-1.));
441:         }

443:         if (j_c*ratioj != j) {
444:           cols[nc] = idx_c[col+m_ghost_c*dof]/dof;
445:           v[nc++]  = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.));
446:         }

448:         if (l_c*ratiok != l) {
449:           cols[nc] = idx_c[col+m_ghost_c*n_ghost_c*dof]/dof;
450:           v[nc++]  = .125*(1. - (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.));
451:         }

453:         if (j_c*ratioj != j && i_c*ratioi != i) {
454:           cols[nc] = idx_c[col+(m_ghost_c+1)*dof]/dof;
455:           v[nc++]  = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. - (2.0*z-1.));
456:         }

458:         if (j_c*ratioj != j && l_c*ratiok != l) {
459:           cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c)*dof]/dof;
460:           v[nc++]  = .125*(1. - (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.));
461:         }

463:         if (i_c*ratioi != i && l_c*ratiok != l) {
464:           cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+1)*dof]/dof;
465:           v[nc++]  = .125*(1. + (2.0*x-1.))*(1. - (2.0*y-1.))*(1. + (2.0*z-1.));
466:         }

468:         if (i_c*ratioi != i && l_c*ratiok != l && j_c*ratioj != j) {
469:           cols[nc] = idx_c[col+(m_ghost_c*n_ghost_c+m_ghost_c+1)*dof]/dof;
470:           v[nc++]  = .125*(1. + (2.0*x-1.))*(1. + (2.0*y-1.))*(1. + (2.0*z-1.));
471:         }
472:         MatSetValues(mat,1,&row,nc,cols,v,INSERT_VALUES);
473:       }
474:     }
475:   }
476:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
477:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);

479:   MatCreateMAIJ(mat,dof,A);
480:   MatDestroy(mat);
481:   PetscLogFlops(13*m_f*n_f);
482:   return(0);
483: }

485: /*@C
486:    DAGetInterpolation - Gets an interpolation matrix that maps between 
487:    grids associated with two DAs.

489:    Collective on DA

491:    Input Parameters:
492: +  dac - the coarse grid DA
493: -  daf - the fine grid DA

495:    Output Parameters:
496: +  A - the interpolation matrix
497: -  scale - a scaling vector used to scale the coarse grid restricted vector before applying the 
498:            grid function or grid Jacobian to it.

500:    Level: intermediate

502: .keywords: interpolation, restriction, multigrid 

504: .seealso: DARefine()
505: @*/
506: int DAGetInterpolation(DA dac,DA daf,Mat *A,Vec *scale)
507: {
508:   int            ierr,dimc,Mc,Nc,Pc,mc,nc,pc,dofc,sc,dimf,Mf,Nf,Pf,mf,nf,pf,doff,sf;
509:   DAPeriodicType wrapc,wrapf;
510:   DAStencilType  stc,stf;


516:   DAGetInfo(dac,&dimc,&Mc,&Nc,&Pc,&mc,&nc,&pc,&dofc,&sc,&wrapc,&stc);
517:   DAGetInfo(daf,&dimf,&Mf,&Nf,&Pf,&mf,&nf,&pf,&doff,&sf,&wrapf,&stf);
518:   if (dimc != dimf) SETERRQ2(1,"Dimensions of DA do not match %d %d",dimc,dimf);
519:   /* if (mc != mf) SETERRQ2(1,"Processor dimensions of DA in X %d %d do not match",mc,mf);
520:      if (nc != nf) SETERRQ2(1,"Processor dimensions of DA in Y %d %d do not match",nc,nf);
521:      if (pc != pf) SETERRQ2(1,"Processor dimensions of DA in Z %d %d do not match",pc,pf); */
522:   if (dofc != doff) SETERRQ2(1,"DOF of DA do not match %d %d",dofc,doff);
523:   if (sc != sf) SETERRQ2(1,"Stencil width of DA do not match %d %d",sc,sf);
524:   if (wrapc != wrapf) SETERRQ(1,"Periodic type different in two DAs");
525:   if (stc != stf) SETERRQ(1,"Stencil type different in two DAs");
526:   if (Mc < 2) SETERRQ(1,"Coarse grid requires at least 2 points in x direction");
527:   if (dimc > 1 && Nc < 2) SETERRQ(1,"Coarse grid requires at least 2 points in y direction");
528:   if (dimc > 2 && Pc < 2) SETERRQ(1,"Coarse grid requires at least 2 points in z direction");

530:   if (dac->interptype == DA_Q1){
531:     if (dimc == 1){
532:       DAGetInterpolation_1D_Q1(dac,daf,A);
533:     } else if (dimc == 2){
534:       DAGetInterpolation_2D_Q1(dac,daf,A);
535:     } else if (dimc == 3){
536:       DAGetInterpolation_3D_Q1(dac,daf,A);
537:     } else {
538:       SETERRQ2(1,"No support for this DA dimension %d for interpolation type %d",dimc,dac->interptype);
539:     }
540:   } else if (dac->interptype == DA_Q0){
541:     if (dimc == 1){
542:       DAGetInterpolation_1D_Q0(dac,daf,A);
543:     } else {
544:       SETERRQ2(1,"No support for this DA dimension %d for interpolation type %d",dimc,dac->interptype);
545:     }
546:   }
547:   if (scale) {
548:     DMGetInterpolationScale((DM)dac,(DM)daf,*A,scale);
549:   }
550:   return(0);
551: }