Actual source code: dmfieldda.c

  1: #include <petsc/private/dmfieldimpl.h>
  2: #include <petsc/private/dmimpl.h>
  3: #include <petscdmda.h>

  5: typedef struct _n_DMField_DA {
  6:   PetscScalar *cornerVals;
  7:   PetscScalar *cornerCoeffs;
  8:   PetscScalar *work;
  9:   PetscReal    coordRange[3][2];
 10: } DMField_DA;

 12: static PetscErrorCode DMFieldDestroy_DA(DMField field)
 13: {
 14:   DMField_DA *dafield;

 16:   PetscFunctionBegin;
 17:   dafield = (DMField_DA *)field->data;
 18:   PetscCall(PetscFree3(dafield->cornerVals, dafield->cornerCoeffs, dafield->work));
 19:   PetscCall(PetscFree(dafield));
 20:   PetscFunctionReturn(PETSC_SUCCESS);
 21: }

 23: static PetscErrorCode DMFieldView_DA(DMField field, PetscViewer viewer)
 24: {
 25:   DMField_DA *dafield = (DMField_DA *)field->data;
 26:   PetscBool   iascii;

 28:   PetscFunctionBegin;
 29:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
 30:   if (iascii) {
 31:     PetscInt i, c, dim;
 32:     PetscInt nc;
 33:     DM       dm = field->dm;

 35:     PetscCall(PetscViewerASCIIPrintf(viewer, "Field corner values:\n"));
 36:     PetscCall(PetscViewerASCIIPushTab(viewer));
 37:     PetscCall(DMGetDimension(dm, &dim));
 38:     nc = field->numComponents;
 39:     for (i = 0, c = 0; i < (1 << dim); i++) {
 40:       PetscInt j;

 42:       for (j = 0; j < nc; j++, c++) {
 43:         PetscScalar val = dafield->cornerVals[nc * i + j];

 45: #if !defined(PETSC_USE_COMPLEX)
 46:         PetscCall(PetscViewerASCIIPrintf(viewer, "%g ", (double)val));
 47: #else
 48:         PetscCall(PetscViewerASCIIPrintf(viewer, "%g+i%g ", (double)PetscRealPart(val), (double)PetscImaginaryPart(val)));
 49: #endif
 50:       }
 51:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
 52:     }
 53:     PetscCall(PetscViewerASCIIPopTab(viewer));
 54:   }
 55:   PetscFunctionReturn(PETSC_SUCCESS);
 56: }

 58: #define MEdot(y, A, x, m, c, cast) \
 59:   do { \
 60:     PetscInt _k, _l; \
 61:     for (_k = 0; _k < (c); _k++) (y)[_k] = 0.; \
 62:     for (_l = 0; _l < (m); _l++) { \
 63:       for (_k = 0; _k < (c); _k++) (y)[_k] += cast((A)[(c) * _l + _k]) * (x)[_l]; \
 64:     } \
 65:   } while (0)

 67: #define MEHess(out, cf, etaB, etaD, dim, nc, cast) \
 68:   do { \
 69:     PetscInt _m, _j, _k; \
 70:     for (_m = 0; _m < (nc) * (dim) * (dim); _m++) (out)[_m] = 0.; \
 71:     for (_j = 0; _j < (dim); _j++) { \
 72:       for (_k = _j + 1; _k < (dim); _k++) { \
 73:         PetscInt _ind = (1 << _j) + (1 << _k); \
 74:         for (_m = 0; _m < (nc); _m++) { \
 75:           PetscScalar c = (cf)[_m] * (etaB)[_ind] * (etaD)[_ind]; \
 76:           (out)[(_m * (dim) + _k) * (dim) + _j] += cast(c); \
 77:           (out)[(_m * (dim) + _j) * (dim) + _k] += cast(c); \
 78:         } \
 79:       } \
 80:     } \
 81:   } while (0)

 83: static void MultilinearEvaluate(PetscInt dim, PetscReal (*coordRange)[2], PetscInt nc, PetscScalar *cf, PetscScalar *cfWork, PetscInt nPoints, const PetscScalar *points, PetscDataType datatype, void *B, void *D, void *H)
 84: {
 85:   PetscInt i, j, k, l, m;
 86:   PetscInt whol = 1 << dim;
 87:   PetscInt half = whol >> 1;

 89:   PetscFunctionBeginHot;
 90:   if (!B && !D && !H) PetscFunctionReturnVoid();
 91:   for (i = 0; i < nPoints; i++) {
 92:     const PetscScalar *point   = &points[dim * i];
 93:     PetscReal          deta[3] = {0.};
 94:     PetscReal          etaB[8] = {1., 1., 1., 1., 1., 1., 1., 1.};
 95:     PetscReal          etaD[8] = {1., 1., 1., 1., 1., 1., 1., 1.};
 96:     PetscReal          work[8];

 98:     for (j = 0; j < dim; j++) {
 99:       PetscReal e, d;

101:       e       = (PetscRealPart(point[j]) - coordRange[j][0]) / coordRange[j][1];
102:       deta[j] = d = 1. / coordRange[j][1];
103:       for (k = 0; k < whol; k++) work[k] = etaB[k];
104:       for (k = 0; k < half; k++) {
105:         etaB[k]        = work[2 * k] * e;
106:         etaB[k + half] = work[2 * k + 1];
107:       }
108:       if (H) {
109:         for (k = 0; k < whol; k++) work[k] = etaD[k];
110:         for (k = 0; k < half; k++) {
111:           etaD[k + half] = work[2 * k];
112:           etaD[k]        = work[2 * k + 1] * d;
113:         }
114:       }
115:     }
116:     if (B) {
117:       if (datatype == PETSC_SCALAR) {
118:         PetscScalar *out = &((PetscScalar *)B)[nc * i];

120:         MEdot(out, cf, etaB, (1 << dim), nc, (PetscScalar));
121:       } else {
122:         PetscReal *out = &((PetscReal *)B)[nc * i];

124:         MEdot(out, cf, etaB, (1 << dim), nc, PetscRealPart);
125:       }
126:     }
127:     if (D) {
128:       if (datatype == PETSC_SCALAR) {
129:         PetscScalar *out = &((PetscScalar *)D)[nc * dim * i];

131:         for (m = 0; m < nc * dim; m++) out[m] = 0.;
132:       } else {
133:         PetscReal *out = &((PetscReal *)D)[nc * dim * i];

135:         for (m = 0; m < nc * dim; m++) out[m] = 0.;
136:       }
137:       for (j = 0; j < dim; j++) {
138:         PetscReal d = deta[j];

140:         for (k = 0; k < whol * nc; k++) cfWork[k] = cf[k];
141:         for (k = 0; k < whol; k++) work[k] = etaB[k];
142:         for (k = 0; k < half; k++) {
143:           PetscReal e;

145:           etaB[k]        = work[2 * k];
146:           etaB[k + half] = e = work[2 * k + 1];

148:           for (l = 0; l < nc; l++) {
149:             cf[k * nc + l]          = cfWork[2 * k * nc + l];
150:             cf[(k + half) * nc + l] = cfWork[(2 * k + 1) * nc + l];
151:           }
152:           if (datatype == PETSC_SCALAR) {
153:             PetscScalar *out = &((PetscScalar *)D)[nc * dim * i];

155:             for (l = 0; l < nc; l++) out[l * dim + j] += d * e * cf[k * nc + l];
156:           } else {
157:             PetscReal *out = &((PetscReal *)D)[nc * dim * i];

159:             for (l = 0; l < nc; l++) out[l * dim + j] += d * e * PetscRealPart(cf[k * nc + l]);
160:           }
161:         }
162:       }
163:     }
164:     if (H) {
165:       if (datatype == PETSC_SCALAR) {
166:         PetscScalar *out = &((PetscScalar *)H)[nc * dim * dim * i];

168:         MEHess(out, cf, etaB, etaD, dim, nc, (PetscScalar));
169:       } else {
170:         PetscReal *out = &((PetscReal *)H)[nc * dim * dim * i];

172:         MEHess(out, cf, etaB, etaD, dim, nc, PetscRealPart);
173:       }
174:     }
175:   }
176:   PetscFunctionReturnVoid();
177: }

