Actual source code: plexfem.c

petsc-master 2017-12-11
Report Typos and Errors
  1:  #include <petsc/private/dmpleximpl.h>
  2:  #include <petscsf.h>

  4:  #include <petsc/private/petscfeimpl.h>
  5:  #include <petsc/private/petscfvimpl.h>

  7: /*@
  8:   DMPlexGetScale - Get the scale for the specified fundamental unit

 10:   Not collective

 12:   Input Arguments:
 13: + dm   - the DM
 14: - unit - The SI unit

 16:   Output Argument:
 17: . scale - The value used to scale all quantities with this unit

 19:   Level: advanced

 21: .seealso: DMPlexSetScale(), PetscUnit
 22: @*/
 23: PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
 24: {
 25:   DM_Plex *mesh = (DM_Plex*) dm->data;

 30:   *scale = mesh->scale[unit];
 31:   return(0);
 32: }

 34: /*@
 35:   DMPlexSetScale - Set the scale for the specified fundamental unit

 37:   Not collective

 39:   Input Arguments:
 40: + dm   - the DM
 41: . unit - The SI unit
 42: - scale - The value used to scale all quantities with this unit

 44:   Level: advanced

 46: .seealso: DMPlexGetScale(), PetscUnit
 47: @*/
 48: PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
 49: {
 50:   DM_Plex *mesh = (DM_Plex*) dm->data;

 54:   mesh->scale[unit] = scale;
 55:   return(0);
 56: }

 58: static PetscErrorCode DMPlexProjectRigidBody_Private(PetscInt dim, PetscReal t, const PetscReal X[], PetscInt Nf, PetscScalar *mode, void *ctx)
 59: {
 60:   const PetscInt eps[3][3][3] = {{{0, 0, 0}, {0, 0, 1}, {0, -1, 0}}, {{0, 0, -1}, {0, 0, 0}, {1, 0, 0}}, {{0, 1, 0}, {-1, 0, 0}, {0, 0, 0}}};
 61:   PetscInt *ctxInt  = (PetscInt *) ctx;
 62:   PetscInt  dim2    = ctxInt[0];
 63:   PetscInt  d       = ctxInt[1];
 64:   PetscInt  i, j, k = dim > 2 ? d - dim : d;

 67:   if (dim != dim2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input dimension %d does not match context dimension %d", dim, dim2);
 68:   for (i = 0; i < dim; i++) mode[i] = 0.;
 69:   if (d < dim) {
 70:     mode[d] = 1.; /* Translation along axis d */
 71:   } else {
 72:     for (i = 0; i < dim; i++) {
 73:       for (j = 0; j < dim; j++) {
 74:         mode[j] += eps[i][j][k]*X[i]; /* Rotation about axis d */
 75:       }
 76:     }
 77:   }
 78:   return(0);
 79: }

 81: /*@C
 82:   DMPlexCreateRigidBody - For the default global section, create rigid body modes by function space interpolation

 84:   Collective on DM

 86:   Input Arguments:
 87: . dm - the DM

 89:   Output Argument:
 90: . sp - the null space

 92:   Note: This is necessary to provide a suitable coarse space for algebraic multigrid

 94:   Level: advanced

 96: .seealso: MatNullSpaceCreate(), PCGAMG
 97: @*/
 98: PetscErrorCode DMPlexCreateRigidBody(DM dm, MatNullSpace *sp)
 99: {
100:   MPI_Comm       comm;
101:   Vec            mode[6];
102:   PetscSection   section, globalSection;
103:   PetscInt       dim, dimEmbed, n, m, d, i, j;

107:   PetscObjectGetComm((PetscObject)dm,&comm);
108:   DMGetDimension(dm, &dim);
109:   DMGetCoordinateDim(dm, &dimEmbed);
110:   if (dim == 1) {
111:     MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);
112:     return(0);
113:   }
114:   DMGetDefaultSection(dm, &section);
115:   DMGetDefaultGlobalSection(dm, &globalSection);
116:   PetscSectionGetConstrainedStorageSize(globalSection, &n);
117:   m    = (dim*(dim+1))/2;
118:   VecCreate(comm, &mode[0]);
119:   VecSetSizes(mode[0], n, PETSC_DETERMINE);
120:   VecSetUp(mode[0]);
121:   for (i = 1; i < m; ++i) {VecDuplicate(mode[0], &mode[i]);}
122:   for (d = 0; d < m; d++) {
123:     PetscInt         ctx[2];
124:     PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *) = DMPlexProjectRigidBody_Private;
125:     void            *voidctx = (void *) (&ctx[0]);

127:     ctx[0] = dimEmbed;
128:     ctx[1] = d;
129:     DMProjectFunction(dm, 0.0, &func, &voidctx, INSERT_VALUES, mode[d]);
130:   }
131:   for (i = 0; i < dim; ++i) {VecNormalize(mode[i], NULL);}
132:   /* Orthonormalize system */
133:   for (i = dim; i < m; ++i) {
134:     PetscScalar dots[6];

136:     VecMDot(mode[i], i, mode, dots);
137:     for (j = 0; j < i; ++j) dots[j] *= -1.0;
138:     VecMAXPY(mode[i], i, dots, mode);
139:     VecNormalize(mode[i], NULL);
140:   }
141:   MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);
142:   for (i = 0; i< m; ++i) {VecDestroy(&mode[i]);}
143:   return(0);
144: }

146: /*@C
147:   DMPlexCreateRigidBodies - For the default global section, create rigid body modes by function space interpolation

149:   Collective on DM

151:   Input Arguments:
152: + dm    - the DM
153: . nb    - The number of bodies
154: . label - The DMLabel marking each domain
155: . nids  - The number of ids per body
156: - ids   - An array of the label ids in sequence for each domain

158:   Output Argument:
159: . sp - the null space

161:   Note: This is necessary to provide a suitable coarse space for algebraic multigrid

163:   Level: advanced

165: .seealso: MatNullSpaceCreate()
166: @*/
167: PetscErrorCode DMPlexCreateRigidBodies(DM dm, PetscInt nb, DMLabel label, const PetscInt nids[], const PetscInt ids[], MatNullSpace *sp)
168: {
169:   MPI_Comm       comm;
170:   PetscSection   section, globalSection;
171:   Vec           *mode;
172:   PetscScalar   *dots;
173:   PetscInt       dim, dimEmbed, n, m, b, d, i, j, off;

177:   PetscObjectGetComm((PetscObject)dm,&comm);
178:   DMGetDimension(dm, &dim);
179:   DMGetCoordinateDim(dm, &dimEmbed);
180:   DMGetDefaultSection(dm, &section);
181:   DMGetDefaultGlobalSection(dm, &globalSection);
182:   PetscSectionGetConstrainedStorageSize(globalSection, &n);
183:   m    = nb * (dim*(dim+1))/2;
184:   PetscMalloc2(m, &mode, m, &dots);
185:   VecCreate(comm, &mode[0]);
186:   VecSetSizes(mode[0], n, PETSC_DETERMINE);
187:   VecSetUp(mode[0]);
188:   for (i = 1; i < m; ++i) {VecDuplicate(mode[0], &mode[i]);}
189:   for (b = 0, off = 0; b < nb; ++b) {
190:     for (d = 0; d < m/nb; ++d) {
191:       PetscInt         ctx[2];
192:       PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *) = DMPlexProjectRigidBody_Private;
193:       void            *voidctx = (void *) (&ctx[0]);

195:       ctx[0] = dimEmbed;
196:       ctx[1] = d;
197:       DMProjectFunctionLabel(dm, 0.0, label, nids[b], &ids[off], 0, NULL, &func, &voidctx, INSERT_VALUES, mode[d]);
198:       off   += nids[b];
199:     }
200:   }
201:   for (i = 0; i < dim; ++i) {VecNormalize(mode[i], NULL);}
202:   /* Orthonormalize system */
203:   for (i = 0; i < m; ++i) {
204:     VecMDot(mode[i], i, mode, dots);
205:     for (j = 0; j < i; ++j) dots[j] *= -1.0;
206:     VecMAXPY(mode[i], i, dots, mode);
207:     VecNormalize(mode[i], NULL);
208:   }
209:   MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);
210:   for (i = 0; i< m; ++i) {VecDestroy(&mode[i]);}
211:   PetscFree2(mode, dots);
212:   return(0);
213: }

215: /*@
216:   DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs
217:   are computed by associating the basis function with one of the mesh points in its transitively-closed support, and
218:   evaluating the dual space basis of that point.  A basis function is associated with the point in its
219:   transitively-closed support whose mesh height is highest (w.r.t. DAG height), but not greater than the maximum
220:   projection height, which is set with this function.  By default, the maximum projection height is zero, which means
221:   that only mesh cells are used to project basis functions.  A height of one, for example, evaluates a cell-interior
222:   basis functions using its cells dual space basis, but all other basis functions with the dual space basis of a face.

224:   Input Parameters:
225: + dm - the DMPlex object
226: - height - the maximum projection height >= 0

228:   Level: advanced

230: .seealso: DMPlexGetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal()
231: @*/
232: PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height)
233: {
234:   DM_Plex *plex = (DM_Plex *) dm->data;

238:   plex->maxProjectionHeight = height;
239:   return(0);
240: }

242: /*@
243:   DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in
244:   DMPlexProjectXXXLocal() functions.

246:   Input Parameters:
247: . dm - the DMPlex object

249:   Output Parameters:
250: . height - the maximum projection height

252:   Level: intermediate

254: .seealso: DMPlexSetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal()
255: @*/
256: PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height)
257: {
258:   DM_Plex *plex = (DM_Plex *) dm->data;

262:   *height = plex->maxProjectionHeight;
263:   return(0);
264: }

266: /*@C
267:   DMPlexInsertBoundaryValuesEssential - Insert boundary values into a local vector

269:   Input Parameters:
270: + dm     - The DM, with a PetscDS that matches the problem being constrained
271: . time   - The time
272: . field  - The field to constrain
273: . Nc     - The number of constrained field components, or 0 for all components
274: . comps  - An array of constrained component numbers, or NULL for all components
275: . label  - The DMLabel defining constrained points
276: . numids - The number of DMLabel ids for constrained points
277: . ids    - An array of ids for constrained points
278: . func   - A pointwise function giving boundary values
279: - ctx    - An optional user context for bcFunc

281:   Output Parameter:
282: . locX   - A local vector to receives the boundary values

284:   Level: developer

286: .seealso: DMPlexInsertBoundaryValuesEssentialField(), DMAddBoundary()
287: @*/
288: PetscErrorCode DMPlexInsertBoundaryValuesEssential(DM dm, PetscReal time, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx, Vec locX)
289: {
290:   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
291:   void            **ctxs;
292:   PetscInt          numFields;
293:   PetscErrorCode    ierr;

296:   DMGetNumFields(dm, &numFields);
297:   PetscCalloc2(numFields,&funcs,numFields,&ctxs);
298:   funcs[field] = func;
299:   ctxs[field]  = ctx;
300:   DMProjectFunctionLabelLocal(dm, time, label, numids, ids, Nc, comps, funcs, ctxs, INSERT_BC_VALUES, locX);
301:   PetscFree2(funcs,ctxs);
302:   return(0);
303: }

305: /*@C
306:   DMPlexInsertBoundaryValuesEssentialField - Insert boundary values into a local vector

308:   Input Parameters:
309: + dm     - The DM, with a PetscDS that matches the problem being constrained
310: . time   - The time
311: . locU   - A local vector with the input solution values
312: . field  - The field to constrain
313: . Nc     - The number of constrained field components, or 0 for all components
314: . comps  - An array of constrained component numbers, or NULL for all components
315: . label  - The DMLabel defining constrained points
316: . numids - The number of DMLabel ids for constrained points
317: . ids    - An array of ids for constrained points
318: . func   - A pointwise function giving boundary values
319: - ctx    - An optional user context for bcFunc

321:   Output Parameter:
322: . locX   - A local vector to receives the boundary values

324:   Level: developer

326: .seealso: DMPlexInsertBoundaryValuesEssential(), DMAddBoundary()
327: @*/
328: PetscErrorCode DMPlexInsertBoundaryValuesEssentialField(DM dm, PetscReal time, Vec locU, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[],
329:                                                         void (*func)(PetscInt, PetscInt, PetscInt,
330:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
331:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
332:                                                                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[],
333:                                                                      PetscScalar[]),
334:                                                         void *ctx, Vec locX)
335: {
336:   void (**funcs)(PetscInt, PetscInt, PetscInt,
337:                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
338:                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
339:                  PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
340:   void            **ctxs;
341:   PetscInt          numFields;
342:   PetscErrorCode    ierr;

345:   DMGetNumFields(dm, &numFields);
346:   PetscCalloc2(numFields,&funcs,numFields,&ctxs);
347:   funcs[field] = func;
348:   ctxs[field]  = ctx;
349:   DMProjectFieldLabelLocal(dm, time, label, numids, ids, Nc, comps, locU, funcs, INSERT_BC_VALUES, locX);
350:   PetscFree2(funcs,ctxs);
351:   return(0);
352: }

354: /*@C
355:   DMPlexInsertBoundaryValuesRiemann - Insert boundary values into a local vector

357:   Input Parameters:
358: + dm     - The DM, with a PetscDS that matches the problem being constrained
359: . time   - The time
360: . faceGeometry - A vector with the FVM face geometry information
361: . cellGeometry - A vector with the FVM cell geometry information
362: . Grad         - A vector with the FVM cell gradient information
363: . field  - The field to constrain
364: . Nc     - The number of constrained field components, or 0 for all components
365: . comps  - An array of constrained component numbers, or NULL for all components
366: . label  - The DMLabel defining constrained points
367: . numids - The number of DMLabel ids for constrained points
368: . ids    - An array of ids for constrained points
369: . func   - A pointwise function giving boundary values
370: - ctx    - An optional user context for bcFunc

372:   Output Parameter:
373: . locX   - A local vector to receives the boundary values

375:   Note: This implementation currently ignores the numcomps/comps argument from DMAddBoundary()

377:   Level: developer

379: .seealso: DMPlexInsertBoundaryValuesEssential(), DMPlexInsertBoundaryValuesEssentialField(), DMAddBoundary()
380: @*/
381: PetscErrorCode DMPlexInsertBoundaryValuesRiemann(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[],
382:                                                  PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *ctx, Vec locX)
383: {
384:   PetscDS            prob;
385:   PetscSF            sf;
386:   DM                 dmFace, dmCell, dmGrad;
387:   const PetscScalar *facegeom, *cellgeom = NULL, *grad;
388:   const PetscInt    *leaves;
389:   PetscScalar       *x, *fx;
390:   PetscInt           dim, nleaves, loc, fStart, fEnd, pdim, i;
391:   PetscErrorCode     ierr, ierru = 0;

394:   DMGetPointSF(dm, &sf);
395:   PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);
396:   nleaves = PetscMax(0, nleaves);
397:   DMGetDimension(dm, &dim);
398:   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
399:   DMGetDS(dm, &prob);
400:   VecGetDM(faceGeometry, &dmFace);
401:   VecGetArrayRead(faceGeometry, &facegeom);
402:   if (cellGeometry) {
403:     VecGetDM(cellGeometry, &dmCell);
404:     VecGetArrayRead(cellGeometry, &cellgeom);
405:   }
406:   if (Grad) {
407:     PetscFV fv;

409:     PetscDSGetDiscretization(prob, field, (PetscObject *) &fv);
410:     VecGetDM(Grad, &dmGrad);
411:     VecGetArrayRead(Grad, &grad);
412:     PetscFVGetNumComponents(fv, &pdim);
413:     DMGetWorkArray(dm, pdim, MPIU_SCALAR, &fx);
414:   }
415:   VecGetArray(locX, &x);
416:   for (i = 0; i < numids; ++i) {
417:     IS              faceIS;
418:     const PetscInt *faces;
419:     PetscInt        numFaces, f;

421:     DMLabelGetStratumIS(label, ids[i], &faceIS);
422:     if (!faceIS) continue; /* No points with that id on this process */
423:     ISGetLocalSize(faceIS, &numFaces);
424:     ISGetIndices(faceIS, &faces);
425:     for (f = 0; f < numFaces; ++f) {
426:       const PetscInt         face = faces[f], *cells;
427:       PetscFVFaceGeom        *fg;

429:       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
430:       PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);
431:       if (loc >= 0) continue;
432:       DMPlexPointLocalRead(dmFace, face, facegeom, &fg);
433:       DMPlexGetSupport(dm, face, &cells);
434:       if (Grad) {
435:         PetscFVCellGeom       *cg;
436:         PetscScalar           *cx, *cgrad;
437:         PetscScalar           *xG;
438:         PetscReal              dx[3];
439:         PetscInt               d;

441:         DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);
442:         DMPlexPointLocalRead(dm, cells[0], x, &cx);
443:         DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);
444:         DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);
445:         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
446:         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
447:         ierru = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);
448:         if (ierru) {
449:           ISRestoreIndices(faceIS, &faces);
450:           ISDestroy(&faceIS);
451:           goto cleanup;
452:         }
453:       } else {
454:         PetscScalar       *xI;
455:         PetscScalar       *xG;

457:         DMPlexPointLocalRead(dm, cells[0], x, &xI);
458:         DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);
459:         ierru = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);
460:         if (ierru) {
461:           ISRestoreIndices(faceIS, &faces);
462:           ISDestroy(&faceIS);
463:           goto cleanup;
464:         }
465:       }
466:     }
467:     ISRestoreIndices(faceIS, &faces);
468:     ISDestroy(&faceIS);
469:   }
470:   cleanup:
471:   VecRestoreArray(locX, &x);
472:   if (Grad) {
473:     DMRestoreWorkArray(dm, pdim, MPIU_SCALAR, &fx);
474:     VecRestoreArrayRead(Grad, &grad);
475:   }
476:   if (cellGeometry) {VecRestoreArrayRead(cellGeometry, &cellgeom);}
477:   VecRestoreArrayRead(faceGeometry, &facegeom);
478:   CHKERRQ(ierru);
479:   return(0);
480: }

482: PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
483: {
484:   PetscInt       numBd, b;

488:   PetscDSGetNumBoundary(dm->prob, &numBd);
489:   for (b = 0; b < numBd; ++b) {
490:     DMBoundaryConditionType type;
491:     const char             *labelname;
492:     DMLabel                 label;
493:     PetscInt                field, Nc;
494:     const PetscInt         *comps;
495:     PetscObject             obj;
496:     PetscClassId            id;
497:     void                    (*func)(void);
498:     PetscInt                numids;
499:     const PetscInt         *ids;
500:     void                   *ctx;

502:     DMGetBoundary(dm, b, &type, NULL, &labelname, &field, &Nc, &comps, &func, &numids, &ids, &ctx);
503:     if (insertEssential != (type & DM_BC_ESSENTIAL)) continue;
504:     DMGetLabel(dm, labelname, &label);
505:     DMGetField(dm, field, &obj);
506:     PetscObjectGetClassId(obj, &id);
507:     if (id == PETSCFE_CLASSID) {
508:       switch (type) {
509:         /* for FEM, there is no insertion to be done for non-essential boundary conditions */
510:       case DM_BC_ESSENTIAL:
511:         DMPlexLabelAddCells(dm,label);
512:         DMPlexInsertBoundaryValuesEssential(dm, time, field, Nc, comps, label, numids, ids, (PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *)) func, ctx, locX);
513:         DMPlexLabelClearCells(dm,label);
514:         break;
515:       case DM_BC_ESSENTIAL_FIELD:
516:         DMPlexLabelAddCells(dm,label);
517:         DMPlexInsertBoundaryValuesEssentialField(dm, time, locX, field, Nc, comps, label, numids, ids,
518:                                                         (void (*)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
519:                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
520:                                                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) func, ctx, locX);
521:         DMPlexLabelClearCells(dm,label);
522:         break;
523:       default: break;
524:       }
525:     } else if (id == PETSCFV_CLASSID) {
526:       if (!faceGeomFVM) continue;
527:       DMPlexInsertBoundaryValuesRiemann(dm, time, faceGeomFVM, cellGeomFVM, gradFVM, field, Nc, comps, label, numids, ids,
528:                                                (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);
529:     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
530:   }
531:   return(0);
532: }

534: /*@
535:   DMPlexInsertBoundaryValues - Puts coefficients which represent boundary values into the local solution vector

537:   Input Parameters:
538: + dm - The DM
539: . insertEssential - Should I insert essential (e.g. Dirichlet) or inessential (e.g. Neumann) boundary conditions
540: . time - The time
541: . faceGeomFVM - Face geometry data for FV discretizations
542: . cellGeomFVM - Cell geometry data for FV discretizations
543: - gradFVM - Gradient reconstruction data for FV discretizations

545:   Output Parameters:
546: . locX - Solution updated with boundary values

548:   Level: developer

550: .seealso: DMProjectFunctionLabelLocal()
551: @*/
552: PetscErrorCode DMPlexInsertBoundaryValues(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
553: {

562:   PetscTryMethod(dm,"DMPlexInsertBoundaryValues_C",(DM,PetscBool,Vec,PetscReal,Vec,Vec,Vec),(dm,insertEssential,locX,time,faceGeomFVM,cellGeomFVM,gradFVM));
563:   return(0);
564: }

566: PetscErrorCode DMComputeL2Diff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
567: {
568:   Vec              localX;
569:   PetscErrorCode   ierr;

572:   DMGetLocalVector(dm, &localX);
573:   DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localX, time, NULL, NULL, NULL);
574:   DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);
575:   DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);
576:   DMPlexComputeL2DiffLocal(dm, time, funcs, ctxs, localX, diff);
577:   DMRestoreLocalVector(dm, &localX);
578:   return(0);
579: }