179: static PetscErrorCode DMFieldEvaluate_DA(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H)
180: {
181:   DM                 dm;
182:   DMField_DA        *dafield;
183:   PetscInt           dim;
184:   PetscInt           N, n, nc;
185:   const PetscScalar *array;
186:   PetscReal(*coordRange)[2];

188:   PetscFunctionBegin;
189:   dm      = field->dm;
190:   nc      = field->numComponents;
191:   dafield = (DMField_DA *)field->data;
192:   PetscCall(DMGetDimension(dm, &dim));
193:   PetscCall(VecGetLocalSize(points, &N));
194:   PetscCheck(N % dim == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Point vector size %" PetscInt_FMT " not divisible by coordinate dimension %" PetscInt_FMT, N, dim);
195:   n          = N / dim;
196:   coordRange = &dafield->coordRange[0];
197:   PetscCall(VecGetArrayRead(points, &array));
198:   MultilinearEvaluate(dim, coordRange, nc, dafield->cornerCoeffs, dafield->work, n, array, datatype, B, D, H);
199:   PetscCall(VecRestoreArrayRead(points, &array));
200:   PetscFunctionReturn(PETSC_SUCCESS);
201: }

203: static PetscErrorCode DMFieldEvaluateFE_DA(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H)
204: {
205:   PetscInt  c, i, j, k, dim, cellsPer[3] = {0}, first[3] = {0}, whol, half;
206:   PetscReal stepPer[3]           = {0.};
207:   PetscReal cellCoordRange[3][2] = {
208:     {0., 1.},
209:     {0., 1.},
210:     {0., 1.}
211:   };
212:   PetscScalar     *cellCoeffs, *work;
213:   DM               dm;
214:   DMDALocalInfo    info;
215:   PetscInt         cStart, cEnd;
216:   PetscInt         nq, nc;
217:   const PetscReal *q;
218: #if defined(PETSC_USE_COMPLEX)
219:   PetscScalar *qs;
220: #else
221:   const PetscScalar *qs;
222: #endif
223:   DMField_DA     *dafield;
224:   PetscBool       isStride;
225:   const PetscInt *cells  = NULL;
226:   PetscInt        sfirst = -1, stride = -1, nCells;

228:   PetscFunctionBegin;
229:   dafield = (DMField_DA *)field->data;
230:   dm      = field->dm;
231:   nc      = field->numComponents;
232:   PetscCall(DMDAGetLocalInfo(dm, &info));
233:   dim         = info.dim;
234:   work        = dafield->work;
235:   stepPer[0]  = 1. / info.mx;
236:   stepPer[1]  = 1. / info.my;
237:   stepPer[2]  = 1. / info.mz;
238:   first[0]    = info.gxs;
239:   first[1]    = info.gys;
240:   first[2]    = info.gzs;
241:   cellsPer[0] = info.gxm;
242:   cellsPer[1] = info.gym;
243:   cellsPer[2] = info.gzm;
244:   /* TODO: probably take components into account */
245:   PetscCall(PetscQuadratureGetData(points, NULL, NULL, &nq, &q, NULL));
246: #if defined(PETSC_USE_COMPLEX)
247:   PetscCall(DMGetWorkArray(dm, nq * dim, MPIU_SCALAR, &qs));
248:   for (i = 0; i < nq * dim; i++) qs[i] = q[i];
249: #else
250:   qs = q;
251: #endif
252:   PetscCall(DMDAGetHeightStratum(dm, 0, &cStart, &cEnd));
253:   PetscCall(DMGetWorkArray(dm, (1 << dim) * nc, MPIU_SCALAR, &cellCoeffs));
254:   whol = (1 << dim);
255:   half = whol >> 1;
256:   PetscCall(ISGetLocalSize(cellIS, &nCells));
257:   PetscCall(PetscObjectTypeCompare((PetscObject)cellIS, ISSTRIDE, &isStride));
258:   if (isStride) {
259:     PetscCall(ISStrideGetInfo(cellIS, &sfirst, &stride));
260:   } else PetscCall(ISGetIndices(cellIS, &cells));
261:   for (c = 0; c < nCells; c++) {
262:     PetscInt cell   = isStride ? (sfirst + c * stride) : cells[c];
263:     PetscInt rem    = cell;
264:     PetscInt ijk[3] = {0};
265:     void    *cB, *cD, *cH;

267:     if (datatype == PETSC_SCALAR) {
268:       cB = PetscSafePointerPlusOffset((PetscScalar *)B, nc * nq * c);
269:       cD = PetscSafePointerPlusOffset((PetscScalar *)D, nc * nq * dim * c);
270:       cH = PetscSafePointerPlusOffset((PetscScalar *)H, nc * nq * dim * dim * c);
271:     } else {
272:       cB = PetscSafePointerPlusOffset(((PetscReal *)B), nc * nq * c);
273:       cD = PetscSafePointerPlusOffset(((PetscReal *)D), nc * nq * dim * c);
274:       cH = PetscSafePointerPlusOffset(((PetscReal *)H), nc * nq * dim * dim * c);
275:     }
276:     PetscCheck(cell >= cStart && cell < cEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " not a cell [%" PetscInt_FMT ",%" PetscInt_FMT "), not implemented yet", cell, cStart, cEnd);
277:     for (i = 0; i < nc * whol; i++) work[i] = dafield->cornerCoeffs[i];
278:     for (j = 0; j < dim; j++) {
279:       PetscReal e, d;
280:       ijk[j] = (rem % cellsPer[j]);
281:       rem /= cellsPer[j];

283:       e = 2. * (ijk[j] + first[j] + 0.5) * stepPer[j] - 1.;
284:       d = stepPer[j];
285:       for (i = 0; i < half; i++) {
286:         for (k = 0; k < nc; k++) {
287:           cellCoeffs[i * nc + k]          = work[2 * i * nc + k] * d;
288:           cellCoeffs[(i + half) * nc + k] = work[2 * i * nc + k] * e + work[(2 * i + 1) * nc + k];
289:         }
290:       }
291:       for (i = 0; i < whol * nc; i++) work[i] = cellCoeffs[i];
292:     }
293:     MultilinearEvaluate(dim, cellCoordRange, nc, cellCoeffs, dafield->work, nq, qs, datatype, cB, cD, cH);
294:   }
295:   if (!isStride) PetscCall(ISRestoreIndices(cellIS, &cells));
296:   PetscCall(DMRestoreWorkArray(dm, (1 << dim) * nc, MPIU_SCALAR, &cellCoeffs));
297: #if defined(PETSC_USE_COMPLEX)
298:   PetscCall(DMRestoreWorkArray(dm, nq * dim, MPIU_SCALAR, &qs));
299: #endif
300:   PetscFunctionReturn(PETSC_SUCCESS);
301: }

303: static PetscErrorCode DMFieldEvaluateFV_DA(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H)
304: {
305:   PetscInt        c, i, dim, cellsPer[3] = {0}, first[3] = {0};
306:   PetscReal       stepPer[3] = {0.};
307:   DM              dm;
308:   DMDALocalInfo   info;
309:   PetscInt        cStart, cEnd, numCells;
310:   PetscInt        nc;
311:   PetscScalar    *points;
312:   DMField_DA     *dafield;
313:   PetscBool       isStride;
314:   const PetscInt *cells  = NULL;
315:   PetscInt        sfirst = -1, stride = -1;

317:   PetscFunctionBegin;
318:   dafield = (DMField_DA *)field->data;
319:   dm      = field->dm;
320:   nc      = field->numComponents;
321:   PetscCall(DMDAGetLocalInfo(dm, &info));
322:   dim         = info.dim;
323:   stepPer[0]  = 1. / info.mx;
324:   stepPer[1]  = 1. / info.my;
325:   stepPer[2]  = 1. / info.mz;
326:   first[0]    = info.gxs;
327:   first[1]    = info.gys;
328:   first[2]    = info.gzs;
329:   cellsPer[0] = info.gxm;
330:   cellsPer[1] = info.gym;
331:   cellsPer[2] = info.gzm;
332:   PetscCall(DMDAGetHeightStratum(dm, 0, &cStart, &cEnd));
333:   PetscCall(ISGetLocalSize(cellIS, &numCells));
334:   PetscCall(DMGetWorkArray(dm, dim * numCells, MPIU_SCALAR, &points));
335:   PetscCall(PetscObjectTypeCompare((PetscObject)cellIS, ISSTRIDE, &isStride));
336:   if (isStride) {
337:     PetscCall(ISStrideGetInfo(cellIS, &sfirst, &stride));
338:   } else PetscCall(ISGetIndices(cellIS, &cells));
339:   for (c = 0; c < numCells; c++) {
340:     PetscInt cell   = isStride ? (sfirst + c * stride) : cells[c];
341:     PetscInt rem    = cell;
342:     PetscInt ijk[3] = {0};

344:     PetscCheck(cell >= cStart && cell < cEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " not a cell [%" PetscInt_FMT ",%" PetscInt_FMT "), not implemented yet", cell, cStart, cEnd);
345:     for (i = 0; i < dim; i++) {
346:       ijk[i] = (rem % cellsPer[i]);
347:       rem /= cellsPer[i];
348:       points[dim * c + i] = (ijk[i] + first[i] + 0.5) * stepPer[i];
349:     }
350:   }
351:   if (!isStride) PetscCall(ISRestoreIndices(cellIS, &cells));
352:   MultilinearEvaluate(dim, dafield->coordRange, nc, dafield->cornerCoeffs, dafield->work, numCells, points, datatype, B, D, H);
353:   PetscCall(DMRestoreWorkArray(dm, dim * numCells, MPIU_SCALAR, &points));
354:   PetscFunctionReturn(PETSC_SUCCESS);
355: }

357: static PetscErrorCode DMFieldGetDegree_DA(DMField field, IS pointIS, PetscInt *minDegree, PetscInt *maxDegree)
358: {
359:   DM       dm;
360:   PetscInt dim, h, imin;

362:   PetscFunctionBegin;
363:   dm = field->dm;
364:   PetscCall(ISGetMinMax(pointIS, &imin, NULL));
365:   PetscCall(DMGetDimension(dm, &dim));
366:   for (h = 0; h <= dim; h++) {
367:     PetscInt hEnd;

369:     PetscCall(DMDAGetHeightStratum(dm, h, NULL, &hEnd));
370:     if (imin < hEnd) break;
371:   }
372:   dim -= h;
373:   if (minDegree) *minDegree = 1;
374:   if (maxDegree) *maxDegree = dim;
375:   PetscFunctionReturn(PETSC_SUCCESS);
376: }

378: static PetscErrorCode DMFieldCreateDefaultQuadrature_DA(DMField field, IS cellIS, PetscQuadrature *quad)
379: {
380:   PetscInt h, dim, imax, imin;
381:   DM       dm;

383:   PetscFunctionBegin;
384:   dm = field->dm;
385:   PetscCall(ISGetMinMax(cellIS, &imin, &imax));
386:   PetscCall(DMGetDimension(dm, &dim));
387:   *quad = NULL;
388:   for (h = 0; h <= dim; h++) {
389:     PetscInt hStart, hEnd;

391:     PetscCall(DMDAGetHeightStratum(dm, h, &hStart, &hEnd));
392:     if (imin >= hStart && imax < hEnd) break;
393:   }
394:   dim -= h;
395:   if (dim > 0) PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 1, -1.0, 1.0, quad));
396:   PetscFunctionReturn(PETSC_SUCCESS);
397: }