581: /*@C
582:   DMComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.

584:   Input Parameters:
585: + dm     - The DM
586: . time   - The time
587: . funcs  - The functions to evaluate for each field component
588: . ctxs   - Optional array of contexts to pass to each function, or NULL.
589: - localX - The coefficient vector u_h, a local vector

591:   Output Parameter:
592: . diff - The diff ||u - u_h||_2

594:   Level: developer

596: .seealso: DMProjectFunction(), DMComputeL2FieldDiff(), DMComputeL2GradientDiff()
597: @*/
598: PetscErrorCode DMPlexComputeL2DiffLocal(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec localX, PetscReal *diff)
599: {
600:   const PetscInt   debug = 0;
601:   PetscSection     section;
602:   PetscQuadrature  quad;
603:   PetscScalar     *funcVal, *interpolant;
604:   PetscReal       *coords, *detJ, *J;
605:   PetscReal        localDiff = 0.0;
606:   const PetscReal *quadWeights;
607:   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cellHeight, cStart, cEnd, cEndInterior, c, field, fieldOffset;
608:   PetscErrorCode   ierr;

611:   DMGetDimension(dm, &dim);
612:   DMGetCoordinateDim(dm, &coordDim);
613:   DMGetDefaultSection(dm, &section);
614:   PetscSectionGetNumFields(section, &numFields);
615:   for (field = 0; field < numFields; ++field) {
616:     PetscObject  obj;
617:     PetscClassId id;
618:     PetscInt     Nc;

620:     DMGetField(dm, field, &obj);
621:     PetscObjectGetClassId(obj, &id);
622:     if (id == PETSCFE_CLASSID) {
623:       PetscFE fe = (PetscFE) obj;

625:       PetscFEGetQuadrature(fe, &quad);
626:       PetscFEGetNumComponents(fe, &Nc);
627:     } else if (id == PETSCFV_CLASSID) {
628:       PetscFV fv = (PetscFV) obj;

630:       PetscFVGetQuadrature(fv, &quad);
631:       PetscFVGetNumComponents(fv, &Nc);
632:     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
633:     numComponents += Nc;
634:   }
635:   PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights);
636:   if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents);
637:   PetscMalloc5(numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);
638:   DMPlexGetVTKCellHeight(dm, &cellHeight);
639:   DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);
640:   DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
641:   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
642:   for (c = cStart; c < cEnd; ++c) {
643:     PetscScalar *x = NULL;
644:     PetscReal    elemDiff = 0.0;
645:     PetscInt     qc = 0;

647:     DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);
648:     DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);

650:     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
651:       PetscObject  obj;
652:       PetscClassId id;
653:       void * const ctx = ctxs ? ctxs[field] : NULL;
654:       PetscInt     Nb, Nc, q, fc;

656:       DMGetField(dm, field, &obj);
657:       PetscObjectGetClassId(obj, &id);
658:       if (id == PETSCFE_CLASSID)      {PetscFEGetNumComponents((PetscFE) obj, &Nc);PetscFEGetDimension((PetscFE) obj, &Nb);}
659:       else if (id == PETSCFV_CLASSID) {PetscFVGetNumComponents((PetscFV) obj, &Nc);Nb = 1;}
660:       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
661:       if (debug) {
662:         char title[1024];
663:         PetscSNPrintf(title, 1023, "Solution for Field %d", field);
664:         DMPrintCellVector(c, title, Nb, &x[fieldOffset]);
665:       }
666:       for (q = 0; q < Nq; ++q) {
667:         if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, point %D", detJ[q], c, q);
668:         (*funcs[field])(coordDim, time, &coords[coordDim * q], Nc, funcVal, ctx);
669:         if (ierr) {
670:           PetscErrorCode ierr2;
671:           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
672:           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
673:           ierr2 = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr2);
674: 
675:         }
676:         if (id == PETSCFE_CLASSID)      {PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);}
677:         else if (id == PETSCFV_CLASSID) {PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);}
678:         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
679:         for (fc = 0; fc < Nc; ++fc) {
680:           const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)];
681:           if (debug) {PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]);}
682:           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q];
683:         }
684:       }
685:       fieldOffset += Nb;
686:       qc += Nc;
687:     }
688:     DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);
689:     if (debug) {PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);}
690:     localDiff += elemDiff;
691:   }
692:   PetscFree5(funcVal,interpolant,coords,detJ,J);
693:   MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));
694:   *diff = PetscSqrtReal(*diff);
695:   return(0);
696: }

698: PetscErrorCode DMComputeL2GradientDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
699: {
700:   const PetscInt   debug = 0;
701:   PetscSection     section;
702:   PetscQuadrature  quad;
703:   Vec              localX;
704:   PetscScalar     *funcVal, *interpolant;
705:   const PetscReal *quadPoints, *quadWeights;
706:   PetscReal       *coords, *realSpaceDer, *J, *invJ, *detJ;
707:   PetscReal        localDiff = 0.0;
708:   PetscInt         dim, coordDim, qNc = 0, Nq = 0, numFields, numComponents = 0, cStart, cEnd, cEndInterior, c, field, fieldOffset;
709:   PetscErrorCode   ierr;

712:   DMGetDimension(dm, &dim);
713:   DMGetCoordinateDim(dm, &coordDim);
714:   DMGetDefaultSection(dm, &section);
715:   PetscSectionGetNumFields(section, &numFields);
716:   DMGetLocalVector(dm, &localX);
717:   DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);
718:   DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);
719:   for (field = 0; field < numFields; ++field) {
720:     PetscFE  fe;
721:     PetscInt Nc;

723:     DMGetField(dm, field, (PetscObject *) &fe);
724:     PetscFEGetQuadrature(fe, &quad);
725:     PetscFEGetNumComponents(fe, &Nc);
726:     numComponents += Nc;
727:   }
728:   PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
729:   if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents);
730:   /* DMProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX); */
731:   PetscMalloc7(numComponents,&funcVal,coordDim*Nq,&coords,coordDim*Nq,&realSpaceDer,coordDim*coordDim*Nq,&J,coordDim*coordDim*Nq,&invJ,numComponents,&interpolant,Nq,&detJ);
732:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
733:   DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
734:   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
735:   for (c = cStart; c < cEnd; ++c) {
736:     PetscScalar *x = NULL;
737:     PetscReal    elemDiff = 0.0;
738:     PetscInt     qc = 0;

740:     DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, invJ, detJ);
741:     DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);

743:     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
744:       PetscFE          fe;
745:       void * const     ctx = ctxs ? ctxs[field] : NULL;
746:       PetscReal       *basisDer;
747:       PetscInt         Nb, Nc, q, fc;

749:       DMGetField(dm, field, (PetscObject *) &fe);
750:       PetscFEGetDimension(fe, &Nb);
751:       PetscFEGetNumComponents(fe, &Nc);
752:       PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);
753:       if (debug) {
754:         char title[1024];
755:         PetscSNPrintf(title, 1023, "Solution for Field %d", field);
756:         DMPrintCellVector(c, title, Nb, &x[fieldOffset]);
757:       }
758:       for (q = 0; q < Nq; ++q) {
759:         if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", detJ[q], c, q);
760:         (*funcs[field])(coordDim, time, &coords[q*coordDim], n, numFields, funcVal, ctx);
761:         if (ierr) {
762:           PetscErrorCode ierr2;
763:           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
764:           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
765:           ierr2 = PetscFree7(funcVal,coords,realSpaceDer,J,invJ,interpolant,detJ);CHKERRQ(ierr2);
766: 
767:         }
768:         PetscFEInterpolateGradient_Static(fe, &x[fieldOffset], coordDim, invJ, n, q, interpolant);
769:         for (fc = 0; fc < Nc; ++fc) {
770:           const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)];
771:           if (debug) {PetscPrintf(PETSC_COMM_SELF, "    elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]);}
772:           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q];
773:         }
774:       }
775:       fieldOffset += Nb;
776:       qc          += Nc;
777:     }
778:     DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);
779:     if (debug) {PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);}
780:     localDiff += elemDiff;
781:   }
782:   PetscFree7(funcVal,coords,realSpaceDer,J,invJ,interpolant,detJ);
783:   DMRestoreLocalVector(dm, &localX);
784:   MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));
785:   *diff = PetscSqrtReal(*diff);
786:   return(0);
787: }

789: PetscErrorCode DMComputeL2FieldDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
790: {
791:   const PetscInt   debug = 0;
792:   PetscSection     section;
793:   PetscQuadrature  quad;
794:   Vec              localX;
795:   PetscScalar     *funcVal, *interpolant;
796:   PetscReal       *coords, *detJ, *J;
797:   PetscReal       *localDiff;
798:   const PetscReal *quadPoints, *quadWeights;
799:   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, cEndInterior, c, field, fieldOffset;
800:   PetscErrorCode   ierr;

803:   DMGetDimension(dm, &dim);
804:   DMGetCoordinateDim(dm, &coordDim);
805:   DMGetDefaultSection(dm, &section);
806:   PetscSectionGetNumFields(section, &numFields);
807:   DMGetLocalVector(dm, &localX);
808:   DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);
809:   DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);
810:   DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);
811:   for (field = 0; field < numFields; ++field) {
812:     PetscObject  obj;
813:     PetscClassId id;
814:     PetscInt     Nc;

816:     DMGetField(dm, field, &obj);
817:     PetscObjectGetClassId(obj, &id);
818:     if (id == PETSCFE_CLASSID) {
819:       PetscFE fe = (PetscFE) obj;

821:       PetscFEGetQuadrature(fe, &quad);
822:       PetscFEGetNumComponents(fe, &Nc);
823:     } else if (id == PETSCFV_CLASSID) {
824:       PetscFV fv = (PetscFV) obj;

826:       PetscFVGetQuadrature(fv, &quad);
827:       PetscFVGetNumComponents(fv, &Nc);
828:     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
829:     numComponents += Nc;
830:   }
831:   PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
832:   if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents);
833:   PetscCalloc6(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);
834:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
835:   DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
836:   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
837:   for (c = cStart; c < cEnd; ++c) {
838:     PetscScalar *x = NULL;
839:     PetscInt     qc = 0;

841:     DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);
842:     DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);

844:     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
845:       PetscObject  obj;
846:       PetscClassId id;
847:       void * const ctx = ctxs ? ctxs[field] : NULL;
848:       PetscInt     Nb, Nc, q, fc;

850:       PetscReal       elemDiff = 0.0;

852:       DMGetField(dm, field, &obj);
853:       PetscObjectGetClassId(obj, &id);
854:       if (id == PETSCFE_CLASSID)      {PetscFEGetNumComponents((PetscFE) obj, &Nc);PetscFEGetDimension((PetscFE) obj, &Nb);}
855:       else if (id == PETSCFV_CLASSID) {PetscFVGetNumComponents((PetscFV) obj, &Nc);Nb = 1;}
856:       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
857:       if (debug) {
858:         char title[1024];
859:         PetscSNPrintf(title, 1023, "Solution for Field %d", field);
860:         DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);
861:       }
862:       for (q = 0; q < Nq; ++q) {
863:         if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature point %D", detJ, c, q);
864:         (*funcs[field])(coordDim, time, &coords[coordDim*q], numFields, funcVal, ctx);
865:         if (ierr) {
866:           PetscErrorCode ierr2;
867:           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
868:           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
869:           ierr2 = PetscFree6(localDiff,funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr2);
870: 
871:         }
872:         if (id == PETSCFE_CLASSID)      {PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);}
873:         else if (id == PETSCFV_CLASSID) {PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);}
874:         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
875:         for (fc = 0; fc < Nc; ++fc) {
876:           const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)];
877:           if (debug) {PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d point %g %g %g diff %g\n", c, field, coordDim > 0 ? coords[0] : 0., coordDim > 1 ? coords[1] : 0., coordDim > 2 ? coords[2] : 0., PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]);}
878:           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q];
879:         }
880:       }
881:       fieldOffset += Nb;
882:       qc          += Nc;
883:       localDiff[field] += elemDiff;
884:     }
885:     DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);
886:   }
887:   DMRestoreLocalVector(dm, &localX);
888:   MPIU_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));
889:   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
890:   PetscFree6(localDiff,funcVal,interpolant,coords,detJ,J);
891:   return(0);
892: }

894: /*@C
895:   DMPlexComputeL2DiffVec - This function computes the cellwise L_2 difference between a function u and an FEM interpolant solution u_h, and stores it in a Vec.

897:   Input Parameters:
898: + dm    - The DM
899: . time  - The time
900: . funcs - The functions to evaluate for each field component: NULL means that component does not contribute to error calculation
901: . ctxs  - Optional array of contexts to pass to each function, or NULL.
902: - X     - The coefficient vector u_h

904:   Output Parameter:
905: . D - A Vec which holds the difference ||u - u_h||_2 for each cell

907:   Level: developer

909: .seealso: DMProjectFunction(), DMComputeL2Diff(), DMPlexComputeL2FieldDiff(), DMComputeL2GradientDiff()
910: @*/
911: PetscErrorCode DMPlexComputeL2DiffVec(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, Vec D)
912: {
913:   PetscSection     section;
914:   PetscQuadrature  quad;
915:   Vec              localX;
916:   PetscScalar     *funcVal, *interpolant;
917:   PetscReal       *coords, *detJ, *J;
918:   const PetscReal *quadPoints, *quadWeights;
919:   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, cEndInterior, c, field, fieldOffset;
920:   PetscErrorCode   ierr;

923:   VecSet(D, 0.0);
924:   DMGetDimension(dm, &dim);
925:   DMGetCoordinateDim(dm, &coordDim);
926:   DMGetDefaultSection(dm, &section);
927:   PetscSectionGetNumFields(section, &numFields);
928:   DMGetLocalVector(dm, &localX);
929:   DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);
930:   DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);
931:   DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);
932:   for (field = 0; field < numFields; ++field) {
933:     PetscObject  obj;
934:     PetscClassId id;
935:     PetscInt     Nc;

937:     DMGetField(dm, field, &obj);
938:     PetscObjectGetClassId(obj, &id);
939:     if (id == PETSCFE_CLASSID) {
940:       PetscFE fe = (PetscFE) obj;

942:       PetscFEGetQuadrature(fe, &quad);
943:       PetscFEGetNumComponents(fe, &Nc);
944:     } else if (id == PETSCFV_CLASSID) {
945:       PetscFV fv = (PetscFV) obj;

947:       PetscFVGetQuadrature(fv, &quad);
948:       PetscFVGetNumComponents(fv, &Nc);
949:     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
950:     numComponents += Nc;
951:   }
952:   PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
953:   if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents);
954:   PetscMalloc5(numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);
955:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
956:   DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
957:   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
958:   for (c = cStart; c < cEnd; ++c) {
959:     PetscScalar *x = NULL;
960:     PetscScalar  elemDiff = 0.0;
961:     PetscInt     qc = 0;

963:     DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);
964:     DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);

966:     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
967:       PetscObject  obj;
968:       PetscClassId id;
969:       void * const ctx = ctxs ? ctxs[field] : NULL;
970:       PetscInt     Nb, Nc, q, fc;

972:       DMGetField(dm, field, &obj);
973:       PetscObjectGetClassId(obj, &id);
974:       if (id == PETSCFE_CLASSID)      {PetscFEGetNumComponents((PetscFE) obj, &Nc);PetscFEGetDimension((PetscFE) obj, &Nb);}
975:       else if (id == PETSCFV_CLASSID) {PetscFVGetNumComponents((PetscFV) obj, &Nc);Nb = 1;}
976:       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
977:       if (funcs[field]) {
978:         for (q = 0; q < Nq; ++q) {
979:           if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", (double)detJ[q], c, q);
980:           (*funcs[field])(dim, time, &coords[q*coordDim], Nc, funcVal, ctx);
981:           if (ierr) {
982:             PetscErrorCode ierr2;
983:             ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
984:             ierr2 = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr2);
985:             ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
986: 
987:           }
988:           if (id == PETSCFE_CLASSID)      {PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);}
989:           else if (id == PETSCFV_CLASSID) {PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);}
990:           else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
991:           for (fc = 0; fc < Nc; ++fc) {
992:             const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)];
993:             elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q];
994:           }
995:         }
996:       }
997:       fieldOffset += Nb;
998:       qc          += Nc;
999:     }
1000:     DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);
1001:     VecSetValue(D, c - cStart, elemDiff, INSERT_VALUES);
1002:   }
1003:   PetscFree5(funcVal,interpolant,coords,detJ,J);
1004:   DMRestoreLocalVector(dm, &localX);
1005:   VecSqrtAbs(D);
1006:   return(0);
1007: }

1009: /*@C
1010:   DMPlexComputeGradientClementInterpolant - This function computes the L2 projection of the cellwise gradient of a function u onto P1, and stores it in a Vec.

1012:   Input Parameters:
1013: + dm - The DM
1014: - LocX  - The coefficient vector u_h

1016:   Output Parameter:
1017: . locC - A Vec which holds the Clement interpolant of the gradient

1019:   Notes: Add citation to (Clement, 1975) and definition of the interpolant
1020:   \nabla u_h(v_i) = \sum_{T_i \in support(v_i)} |T_i| \nabla u_h(T_i) / \sum_{T_i \in support(v_i)} |T_i| where |T_i| is the cell volume

1022:   Level: developer

1024: .seealso: DMProjectFunction(), DMComputeL2Diff(), DMPlexComputeL2FieldDiff(), DMComputeL2GradientDiff()
1025: @*/
1026: PetscErrorCode DMPlexComputeGradientClementInterpolant(DM dm, Vec locX, Vec locC)
1027: {
1028:   DM_Plex         *mesh  = (DM_Plex *) dm->data;
1029:   PetscInt         debug = mesh->printFEM;
1030:   DM               dmC;
1031:   PetscSection     section;
1032:   PetscQuadrature  quad;
1033:   PetscScalar     *interpolant, *gradsum;
1034:   PetscReal       *coords, *detJ, *J, *invJ;
1035:   const PetscReal *quadPoints, *quadWeights;
1036:   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, cEndInterior, vStart, vEnd, v, field, fieldOffset;
1037:   PetscErrorCode   ierr;

1040:   VecGetDM(locC, &dmC);
1041:   VecSet(locC, 0.0);
1042:   DMGetDimension(dm, &dim);
1043:   DMGetCoordinateDim(dm, &coordDim);
1044:   DMGetDefaultSection(dm, &section);
1045:   PetscSectionGetNumFields(section, &numFields);
1046:   for (field = 0; field < numFields; ++field) {
1047:     PetscObject  obj;
1048:     PetscClassId id;
1049:     PetscInt     Nc;

1051:     DMGetField(dm, field, &obj);
1052:     PetscObjectGetClassId(obj, &id);
1053:     if (id == PETSCFE_CLASSID) {
1054:       PetscFE fe = (PetscFE) obj;

1056:       PetscFEGetQuadrature(fe, &quad);
1057:       PetscFEGetNumComponents(fe, &Nc);
1058:     } else if (id == PETSCFV_CLASSID) {
1059:       PetscFV fv = (PetscFV) obj;

1061:       PetscFVGetQuadrature(fv, &quad);
1062:       PetscFVGetNumComponents(fv, &Nc);
1063:     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
1064:     numComponents += Nc;
1065:   }
1066:   PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
1067:   if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents);
1068:   PetscMalloc6(coordDim*numComponents*2,&gradsum,coordDim*numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J,coordDim*coordDim*Nq,&invJ);
1069:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1070:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
1071:   DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
1072:   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1073:   for (v = vStart; v < vEnd; ++v) {
1074:     PetscScalar volsum = 0.0;
1075:     PetscInt   *star = NULL;
1076:     PetscInt    starSize, st, d, fc;

1078:     PetscMemzero(gradsum, coordDim*numComponents * sizeof(PetscScalar));
1079:     DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);
1080:     for (st = 0; st < starSize*2; st += 2) {
1081:       const PetscInt cell = star[st];
1082:       PetscScalar   *grad = &gradsum[coordDim*numComponents];
1083:       PetscScalar   *x    = NULL;
1084:       PetscReal      vol  = 0.0;

1086:       if ((cell < cStart) || (cell >= cEnd)) continue;
1087:       DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, J, invJ, detJ);
1088:       DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x);
1089:       for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1090:         PetscObject  obj;
1091:         PetscClassId id;
1092:         PetscInt     Nb, Nc, q, qc = 0;

1094:         PetscMemzero(grad, coordDim*numComponents * sizeof(PetscScalar));
1095:         DMGetField(dm, field, &obj);
1096:         PetscObjectGetClassId(obj, &id);
1097:         if (id == PETSCFE_CLASSID)      {PetscFEGetNumComponents((PetscFE) obj, &Nc);PetscFEGetDimension((PetscFE) obj, &Nb);}
1098:         else if (id == PETSCFV_CLASSID) {PetscFVGetNumComponents((PetscFV) obj, &Nc);Nb = 1;}
1099:         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
1100:         for (q = 0; q < Nq; ++q) {
1101:           if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", (double)detJ[q], cell, q);
1102:           if (ierr) {
1103:             PetscErrorCode ierr2;
1104:             ierr2 = DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x);CHKERRQ(ierr2);
1105:             ierr2 = DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr2);
1106:             ierr2 = PetscFree4(interpolant,coords,detJ,J);CHKERRQ(ierr2);
1107: 
1108:           }
1109:           if (id == PETSCFE_CLASSID)      {PetscFEInterpolateGradient_Static((PetscFE) obj, &x[fieldOffset], coordDim, invJ, NULL, q, interpolant);}
1110:           else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
1111:           for (fc = 0; fc < Nc; ++fc) {
1112:             const PetscReal wt = quadWeights[q*qNc+qc+fc];

1114:             for (d = 0; d < coordDim; ++d) grad[fc*coordDim+d] += interpolant[fc*dim+d]*wt*detJ[q];
1115:           }
1116:           vol += quadWeights[q*qNc]*detJ[q];
1117:         }
1118:         fieldOffset += Nb;
1119:         qc          += Nc;
1120:       }
1121:       DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x);
1122:       for (fc = 0; fc < numComponents; ++fc) {
1123:         for (d = 0; d < coordDim; ++d) {
1124:           gradsum[fc*coordDim+d] += grad[fc*coordDim+d];
1125:         }
1126:       }
1127:       volsum += vol;
1128:       if (debug) {
1129:         PetscPrintf(PETSC_COMM_SELF, "Cell %D gradient: [", cell);
1130:         for (fc = 0; fc < numComponents; ++fc) {
1131:           for (d = 0; d < coordDim; ++d) {
1132:             if (fc || d > 0) {PetscPrintf(PETSC_COMM_SELF, ", ");}
1133:             PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc*coordDim+d]));
1134:           }
1135:         }
1136:         PetscPrintf(PETSC_COMM_SELF, "]\n");
1137:       }
1138:     }
1139:     for (fc = 0; fc < numComponents; ++fc) {
1140:       for (d = 0; d < coordDim; ++d) gradsum[fc*coordDim+d] /= volsum;
1141:     }
1142:     DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);
1143:     DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES);
1144:   }
1145:   PetscFree6(gradsum,interpolant,coords,detJ,J,invJ);
1146:   return(0);
1147: }