399: static PetscErrorCode DMFieldInitialize_DA(DMField field)
400: {
401:   DM          dm;
402:   Vec         coords = NULL;
403:   PetscInt    dim, i, j, k;
404:   DMField_DA *dafield = (DMField_DA *)field->data;

406:   PetscFunctionBegin;
407:   field->ops->destroy                 = DMFieldDestroy_DA;
408:   field->ops->evaluate                = DMFieldEvaluate_DA;
409:   field->ops->evaluateFE              = DMFieldEvaluateFE_DA;
410:   field->ops->evaluateFV              = DMFieldEvaluateFV_DA;
411:   field->ops->getDegree               = DMFieldGetDegree_DA;
412:   field->ops->createDefaultQuadrature = DMFieldCreateDefaultQuadrature_DA;
413:   field->ops->view                    = DMFieldView_DA;
414:   dm                                  = field->dm;
415:   PetscCall(DMGetDimension(dm, &dim));
416:   PetscCall(DMGetCoordinates(dm, &coords));
417:   if (coords) {
418:     PetscInt           n;
419:     const PetscScalar *array;
420:     PetscReal          mins[3][2] = {
421:       {PETSC_MAX_REAL, PETSC_MAX_REAL},
422:       {PETSC_MAX_REAL, PETSC_MAX_REAL},
423:       {PETSC_MAX_REAL, PETSC_MAX_REAL}
424:     };

426:     PetscCall(VecGetLocalSize(coords, &n));
427:     n /= dim;
428:     PetscCall(VecGetArrayRead(coords, &array));
429:     for (i = 0, k = 0; i < n; i++) {
430:       for (j = 0; j < dim; j++, k++) {
431:         PetscReal val = PetscRealPart(array[k]);

433:         mins[j][0] = PetscMin(mins[j][0], val);
434:         mins[j][1] = PetscMin(mins[j][1], -val);
435:       }
436:     }
437:     PetscCall(VecRestoreArrayRead(coords, &array));
438:     PetscCall(MPIU_Allreduce((PetscReal *)mins, &dafield->coordRange[0][0], 2 * dim, MPIU_REAL, MPI_MIN, PetscObjectComm((PetscObject)dm)));
439:     for (j = 0; j < dim; j++) dafield->coordRange[j][1] = -dafield->coordRange[j][1];
440:   } else {
441:     for (j = 0; j < dim; j++) {
442:       dafield->coordRange[j][0] = 0.;
443:       dafield->coordRange[j][1] = 1.;
444:     }
445:   }
446:   for (j = 0; j < dim; j++) {
447:     PetscReal avg = 0.5 * (dafield->coordRange[j][1] + dafield->coordRange[j][0]);
448:     PetscReal dif = 0.5 * (dafield->coordRange[j][1] - dafield->coordRange[j][0]);

450:     dafield->coordRange[j][0] = avg;
451:     dafield->coordRange[j][1] = dif;
452:   }
453:   PetscFunctionReturn(PETSC_SUCCESS);
454: }