1149: /*@
1150:   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user

1152:   Input Parameters:
1153: + dm - The mesh
1154: . X  - Local input vector
1155: - user - The user context

1157:   Output Parameter:
1158: . integral - Local integral for each field

1160:   Level: developer

1162: .seealso: DMPlexComputeResidualFEM()
1163: @*/
1164: PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
1165: {
1166:   DM_Plex           *mesh  = (DM_Plex *) dm->data;
1167:   DM                 dmAux, dmGrad;
1168:   Vec                localX, A, cellGeometryFVM = NULL, faceGeometryFVM = NULL, locGrad = NULL;
1169:   PetscDS            prob, probAux = NULL;
1170:   PetscSection       section, sectionAux;
1171:   PetscFV            fvm = NULL;
1172:   PetscFECellGeom   *cgeomFEM;
1173:   PetscFVFaceGeom   *fgeomFVM;
1174:   PetscFVCellGeom   *cgeomFVM;
1175:   PetscScalar       *u, *a = NULL;
1176:   const PetscScalar *constants, *lgrad;
1177:   PetscReal         *lintegral;
1178:   PetscInt          *uOff, *uOff_x, *aOff = NULL;
1179:   PetscInt           dim, numConstants, Nf, NfAux = 0, f, numCells, cStart, cEnd, cEndInterior, c;
1180:   PetscInt           totDim, totDimAux;
1181:   PetscBool          useFVM = PETSC_FALSE;
1182:   PetscErrorCode     ierr;

1185:   PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);
1186:   DMGetLocalVector(dm, &localX);
1187:   DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localX, 0.0, NULL, NULL, NULL);
1188:   DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);
1189:   DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);
1190:   DMGetDimension(dm, &dim);
1191:   DMGetDefaultSection(dm, &section);
1192:   DMGetDS(dm, &prob);
1193:   PetscDSGetTotalDimension(prob, &totDim);
1194:   PetscDSGetComponentOffsets(prob, &uOff);
1195:   PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);
1196:   PetscDSGetConstants(prob, &numConstants, &constants);
1197:   PetscSectionGetNumFields(section, &Nf);
1198:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
1199:   DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
1200:   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1201:   numCells = cEnd - cStart;
1202:   PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);
1203:   PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);
1204:   if (dmAux) {
1205:     DMGetDS(dmAux, &probAux);
1206:     PetscDSGetNumFields(probAux, &NfAux);
1207:     DMGetDefaultSection(dmAux, &sectionAux);
1208:     PetscDSGetTotalDimension(probAux, &totDimAux);
1209:     PetscDSGetComponentOffsets(probAux, &aOff);
1210:     PetscDSGetEvaluationArrays(probAux, &a, NULL, NULL);
1211:   }
1212:   for (f = 0; f < Nf; ++f) {
1213:     PetscObject  obj;
1214:     PetscClassId id;

1216:     PetscDSGetDiscretization(prob, f, &obj);
1217:     PetscObjectGetClassId(obj, &id);
1218:     if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;}
1219:   }
1220:   if (useFVM) {
1221:     Vec       grad;
1222:     PetscInt  fStart, fEnd;
1223:     PetscBool compGrad;

1225:     PetscFVGetComputeGradients(fvm, &compGrad);
1226:     PetscFVSetComputeGradients(fvm, PETSC_TRUE);
1227:     DMPlexComputeGeometryFVM(dm, &cellGeometryFVM, &faceGeometryFVM);
1228:     DMPlexComputeGradientFVM(dm, fvm, faceGeometryFVM, cellGeometryFVM, &dmGrad);
1229:     PetscFVSetComputeGradients(fvm, compGrad);
1230:     VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);
1231:     VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);
1232:     /* Reconstruct and limit cell gradients */
1233:     DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
1234:     DMGetGlobalVector(dmGrad, &grad);
1235:     DMPlexReconstructGradients_Internal(dm, fvm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, localX, grad);
1236:     /* Communicate gradient values */
1237:     DMGetLocalVector(dmGrad, &locGrad);
1238:     DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);
1239:     DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);
1240:     DMRestoreGlobalVector(dmGrad, &grad);
1241:     /* Handle non-essential (e.g. outflow) boundary values */
1242:     DMPlexInsertBoundaryValues(dm, PETSC_FALSE, localX, 0.0, faceGeometryFVM, cellGeometryFVM, locGrad);
1243:     VecGetArrayRead(locGrad, &lgrad);
1244:   }
1245:   PetscMalloc3(Nf,&lintegral,numCells*totDim,&u,numCells,&cgeomFEM);
1246:   if (dmAux) {PetscMalloc1(numCells*totDimAux, &a);}
1247:   for (f = 0; f < Nf; ++f) {lintegral[f] = 0.0;}
1248:   for (c = cStart; c < cEnd; ++c) {
1249:     PetscScalar *x = NULL;
1250:     PetscInt     i;

1252:     DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeomFEM[c].v0, cgeomFEM[c].J, cgeomFEM[c].invJ, &cgeomFEM[c].detJ);
1253:     if (cgeomFEM[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeomFEM[c].detJ, c);
1254:     DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);
1255:     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
1256:     DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);
1257:     if (dmAux) {
1258:       DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);
1259:       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
1260:       DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);
1261:     }
1262:   }
1263:   for (f = 0; f < Nf; ++f) {
1264:     PetscObject  obj;
1265:     PetscClassId id;
1266:     PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;

1268:     PetscDSGetDiscretization(prob, f, &obj);
1269:     PetscObjectGetClassId(obj, &id);
1270:     if (id == PETSCFE_CLASSID) {
1271:       PetscFE         fe = (PetscFE) obj;
1272:       PetscQuadrature q;
1273:       PetscInt        Nq, Nb;

1275:       PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
1276:       PetscFEGetQuadrature(fe, &q);
1277:       PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL);
1278:       PetscFEGetDimension(fe, &Nb);
1279:       blockSize = Nb*Nq;
1280:       batchSize = numBlocks * blockSize;
1281:       PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
1282:       numChunks = numCells / (numBatches*batchSize);
1283:       Ne        = numChunks*numBatches*batchSize;
1284:       Nr        = numCells % (numBatches*batchSize);
1285:       offset    = numCells - Nr;
1286:       PetscFEIntegrate(fe, prob, f, Ne, cgeomFEM, u, probAux, a, lintegral);
1287:       PetscFEIntegrate(fe, prob, f, Nr, &cgeomFEM[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], lintegral);
1288:     } else if (id == PETSCFV_CLASSID) {
1289:       /* PetscFV  fv = (PetscFV) obj; */
1290:       PetscInt       foff;
1291:       PetscPointFunc obj_func;
1292:       PetscScalar    lint;

1294:       PetscDSGetObjective(prob, f, &obj_func);
1295:       PetscDSGetFieldOffset(prob, f, &foff);
1296:       if (obj_func) {
1297:         for (c = 0; c < numCells; ++c) {
1298:           PetscScalar *u_x;

1300:           DMPlexPointLocalRead(dmGrad, c, lgrad, &u_x);
1301:           obj_func(dim, Nf, NfAux, uOff, uOff_x, &u[totDim*c+foff], NULL, u_x, aOff, NULL, a, NULL, NULL, 0.0, cgeomFVM[c].centroid, numConstants, constants, &lint);
1302:           lintegral[f] += PetscRealPart(lint)*cgeomFVM[c].volume;
1303:         }
1304:       }
1305:     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
1306:   }
1307:   if (useFVM) {
1308:     VecRestoreArrayRead(locGrad, &lgrad);
1309:     VecRestoreArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);
1310:     VecRestoreArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);
1311:     DMRestoreLocalVector(dmGrad, &locGrad);
1312:     VecDestroy(&faceGeometryFVM);
1313:     VecDestroy(&cellGeometryFVM);
1314:     DMDestroy(&dmGrad);
1315:   }
1316:   if (dmAux) {PetscFree(a);}
1317:   if (mesh->printFEM) {
1318:     PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");
1319:     for (f = 0; f < Nf; ++f) {PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", lintegral[f]);}
1320:     PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");
1321:   }
1322:   DMRestoreLocalVector(dm, &localX);
1323:   MPIU_Allreduce(lintegral, integral, Nf, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));
1324:   PetscFree3(lintegral,u,cgeomFEM);
1325:   PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);
1326:   return(0);
1327: }

1329: /*@
1330:   DMPlexComputeInterpolatorNested - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.

1332:   Input Parameters:
1333: + dmf  - The fine mesh
1334: . dmc  - The coarse mesh
1335: - user - The user context

1337:   Output Parameter:
1338: . In  - The interpolation matrix

1340:   Level: developer

1342: .seealso: DMPlexComputeInterpolatorGeneral(), DMPlexComputeJacobianFEM()
1343: @*/
1344: PetscErrorCode DMPlexComputeInterpolatorNested(DM dmc, DM dmf, Mat In, void *user)
1345: {
1346:   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1347:   const char       *name  = "Interpolator";
1348:   PetscDS           prob;
1349:   PetscFE          *feRef;
1350:   PetscFV          *fvRef;
1351:   PetscSection      fsection, fglobalSection;
1352:   PetscSection      csection, cglobalSection;
1353:   PetscScalar      *elemMat;
1354:   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, cEndInterior, c;
1355:   PetscInt          cTotDim, rTotDim = 0;
1356:   PetscErrorCode    ierr;

1359:   PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);
1360:   DMGetDimension(dmf, &dim);
1361:   DMGetDefaultSection(dmf, &fsection);
1362:   DMGetDefaultGlobalSection(dmf, &fglobalSection);
1363:   DMGetDefaultSection(dmc, &csection);
1364:   DMGetDefaultGlobalSection(dmc, &cglobalSection);
1365:   PetscSectionGetNumFields(fsection, &Nf);
1366:   DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);
1367:   DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);
1368:   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1369:   DMGetDS(dmf, &prob);
1370:   PetscCalloc2(Nf,&feRef,Nf,&fvRef);
1371:   for (f = 0; f < Nf; ++f) {
1372:     PetscObject  obj;
1373:     PetscClassId id;
1374:     PetscInt     rNb = 0, Nc = 0;

1376:     PetscDSGetDiscretization(prob, f, &obj);
1377:     PetscObjectGetClassId(obj, &id);
1378:     if (id == PETSCFE_CLASSID) {
1379:       PetscFE fe = (PetscFE) obj;

1381:       PetscFERefine(fe, &feRef[f]);
1382:       PetscFEGetDimension(feRef[f], &rNb);
1383:       PetscFEGetNumComponents(fe, &Nc);
1384:     } else if (id == PETSCFV_CLASSID) {
1385:       PetscFV        fv = (PetscFV) obj;
1386:       PetscDualSpace Q;

1388:       PetscFVRefine(fv, &fvRef[f]);
1389:       PetscFVGetDualSpace(fvRef[f], &Q);
1390:       PetscDualSpaceGetDimension(Q, &rNb);
1391:       PetscFVGetNumComponents(fv, &Nc);
1392:     }
1393:     rTotDim += rNb;
1394:   }
1395:   PetscDSGetTotalDimension(prob, &cTotDim);
1396:   PetscMalloc1(rTotDim*cTotDim,&elemMat);
1397:   PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));
1398:   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1399:     PetscDualSpace   Qref;
1400:     PetscQuadrature  f;
1401:     const PetscReal *qpoints, *qweights;
1402:     PetscReal       *points;
1403:     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;

1405:     /* Compose points from all dual basis functionals */
1406:     if (feRef[fieldI]) {
1407:       PetscFEGetDualSpace(feRef[fieldI], &Qref);
1408:       PetscFEGetNumComponents(feRef[fieldI], &Nc);
1409:     } else {
1410:       PetscFVGetDualSpace(fvRef[fieldI], &Qref);
1411:       PetscFVGetNumComponents(fvRef[fieldI], &Nc);
1412:     }
1413:     PetscDualSpaceGetDimension(Qref, &fpdim);
1414:     for (i = 0; i < fpdim; ++i) {
1415:       PetscDualSpaceGetFunctional(Qref, i, &f);
1416:       PetscQuadratureGetData(f, NULL, NULL, &Np, NULL, NULL);
1417:       npoints += Np;
1418:     }
1419:     PetscMalloc1(npoints*dim,&points);
1420:     for (i = 0, k = 0; i < fpdim; ++i) {
1421:       PetscDualSpaceGetFunctional(Qref, i, &f);
1422:       PetscQuadratureGetData(f, NULL, NULL, &Np, &qpoints, NULL);
1423:       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1424:     }

1426:     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
1427:       PetscObject  obj;
1428:       PetscClassId id;
1429:       PetscReal   *B;
1430:       PetscInt     NcJ = 0, cpdim = 0, j, qNc;

1432:       PetscDSGetDiscretization(prob, fieldJ, &obj);
1433:       PetscObjectGetClassId(obj, &id);
1434:       if (id == PETSCFE_CLASSID) {
1435:         PetscFE fe = (PetscFE) obj;

1437:         /* Evaluate basis at points */
1438:         PetscFEGetNumComponents(fe, &NcJ);
1439:         PetscFEGetDimension(fe, &cpdim);
1440:         /* For now, fields only interpolate themselves */
1441:         if (fieldI == fieldJ) {
1442:           if (Nc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %D does not match coarse field %D", Nc, NcJ);
1443:           PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);
1444:           for (i = 0, k = 0; i < fpdim; ++i) {
1445:             PetscDualSpaceGetFunctional(Qref, i, &f);
1446:             PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights);
1447:             if (qNc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %D does not match coarse field %D", qNc, NcJ);
1448:             for (p = 0; p < Np; ++p, ++k) {
1449:               for (j = 0; j < cpdim; ++j) {
1450:                 /* NOTE: This is not quite right, unless fpdim == number of fine grid functional quad points */
1451:                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i)*cTotDim + offsetJ + j] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p*qNc+c];
1452:               }
1453:             }
1454:           }
1455:           PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);
1456:         }
1457:       } else if (id == PETSCFV_CLASSID) {
1458:         PetscFV        fv = (PetscFV) obj;

1460:         /* Evaluate constant function at points */
1461:         PetscFVGetNumComponents(fv, &NcJ);
1462:         cpdim = 1;
1463:         /* For now, fields only interpolate themselves */
1464:         if (fieldI == fieldJ) {
1465:           if (Nc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %d does not match coarse field %d", Nc, NcJ);
1466:           for (i = 0, k = 0; i < fpdim; ++i) {
1467:             PetscDualSpaceGetFunctional(Qref, i, &f);
1468:             PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights);
1469:             if (qNc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %D does not match coarse field %D", qNc, NcJ);
1470:             for (p = 0; p < Np; ++p, ++k) {
1471:               for (j = 0; j < cpdim; ++j) {
1472:                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += 1.0*qweights[p*qNc+c];
1473:               }
1474:             }
1475:           }
1476:         }
1477:       }
1478:       offsetJ += cpdim*NcJ;
1479:     }
1480:     offsetI += fpdim*Nc;
1481:     PetscFree(points);
1482:   }
1483:   if (mesh->printFEM > 1) {DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);}
1484:   /* Preallocate matrix */
1485:   {
1486:     Mat          preallocator;
1487:     PetscScalar *vals;
1488:     PetscInt    *cellCIndices, *cellFIndices;
1489:     PetscInt     locRows, locCols, cell;

1491:     MatGetLocalSize(In, &locRows, &locCols);
1492:     MatCreate(PetscObjectComm((PetscObject) In), &preallocator);
1493:     MatSetType(preallocator, MATPREALLOCATOR);
1494:     MatSetSizes(preallocator, locRows, locCols, PETSC_DETERMINE, PETSC_DETERMINE);
1495:     MatSetUp(preallocator);
1496:     PetscCalloc3(rTotDim*cTotDim, &vals,cTotDim,&cellCIndices,rTotDim,&cellFIndices);
1497:     for (cell = cStart; cell < cEnd; ++cell) {
1498:       DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);
1499:       MatSetValues(preallocator, rTotDim, cellFIndices, cTotDim, cellCIndices, vals, INSERT_VALUES);
1500:     }
1501:     PetscFree3(vals,cellCIndices,cellFIndices);
1502:     MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY);
1503:     MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY);
1504:     MatPreallocatorPreallocate(preallocator, PETSC_TRUE, In);
1505:     MatDestroy(&preallocator);
1506:   }
1507:   /* Fill matrix */
1508:   MatZeroEntries(In);
1509:   for (c = cStart; c < cEnd; ++c) {
1510:     DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);
1511:   }
1512:   for (f = 0; f < Nf; ++f) {PetscFEDestroy(&feRef[f]);}
1513:   PetscFree2(feRef,fvRef);
1514:   PetscFree(elemMat);
1515:   MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);
1516:   MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);
1517:   if (mesh->printFEM) {
1518:     PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);
1519:     MatChop(In, 1.0e-10);
1520:     MatView(In, PETSC_VIEWER_STDOUT_WORLD);
1521:   }
1522:   PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);
1523:   return(0);
1524: }

1526: PetscErrorCode DMPlexComputeMassMatrixNested(DM dmc, DM dmf, Mat mass, void *user)
1527: {
1528:   SETERRQ(PetscObjectComm((PetscObject) dmc), PETSC_ERR_SUP, "Laziness");
1529: }

1531: /*@
1532:   DMPlexComputeInterpolatorGeneral - Form the local portion of the interpolation matrix I from the coarse DM to a non-nested fine DM.

1534:   Input Parameters:
1535: + dmf  - The fine mesh
1536: . dmc  - The coarse mesh
1537: - user - The user context

1539:   Output Parameter:
1540: . In  - The interpolation matrix

1542:   Level: developer

1544: .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM()
1545: @*/
1546: PetscErrorCode DMPlexComputeInterpolatorGeneral(DM dmc, DM dmf, Mat In, void *user)
1547: {
1548:   DM_Plex       *mesh = (DM_Plex *) dmf->data;
1549:   const char    *name = "Interpolator";
1550:   PetscDS        prob;
1551:   PetscSection   fsection, csection, globalFSection, globalCSection;
1552:   PetscHashJK    ht;
1553:   PetscLayout    rLayout;
1554:   PetscInt      *dnz, *onz;
1555:   PetscInt       locRows, rStart, rEnd;
1556:   PetscReal     *x, *v0, *J, *invJ, detJ;
1557:   PetscReal     *v0c, *Jc, *invJc, detJc;
1558:   PetscScalar   *elemMat;
1559:   PetscInt       dim, Nf, field, totDim, cStart, cEnd, cell, ccell;

1563:   PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);
1564:   DMGetCoordinateDim(dmc, &dim);
1565:   DMGetDS(dmc, &prob);
1566:   PetscDSGetRefCoordArrays(prob, &x, NULL);
1567:   PetscDSGetNumFields(prob, &Nf);
1568:   PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);
1569:   PetscMalloc3(dim,&v0c,dim*dim,&Jc,dim*dim,&invJc);
1570:   DMGetDefaultSection(dmf, &fsection);
1571:   DMGetDefaultGlobalSection(dmf, &globalFSection);
1572:   DMGetDefaultSection(dmc, &csection);
1573:   DMGetDefaultGlobalSection(dmc, &globalCSection);
1574:   DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);
1575:   PetscDSGetTotalDimension(prob, &totDim);
1576:   PetscMalloc1(totDim, &elemMat);

1578:   MatGetLocalSize(In, &locRows, NULL);
1579:   PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);
1580:   PetscLayoutSetLocalSize(rLayout, locRows);
1581:   PetscLayoutSetBlockSize(rLayout, 1);
1582:   PetscLayoutSetUp(rLayout);
1583:   PetscLayoutGetRange(rLayout, &rStart, &rEnd);
1584:   PetscLayoutDestroy(&rLayout);
1585:   PetscCalloc2(locRows,&dnz,locRows,&onz);
1586:   PetscHashJKCreate(&ht);
1587:   for (field = 0; field < Nf; ++field) {
1588:     PetscObject      obj;
1589:     PetscClassId     id;
1590:     PetscDualSpace   Q = NULL;
1591:     PetscQuadrature  f;
1592:     const PetscReal *qpoints;
1593:     PetscInt         Nc, Np, fpdim, i, d;

1595:     PetscDSGetDiscretization(prob, field, &obj);
1596:     PetscObjectGetClassId(obj, &id);
1597:     if (id == PETSCFE_CLASSID) {
1598:       PetscFE fe = (PetscFE) obj;

1600:       PetscFEGetDualSpace(fe, &Q);
1601:       PetscFEGetNumComponents(fe, &Nc);
1602:     } else if (id == PETSCFV_CLASSID) {
1603:       PetscFV fv = (PetscFV) obj;

1605:       PetscFVGetDualSpace(fv, &Q);
1606:       Nc   = 1;
1607:     }
1608:     PetscDualSpaceGetDimension(Q, &fpdim);
1609:     /* For each fine grid cell */
1610:     for (cell = cStart; cell < cEnd; ++cell) {
1611:       PetscInt *findices,   *cindices;
1612:       PetscInt  numFIndices, numCIndices;

1614:       DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);
1615:       DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);
1616:       if (numFIndices != fpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim);
1617:       for (i = 0; i < fpdim; ++i) {
1618:         Vec             pointVec;
1619:         PetscScalar    *pV;
1620:         PetscSF         coarseCellSF = NULL;
1621:         const PetscSFNode *coarseCells;
1622:         PetscInt        numCoarseCells, q, c;

1624:         /* Get points from the dual basis functional quadrature */
1625:         PetscDualSpaceGetFunctional(Q, i, &f);
1626:         PetscQuadratureGetData(f, NULL, NULL, &Np, &qpoints, NULL);
1627:         VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);
1628:         VecSetBlockSize(pointVec, dim);
1629:         VecGetArray(pointVec, &pV);
1630:         for (q = 0; q < Np; ++q) {
1631:           /* Transform point to real space */
1632:           CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x);
1633:           for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
1634:         }
1635:         VecRestoreArray(pointVec, &pV);
1636:         /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
1637:         /* OPT: Pack all quad points from fine cell */
1638:         DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);
1639:         PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view");
1640:         /* Update preallocation info */
1641:         PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);
1642:         if (numCoarseCells != Np) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
1643:         {
1644:           PetscHashJKKey  key;
1645:           PetscHashJKIter missing, iter;

1647:           key.j = findices[i];
1648:           if (key.j >= 0) {
1649:             /* Get indices for coarse elements */
1650:             for (ccell = 0; ccell < numCoarseCells; ++ccell) {
1651:               DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);
1652:               for (c = 0; c < numCIndices; ++c) {
1653:                 key.k = cindices[c];
1654:                 if (key.k < 0) continue;
1655:                 PetscHashJKPut(ht, key, &missing, &iter);
1656:                 if (missing) {
1657:                   PetscHashJKSet(ht, iter, 1);
1658:                   if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
1659:                   else                                     ++onz[key.j-rStart];
1660:                 }
1661:               }
1662:               DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);
1663:             }
1664:           }
1665:         }
1666:         PetscSFDestroy(&coarseCellSF);
1667:         VecDestroy(&pointVec);
1668:       }
1669:       DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);
1670:     }
1671:   }
1672:   PetscHashJKDestroy(&ht);
1673:   MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);
1674:   MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
1675:   PetscFree2(dnz,onz);
1676:   for (field = 0; field < Nf; ++field) {
1677:     PetscObject      obj;
1678:     PetscClassId     id;
1679:     PetscDualSpace   Q = NULL;
1680:     PetscQuadrature  f;
1681:     const PetscReal *qpoints, *qweights;
1682:     PetscInt         Nc, qNc, Np, fpdim, i, d;

1684:     PetscDSGetDiscretization(prob, field, &obj);
1685:     PetscObjectGetClassId(obj, &id);
1686:     if (id == PETSCFE_CLASSID) {
1687:       PetscFE fe = (PetscFE) obj;

1689:       PetscFEGetDualSpace(fe, &Q);
1690:       PetscFEGetNumComponents(fe, &Nc);
1691:     } else if (id == PETSCFV_CLASSID) {
1692:       PetscFV fv = (PetscFV) obj;

1694:       PetscFVGetDualSpace(fv, &Q);
1695:       Nc   = 1;
1696:     }
1697:     PetscDualSpaceGetDimension(Q, &fpdim);
1698:     /* For each fine grid cell */
1699:     for (cell = cStart; cell < cEnd; ++cell) {
1700:       PetscInt *findices,   *cindices;
1701:       PetscInt  numFIndices, numCIndices;

1703:       DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);
1704:       DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);
1705:       if (numFIndices != fpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim);
1706:       for (i = 0; i < fpdim; ++i) {
1707:         Vec             pointVec;
1708:         PetscScalar    *pV;
1709:         PetscSF         coarseCellSF = NULL;
1710:         const PetscSFNode *coarseCells;
1711:         PetscInt        numCoarseCells, cpdim, q, c, j;

1713:         /* Get points from the dual basis functional quadrature */
1714:         PetscDualSpaceGetFunctional(Q, i, &f);
1715:         PetscQuadratureGetData(f, NULL, &qNc, &Np, &qpoints, &qweights);
1716:         if (qNc != Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %D does not match coarse field %D", qNc, Nc);
1717:         VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);
1718:         VecSetBlockSize(pointVec, dim);
1719:         VecGetArray(pointVec, &pV);
1720:         for (q = 0; q < Np; ++q) {
1721:           /* Transform point to real space */
1722:           CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x);
1723:           for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
1724:         }
1725:         VecRestoreArray(pointVec, &pV);
1726:         /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
1727:         /* OPT: Read this out from preallocation information */
1728:         DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);
1729:         /* Update preallocation info */
1730:         PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);
1731:         if (numCoarseCells != Np) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
1732:         VecGetArray(pointVec, &pV);
1733:         for (ccell = 0; ccell < numCoarseCells; ++ccell) {
1734:           PetscReal pVReal[3];

1736:           DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);
1737:           /* Transform points from real space to coarse reference space */
1738:           DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc);
1739:           for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell*dim+d]);
1740:           CoordinatesRealToRef(dim, dim, v0c, invJc, pVReal, x);