456: PETSC_INTERN PetscErrorCode DMFieldCreate_DA(DMField field)
457: {
458:   DMField_DA *dafield;

460:   PetscFunctionBegin;
461:   PetscCall(PetscNew(&dafield));
462:   field->data = dafield;
463:   PetscCall(DMFieldInitialize_DA(field));
464:   PetscFunctionReturn(PETSC_SUCCESS);
465: }

467: PetscErrorCode DMFieldCreateDA(DM dm, PetscInt nc, const PetscScalar *cornerValues, DMField *field)
468: {
469:   DMField      b;
470:   DMField_DA  *dafield;
471:   PetscInt     dim, nv, i, j, k;
472:   PetscInt     half;
473:   PetscScalar *cv, *cf, *work;

475:   PetscFunctionBegin;
476:   PetscCall(DMFieldCreate(dm, nc, DMFIELD_VERTEX, &b));
477:   PetscCall(DMFieldSetType(b, DMFIELDDA));
478:   dafield = (DMField_DA *)b->data;
479:   PetscCall(DMGetDimension(dm, &dim));
480:   nv = (1 << dim) * nc;
481:   PetscCall(PetscMalloc3(nv, &cv, nv, &cf, nv, &work));
482:   for (i = 0; i < nv; i++) cv[i] = cornerValues[i];
483:   for (i = 0; i < nv; i++) cf[i] = cv[i];
484:   dafield->cornerVals   = cv;
485:   dafield->cornerCoeffs = cf;
486:   dafield->work         = work;
487:   half                  = (1 << (dim - 1));
488:   for (i = 0; i < dim; i++) {
489:     PetscScalar *w;

491:     w = work;
492:     for (j = 0; j < half; j++) {
493:       for (k = 0; k < nc; k++) w[j * nc + k] = 0.5 * (cf[(2 * j + 1) * nc + k] - cf[(2 * j) * nc + k]);
494:     }
495:     w = &work[j * nc];
496:     for (j = 0; j < half; j++) {
497:       for (k = 0; k < nc; k++) w[j * nc + k] = 0.5 * (cf[(2 * j) * nc + k] + cf[(2 * j + 1) * nc + k]);
498:     }
499:     for (j = 0; j < nv; j++) cf[j] = work[j];
500:   }
501:   *field = b;
502:   PetscFunctionReturn(PETSC_SUCCESS);
503: }