1742:           if (id == PETSCFE_CLASSID) {
1743:             PetscFE    fe = (PetscFE) obj;
1744:             PetscReal *B;

1746:             /* Evaluate coarse basis on contained point */
1747:             PetscFEGetDimension(fe, &cpdim);
1748:             PetscFEGetTabulation(fe, 1, x, &B, NULL, NULL);
1749:             PetscMemzero(elemMat, cpdim * sizeof(PetscScalar));
1750:             /* Get elemMat entries by multiplying by weight */
1751:             for (j = 0; j < cpdim; ++j) {
1752:               for (c = 0; c < Nc; ++c) elemMat[j] += B[j*Nc + c]*qweights[ccell*qNc + c];
1753:             }
1754:             PetscFERestoreTabulation(fe, 1, x, &B, NULL, NULL);
1755:           } else {
1756:             cpdim = 1;
1757:             for (j = 0; j < cpdim; ++j) {
1758:               for (c = 0; c < Nc; ++c) elemMat[j] += 1.0*qweights[ccell*qNc + c];
1759:             }
1760:           }
1761:           /* Update interpolator */
1762:           if (mesh->printFEM > 1) {DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);}
1763:           if (numCIndices != cpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D", numCIndices, cpdim);
1764:           MatSetValues(In, 1, &findices[i], numCIndices, cindices, elemMat, INSERT_VALUES);
1765:           DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);
1766:         }
1767:         VecRestoreArray(pointVec, &pV);
1768:         PetscSFDestroy(&coarseCellSF);
1769:         VecDestroy(&pointVec);
1770:       }
1771:       DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);
1772:     }
1773:   }
1774:   PetscFree3(v0,J,invJ);
1775:   PetscFree3(v0c,Jc,invJc);
1776:   PetscFree(elemMat);
1777:   MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);
1778:   MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);
1779:   PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);
1780:   return(0);
1781: }

1783: /*@
1784:   DMPlexComputeMassMatrixGeneral - Form the local portion of the mass matrix M from the coarse DM to a non-nested fine DM.

1786:   Input Parameters:
1787: + dmf  - The fine mesh
1788: . dmc  - The coarse mesh
1789: - user - The user context

1791:   Output Parameter:
1792: . mass  - The mass matrix

1794:   Level: developer

1796: .seealso: DMPlexComputeMassMatrixNested(), DMPlexComputeInterpolatorNested(), DMPlexComputeInterpolatorGeneral(), DMPlexComputeJacobianFEM()
1797: @*/
1798: PetscErrorCode DMPlexComputeMassMatrixGeneral(DM dmc, DM dmf, Mat mass, void *user)
1799: {
1800:   DM_Plex       *mesh = (DM_Plex *) dmf->data;
1801:   const char    *name = "Mass Matrix";
1802:   PetscDS        prob;
1803:   PetscSection   fsection, csection, globalFSection, globalCSection;
1804:   PetscHashJK    ht;
1805:   PetscLayout    rLayout;
1806:   PetscInt      *dnz, *onz;
1807:   PetscInt       locRows, rStart, rEnd;
1808:   PetscReal     *x, *v0, *J, *invJ, detJ;
1809:   PetscReal     *v0c, *Jc, *invJc, detJc;
1810:   PetscScalar   *elemMat;
1811:   PetscInt       dim, Nf, field, totDim, cStart, cEnd, cell, ccell;

1815:   DMGetCoordinateDim(dmc, &dim);
1816:   DMGetDS(dmc, &prob);
1817:   PetscDSGetRefCoordArrays(prob, &x, NULL);
1818:   PetscDSGetNumFields(prob, &Nf);
1819:   PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);
1820:   PetscMalloc3(dim,&v0c,dim*dim,&Jc,dim*dim,&invJc);
1821:   DMGetDefaultSection(dmf, &fsection);
1822:   DMGetDefaultGlobalSection(dmf, &globalFSection);
1823:   DMGetDefaultSection(dmc, &csection);
1824:   DMGetDefaultGlobalSection(dmc, &globalCSection);
1825:   DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);
1826:   PetscDSGetTotalDimension(prob, &totDim);
1827:   PetscMalloc1(totDim, &elemMat);

1829:   MatGetLocalSize(mass, &locRows, NULL);
1830:   PetscLayoutCreate(PetscObjectComm((PetscObject) mass), &rLayout);
1831:   PetscLayoutSetLocalSize(rLayout, locRows);
1832:   PetscLayoutSetBlockSize(rLayout, 1);
1833:   PetscLayoutSetUp(rLayout);
1834:   PetscLayoutGetRange(rLayout, &rStart, &rEnd);
1835:   PetscLayoutDestroy(&rLayout);
1836:   PetscCalloc2(locRows,&dnz,locRows,&onz);
1837:   PetscHashJKCreate(&ht);
1838:   for (field = 0; field < Nf; ++field) {
1839:     PetscObject      obj;
1840:     PetscClassId     id;
1841:     PetscQuadrature  quad;
1842:     const PetscReal *qpoints;
1843:     PetscInt         Nq, Nc, i, d;

1845:     PetscDSGetDiscretization(prob, field, &obj);
1846:     PetscObjectGetClassId(obj, &id);
1847:     if (id == PETSCFE_CLASSID) {PetscFEGetQuadrature((PetscFE) obj, &quad);}
1848:     else                       {PetscFVGetQuadrature((PetscFV) obj, &quad);}
1849:     PetscQuadratureGetData(quad, NULL, &Nc, &Nq, &qpoints, NULL);
1850:     /* For each fine grid cell */
1851:     for (cell = cStart; cell < cEnd; ++cell) {
1852:       Vec                pointVec;
1853:       PetscScalar       *pV;
1854:       PetscSF            coarseCellSF = NULL;
1855:       const PetscSFNode *coarseCells;
1856:       PetscInt           numCoarseCells, q, c;
1857:       PetscInt          *findices,   *cindices;
1858:       PetscInt           numFIndices, numCIndices;

1860:       DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);
1861:       DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);
1862:       /* Get points from the quadrature */
1863:       VecCreateSeq(PETSC_COMM_SELF, Nq*dim, &pointVec);
1864:       VecSetBlockSize(pointVec, dim);
1865:       VecGetArray(pointVec, &pV);
1866:       for (q = 0; q < Nq; ++q) {
1867:         /* Transform point to real space */
1868:         CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x);
1869:         for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
1870:       }
1871:       VecRestoreArray(pointVec, &pV);
1872:       /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
1873:       DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);
1874:       PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view");
1875:       /* Update preallocation info */
1876:       PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);
1877:       if (numCoarseCells != Nq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
1878:       {
1879:         PetscHashJKKey  key;
1880:         PetscHashJKIter missing, iter;

1882:         for (i = 0; i < numFIndices; ++i) {
1883:           key.j = findices[i];
1884:           if (key.j >= 0) {
1885:             /* Get indices for coarse elements */
1886:             for (ccell = 0; ccell < numCoarseCells; ++ccell) {
1887:               DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);
1888:               for (c = 0; c < numCIndices; ++c) {
1889:                 key.k = cindices[c];
1890:                 if (key.k < 0) continue;
1891:                 PetscHashJKPut(ht, key, &missing, &iter);
1892:                 if (missing) {
1893:                   PetscHashJKSet(ht, iter, 1);
1894:                   if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
1895:                   else                                     ++onz[key.j-rStart];
1896:                 }
1897:               }
1898:               DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);
1899:             }
1900:           }
1901:         }
1902:       }
1903:       PetscSFDestroy(&coarseCellSF);
1904:       VecDestroy(&pointVec);
1905:       DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);
1906:     }
1907:   }
1908:   PetscHashJKDestroy(&ht);
1909:   MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL);
1910:   MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
1911:   PetscFree2(dnz,onz);
1912:   for (field = 0; field < Nf; ++field) {
1913:     PetscObject      obj;
1914:     PetscClassId     id;
1915:     PetscQuadrature  quad;
1916:     PetscReal       *Bfine;
1917:     const PetscReal *qpoints, *qweights;
1918:     PetscInt         Nq, Nc, i, d;

1920:     PetscDSGetDiscretization(prob, field, &obj);
1921:     PetscObjectGetClassId(obj, &id);
1922:     if (id == PETSCFE_CLASSID) {PetscFEGetQuadrature((PetscFE) obj, &quad);PetscFEGetDefaultTabulation((PetscFE) obj, &Bfine, NULL, NULL);}
1923:     else                       {PetscFVGetQuadrature((PetscFV) obj, &quad);}
1924:     PetscQuadratureGetData(quad, NULL, &Nc, &Nq, &qpoints, &qweights);
1925:     /* For each fine grid cell */
1926:     for (cell = cStart; cell < cEnd; ++cell) {
1927:       Vec                pointVec;
1928:       PetscScalar       *pV;
1929:       PetscSF            coarseCellSF = NULL;
1930:       const PetscSFNode *coarseCells;
1931:       PetscInt           numCoarseCells, cpdim, q, c, j;
1932:       PetscInt          *findices,   *cindices;
1933:       PetscInt           numFIndices, numCIndices;

1935:       DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);
1936:       DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);
1937:       /* Get points from the quadrature */
1938:       VecCreateSeq(PETSC_COMM_SELF, Nq*dim, &pointVec);
1939:       VecSetBlockSize(pointVec, dim);
1940:       VecGetArray(pointVec, &pV);
1941:       for (q = 0; q < Nq; ++q) {
1942:         /* Transform point to real space */
1943:         CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x);
1944:         for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
1945:       }
1946:       VecRestoreArray(pointVec, &pV);
1947:       /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
1948:       DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);
1949:       /* Update matrix */
1950:       PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);
1951:       if (numCoarseCells != Nq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
1952:       VecGetArray(pointVec, &pV);
1953:       for (ccell = 0; ccell < numCoarseCells; ++ccell) {
1954:         PetscReal pVReal[3];

1956:         DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);
1957:         /* Transform points from real space to coarse reference space */
1958:         DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc);
1959:         for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell*dim+d]);
1960:         CoordinatesRealToRef(dim, dim, v0c, invJc, pVReal, x);

1962:         if (id == PETSCFE_CLASSID) {
1963:           PetscFE    fe = (PetscFE) obj;
1964:           PetscReal *B;

1966:           /* Evaluate coarse basis on contained point */
1967:           PetscFEGetDimension(fe, &cpdim);
1968:           PetscFEGetTabulation(fe, 1, x, &B, NULL, NULL);
1969:           /* Get elemMat entries by multiplying by weight */
1970:           for (i = 0; i < numFIndices; ++i) {
1971:             PetscMemzero(elemMat, cpdim * sizeof(PetscScalar));
1972:             for (j = 0; j < cpdim; ++j) {
1973:               for (c = 0; c < Nc; ++c) elemMat[j] += B[j*Nc + c]*Bfine[(ccell*numFIndices + i)*Nc + c]*qweights[ccell*Nc + c]*detJ;
1974:             }
1975:             /* Update interpolator */
1976:             if (mesh->printFEM > 1) {DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);}
1977:             if (numCIndices != cpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D", numCIndices, cpdim);
1978:             MatSetValues(mass, 1, &findices[i], numCIndices, cindices, elemMat, ADD_VALUES);
1979:           }
1980:           PetscFERestoreTabulation(fe, 1, x, &B, NULL, NULL);
1981:         } else {
1982:           cpdim = 1;
1983:           for (i = 0; i < numFIndices; ++i) {
1984:             PetscMemzero(elemMat, cpdim * sizeof(PetscScalar));
1985:             for (j = 0; j < cpdim; ++j) {
1986:               for (c = 0; c < Nc; ++c) elemMat[j] += 1.0*1.0*qweights[ccell*Nc + c]*detJ;
1987:             }
1988:             /* Update interpolator */
1989:             if (mesh->printFEM > 1) {DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);}
1990:             PetscPrintf(PETSC_COMM_SELF, "Nq: %d %d Nf: %d %d Nc: %d %d\n", ccell, Nq, i, numFIndices, j, numCIndices);
1991:             if (numCIndices != cpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D", numCIndices, cpdim);
1992:             MatSetValues(mass, 1, &findices[i], numCIndices, cindices, elemMat, ADD_VALUES);
1993:           }
1994:         }
1995:         DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);
1996:       }
1997:       VecRestoreArray(pointVec, &pV);
1998:       PetscSFDestroy(&coarseCellSF);
1999:       VecDestroy(&pointVec);
2000:       DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);
2001:     }
2002:   }
2003:   PetscFree3(v0,J,invJ);
2004:   PetscFree3(v0c,Jc,invJc);
2005:   PetscFree(elemMat);
2006:   MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY);
2007:   MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY);
2008:   return(0);
2009: }

2011: /*@
2012:   DMPlexComputeInjectorFEM - Compute a mapping from coarse unknowns to fine unknowns

2014:   Input Parameters:
2015: + dmc  - The coarse mesh
2016: - dmf  - The fine mesh
2017: - user - The user context

2019:   Output Parameter:
2020: . sc   - The mapping

2022:   Level: developer

2024: .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM()
2025: @*/
2026: PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
2027: {
2028:   PetscDS        prob;
2029:   PetscFE       *feRef;
2030:   PetscFV       *fvRef;
2031:   Vec            fv, cv;
2032:   IS             fis, cis;
2033:   PetscSection   fsection, fglobalSection, csection, cglobalSection;
2034:   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
2035:   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, cEndInterior, c, dim, d, startC, endC, offsetC, offsetF, m;

2039:   PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);
2040:   DMGetDimension(dmf, &dim);
2041:   DMGetDefaultSection(dmf, &fsection);
2042:   DMGetDefaultGlobalSection(dmf, &fglobalSection);
2043:   DMGetDefaultSection(dmc, &csection);
2044:   DMGetDefaultGlobalSection(dmc, &cglobalSection);
2045:   PetscSectionGetNumFields(fsection, &Nf);
2046:   DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);
2047:   DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);
2048:   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
2049:   DMGetDS(dmc, &prob);
2050:   PetscCalloc2(Nf,&feRef,Nf,&fvRef);
2051:   for (f = 0; f < Nf; ++f) {
2052:     PetscObject  obj;
2053:     PetscClassId id;
2054:     PetscInt     fNb = 0, Nc = 0;

2056:     PetscDSGetDiscretization(prob, f, &obj);
2057:     PetscObjectGetClassId(obj, &id);
2058:     if (id == PETSCFE_CLASSID) {
2059:       PetscFE fe = (PetscFE) obj;

2061:       PetscFERefine(fe, &feRef[f]);
2062:       PetscFEGetDimension(feRef[f], &fNb);
2063:       PetscFEGetNumComponents(fe, &Nc);
2064:     } else if (id == PETSCFV_CLASSID) {
2065:       PetscFV        fv = (PetscFV) obj;
2066:       PetscDualSpace Q;

2068:       PetscFVRefine(fv, &fvRef[f]);
2069:       PetscFVGetDualSpace(fvRef[f], &Q);
2070:       PetscDualSpaceGetDimension(Q, &fNb);
2071:       PetscFVGetNumComponents(fv, &Nc);
2072:     }
2073:     fTotDim += fNb*Nc;
2074:   }
2075:   PetscDSGetTotalDimension(prob, &cTotDim);
2076:   PetscMalloc1(cTotDim,&cmap);
2077:   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
2078:     PetscFE        feC;
2079:     PetscFV        fvC;
2080:     PetscDualSpace QF, QC;
2081:     PetscInt       NcF, NcC, fpdim, cpdim;

2083:     if (feRef[field]) {
2084:       PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);
2085:       PetscFEGetNumComponents(feC, &NcC);
2086:       PetscFEGetNumComponents(feRef[field], &NcF);
2087:       PetscFEGetDualSpace(feRef[field], &QF);
2088:       PetscDualSpaceGetDimension(QF, &fpdim);
2089:       PetscFEGetDualSpace(feC, &QC);
2090:       PetscDualSpaceGetDimension(QC, &cpdim);
2091:     } else {
2092:       PetscDSGetDiscretization(prob, field, (PetscObject *) &fvC);
2093:       PetscFVGetNumComponents(fvC, &NcC);
2094:       PetscFVGetNumComponents(fvRef[field], &NcF);
2095:       PetscFVGetDualSpace(fvRef[field], &QF);
2096:       PetscDualSpaceGetDimension(QF, &fpdim);
2097:       PetscFVGetDualSpace(fvC, &QC);
2098:       PetscDualSpaceGetDimension(QC, &cpdim);
2099:     }
2100:     if (NcF != NcC) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %d does not match coarse field %d", NcF, NcC);
2101:     for (c = 0; c < cpdim; ++c) {
2102:       PetscQuadrature  cfunc;
2103:       const PetscReal *cqpoints;
2104:       PetscInt         NpC;
2105:       PetscBool        found = PETSC_FALSE;

2107:       PetscDualSpaceGetFunctional(QC, c, &cfunc);
2108:       PetscQuadratureGetData(cfunc, NULL, NULL, &NpC, &cqpoints, NULL);
2109:       if (NpC != 1 && feRef[field]) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
2110:       for (f = 0; f < fpdim; ++f) {
2111:         PetscQuadrature  ffunc;
2112:         const PetscReal *fqpoints;
2113:         PetscReal        sum = 0.0;
2114:         PetscInt         NpF, comp;

2116:         PetscDualSpaceGetFunctional(QF, f, &ffunc);
2117:         PetscQuadratureGetData(ffunc, NULL, NULL, &NpF, &fqpoints, NULL);
2118:         if (NpC != NpF) continue;
2119:         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
2120:         if (sum > 1.0e-9) continue;
2121:         for (comp = 0; comp < NcC; ++comp) {
2122:           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
2123:         }
2124:         found = PETSC_TRUE;
2125:         break;
2126:       }
2127:       if (!found) {
2128:         /* TODO We really want the average here, but some asshole put VecScatter in the interface */
2129:         if (fvRef[field]) {
2130:           PetscInt comp;
2131:           for (comp = 0; comp < NcC; ++comp) {
2132:             cmap[(offsetC+c)*NcC+comp] = (offsetF+0)*NcF+comp;
2133:           }
2134:         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate matching functional for injection");
2135:       }
2136:     }
2137:     offsetC += cpdim*NcC;
2138:     offsetF += fpdim*NcF;
2139:   }
2140:   for (f = 0; f < Nf; ++f) {PetscFEDestroy(&feRef[f]);PetscFVDestroy(&fvRef[f]);}
2141:   PetscFree2(feRef,fvRef);

2143:   DMGetGlobalVector(dmf, &fv);
2144:   DMGetGlobalVector(dmc, &cv);
2145:   VecGetOwnershipRange(cv, &startC, &endC);
2146:   PetscSectionGetConstrainedStorageSize(cglobalSection, &m);
2147:   PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);
2148:   PetscMalloc1(m,&cindices);
2149:   PetscMalloc1(m,&findices);
2150:   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
2151:   for (c = cStart; c < cEnd; ++c) {
2152:     DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);
2153:     for (d = 0; d < cTotDim; ++d) {
2154:       if ((cellCIndices[d] < startC) || (cellCIndices[d] >= endC)) continue;
2155:       if ((findices[cellCIndices[d]-startC] >= 0) && (findices[cellCIndices[d]-startC] != cellFIndices[cmap[d]])) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coarse dof %d maps to both %d and %d", cindices[cellCIndices[d]-startC], findices[cellCIndices[d]-startC], cellFIndices[cmap[d]]);
2156:       cindices[cellCIndices[d]-startC] = cellCIndices[d];
2157:       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
2158:     }
2159:   }
2160:   PetscFree(cmap);
2161:   PetscFree2(cellCIndices,cellFIndices);

2163:   ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);
2164:   ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);
2165:   VecScatterCreate(cv, cis, fv, fis, sc);
2166:   ISDestroy(&cis);
2167:   ISDestroy(&fis);
2168:   DMRestoreGlobalVector(dmf, &fv);
2169:   DMRestoreGlobalVector(dmc, &cv);
2170:   PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);
2171:   return(0);
2172: }