Actual source code: plexfem.c

  1: #include <petsc/private/dmpleximpl.h>
  2: #include <petscsf.h>

  4: #include <petscblaslapack.h>
  5: #include <petsc/private/hashsetij.h>
  6: #include <petsc/private/petscfeimpl.h>
  7: #include <petsc/private/petscfvimpl.h>

  9: PetscBool Clementcite = PETSC_FALSE;
 10: const char ClementCitation[] = "@article{clement1975approximation,\n"
 11:                                "  title   = {Approximation by finite element functions using local regularization},\n"
 12:                                "  author  = {Philippe Cl{\\'e}ment},\n"
 13:                                "  journal = {Revue fran{\\c{c}}aise d'automatique, informatique, recherche op{\\'e}rationnelle. Analyse num{\\'e}rique},\n"
 14:                                "  volume  = {9},\n"
 15:                                "  number  = {R2},\n"
 16:                                "  pages   = {77--84},\n"
 17:                                "  year    = {1975}\n}\n";

 19: static PetscErrorCode DMPlexConvertPlex(DM dm, DM *plex, PetscBool copy)
 20: {
 21:   PetscBool      isPlex;

 23:   PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);
 24:   if (isPlex) {
 25:     *plex = dm;
 26:     PetscObjectReference((PetscObject) dm);
 27:   } else {
 28:     PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);
 29:     if (!*plex) {
 30:       DMConvert(dm, DMPLEX, plex);
 31:       PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);
 32:       if (copy) {
 33:         DMSubDomainHookLink link;

 35:         DMCopyAuxiliaryVec(dm, *plex);
 36:         /* Run the subdomain hook (this will copy the DMSNES/DMTS) */
 37:         for (link = dm->subdomainhook; link; link = link->next) {
 38:           if (link->ddhook) (*link->ddhook)(dm, *plex, link->ctx);
 39:         }
 40:       }
 41:     } else {
 42:       PetscObjectReference((PetscObject) *plex);
 43:     }
 44:   }
 45:   return 0;
 46: }

 48: static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom (void *ctx)
 49: {
 50:   PetscFEGeom *geom = (PetscFEGeom *) ctx;

 52:   PetscFEGeomDestroy(&geom);
 53:   return 0;
 54: }

 56: static PetscErrorCode DMPlexGetFEGeom(DMField coordField, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
 57: {
 58:   char            composeStr[33] = {0};
 59:   PetscObjectId   id;
 60:   PetscContainer  container;

 62:   PetscObjectGetId((PetscObject)quad,&id);
 63:   PetscSNPrintf(composeStr, 32, "DMPlexGetFEGeom_%x\n", id);
 64:   PetscObjectQuery((PetscObject) pointIS, composeStr, (PetscObject *) &container);
 65:   if (container) {
 66:     PetscContainerGetPointer(container, (void **) geom);
 67:   } else {
 68:     DMFieldCreateFEGeom(coordField, pointIS, quad, faceData, geom);
 69:     PetscContainerCreate(PETSC_COMM_SELF,&container);
 70:     PetscContainerSetPointer(container, (void *) *geom);
 71:     PetscContainerSetUserDestroy(container, PetscContainerUserDestroy_PetscFEGeom);
 72:     PetscObjectCompose((PetscObject) pointIS, composeStr, (PetscObject) container);
 73:     PetscContainerDestroy(&container);
 74:   }
 75:   return 0;
 76: }

 78: static PetscErrorCode DMPlexRestoreFEGeom(DMField coordField, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
 79: {
 80:   *geom = NULL;
 81:   return 0;
 82: }

 84: /*@
 85:   DMPlexGetScale - Get the scale for the specified fundamental unit

 87:   Not collective

 89:   Input Parameters:
 90: + dm   - the DM
 91: - unit - The SI unit

 93:   Output Parameter:
 94: . scale - The value used to scale all quantities with this unit

 96:   Level: advanced

 98: .seealso: DMPlexSetScale(), PetscUnit
 99: @*/
100: PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
101: {
102:   DM_Plex *mesh = (DM_Plex*) dm->data;

106:   *scale = mesh->scale[unit];
107:   return 0;
108: }

110: /*@
111:   DMPlexSetScale - Set the scale for the specified fundamental unit

113:   Not collective

115:   Input Parameters:
116: + dm   - the DM
117: . unit - The SI unit
118: - scale - The value used to scale all quantities with this unit

120:   Level: advanced

122: .seealso: DMPlexGetScale(), PetscUnit
123: @*/
124: PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
125: {
126:   DM_Plex *mesh = (DM_Plex*) dm->data;

129:   mesh->scale[unit] = scale;
130:   return 0;
131: }

133: static PetscErrorCode DMPlexProjectRigidBody_Private(PetscInt dim, PetscReal t, const PetscReal X[], PetscInt Nc, PetscScalar *mode, void *ctx)
134: {
135:   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}}};
136:   PetscInt *ctxInt  = (PetscInt *) ctx;
137:   PetscInt  dim2    = ctxInt[0];
138:   PetscInt  d       = ctxInt[1];
139:   PetscInt  i, j, k = dim > 2 ? d - dim : d;

142:   for (i = 0; i < dim; i++) mode[i] = 0.;
143:   if (d < dim) {
144:     mode[d] = 1.; /* Translation along axis d */
145:   } else {
146:     for (i = 0; i < dim; i++) {
147:       for (j = 0; j < dim; j++) {
148:         mode[j] += eps[i][j][k]*X[i]; /* Rotation about axis d */
149:       }
150:     }
151:   }
152:   return 0;
153: }

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

158:   Collective on dm

160:   Input Parameters:
161: + dm - the DM
162: - field - The field number for the rigid body space, or 0 for the default

164:   Output Parameter:
165: . sp - the null space

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

169:   Level: advanced

171: .seealso: MatNullSpaceCreate(), PCGAMG
172: @*/
173: PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscInt field, MatNullSpace *sp)
174: {
175:   PetscErrorCode (**func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *);
176:   MPI_Comm          comm;
177:   Vec               mode[6];
178:   PetscSection      section, globalSection;
179:   PetscInt          dim, dimEmbed, Nf, n, m, mmin, d, i, j;

181:   PetscObjectGetComm((PetscObject) dm, &comm);
182:   DMGetDimension(dm, &dim);
183:   DMGetCoordinateDim(dm, &dimEmbed);
184:   DMGetNumFields(dm, &Nf);
186:   if (dim == 1 && Nf < 2) {
187:     MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);
188:     return 0;
189:   }
190:   DMGetLocalSection(dm, &section);
191:   DMGetGlobalSection(dm, &globalSection);
192:   PetscSectionGetConstrainedStorageSize(globalSection, &n);
193:   PetscCalloc1(Nf, &func);
194:   m    = (dim*(dim+1))/2;
195:   VecCreate(comm, &mode[0]);
196:   VecSetType(mode[0], dm->vectype);
197:   VecSetSizes(mode[0], n, PETSC_DETERMINE);
198:   VecSetUp(mode[0]);
199:   VecGetSize(mode[0], &n);
200:   mmin = PetscMin(m, n);
201:   func[field] = DMPlexProjectRigidBody_Private;
202:   for (i = 1; i < m; ++i) VecDuplicate(mode[0], &mode[i]);
203:   for (d = 0; d < m; d++) {
204:     PetscInt ctx[2];
205:     void    *voidctx = (void *) (&ctx[0]);

207:     ctx[0] = dimEmbed;
208:     ctx[1] = d;
209:     DMProjectFunction(dm, 0.0, func, &voidctx, INSERT_VALUES, mode[d]);
210:   }
211:   /* Orthonormalize system */
212:   for (i = 0; i < mmin; ++i) {
213:     PetscScalar dots[6];

215:     VecNormalize(mode[i], NULL);
216:     VecMDot(mode[i], mmin-i-1, mode+i+1, dots+i+1);
217:     for (j = i+1; j < mmin; ++j) {
218:       dots[j] *= -1.0;
219:       VecAXPY(mode[j], dots[j], mode[i]);
220:     }
221:   }
222:   MatNullSpaceCreate(comm, PETSC_FALSE, mmin, mode, sp);
223:   for (i = 0; i < m; ++i) VecDestroy(&mode[i]);
224:   PetscFree(func);
225:   return 0;
226: }

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

231:   Collective on dm

233:   Input Parameters:
234: + dm    - the DM
235: . nb    - The number of bodies
236: . label - The DMLabel marking each domain
237: . nids  - The number of ids per body
238: - ids   - An array of the label ids in sequence for each domain

240:   Output Parameter:
241: . sp - the null space

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

245:   Level: advanced

247: .seealso: MatNullSpaceCreate()
248: @*/
249: PetscErrorCode DMPlexCreateRigidBodies(DM dm, PetscInt nb, DMLabel label, const PetscInt nids[], const PetscInt ids[], MatNullSpace *sp)
250: {
251:   MPI_Comm       comm;
252:   PetscSection   section, globalSection;
253:   Vec           *mode;
254:   PetscScalar   *dots;
255:   PetscInt       dim, dimEmbed, n, m, b, d, i, j, off;

257:   PetscObjectGetComm((PetscObject)dm,&comm);
258:   DMGetDimension(dm, &dim);
259:   DMGetCoordinateDim(dm, &dimEmbed);
260:   DMGetLocalSection(dm, &section);
261:   DMGetGlobalSection(dm, &globalSection);
262:   PetscSectionGetConstrainedStorageSize(globalSection, &n);
263:   m    = nb * (dim*(dim+1))/2;
264:   PetscMalloc2(m, &mode, m, &dots);
265:   VecCreate(comm, &mode[0]);
266:   VecSetSizes(mode[0], n, PETSC_DETERMINE);
267:   VecSetUp(mode[0]);
268:   for (i = 1; i < m; ++i) VecDuplicate(mode[0], &mode[i]);
269:   for (b = 0, off = 0; b < nb; ++b) {
270:     for (d = 0; d < m/nb; ++d) {
271:       PetscInt         ctx[2];
272:       PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *) = DMPlexProjectRigidBody_Private;
273:       void            *voidctx = (void *) (&ctx[0]);

275:       ctx[0] = dimEmbed;
276:       ctx[1] = d;
277:       DMProjectFunctionLabel(dm, 0.0, label, nids[b], &ids[off], 0, NULL, &func, &voidctx, INSERT_VALUES, mode[d]);
278:       off   += nids[b];
279:     }
280:   }
281:   /* Orthonormalize system */
282:   for (i = 0; i < m; ++i) {
283:     PetscScalar dots[6];

285:     VecNormalize(mode[i], NULL);
286:     VecMDot(mode[i], m-i-1, mode+i+1, dots+i+1);
287:     for (j = i+1; j < m; ++j) {
288:       dots[j] *= -1.0;
289:       VecAXPY(mode[j], dots[j], mode[i]);
290:     }
291:   }
292:   MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);
293:   for (i = 0; i< m; ++i) VecDestroy(&mode[i]);
294:   PetscFree2(mode, dots);
295:   return 0;
296: }

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

307:   Input Parameters:
308: + dm - the DMPlex object
309: - height - the maximum projection height >= 0

311:   Level: advanced

313: .seealso: DMPlexGetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal()
314: @*/
315: PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height)
316: {
317:   DM_Plex *plex = (DM_Plex *) dm->data;

320:   plex->maxProjectionHeight = height;
321:   return 0;
322: }

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

328:   Input Parameters:
329: . dm - the DMPlex object

331:   Output Parameters:
332: . height - the maximum projection height

334:   Level: intermediate

336: .seealso: DMPlexSetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal()
337: @*/
338: PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height)
339: {
340:   DM_Plex *plex = (DM_Plex *) dm->data;

343:   *height = plex->maxProjectionHeight;
344:   return 0;
345: }

347: typedef struct {
348:   PetscReal    alpha; /* The first Euler angle, and in 2D the only one */
349:   PetscReal    beta;  /* The second Euler angle */
350:   PetscReal    gamma; /* The third Euler angle */
351:   PetscInt     dim;   /* The dimension of R */
352:   PetscScalar *R;     /* The rotation matrix, transforming a vector in the local basis to the global basis */
353:   PetscScalar *RT;    /* The transposed rotation matrix, transforming a vector in the global basis to the local basis */
354: } RotCtx;

356: /*
357:   Note: Following https://en.wikipedia.org/wiki/Euler_angles, we will specify Euler angles by extrinsic rotations, meaning that
358:   we rotate with respect to a fixed initial coordinate system, the local basis (x-y-z). The global basis (X-Y-Z) is reached as follows:
359:   $ The XYZ system rotates about the z axis by alpha. The X axis is now at angle alpha with respect to the x axis.
360:   $ The XYZ system rotates again about the x axis by beta. The Z axis is now at angle beta with respect to the z axis.
361:   $ The XYZ system rotates a third time about the z axis by gamma.
362: */
363: static PetscErrorCode DMPlexBasisTransformSetUp_Rotation_Internal(DM dm, void *ctx)
364: {
365:   RotCtx        *rc  = (RotCtx *) ctx;
366:   PetscInt       dim = rc->dim;
367:   PetscReal      c1, s1, c2, s2, c3, s3;

369:   PetscMalloc2(PetscSqr(dim), &rc->R, PetscSqr(dim), &rc->RT);
370:   switch (dim) {
371:   case 2:
372:     c1 = PetscCosReal(rc->alpha);s1 = PetscSinReal(rc->alpha);
373:     rc->R[0] =  c1;rc->R[1] = s1;
374:     rc->R[2] = -s1;rc->R[3] = c1;
375:     PetscArraycpy(rc->RT, rc->R, PetscSqr(dim));
376:     DMPlex_Transpose2D_Internal(rc->RT);
377:     break;
378:   case 3:
379:     c1 = PetscCosReal(rc->alpha);s1 = PetscSinReal(rc->alpha);
380:     c2 = PetscCosReal(rc->beta); s2 = PetscSinReal(rc->beta);
381:     c3 = PetscCosReal(rc->gamma);s3 = PetscSinReal(rc->gamma);
382:     rc->R[0] =  c1*c3 - c2*s1*s3;rc->R[1] =  c3*s1    + c1*c2*s3;rc->R[2] = s2*s3;
383:     rc->R[3] = -c1*s3 - c2*c3*s1;rc->R[4] =  c1*c2*c3 - s1*s3;   rc->R[5] = c3*s2;
384:     rc->R[6] =  s1*s2;           rc->R[7] = -c1*s2;              rc->R[8] = c2;
385:     PetscArraycpy(rc->RT, rc->R, PetscSqr(dim));
386:     DMPlex_Transpose3D_Internal(rc->RT);
387:     break;
388:   default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not supported", dim);
389:   }
390:   return 0;
391: }

393: static PetscErrorCode DMPlexBasisTransformDestroy_Rotation_Internal(DM dm, void *ctx)
394: {
395:   RotCtx        *rc = (RotCtx *) ctx;

397:   PetscFree2(rc->R, rc->RT);
398:   PetscFree(rc);
399:   return 0;
400: }

402: static PetscErrorCode DMPlexBasisTransformGetMatrix_Rotation_Internal(DM dm, const PetscReal x[], PetscBool l2g, const PetscScalar **A, void *ctx)
403: {
404:   RotCtx *rc = (RotCtx *) ctx;

408:   if (l2g) {*A = rc->R;}
409:   else     {*A = rc->RT;}
410:   return 0;
411: }

413: PetscErrorCode DMPlexBasisTransformApplyReal_Internal(DM dm, const PetscReal x[], PetscBool l2g, PetscInt dim, const PetscReal *y, PetscReal *z, void *ctx)
414: {
415:   #if defined(PETSC_USE_COMPLEX)
416:   switch (dim) {
417:     case 2:
418:     {
419:       PetscScalar yt[2] = {y[0], y[1]}, zt[2] = {0.0,0.0};

421:       DMPlexBasisTransformApply_Internal(dm, x, l2g, dim, yt, zt, ctx);
422:       z[0] = PetscRealPart(zt[0]); z[1] = PetscRealPart(zt[1]);
423:     }
424:     break;
425:     case 3:
426:     {
427:       PetscScalar yt[3] = {y[0], y[1], y[2]}, zt[3] = {0.0,0.0,0.0};

429:       DMPlexBasisTransformApply_Internal(dm, x, l2g, dim, yt, zt, ctx);
430:       z[0] = PetscRealPart(zt[0]); z[1] = PetscRealPart(zt[1]); z[2] = PetscRealPart(zt[2]);
431:     }
432:     break;
433:   }
434:   #else
435:   DMPlexBasisTransformApply_Internal(dm, x, l2g, dim, y, z, ctx);
436:   #endif
437:   return 0;
438: }

440: PetscErrorCode DMPlexBasisTransformApply_Internal(DM dm, const PetscReal x[], PetscBool l2g, PetscInt dim, const PetscScalar *y, PetscScalar *z, void *ctx)
441: {
442:   const PetscScalar *A;

445:   (*dm->transformGetMatrix)(dm, x, l2g, &A, ctx);
446:   switch (dim) {
447:   case 2: DMPlex_Mult2D_Internal(A, 1, y, z);break;
448:   case 3: DMPlex_Mult3D_Internal(A, 1, y, z);break;
449:   }
450:   return 0;
451: }

453: static PetscErrorCode DMPlexBasisTransformField_Internal(DM dm, DM tdm, Vec tv, PetscInt p, PetscInt f, PetscBool l2g, PetscScalar *a)
454: {
455:   PetscSection       ts;
456:   const PetscScalar *ta, *tva;
457:   PetscInt           dof;

460:   DMGetLocalSection(tdm, &ts);
461:   PetscSectionGetFieldDof(ts, p, f, &dof);
462:   VecGetArrayRead(tv, &ta);
463:   DMPlexPointLocalFieldRead(tdm, p, f, ta, &tva);
464:   if (l2g) {
465:     switch (dof) {
466:     case 4: DMPlex_Mult2D_Internal(tva, 1, a, a);break;
467:     case 9: DMPlex_Mult3D_Internal(tva, 1, a, a);break;
468:     }
469:   } else {
470:     switch (dof) {
471:     case 4: DMPlex_MultTranspose2D_Internal(tva, 1, a, a);break;
472:     case 9: DMPlex_MultTranspose3D_Internal(tva, 1, a, a);break;
473:     }
474:   }
475:   VecRestoreArrayRead(tv, &ta);
476:   return 0;
477: }

479: static PetscErrorCode DMPlexBasisTransformFieldTensor_Internal(DM dm, DM tdm, Vec tv, PetscInt pf, PetscInt f, PetscInt pg, PetscInt g, PetscBool l2g, PetscInt lda, PetscScalar *a)
480: {
481:   PetscSection       s, ts;
482:   const PetscScalar *ta, *tvaf, *tvag;
483:   PetscInt           fdof, gdof, fpdof, gpdof;

486:   DMGetLocalSection(dm, &s);
487:   DMGetLocalSection(tdm, &ts);
488:   PetscSectionGetFieldDof(s, pf, f, &fpdof);
489:   PetscSectionGetFieldDof(s, pg, g, &gpdof);
490:   PetscSectionGetFieldDof(ts, pf, f, &fdof);
491:   PetscSectionGetFieldDof(ts, pg, g, &gdof);
492:   VecGetArrayRead(tv, &ta);
493:   DMPlexPointLocalFieldRead(tdm, pf, f, ta, &tvaf);
494:   DMPlexPointLocalFieldRead(tdm, pg, g, ta, &tvag);
495:   if (l2g) {
496:     switch (fdof) {
497:     case 4: DMPlex_MatMult2D_Internal(tvaf, gpdof, lda, a, a);break;
498:     case 9: DMPlex_MatMult3D_Internal(tvaf, gpdof, lda, a, a);break;
499:     }
500:     switch (gdof) {
501:     case 4: DMPlex_MatMultTransposeLeft2D_Internal(tvag, fpdof, lda, a, a);break;
502:     case 9: DMPlex_MatMultTransposeLeft3D_Internal(tvag, fpdof, lda, a, a);break;
503:     }
504:   } else {
505:     switch (fdof) {
506:     case 4: DMPlex_MatMultTranspose2D_Internal(tvaf, gpdof, lda, a, a);break;
507:     case 9: DMPlex_MatMultTranspose3D_Internal(tvaf, gpdof, lda, a, a);break;
508:     }
509:     switch (gdof) {
510:     case 4: DMPlex_MatMultLeft2D_Internal(tvag, fpdof, lda, a, a);break;
511:     case 9: DMPlex_MatMultLeft3D_Internal(tvag, fpdof, lda, a, a);break;
512:     }
513:   }
514:   VecRestoreArrayRead(tv, &ta);
515:   return 0;
516: }

518: PetscErrorCode DMPlexBasisTransformPoint_Internal(DM dm, DM tdm, Vec tv, PetscInt p, PetscBool fieldActive[], PetscBool l2g, PetscScalar *a)
519: {
520:   PetscSection    s;
521:   PetscSection    clSection;
522:   IS              clPoints;
523:   const PetscInt *clp;
524:   PetscInt       *points = NULL;
525:   PetscInt        Nf, f, Np, cp, dof, d = 0;

527:   DMGetLocalSection(dm, &s);
528:   PetscSectionGetNumFields(s, &Nf);
529:   DMPlexGetCompressedClosure(dm, s, p, &Np, &points, &clSection, &clPoints, &clp);
530:   for (f = 0; f < Nf; ++f) {
531:     for (cp = 0; cp < Np*2; cp += 2) {
532:       PetscSectionGetFieldDof(s, points[cp], f, &dof);
533:       if (!dof) continue;
534:       if (fieldActive[f]) DMPlexBasisTransformField_Internal(dm, tdm, tv, points[cp], f, l2g, &a[d]);
535:       d += dof;
536:     }
537:   }
538:   DMPlexRestoreCompressedClosure(dm, s, p, &Np, &points, &clSection, &clPoints, &clp);
539:   return 0;
540: }

542: PetscErrorCode DMPlexBasisTransformPointTensor_Internal(DM dm, DM tdm, Vec tv, PetscInt p, PetscBool l2g, PetscInt lda, PetscScalar *a)
543: {
544:   PetscSection    s;
545:   PetscSection    clSection;
546:   IS              clPoints;
547:   const PetscInt *clp;
548:   PetscInt       *points = NULL;
549:   PetscInt        Nf, f, g, Np, cpf, cpg, fdof, gdof, r, c = 0;

551:   DMGetLocalSection(dm, &s);
552:   PetscSectionGetNumFields(s, &Nf);
553:   DMPlexGetCompressedClosure(dm, s, p, &Np, &points, &clSection, &clPoints, &clp);
554:   for (f = 0, r = 0; f < Nf; ++f) {
555:     for (cpf = 0; cpf < Np*2; cpf += 2) {
556:       PetscSectionGetFieldDof(s, points[cpf], f, &fdof);
557:       for (g = 0, c = 0; g < Nf; ++g) {
558:         for (cpg = 0; cpg < Np*2; cpg += 2) {
559:           PetscSectionGetFieldDof(s, points[cpg], g, &gdof);
560:           DMPlexBasisTransformFieldTensor_Internal(dm, tdm, tv, points[cpf], f, points[cpg], g, l2g, lda, &a[r*lda+c]);
561:           c += gdof;
562:         }
563:       }
565:       r += fdof;
566:     }
567:   }
569:   DMPlexRestoreCompressedClosure(dm, s, p, &Np, &points, &clSection, &clPoints, &clp);
570:   return 0;
571: }

573: static PetscErrorCode DMPlexBasisTransform_Internal(DM dm, Vec lv, PetscBool l2g)
574: {
575:   DM                 tdm;
576:   Vec                tv;
577:   PetscSection       ts, s;
578:   const PetscScalar *ta;
579:   PetscScalar       *a, *va;
580:   PetscInt           pStart, pEnd, p, Nf, f;

582:   DMGetBasisTransformDM_Internal(dm, &tdm);
583:   DMGetBasisTransformVec_Internal(dm, &tv);
584:   DMGetLocalSection(tdm, &ts);
585:   DMGetLocalSection(dm, &s);
586:   PetscSectionGetChart(s, &pStart, &pEnd);
587:   PetscSectionGetNumFields(s, &Nf);
588:   VecGetArray(lv, &a);
589:   VecGetArrayRead(tv, &ta);
590:   for (p = pStart; p < pEnd; ++p) {
591:     for (f = 0; f < Nf; ++f) {
592:       DMPlexPointLocalFieldRef(dm, p, f, a, &va);
593:       DMPlexBasisTransformField_Internal(dm, tdm, tv, p, f, l2g, va);
594:     }
595:   }
596:   VecRestoreArray(lv, &a);
597:   VecRestoreArrayRead(tv, &ta);
598:   return 0;
599: }

601: /*@
602:   DMPlexGlobalToLocalBasis - Transform the values in the given local vector from the global basis to the local basis

604:   Input Parameters:
605: + dm - The DM
606: - lv - A local vector with values in the global basis

608:   Output Parameters:
609: . lv - A local vector with values in the local basis

611:   Note: This method is only intended to be called inside DMGlobalToLocal(). It is unlikely that a user will have a local vector full of coefficients for the global basis unless they are reimplementing GlobalToLocal.

613:   Level: developer

615: .seealso: DMPlexLocalToGlobalBasis(), DMGetLocalSection(), DMPlexCreateBasisRotation()
616: @*/
617: PetscErrorCode DMPlexGlobalToLocalBasis(DM dm, Vec lv)
618: {
621:   DMPlexBasisTransform_Internal(dm, lv, PETSC_FALSE);
622:   return 0;
623: }

625: /*@
626:   DMPlexLocalToGlobalBasis - Transform the values in the given local vector from the local basis to the global basis

628:   Input Parameters:
629: + dm - The DM
630: - lv - A local vector with values in the local basis

632:   Output Parameters:
633: . lv - A local vector with values in the global basis

635:   Note: This method is only intended to be called inside DMGlobalToLocal(). It is unlikely that a user would want a local vector full of coefficients for the global basis unless they are reimplementing GlobalToLocal.

637:   Level: developer

639: .seealso: DMPlexGlobalToLocalBasis(), DMGetLocalSection(), DMPlexCreateBasisRotation()
640: @*/
641: PetscErrorCode DMPlexLocalToGlobalBasis(DM dm, Vec lv)
642: {
645:   DMPlexBasisTransform_Internal(dm, lv, PETSC_TRUE);
646:   return 0;
647: }

649: /*@
650:   DMPlexCreateBasisRotation - Create an internal transformation from the global basis, used to specify boundary conditions
651:     and global solutions, to a local basis, appropriate for discretization integrals and assembly.

653:   Input Parameters:
654: + dm    - The DM
655: . alpha - The first Euler angle, and in 2D the only one
656: . beta  - The second Euler angle
657: - gamma - The third Euler angle

659:   Note: Following https://en.wikipedia.org/wiki/Euler_angles, we will specify Euler angles by extrinsic rotations, meaning that
660:   we rotate with respect to a fixed initial coordinate system, the local basis (x-y-z). The global basis (X-Y-Z) is reached as follows:
661:   $ The XYZ system rotates about the z axis by alpha. The X axis is now at angle alpha with respect to the x axis.
662:   $ The XYZ system rotates again about the x axis by beta. The Z axis is now at angle beta with respect to the z axis.
663:   $ The XYZ system rotates a third time about the z axis by gamma.

665:   Level: developer

667: .seealso: DMPlexGlobalToLocalBasis(), DMPlexLocalToGlobalBasis()
668: @*/
669: PetscErrorCode DMPlexCreateBasisRotation(DM dm, PetscReal alpha, PetscReal beta, PetscReal gamma)
670: {
671:   RotCtx        *rc;
672:   PetscInt       cdim;

674:   DMGetCoordinateDim(dm, &cdim);
675:   PetscMalloc1(1, &rc);
676:   dm->transformCtx       = rc;
677:   dm->transformSetUp     = DMPlexBasisTransformSetUp_Rotation_Internal;
678:   dm->transformDestroy   = DMPlexBasisTransformDestroy_Rotation_Internal;
679:   dm->transformGetMatrix = DMPlexBasisTransformGetMatrix_Rotation_Internal;
680:   rc->dim   = cdim;
681:   rc->alpha = alpha;
682:   rc->beta  = beta;
683:   rc->gamma = gamma;
684:   (*dm->transformSetUp)(dm, dm->transformCtx);
685:   DMConstructBasisTransform_Internal(dm);
686:   return 0;
687: }

689: /*@C
690:   DMPlexInsertBoundaryValuesEssential - Insert boundary values into a local vector using a function of the coordinates

692:   Input Parameters:
693: + dm     - The DM, with a PetscDS that matches the problem being constrained
694: . time   - The time
695: . field  - The field to constrain
696: . Nc     - The number of constrained field components, or 0 for all components
697: . comps  - An array of constrained component numbers, or NULL for all components
698: . label  - The DMLabel defining constrained points
699: . numids - The number of DMLabel ids for constrained points
700: . ids    - An array of ids for constrained points
701: . func   - A pointwise function giving boundary values
702: - ctx    - An optional user context for bcFunc

704:   Output Parameter:
705: . locX   - A local vector to receives the boundary values

707:   Level: developer

709: .seealso: DMPlexInsertBoundaryValuesEssentialField(), DMPlexInsertBoundaryValuesEssentialBdField(), DMAddBoundary()
710: @*/
711: 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)
712: {
713:   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
714:   void            **ctxs;
715:   PetscInt          numFields;

717:   DMGetNumFields(dm, &numFields);
718:   PetscCalloc2(numFields,&funcs,numFields,&ctxs);
719:   funcs[field] = func;
720:   ctxs[field]  = ctx;
721:   DMProjectFunctionLabelLocal(dm, time, label, numids, ids, Nc, comps, funcs, ctxs, INSERT_BC_VALUES, locX);
722:   PetscFree2(funcs,ctxs);
723:   return 0;
724: }

726: /*@C
727:   DMPlexInsertBoundaryValuesEssentialField - Insert boundary values into a local vector using a function of the coordinates and field data

729:   Input Parameters:
730: + dm     - The DM, with a PetscDS that matches the problem being constrained
731: . time   - The time
732: . locU   - A local vector with the input solution values
733: . field  - The field to constrain
734: . Nc     - The number of constrained field components, or 0 for all components
735: . comps  - An array of constrained component numbers, or NULL for all components
736: . label  - The DMLabel defining constrained points
737: . numids - The number of DMLabel ids for constrained points
738: . ids    - An array of ids for constrained points
739: . func   - A pointwise function giving boundary values
740: - ctx    - An optional user context for bcFunc

742:   Output Parameter:
743: . locX   - A local vector to receives the boundary values

745:   Level: developer

747: .seealso: DMPlexInsertBoundaryValuesEssential(), DMPlexInsertBoundaryValuesEssentialBdField(), DMAddBoundary()
748: @*/
749: PetscErrorCode DMPlexInsertBoundaryValuesEssentialField(DM dm, PetscReal time, Vec locU, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[],
750:                                                         void (*func)(PetscInt, PetscInt, PetscInt,
751:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
752:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
753:                                                                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[],
754:                                                                      PetscScalar[]),
755:                                                         void *ctx, Vec locX)
756: {
757:   void (**funcs)(PetscInt, PetscInt, PetscInt,
758:                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
759:                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
760:                  PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
761:   void            **ctxs;
762:   PetscInt          numFields;

764:   DMGetNumFields(dm, &numFields);
765:   PetscCalloc2(numFields,&funcs,numFields,&ctxs);
766:   funcs[field] = func;
767:   ctxs[field]  = ctx;
768:   DMProjectFieldLabelLocal(dm, time, label, numids, ids, Nc, comps, locU, funcs, INSERT_BC_VALUES, locX);
769:   PetscFree2(funcs,ctxs);
770:   return 0;
771: }

773: /*@C
774:   DMPlexInsertBoundaryValuesEssentialBdField - Insert boundary values into a local vector using a function of the coodinates and boundary field data

776:   Collective on dm

778:   Input Parameters:
779: + dm     - The DM, with a PetscDS that matches the problem being constrained
780: . time   - The time
781: . locU   - A local vector with the input solution values
782: . field  - The field to constrain
783: . Nc     - The number of constrained field components, or 0 for all components
784: . comps  - An array of constrained component numbers, or NULL for all components
785: . label  - The DMLabel defining constrained points
786: . numids - The number of DMLabel ids for constrained points
787: . ids    - An array of ids for constrained points
788: . func   - A pointwise function giving boundary values, the calling sequence is given in DMProjectBdFieldLabelLocal()
789: - ctx    - An optional user context for bcFunc

791:   Output Parameter:
792: . locX   - A local vector to receive the boundary values

794:   Level: developer

796: .seealso: DMProjectBdFieldLabelLocal(), DMPlexInsertBoundaryValuesEssential(), DMPlexInsertBoundaryValuesEssentialField(), DMAddBoundary()
797: @*/
798: PetscErrorCode DMPlexInsertBoundaryValuesEssentialBdField(DM dm, PetscReal time, Vec locU, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[],
799:                                                           void (*func)(PetscInt, PetscInt, PetscInt,
800:                                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
801:                                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
802:                                                                        PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[],
803:                                                                        PetscScalar[]),
804:                                                           void *ctx, Vec locX)
805: {
806:   void (**funcs)(PetscInt, PetscInt, PetscInt,
807:                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
808:                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
809:                  PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
810:   void            **ctxs;
811:   PetscInt          numFields;

813:   DMGetNumFields(dm, &numFields);
814:   PetscCalloc2(numFields,&funcs,numFields,&ctxs);
815:   funcs[field] = func;
816:   ctxs[field]  = ctx;
817:   DMProjectBdFieldLabelLocal(dm, time, label, numids, ids, Nc, comps, locU, funcs, INSERT_BC_VALUES, locX);
818:   PetscFree2(funcs,ctxs);
819:   return 0;
820: }

822: /*@C
823:   DMPlexInsertBoundaryValuesRiemann - Insert boundary values into a local vector

825:   Input Parameters:
826: + dm     - The DM, with a PetscDS that matches the problem being constrained
827: . time   - The time
828: . faceGeometry - A vector with the FVM face geometry information
829: . cellGeometry - A vector with the FVM cell geometry information
830: . Grad         - A vector with the FVM cell gradient information
831: . field  - The field to constrain
832: . Nc     - The number of constrained field components, or 0 for all components
833: . comps  - An array of constrained component numbers, or NULL for all components
834: . label  - The DMLabel defining constrained points
835: . numids - The number of DMLabel ids for constrained points
836: . ids    - An array of ids for constrained points
837: . func   - A pointwise function giving boundary values
838: - ctx    - An optional user context for bcFunc

840:   Output Parameter:
841: . locX   - A local vector to receives the boundary values

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

845:   Level: developer

847: .seealso: DMPlexInsertBoundaryValuesEssential(), DMPlexInsertBoundaryValuesEssentialField(), DMAddBoundary()
848: @*/
849: 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[],
850:                                                  PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *ctx, Vec locX)
851: {
852:   PetscDS            prob;
853:   PetscSF            sf;
854:   DM                 dmFace, dmCell, dmGrad;
855:   const PetscScalar *facegeom, *cellgeom = NULL, *grad;
856:   const PetscInt    *leaves;
857:   PetscScalar       *x, *fx;
858:   PetscInt           dim, nleaves, loc, fStart, fEnd, pdim, i;
859:   PetscErrorCode     ierru = 0;

861:   DMGetPointSF(dm, &sf);
862:   PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);
863:   nleaves = PetscMax(0, nleaves);
864:   DMGetDimension(dm, &dim);
865:   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
866:   DMGetDS(dm, &prob);
867:   VecGetDM(faceGeometry, &dmFace);
868:   VecGetArrayRead(faceGeometry, &facegeom);
869:   if (cellGeometry) {
870:     VecGetDM(cellGeometry, &dmCell);
871:     VecGetArrayRead(cellGeometry, &cellgeom);
872:   }
873:   if (Grad) {
874:     PetscFV fv;

876:     PetscDSGetDiscretization(prob, field, (PetscObject *) &fv);
877:     VecGetDM(Grad, &dmGrad);
878:     VecGetArrayRead(Grad, &grad);
879:     PetscFVGetNumComponents(fv, &pdim);
880:     DMGetWorkArray(dm, pdim, MPIU_SCALAR, &fx);
881:   }
882:   VecGetArray(locX, &x);
883:   for (i = 0; i < numids; ++i) {
884:     IS              faceIS;
885:     const PetscInt *faces;
886:     PetscInt        numFaces, f;

888:     DMLabelGetStratumIS(label, ids[i], &faceIS);
889:     if (!faceIS) continue; /* No points with that id on this process */
890:     ISGetLocalSize(faceIS, &numFaces);
891:     ISGetIndices(faceIS, &faces);
892:     for (f = 0; f < numFaces; ++f) {
893:       const PetscInt         face = faces[f], *cells;
894:       PetscFVFaceGeom        *fg;

896:       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
897:       PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);
898:       if (loc >= 0) continue;
899:       DMPlexPointLocalRead(dmFace, face, facegeom, &fg);
900:       DMPlexGetSupport(dm, face, &cells);
901:       if (Grad) {
902:         PetscFVCellGeom       *cg;
903:         PetscScalar           *cx, *cgrad;
904:         PetscScalar           *xG;
905:         PetscReal              dx[3];
906:         PetscInt               d;

908:         DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);
909:         DMPlexPointLocalRead(dm, cells[0], x, &cx);
910:         DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);
911:         DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);
912:         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
913:         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
914:         (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);
915:       } else {
916:         PetscScalar *xI;
917:         PetscScalar *xG;

919:         DMPlexPointLocalRead(dm, cells[0], x, &xI);
920:         DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);
921:         ierru = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);
922:         if (ierru) {
923:           ISRestoreIndices(faceIS, &faces);
924:           ISDestroy(&faceIS);
925:           goto cleanup;
926:         }
927:       }
928:     }
929:     ISRestoreIndices(faceIS, &faces);
930:     ISDestroy(&faceIS);
931:   }
932:   cleanup:
933:   VecRestoreArray(locX, &x);
934:   if (Grad) {
935:     DMRestoreWorkArray(dm, pdim, MPIU_SCALAR, &fx);
936:     VecRestoreArrayRead(Grad, &grad);
937:   }
938:   if (cellGeometry) VecRestoreArrayRead(cellGeometry, &cellgeom);
939:   VecRestoreArrayRead(faceGeometry, &facegeom);
940:   ierru;
941:   return 0;
942: }

944: static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
945: {
946:   PetscInt c;
947:   for (c = 0; c < Nc; ++c) u[c] = 0.0;
948:   return 0;
949: }

951: PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
952: {
953:   PetscObject    isZero;
954:   PetscDS        prob;
955:   PetscInt       numBd, b;

957:   DMGetDS(dm, &prob);
958:   PetscDSGetNumBoundary(prob, &numBd);
959:   PetscObjectQuery((PetscObject) locX, "__Vec_bc_zero__", &isZero);
960:   for (b = 0; b < numBd; ++b) {
961:     PetscWeakForm           wf;
962:     DMBoundaryConditionType type;
963:     const char             *name;
964:     DMLabel                 label;
965:     PetscInt                field, Nc;
966:     const PetscInt         *comps;
967:     PetscObject             obj;
968:     PetscClassId            id;
969:     void                  (*bvfunc)(void);
970:     PetscInt                numids;
971:     const PetscInt         *ids;
972:     void                   *ctx;

974:     PetscDSGetBoundary(prob, b, &wf, &type, &name, &label, &numids, &ids, &field, &Nc, &comps, &bvfunc, NULL, &ctx);
975:     if (insertEssential != (type & DM_BC_ESSENTIAL)) continue;
976:     DMGetField(dm, field, NULL, &obj);
977:     PetscObjectGetClassId(obj, &id);
978:     if (id == PETSCFE_CLASSID) {
979:       switch (type) {
980:         /* for FEM, there is no insertion to be done for non-essential boundary conditions */
981:       case DM_BC_ESSENTIAL:
982:         {
983:           PetscSimplePointFunc func = (PetscSimplePointFunc) bvfunc;

985:           if (isZero) func = zero;
986:           DMPlexLabelAddCells(dm,label);
987:           DMPlexInsertBoundaryValuesEssential(dm, time, field, Nc, comps, label, numids, ids, func, ctx, locX);
988:           DMPlexLabelClearCells(dm,label);
989:         }
990:         break;
991:       case DM_BC_ESSENTIAL_FIELD:
992:         {
993:           PetscPointFunc func = (PetscPointFunc) bvfunc;

995:           DMPlexLabelAddCells(dm,label);
996:           DMPlexInsertBoundaryValuesEssentialField(dm, time, locX, field, Nc, comps, label, numids, ids, func, ctx, locX);
997:           DMPlexLabelClearCells(dm,label);
998:         }
999:         break;
1000:       default: break;
1001:       }
1002:     } else if (id == PETSCFV_CLASSID) {
1003:       {
1004:         PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*) = (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) bvfunc;

1006:         if (!faceGeomFVM) continue;
1007:         DMPlexInsertBoundaryValuesRiemann(dm, time, faceGeomFVM, cellGeomFVM, gradFVM, field, Nc, comps, label, numids, ids, func, ctx, locX);
1008:       }
1009:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", field);
1010:   }
1011:   return 0;
1012: }

1014: PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues_Plex(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
1015: {
1016:   PetscObject    isZero;
1017:   PetscDS        prob;
1018:   PetscInt       numBd, b;

1020:   if (!locX) return 0;
1021:   DMGetDS(dm, &prob);
1022:   PetscDSGetNumBoundary(prob, &numBd);
1023:   PetscObjectQuery((PetscObject) locX, "__Vec_bc_zero__", &isZero);
1024:   for (b = 0; b < numBd; ++b) {
1025:     PetscWeakForm           wf;
1026:     DMBoundaryConditionType type;
1027:     const char             *name;
1028:     DMLabel                 label;
1029:     PetscInt                field, Nc;
1030:     const PetscInt         *comps;
1031:     PetscObject             obj;
1032:     PetscClassId            id;
1033:     PetscInt                numids;
1034:     const PetscInt         *ids;
1035:     void                  (*bvfunc)(void);
1036:     void                   *ctx;

1038:     PetscDSGetBoundary(prob, b, &wf, &type, &name, &label, &numids, &ids, &field, &Nc, &comps, NULL, &bvfunc, &ctx);
1039:     if (insertEssential != (type & DM_BC_ESSENTIAL)) continue;
1040:     DMGetField(dm, field, NULL, &obj);
1041:     PetscObjectGetClassId(obj, &id);
1042:     if (id == PETSCFE_CLASSID) {
1043:       switch (type) {
1044:         /* for FEM, there is no insertion to be done for non-essential boundary conditions */
1045:       case DM_BC_ESSENTIAL:
1046:         {
1047:           PetscSimplePointFunc func_t = (PetscSimplePointFunc) bvfunc;

1049:           if (isZero) func_t = zero;
1050:           DMPlexLabelAddCells(dm,label);
1051:           DMPlexInsertBoundaryValuesEssential(dm, time, field, Nc, comps, label, numids, ids, func_t, ctx, locX);
1052:           DMPlexLabelClearCells(dm,label);
1053:         }
1054:         break;
1055:       case DM_BC_ESSENTIAL_FIELD:
1056:         {
1057:           PetscPointFunc func_t = (PetscPointFunc) bvfunc;

1059:           DMPlexLabelAddCells(dm,label);
1060:           DMPlexInsertBoundaryValuesEssentialField(dm, time, locX, field, Nc, comps, label, numids, ids, func_t, ctx, locX);
1061:           DMPlexLabelClearCells(dm,label);
1062:         }
1063:         break;
1064:       default: break;
1065:       }
1066:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", field);
1067:   }
1068:   return 0;
1069: }

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

1074:   Not Collective

1076:   Input Parameters:
1077: + dm - The DM
1078: . insertEssential - Should I insert essential (e.g. Dirichlet) or inessential (e.g. Neumann) boundary conditions
1079: . time - The time
1080: . faceGeomFVM - Face geometry data for FV discretizations
1081: . cellGeomFVM - Cell geometry data for FV discretizations
1082: - gradFVM - Gradient reconstruction data for FV discretizations

1084:   Output Parameters:
1085: . locX - Solution updated with boundary values

1087:   Level: intermediate

1089: .seealso: DMProjectFunctionLabelLocal(), DMAddBoundary()
1090: @*/
1091: PetscErrorCode DMPlexInsertBoundaryValues(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
1092: {
1098:   PetscTryMethod(dm,"DMPlexInsertBoundaryValues_C",(DM,PetscBool,Vec,PetscReal,Vec,Vec,Vec),(dm,insertEssential,locX,time,faceGeomFVM,cellGeomFVM,gradFVM));
1099:   return 0;
1100: }

1102: /*@
1103:   DMPlexInsertTimeDerivativeBoundaryValues - Puts coefficients which represent boundary values of the time derivative into the local solution vector

1105:   Input Parameters:
1106: + dm - The DM
1107: . insertEssential - Should I insert essential (e.g. Dirichlet) or inessential (e.g. Neumann) boundary conditions
1108: . time - The time
1109: . faceGeomFVM - Face geometry data for FV discretizations
1110: . cellGeomFVM - Cell geometry data for FV discretizations
1111: - gradFVM - Gradient reconstruction data for FV discretizations

1113:   Output Parameters:
1114: . locX_t - Solution updated with boundary values

1116:   Level: developer

1118: .seealso: DMProjectFunctionLabelLocal()
1119: @*/
1120: PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues(DM dm, PetscBool insertEssential, Vec locX_t, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
1121: {
1127:   PetscTryMethod(dm,"DMPlexInsertTimeDerviativeBoundaryValues_C",(DM,PetscBool,Vec,PetscReal,Vec,Vec,Vec),(dm,insertEssential,locX_t,time,faceGeomFVM,cellGeomFVM,gradFVM));
1128:   return 0;
1129: }

1131: PetscErrorCode DMComputeL2Diff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
1132: {
1133:   Vec              localX;

1135:   DMGetLocalVector(dm, &localX);
1136:   DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localX, time, NULL, NULL, NULL);
1137:   DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);
1138:   DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);
1139:   DMPlexComputeL2DiffLocal(dm, time, funcs, ctxs, localX, diff);
1140:   DMRestoreLocalVector(dm, &localX);
1141:   return 0;
1142: }

1144: /*@C
1145:   DMComputeL2DiffLocal - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.

1147:   Collective on dm

1149:   Input Parameters:
1150: + dm     - The DM
1151: . time   - The time
1152: . funcs  - The functions to evaluate for each field component
1153: . ctxs   - Optional array of contexts to pass to each function, or NULL.
1154: - localX - The coefficient vector u_h, a local vector

1156:   Output Parameter:
1157: . diff - The diff ||u - u_h||_2

1159:   Level: developer

1161: .seealso: DMProjectFunction(), DMComputeL2FieldDiff(), DMComputeL2GradientDiff()
1162: @*/
1163: PetscErrorCode DMPlexComputeL2DiffLocal(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec localX, PetscReal *diff)
1164: {
1165:   const PetscInt   debug = ((DM_Plex*)dm->data)->printL2;
1166:   DM               tdm;
1167:   Vec              tv;
1168:   PetscSection     section;
1169:   PetscQuadrature  quad;
1170:   PetscFEGeom      fegeom;
1171:   PetscScalar     *funcVal, *interpolant;
1172:   PetscReal       *coords, *gcoords;
1173:   PetscReal        localDiff = 0.0;
1174:   const PetscReal *quadWeights;
1175:   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cellHeight, cStart, cEnd, c, field, fieldOffset;
1176:   PetscBool        transform;
1177:   PetscErrorCode   ierr;

1179:   DMGetDimension(dm, &dim);
1180:   DMGetCoordinateDim(dm, &coordDim);
1181:   fegeom.dimEmbed = coordDim;
1182:   DMGetLocalSection(dm, &section);
1183:   PetscSectionGetNumFields(section, &numFields);
1184:   DMGetBasisTransformDM_Internal(dm, &tdm);
1185:   DMGetBasisTransformVec_Internal(dm, &tv);
1186:   DMHasBasisTransform(dm, &transform);
1187:   for (field = 0; field < numFields; ++field) {
1188:     PetscObject  obj;
1189:     PetscClassId id;
1190:     PetscInt     Nc;

1192:     DMGetField(dm, field, NULL, &obj);
1193:     PetscObjectGetClassId(obj, &id);
1194:     if (id == PETSCFE_CLASSID) {
1195:       PetscFE fe = (PetscFE) obj;

1197:       PetscFEGetQuadrature(fe, &quad);
1198:       PetscFEGetNumComponents(fe, &Nc);
1199:     } else if (id == PETSCFV_CLASSID) {
1200:       PetscFV fv = (PetscFV) obj;

1202:       PetscFVGetQuadrature(fv, &quad);
1203:       PetscFVGetNumComponents(fv, &Nc);
1204:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", field);
1205:     numComponents += Nc;
1206:   }
1207:   PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights);
1209:   PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&fegeom.detJ,coordDim*coordDim*Nq,&fegeom.J,coordDim*coordDim*Nq,&fegeom.invJ);
1210:   DMPlexGetVTKCellHeight(dm, &cellHeight);
1211:   DMPlexGetSimplexOrBoxCells(dm, cellHeight, &cStart, &cEnd);
1212:   for (c = cStart; c < cEnd; ++c) {
1213:     PetscScalar *x = NULL;
1214:     PetscReal    elemDiff = 0.0;
1215:     PetscInt     qc = 0;

1217:     DMPlexComputeCellGeometryFEM(dm, c, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ);
1218:     DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);

1220:     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1221:       PetscObject  obj;
1222:       PetscClassId id;
1223:       void * const ctx = ctxs ? ctxs[field] : NULL;
1224:       PetscInt     Nb, Nc, q, fc;

1226:       DMGetField(dm, field, NULL, &obj);
1227:       PetscObjectGetClassId(obj, &id);
1228:       if (id == PETSCFE_CLASSID)      {PetscFEGetNumComponents((PetscFE) obj, &Nc));PetscCall(PetscFEGetDimension((PetscFE) obj, &Nb);}
1229:       else if (id == PETSCFV_CLASSID) {PetscFVGetNumComponents((PetscFV) obj, &Nc);Nb = 1;}
1230:       else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", field);
1231:       if (debug) {
1232:         char title[1024];
1233:         PetscSNPrintf(title, 1023, "Solution for Field %D", field);
1234:         DMPrintCellVector(c, title, Nb, &x[fieldOffset]);
1235:       }
1236:       for (q = 0; q < Nq; ++q) {
1237:         PetscFEGeom qgeom;

1239:         qgeom.dimEmbed = fegeom.dimEmbed;
1240:         qgeom.J        = &fegeom.J[q*coordDim*coordDim];
1241:         qgeom.invJ     = &fegeom.invJ[q*coordDim*coordDim];
1242:         qgeom.detJ     = &fegeom.detJ[q];
1244:         if (transform) {
1245:           gcoords = &coords[coordDim*Nq];
1246:           DMPlexBasisTransformApplyReal_Internal(dm, &coords[coordDim*q], PETSC_TRUE, coordDim, &coords[coordDim*q], gcoords, dm->transformCtx);
1247:         } else {
1248:           gcoords = &coords[coordDim*q];
1249:         }
1250:         PetscArrayzero(funcVal,Nc);
1251:         (*funcs[field])(coordDim, time, gcoords, Nc, funcVal, ctx);
1252:         if (ierr) {
1253:           DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);
1254:           DMRestoreLocalVector(dm, &localX);
1255:           PetscFree6(funcVal,interpolant,coords,fegeom.detJ,fegeom.J,fegeom.invJ);
1256:           ierr;
1257:         }
1258:         if (transform) DMPlexBasisTransformApply_Internal(dm, &coords[coordDim*q], PETSC_FALSE, Nc, funcVal, funcVal, dm->transformCtx);
1259:         if (id == PETSCFE_CLASSID)      PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], &qgeom, q, interpolant);
1260:         else if (id == PETSCFV_CLASSID) PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);
1261:         else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1262:         for (fc = 0; fc < Nc; ++fc) {
1263:           const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)];
1264:           if (debug) PetscPrintf(PETSC_COMM_SELF, "    elem %" PetscInt_FMT " field %" PetscInt_FMT ",%" PetscInt_FMT " point %g %g %g diff %g (%g, %g)\n", c, field, fc, (double)(coordDim > 0 ? coords[coordDim*q] : 0.), (double)(coordDim > 1 ? coords[coordDim*q+1] : 0.),(double)(coordDim > 2 ? coords[coordDim*q+2] : 0.), (double)(PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*fegeom.detJ[q]), PetscRealPart(interpolant[fc]), PetscRealPart(funcVal[fc]));
1265:           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*fegeom.detJ[q];
1266:         }
1267:       }
1268:       fieldOffset += Nb;
1269:       qc += Nc;
1270:     }
1271:     DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);
1272:     if (debug) PetscPrintf(PETSC_COMM_SELF, "  elem %" PetscInt_FMT " diff %g\n", c, (double)elemDiff);
1273:     localDiff += elemDiff;
1274:   }
1275:   PetscFree6(funcVal,interpolant,coords,fegeom.detJ,fegeom.J,fegeom.invJ);
1276:   MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));
1277:   *diff = PetscSqrtReal(*diff);
1278:   return 0;
1279: }

1281: 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)
1282: {
1283:   const PetscInt   debug = ((DM_Plex*)dm->data)->printL2;
1284:   DM               tdm;
1285:   PetscSection     section;
1286:   PetscQuadrature  quad;
1287:   Vec              localX, tv;
1288:   PetscScalar     *funcVal, *interpolant;
1289:   const PetscReal *quadWeights;
1290:   PetscFEGeom      fegeom;
1291:   PetscReal       *coords, *gcoords;
1292:   PetscReal        localDiff = 0.0;
1293:   PetscInt         dim, coordDim, qNc = 0, Nq = 0, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset;
1294:   PetscBool        transform;
1295:   PetscErrorCode   ierr;

1297:   DMGetDimension(dm, &dim);
1298:   DMGetCoordinateDim(dm, &coordDim);
1299:   fegeom.dimEmbed = coordDim;
1300:   DMGetLocalSection(dm, &section);
1301:   PetscSectionGetNumFields(section, &numFields);
1302:   DMGetLocalVector(dm, &localX);
1303:   DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);
1304:   DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);
1305:   DMGetBasisTransformDM_Internal(dm, &tdm);
1306:   DMGetBasisTransformVec_Internal(dm, &tv);
1307:   DMHasBasisTransform(dm, &transform);
1308:   for (field = 0; field < numFields; ++field) {
1309:     PetscFE  fe;
1310:     PetscInt Nc;

1312:     DMGetField(dm, field, NULL, (PetscObject *) &fe);
1313:     PetscFEGetQuadrature(fe, &quad);
1314:     PetscFEGetNumComponents(fe, &Nc);
1315:     numComponents += Nc;
1316:   }
1317:   PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights);
1319:   /* DMProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX); */
1320:   PetscMalloc6(numComponents,&funcVal,coordDim*Nq,&coords,coordDim*coordDim*Nq,&fegeom.J,coordDim*coordDim*Nq,&fegeom.invJ,numComponents*coordDim,&interpolant,Nq,&fegeom.detJ);
1321:   DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);
1322:   for (c = cStart; c < cEnd; ++c) {
1323:     PetscScalar *x = NULL;
1324:     PetscReal    elemDiff = 0.0;
1325:     PetscInt     qc = 0;

1327:     DMPlexComputeCellGeometryFEM(dm, c, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ);
1328:     DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);

1330:     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1331:       PetscFE          fe;
1332:       void * const     ctx = ctxs ? ctxs[field] : NULL;
1333:       PetscInt         Nb, Nc, q, fc;

1335:       DMGetField(dm, field, NULL, (PetscObject *) &fe);
1336:       PetscFEGetDimension(fe, &Nb);
1337:       PetscFEGetNumComponents(fe, &Nc);
1338:       if (debug) {
1339:         char title[1024];
1340:         PetscSNPrintf(title, 1023, "Solution for Field %" PetscInt_FMT, field);
1341:         DMPrintCellVector(c, title, Nb, &x[fieldOffset]);
1342:       }
1343:       for (q = 0; q < Nq; ++q) {
1344:         PetscFEGeom qgeom;

1346:         qgeom.dimEmbed = fegeom.dimEmbed;
1347:         qgeom.J        = &fegeom.J[q*coordDim*coordDim];
1348:         qgeom.invJ     = &fegeom.invJ[q*coordDim*coordDim];
1349:         qgeom.detJ     = &fegeom.detJ[q];
1351:         if (transform) {
1352:           gcoords = &coords[coordDim*Nq];
1353:           DMPlexBasisTransformApplyReal_Internal(dm, &coords[coordDim*q], PETSC_TRUE, coordDim, &coords[coordDim*q], gcoords, dm->transformCtx);
1354:         } else {
1355:           gcoords = &coords[coordDim*q];
1356:         }
1357:         PetscArrayzero(funcVal,Nc);
1358:         (*funcs[field])(coordDim, time, gcoords, n, Nc, funcVal, ctx);
1359:         if (ierr) {
1360:           DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);
1361:           DMRestoreLocalVector(dm, &localX);
1362:           PetscFree6(funcVal,coords,fegeom.J,fegeom.invJ,interpolant,fegeom.detJ);
1363:           ierr;
1364:         }
1365:         if (transform) DMPlexBasisTransformApply_Internal(dm, &coords[coordDim*q], PETSC_FALSE, Nc, funcVal, funcVal, dm->transformCtx);
1366:         PetscFEInterpolateGradient_Static(fe, 1, &x[fieldOffset], &qgeom, q, interpolant);
1367:         /* Overwrite with the dot product if the normal is given */
1368:         if (n) {
1369:           for (fc = 0; fc < Nc; ++fc) {
1370:             PetscScalar sum = 0.0;
1371:             PetscInt    d;
1372:             for (d = 0; d < dim; ++d) sum += interpolant[fc*dim+d]*n[d];
1373:             interpolant[fc] = sum;
1374:           }
1375:         }
1376:         for (fc = 0; fc < Nc; ++fc) {
1377:           const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)];
1378:           if (debug) PetscPrintf(PETSC_COMM_SELF, "    elem %D fieldDer %D,%D diff %g\n", c, field, fc, (double)(PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*fegeom.detJ[q]));
1379:           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*fegeom.detJ[q];
1380:         }
1381:       }
1382:       fieldOffset += Nb;
1383:       qc          += Nc;
1384:     }
1385:     DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);
1386:     if (debug) PetscPrintf(PETSC_COMM_SELF, "  elem %D diff %g\n", c, (double)elemDiff);
1387:     localDiff += elemDiff;
1388:   }
1389:   PetscFree6(funcVal,coords,fegeom.J,fegeom.invJ,interpolant,fegeom.detJ);
1390:   DMRestoreLocalVector(dm, &localX);
1391:   MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));
1392:   *diff = PetscSqrtReal(*diff);
1393:   return 0;
1394: }

1396: PetscErrorCode DMComputeL2FieldDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
1397: {
1398:   const PetscInt   debug = ((DM_Plex*)dm->data)->printL2;
1399:   DM               tdm;
1400:   DMLabel          depthLabel;
1401:   PetscSection     section;
1402:   Vec              localX, tv;
1403:   PetscReal       *localDiff;
1404:   PetscInt         dim, depth, dE, Nf, f, Nds, s;
1405:   PetscBool        transform;
1406:   PetscErrorCode   ierr;

1408:   DMGetDimension(dm, &dim);
1409:   DMGetCoordinateDim(dm, &dE);
1410:   DMGetLocalSection(dm, &section);
1411:   DMGetLocalVector(dm, &localX);
1412:   DMGetBasisTransformDM_Internal(dm, &tdm);
1413:   DMGetBasisTransformVec_Internal(dm, &tv);
1414:   DMHasBasisTransform(dm, &transform);
1415:   DMGetNumFields(dm, &Nf);
1416:   DMPlexGetDepthLabel(dm, &depthLabel);
1417:   DMLabelGetNumValues(depthLabel, &depth);

1419:   VecSet(localX, 0.0);
1420:   DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);
1421:   DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);
1422:   DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);
1423:   DMGetNumDS(dm, &Nds);
1424:   PetscCalloc1(Nf, &localDiff);
1425:   for (s = 0; s < Nds; ++s) {
1426:     PetscDS          ds;
1427:     DMLabel          label;
1428:     IS               fieldIS, pointIS;
1429:     const PetscInt  *fields, *points = NULL;
1430:     PetscQuadrature  quad;
1431:     const PetscReal *quadPoints, *quadWeights;
1432:     PetscFEGeom      fegeom;
1433:     PetscReal       *coords, *gcoords;
1434:     PetscScalar     *funcVal, *interpolant;
1435:     PetscBool        isCohesive;
1436:     PetscInt         qNc, Nq, totNc, cStart = 0, cEnd, c, dsNf;

1438:     DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds);
1439:     ISGetIndices(fieldIS, &fields);
1440:     PetscDSIsCohesive(ds, &isCohesive);
1441:     PetscDSGetNumFields(ds, &dsNf);
1442:     PetscDSGetTotalComponents(ds, &totNc);
1443:     PetscDSGetQuadrature(ds, &quad);
1444:     PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
1446:     PetscCalloc6(totNc, &funcVal, totNc, &interpolant, dE*(Nq+1), &coords,Nq, &fegeom.detJ, dE*dE*Nq, &fegeom.J, dE*dE*Nq, &fegeom.invJ);
1447:     if (!label) {
1448:       DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);
1449:     } else {
1450:       DMLabelGetStratumIS(label, 1, &pointIS);
1451:       ISGetLocalSize(pointIS, &cEnd);
1452:       ISGetIndices(pointIS, &points);
1453:     }
1454:     for (c = cStart; c < cEnd; ++c) {
1455:       const PetscInt cell  = points ? points[c] : c;
1456:       PetscScalar    *x    = NULL;
1457:       const PetscInt *cone;
1458:       PetscInt        qc   = 0, fOff = 0, dep;

1460:       DMLabelGetValue(depthLabel, cell, &dep);
1461:       if (dep != depth-1) continue;
1462:       if (isCohesive) {
1463:         DMPlexGetCone(dm, cell, &cone);
1464:         DMPlexComputeCellGeometryFEM(dm, cone[0], quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ);
1465:       } else {
1466:         DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ);
1467:       }
1468:       DMPlexVecGetClosure(dm, NULL, localX, cell, NULL, &x);
1469:       for (f = 0; f < dsNf; ++f) {
1470:         PetscObject  obj;
1471:         PetscClassId id;
1472:         void * const ctx = ctxs ? ctxs[fields[f]] : NULL;
1473:         PetscInt     Nb, Nc, q, fc;
1474:         PetscReal    elemDiff = 0.0;
1475:         PetscBool    cohesive;

1477:         PetscDSGetCohesive(ds, f, &cohesive);
1478:         if (isCohesive && !cohesive) continue;
1479:         PetscDSGetDiscretization(ds, f, &obj);
1480:         PetscObjectGetClassId(obj, &id);
1481:         if (id == PETSCFE_CLASSID)      {PetscFEGetNumComponents((PetscFE) obj, &Nc));PetscCall(PetscFEGetDimension((PetscFE) obj, &Nb);}
1482:         else if (id == PETSCFV_CLASSID) {PetscFVGetNumComponents((PetscFV) obj, &Nc);Nb = 1;}
1483:         else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", fields[f]);
1484:         if (debug) {
1485:           char title[1024];
1486:           PetscSNPrintf(title, 1023, "Solution for Field %D", fields[f]);
1487:           DMPrintCellVector(cell, title, Nb, &x[fOff]);
1488:         }
1489:         for (q = 0; q < Nq; ++q) {
1490:           PetscFEGeom qgeom;

1492:           qgeom.dimEmbed = fegeom.dimEmbed;
1493:           qgeom.J        = &fegeom.J[q*dE*dE];
1494:           qgeom.invJ     = &fegeom.invJ[q*dE*dE];
1495:           qgeom.detJ     = &fegeom.detJ[q];
1497:           if (transform) {
1498:             gcoords = &coords[dE*Nq];
1499:             DMPlexBasisTransformApplyReal_Internal(dm, &coords[dE*q], PETSC_TRUE, dE, &coords[dE*q], gcoords, dm->transformCtx);
1500:           } else {
1501:             gcoords = &coords[dE*q];
1502:           }
1503:           for (fc = 0; fc < Nc; ++fc) funcVal[fc] = 0.;
1504:           (*funcs[fields[f]])(dE, time, gcoords, Nc, funcVal, ctx);
1505:           if (ierr) {
1506:             DMPlexVecRestoreClosure(dm, NULL, localX, cell, NULL, &x);
1507:             DMRestoreLocalVector(dm, &localX);
1508:             PetscFree6(funcVal,interpolant,coords,fegeom.detJ,fegeom.J,fegeom.invJ);
1509:             ierr;
1510:           }
1511:           if (transform) DMPlexBasisTransformApply_Internal(dm, &coords[dE*q], PETSC_FALSE, Nc, funcVal, funcVal, dm->transformCtx);
1512:           /* Call once for each face, except for lagrange field */
1513:           if (id == PETSCFE_CLASSID)      PetscFEInterpolate_Static((PetscFE) obj, &x[fOff], &qgeom, q, interpolant);
1514:           else if (id == PETSCFV_CLASSID) PetscFVInterpolate_Static((PetscFV) obj, &x[fOff], q, interpolant);
1515:           else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", fields[f]);
1516:           for (fc = 0; fc < Nc; ++fc) {
1517:             const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)];
1518:             if (debug) PetscPrintf(PETSC_COMM_SELF, "    cell %D field %D,%D point %g %g %g diff %g\n", cell, fields[f], fc, (double)(dE > 0 ? coords[dE*q] : 0.), (double)(dE > 1 ? coords[dE*q+1] : 0.),(double)(dE > 2 ? coords[dE*q+2] : 0.), (double)(PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*fegeom.detJ[q]));
1519:             elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*fegeom.detJ[q];
1520:           }
1521:         }
1522:         fOff += Nb;
1523:         qc   += Nc;
1524:         localDiff[fields[f]] += elemDiff;
1525:         if (debug) PetscPrintf(PETSC_COMM_SELF, "  cell %D field %D cum diff %g\n", cell, fields[f], (double)localDiff[fields[f]]);
1526:       }
1527:       DMPlexVecRestoreClosure(dm, NULL, localX, cell, NULL, &x);
1528:     }
1529:     if (label) {
1530:       ISRestoreIndices(pointIS, &points);
1531:       ISDestroy(&pointIS);
1532:     }
1533:     ISRestoreIndices(fieldIS, &fields);
1534:     PetscFree6(funcVal, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ);
1535:   }
1536:   DMRestoreLocalVector(dm, &localX);
1537:   MPIU_Allreduce(localDiff, diff, Nf, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));
1538:   PetscFree(localDiff);
1539:   for (f = 0; f < Nf; ++f) diff[f] = PetscSqrtReal(diff[f]);
1540:   return 0;
1541: }

1543: /*@C
1544:   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.

1546:   Collective on dm

1548:   Input Parameters:
1549: + dm    - The DM
1550: . time  - The time
1551: . funcs - The functions to evaluate for each field component: NULL means that component does not contribute to error calculation
1552: . ctxs  - Optional array of contexts to pass to each function, or NULL.
1553: - X     - The coefficient vector u_h

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

1558:   Level: developer

1560: .seealso: DMProjectFunction(), DMComputeL2Diff(), DMPlexComputeL2FieldDiff(), DMComputeL2GradientDiff()
1561: @*/
1562: PetscErrorCode DMPlexComputeL2DiffVec(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, Vec D)
1563: {
1564:   PetscSection     section;
1565:   PetscQuadrature  quad;
1566:   Vec              localX;
1567:   PetscFEGeom      fegeom;
1568:   PetscScalar     *funcVal, *interpolant;
1569:   PetscReal       *coords;
1570:   const PetscReal *quadPoints, *quadWeights;
1571:   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, c, field, fieldOffset;

1573:   VecSet(D, 0.0);
1574:   DMGetDimension(dm, &dim);
1575:   DMGetCoordinateDim(dm, &coordDim);
1576:   DMGetLocalSection(dm, &section);
1577:   PetscSectionGetNumFields(section, &numFields);
1578:   DMGetLocalVector(dm, &localX);
1579:   DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);
1580:   DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);
1581:   DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);
1582:   for (field = 0; field < numFields; ++field) {
1583:     PetscObject  obj;
1584:     PetscClassId id;
1585:     PetscInt     Nc;

1587:     DMGetField(dm, field, NULL, &obj);
1588:     PetscObjectGetClassId(obj, &id);
1589:     if (id == PETSCFE_CLASSID) {
1590:       PetscFE fe = (PetscFE) obj;

1592:       PetscFEGetQuadrature(fe, &quad);
1593:       PetscFEGetNumComponents(fe, &Nc);
1594:     } else if (id == PETSCFV_CLASSID) {
1595:       PetscFV fv = (PetscFV) obj;

1597:       PetscFVGetQuadrature(fv, &quad);
1598:       PetscFVGetNumComponents(fv, &Nc);
1599:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", field);
1600:     numComponents += Nc;
1601:   }
1602:   PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
1604:   PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&fegeom.detJ,coordDim*coordDim*Nq,&fegeom.J,coordDim*coordDim*Nq,&fegeom.invJ);
1605:   DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);
1606:   for (c = cStart; c < cEnd; ++c) {
1607:     PetscScalar *x = NULL;
1608:     PetscScalar  elemDiff = 0.0;
1609:     PetscInt     qc = 0;

1611:     DMPlexComputeCellGeometryFEM(dm, c, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ);
1612:     DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);

1614:     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1615:       PetscObject  obj;
1616:       PetscClassId id;
1617:       void * const ctx = ctxs ? ctxs[field] : NULL;
1618:       PetscInt     Nb, Nc, q, fc;

1620:       DMGetField(dm, field, NULL, &obj);
1621:       PetscObjectGetClassId(obj, &id);
1622:       if (id == PETSCFE_CLASSID)      {PetscFEGetNumComponents((PetscFE) obj, &Nc));PetscCall(PetscFEGetDimension((PetscFE) obj, &Nb);}
1623:       else if (id == PETSCFV_CLASSID) {PetscFVGetNumComponents((PetscFV) obj, &Nc);Nb = 1;}
1624:       else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", field);
1625:       if (funcs[field]) {
1626:         for (q = 0; q < Nq; ++q) {
1627:           PetscFEGeom qgeom;

1629:           qgeom.dimEmbed = fegeom.dimEmbed;
1630:           qgeom.J        = &fegeom.J[q*coordDim*coordDim];
1631:           qgeom.invJ     = &fegeom.invJ[q*coordDim*coordDim];
1632:           qgeom.detJ     = &fegeom.detJ[q];
1634:           (*funcs[field])(coordDim, time, &coords[q*coordDim], Nc, funcVal, ctx);
1635: #if defined(needs_fix_with_return_code_argument)
1636:           if (ierr) {
1637:             DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);
1638:             PetscFree6(funcVal,interpolant,coords,fegeom.detJ,fegeom.J,fegeom.invJ);
1639:             DMRestoreLocalVector(dm, &localX);
1640:           }
1641: #endif
1642:           if (id == PETSCFE_CLASSID)      PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], &qgeom, q, interpolant);
1643:           else if (id == PETSCFV_CLASSID) PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);
1644:           else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", field);
1645:           for (fc = 0; fc < Nc; ++fc) {
1646:             const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)];
1647:             elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*fegeom.detJ[q];
1648:           }
1649:         }
1650:       }
1651:       fieldOffset += Nb;
1652:       qc          += Nc;
1653:     }
1654:     DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);
1655:     VecSetValue(D, c - cStart, elemDiff, INSERT_VALUES);
1656:   }
1657:   PetscFree6(funcVal,interpolant,coords,fegeom.detJ,fegeom.J,fegeom.invJ);
1658:   DMRestoreLocalVector(dm, &localX);
1659:   VecSqrtAbs(D);
1660:   return 0;
1661: }

1663: /*@
1664:   DMPlexComputeClementInterpolant - This function computes the L2 projection of the cellwise values of a function u onto P1, and stores it in a Vec.

1666:   Collective on dm

1668:   Input Parameters:
1669: + dm - The DM
1670: - locX  - The coefficient vector u_h

1672:   Output Parameter:
1673: . locC - A Vec which holds the Clement interpolant of the function

1675:   Notes:
1676:   u_h(v_i) = \sum_{T_i \in support(v_i)} |T_i| u_h(T_i) / \sum_{T_i \in support(v_i)} |T_i| where |T_i| is the cell volume

1678:   Level: developer

1680: .seealso: DMProjectFunction(), DMComputeL2Diff(), DMPlexComputeL2FieldDiff(), DMComputeL2GradientDiff()
1681: @*/
1682: PetscErrorCode DMPlexComputeClementInterpolant(DM dm, Vec locX, Vec locC)
1683: {
1684:   PetscInt         debug = ((DM_Plex *) dm->data)->printFEM;
1685:   DM               dmc;
1686:   PetscQuadrature  quad;
1687:   PetscScalar     *interpolant, *valsum;
1688:   PetscFEGeom      fegeom;
1689:   PetscReal       *coords;
1690:   const PetscReal *quadPoints, *quadWeights;
1691:   PetscInt         dim, cdim, Nf, f, Nc = 0, Nq, qNc, cStart, cEnd, vStart, vEnd, v;

1693:   PetscCitationsRegister(ClementCitation, &Clementcite);
1694:   VecGetDM(locC, &dmc);
1695:   VecSet(locC, 0.0);
1696:   DMGetDimension(dm, &dim);
1697:   DMGetCoordinateDim(dm, &cdim);
1698:   fegeom.dimEmbed = cdim;
1699:   DMGetNumFields(dm, &Nf);
1701:   for (f = 0; f < Nf; ++f) {
1702:     PetscObject  obj;
1703:     PetscClassId id;
1704:     PetscInt     fNc;

1706:     DMGetField(dm, f, NULL, &obj);
1707:     PetscObjectGetClassId(obj, &id);
1708:     if (id == PETSCFE_CLASSID) {
1709:       PetscFE fe = (PetscFE) obj;

1711:       PetscFEGetQuadrature(fe, &quad);
1712:       PetscFEGetNumComponents(fe, &fNc);
1713:     } else if (id == PETSCFV_CLASSID) {
1714:       PetscFV fv = (PetscFV) obj;

1716:       PetscFVGetQuadrature(fv, &quad);
1717:       PetscFVGetNumComponents(fv, &fNc);
1718:     } else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f);
1719:     Nc += fNc;
1720:   }
1721:   PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
1723:   PetscMalloc6(Nc*2, &valsum, Nc, &interpolant, cdim*Nq, &coords, Nq, &fegeom.detJ, cdim*cdim*Nq, &fegeom.J, cdim*cdim*Nq, &fegeom.invJ);
1724:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1725:   DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);
1726:   for (v = vStart; v < vEnd; ++v) {
1727:     PetscScalar volsum = 0.0;
1728:     PetscInt   *star   = NULL;
1729:     PetscInt    starSize, st, fc;

1731:     PetscArrayzero(valsum, Nc);
1732:     DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);
1733:     for (st = 0; st < starSize*2; st += 2) {
1734:       const PetscInt cell = star[st];
1735:       PetscScalar   *val  = &valsum[Nc];
1736:       PetscScalar   *x    = NULL;
1737:       PetscReal      vol  = 0.0;
1738:       PetscInt       foff = 0;

1740:       if ((cell < cStart) || (cell >= cEnd)) continue;
1741:       DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ);
1742:       DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x);
1743:       for (f = 0; f < Nf; ++f) {
1744:         PetscObject  obj;
1745:         PetscClassId id;
1746:         PetscInt     Nb, fNc, q;

1748:         PetscArrayzero(val, Nc);
1749:         DMGetField(dm, f, NULL, &obj);
1750:         PetscObjectGetClassId(obj, &id);
1751:         if (id == PETSCFE_CLASSID)      {PetscFEGetNumComponents((PetscFE) obj, &fNc));PetscCall(PetscFEGetDimension((PetscFE) obj, &Nb);}
1752:         else if (id == PETSCFV_CLASSID) {PetscFVGetNumComponents((PetscFV) obj, &fNc);Nb = 1;}
1753:         else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f);
1754:         for (q = 0; q < Nq; ++q) {
1755:           const PetscReal wt = quadWeights[q]*fegeom.detJ[q];
1756:           PetscFEGeom     qgeom;

1758:           qgeom.dimEmbed = fegeom.dimEmbed;
1759:           qgeom.J        = &fegeom.J[q*cdim*cdim];
1760:           qgeom.invJ     = &fegeom.invJ[q*cdim*cdim];
1761:           qgeom.detJ     = &fegeom.detJ[q];
1763:           if (id == PETSCFE_CLASSID) PetscFEInterpolate_Static((PetscFE) obj, &x[foff], &qgeom, q, interpolant);
1764:           else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f);
1765:           for (fc = 0; fc < fNc; ++fc) val[foff+fc] += interpolant[fc]*wt;
1766:           vol += wt;
1767:         }
1768:         foff += Nb;
1769:       }
1770:       DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x);
1771:       for (fc = 0; fc < Nc; ++fc) valsum[fc] += val[fc];
1772:       volsum += vol;
1773:       if (debug) {
1774:         PetscPrintf(PETSC_COMM_SELF, "Vertex %" PetscInt_FMT " Cell %" PetscInt_FMT " value: [", v, cell);
1775:         for (fc = 0; fc < Nc; ++fc) {
1776:           if (fc) PetscPrintf(PETSC_COMM_SELF, ", ");
1777:           PetscPrintf(PETSC_COMM_SELF, "%g", (double) PetscRealPart(val[fc]));
1778:         }
1779:         PetscPrintf(PETSC_COMM_SELF, "]\n");
1780:       }
1781:     }
1782:     for (fc = 0; fc < Nc; ++fc) valsum[fc] /= volsum;
1783:     DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);
1784:     DMPlexVecSetClosure(dmc, NULL, locC, v, valsum, INSERT_VALUES);
1785:   }
1786:   PetscFree6(valsum, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ);
1787:   return 0;
1788: }

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

1793:   Collective on dm

1795:   Input Parameters:
1796: + dm - The DM
1797: - locX  - The coefficient vector u_h

1799:   Output Parameter:
1800: . locC - A Vec which holds the Clement interpolant of the gradient

1802:   Notes:
1803:   \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

1805:   Level: developer

1807: .seealso: DMProjectFunction(), DMComputeL2Diff(), DMPlexComputeL2FieldDiff(), DMComputeL2GradientDiff()
1808: @*/
1809: PetscErrorCode DMPlexComputeGradientClementInterpolant(DM dm, Vec locX, Vec locC)
1810: {
1811:   DM_Plex         *mesh  = (DM_Plex *) dm->data;
1812:   PetscInt         debug = mesh->printFEM;
1813:   DM               dmC;
1814:   PetscQuadrature  quad;
1815:   PetscScalar     *interpolant, *gradsum;
1816:   PetscFEGeom      fegeom;
1817:   PetscReal       *coords;
1818:   const PetscReal *quadPoints, *quadWeights;
1819:   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, vStart, vEnd, v, field, fieldOffset;

1821:   PetscCitationsRegister(ClementCitation, &Clementcite);
1822:   VecGetDM(locC, &dmC);
1823:   VecSet(locC, 0.0);
1824:   DMGetDimension(dm, &dim);
1825:   DMGetCoordinateDim(dm, &coordDim);
1826:   fegeom.dimEmbed = coordDim;
1827:   DMGetNumFields(dm, &numFields);
1829:   for (field = 0; field < numFields; ++field) {
1830:     PetscObject  obj;
1831:     PetscClassId id;
1832:     PetscInt     Nc;

1834:     DMGetField(dm, field, NULL, &obj);
1835:     PetscObjectGetClassId(obj, &id);
1836:     if (id == PETSCFE_CLASSID) {
1837:       PetscFE fe = (PetscFE) obj;

1839:       PetscFEGetQuadrature(fe, &quad);
1840:       PetscFEGetNumComponents(fe, &Nc);
1841:     } else if (id == PETSCFV_CLASSID) {
1842:       PetscFV fv = (PetscFV) obj;

1844:       PetscFVGetQuadrature(fv, &quad);
1845:       PetscFVGetNumComponents(fv, &Nc);
1846:     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", field);
1847:     numComponents += Nc;
1848:   }
1849:   PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
1851:   PetscMalloc6(coordDim*numComponents*2,&gradsum,coordDim*numComponents,&interpolant,coordDim*Nq,&coords,Nq,&fegeom.detJ,coordDim*coordDim*Nq,&fegeom.J,coordDim*coordDim*Nq,&fegeom.invJ);
1852:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1853:   DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);
1854:   for (v = vStart; v < vEnd; ++v) {
1855:     PetscScalar volsum = 0.0;
1856:     PetscInt   *star   = NULL;
1857:     PetscInt    starSize, st, d, fc;

1859:     PetscArrayzero(gradsum, coordDim*numComponents);
1860:     DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);
1861:     for (st = 0; st < starSize*2; st += 2) {
1862:       const PetscInt cell = star[st];
1863:       PetscScalar   *grad = &gradsum[coordDim*numComponents];
1864:       PetscScalar   *x    = NULL;
1865:       PetscReal      vol  = 0.0;

1867:       if ((cell < cStart) || (cell >= cEnd)) continue;
1868:       DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ);
1869:       DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x);
1870:       for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1871:         PetscObject  obj;
1872:         PetscClassId id;
1873:         PetscInt     Nb, Nc, q, qc = 0;

1875:         PetscArrayzero(grad, coordDim*numComponents);
1876:         DMGetField(dm, field, NULL, &obj);
1877:         PetscObjectGetClassId(obj, &id);
1878:         if (id == PETSCFE_CLASSID)      {PetscFEGetNumComponents((PetscFE) obj, &Nc));PetscCall(PetscFEGetDimension((PetscFE) obj, &Nb);}
1879:         else if (id == PETSCFV_CLASSID) {PetscFVGetNumComponents((PetscFV) obj, &Nc);Nb = 1;}
1880:         else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", field);
1881:         for (q = 0; q < Nq; ++q) {
1882:           PetscFEGeom qgeom;

1884:           qgeom.dimEmbed = fegeom.dimEmbed;
1885:           qgeom.J        = &fegeom.J[q*coordDim*coordDim];
1886:           qgeom.invJ     = &fegeom.invJ[q*coordDim*coordDim];
1887:           qgeom.detJ     = &fegeom.detJ[q];
1889:           if (id == PETSCFE_CLASSID)      PetscFEInterpolateGradient_Static((PetscFE) obj, 1, &x[fieldOffset], &qgeom, q, interpolant);
1890:           else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", field);
1891:           for (fc = 0; fc < Nc; ++fc) {
1892:             const PetscReal wt = quadWeights[q*qNc+qc];

1894:             for (d = 0; d < coordDim; ++d) grad[fc*coordDim+d] += interpolant[fc*dim+d]*wt*fegeom.detJ[q];
1895:           }
1896:           vol += quadWeights[q*qNc]*fegeom.detJ[q];
1897:         }
1898:         fieldOffset += Nb;
1899:         qc          += Nc;
1900:       }
1901:       DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x);
1902:       for (fc = 0; fc < numComponents; ++fc) {
1903:         for (d = 0; d < coordDim; ++d) {
1904:           gradsum[fc*coordDim+d] += grad[fc*coordDim+d];
1905:         }
1906:       }
1907:       volsum += vol;
1908:       if (debug) {
1909:         PetscPrintf(PETSC_COMM_SELF, "Vertex %" PetscInt_FMT " Cell %" PetscInt_FMT " gradient: [", v, cell);
1910:         for (fc = 0; fc < numComponents; ++fc) {
1911:           for (d = 0; d < coordDim; ++d) {
1912:             if (fc || d > 0) PetscPrintf(PETSC_COMM_SELF, ", ");
1913:             PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc*coordDim+d]));
1914:           }
1915:         }
1916:         PetscPrintf(PETSC_COMM_SELF, "]\n");
1917:       }
1918:     }
1919:     for (fc = 0; fc < numComponents; ++fc) {
1920:       for (d = 0; d < coordDim; ++d) gradsum[fc*coordDim+d] /= volsum;
1921:     }
1922:     DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);
1923:     DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES);
1924:   }
1925:   PetscFree6(gradsum,interpolant,coords,fegeom.detJ,fegeom.J,fegeom.invJ);
1926:   return 0;
1927: }

1929: static PetscErrorCode DMPlexComputeIntegral_Internal(DM dm, Vec X, PetscInt cStart, PetscInt cEnd, PetscScalar *cintegral, void *user)
1930: {
1931:   DM                 dmAux = NULL;
1932:   PetscDS            prob,    probAux = NULL;
1933:   PetscSection       section, sectionAux;
1934:   Vec                locX,    locA;
1935:   PetscInt           dim, numCells = cEnd - cStart, c, f;
1936:   PetscBool          useFVM = PETSC_FALSE;
1937:   /* DS */
1938:   PetscInt           Nf,    totDim,    *uOff, *uOff_x, numConstants;
1939:   PetscInt           NfAux, totDimAux, *aOff;
1940:   PetscScalar       *u, *a;
1941:   const PetscScalar *constants;
1942:   /* Geometry */
1943:   PetscFEGeom       *cgeomFEM;
1944:   DM                 dmGrad;
1945:   PetscQuadrature    affineQuad = NULL;
1946:   Vec                cellGeometryFVM = NULL, faceGeometryFVM = NULL, locGrad = NULL;
1947:   PetscFVCellGeom   *cgeomFVM;
1948:   const PetscScalar *lgrad;
1949:   PetscInt           maxDegree;
1950:   DMField            coordField;
1951:   IS                 cellIS;

1953:   DMGetDS(dm, &prob);
1954:   DMGetDimension(dm, &dim);
1955:   DMGetLocalSection(dm, &section);
1956:   DMGetNumFields(dm, &Nf);
1957:   /* Determine which discretizations we have */
1958:   for (f = 0; f < Nf; ++f) {
1959:     PetscObject  obj;
1960:     PetscClassId id;

1962:     PetscDSGetDiscretization(prob, f, &obj);
1963:     PetscObjectGetClassId(obj, &id);
1964:     if (id == PETSCFV_CLASSID) useFVM = PETSC_TRUE;
1965:   }
1966:   /* Get local solution with boundary values */
1967:   DMGetLocalVector(dm, &locX);
1968:   DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locX, 0.0, NULL, NULL, NULL);
1969:   DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);
1970:   DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);
1971:   /* Read DS information */
1972:   PetscDSGetTotalDimension(prob, &totDim);
1973:   PetscDSGetComponentOffsets(prob, &uOff);
1974:   PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);
1975:   ISCreateStride(PETSC_COMM_SELF,numCells,cStart,1,&cellIS);
1976:   PetscDSGetConstants(prob, &numConstants, &constants);
1977:   /* Read Auxiliary DS information */
1978:   DMGetAuxiliaryVec(dm, NULL, 0, 0, &locA);
1979:   if (locA) {
1980:     VecGetDM(locA, &dmAux);
1981:     DMGetDS(dmAux, &probAux);
1982:     PetscDSGetNumFields(probAux, &NfAux);
1983:     DMGetLocalSection(dmAux, &sectionAux);
1984:     PetscDSGetTotalDimension(probAux, &totDimAux);
1985:     PetscDSGetComponentOffsets(probAux, &aOff);
1986:   }
1987:   /* Allocate data  arrays */
1988:   PetscCalloc1(numCells*totDim, &u);
1989:   if (dmAux) PetscMalloc1(numCells*totDimAux, &a);
1990:   /* Read out geometry */
1991:   DMGetCoordinateField(dm,&coordField);
1992:   DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);
1993:   if (maxDegree <= 1) {
1994:     DMFieldCreateDefaultQuadrature(coordField,cellIS,&affineQuad);
1995:     if (affineQuad) {
1996:       DMFieldCreateFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&cgeomFEM);
1997:     }
1998:   }
1999:   if (useFVM) {
2000:     PetscFV   fv = NULL;
2001:     Vec       grad;
2002:     PetscInt  fStart, fEnd;
2003:     PetscBool compGrad;

2005:     for (f = 0; f < Nf; ++f) {
2006:       PetscObject  obj;
2007:       PetscClassId id;

2009:       PetscDSGetDiscretization(prob, f, &obj);
2010:       PetscObjectGetClassId(obj, &id);
2011:       if (id == PETSCFV_CLASSID) {fv = (PetscFV) obj; break;}
2012:     }
2013:     PetscFVGetComputeGradients(fv, &compGrad);
2014:     PetscFVSetComputeGradients(fv, PETSC_TRUE);
2015:     DMPlexComputeGeometryFVM(dm, &cellGeometryFVM, &faceGeometryFVM);
2016:     DMPlexComputeGradientFVM(dm, fv, faceGeometryFVM, cellGeometryFVM, &dmGrad);
2017:     PetscFVSetComputeGradients(fv, compGrad);
2018:     VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);
2019:     /* Reconstruct and limit cell gradients */
2020:     DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
2021:     DMGetGlobalVector(dmGrad, &grad);
2022:     DMPlexReconstructGradients_Internal(dm, fv, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad);
2023:     /* Communicate gradient values */
2024:     DMGetLocalVector(dmGrad, &locGrad);
2025:     DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);
2026:     DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);
2027:     DMRestoreGlobalVector(dmGrad, &grad);
2028:     /* Handle non-essential (e.g. outflow) boundary values */
2029:     DMPlexInsertBoundaryValues(dm, PETSC_FALSE, locX, 0.0, faceGeometryFVM, cellGeometryFVM, locGrad);
2030:     VecGetArrayRead(locGrad, &lgrad);
2031:   }
2032:   /* Read out data from inputs */
2033:   for (c = cStart; c < cEnd; ++c) {
2034:     PetscScalar *x = NULL;
2035:     PetscInt     i;

2037:     DMPlexVecGetClosure(dm, section, locX, c, NULL, &x);
2038:     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
2039:     DMPlexVecRestoreClosure(dm, section, locX, c, NULL, &x);
2040:     if (dmAux) {
2041:       DMPlexVecGetClosure(dmAux, sectionAux, locA, c, NULL, &x);
2042:       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
2043:       DMPlexVecRestoreClosure(dmAux, sectionAux, locA, c, NULL, &x);
2044:     }
2045:   }
2046:   /* Do integration for each field */
2047:   for (f = 0; f < Nf; ++f) {
2048:     PetscObject  obj;
2049:     PetscClassId id;
2050:     PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;

2052:     PetscDSGetDiscretization(prob, f, &obj);
2053:     PetscObjectGetClassId(obj, &id);
2054:     if (id == PETSCFE_CLASSID) {
2055:       PetscFE         fe = (PetscFE) obj;
2056:       PetscQuadrature q;
2057:       PetscFEGeom     *chunkGeom = NULL;
2058:       PetscInt        Nq, Nb;

2060:       PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
2061:       PetscFEGetQuadrature(fe, &q);
2062:       PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL);
2063:       PetscFEGetDimension(fe, &Nb);
2064:       blockSize = Nb*Nq;
2065:       batchSize = numBlocks * blockSize;
2066:       PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
2067:       numChunks = numCells / (numBatches*batchSize);
2068:       Ne        = numChunks*numBatches*batchSize;
2069:       Nr        = numCells % (numBatches*batchSize);
2070:       offset    = numCells - Nr;
2071:       if (!affineQuad) {
2072:         DMFieldCreateFEGeom(coordField,cellIS,q,PETSC_FALSE,&cgeomFEM);
2073:       }
2074:       PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);
2075:       PetscFEIntegrate(prob, f, Ne, chunkGeom, u, probAux, a, cintegral);
2076:       PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&chunkGeom);
2077:       PetscFEIntegrate(prob, f, Nr, chunkGeom, &u[offset*totDim], probAux, &a[offset*totDimAux], &cintegral[offset*Nf]);
2078:       PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&chunkGeom);
2079:       if (!affineQuad) {
2080:         PetscFEGeomDestroy(&cgeomFEM);
2081:       }
2082:     } else if (id == PETSCFV_CLASSID) {
2083:       PetscInt       foff;
2084:       PetscPointFunc obj_func;
2085:       PetscScalar    lint;

2087:       PetscDSGetObjective(prob, f, &obj_func);
2088:       PetscDSGetFieldOffset(prob, f, &foff);
2089:       if (obj_func) {
2090:         for (c = 0; c < numCells; ++c) {
2091:           PetscScalar *u_x;

2093:           DMPlexPointLocalRead(dmGrad, c, lgrad, &u_x);
2094:           obj_func(dim, Nf, NfAux, uOff, uOff_x, &u[totDim*c+foff], NULL, u_x, aOff, NULL, &a[totDimAux*c], NULL, NULL, 0.0, cgeomFVM[c].centroid, numConstants, constants, &lint);
2095:           cintegral[c*Nf+f] += PetscRealPart(lint)*cgeomFVM[c].volume;
2096:         }
2097:       }
2098:     } else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f);
2099:   }
2100:   /* Cleanup data arrays */
2101:   if (useFVM) {
2102:     VecRestoreArrayRead(locGrad, &lgrad);
2103:     VecRestoreArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);
2104:     DMRestoreLocalVector(dmGrad, &locGrad);
2105:     VecDestroy(&faceGeometryFVM);
2106:     VecDestroy(&cellGeometryFVM);
2107:     DMDestroy(&dmGrad);
2108:   }
2109:   if (dmAux) PetscFree(a);
2110:   PetscFree(u);
2111:   /* Cleanup */
2112:   if (affineQuad) {
2113:     PetscFEGeomDestroy(&cgeomFEM);
2114:   }
2115:   PetscQuadratureDestroy(&affineQuad);
2116:   ISDestroy(&cellIS);
2117:   DMRestoreLocalVector(dm, &locX);
2118:   return 0;
2119: }

2121: /*@
2122:   DMPlexComputeIntegralFEM - Form the integral over the domain from the global input X using pointwise functions specified by the user

2124:   Input Parameters:
2125: + dm - The mesh
2126: . X  - Global input vector
2127: - user - The user context

2129:   Output Parameter:
2130: . integral - Integral for each field

2132:   Level: developer

2134: .seealso: DMPlexSNESComputeResidualFEM()
2135: @*/
2136: PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscScalar *integral, void *user)
2137: {
2138:   DM_Plex       *mesh = (DM_Plex *) dm->data;
2139:   PetscScalar   *cintegral, *lintegral;
2140:   PetscInt       Nf, f, cellHeight, cStart, cEnd, cell;

2145:   PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);
2146:   DMGetNumFields(dm, &Nf);
2147:   DMPlexGetVTKCellHeight(dm, &cellHeight);
2148:   DMPlexGetSimplexOrBoxCells(dm, cellHeight, &cStart, &cEnd);
2149:   /* TODO Introduce a loop over large chunks (right now this is a single chunk) */
2150:   PetscCalloc2(Nf, &lintegral, (cEnd-cStart)*Nf, &cintegral);
2151:   DMPlexComputeIntegral_Internal(dm, X, cStart, cEnd, cintegral, user);
2152:   /* Sum up values */
2153:   for (cell = cStart; cell < cEnd; ++cell) {
2154:     const PetscInt c = cell - cStart;

2156:     if (mesh->printFEM > 1) DMPrintCellVector(cell, "Cell Integral", Nf, &cintegral[c*Nf]);
2157:     for (f = 0; f < Nf; ++f) lintegral[f] += cintegral[c*Nf+f];
2158:   }
2159:   MPIU_Allreduce(lintegral, integral, Nf, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject) dm));
2160:   if (mesh->printFEM) {
2161:     PetscPrintf(PetscObjectComm((PetscObject) dm), "Integral:");
2162:     for (f = 0; f < Nf; ++f) PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", (double) PetscRealPart(integral[f]));
2163:     PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");
2164:   }
2165:   PetscFree2(lintegral, cintegral);
2166:   PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);
2167:   return 0;
2168: }

2170: /*@
2171:   DMPlexComputeCellwiseIntegralFEM - Form the vector of cellwise integrals F from the global input X using pointwise functions specified by the user

2173:   Input Parameters:
2174: + dm - The mesh
2175: . X  - Global input vector
2176: - user - The user context

2178:   Output Parameter:
2179: . integral - Cellwise integrals for each field

2181:   Level: developer

2183: .seealso: DMPlexSNESComputeResidualFEM()
2184: @*/
2185: PetscErrorCode DMPlexComputeCellwiseIntegralFEM(DM dm, Vec X, Vec F, void *user)
2186: {
2187:   DM_Plex       *mesh = (DM_Plex *) dm->data;
2188:   DM             dmF;
2189:   PetscSection   sectionF;
2190:   PetscScalar   *cintegral, *af;
2191:   PetscInt       Nf, f, cellHeight, cStart, cEnd, cell;

2196:   PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);
2197:   DMGetNumFields(dm, &Nf);
2198:   DMPlexGetVTKCellHeight(dm, &cellHeight);
2199:   DMPlexGetSimplexOrBoxCells(dm, cellHeight, &cStart, &cEnd);
2200:   /* TODO Introduce a loop over large chunks (right now this is a single chunk) */
2201:   PetscCalloc1((cEnd-cStart)*Nf, &cintegral);
2202:   DMPlexComputeIntegral_Internal(dm, X, cStart, cEnd, cintegral, user);
2203:   /* Put values in F*/
2204:   VecGetDM(F, &dmF);
2205:   DMGetLocalSection(dmF, &sectionF);
2206:   VecGetArray(F, &af);
2207:   for (cell = cStart; cell < cEnd; ++cell) {
2208:     const PetscInt c = cell - cStart;
2209:     PetscInt       dof, off;

2211:     if (mesh->printFEM > 1) DMPrintCellVector(cell, "Cell Integral", Nf, &cintegral[c*Nf]);
2212:     PetscSectionGetDof(sectionF, cell, &dof);
2213:     PetscSectionGetOffset(sectionF, cell, &off);
2215:     for (f = 0; f < Nf; ++f) af[off+f] = cintegral[c*Nf+f];
2216:   }
2217:   VecRestoreArray(F, &af);
2218:   PetscFree(cintegral);
2219:   PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);
2220:   return 0;
2221: }

2223: static PetscErrorCode DMPlexComputeBdIntegral_Internal(DM dm, Vec locX, IS pointIS,
2224:                                                        void (*func)(PetscInt, PetscInt, PetscInt,
2225:                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
2226:                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
2227:                                                                     PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
2228:                                                        PetscScalar *fintegral, void *user)
2229: {
2230:   DM                 plex = NULL, plexA = NULL;
2231:   DMEnclosureType    encAux;
2232:   PetscDS            prob, probAux = NULL;
2233:   PetscSection       section, sectionAux = NULL;
2234:   Vec                locA = NULL;
2235:   DMField            coordField;
2236:   PetscInt           Nf,        totDim,        *uOff, *uOff_x;
2237:   PetscInt           NfAux = 0, totDimAux = 0, *aOff = NULL;
2238:   PetscScalar       *u, *a = NULL;
2239:   const PetscScalar *constants;
2240:   PetscInt           numConstants, f;

2242:   DMGetCoordinateField(dm, &coordField);
2243:   DMConvert(dm, DMPLEX, &plex);
2244:   DMGetDS(dm, &prob);
2245:   DMGetLocalSection(dm, &section);
2246:   PetscSectionGetNumFields(section, &Nf);
2247:   /* Determine which discretizations we have */
2248:   for (f = 0; f < Nf; ++f) {
2249:     PetscObject  obj;
2250:     PetscClassId id;

2252:     PetscDSGetDiscretization(prob, f, &obj);
2253:     PetscObjectGetClassId(obj, &id);
2255:   }
2256:   /* Read DS information */
2257:   PetscDSGetTotalDimension(prob, &totDim);
2258:   PetscDSGetComponentOffsets(prob, &uOff);
2259:   PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);
2260:   PetscDSGetConstants(prob, &numConstants, &constants);
2261:   /* Read Auxiliary DS information */
2262:   DMGetAuxiliaryVec(dm, NULL, 0, 0, &locA);
2263:   if (locA) {
2264:     DM dmAux;

2266:     VecGetDM(locA, &dmAux);
2267:     DMGetEnclosureRelation(dmAux, dm, &encAux);
2268:     DMConvert(dmAux, DMPLEX, &plexA);
2269:     DMGetDS(dmAux, &probAux);
2270:     PetscDSGetNumFields(probAux, &NfAux);
2271:     DMGetLocalSection(dmAux, &sectionAux);
2272:     PetscDSGetTotalDimension(probAux, &totDimAux);
2273:     PetscDSGetComponentOffsets(probAux, &aOff);
2274:   }
2275:   /* Integrate over points */
2276:   {
2277:     PetscFEGeom    *fgeom, *chunkGeom = NULL;
2278:     PetscInt        maxDegree;
2279:     PetscQuadrature qGeom = NULL;
2280:     const PetscInt *points;
2281:     PetscInt        numFaces, face, Nq, field;
2282:     PetscInt        numChunks, chunkSize, chunk, Nr, offset;

2284:     ISGetLocalSize(pointIS, &numFaces);
2285:     ISGetIndices(pointIS, &points);
2286:     PetscCalloc2(numFaces*totDim, &u, locA ? numFaces*totDimAux : 0, &a);
2287:     DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree);
2288:     for (field = 0; field < Nf; ++field) {
2289:       PetscFE fe;

2291:       PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
2292:       if (maxDegree <= 1) DMFieldCreateDefaultQuadrature(coordField, pointIS, &qGeom);
2293:       if (!qGeom) {
2294:         PetscFEGetFaceQuadrature(fe, &qGeom);
2295:         PetscObjectReference((PetscObject) qGeom);
2296:       }
2297:       PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
2298:       DMPlexGetFEGeom(coordField, pointIS, qGeom, PETSC_TRUE, &fgeom);
2299:       for (face = 0; face < numFaces; ++face) {
2300:         const PetscInt point = points[face], *support;
2301:         PetscScalar    *x    = NULL;
2302:         PetscInt       i;

2304:         DMPlexGetSupport(dm, point, &support);
2305:         DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x);
2306:         for (i = 0; i < totDim; ++i) u[face*totDim+i] = x[i];
2307:         DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x);
2308:         if (locA) {
2309:           PetscInt subp;
2310:           DMGetEnclosurePoint(plexA, dm, encAux, support[0], &subp);
2311:           DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x);
2312:           for (i = 0; i < totDimAux; ++i) a[f*totDimAux+i] = x[i];
2313:           DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x);
2314:         }
2315:       }
2316:       /* Get blocking */
2317:       {
2318:         PetscQuadrature q;
2319:         PetscInt        numBatches, batchSize, numBlocks, blockSize;
2320:         PetscInt        Nq, Nb;

2322:         PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
2323:         PetscFEGetQuadrature(fe, &q);
2324:         PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL);
2325:         PetscFEGetDimension(fe, &Nb);
2326:         blockSize = Nb*Nq;
2327:         batchSize = numBlocks * blockSize;
2328:         chunkSize = numBatches*batchSize;
2329:         PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
2330:         numChunks = numFaces / chunkSize;
2331:         Nr        = numFaces % chunkSize;
2332:         offset    = numFaces - Nr;
2333:       }
2334:       /* Do integration for each field */
2335:       for (chunk = 0; chunk < numChunks; ++chunk) {
2336:         PetscFEGeomGetChunk(fgeom, chunk*chunkSize, (chunk+1)*chunkSize, &chunkGeom);
2337:         PetscFEIntegrateBd(prob, field, func, chunkSize, chunkGeom, u, probAux, a, fintegral);
2338:         PetscFEGeomRestoreChunk(fgeom, 0, offset, &chunkGeom);
2339:       }
2340:       PetscFEGeomGetChunk(fgeom, offset, numFaces, &chunkGeom);
2341:       PetscFEIntegrateBd(prob, field, func, Nr, chunkGeom, &u[offset*totDim], probAux, a ? &a[offset*totDimAux] : NULL, &fintegral[offset*Nf]);
2342:       PetscFEGeomRestoreChunk(fgeom, offset, numFaces, &chunkGeom);
2343:       /* Cleanup data arrays */
2344:       DMPlexRestoreFEGeom(coordField, pointIS, qGeom, PETSC_TRUE, &fgeom);
2345:       PetscQuadratureDestroy(&qGeom);
2346:       PetscFree2(u, a);
2347:       ISRestoreIndices(pointIS, &points);
2348:     }
2349:   }
2350:   if (plex)  DMDestroy(&plex);
2351:   if (plexA) DMDestroy(&plexA);
2352:   return 0;
2353: }

2355: /*@
2356:   DMPlexComputeBdIntegral - Form the integral over the specified boundary from the global input X using pointwise functions specified by the user

2358:   Input Parameters:
2359: + dm      - The mesh
2360: . X       - Global input vector
2361: . label   - The boundary DMLabel
2362: . numVals - The number of label values to use, or PETSC_DETERMINE for all values
2363: . vals    - The label values to use, or PETSC_NULL for all values
2364: . func    - The function to integrate along the boundary
2365: - user    - The user context

2367:   Output Parameter:
2368: . integral - Integral for each field

2370:   Level: developer

2372: .seealso: DMPlexComputeIntegralFEM(), DMPlexComputeBdResidualFEM()
2373: @*/
2374: PetscErrorCode DMPlexComputeBdIntegral(DM dm, Vec X, DMLabel label, PetscInt numVals, const PetscInt vals[],
2375:                                        void (*func)(PetscInt, PetscInt, PetscInt,
2376:                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
2377:                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
2378:                                                     PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
2379:                                        PetscScalar *integral, void *user)
2380: {
2381:   Vec            locX;
2382:   PetscSection   section;
2383:   DMLabel        depthLabel;
2384:   IS             facetIS;
2385:   PetscInt       dim, Nf, f, v;

2392:   PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);
2393:   DMPlexGetDepthLabel(dm, &depthLabel);
2394:   DMGetDimension(dm, &dim);
2395:   DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);
2396:   DMGetLocalSection(dm, &section);
2397:   PetscSectionGetNumFields(section, &Nf);
2398:   /* Get local solution with boundary values */
2399:   DMGetLocalVector(dm, &locX);
2400:   DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locX, 0.0, NULL, NULL, NULL);
2401:   DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);
2402:   DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);
2403:   /* Loop over label values */
2404:   PetscArrayzero(integral, Nf);
2405:   for (v = 0; v < numVals; ++v) {
2406:     IS           pointIS;
2407:     PetscInt     numFaces, face;
2408:     PetscScalar *fintegral;

2410:     DMLabelGetStratumIS(label, vals[v], &pointIS);
2411:     if (!pointIS) continue; /* No points with that id on this process */
2412:     {
2413:       IS isectIS;

2415:       /* TODO: Special cases of ISIntersect where it is quick to check a priori if one is a superset of the other */
2416:       ISIntersect_Caching_Internal(facetIS, pointIS, &isectIS);
2417:       ISDestroy(&pointIS);
2418:       pointIS = isectIS;
2419:     }
2420:     ISGetLocalSize(pointIS, &numFaces);
2421:     PetscCalloc1(numFaces*Nf, &fintegral);
2422:     DMPlexComputeBdIntegral_Internal(dm, locX, pointIS, func, fintegral, user);
2423:     /* Sum point contributions into integral */
2424:     for (f = 0; f < Nf; ++f) for (face = 0; face < numFaces; ++face) integral[f] += fintegral[face*Nf+f];
2425:     PetscFree(fintegral);
2426:     ISDestroy(&pointIS);
2427:   }
2428:   DMRestoreLocalVector(dm, &locX);
2429:   ISDestroy(&facetIS);
2430:   PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);
2431:   return 0;
2432: }

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

2437:   Input Parameters:
2438: + dmc  - The coarse mesh
2439: . dmf  - The fine mesh
2440: . isRefined - Flag indicating regular refinement, rather than the same topology
2441: - user - The user context

2443:   Output Parameter:
2444: . In  - The interpolation matrix

2446:   Level: developer

2448: .seealso: DMPlexComputeInterpolatorGeneral(), DMPlexComputeJacobianFEM()
2449: @*/
2450: PetscErrorCode DMPlexComputeInterpolatorNested(DM dmc, DM dmf, PetscBool isRefined, Mat In, void *user)
2451: {
2452:   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
2453:   const char       *name  = "Interpolator";
2454:   PetscFE          *feRef;
2455:   PetscFV          *fvRef;
2456:   PetscSection      fsection, fglobalSection;
2457:   PetscSection      csection, cglobalSection;
2458:   PetscScalar      *elemMat;
2459:   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
2460:   PetscInt          cTotDim=0, rTotDim = 0;

2462:   PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);
2463:   DMGetDimension(dmf, &dim);
2464:   DMGetLocalSection(dmf, &fsection);
2465:   DMGetGlobalSection(dmf, &fglobalSection);
2466:   DMGetLocalSection(dmc, &csection);
2467:   DMGetGlobalSection(dmc, &cglobalSection);
2468:   PetscSectionGetNumFields(fsection, &Nf);
2469:   DMPlexGetSimplexOrBoxCells(dmc, 0, &cStart, &cEnd);
2470:   PetscCalloc2(Nf, &feRef, Nf, &fvRef);
2471:   for (f = 0; f < Nf; ++f) {
2472:     PetscObject  obj, objc;
2473:     PetscClassId id, idc;
2474:     PetscInt     rNb = 0, Nc = 0, cNb = 0;

2476:     DMGetField(dmf, f, NULL, &obj);
2477:     PetscObjectGetClassId(obj, &id);
2478:     if (id == PETSCFE_CLASSID) {
2479:       PetscFE fe = (PetscFE) obj;

2481:       if (isRefined) {
2482:         PetscFERefine(fe, &feRef[f]);
2483:       } else {
2484:         PetscObjectReference((PetscObject) fe);
2485:         feRef[f] = fe;
2486:       }
2487:       PetscFEGetDimension(feRef[f], &rNb);
2488:       PetscFEGetNumComponents(fe, &Nc);
2489:     } else if (id == PETSCFV_CLASSID) {
2490:       PetscFV        fv = (PetscFV) obj;
2491:       PetscDualSpace Q;

2493:       if (isRefined) {
2494:         PetscFVRefine(fv, &fvRef[f]);
2495:       } else {
2496:         PetscObjectReference((PetscObject) fv);
2497:         fvRef[f] = fv;
2498:       }
2499:       PetscFVGetDualSpace(fvRef[f], &Q);
2500:       PetscDualSpaceGetDimension(Q, &rNb);
2501:       PetscFVGetDualSpace(fv, &Q);
2502:       PetscFVGetNumComponents(fv, &Nc);
2503:     }
2504:     DMGetField(dmc, f, NULL, &objc);
2505:     PetscObjectGetClassId(objc, &idc);
2506:     if (idc == PETSCFE_CLASSID) {
2507:       PetscFE fe = (PetscFE) objc;

2509:       PetscFEGetDimension(fe, &cNb);
2510:     } else if (id == PETSCFV_CLASSID) {
2511:       PetscFV        fv = (PetscFV) obj;
2512:       PetscDualSpace Q;

2514:       PetscFVGetDualSpace(fv, &Q);
2515:       PetscDualSpaceGetDimension(Q, &cNb);
2516:     }
2517:     rTotDim += rNb;
2518:     cTotDim += cNb;
2519:   }
2520:   PetscMalloc1(rTotDim*cTotDim,&elemMat);
2521:   PetscArrayzero(elemMat, rTotDim*cTotDim);
2522:   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
2523:     PetscDualSpace   Qref;
2524:     PetscQuadrature  f;
2525:     const PetscReal *qpoints, *qweights;
2526:     PetscReal       *points;
2527:     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;

2529:     /* Compose points from all dual basis functionals */
2530:     if (feRef[fieldI]) {
2531:       PetscFEGetDualSpace(feRef[fieldI], &Qref);
2532:       PetscFEGetNumComponents(feRef[fieldI], &Nc);
2533:     } else {
2534:       PetscFVGetDualSpace(fvRef[fieldI], &Qref);
2535:       PetscFVGetNumComponents(fvRef[fieldI], &Nc);
2536:     }
2537:     PetscDualSpaceGetDimension(Qref, &fpdim);
2538:     for (i = 0; i < fpdim; ++i) {
2539:       PetscDualSpaceGetFunctional(Qref, i, &f);
2540:       PetscQuadratureGetData(f, NULL, NULL, &Np, NULL, NULL);
2541:       npoints += Np;
2542:     }
2543:     PetscMalloc1(npoints*dim,&points);
2544:     for (i = 0, k = 0; i < fpdim; ++i) {
2545:       PetscDualSpaceGetFunctional(Qref, i, &f);
2546:       PetscQuadratureGetData(f, NULL, NULL, &Np, &qpoints, NULL);
2547:       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
2548:     }

2550:     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
2551:       PetscObject  obj;
2552:       PetscClassId id;
2553:       PetscInt     NcJ = 0, cpdim = 0, j, qNc;

2555:       DMGetField(dmc, fieldJ, NULL, &obj);
2556:       PetscObjectGetClassId(obj, &id);
2557:       if (id == PETSCFE_CLASSID) {
2558:         PetscFE           fe = (PetscFE) obj;
2559:         PetscTabulation T  = NULL;

2561:         /* Evaluate basis at points */
2562:         PetscFEGetNumComponents(fe, &NcJ);
2563:         PetscFEGetDimension(fe, &cpdim);
2564:         /* For now, fields only interpolate themselves */
2565:         if (fieldI == fieldJ) {
2567:           PetscFECreateTabulation(fe, 1, npoints, points, 0, &T);
2568:           for (i = 0, k = 0; i < fpdim; ++i) {
2569:             PetscDualSpaceGetFunctional(Qref, i, &f);
2570:             PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights);
2572:             for (p = 0; p < Np; ++p, ++k) {
2573:               for (j = 0; j < cpdim; ++j) {
2574:                 /*
2575:                    cTotDim:            Total columns in element interpolation matrix, sum of number of dual basis functionals in each field
2576:                    offsetI, offsetJ:   Offsets into the larger element interpolation matrix for different fields
2577:                    fpdim, i, cpdim, j: Dofs for fine and coarse grids, correspond to dual space basis functionals
2578:                    qNC, Nc, Ncj, c:    Number of components in this field
2579:                    Np, p:              Number of quad points in the fine grid functional i
2580:                    k:                  i*Np + p, overall point number for the interpolation
2581:                 */
2582:                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i)*cTotDim + offsetJ + j] += T->T[0][k*cpdim*NcJ+j*Nc+c]*qweights[p*qNc+c];
2583:               }
2584:             }
2585:           }
2586:           PetscTabulationDestroy(&T);
2587:         }
2588:       } else if (id == PETSCFV_CLASSID) {
2589:         PetscFV        fv = (PetscFV) obj;

2591:         /* Evaluate constant function at points */
2592:         PetscFVGetNumComponents(fv, &NcJ);
2593:         cpdim = 1;
2594:         /* For now, fields only interpolate themselves */
2595:         if (fieldI == fieldJ) {
2597:           for (i = 0, k = 0; i < fpdim; ++i) {
2598:             PetscDualSpaceGetFunctional(Qref, i, &f);
2599:             PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights);
2601:             for (p = 0; p < Np; ++p, ++k) {
2602:               for (j = 0; j < cpdim; ++j) {
2603:                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i)*cTotDim + offsetJ + j] += 1.0*qweights[p*qNc+c];
2604:               }
2605:             }
2606:           }
2607:         }
2608:       }
2609:       offsetJ += cpdim;
2610:     }
2611:     offsetI += fpdim;
2612:     PetscFree(points);
2613:   }
2614:   if (mesh->printFEM > 1) DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);
2615:   /* Preallocate matrix */
2616:   {
2617:     Mat          preallocator;
2618:     PetscScalar *vals;
2619:     PetscInt    *cellCIndices, *cellFIndices;
2620:     PetscInt     locRows, locCols, cell;

2622:     MatGetLocalSize(In, &locRows, &locCols);
2623:     MatCreate(PetscObjectComm((PetscObject) In), &preallocator);
2624:     MatSetType(preallocator, MATPREALLOCATOR);
2625:     MatSetSizes(preallocator, locRows, locCols, PETSC_DETERMINE, PETSC_DETERMINE);
2626:     MatSetUp(preallocator);
2627:     PetscCalloc3(rTotDim*cTotDim, &vals,cTotDim,&cellCIndices,rTotDim,&cellFIndices);
2628:     for (cell = cStart; cell < cEnd; ++cell) {
2629:       if (isRefined) {
2630:         DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);
2631:         MatSetValues(preallocator, rTotDim, cellFIndices, cTotDim, cellCIndices, vals, INSERT_VALUES);
2632:       } else {
2633:         DMPlexMatSetClosureGeneral(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, preallocator, cell, vals, INSERT_VALUES);
2634:       }
2635:     }
2636:     PetscFree3(vals,cellCIndices,cellFIndices);
2637:     MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY);
2638:     MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY);
2639:     MatPreallocatorPreallocate(preallocator, PETSC_TRUE, In);
2640:     MatDestroy(&preallocator);
2641:   }
2642:   /* Fill matrix */
2643:   MatZeroEntries(In);
2644:   for (c = cStart; c < cEnd; ++c) {
2645:     if (isRefined) {
2646:       DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);
2647:     } else {
2648:       DMPlexMatSetClosureGeneral(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);
2649:     }
2650:   }
2651:   for (f = 0; f < Nf; ++f) PetscFEDestroy(&feRef[f]);
2652:   PetscFree2(feRef,fvRef);
2653:   PetscFree(elemMat);
2654:   MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);
2655:   MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);
2656:   if (mesh->printFEM > 1) {
2657:     PetscPrintf(PetscObjectComm((PetscObject)In), "%s:\n", name);
2658:     MatChop(In, 1.0e-10);
2659:     MatView(In, NULL);
2660:   }
2661:   PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);
2662:   return 0;
2663: }

2665: PetscErrorCode DMPlexComputeMassMatrixNested(DM dmc, DM dmf, Mat mass, void *user)
2666: {
2667:   SETERRQ(PetscObjectComm((PetscObject) dmc), PETSC_ERR_SUP, "Laziness");
2668: }

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

2673:   Input Parameters:
2674: + dmf  - The fine mesh
2675: . dmc  - The coarse mesh
2676: - user - The user context

2678:   Output Parameter:
2679: . In  - The interpolation matrix

2681:   Level: developer

2683: .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM()
2684: @*/
2685: PetscErrorCode DMPlexComputeInterpolatorGeneral(DM dmc, DM dmf, Mat In, void *user)
2686: {
2687:   DM_Plex       *mesh = (DM_Plex *) dmf->data;
2688:   const char    *name = "Interpolator";
2689:   PetscDS        prob;
2690:   PetscSection   fsection, csection, globalFSection, globalCSection;
2691:   PetscHSetIJ    ht;
2692:   PetscLayout    rLayout;
2693:   PetscInt      *dnz, *onz;
2694:   PetscInt       locRows, rStart, rEnd;
2695:   PetscReal     *x, *v0, *J, *invJ, detJ;
2696:   PetscReal     *v0c, *Jc, *invJc, detJc;
2697:   PetscScalar   *elemMat;
2698:   PetscInt       dim, Nf, field, totDim, cStart, cEnd, cell, ccell;

2700:   PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);
2701:   DMGetCoordinateDim(dmc, &dim);
2702:   DMGetDS(dmc, &prob);
2703:   PetscDSGetWorkspace(prob, &x, NULL, NULL, NULL, NULL);
2704:   PetscDSGetNumFields(prob, &Nf);
2705:   PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);
2706:   PetscMalloc3(dim,&v0c,dim*dim,&Jc,dim*dim,&invJc);
2707:   DMGetLocalSection(dmf, &fsection);
2708:   DMGetGlobalSection(dmf, &globalFSection);
2709:   DMGetLocalSection(dmc, &csection);
2710:   DMGetGlobalSection(dmc, &globalCSection);
2711:   DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);
2712:   PetscDSGetTotalDimension(prob, &totDim);
2713:   PetscMalloc1(totDim, &elemMat);

2715:   MatGetLocalSize(In, &locRows, NULL);
2716:   PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);
2717:   PetscLayoutSetLocalSize(rLayout, locRows);
2718:   PetscLayoutSetBlockSize(rLayout, 1);
2719:   PetscLayoutSetUp(rLayout);
2720:   PetscLayoutGetRange(rLayout, &rStart, &rEnd);
2721:   PetscLayoutDestroy(&rLayout);
2722:   PetscCalloc2(locRows,&dnz,locRows,&onz);
2723:   PetscHSetIJCreate(&ht);
2724:   for (field = 0; field < Nf; ++field) {
2725:     PetscObject      obj;
2726:     PetscClassId     id;
2727:     PetscDualSpace   Q = NULL;
2728:     PetscQuadrature  f;
2729:     const PetscReal *qpoints;
2730:     PetscInt         Nc, Np, fpdim, i, d;

2732:     PetscDSGetDiscretization(prob, field, &obj);
2733:     PetscObjectGetClassId(obj, &id);
2734:     if (id == PETSCFE_CLASSID) {
2735:       PetscFE fe = (PetscFE) obj;

2737:       PetscFEGetDualSpace(fe, &Q);
2738:       PetscFEGetNumComponents(fe, &Nc);
2739:     } else if (id == PETSCFV_CLASSID) {
2740:       PetscFV fv = (PetscFV) obj;

2742:       PetscFVGetDualSpace(fv, &Q);
2743:       Nc   = 1;
2744:     }
2745:     PetscDualSpaceGetDimension(Q, &fpdim);
2746:     /* For each fine grid cell */
2747:     for (cell = cStart; cell < cEnd; ++cell) {
2748:       PetscInt *findices,   *cindices;
2749:       PetscInt  numFIndices, numCIndices;

2751:       DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);
2752:       DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);
2754:       for (i = 0; i < fpdim; ++i) {
2755:         Vec             pointVec;
2756:         PetscScalar    *pV;
2757:         PetscSF         coarseCellSF = NULL;
2758:         const PetscSFNode *coarseCells;
2759:         PetscInt        numCoarseCells, q, c;

2761:         /* Get points from the dual basis functional quadrature */
2762:         PetscDualSpaceGetFunctional(Q, i, &f);
2763:         PetscQuadratureGetData(f, NULL, NULL, &Np, &qpoints, NULL);
2764:         VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);
2765:         VecSetBlockSize(pointVec, dim);
2766:         VecGetArray(pointVec, &pV);
2767:         for (q = 0; q < Np; ++q) {
2768:           const PetscReal xi0[3] = {-1., -1., -1.};

2770:           /* Transform point to real space */
2771:           CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q*dim], x);
2772:           for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
2773:         }
2774:         VecRestoreArray(pointVec, &pV);
2775:         /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
2776:         /* OPT: Pack all quad points from fine cell */
2777:         DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);
2778:         PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view");
2779:         /* Update preallocation info */
2780:         PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);
2782:         {
2783:           PetscHashIJKey key;
2784:           PetscBool      missing;

2786:           key.i = findices[i];
2787:           if (key.i >= 0) {
2788:             /* Get indices for coarse elements */
2789:             for (ccell = 0; ccell < numCoarseCells; ++ccell) {
2790:               DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL);
2791:               for (c = 0; c < numCIndices; ++c) {
2792:                 key.j = cindices[c];
2793:                 if (key.j < 0) continue;
2794:                 PetscHSetIJQueryAdd(ht, key, &missing);
2795:                 if (missing) {
2796:                   if ((key.j >= rStart) && (key.j < rEnd)) ++dnz[key.i-rStart];
2797:                   else                                     ++onz[key.i-rStart];
2798:                 }
2799:               }
2800:               DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL);
2801:             }
2802:           }
2803:         }
2804:         PetscSFDestroy(&coarseCellSF);
2805:         VecDestroy(&pointVec);
2806:       }
2807:       DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);
2808:     }
2809:   }
2810:   PetscHSetIJDestroy(&ht);
2811:   MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);
2812:   MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
2813:   PetscFree2(dnz,onz);
2814:   for (field = 0; field < Nf; ++field) {
2815:     PetscObject       obj;
2816:     PetscClassId      id;
2817:     PetscDualSpace    Q = NULL;
2818:     PetscTabulation T = NULL;
2819:     PetscQuadrature   f;
2820:     const PetscReal  *qpoints, *qweights;
2821:     PetscInt          Nc, qNc, Np, fpdim, i, d;

2823:     PetscDSGetDiscretization(prob, field, &obj);
2824:     PetscObjectGetClassId(obj, &id);
2825:     if (id == PETSCFE_CLASSID) {
2826:       PetscFE fe = (PetscFE) obj;

2828:       PetscFEGetDualSpace(fe, &Q);
2829:       PetscFEGetNumComponents(fe, &Nc);
2830:       PetscFECreateTabulation(fe, 1, 1, x, 0, &T);
2831:     } else if (id == PETSCFV_CLASSID) {
2832:       PetscFV fv = (PetscFV) obj;

2834:       PetscFVGetDualSpace(fv, &Q);
2835:       Nc   = 1;
2836:     } else SETERRQ(PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %D",field);
2837:     PetscDualSpaceGetDimension(Q, &fpdim);
2838:     /* For each fine grid cell */
2839:     for (cell = cStart; cell < cEnd; ++cell) {
2840:       PetscInt *findices,   *cindices;
2841:       PetscInt  numFIndices, numCIndices;

2843:       DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);
2844:       DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);
2846:       for (i = 0; i < fpdim; ++i) {
2847:         Vec             pointVec;
2848:         PetscScalar    *pV;
2849:         PetscSF         coarseCellSF = NULL;
2850:         const PetscSFNode *coarseCells;
2851:         PetscInt        numCoarseCells, cpdim, q, c, j;

2853:         /* Get points from the dual basis functional quadrature */
2854:         PetscDualSpaceGetFunctional(Q, i, &f);
2855:         PetscQuadratureGetData(f, NULL, &qNc, &Np, &qpoints, &qweights);
2857:         VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);
2858:         VecSetBlockSize(pointVec, dim);
2859:         VecGetArray(pointVec, &pV);
2860:         for (q = 0; q < Np; ++q) {
2861:           const PetscReal xi0[3] = {-1., -1., -1.};

2863:           /* Transform point to real space */
2864:           CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q*dim], x);
2865:           for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
2866:         }
2867:         VecRestoreArray(pointVec, &pV);
2868:         /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
2869:         /* OPT: Read this out from preallocation information */
2870:         DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);
2871:         /* Update preallocation info */
2872:         PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);
2874:         VecGetArray(pointVec, &pV);
2875:         for (ccell = 0; ccell < numCoarseCells; ++ccell) {
2876:           PetscReal pVReal[3];
2877:           const PetscReal xi0[3] = {-1., -1., -1.};

2879:           DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL);
2880:           /* Transform points from real space to coarse reference space */
2881:           DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc);
2882:           for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell*dim+d]);
2883:           CoordinatesRealToRef(dim, dim, xi0, v0c, invJc, pVReal, x);

2885:           if (id == PETSCFE_CLASSID) {
2886:             PetscFE fe = (PetscFE) obj;

2888:             /* Evaluate coarse basis on contained point */
2889:             PetscFEGetDimension(fe, &cpdim);
2890:             PetscFEComputeTabulation(fe, 1, x, 0, T);
2891:             PetscArrayzero(elemMat, cpdim);
2892:             /* Get elemMat entries by multiplying by weight */
2893:             for (j = 0; j < cpdim; ++j) {
2894:               for (c = 0; c < Nc; ++c) elemMat[j] += T->T[0][j*Nc + c]*qweights[ccell*qNc + c];
2895:             }
2896:           } else {
2897:             cpdim = 1;
2898:             for (j = 0; j < cpdim; ++j) {
2899:               for (c = 0; c < Nc; ++c) elemMat[j] += 1.0*qweights[ccell*qNc + c];
2900:             }
2901:           }
2902:           /* Update interpolator */
2903:           if (mesh->printFEM > 1) DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);
2905:           MatSetValues(In, 1, &findices[i], numCIndices, cindices, elemMat, INSERT_VALUES);
2906:           DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL);
2907:         }
2908:         VecRestoreArray(pointVec, &pV);
2909:         PetscSFDestroy(&coarseCellSF);
2910:         VecDestroy(&pointVec);
2911:       }
2912:       DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);
2913:     }
2914:     if (id == PETSCFE_CLASSID) PetscTabulationDestroy(&T);
2915:   }
2916:   PetscFree3(v0,J,invJ);
2917:   PetscFree3(v0c,Jc,invJc);
2918:   PetscFree(elemMat);
2919:   MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);
2920:   MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);
2921:   PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);
2922:   return 0;
2923: }

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

2928:   Input Parameters:
2929: + dmf  - The fine mesh
2930: . dmc  - The coarse mesh
2931: - user - The user context

2933:   Output Parameter:
2934: . mass  - The mass matrix

2936:   Level: developer

2938: .seealso: DMPlexComputeMassMatrixNested(), DMPlexComputeInterpolatorNested(), DMPlexComputeInterpolatorGeneral(), DMPlexComputeJacobianFEM()
2939: @*/
2940: PetscErrorCode DMPlexComputeMassMatrixGeneral(DM dmc, DM dmf, Mat mass, void *user)
2941: {
2942:   DM_Plex       *mesh = (DM_Plex *) dmf->data;
2943:   const char    *name = "Mass Matrix";
2944:   PetscDS        prob;
2945:   PetscSection   fsection, csection, globalFSection, globalCSection;
2946:   PetscHSetIJ    ht;
2947:   PetscLayout    rLayout;
2948:   PetscInt      *dnz, *onz;
2949:   PetscInt       locRows, rStart, rEnd;
2950:   PetscReal     *x, *v0, *J, *invJ, detJ;
2951:   PetscReal     *v0c, *Jc, *invJc, detJc;
2952:   PetscScalar   *elemMat;
2953:   PetscInt       dim, Nf, field, totDim, cStart, cEnd, cell, ccell;

2955:   DMGetCoordinateDim(dmc, &dim);
2956:   DMGetDS(dmc, &prob);
2957:   PetscDSGetWorkspace(prob, &x, NULL, NULL, NULL, NULL);
2958:   PetscDSGetNumFields(prob, &Nf);
2959:   PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);
2960:   PetscMalloc3(dim,&v0c,dim*dim,&Jc,dim*dim,&invJc);
2961:   DMGetLocalSection(dmf, &fsection);
2962:   DMGetGlobalSection(dmf, &globalFSection);
2963:   DMGetLocalSection(dmc, &csection);
2964:   DMGetGlobalSection(dmc, &globalCSection);
2965:   DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);
2966:   PetscDSGetTotalDimension(prob, &totDim);
2967:   PetscMalloc1(totDim, &elemMat);

2969:   MatGetLocalSize(mass, &locRows, NULL);
2970:   PetscLayoutCreate(PetscObjectComm((PetscObject) mass), &rLayout);
2971:   PetscLayoutSetLocalSize(rLayout, locRows);
2972:   PetscLayoutSetBlockSize(rLayout, 1);
2973:   PetscLayoutSetUp(rLayout);
2974:   PetscLayoutGetRange(rLayout, &rStart, &rEnd);
2975:   PetscLayoutDestroy(&rLayout);
2976:   PetscCalloc2(locRows,&dnz,locRows,&onz);
2977:   PetscHSetIJCreate(&ht);
2978:   for (field = 0; field < Nf; ++field) {
2979:     PetscObject      obj;
2980:     PetscClassId     id;
2981:     PetscQuadrature  quad;
2982:     const PetscReal *qpoints;
2983:     PetscInt         Nq, Nc, i, d;

2985:     PetscDSGetDiscretization(prob, field, &obj);
2986:     PetscObjectGetClassId(obj, &id);
2987:     if (id == PETSCFE_CLASSID) PetscFEGetQuadrature((PetscFE) obj, &quad);
2988:     else                       PetscFVGetQuadrature((PetscFV) obj, &quad);
2989:     PetscQuadratureGetData(quad, NULL, &Nc, &Nq, &qpoints, NULL);
2990:     /* For each fine grid cell */
2991:     for (cell = cStart; cell < cEnd; ++cell) {
2992:       Vec                pointVec;
2993:       PetscScalar       *pV;
2994:       PetscSF            coarseCellSF = NULL;
2995:       const PetscSFNode *coarseCells;
2996:       PetscInt           numCoarseCells, q, c;
2997:       PetscInt          *findices,   *cindices;
2998:       PetscInt           numFIndices, numCIndices;

3000:       DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);
3001:       DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);
3002:       /* Get points from the quadrature */
3003:       VecCreateSeq(PETSC_COMM_SELF, Nq*dim, &pointVec);
3004:       VecSetBlockSize(pointVec, dim);
3005:       VecGetArray(pointVec, &pV);
3006:       for (q = 0; q < Nq; ++q) {
3007:         const PetscReal xi0[3] = {-1., -1., -1.};

3009:         /* Transform point to real space */
3010:         CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q*dim], x);
3011:         for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
3012:       }
3013:       VecRestoreArray(pointVec, &pV);
3014:       /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
3015:       DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);
3016:       PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view");
3017:       /* Update preallocation info */
3018:       PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);
3020:       {
3021:         PetscHashIJKey key;
3022:         PetscBool      missing;

3024:         for (i = 0; i < numFIndices; ++i) {
3025:           key.i = findices[i];
3026:           if (key.i >= 0) {
3027:             /* Get indices for coarse elements */
3028:             for (ccell = 0; ccell < numCoarseCells; ++ccell) {
3029:               DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL);
3030:               for (c = 0; c < numCIndices; ++c) {
3031:                 key.j = cindices[c];
3032:                 if (key.j < 0) continue;
3033:                 PetscHSetIJQueryAdd(ht, key, &missing);
3034:                 if (missing) {
3035:                   if ((key.j >= rStart) && (key.j < rEnd)) ++dnz[key.i-rStart];
3036:                   else                                     ++onz[key.i-rStart];
3037:                 }
3038:               }
3039:               DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL);
3040:             }
3041:           }
3042:         }
3043:       }
3044:       PetscSFDestroy(&coarseCellSF);
3045:       VecDestroy(&pointVec);
3046:       DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);
3047:     }
3048:   }
3049:   PetscHSetIJDestroy(&ht);
3050:   MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL);
3051:   MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
3052:   PetscFree2(dnz,onz);
3053:   for (field = 0; field < Nf; ++field) {
3054:     PetscObject       obj;
3055:     PetscClassId      id;
3056:     PetscTabulation T, Tfine;
3057:     PetscQuadrature   quad;
3058:     const PetscReal  *qpoints, *qweights;
3059:     PetscInt          Nq, Nc, i, d;

3061:     PetscDSGetDiscretization(prob, field, &obj);
3062:     PetscObjectGetClassId(obj, &id);
3063:     if (id == PETSCFE_CLASSID) {
3064:       PetscFEGetQuadrature((PetscFE) obj, &quad);
3065:       PetscFEGetCellTabulation((PetscFE) obj, 1, &Tfine);
3066:       PetscFECreateTabulation((PetscFE) obj, 1, 1, x, 0, &T);
3067:     } else {
3068:       PetscFVGetQuadrature((PetscFV) obj, &quad);
3069:     }
3070:     PetscQuadratureGetData(quad, NULL, &Nc, &Nq, &qpoints, &qweights);
3071:     /* For each fine grid cell */
3072:     for (cell = cStart; cell < cEnd; ++cell) {
3073:       Vec                pointVec;
3074:       PetscScalar       *pV;
3075:       PetscSF            coarseCellSF = NULL;
3076:       const PetscSFNode *coarseCells;
3077:       PetscInt           numCoarseCells, cpdim, q, c, j;
3078:       PetscInt          *findices,   *cindices;
3079:       PetscInt           numFIndices, numCIndices;

3081:       DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);
3082:       DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);
3083:       /* Get points from the quadrature */
3084:       VecCreateSeq(PETSC_COMM_SELF, Nq*dim, &pointVec);
3085:       VecSetBlockSize(pointVec, dim);
3086:       VecGetArray(pointVec, &pV);
3087:       for (q = 0; q < Nq; ++q) {
3088:         const PetscReal xi0[3] = {-1., -1., -1.};

3090:         /* Transform point to real space */
3091:         CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q*dim], x);
3092:         for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
3093:       }
3094:       VecRestoreArray(pointVec, &pV);
3095:       /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
3096:       DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);
3097:       /* Update matrix */
3098:       PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);
3100:       VecGetArray(pointVec, &pV);
3101:       for (ccell = 0; ccell < numCoarseCells; ++ccell) {
3102:         PetscReal pVReal[3];
3103:         const PetscReal xi0[3] = {-1., -1., -1.};

3105:         DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL);
3106:         /* Transform points from real space to coarse reference space */
3107:         DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc);
3108:         for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell*dim+d]);
3109:         CoordinatesRealToRef(dim, dim, xi0, v0c, invJc, pVReal, x);

3111:         if (id == PETSCFE_CLASSID) {
3112:           PetscFE fe = (PetscFE) obj;

3114:           /* Evaluate coarse basis on contained point */
3115:           PetscFEGetDimension(fe, &cpdim);
3116:           PetscFEComputeTabulation(fe, 1, x, 0, T);
3117:           /* Get elemMat entries by multiplying by weight */
3118:           for (i = 0; i < numFIndices; ++i) {
3119:             PetscArrayzero(elemMat, cpdim);
3120:             for (j = 0; j < cpdim; ++j) {
3121:               for (c = 0; c < Nc; ++c) elemMat[j] += T->T[0][j*Nc + c]*Tfine->T[0][(ccell*numFIndices + i)*Nc + c]*qweights[ccell*Nc + c]*detJ;
3122:             }
3123:             /* Update interpolator */
3124:             if (mesh->printFEM > 1) DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);
3126:             MatSetValues(mass, 1, &findices[i], numCIndices, cindices, elemMat, ADD_VALUES);
3127:           }
3128:         } else {
3129:           cpdim = 1;
3130:           for (i = 0; i < numFIndices; ++i) {
3131:             PetscArrayzero(elemMat, cpdim);
3132:             for (j = 0; j < cpdim; ++j) {
3133:               for (c = 0; c < Nc; ++c) elemMat[j] += 1.0*1.0*qweights[ccell*Nc + c]*detJ;
3134:             }
3135:             /* Update interpolator */
3136:             if (mesh->printFEM > 1) DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);
3137:             PetscPrintf(PETSC_COMM_SELF, "Nq: %D %D Nf: %D %D Nc: %D %D\n", ccell, Nq, i, numFIndices, j, numCIndices);
3139:             MatSetValues(mass, 1, &findices[i], numCIndices, cindices, elemMat, ADD_VALUES);
3140:           }
3141:         }
3142:         DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL);
3143:       }
3144:       VecRestoreArray(pointVec, &pV);
3145:       PetscSFDestroy(&coarseCellSF);
3146:       VecDestroy(&pointVec);
3147:       DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);
3148:     }
3149:     if (id == PETSCFE_CLASSID) PetscTabulationDestroy(&T);
3150:   }
3151:   PetscFree3(v0,J,invJ);
3152:   PetscFree3(v0c,Jc,invJc);
3153:   PetscFree(elemMat);
3154:   MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY);
3155:   MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY);
3156:   return 0;
3157: }

3159: /*@
3160:   DMPlexComputeInjectorFEM - Compute a mapping from coarse unknowns to fine unknowns

3162:   Input Parameters:
3163: + dmc  - The coarse mesh
3164: - dmf  - The fine mesh
3165: - user - The user context

3167:   Output Parameter:
3168: . sc   - The mapping

3170:   Level: developer

3172: .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM()
3173: @*/
3174: PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
3175: {
3176:   PetscDS        prob;
3177:   PetscFE       *feRef;
3178:   PetscFV       *fvRef;
3179:   Vec            fv, cv;
3180:   IS             fis, cis;
3181:   PetscSection   fsection, fglobalSection, csection, cglobalSection;
3182:   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
3183:   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, c, dim, d, startC, endC, offsetC, offsetF, m;
3184:   PetscBool     *needAvg;

3186:   PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);
3187:   DMGetDimension(dmf, &dim);
3188:   DMGetLocalSection(dmf, &fsection);
3189:   DMGetGlobalSection(dmf, &fglobalSection);
3190:   DMGetLocalSection(dmc, &csection);
3191:   DMGetGlobalSection(dmc, &cglobalSection);
3192:   PetscSectionGetNumFields(fsection, &Nf);
3193:   DMPlexGetSimplexOrBoxCells(dmc, 0, &cStart, &cEnd);
3194:   DMGetDS(dmc, &prob);
3195:   PetscCalloc3(Nf,&feRef,Nf,&fvRef,Nf,&needAvg);
3196:   for (f = 0; f < Nf; ++f) {
3197:     PetscObject  obj;
3198:     PetscClassId id;
3199:     PetscInt     fNb = 0, Nc = 0;

3201:     PetscDSGetDiscretization(prob, f, &obj);
3202:     PetscObjectGetClassId(obj, &id);
3203:     if (id == PETSCFE_CLASSID) {
3204:       PetscFE    fe = (PetscFE) obj;
3205:       PetscSpace sp;
3206:       PetscInt   maxDegree;

3208:       PetscFERefine(fe, &feRef[f]);
3209:       PetscFEGetDimension(feRef[f], &fNb);
3210:       PetscFEGetNumComponents(fe, &Nc);
3211:       PetscFEGetBasisSpace(fe, &sp);
3212:       PetscSpaceGetDegree(sp, NULL, &maxDegree);
3213:       if (!maxDegree) needAvg[f] = PETSC_TRUE;
3214:     } else if (id == PETSCFV_CLASSID) {
3215:       PetscFV        fv = (PetscFV) obj;
3216:       PetscDualSpace Q;

3218:       PetscFVRefine(fv, &fvRef[f]);
3219:       PetscFVGetDualSpace(fvRef[f], &Q);
3220:       PetscDualSpaceGetDimension(Q, &fNb);
3221:       PetscFVGetNumComponents(fv, &Nc);
3222:       needAvg[f] = PETSC_TRUE;
3223:     }
3224:     fTotDim += fNb;
3225:   }
3226:   PetscDSGetTotalDimension(prob, &cTotDim);
3227:   PetscMalloc1(cTotDim,&cmap);
3228:   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
3229:     PetscFE        feC;
3230:     PetscFV        fvC;
3231:     PetscDualSpace QF, QC;
3232:     PetscInt       order = -1, NcF, NcC, fpdim, cpdim;

3234:     if (feRef[field]) {
3235:       PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);
3236:       PetscFEGetNumComponents(feC, &NcC);
3237:       PetscFEGetNumComponents(feRef[field], &NcF);
3238:       PetscFEGetDualSpace(feRef[field], &QF);
3239:       PetscDualSpaceGetOrder(QF, &order);
3240:       PetscDualSpaceGetDimension(QF, &fpdim);
3241:       PetscFEGetDualSpace(feC, &QC);
3242:       PetscDualSpaceGetDimension(QC, &cpdim);
3243:     } else {
3244:       PetscDSGetDiscretization(prob, field, (PetscObject *) &fvC);
3245:       PetscFVGetNumComponents(fvC, &NcC);
3246:       PetscFVGetNumComponents(fvRef[field], &NcF);
3247:       PetscFVGetDualSpace(fvRef[field], &QF);
3248:       PetscDualSpaceGetDimension(QF, &fpdim);
3249:       PetscFVGetDualSpace(fvC, &QC);
3250:       PetscDualSpaceGetDimension(QC, &cpdim);
3251:     }
3253:     for (c = 0; c < cpdim; ++c) {
3254:       PetscQuadrature  cfunc;
3255:       const PetscReal *cqpoints, *cqweights;
3256:       PetscInt         NqcC, NpC;
3257:       PetscBool        found = PETSC_FALSE;

3259:       PetscDualSpaceGetFunctional(QC, c, &cfunc);
3260:       PetscQuadratureGetData(cfunc, NULL, &NqcC, &NpC, &cqpoints, &cqweights);
3263:       for (f = 0; f < fpdim; ++f) {
3264:         PetscQuadrature  ffunc;
3265:         const PetscReal *fqpoints, *fqweights;
3266:         PetscReal        sum = 0.0;
3267:         PetscInt         NqcF, NpF;

3269:         PetscDualSpaceGetFunctional(QF, f, &ffunc);
3270:         PetscQuadratureGetData(ffunc, NULL, &NqcF, &NpF, &fqpoints, &fqweights);
3272:         if (NpC != NpF) continue;
3273:         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
3274:         if (sum > 1.0e-9) continue;
3275:         for (d = 0; d < NcC; ++d) sum += PetscAbsReal(cqweights[d]*fqweights[d]);
3276:         if (sum < 1.0e-9) continue;
3277:         cmap[offsetC+c] = offsetF+f;
3278:         found = PETSC_TRUE;
3279:         break;
3280:       }
3281:       if (!found) {
3282:         /* TODO We really want the average here, but some asshole put VecScatter in the interface */
3283:         if (fvRef[field] || (feRef[field] && order == 0)) {
3284:           cmap[offsetC+c] = offsetF+0;
3285:         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate matching functional for injection");
3286:       }
3287:     }
3288:     offsetC += cpdim;
3289:     offsetF += fpdim;
3290:   }
3291:   for (f = 0; f < Nf; ++f) {PetscFEDestroy(&feRef[f]));PetscCall(PetscFVDestroy(&fvRef[f]);}
3292:   PetscFree3(feRef,fvRef,needAvg);

3294:   DMGetGlobalVector(dmf, &fv);
3295:   DMGetGlobalVector(dmc, &cv);
3296:   VecGetOwnershipRange(cv, &startC, &endC);
3297:   PetscSectionGetConstrainedStorageSize(cglobalSection, &m);
3298:   PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);
3299:   PetscMalloc1(m,&cindices);
3300:   PetscMalloc1(m,&findices);
3301:   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
3302:   for (c = cStart; c < cEnd; ++c) {
3303:     DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);
3304:     for (d = 0; d < cTotDim; ++d) {
3305:       if ((cellCIndices[d] < startC) || (cellCIndices[d] >= endC)) continue;
3307:       cindices[cellCIndices[d]-startC] = cellCIndices[d];
3308:       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
3309:     }
3310:   }
3311:   PetscFree(cmap);
3312:   PetscFree2(cellCIndices,cellFIndices);

3314:   ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);
3315:   ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);
3316:   VecScatterCreate(cv, cis, fv, fis, sc);
3317:   ISDestroy(&cis);
3318:   ISDestroy(&fis);
3319:   DMRestoreGlobalVector(dmf, &fv);
3320:   DMRestoreGlobalVector(dmc, &cv);
3321:   PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);
3322:   return 0;
3323: }

3325: /*@C
3326:   DMPlexGetCellFields - Retrieve the field values values for a chunk of cells

3328:   Input Parameters:
3329: + dm     - The DM
3330: . cellIS - The cells to include
3331: . locX   - A local vector with the solution fields
3332: . locX_t - A local vector with solution field time derivatives, or NULL
3333: - locA   - A local vector with auxiliary fields, or NULL

3335:   Output Parameters:
3336: + u   - The field coefficients
3337: . u_t - The fields derivative coefficients
3338: - a   - The auxiliary field coefficients

3340:   Level: developer

3342: .seealso: DMPlexGetFaceFields()
3343: @*/
3344: PetscErrorCode DMPlexGetCellFields(DM dm, IS cellIS, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a)
3345: {
3346:   DM              plex, plexA = NULL;
3347:   DMEnclosureType encAux;
3348:   PetscSection    section, sectionAux;
3349:   PetscDS         prob;
3350:   const PetscInt *cells;
3351:   PetscInt        cStart, cEnd, numCells, totDim, totDimAux, c;

3360:   DMPlexConvertPlex(dm, &plex, PETSC_FALSE);
3361:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
3362:   DMGetLocalSection(dm, &section);
3363:   DMGetCellDS(dm, cells ? cells[cStart] : cStart, &prob);
3364:   PetscDSGetTotalDimension(prob, &totDim);
3365:   if (locA) {
3366:     DM      dmAux;
3367:     PetscDS probAux;

3369:     VecGetDM(locA, &dmAux);
3370:     DMGetEnclosureRelation(dmAux, dm, &encAux);
3371:     DMPlexConvertPlex(dmAux, &plexA, PETSC_FALSE);
3372:     DMGetLocalSection(dmAux, &sectionAux);
3373:     DMGetDS(dmAux, &probAux);
3374:     PetscDSGetTotalDimension(probAux, &totDimAux);
3375:   }
3376:   numCells = cEnd - cStart;
3377:   DMGetWorkArray(dm, numCells*totDim, MPIU_SCALAR, u);
3378:   if (locX_t) DMGetWorkArray(dm, numCells*totDim, MPIU_SCALAR, u_t); else {*u_t = NULL;}
3379:   if (locA)   DMGetWorkArray(dm, numCells*totDimAux, MPIU_SCALAR, a); else {*a = NULL;}
3380:   for (c = cStart; c < cEnd; ++c) {
3381:     const PetscInt cell = cells ? cells[c] : c;
3382:     const PetscInt cind = c - cStart;
3383:     PetscScalar   *x = NULL, *x_t = NULL, *ul = *u, *ul_t = *u_t, *al = *a;
3384:     PetscInt       i;

3386:     DMPlexVecGetClosure(plex, section, locX, cell, NULL, &x);
3387:     for (i = 0; i < totDim; ++i) ul[cind*totDim+i] = x[i];
3388:     DMPlexVecRestoreClosure(plex, section, locX, cell, NULL, &x);
3389:     if (locX_t) {
3390:       DMPlexVecGetClosure(plex, section, locX_t, cell, NULL, &x_t);
3391:       for (i = 0; i < totDim; ++i) ul_t[cind*totDim+i] = x_t[i];
3392:       DMPlexVecRestoreClosure(plex, section, locX_t, cell, NULL, &x_t);
3393:     }
3394:     if (locA) {
3395:       PetscInt subcell;
3396:       DMGetEnclosurePoint(plexA, dm, encAux, cell, &subcell);
3397:       DMPlexVecGetClosure(plexA, sectionAux, locA, subcell, NULL, &x);
3398:       for (i = 0; i < totDimAux; ++i) al[cind*totDimAux+i] = x[i];
3399:       DMPlexVecRestoreClosure(plexA, sectionAux, locA, subcell, NULL, &x);
3400:     }
3401:   }
3402:   DMDestroy(&plex);
3403:   if (locA) DMDestroy(&plexA);
3404:   ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);
3405:   return 0;
3406: }

3408: /*@C
3409:   DMPlexRestoreCellFields - Restore the field values values for a chunk of cells

3411:   Input Parameters:
3412: + dm     - The DM
3413: . cellIS - The cells to include
3414: . locX   - A local vector with the solution fields
3415: . locX_t - A local vector with solution field time derivatives, or NULL
3416: - locA   - A local vector with auxiliary fields, or NULL

3418:   Output Parameters:
3419: + u   - The field coefficients
3420: . u_t - The fields derivative coefficients
3421: - a   - The auxiliary field coefficients

3423:   Level: developer

3425: .seealso: DMPlexGetFaceFields()
3426: @*/
3427: PetscErrorCode DMPlexRestoreCellFields(DM dm, IS cellIS, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a)
3428: {
3429:   DMRestoreWorkArray(dm, 0, MPIU_SCALAR, u);
3430:   if (locX_t) DMRestoreWorkArray(dm, 0, MPIU_SCALAR, u_t);
3431:   if (locA)   DMRestoreWorkArray(dm, 0, MPIU_SCALAR, a);
3432:   return 0;
3433: }

3435: /*
3436:   Get the auxiliary field vectors for the negative side (s = 0) and positive side (s = 1) of the interfaace
3437: */
3438: static PetscErrorCode DMPlexGetHybridAuxFields(DM dm, DM dmAux[], PetscDS dsAux[], IS cellIS, Vec locA[], PetscScalar *a[])
3439: {
3440:   DM              plexA[2];
3441:   DMEnclosureType encAux[2];
3442:   PetscSection    sectionAux[2];
3443:   const PetscInt *cells;
3444:   PetscInt        cStart, cEnd, numCells, c, s, totDimAux[2];

3447:   if (!locA[0] || !locA[1]) return 0;
3451:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
3452:   numCells = cEnd - cStart;
3453:   for (s = 0; s < 2; ++s) {
3457:     DMPlexConvertPlex(dmAux[s], &plexA[s], PETSC_FALSE);
3458:     DMGetEnclosureRelation(dmAux[s], dm, &encAux[s]);
3459:     DMGetLocalSection(dmAux[s], &sectionAux[s]);
3460:     PetscDSGetTotalDimension(dsAux[s], &totDimAux[s]);
3461:     DMGetWorkArray(dmAux[s], numCells*totDimAux[s], MPIU_SCALAR, &a[s]);
3462:   }
3463:   for (c = cStart; c < cEnd; ++c) {
3464:     const PetscInt  cell = cells ? cells[c] : c;
3465:     const PetscInt  cind = c - cStart;
3466:     const PetscInt *cone, *ornt;

3468:     DMPlexGetCone(dm, cell, &cone);
3469:     DMPlexGetConeOrientation(dm, cell, &ornt);
3471:     for (s = 0; s < 2; ++s) {
3472:       const PetscInt *support;
3473:       PetscScalar    *x = NULL, *al = a[s];
3474:       const PetscInt  tdA = totDimAux[s];
3475:       PetscInt        ssize, scell;
3476:       PetscInt        subface, Na, i;

3478:       DMPlexGetSupport(dm, cone[s], &support);
3479:       DMPlexGetSupportSize(dm, cone[s], &ssize);
3481:       if      (support[0] == cell) scell = support[1];
3482:       else if (support[1] == cell) scell = support[0];
3483:       else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %D does not have cell %D in its support", cone[s], cell);

3485:       DMGetEnclosurePoint(plexA[s], dm, encAux[s], scell, &subface);
3486:       DMPlexVecGetClosure(plexA[s], sectionAux[s], locA[s], subface, &Na, &x);
3487:       for (i = 0; i < Na; ++i) al[cind*tdA+i] = x[i];
3488:       DMPlexVecRestoreClosure(plexA[s], sectionAux[s], locA[s], subface, &Na, &x);
3489:     }
3490:   }
3491:   for (s = 0; s < 2; ++s) DMDestroy(&plexA[s]);
3492:   ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);
3493:   return 0;
3494: }

3496: static PetscErrorCode DMPlexRestoreHybridAuxFields(DM dmAux[], PetscDS dsAux[], IS cellIS, Vec locA[], PetscScalar *a[])
3497: {
3498:   if (!locA[0] || !locA[1]) return 0;
3499:   DMRestoreWorkArray(dmAux[0], 0, MPIU_SCALAR, &a[0]);
3500:   DMRestoreWorkArray(dmAux[1], 0, MPIU_SCALAR, &a[1]);
3501:   return 0;
3502: }

3504: /*@C
3505:   DMPlexGetFaceFields - Retrieve the field values values for a chunk of faces

3507:   Input Parameters:
3508: + dm     - The DM
3509: . fStart - The first face to include
3510: . fEnd   - The first face to exclude
3511: . locX   - A local vector with the solution fields
3512: . locX_t - A local vector with solution field time derivatives, or NULL
3513: . faceGeometry - A local vector with face geometry
3514: . cellGeometry - A local vector with cell geometry
3515: - locaGrad - A local vector with field gradients, or NULL

3517:   Output Parameters:
3518: + Nface - The number of faces with field values
3519: . uL - The field values at the left side of the face
3520: - uR - The field values at the right side of the face

3522:   Level: developer

3524: .seealso: DMPlexGetCellFields()
3525: @*/
3526: PetscErrorCode DMPlexGetFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, Vec locX_t, Vec faceGeometry, Vec cellGeometry, Vec locGrad, PetscInt *Nface, PetscScalar **uL, PetscScalar **uR)
3527: {
3528:   DM                 dmFace, dmCell, dmGrad = NULL;
3529:   PetscSection       section;
3530:   PetscDS            prob;
3531:   DMLabel            ghostLabel;
3532:   const PetscScalar *facegeom, *cellgeom, *x, *lgrad;
3533:   PetscBool         *isFE;
3534:   PetscInt           dim, Nf, f, Nc, numFaces = fEnd - fStart, iface, face;

3544:   DMGetDimension(dm, &dim);
3545:   DMGetDS(dm, &prob);
3546:   DMGetLocalSection(dm, &section);
3547:   PetscDSGetNumFields(prob, &Nf);
3548:   PetscDSGetTotalComponents(prob, &Nc);
3549:   PetscMalloc1(Nf, &isFE);
3550:   for (f = 0; f < Nf; ++f) {
3551:     PetscObject  obj;
3552:     PetscClassId id;

3554:     PetscDSGetDiscretization(prob, f, &obj);
3555:     PetscObjectGetClassId(obj, &id);
3556:     if (id == PETSCFE_CLASSID)      {isFE[f] = PETSC_TRUE;}
3557:     else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
3558:     else                            {isFE[f] = PETSC_FALSE;}
3559:   }
3560:   DMGetLabel(dm, "ghost", &ghostLabel);
3561:   VecGetArrayRead(locX, &x);
3562:   VecGetDM(faceGeometry, &dmFace);
3563:   VecGetArrayRead(faceGeometry, &facegeom);
3564:   VecGetDM(cellGeometry, &dmCell);
3565:   VecGetArrayRead(cellGeometry, &cellgeom);
3566:   if (locGrad) {
3567:     VecGetDM(locGrad, &dmGrad);
3568:     VecGetArrayRead(locGrad, &lgrad);
3569:   }
3570:   DMGetWorkArray(dm, numFaces*Nc, MPIU_SCALAR, uL);
3571:   DMGetWorkArray(dm, numFaces*Nc, MPIU_SCALAR, uR);
3572:   /* Right now just eat the extra work for FE (could make a cell loop) */
3573:   for (face = fStart, iface = 0; face < fEnd; ++face) {
3574:     const PetscInt        *cells;
3575:     PetscFVFaceGeom       *fg;
3576:     PetscFVCellGeom       *cgL, *cgR;
3577:     PetscScalar           *xL, *xR, *gL, *gR;
3578:     PetscScalar           *uLl = *uL, *uRl = *uR;
3579:     PetscInt               ghost, nsupp, nchild;

3581:     DMLabelGetValue(ghostLabel, face, &ghost);
3582:     DMPlexGetSupportSize(dm, face, &nsupp);
3583:     DMPlexGetTreeChildren(dm, face, &nchild, NULL);
3584:     if (ghost >= 0 || nsupp > 2 || nchild > 0) continue;
3585:     DMPlexPointLocalRead(dmFace, face, facegeom, &fg);
3586:     DMPlexGetSupport(dm, face, &cells);
3587:     DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);
3588:     DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);
3589:     for (f = 0; f < Nf; ++f) {
3590:       PetscInt off;

3592:       PetscDSGetComponentOffset(prob, f, &off);
3593:       if (isFE[f]) {
3594:         const PetscInt *cone;
3595:         PetscInt        comp, coneSizeL, coneSizeR, faceLocL, faceLocR, ldof, rdof, d;

3597:         xL = xR = NULL;
3598:         PetscSectionGetFieldComponents(section, f, &comp);
3599:         DMPlexVecGetClosure(dm, section, locX, cells[0], &ldof, (PetscScalar **) &xL);
3600:         DMPlexVecGetClosure(dm, section, locX, cells[1], &rdof, (PetscScalar **) &xR);
3601:         DMPlexGetCone(dm, cells[0], &cone);
3602:         DMPlexGetConeSize(dm, cells[0], &coneSizeL);
3603:         for (faceLocL = 0; faceLocL < coneSizeL; ++faceLocL) if (cone[faceLocL] == face) break;
3604:         DMPlexGetCone(dm, cells[1], &cone);
3605:         DMPlexGetConeSize(dm, cells[1], &coneSizeR);
3606:         for (faceLocR = 0; faceLocR < coneSizeR; ++faceLocR) if (cone[faceLocR] == face) break;
3608:         /* Check that FEM field has values in the right cell (sometimes its an FV ghost cell) */
3609:         /* TODO: this is a hack that might not be right for nonconforming */
3610:         if (faceLocL < coneSizeL) {
3611:           PetscFEEvaluateFaceFields_Internal(prob, f, faceLocL, xL, &uLl[iface*Nc+off]);
3612:           if (rdof == ldof && faceLocR < coneSizeR) PetscFEEvaluateFaceFields_Internal(prob, f, faceLocR, xR, &uRl[iface*Nc+off]);
3613:           else              {for (d = 0; d < comp; ++d) uRl[iface*Nc+off+d] = uLl[iface*Nc+off+d];}
3614:         }
3615:         else {
3616:           PetscFEEvaluateFaceFields_Internal(prob, f, faceLocR, xR, &uRl[iface*Nc+off]);
3617:           PetscSectionGetFieldComponents(section, f, &comp);
3618:           for (d = 0; d < comp; ++d) uLl[iface*Nc+off+d] = uRl[iface*Nc+off+d];
3619:         }
3620:         DMPlexVecRestoreClosure(dm, section, locX, cells[0], &ldof, (PetscScalar **) &xL);
3621:         DMPlexVecRestoreClosure(dm, section, locX, cells[1], &rdof, (PetscScalar **) &xR);
3622:       } else {
3623:         PetscFV  fv;
3624:         PetscInt numComp, c;

3626:         PetscDSGetDiscretization(prob, f, (PetscObject *) &fv);
3627:         PetscFVGetNumComponents(fv, &numComp);
3628:         DMPlexPointLocalFieldRead(dm, cells[0], f, x, &xL);
3629:         DMPlexPointLocalFieldRead(dm, cells[1], f, x, &xR);
3630:         if (dmGrad) {
3631:           PetscReal dxL[3], dxR[3];

3633:           DMPlexPointLocalRead(dmGrad, cells[0], lgrad, &gL);
3634:           DMPlexPointLocalRead(dmGrad, cells[1], lgrad, &gR);
3635:           DMPlex_WaxpyD_Internal(dim, -1, cgL->centroid, fg->centroid, dxL);
3636:           DMPlex_WaxpyD_Internal(dim, -1, cgR->centroid, fg->centroid, dxR);
3637:           for (c = 0; c < numComp; ++c) {
3638:             uLl[iface*Nc+off+c] = xL[c] + DMPlex_DotD_Internal(dim, &gL[c*dim], dxL);
3639:             uRl[iface*Nc+off+c] = xR[c] + DMPlex_DotD_Internal(dim, &gR[c*dim], dxR);
3640:           }
3641:         } else {
3642:           for (c = 0; c < numComp; ++c) {
3643:             uLl[iface*Nc+off+c] = xL[c];
3644:             uRl[iface*Nc+off+c] = xR[c];
3645:           }
3646:         }
3647:       }
3648:     }
3649:     ++iface;
3650:   }
3651:   *Nface = iface;
3652:   VecRestoreArrayRead(locX, &x);
3653:   VecRestoreArrayRead(faceGeometry, &facegeom);
3654:   VecRestoreArrayRead(cellGeometry, &cellgeom);
3655:   if (locGrad) {
3656:     VecRestoreArrayRead(locGrad, &lgrad);
3657:   }
3658:   PetscFree(isFE);
3659:   return 0;
3660: }

3662: /*@C
3663:   DMPlexRestoreFaceFields - Restore the field values values for a chunk of faces

3665:   Input Parameters:
3666: + dm     - The DM
3667: . fStart - The first face to include
3668: . fEnd   - The first face to exclude
3669: . locX   - A local vector with the solution fields
3670: . locX_t - A local vector with solution field time derivatives, or NULL
3671: . faceGeometry - A local vector with face geometry
3672: . cellGeometry - A local vector with cell geometry
3673: - locaGrad - A local vector with field gradients, or NULL

3675:   Output Parameters:
3676: + Nface - The number of faces with field values
3677: . uL - The field values at the left side of the face
3678: - uR - The field values at the right side of the face

3680:   Level: developer

3682: .seealso: DMPlexGetFaceFields()
3683: @*/
3684: PetscErrorCode DMPlexRestoreFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, Vec locX_t, Vec faceGeometry, Vec cellGeometry, Vec locGrad, PetscInt *Nface, PetscScalar **uL, PetscScalar **uR)
3685: {
3686:   DMRestoreWorkArray(dm, 0, MPIU_SCALAR, uL);
3687:   DMRestoreWorkArray(dm, 0, MPIU_SCALAR, uR);
3688:   return 0;
3689: }

3691: /*@C
3692:   DMPlexGetFaceGeometry - Retrieve the geometric values for a chunk of faces

3694:   Input Parameters:
3695: + dm     - The DM
3696: . fStart - The first face to include
3697: . fEnd   - The first face to exclude
3698: . faceGeometry - A local vector with face geometry
3699: - cellGeometry - A local vector with cell geometry

3701:   Output Parameters:
3702: + Nface - The number of faces with field values
3703: . fgeom - The extract the face centroid and normal
3704: - vol   - The cell volume

3706:   Level: developer

3708: .seealso: DMPlexGetCellFields()
3709: @*/
3710: PetscErrorCode DMPlexGetFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscInt *Nface, PetscFVFaceGeom **fgeom, PetscReal **vol)
3711: {
3712:   DM                 dmFace, dmCell;
3713:   DMLabel            ghostLabel;
3714:   const PetscScalar *facegeom, *cellgeom;
3715:   PetscInt           dim, numFaces = fEnd - fStart, iface, face;

3722:   DMGetDimension(dm, &dim);
3723:   DMGetLabel(dm, "ghost", &ghostLabel);
3724:   VecGetDM(faceGeometry, &dmFace);
3725:   VecGetArrayRead(faceGeometry, &facegeom);
3726:   VecGetDM(cellGeometry, &dmCell);
3727:   VecGetArrayRead(cellGeometry, &cellgeom);
3728:   PetscMalloc1(numFaces, fgeom);
3729:   DMGetWorkArray(dm, numFaces*2, MPIU_SCALAR, vol);
3730:   for (face = fStart, iface = 0; face < fEnd; ++face) {
3731:     const PetscInt        *cells;
3732:     PetscFVFaceGeom       *fg;
3733:     PetscFVCellGeom       *cgL, *cgR;
3734:     PetscFVFaceGeom       *fgeoml = *fgeom;
3735:     PetscReal             *voll   = *vol;
3736:     PetscInt               ghost, d, nchild, nsupp;

3738:     DMLabelGetValue(ghostLabel, face, &ghost);
3739:     DMPlexGetSupportSize(dm, face, &nsupp);
3740:     DMPlexGetTreeChildren(dm, face, &nchild, NULL);
3741:     if (ghost >= 0 || nsupp > 2 || nchild > 0) continue;
3742:     DMPlexPointLocalRead(dmFace, face, facegeom, &fg);
3743:     DMPlexGetSupport(dm, face, &cells);
3744:     DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);
3745:     DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);
3746:     for (d = 0; d < dim; ++d) {
3747:       fgeoml[iface].centroid[d] = fg->centroid[d];
3748:       fgeoml[iface].normal[d]   = fg->normal[d];
3749:     }
3750:     voll[iface*2+0] = cgL->volume;
3751:     voll[iface*2+1] = cgR->volume;
3752:     ++iface;
3753:   }
3754:   *Nface = iface;
3755:   VecRestoreArrayRead(faceGeometry, &facegeom);
3756:   VecRestoreArrayRead(cellGeometry, &cellgeom);
3757:   return 0;
3758: }

3760: /*@C
3761:   DMPlexRestoreFaceGeometry - Restore the field values values for a chunk of faces

3763:   Input Parameters:
3764: + dm     - The DM
3765: . fStart - The first face to include
3766: . fEnd   - The first face to exclude
3767: . faceGeometry - A local vector with face geometry
3768: - cellGeometry - A local vector with cell geometry

3770:   Output Parameters:
3771: + Nface - The number of faces with field values
3772: . fgeom - The extract the face centroid and normal
3773: - vol   - The cell volume

3775:   Level: developer

3777: .seealso: DMPlexGetFaceFields()
3778: @*/
3779: PetscErrorCode DMPlexRestoreFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscInt *Nface, PetscFVFaceGeom **fgeom, PetscReal **vol)
3780: {
3781:   PetscFree(*fgeom);
3782:   DMRestoreWorkArray(dm, 0, MPIU_REAL, vol);
3783:   return 0;
3784: }

3786: PetscErrorCode DMSNESGetFEGeom(DMField coordField, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
3787: {
3788:   char            composeStr[33] = {0};
3789:   PetscObjectId   id;
3790:   PetscContainer  container;

3792:   PetscObjectGetId((PetscObject)quad,&id);
3793:   PetscSNPrintf(composeStr, 32, "DMSNESGetFEGeom_%x\n", id);
3794:   PetscObjectQuery((PetscObject) pointIS, composeStr, (PetscObject *) &container);
3795:   if (container) {
3796:     PetscContainerGetPointer(container, (void **) geom);
3797:   } else {
3798:     DMFieldCreateFEGeom(coordField, pointIS, quad, faceData, geom);
3799:     PetscContainerCreate(PETSC_COMM_SELF,&container);
3800:     PetscContainerSetPointer(container, (void *) *geom);
3801:     PetscContainerSetUserDestroy(container, PetscContainerUserDestroy_PetscFEGeom);
3802:     PetscObjectCompose((PetscObject) pointIS, composeStr, (PetscObject) container);
3803:     PetscContainerDestroy(&container);
3804:   }
3805:   return 0;
3806: }

3808: PetscErrorCode DMSNESRestoreFEGeom(DMField coordField, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
3809: {
3810:   *geom = NULL;
3811:   return 0;
3812: }

3814: PetscErrorCode DMPlexComputeResidual_Patch_Internal(DM dm, PetscSection section, IS cellIS, PetscReal t, Vec locX, Vec locX_t, Vec locF, void *user)
3815: {
3816:   DM_Plex         *mesh       = (DM_Plex *) dm->data;
3817:   const char      *name       = "Residual";
3818:   DM               dmAux      = NULL;
3819:   DMLabel          ghostLabel = NULL;
3820:   PetscDS          prob       = NULL;
3821:   PetscDS          probAux    = NULL;
3822:   PetscBool        useFEM     = PETSC_FALSE;
3823:   PetscBool        isImplicit = (locX_t || t == PETSC_MIN_REAL) ? PETSC_TRUE : PETSC_FALSE;
3824:   DMField          coordField = NULL;
3825:   Vec              locA;
3826:   PetscScalar     *u = NULL, *u_t, *a, *uL = NULL, *uR = NULL;
3827:   IS               chunkIS;
3828:   const PetscInt  *cells;
3829:   PetscInt         cStart, cEnd, numCells;
3830:   PetscInt         Nf, f, totDim, totDimAux, numChunks, cellChunkSize, chunk, fStart, fEnd;
3831:   PetscInt         maxDegree = PETSC_MAX_INT;
3832:   PetscFormKey key;
3833:   PetscQuadrature  affineQuad = NULL, *quads = NULL;
3834:   PetscFEGeom     *affineGeom = NULL, **geoms = NULL;

3836:   PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);
3837:   /* FEM+FVM */
3838:   /* 1: Get sizes from dm and dmAux */
3839:   DMGetLabel(dm, "ghost", &ghostLabel);
3840:   DMGetDS(dm, &prob);
3841:   PetscDSGetNumFields(prob, &Nf);
3842:   PetscDSGetTotalDimension(prob, &totDim);
3843:   DMGetAuxiliaryVec(dm, NULL, 0, 0, &locA);
3844:   if (locA) {
3845:     VecGetDM(locA, &dmAux);
3846:     DMGetDS(dmAux, &probAux);
3847:     PetscDSGetTotalDimension(probAux, &totDimAux);
3848:   }
3849:   /* 2: Get geometric data */
3850:   for (f = 0; f < Nf; ++f) {
3851:     PetscObject  obj;
3852:     PetscClassId id;
3853:     PetscBool    fimp;

3855:     PetscDSGetImplicit(prob, f, &fimp);
3856:     if (isImplicit != fimp) continue;
3857:     PetscDSGetDiscretization(prob, f, &obj);
3858:     PetscObjectGetClassId(obj, &id);
3859:     if (id == PETSCFE_CLASSID) {useFEM = PETSC_TRUE;}
3861:   }
3862:   if (useFEM) {
3863:     DMGetCoordinateField(dm, &coordField);
3864:     DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);
3865:     if (maxDegree <= 1) {
3866:       DMFieldCreateDefaultQuadrature(coordField,cellIS,&affineQuad);
3867:       if (affineQuad) {
3868:         DMSNESGetFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);
3869:       }
3870:     } else {
3871:       PetscCalloc2(Nf,&quads,Nf,&geoms);
3872:       for (f = 0; f < Nf; ++f) {
3873:         PetscObject  obj;
3874:         PetscClassId id;
3875:         PetscBool    fimp;

3877:         PetscDSGetImplicit(prob, f, &fimp);
3878:         if (isImplicit != fimp) continue;
3879:         PetscDSGetDiscretization(prob, f, &obj);
3880:         PetscObjectGetClassId(obj, &id);
3881:         if (id == PETSCFE_CLASSID) {
3882:           PetscFE fe = (PetscFE) obj;

3884:           PetscFEGetQuadrature(fe, &quads[f]);
3885:           PetscObjectReference((PetscObject)quads[f]);
3886:           DMSNESGetFEGeom(coordField,cellIS,quads[f],PETSC_FALSE,&geoms[f]);
3887:         }
3888:       }
3889:     }
3890:   }
3891:   /* Loop over chunks */
3892:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
3893:   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
3894:   if (useFEM) ISCreate(PETSC_COMM_SELF, &chunkIS);
3895:   numCells      = cEnd - cStart;
3896:   numChunks     = 1;
3897:   cellChunkSize = numCells/numChunks;
3898:   numChunks     = PetscMin(1,numCells);
3899:   key.label     = NULL;
3900:   key.value     = 0;
3901:   key.part      = 0;
3902:   for (chunk = 0; chunk < numChunks; ++chunk) {
3903:     PetscScalar     *elemVec, *fluxL = NULL, *fluxR = NULL;
3904:     PetscReal       *vol = NULL;
3905:     PetscFVFaceGeom *fgeom = NULL;
3906:     PetscInt         cS = cStart+chunk*cellChunkSize, cE = PetscMin(cS+cellChunkSize, cEnd), numCells = cE - cS, c;
3907:     PetscInt         numFaces = 0;

3909:     /* Extract field coefficients */
3910:     if (useFEM) {
3911:       ISGetPointSubrange(chunkIS, cS, cE, cells);
3912:       DMPlexGetCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a);
3913:       DMGetWorkArray(dm, numCells*totDim, MPIU_SCALAR, &elemVec);
3914:       PetscArrayzero(elemVec, numCells*totDim);
3915:     }
3916:     /* TODO We will interlace both our field coefficients (u, u_t, uL, uR, etc.) and our output (elemVec, fL, fR). I think this works */
3917:     /* Loop over fields */
3918:     for (f = 0; f < Nf; ++f) {
3919:       PetscObject  obj;
3920:       PetscClassId id;
3921:       PetscBool    fimp;
3922:       PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;

3924:       key.field = f;
3925:       PetscDSGetImplicit(prob, f, &fimp);
3926:       if (isImplicit != fimp) continue;
3927:       PetscDSGetDiscretization(prob, f, &obj);
3928:       PetscObjectGetClassId(obj, &id);
3929:       if (id == PETSCFE_CLASSID) {
3930:         PetscFE         fe = (PetscFE) obj;
3931:         PetscFEGeom    *geom = affineGeom ? affineGeom : geoms[f];
3932:         PetscFEGeom    *chunkGeom = NULL;
3933:         PetscQuadrature quad = affineQuad ? affineQuad : quads[f];
3934:         PetscInt        Nq, Nb;

3936:         PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
3937:         PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);
3938:         PetscFEGetDimension(fe, &Nb);
3939:         blockSize = Nb;
3940:         batchSize = numBlocks * blockSize;
3941:         PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
3942:         numChunks = numCells / (numBatches*batchSize);
3943:         Ne        = numChunks*numBatches*batchSize;
3944:         Nr        = numCells % (numBatches*batchSize);
3945:         offset    = numCells - Nr;
3946:         /* Integrate FE residual to get elemVec (need fields at quadrature points) */
3947:         /*   For FV, I think we use a P0 basis and the cell coefficients (for subdivided cells, we can tweak the basis tabulation to be the indicator function) */
3948:         PetscFEGeomGetChunk(geom,0,offset,&chunkGeom);
3949:         PetscFEIntegrateResidual(prob, key, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec);
3950:         PetscFEGeomGetChunk(geom,offset,numCells,&chunkGeom);
3951:         PetscFEIntegrateResidual(prob, key, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, &elemVec[offset*totDim]);
3952:         PetscFEGeomRestoreChunk(geom,offset,numCells,&chunkGeom);
3953:       } else if (id == PETSCFV_CLASSID) {
3954:         PetscFV fv = (PetscFV) obj;

3956:         Ne = numFaces;
3957:         /* Riemann solve over faces (need fields at face centroids) */
3958:         /*   We need to evaluate FE fields at those coordinates */
3959:         PetscFVIntegrateRHSFunction(fv, prob, f, Ne, fgeom, vol, uL, uR, fluxL, fluxR);
3960:       } else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f);
3961:     }
3962:     /* Loop over domain */
3963:     if (useFEM) {
3964:       /* Add elemVec to locX */
3965:       for (c = cS; c < cE; ++c) {
3966:         const PetscInt cell = cells ? cells[c] : c;
3967:         const PetscInt cind = c - cStart;

3969:         if (mesh->printFEM > 1) DMPrintCellVector(cell, name, totDim, &elemVec[cind*totDim]);
3970:         if (ghostLabel) {
3971:           PetscInt ghostVal;

3973:           DMLabelGetValue(ghostLabel,cell,&ghostVal);
3974:           if (ghostVal > 0) continue;
3975:         }
3976:         DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cind*totDim], ADD_ALL_VALUES);
3977:       }
3978:     }
3979:     /* Handle time derivative */
3980:     if (locX_t) {
3981:       PetscScalar *x_t, *fa;

3983:       VecGetArray(locF, &fa);
3984:       VecGetArray(locX_t, &x_t);
3985:       for (f = 0; f < Nf; ++f) {
3986:         PetscFV      fv;
3987:         PetscObject  obj;
3988:         PetscClassId id;
3989:         PetscInt     pdim, d;

3991:         PetscDSGetDiscretization(prob, f, &obj);
3992:         PetscObjectGetClassId(obj, &id);
3993:         if (id != PETSCFV_CLASSID) continue;
3994:         fv   = (PetscFV) obj;
3995:         PetscFVGetNumComponents(fv, &pdim);
3996:         for (c = cS; c < cE; ++c) {
3997:           const PetscInt cell = cells ? cells[c] : c;
3998:           PetscScalar   *u_t, *r;

4000:           if (ghostLabel) {
4001:             PetscInt ghostVal;

4003:             DMLabelGetValue(ghostLabel, cell, &ghostVal);
4004:             if (ghostVal > 0) continue;
4005:           }
4006:           DMPlexPointLocalFieldRead(dm, cell, f, x_t, &u_t);
4007:           DMPlexPointLocalFieldRef(dm, cell, f, fa, &r);
4008:           for (d = 0; d < pdim; ++d) r[d] += u_t[d];
4009:         }
4010:       }
4011:       VecRestoreArray(locX_t, &x_t);
4012:       VecRestoreArray(locF, &fa);
4013:     }
4014:     if (useFEM) {
4015:       DMPlexRestoreCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a);
4016:       DMRestoreWorkArray(dm, numCells*totDim, MPIU_SCALAR, &elemVec);
4017:     }
4018:   }
4019:   if (useFEM) ISDestroy(&chunkIS);
4020:   ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);
4021:   /* TODO Could include boundary residual here (see DMPlexComputeResidual_Internal) */
4022:   if (useFEM) {
4023:     if (maxDegree <= 1) {
4024:       DMSNESRestoreFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);
4025:       PetscQuadratureDestroy(&affineQuad);
4026:     } else {
4027:       for (f = 0; f < Nf; ++f) {
4028:         DMSNESRestoreFEGeom(coordField,cellIS,quads[f],PETSC_FALSE,&geoms[f]);
4029:         PetscQuadratureDestroy(&quads[f]);
4030:       }
4031:       PetscFree2(quads,geoms);
4032:     }
4033:   }
4034:   PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);
4035:   return 0;
4036: }

4038: /*
4039:   We always assemble JacP, and if the matrix is different from Jac and two different sets of point functions are provided, we also assemble Jac

4041:   X   - The local solution vector
4042:   X_t - The local solution time derivative vector, or NULL
4043: */
4044: PetscErrorCode DMPlexComputeJacobian_Patch_Internal(DM dm, PetscSection section, PetscSection globalSection, IS cellIS,
4045:                                                     PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Mat Jac, Mat JacP, void *ctx)
4046: {
4047:   DM_Plex         *mesh  = (DM_Plex *) dm->data;
4048:   const char      *name = "Jacobian", *nameP = "JacobianPre";
4049:   DM               dmAux = NULL;
4050:   PetscDS          prob,   probAux = NULL;
4051:   PetscSection     sectionAux = NULL;
4052:   Vec              A;
4053:   DMField          coordField;
4054:   PetscFEGeom     *cgeomFEM;
4055:   PetscQuadrature  qGeom = NULL;
4056:   Mat              J = Jac, JP = JacP;
4057:   PetscScalar     *work, *u = NULL, *u_t = NULL, *a = NULL, *elemMat = NULL, *elemMatP = NULL, *elemMatD = NULL;
4058:   PetscBool        hasJac, hasPrec, hasDyn, assembleJac, *isFE, hasFV = PETSC_FALSE;
4059:   const PetscInt  *cells;
4060:   PetscFormKey key;
4061:   PetscInt         Nf, fieldI, fieldJ, maxDegree, numCells, cStart, cEnd, numChunks, chunkSize, chunk, totDim, totDimAux = 0, sz, wsz, off = 0, offCell = 0;

4063:   ISGetLocalSize(cellIS, &numCells);
4064:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
4065:   PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);
4066:   DMGetDS(dm, &prob);
4067:   DMGetAuxiliaryVec(dm, NULL, 0, 0, &A);
4068:   if (A) {
4069:     VecGetDM(A, &dmAux);
4070:     DMGetLocalSection(dmAux, &sectionAux);
4071:     DMGetDS(dmAux, &probAux);
4072:   }
4073:   /* Get flags */
4074:   PetscDSGetNumFields(prob, &Nf);
4075:   DMGetWorkArray(dm, Nf, MPIU_BOOL, &isFE);
4076:   for (fieldI = 0; fieldI < Nf; ++fieldI) {
4077:     PetscObject  disc;
4078:     PetscClassId id;
4079:     PetscDSGetDiscretization(prob, fieldI, &disc);
4080:     PetscObjectGetClassId(disc, &id);
4081:     if (id == PETSCFE_CLASSID)      {isFE[fieldI] = PETSC_TRUE;}
4082:     else if (id == PETSCFV_CLASSID) {hasFV = PETSC_TRUE; isFE[fieldI] = PETSC_FALSE;}
4083:   }
4084:   PetscDSHasJacobian(prob, &hasJac);
4085:   PetscDSHasJacobianPreconditioner(prob, &hasPrec);
4086:   PetscDSHasDynamicJacobian(prob, &hasDyn);
4087:   assembleJac = hasJac && hasPrec && (Jac != JacP) ? PETSC_TRUE : PETSC_FALSE;
4088:   hasDyn      = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
4089:   if (hasFV) MatSetOption(JP, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE); /* No allocated space for FV stuff, so ignore the zero entries */
4090:   PetscDSGetTotalDimension(prob, &totDim);
4091:   if (probAux) PetscDSGetTotalDimension(probAux, &totDimAux);
4092:   /* Compute batch sizes */
4093:   if (isFE[0]) {
4094:     PetscFE         fe;
4095:     PetscQuadrature q;
4096:     PetscInt        numQuadPoints, numBatches, batchSize, numBlocks, blockSize, Nb;

4098:     PetscDSGetDiscretization(prob, 0, (PetscObject *) &fe);
4099:     PetscFEGetQuadrature(fe, &q);
4100:     PetscQuadratureGetData(q, NULL, NULL, &numQuadPoints, NULL, NULL);
4101:     PetscFEGetDimension(fe, &Nb);
4102:     PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
4103:     blockSize = Nb*numQuadPoints;
4104:     batchSize = numBlocks  * blockSize;
4105:     chunkSize = numBatches * batchSize;
4106:     numChunks = numCells / chunkSize + numCells % chunkSize;
4107:     PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
4108:   } else {
4109:     chunkSize = numCells;
4110:     numChunks = 1;
4111:   }
4112:   /* Get work space */
4113:   wsz  = (((X?1:0) + (X_t?1:0) + (dmAux?1:0))*totDim + ((hasJac?1:0) + (hasPrec?1:0) + (hasDyn?1:0))*totDim*totDim)*chunkSize;
4114:   DMGetWorkArray(dm, wsz, MPIU_SCALAR, &work);
4115:   PetscArrayzero(work, wsz);
4116:   off      = 0;
4117:   u        = X       ? (sz = chunkSize*totDim,        off += sz, work+off-sz) : NULL;
4118:   u_t      = X_t     ? (sz = chunkSize*totDim,        off += sz, work+off-sz) : NULL;
4119:   a        = dmAux   ? (sz = chunkSize*totDimAux,     off += sz, work+off-sz) : NULL;
4120:   elemMat  = hasJac  ? (sz = chunkSize*totDim*totDim, off += sz, work+off-sz) : NULL;
4121:   elemMatP = hasPrec ? (sz = chunkSize*totDim*totDim, off += sz, work+off-sz) : NULL;
4122:   elemMatD = hasDyn  ? (sz = chunkSize*totDim*totDim, off += sz, work+off-sz) : NULL;
4124:   /* Setup geometry */
4125:   DMGetCoordinateField(dm, &coordField);
4126:   DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree);
4127:   if (maxDegree <= 1) DMFieldCreateDefaultQuadrature(coordField, cellIS, &qGeom);
4128:   if (!qGeom) {
4129:     PetscFE fe;

4131:     PetscDSGetDiscretization(prob, 0, (PetscObject *) &fe);
4132:     PetscFEGetQuadrature(fe, &qGeom);
4133:     PetscObjectReference((PetscObject) qGeom);
4134:   }
4135:   DMSNESGetFEGeom(coordField, cellIS, qGeom, PETSC_FALSE, &cgeomFEM);
4136:   /* Compute volume integrals */
4137:   if (assembleJac) MatZeroEntries(J);
4138:   MatZeroEntries(JP);
4139:   key.label = NULL;
4140:   key.value = 0;
4141:   key.part  = 0;
4142:   for (chunk = 0; chunk < numChunks; ++chunk, offCell += chunkSize) {
4143:     const PetscInt   Ncell = PetscMin(chunkSize, numCells - offCell);
4144:     PetscInt         c;

4146:     /* Extract values */
4147:     for (c = 0; c < Ncell; ++c) {
4148:       const PetscInt cell = cells ? cells[c+offCell] : c+offCell;
4149:       PetscScalar   *x = NULL,  *x_t = NULL;
4150:       PetscInt       i;

4152:       if (X) {
4153:         DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);
4154:         for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
4155:         DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);
4156:       }
4157:       if (X_t) {
4158:         DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);
4159:         for (i = 0; i < totDim; ++i) u_t[c*totDim+i] = x_t[i];
4160:         DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);
4161:       }
4162:       if (dmAux) {
4163:         DMPlexVecGetClosure(dmAux, sectionAux, A, cell, NULL, &x);
4164:         for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
4165:         DMPlexVecRestoreClosure(dmAux, sectionAux, A, cell, NULL, &x);
4166:       }
4167:     }
4168:     for (fieldI = 0; fieldI < Nf; ++fieldI) {
4169:       PetscFE fe;
4170:       PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
4171:       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
4172:         key.field = fieldI*Nf + fieldJ;
4173:         if (hasJac)  PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN,     key, Ncell, cgeomFEM, u, u_t, probAux, a, t, X_tShift, elemMat);
4174:         if (hasPrec) PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_PRE, key, Ncell, cgeomFEM, u, u_t, probAux, a, t, X_tShift, elemMatP);
4175:         if (hasDyn)  PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, key, Ncell, cgeomFEM, u, u_t, probAux, a, t, X_tShift, elemMatD);
4176:       }
4177:       /* For finite volume, add the identity */
4178:       if (!isFE[fieldI]) {
4179:         PetscFV  fv;
4180:         PetscInt eOffset = 0, Nc, fc, foff;

4182:         PetscDSGetFieldOffset(prob, fieldI, &foff);
4183:         PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fv);
4184:         PetscFVGetNumComponents(fv, &Nc);
4185:         for (c = 0; c < chunkSize; ++c, eOffset += totDim*totDim) {
4186:           for (fc = 0; fc < Nc; ++fc) {
4187:             const PetscInt i = foff + fc;
4188:             if (hasJac)  {elemMat [eOffset+i*totDim+i] = 1.0;}
4189:             if (hasPrec) {elemMatP[eOffset+i*totDim+i] = 1.0;}
4190:           }
4191:         }
4192:       }
4193:     }
4194:     /*   Add contribution from X_t */
4195:     if (hasDyn) {for (c = 0; c < chunkSize*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];}
4196:     /* Insert values into matrix */
4197:     for (c = 0; c < Ncell; ++c) {
4198:       const PetscInt cell = cells ? cells[c+offCell] : c+offCell;
4199:       if (mesh->printFEM > 1) {
4200:         if (hasJac)  DMPrintCellMatrix(cell, name,  totDim, totDim, &elemMat[(c-cStart)*totDim*totDim]);
4201:         if (hasPrec) DMPrintCellMatrix(cell, nameP, totDim, totDim, &elemMatP[(c-cStart)*totDim*totDim]);
4202:       }
4203:       if (assembleJac) DMPlexMatSetClosure(dm, section, globalSection, Jac, cell, &elemMat[(c-cStart)*totDim*totDim], ADD_VALUES);
4204:       DMPlexMatSetClosure(dm, section, globalSection, JP, cell, &elemMat[(c-cStart)*totDim*totDim], ADD_VALUES);
4205:     }
4206:   }
4207:   /* Cleanup */
4208:   DMSNESRestoreFEGeom(coordField, cellIS, qGeom, PETSC_FALSE, &cgeomFEM);
4209:   PetscQuadratureDestroy(&qGeom);
4210:   if (hasFV) MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE);
4211:   DMRestoreWorkArray(dm, Nf, MPIU_BOOL, &isFE);
4212:   DMRestoreWorkArray(dm, ((1 + (X_t?1:0) + (dmAux?1:0))*totDim + ((hasJac?1:0) + (hasPrec?1:0) + (hasDyn?1:0))*totDim*totDim)*chunkSize, MPIU_SCALAR, &work);
4213:   /* Compute boundary integrals */
4214:   /* DMPlexComputeBdJacobian_Internal(dm, X, X_t, t, X_tShift, Jac, JacP, ctx); */
4215:   /* Assemble matrix */
4216:   if (assembleJac) {MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY));PetscCall(MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY);}
4217:   MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY));PetscCall(MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);
4218:   PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);
4219:   return 0;
4220: }

4222: /******** FEM Assembly Function ********/

4224: static PetscErrorCode DMConvertPlex_Internal(DM dm, DM *plex, PetscBool copy)
4225: {
4226:   PetscBool      isPlex;

4228:   PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);
4229:   if (isPlex) {
4230:     *plex = dm;
4231:     PetscObjectReference((PetscObject) dm);
4232:   } else {
4233:     PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);
4234:     if (!*plex) {
4235:       DMConvert(dm,DMPLEX,plex);
4236:       PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);
4237:       if (copy) {
4238:         DMCopyAuxiliaryVec(dm, *plex);
4239:       }
4240:     } else {
4241:       PetscObjectReference((PetscObject) *plex);
4242:     }
4243:   }
4244:   return 0;
4245: }

4247: /*@
4248:   DMPlexGetGeometryFVM - Return precomputed geometric data

4250:   Collective on DM

4252:   Input Parameter:
4253: . dm - The DM

4255:   Output Parameters:
4256: + facegeom - The values precomputed from face geometry
4257: . cellgeom - The values precomputed from cell geometry
4258: - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell

4260:   Level: developer

4262: .seealso: DMTSSetRHSFunctionLocal()
4263: @*/
4264: PetscErrorCode DMPlexGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
4265: {
4266:   DM             plex;

4269:   DMConvertPlex_Internal(dm,&plex,PETSC_TRUE);
4270:   DMPlexGetDataFVM(plex, NULL, cellgeom, facegeom, NULL);
4271:   if (minRadius) DMPlexGetMinRadius(plex, minRadius);
4272:   DMDestroy(&plex);
4273:   return 0;
4274: }

4276: /*@
4277:   DMPlexGetGradientDM - Return gradient data layout

4279:   Collective on DM

4281:   Input Parameters:
4282: + dm - The DM
4283: - fv - The PetscFV

4285:   Output Parameter:
4286: . dmGrad - The layout for gradient values

4288:   Level: developer

4290: .seealso: DMPlexGetGeometryFVM()
4291: @*/
4292: PetscErrorCode DMPlexGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
4293: {
4294:   DM             plex;
4295:   PetscBool      computeGradients;

4300:   PetscFVGetComputeGradients(fv, &computeGradients);
4301:   if (!computeGradients) {*dmGrad = NULL; return 0;}
4302:   DMConvertPlex_Internal(dm,&plex,PETSC_TRUE);
4303:   DMPlexGetDataFVM(plex, fv, NULL, NULL, dmGrad);
4304:   DMDestroy(&plex);
4305:   return 0;
4306: }

4308: static PetscErrorCode DMPlexComputeBdResidual_Single_Internal(DM dm, PetscReal t, PetscWeakForm wf, PetscFormKey key, Vec locX, Vec locX_t, Vec locF, DMField coordField, IS facetIS)
4309: {
4310:   DM_Plex         *mesh = (DM_Plex *) dm->data;
4311:   DM               plex = NULL, plexA = NULL;
4312:   DMEnclosureType  encAux;
4313:   PetscDS          prob, probAux = NULL;
4314:   PetscSection     section, sectionAux = NULL;
4315:   Vec              locA = NULL;
4316:   PetscScalar     *u = NULL, *u_t = NULL, *a = NULL, *elemVec = NULL;
4317:   PetscInt         totDim, totDimAux = 0;

4319:   DMConvert(dm, DMPLEX, &plex);
4320:   DMGetLocalSection(dm, &section);
4321:   DMGetDS(dm, &prob);
4322:   PetscDSGetTotalDimension(prob, &totDim);
4323:   DMGetAuxiliaryVec(dm, key.label, key.value, key.part, &locA);
4324:   if (locA) {
4325:     DM dmAux;

4327:     VecGetDM(locA, &dmAux);
4328:     DMGetEnclosureRelation(dmAux, dm, &encAux);
4329:     DMConvert(dmAux, DMPLEX, &plexA);
4330:     DMGetDS(plexA, &probAux);
4331:     PetscDSGetTotalDimension(probAux, &totDimAux);
4332:     DMGetLocalSection(plexA, &sectionAux);
4333:   }
4334:   {
4335:     PetscFEGeom     *fgeom;
4336:     PetscInt         maxDegree;
4337:     PetscQuadrature  qGeom = NULL;
4338:     IS               pointIS;
4339:     const PetscInt  *points;
4340:     PetscInt         numFaces, face, Nq;

4342:     DMLabelGetStratumIS(key.label, key.value, &pointIS);
4343:     if (!pointIS) goto end; /* No points with that id on this process */
4344:     {
4345:       IS isectIS;

4347:       /* TODO: Special cases of ISIntersect where it is quick to check a priori if one is a superset of the other */
4348:       ISIntersect_Caching_Internal(facetIS,pointIS,&isectIS);
4349:       ISDestroy(&pointIS);
4350:       pointIS = isectIS;
4351:     }
4352:     ISGetLocalSize(pointIS,&numFaces);
4353:     ISGetIndices(pointIS,&points);
4354:     PetscMalloc4(numFaces*totDim, &u, locX_t ? numFaces*totDim : 0, &u_t, numFaces*totDim, &elemVec, locA ? numFaces*totDimAux : 0, &a);
4355:     DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);
4356:     if (maxDegree <= 1) {
4357:       DMFieldCreateDefaultQuadrature(coordField,pointIS,&qGeom);
4358:     }
4359:     if (!qGeom) {
4360:       PetscFE fe;

4362:       PetscDSGetDiscretization(prob, key.field, (PetscObject *) &fe);
4363:       PetscFEGetFaceQuadrature(fe, &qGeom);
4364:       PetscObjectReference((PetscObject)qGeom);
4365:     }
4366:     PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
4367:     DMSNESGetFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);
4368:     for (face = 0; face < numFaces; ++face) {
4369:       const PetscInt point = points[face], *support;
4370:       PetscScalar   *x     = NULL;
4371:       PetscInt       i;

4373:       DMPlexGetSupport(dm, point, &support);
4374:       DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x);
4375:       for (i = 0; i < totDim; ++i) u[face*totDim+i] = x[i];
4376:       DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x);
4377:       if (locX_t) {
4378:         DMPlexVecGetClosure(plex, section, locX_t, support[0], NULL, &x);
4379:         for (i = 0; i < totDim; ++i) u_t[face*totDim+i] = x[i];
4380:         DMPlexVecRestoreClosure(plex, section, locX_t, support[0], NULL, &x);
4381:       }
4382:       if (locA) {
4383:         PetscInt subp;

4385:         DMGetEnclosurePoint(plexA, dm, encAux, support[0], &subp);
4386:         DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x);
4387:         for (i = 0; i < totDimAux; ++i) a[face*totDimAux+i] = x[i];
4388:         DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x);
4389:       }
4390:     }
4391:     PetscArrayzero(elemVec, numFaces*totDim);
4392:     {
4393:       PetscFE         fe;
4394:       PetscInt        Nb;
4395:       PetscFEGeom     *chunkGeom = NULL;
4396:       /* Conforming batches */
4397:       PetscInt        numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
4398:       /* Remainder */
4399:       PetscInt        Nr, offset;

4401:       PetscDSGetDiscretization(prob, key.field, (PetscObject *) &fe);
4402:       PetscFEGetDimension(fe, &Nb);
4403:       PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
4404:       /* TODO: documentation is unclear about what is going on with these numbers: how should Nb / Nq factor in ? */
4405:       blockSize = Nb;
4406:       batchSize = numBlocks * blockSize;
4407:       PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
4408:       numChunks = numFaces / (numBatches*batchSize);
4409:       Ne        = numChunks*numBatches*batchSize;
4410:       Nr        = numFaces % (numBatches*batchSize);
4411:       offset    = numFaces - Nr;
4412:       PetscFEGeomGetChunk(fgeom,0,offset,&chunkGeom);
4413:       PetscFEIntegrateBdResidual(prob, wf, key, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec);
4414:       PetscFEGeomRestoreChunk(fgeom, 0, offset, &chunkGeom);
4415:       PetscFEGeomGetChunk(fgeom,offset,numFaces,&chunkGeom);
4416:       PetscFEIntegrateBdResidual(prob, wf, key, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, a ? &a[offset*totDimAux] : NULL, t, &elemVec[offset*totDim]);
4417:       PetscFEGeomRestoreChunk(fgeom,offset,numFaces,&chunkGeom);
4418:     }
4419:     for (face = 0; face < numFaces; ++face) {
4420:       const PetscInt point = points[face], *support;

4422:       if (mesh->printFEM > 1) DMPrintCellVector(point, "BdResidual", totDim, &elemVec[face*totDim]);
4423:       DMPlexGetSupport(plex, point, &support);
4424:       DMPlexVecSetClosure(plex, NULL, locF, support[0], &elemVec[face*totDim], ADD_ALL_VALUES);
4425:     }
4426:     DMSNESRestoreFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);
4427:     PetscQuadratureDestroy(&qGeom);
4428:     ISRestoreIndices(pointIS, &points);
4429:     ISDestroy(&pointIS);
4430:     PetscFree4(u, u_t, elemVec, a);
4431:   }
4432:   end:
4433:   DMDestroy(&plex);
4434:   DMDestroy(&plexA);
4435:   return 0;
4436: }

4438: PetscErrorCode DMPlexComputeBdResidualSingle(DM dm, PetscReal t, PetscWeakForm wf, PetscFormKey key, Vec locX, Vec locX_t, Vec locF)
4439: {
4440:   DMField        coordField;
4441:   DMLabel        depthLabel;
4442:   IS             facetIS;
4443:   PetscInt       dim;

4445:   DMGetDimension(dm, &dim);
4446:   DMPlexGetDepthLabel(dm, &depthLabel);
4447:   DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);
4448:   DMGetCoordinateField(dm, &coordField);
4449:   DMPlexComputeBdResidual_Single_Internal(dm, t, wf, key, locX, locX_t, locF, coordField, facetIS);
4450:   ISDestroy(&facetIS);
4451:   return 0;
4452: }

4454: PetscErrorCode DMPlexComputeBdResidual_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
4455: {
4456:   PetscDS        prob;
4457:   PetscInt       numBd, bd;
4458:   DMField        coordField = NULL;
4459:   IS             facetIS    = NULL;
4460:   DMLabel        depthLabel;
4461:   PetscInt       dim;

4463:   DMGetDS(dm, &prob);
4464:   DMPlexGetDepthLabel(dm, &depthLabel);
4465:   DMGetDimension(dm, &dim);
4466:   DMLabelGetStratumIS(depthLabel,dim - 1,&facetIS);
4467:   PetscDSGetNumBoundary(prob, &numBd);
4468:   for (bd = 0; bd < numBd; ++bd) {
4469:     PetscWeakForm           wf;
4470:     DMBoundaryConditionType type;
4471:     DMLabel                 label;
4472:     const PetscInt         *values;
4473:     PetscInt                field, numValues, v;
4474:     PetscObject             obj;
4475:     PetscClassId            id;
4476:     PetscFormKey            key;

4478:     PetscDSGetBoundary(prob, bd, &wf, &type, NULL, &label, &numValues, &values, &field, NULL, NULL, NULL, NULL, NULL);
4479:     PetscDSGetDiscretization(prob, field, &obj);
4480:     PetscObjectGetClassId(obj, &id);
4481:     if ((id != PETSCFE_CLASSID) || (type & DM_BC_ESSENTIAL)) continue;
4482:     if (!facetIS) {
4483:       DMLabel  depthLabel;
4484:       PetscInt dim;

4486:       DMPlexGetDepthLabel(dm, &depthLabel);
4487:       DMGetDimension(dm, &dim);
4488:       DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS);
4489:     }
4490:     DMGetCoordinateField(dm, &coordField);
4491:     for (v = 0; v < numValues; ++v) {
4492:       key.label = label;
4493:       key.value = values[v];
4494:       key.field = field;
4495:       key.part  = 0;
4496:       DMPlexComputeBdResidual_Single_Internal(dm, t, wf, key, locX, locX_t, locF, coordField, facetIS);
4497:     }
4498:   }
4499:   ISDestroy(&facetIS);
4500:   return 0;
4501: }

4503: PetscErrorCode DMPlexComputeResidual_Internal(DM dm, PetscFormKey key, IS cellIS, PetscReal time, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
4504: {
4505:   DM_Plex         *mesh       = (DM_Plex *) dm->data;
4506:   const char      *name       = "Residual";
4507:   DM               dmAux      = NULL;
4508:   DM               dmGrad     = NULL;
4509:   DMLabel          ghostLabel = NULL;
4510:   PetscDS          ds         = NULL;
4511:   PetscDS          dsAux      = NULL;
4512:   PetscSection     section    = NULL;
4513:   PetscBool        useFEM     = PETSC_FALSE;
4514:   PetscBool        useFVM     = PETSC_FALSE;
4515:   PetscBool        isImplicit = (locX_t || time == PETSC_MIN_REAL) ? PETSC_TRUE : PETSC_FALSE;
4516:   PetscFV          fvm        = NULL;
4517:   PetscFVCellGeom *cgeomFVM   = NULL;
4518:   PetscFVFaceGeom *fgeomFVM   = NULL;
4519:   DMField          coordField = NULL;
4520:   Vec              locA, cellGeometryFVM = NULL, faceGeometryFVM = NULL, grad, locGrad = NULL;
4521:   PetscScalar     *u = NULL, *u_t, *a, *uL, *uR;
4522:   IS               chunkIS;
4523:   const PetscInt  *cells;
4524:   PetscInt         cStart, cEnd, numCells;
4525:   PetscInt         Nf, f, totDim, totDimAux, numChunks, cellChunkSize, faceChunkSize, chunk, fStart, fEnd;
4526:   PetscInt         maxDegree = PETSC_MAX_INT;
4527:   PetscQuadrature  affineQuad = NULL, *quads = NULL;
4528:   PetscFEGeom     *affineGeom = NULL, **geoms = NULL;

4530:   PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);
4531:   /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */
4532:   /* TODO The FVM geometry is over-manipulated. Make the precalc functions return exactly what we need */
4533:   /* FEM+FVM */
4534:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
4535:   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
4536:   /* 1: Get sizes from dm and dmAux */
4537:   DMGetLocalSection(dm, &section);
4538:   DMGetLabel(dm, "ghost", &ghostLabel);
4539:   DMGetCellDS(dm, cells ? cells[cStart] : cStart, &ds);
4540:   PetscDSGetNumFields(ds, &Nf);
4541:   PetscDSGetTotalDimension(ds, &totDim);
4542:   DMGetAuxiliaryVec(dm, key.label, key.value, key.part, &locA);
4543:   if (locA) {
4544:     PetscInt subcell;
4545:     VecGetDM(locA, &dmAux);
4546:     DMGetEnclosurePoint(dmAux, dm, DM_ENC_UNKNOWN, cStart, &subcell);
4547:     DMGetCellDS(dmAux, subcell, &dsAux);
4548:     PetscDSGetTotalDimension(dsAux, &totDimAux);
4549:   }
4550:   /* 2: Get geometric data */
4551:   for (f = 0; f < Nf; ++f) {
4552:     PetscObject  obj;
4553:     PetscClassId id;
4554:     PetscBool    fimp;

4556:     PetscDSGetImplicit(ds, f, &fimp);
4557:     if (isImplicit != fimp) continue;
4558:     PetscDSGetDiscretization(ds, f, &obj);
4559:     PetscObjectGetClassId(obj, &id);
4560:     if (id == PETSCFE_CLASSID) {useFEM = PETSC_TRUE;}
4561:     if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;}
4562:   }
4563:   if (useFEM) {
4564:     DMGetCoordinateField(dm, &coordField);
4565:     DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);
4566:     if (maxDegree <= 1) {
4567:       DMFieldCreateDefaultQuadrature(coordField,cellIS,&affineQuad);
4568:       if (affineQuad) {
4569:         DMSNESGetFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);
4570:       }
4571:     } else {
4572:       PetscCalloc2(Nf,&quads,Nf,&geoms);
4573:       for (f = 0; f < Nf; ++f) {
4574:         PetscObject  obj;
4575:         PetscClassId id;
4576:         PetscBool    fimp;

4578:         PetscDSGetImplicit(ds, f, &fimp);
4579:         if (isImplicit != fimp) continue;
4580:         PetscDSGetDiscretization(ds, f, &obj);
4581:         PetscObjectGetClassId(obj, &id);
4582:         if (id == PETSCFE_CLASSID) {
4583:           PetscFE fe = (PetscFE) obj;

4585:           PetscFEGetQuadrature(fe, &quads[f]);
4586:           PetscObjectReference((PetscObject)quads[f]);
4587:           DMSNESGetFEGeom(coordField,cellIS,quads[f],PETSC_FALSE,&geoms[f]);
4588:         }
4589:       }
4590:     }
4591:   }
4592:   if (useFVM) {
4593:     DMPlexGetGeometryFVM(dm, &faceGeometryFVM, &cellGeometryFVM, NULL);
4594:     VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);
4595:     VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);
4596:     /* Reconstruct and limit cell gradients */
4597:     DMPlexGetGradientDM(dm, fvm, &dmGrad);
4598:     if (dmGrad) {
4599:       DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
4600:       DMGetGlobalVector(dmGrad, &grad);
4601:       DMPlexReconstructGradients_Internal(dm, fvm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad);
4602:       /* Communicate gradient values */
4603:       DMGetLocalVector(dmGrad, &locGrad);
4604:       DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);
4605:       DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);
4606:       DMRestoreGlobalVector(dmGrad, &grad);
4607:     }
4608:     /* Handle non-essential (e.g. outflow) boundary values */
4609:     DMPlexInsertBoundaryValues(dm, PETSC_FALSE, locX, time, faceGeometryFVM, cellGeometryFVM, locGrad);
4610:   }
4611:   /* Loop over chunks */
4612:   if (useFEM) ISCreate(PETSC_COMM_SELF, &chunkIS);
4613:   numCells      = cEnd - cStart;
4614:   numChunks     = 1;
4615:   cellChunkSize = numCells/numChunks;
4616:   faceChunkSize = (fEnd - fStart)/numChunks;
4617:   numChunks     = PetscMin(1,numCells);
4618:   for (chunk = 0; chunk < numChunks; ++chunk) {
4619:     PetscScalar     *elemVec, *fluxL, *fluxR;
4620:     PetscReal       *vol;
4621:     PetscFVFaceGeom *fgeom;
4622:     PetscInt         cS = cStart+chunk*cellChunkSize, cE = PetscMin(cS+cellChunkSize, cEnd), numCells = cE - cS, c;
4623:     PetscInt         fS = fStart+chunk*faceChunkSize, fE = PetscMin(fS+faceChunkSize, fEnd), numFaces = 0, face;

4625:     /* Extract field coefficients */
4626:     if (useFEM) {
4627:       ISGetPointSubrange(chunkIS, cS, cE, cells);
4628:       DMPlexGetCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a);
4629:       DMGetWorkArray(dm, numCells*totDim, MPIU_SCALAR, &elemVec);
4630:       PetscArrayzero(elemVec, numCells*totDim);
4631:     }
4632:     if (useFVM) {
4633:       DMPlexGetFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &numFaces, &uL, &uR);
4634:       DMPlexGetFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &numFaces, &fgeom, &vol);
4635:       DMGetWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxL);
4636:       DMGetWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxR);
4637:       PetscArrayzero(fluxL, numFaces*totDim);
4638:       PetscArrayzero(fluxR, numFaces*totDim);
4639:     }
4640:     /* TODO We will interlace both our field coefficients (u, u_t, uL, uR, etc.) and our output (elemVec, fL, fR). I think this works */
4641:     /* Loop over fields */
4642:     for (f = 0; f < Nf; ++f) {
4643:       PetscObject  obj;
4644:       PetscClassId id;
4645:       PetscBool    fimp;
4646:       PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;

4648:       key.field = f;
4649:       PetscDSGetImplicit(ds, f, &fimp);
4650:       if (isImplicit != fimp) continue;
4651:       PetscDSGetDiscretization(ds, f, &obj);
4652:       PetscObjectGetClassId(obj, &id);
4653:       if (id == PETSCFE_CLASSID) {
4654:         PetscFE         fe = (PetscFE) obj;
4655:         PetscFEGeom    *geom = affineGeom ? affineGeom : geoms[f];
4656:         PetscFEGeom    *chunkGeom = NULL;
4657:         PetscQuadrature quad = affineQuad ? affineQuad : quads[f];
4658:         PetscInt        Nq, Nb;

4660:         PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
4661:         PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);
4662:         PetscFEGetDimension(fe, &Nb);
4663:         blockSize = Nb;
4664:         batchSize = numBlocks * blockSize;
4665:         PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
4666:         numChunks = numCells / (numBatches*batchSize);
4667:         Ne        = numChunks*numBatches*batchSize;
4668:         Nr        = numCells % (numBatches*batchSize);
4669:         offset    = numCells - Nr;
4670:         /* Integrate FE residual to get elemVec (need fields at quadrature points) */
4671:         /*   For FV, I think we use a P0 basis and the cell coefficients (for subdivided cells, we can tweak the basis tabulation to be the indicator function) */
4672:         PetscFEGeomGetChunk(geom,0,offset,&chunkGeom);
4673:         PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, u_t, dsAux, a, t, elemVec);
4674:         PetscFEGeomGetChunk(geom,offset,numCells,&chunkGeom);
4675:         PetscFEIntegrateResidual(ds, key, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, dsAux, &a[offset*totDimAux], t, &elemVec[offset*totDim]);
4676:         PetscFEGeomRestoreChunk(geom,offset,numCells,&chunkGeom);
4677:       } else if (id == PETSCFV_CLASSID) {
4678:         PetscFV fv = (PetscFV) obj;

4680:         Ne = numFaces;
4681:         /* Riemann solve over faces (need fields at face centroids) */
4682:         /*   We need to evaluate FE fields at those coordinates */
4683:         PetscFVIntegrateRHSFunction(fv, ds, f, Ne, fgeom, vol, uL, uR, fluxL, fluxR);
4684:       } else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f);
4685:     }
4686:     /* Loop over domain */
4687:     if (useFEM) {
4688:       /* Add elemVec to locX */
4689:       for (c = cS; c < cE; ++c) {
4690:         const PetscInt cell = cells ? cells[c] : c;
4691:         const PetscInt cind = c - cStart;

4693:         if (mesh->printFEM > 1) DMPrintCellVector(cell, name, totDim, &elemVec[cind*totDim]);
4694:         if (ghostLabel) {
4695:           PetscInt ghostVal;

4697:           DMLabelGetValue(ghostLabel,cell,&ghostVal);
4698:           if (ghostVal > 0) continue;
4699:         }
4700:         DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cind*totDim], ADD_ALL_VALUES);
4701:       }
4702:     }
4703:     if (useFVM) {
4704:       PetscScalar *fa;
4705:       PetscInt     iface;

4707:       VecGetArray(locF, &fa);
4708:       for (f = 0; f < Nf; ++f) {
4709:         PetscFV      fv;
4710:         PetscObject  obj;
4711:         PetscClassId id;
4712:         PetscInt     foff, pdim;

4714:         PetscDSGetDiscretization(ds, f, &obj);
4715:         PetscDSGetFieldOffset(ds, f, &foff);
4716:         PetscObjectGetClassId(obj, &id);
4717:         if (id != PETSCFV_CLASSID) continue;
4718:         fv   = (PetscFV) obj;
4719:         PetscFVGetNumComponents(fv, &pdim);
4720:         /* Accumulate fluxes to cells */
4721:         for (face = fS, iface = 0; face < fE; ++face) {
4722:           const PetscInt *scells;
4723:           PetscScalar    *fL = NULL, *fR = NULL;
4724:           PetscInt        ghost, d, nsupp, nchild;

4726:           DMLabelGetValue(ghostLabel, face, &ghost);
4727:           DMPlexGetSupportSize(dm, face, &nsupp);
4728:           DMPlexGetTreeChildren(dm, face, &nchild, NULL);
4729:           if (ghost >= 0 || nsupp > 2 || nchild > 0) continue;
4730:           DMPlexGetSupport(dm, face, &scells);
4731:           DMLabelGetValue(ghostLabel,scells[0],&ghost);
4732:           if (ghost <= 0) DMPlexPointLocalFieldRef(dm, scells[0], f, fa, &fL);
4733:           DMLabelGetValue(ghostLabel,scells[1],&ghost);
4734:           if (ghost <= 0) DMPlexPointLocalFieldRef(dm, scells[1], f, fa, &fR);
4735:           for (d = 0; d < pdim; ++d) {
4736:             if (fL) fL[d] -= fluxL[iface*totDim+foff+d];
4737:             if (fR) fR[d] += fluxR[iface*totDim+foff+d];
4738:           }
4739:           ++iface;
4740:         }
4741:       }
4742:       VecRestoreArray(locF, &fa);
4743:     }
4744:     /* Handle time derivative */
4745:     if (locX_t) {
4746:       PetscScalar *x_t, *fa;

4748:       VecGetArray(locF, &fa);
4749:       VecGetArray(locX_t, &x_t);
4750:       for (f = 0; f < Nf; ++f) {
4751:         PetscFV      fv;
4752:         PetscObject  obj;
4753:         PetscClassId id;
4754:         PetscInt     pdim, d;

4756:         PetscDSGetDiscretization(ds, f, &obj);
4757:         PetscObjectGetClassId(obj, &id);
4758:         if (id != PETSCFV_CLASSID) continue;
4759:         fv   = (PetscFV) obj;
4760:         PetscFVGetNumComponents(fv, &pdim);
4761:         for (c = cS; c < cE; ++c) {
4762:           const PetscInt cell = cells ? cells[c] : c;
4763:           PetscScalar   *u_t, *r;

4765:           if (ghostLabel) {
4766:             PetscInt ghostVal;

4768:             DMLabelGetValue(ghostLabel, cell, &ghostVal);
4769:             if (ghostVal > 0) continue;
4770:           }
4771:           DMPlexPointLocalFieldRead(dm, cell, f, x_t, &u_t);
4772:           DMPlexPointLocalFieldRef(dm, cell, f, fa, &r);
4773:           for (d = 0; d < pdim; ++d) r[d] += u_t[d];
4774:         }
4775:       }
4776:       VecRestoreArray(locX_t, &x_t);
4777:       VecRestoreArray(locF, &fa);
4778:     }
4779:     if (useFEM) {
4780:       DMPlexRestoreCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a);
4781:       DMRestoreWorkArray(dm, numCells*totDim, MPIU_SCALAR, &elemVec);
4782:     }
4783:     if (useFVM) {
4784:       DMPlexRestoreFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &numFaces, &uL, &uR);
4785:       DMPlexRestoreFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &numFaces, &fgeom, &vol);
4786:       DMRestoreWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxL);
4787:       DMRestoreWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxR);
4788:       if (dmGrad) DMRestoreLocalVector(dmGrad, &locGrad);
4789:     }
4790:   }
4791:   if (useFEM) ISDestroy(&chunkIS);
4792:   ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);

4794:   if (useFEM) {
4795:     DMPlexComputeBdResidual_Internal(dm, locX, locX_t, t, locF, user);

4797:     if (maxDegree <= 1) {
4798:       DMSNESRestoreFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);
4799:       PetscQuadratureDestroy(&affineQuad);
4800:     } else {
4801:       for (f = 0; f < Nf; ++f) {
4802:         DMSNESRestoreFEGeom(coordField,cellIS,quads[f],PETSC_FALSE,&geoms[f]);
4803:         PetscQuadratureDestroy(&quads[f]);
4804:       }
4805:       PetscFree2(quads,geoms);
4806:     }
4807:   }

4809:   /* FEM */
4810:   /* 1: Get sizes from dm and dmAux */
4811:   /* 2: Get geometric data */
4812:   /* 3: Handle boundary values */
4813:   /* 4: Loop over domain */
4814:   /*   Extract coefficients */
4815:   /* Loop over fields */
4816:   /*   Set tiling for FE*/
4817:   /*   Integrate FE residual to get elemVec */
4818:   /*     Loop over subdomain */
4819:   /*       Loop over quad points */
4820:   /*         Transform coords to real space */
4821:   /*         Evaluate field and aux fields at point */
4822:   /*         Evaluate residual at point */
4823:   /*         Transform residual to real space */
4824:   /*       Add residual to elemVec */
4825:   /* Loop over domain */
4826:   /*   Add elemVec to locX */

4828:   /* FVM */
4829:   /* Get geometric data */
4830:   /* If using gradients */
4831:   /*   Compute gradient data */
4832:   /*   Loop over domain faces */
4833:   /*     Count computational faces */
4834:   /*     Reconstruct cell gradient */
4835:   /*   Loop over domain cells */
4836:   /*     Limit cell gradients */
4837:   /* Handle boundary values */
4838:   /* Loop over domain faces */
4839:   /*   Read out field, centroid, normal, volume for each side of face */
4840:   /* Riemann solve over faces */
4841:   /* Loop over domain faces */
4842:   /*   Accumulate fluxes to cells */
4843:   /* TODO Change printFEM to printDisc here */
4844:   if (mesh->printFEM) {
4845:     Vec         locFbc;
4846:     PetscInt    pStart, pEnd, p, maxDof;
4847:     PetscScalar *zeroes;

4849:     VecDuplicate(locF,&locFbc);
4850:     VecCopy(locF,locFbc);
4851:     PetscSectionGetChart(section,&pStart,&pEnd);
4852:     PetscSectionGetMaxDof(section,&maxDof);
4853:     PetscCalloc1(maxDof,&zeroes);
4854:     for (p = pStart; p < pEnd; p++) {
4855:       VecSetValuesSection(locFbc,section,p,zeroes,INSERT_BC_VALUES);
4856:     }
4857:     PetscFree(zeroes);
4858:     DMPrintLocalVec(dm, name, mesh->printTol, locFbc);
4859:     VecDestroy(&locFbc);
4860:   }
4861:   PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);
4862:   return 0;
4863: }

4865: /*
4866:   1) Allow multiple kernels for BdResidual for hybrid DS

4868:   DONE 2) Get out dsAux for either side at the same time as cohesive cell dsAux

4870:   DONE 3) Change DMGetCellFields() to get different aux data a[] for each side
4871:      - I think I just need to replace a[] with the closure from each face

4873:   4) Run both kernels for each non-hybrid field with correct dsAux, and then hybrid field as before
4874: */
4875: PetscErrorCode DMPlexComputeResidual_Hybrid_Internal(DM dm, PetscFormKey key[], IS cellIS, PetscReal time, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
4876: {
4877:   DM_Plex         *mesh       = (DM_Plex *) dm->data;
4878:   const char      *name       = "Hybrid Residual";
4879:   DM               dmAux[3]   = {NULL, NULL, NULL};
4880:   DMLabel          ghostLabel = NULL;
4881:   PetscDS          ds         = NULL;
4882:   PetscDS          dsAux[3]   = {NULL, NULL, NULL};
4883:   Vec              locA[3]    = {NULL, NULL, NULL};
4884:   PetscSection     section    = NULL;
4885:   DMField          coordField = NULL;
4886:   PetscScalar     *u = NULL, *u_t, *a[3];
4887:   PetscScalar     *elemVec;
4888:   IS               chunkIS;
4889:   const PetscInt  *cells;
4890:   PetscInt        *faces;
4891:   PetscInt         cStart, cEnd, numCells;
4892:   PetscInt         Nf, f, totDim, totDimAux[3], numChunks, cellChunkSize, chunk;
4893:   PetscInt         maxDegree = PETSC_MAX_INT;
4894:   PetscQuadrature  affineQuad = NULL, *quads = NULL;
4895:   PetscFEGeom     *affineGeom = NULL, **geoms = NULL;

4897:   if ((key[0].label == key[1].label) && (key[0].value == key[1].value) && (key[0].part == key[1].part)) {
4898:     const char *name;
4899:     PetscObjectGetName((PetscObject) key[0].label, &name);
4900:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Form keys for each side of a cohesive surface must be different (%s, %D, %D)", name, key[0].value, key[0].part);
4901:   }
4902:   PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);
4903:   /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */
4904:   /* FEM */
4905:   ISGetLocalSize(cellIS, &numCells);
4906:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
4907:   /* 1: Get sizes from dm and dmAux */
4908:   DMGetSection(dm, &section);
4909:   DMGetLabel(dm, "ghost", &ghostLabel);
4910:   DMGetCellDS(dm, cStart, &ds);
4911:   PetscDSGetNumFields(ds, &Nf);
4912:   PetscDSGetTotalDimension(ds, &totDim);
4913:   DMGetAuxiliaryVec(dm, key[2].label, key[2].value, key[2].part, &locA[2]);
4914:   if (locA[2]) {
4915:     VecGetDM(locA[2], &dmAux[2]);
4916:     DMGetCellDS(dmAux[2], cStart, &dsAux[2]);
4917:     PetscDSGetTotalDimension(dsAux[2], &totDimAux[2]);
4918:     {
4919:       const PetscInt *cone;
4920:       PetscInt        c;

4922:       DMPlexGetCone(dm, cStart, &cone);
4923:       for (c = 0; c < 2; ++c) {
4924:         const PetscInt *support;
4925:         PetscInt ssize, s;

4927:         DMPlexGetSupport(dm, cone[c], &support);
4928:         DMPlexGetSupportSize(dm, cone[c], &ssize);
4930:         if      (support[0] == cStart) s = 1;
4931:         else if (support[1] == cStart) s = 0;
4932:         else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %D does not have cell %D in its support", cone[c], cStart);
4933:         DMGetAuxiliaryVec(dm, key[c].label, key[c].value, key[c].part, &locA[c]);
4934:         if (locA[c]) VecGetDM(locA[c], &dmAux[c]);
4935:         else         {dmAux[c] = dmAux[2];}
4936:         DMGetCellDS(dmAux[c], support[s], &dsAux[c]);
4937:         PetscDSGetTotalDimension(dsAux[c], &totDimAux[c]);
4938:       }
4939:     }
4940:   }
4941:   /* 2: Setup geometric data */
4942:   DMGetCoordinateField(dm, &coordField);
4943:   DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree);
4944:   if (maxDegree > 1) {
4945:     PetscCalloc2(Nf, &quads, Nf, &geoms);
4946:     for (f = 0; f < Nf; ++f) {
4947:       PetscFE fe;

4949:       PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);
4950:       if (fe) {
4951:         PetscFEGetQuadrature(fe, &quads[f]);
4952:         PetscObjectReference((PetscObject) quads[f]);
4953:       }
4954:     }
4955:   }
4956:   /* Loop over chunks */
4957:   cellChunkSize = numCells;
4958:   numChunks     = !numCells ? 0 : PetscCeilReal(((PetscReal) numCells)/cellChunkSize);
4959:   PetscCalloc1(2*cellChunkSize, &faces);
4960:   ISCreateGeneral(PETSC_COMM_SELF, cellChunkSize, faces, PETSC_USE_POINTER, &chunkIS);
4961:   /* Extract field coefficients */
4962:   /* NOTE This needs the end cap faces to have identical orientations */
4963:   DMPlexGetCellFields(dm, cellIS, locX, locX_t, locA[2], &u, &u_t, &a[2]);
4964:   DMPlexGetHybridAuxFields(dm, dmAux, dsAux, cellIS, locA, a);
4965:   DMGetWorkArray(dm, cellChunkSize*totDim, MPIU_SCALAR, &elemVec);
4966:   for (chunk = 0; chunk < numChunks; ++chunk) {
4967:     PetscInt cS = cStart+chunk*cellChunkSize, cE = PetscMin(cS+cellChunkSize, cEnd), numCells = cE - cS, c;

4969:     PetscMemzero(elemVec, cellChunkSize*totDim * sizeof(PetscScalar));
4970:     /* Get faces */
4971:     for (c = cS; c < cE; ++c) {
4972:       const PetscInt  cell = cells ? cells[c] : c;
4973:       const PetscInt *cone;
4974:       DMPlexGetCone(dm, cell, &cone);
4975:       faces[(c-cS)*2+0] = cone[0];
4976:       faces[(c-cS)*2+1] = cone[1];
4977:     }
4978:     ISGeneralSetIndices(chunkIS, cellChunkSize, faces, PETSC_USE_POINTER);
4979:     /* Get geometric data */
4980:     if (maxDegree <= 1) {
4981:       if (!affineQuad) DMFieldCreateDefaultQuadrature(coordField, chunkIS, &affineQuad);
4982:       if (affineQuad)  DMSNESGetFEGeom(coordField, chunkIS, affineQuad, PETSC_TRUE, &affineGeom);
4983:     } else {
4984:       for (f = 0; f < Nf; ++f) {
4985:         if (quads[f]) DMSNESGetFEGeom(coordField, chunkIS, quads[f], PETSC_TRUE, &geoms[f]);
4986:       }
4987:     }
4988:     /* Loop over fields */
4989:     for (f = 0; f < Nf; ++f) {
4990:       PetscFE         fe;
4991:       PetscFEGeom    *geom = affineGeom ? affineGeom : geoms[f];
4992:       PetscFEGeom    *chunkGeom = NULL, *remGeom = NULL;
4993:       PetscQuadrature quad = affineQuad ? affineQuad : quads[f];
4994:       PetscInt        numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset, Nq, Nb;
4995:       PetscBool       isCohesiveField;

4997:       PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);
4998:       if (!fe) continue;
4999:       PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
5000:       PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);
5001:       PetscFEGetDimension(fe, &Nb);
5002:       blockSize = Nb;
5003:       batchSize = numBlocks * blockSize;
5004:       PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
5005:       numChunks = numCells / (numBatches*batchSize);
5006:       Ne        = numChunks*numBatches*batchSize;
5007:       Nr        = numCells % (numBatches*batchSize);
5008:       offset    = numCells - Nr;
5009:       PetscFEGeomGetChunk(geom,0,offset,&chunkGeom);
5010:       PetscFEGeomGetChunk(geom,offset,numCells,&remGeom);
5011:       PetscDSGetCohesive(ds, f, &isCohesiveField);
5012:       chunkGeom->isCohesive = remGeom->isCohesive = PETSC_TRUE;
5013:       key[0].field = f;
5014:       key[1].field = f;
5015:       key[2].field = f;
5016:       PetscFEIntegrateHybridResidual(ds, key[0], 0, Ne, chunkGeom, u, u_t, dsAux[0], a[0], t, elemVec);
5017:       PetscFEIntegrateHybridResidual(ds, key[0], 0, Nr, remGeom,  &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, dsAux[0], &a[0][offset*totDimAux[0]], t, &elemVec[offset*totDim]);
5018:       PetscFEIntegrateHybridResidual(ds, key[1], 1, Ne, chunkGeom, u, u_t, dsAux[1], a[1], t, elemVec);
5019:       PetscFEIntegrateHybridResidual(ds, key[1], 1, Nr, remGeom,  &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, dsAux[1], &a[1][offset*totDimAux[1]], t, &elemVec[offset*totDim]);
5020:       PetscFEIntegrateHybridResidual(ds, key[2], 2, Ne, chunkGeom, u, u_t, dsAux[2], a[2], t, elemVec);
5021:       PetscFEIntegrateHybridResidual(ds, key[2], 2, Nr, remGeom,  &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, dsAux[2], &a[2][offset*totDimAux[2]], t, &elemVec[offset*totDim]);
5022:       PetscFEGeomRestoreChunk(geom,offset,numCells,&remGeom);
5023:       PetscFEGeomRestoreChunk(geom,0,offset,&chunkGeom);
5024:     }
5025:     /* Add elemVec to locX */
5026:     for (c = cS; c < cE; ++c) {
5027:       const PetscInt cell = cells ? cells[c] : c;
5028:       const PetscInt cind = c - cStart;

5030:       if (mesh->printFEM > 1) DMPrintCellVector(cell, name, totDim, &elemVec[cind*totDim]);
5031:       if (ghostLabel) {
5032:         PetscInt ghostVal;

5034:         DMLabelGetValue(ghostLabel,cell,&ghostVal);
5035:         if (ghostVal > 0) continue;
5036:       }
5037:       DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cind*totDim], ADD_ALL_VALUES);
5038:     }
5039:   }
5040:   DMPlexRestoreCellFields(dm, cellIS, locX, locX_t, locA[2], &u, &u_t, &a[2]);
5041:   DMPlexRestoreHybridAuxFields(dmAux, dsAux, cellIS, locA, a);
5042:   DMRestoreWorkArray(dm, numCells*totDim, MPIU_SCALAR, &elemVec);
5043:   PetscFree(faces);
5044:   ISDestroy(&chunkIS);
5045:   ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);
5046:   if (maxDegree <= 1) {
5047:     DMSNESRestoreFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);
5048:     PetscQuadratureDestroy(&affineQuad);
5049:   } else {
5050:     for (f = 0; f < Nf; ++f) {
5051:       if (geoms) DMSNESRestoreFEGeom(coordField,cellIS,quads[f],PETSC_FALSE,&geoms[f]);
5052:       if (quads) PetscQuadratureDestroy(&quads[f]);
5053:     }
5054:     PetscFree2(quads,geoms);
5055:   }
5056:   PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);
5057:   return 0;
5058: }

5060: PetscErrorCode DMPlexComputeBdJacobian_Single_Internal(DM dm, PetscReal t, PetscWeakForm wf, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt fieldI, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, DMField coordField, IS facetIS)
5061: {
5062:   DM_Plex        *mesh = (DM_Plex *) dm->data;
5063:   DM              plex = NULL, plexA = NULL, tdm;
5064:   DMEnclosureType encAux;
5065:   PetscDS         prob, probAux = NULL;
5066:   PetscSection    section, sectionAux = NULL;
5067:   PetscSection    globalSection;
5068:   Vec             locA = NULL, tv;
5069:   PetscScalar    *u = NULL, *u_t = NULL, *a = NULL, *elemMat = NULL;
5070:   PetscInt        v;
5071:   PetscInt        Nf, totDim, totDimAux = 0;
5072:   PetscBool       transform;

5074:   DMConvert(dm, DMPLEX, &plex);
5075:   DMHasBasisTransform(dm, &transform);
5076:   DMGetBasisTransformDM_Internal(dm, &tdm);
5077:   DMGetBasisTransformVec_Internal(dm, &tv);
5078:   DMGetLocalSection(dm, &section);
5079:   DMGetDS(dm, &prob);
5080:   PetscDSGetNumFields(prob, &Nf);
5081:   PetscDSGetTotalDimension(prob, &totDim);
5082:   DMGetAuxiliaryVec(dm, label, values[0], 0, &locA);
5083:   if (locA) {
5084:     DM dmAux;

5086:     VecGetDM(locA, &dmAux);
5087:     DMGetEnclosureRelation(dmAux, dm, &encAux);
5088:     DMConvert(dmAux, DMPLEX, &plexA);
5089:     DMGetDS(plexA, &probAux);
5090:     PetscDSGetTotalDimension(probAux, &totDimAux);
5091:     DMGetLocalSection(plexA, &sectionAux);
5092:   }

5094:   DMGetGlobalSection(dm, &globalSection);
5095:   for (v = 0; v < numValues; ++v) {
5096:     PetscFEGeom     *fgeom;
5097:     PetscInt         maxDegree;
5098:     PetscQuadrature  qGeom = NULL;
5099:     IS               pointIS;
5100:     const PetscInt  *points;
5101:     PetscFormKey key;
5102:     PetscInt         numFaces, face, Nq;

5104:     key.label = label;
5105:     key.value = values[v];
5106:     key.part  = 0;
5107:     DMLabelGetStratumIS(label, values[v], &pointIS);
5108:     if (!pointIS) continue; /* No points with that id on this process */
5109:     {
5110:       IS isectIS;

5112:       /* TODO: Special cases of ISIntersect where it is quick to check a prior if one is a superset of the other */
5113:       ISIntersect_Caching_Internal(facetIS,pointIS,&isectIS);
5114:       ISDestroy(&pointIS);
5115:       pointIS = isectIS;
5116:     }
5117:     ISGetLocalSize(pointIS, &numFaces);
5118:     ISGetIndices(pointIS, &points);
5119:     PetscMalloc4(numFaces*totDim, &u, locX_t ? numFaces*totDim : 0, &u_t, numFaces*totDim*totDim, &elemMat, locA ? numFaces*totDimAux : 0, &a);
5120:     DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);
5121:     if (maxDegree <= 1) {
5122:       DMFieldCreateDefaultQuadrature(coordField,pointIS,&qGeom);
5123:     }
5124:     if (!qGeom) {
5125:       PetscFE fe;

5127:       PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
5128:       PetscFEGetFaceQuadrature(fe, &qGeom);
5129:       PetscObjectReference((PetscObject)qGeom);
5130:     }
5131:     PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
5132:     DMSNESGetFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);
5133:     for (face = 0; face < numFaces; ++face) {
5134:       const PetscInt point = points[face], *support;
5135:       PetscScalar   *x     = NULL;
5136:       PetscInt       i;

5138:       DMPlexGetSupport(dm, point, &support);
5139:       DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x);
5140:       for (i = 0; i < totDim; ++i) u[face*totDim+i] = x[i];
5141:       DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x);
5142:       if (locX_t) {
5143:         DMPlexVecGetClosure(plex, section, locX_t, support[0], NULL, &x);
5144:         for (i = 0; i < totDim; ++i) u_t[face*totDim+i] = x[i];
5145:         DMPlexVecRestoreClosure(plex, section, locX_t, support[0], NULL, &x);
5146:       }
5147:       if (locA) {
5148:         PetscInt subp;
5149:         DMGetEnclosurePoint(plexA, dm, encAux, support[0], &subp);
5150:         DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x);
5151:         for (i = 0; i < totDimAux; ++i) a[face*totDimAux+i] = x[i];
5152:         DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x);
5153:       }
5154:     }
5155:     PetscArrayzero(elemMat, numFaces*totDim*totDim);
5156:     {
5157:       PetscFE         fe;
5158:       PetscInt        Nb;
5159:       /* Conforming batches */
5160:       PetscInt        numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
5161:       /* Remainder */
5162:       PetscFEGeom    *chunkGeom = NULL;
5163:       PetscInt        fieldJ, Nr, offset;

5165:       PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
5166:       PetscFEGetDimension(fe, &Nb);
5167:       PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
5168:       blockSize = Nb;
5169:       batchSize = numBlocks * blockSize;
5170:       PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
5171:       numChunks = numFaces / (numBatches*batchSize);
5172:       Ne        = numChunks*numBatches*batchSize;
5173:       Nr        = numFaces % (numBatches*batchSize);
5174:       offset    = numFaces - Nr;
5175:       PetscFEGeomGetChunk(fgeom,0,offset,&chunkGeom);
5176:       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
5177:         key.field = fieldI*Nf+fieldJ;
5178:         PetscFEIntegrateBdJacobian(prob, wf, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);
5179:       }
5180:       PetscFEGeomGetChunk(fgeom,offset,numFaces,&chunkGeom);
5181:       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
5182:         key.field = fieldI*Nf+fieldJ;
5183:         PetscFEIntegrateBdJacobian(prob, wf, key, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, a ? &a[offset*totDimAux] : NULL, t, X_tShift, &elemMat[offset*totDim*totDim]);
5184:       }
5185:       PetscFEGeomRestoreChunk(fgeom,offset,numFaces,&chunkGeom);
5186:     }
5187:     for (face = 0; face < numFaces; ++face) {
5188:       const PetscInt point = points[face], *support;

5190:       /* Transform to global basis before insertion in Jacobian */
5191:       DMPlexGetSupport(plex, point, &support);
5192:       if (transform) DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, support[0], PETSC_TRUE, totDim, &elemMat[face*totDim*totDim]);
5193:       if (mesh->printFEM > 1) DMPrintCellMatrix(point, "BdJacobian", totDim, totDim, &elemMat[face*totDim*totDim]);
5194:       DMPlexMatSetClosure(plex, section, globalSection, JacP, support[0], &elemMat[face*totDim*totDim], ADD_VALUES);
5195:     }
5196:     DMSNESRestoreFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);
5197:     PetscQuadratureDestroy(&qGeom);
5198:     ISRestoreIndices(pointIS, &points);
5199:     ISDestroy(&pointIS);
5200:     PetscFree4(u, u_t, elemMat, a);
5201:   }
5202:   if (plex)  DMDestroy(&plex);
5203:   if (plexA) DMDestroy(&plexA);
5204:   return 0;
5205: }

5207: PetscErrorCode DMPlexComputeBdJacobianSingle(DM dm, PetscReal t, PetscWeakForm wf, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt field, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP)
5208: {
5209:   DMField        coordField;
5210:   DMLabel        depthLabel;
5211:   IS             facetIS;
5212:   PetscInt       dim;

5214:   DMGetDimension(dm, &dim);
5215:   DMPlexGetDepthLabel(dm, &depthLabel);
5216:   DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);
5217:   DMGetCoordinateField(dm, &coordField);
5218:   DMPlexComputeBdJacobian_Single_Internal(dm, t, wf, label, numValues, values, field, locX, locX_t, X_tShift, Jac, JacP, coordField, facetIS);
5219:   ISDestroy(&facetIS);
5220:   return 0;
5221: }

5223: PetscErrorCode DMPlexComputeBdJacobian_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
5224: {
5225:   PetscDS          prob;
5226:   PetscInt         dim, numBd, bd;
5227:   DMLabel          depthLabel;
5228:   DMField          coordField = NULL;
5229:   IS               facetIS;

5231:   DMGetDS(dm, &prob);
5232:   DMPlexGetDepthLabel(dm, &depthLabel);
5233:   DMGetDimension(dm, &dim);
5234:   DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);
5235:   PetscDSGetNumBoundary(prob, &numBd);
5236:   DMGetCoordinateField(dm, &coordField);
5237:   for (bd = 0; bd < numBd; ++bd) {
5238:     PetscWeakForm           wf;
5239:     DMBoundaryConditionType type;
5240:     DMLabel                 label;
5241:     const PetscInt         *values;
5242:     PetscInt                fieldI, numValues;
5243:     PetscObject             obj;
5244:     PetscClassId            id;

5246:     PetscDSGetBoundary(prob, bd, &wf, &type, NULL, &label, &numValues, &values, &fieldI, NULL, NULL, NULL, NULL, NULL);
5247:     PetscDSGetDiscretization(prob, fieldI, &obj);
5248:     PetscObjectGetClassId(obj, &id);
5249:     if ((id != PETSCFE_CLASSID) || (type & DM_BC_ESSENTIAL)) continue;
5250:     DMPlexComputeBdJacobian_Single_Internal(dm, t, wf, label, numValues, values, fieldI, locX, locX_t, X_tShift, Jac, JacP, coordField, facetIS);
5251:   }
5252:   ISDestroy(&facetIS);
5253:   return 0;
5254: }

5256: PetscErrorCode DMPlexComputeJacobian_Internal(DM dm, PetscFormKey key, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Mat Jac, Mat JacP,void *user)
5257: {
5258:   DM_Plex        *mesh  = (DM_Plex *) dm->data;
5259:   const char     *name  = "Jacobian";
5260:   DM              dmAux = NULL, plex, tdm;
5261:   DMEnclosureType encAux;
5262:   Vec             A, tv;
5263:   DMField         coordField;
5264:   PetscDS         prob, probAux = NULL;
5265:   PetscSection    section, globalSection, sectionAux;
5266:   PetscScalar    *elemMat, *elemMatP, *elemMatD, *u, *u_t, *a = NULL;
5267:   const PetscInt *cells;
5268:   PetscInt        Nf, fieldI, fieldJ;
5269:   PetscInt        totDim, totDimAux, cStart, cEnd, numCells, c;
5270:   PetscBool       hasJac, hasPrec, hasDyn, hasFV = PETSC_FALSE, transform;

5272:   PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);
5273:   ISGetLocalSize(cellIS, &numCells);
5274:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
5275:   DMHasBasisTransform(dm, &transform);
5276:   DMGetBasisTransformDM_Internal(dm, &tdm);
5277:   DMGetBasisTransformVec_Internal(dm, &tv);
5278:   DMGetLocalSection(dm, &section);
5279:   DMGetGlobalSection(dm, &globalSection);
5280:   DMGetCellDS(dm, cells ? cells[cStart] : cStart, &prob);
5281:   PetscDSGetNumFields(prob, &Nf);
5282:   PetscDSGetTotalDimension(prob, &totDim);
5283:   PetscDSHasJacobian(prob, &hasJac);
5284:   PetscDSHasJacobianPreconditioner(prob, &hasPrec);
5285:   /* user passed in the same matrix, avoid double contributions and
5286:      only assemble the Jacobian */
5287:   if (hasJac && Jac == JacP) hasPrec = PETSC_FALSE;
5288:   PetscDSHasDynamicJacobian(prob, &hasDyn);
5289:   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
5290:   DMGetAuxiliaryVec(dm, key.label, key.value, key.part, &A);
5291:   if (A) {
5292:     VecGetDM(A, &dmAux);
5293:     DMGetEnclosureRelation(dmAux, dm, &encAux);
5294:     DMConvert(dmAux, DMPLEX, &plex);
5295:     DMGetLocalSection(plex, &sectionAux);
5296:     DMGetDS(dmAux, &probAux);
5297:     PetscDSGetTotalDimension(probAux, &totDimAux);
5298:   }
5299:   PetscMalloc5(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,hasJac ? numCells*totDim*totDim : 0,&elemMat,hasPrec ? numCells*totDim*totDim : 0, &elemMatP,hasDyn ? numCells*totDim*totDim : 0, &elemMatD);
5300:   if (dmAux) PetscMalloc1(numCells*totDimAux, &a);
5301:   DMGetCoordinateField(dm, &coordField);
5302:   for (c = cStart; c < cEnd; ++c) {
5303:     const PetscInt cell = cells ? cells[c] : c;
5304:     const PetscInt cind = c - cStart;
5305:     PetscScalar   *x = NULL,  *x_t = NULL;
5306:     PetscInt       i;

5308:     DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);
5309:     for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i];
5310:     DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);
5311:     if (X_t) {
5312:       DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);
5313:       for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i];
5314:       DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);
5315:     }
5316:     if (dmAux) {
5317:       PetscInt subcell;
5318:       DMGetEnclosurePoint(dmAux, dm, encAux, cell, &subcell);
5319:       DMPlexVecGetClosure(plex, sectionAux, A, subcell, NULL, &x);
5320:       for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i];
5321:       DMPlexVecRestoreClosure(plex, sectionAux, A, subcell, NULL, &x);
5322:     }
5323:   }
5324:   if (hasJac)  PetscArrayzero(elemMat,  numCells*totDim*totDim);
5325:   if (hasPrec) PetscArrayzero(elemMatP, numCells*totDim*totDim);
5326:   if (hasDyn)  PetscArrayzero(elemMatD, numCells*totDim*totDim);
5327:   for (fieldI = 0; fieldI < Nf; ++fieldI) {
5328:     PetscClassId    id;
5329:     PetscFE         fe;
5330:     PetscQuadrature qGeom = NULL;
5331:     PetscInt        Nb;
5332:     /* Conforming batches */
5333:     PetscInt        numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
5334:     /* Remainder */
5335:     PetscInt        Nr, offset, Nq;
5336:     PetscInt        maxDegree;
5337:     PetscFEGeom     *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;

5339:     PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
5340:     PetscObjectGetClassId((PetscObject) fe, &id);
5341:     if (id == PETSCFV_CLASSID) {hasFV = PETSC_TRUE; continue;}
5342:     PetscFEGetDimension(fe, &Nb);
5343:     PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
5344:     DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);
5345:     if (maxDegree <= 1) {
5346:       DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);
5347:     }
5348:     if (!qGeom) {
5349:       PetscFEGetQuadrature(fe,&qGeom);
5350:       PetscObjectReference((PetscObject)qGeom);
5351:     }
5352:     PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
5353:     DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
5354:     blockSize = Nb;
5355:     batchSize = numBlocks * blockSize;
5356:     PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
5357:     numChunks = numCells / (numBatches*batchSize);
5358:     Ne        = numChunks*numBatches*batchSize;
5359:     Nr        = numCells % (numBatches*batchSize);
5360:     offset    = numCells - Nr;
5361:     PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);
5362:     PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&remGeom);
5363:     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
5364:       key.field = fieldI*Nf+fieldJ;
5365:       if (hasJac) {
5366:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);
5367:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, key, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMat[offset*totDim*totDim]);
5368:       }
5369:       if (hasPrec) {
5370:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_PRE, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatP);
5371:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_PRE, key, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatP[offset*totDim*totDim]);
5372:       }
5373:       if (hasDyn) {
5374:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD);
5375:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, key, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatD[offset*totDim*totDim]);
5376:       }
5377:     }
5378:     PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&remGeom);
5379:     PetscFEGeomRestoreChunk(cgeomFEM,0,offset,&chunkGeom);
5380:     DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
5381:     PetscQuadratureDestroy(&qGeom);
5382:   }
5383:   /*   Add contribution from X_t */
5384:   if (hasDyn) {for (c = 0; c < numCells*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];}
5385:   if (hasFV) {
5386:     PetscClassId id;
5387:     PetscFV      fv;
5388:     PetscInt     offsetI, NcI, NbI = 1, fc, f;

5390:     for (fieldI = 0; fieldI < Nf; ++fieldI) {
5391:       PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fv);
5392:       PetscDSGetFieldOffset(prob, fieldI, &offsetI);
5393:       PetscObjectGetClassId((PetscObject) fv, &id);
5394:       if (id != PETSCFV_CLASSID) continue;
5395:       /* Put in the identity */
5396:       PetscFVGetNumComponents(fv, &NcI);
5397:       for (c = cStart; c < cEnd; ++c) {
5398:         const PetscInt cind    = c - cStart;
5399:         const PetscInt eOffset = cind*totDim*totDim;
5400:         for (fc = 0; fc < NcI; ++fc) {
5401:           for (f = 0; f < NbI; ++f) {
5402:             const PetscInt i = offsetI + f*NcI+fc;
5403:             if (hasPrec) {
5404:               if (hasJac) {elemMat[eOffset+i*totDim+i] = 1.0;}
5405:               elemMatP[eOffset+i*totDim+i] = 1.0;
5406:             } else {elemMat[eOffset+i*totDim+i] = 1.0;}
5407:           }
5408:         }
5409:       }
5410:     }
5411:     /* No allocated space for FV stuff, so ignore the zero entries */
5412:     MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);
5413:   }
5414:   /* Insert values into matrix */
5415:   for (c = cStart; c < cEnd; ++c) {
5416:     const PetscInt cell = cells ? cells[c] : c;
5417:     const PetscInt cind = c - cStart;

5419:     /* Transform to global basis before insertion in Jacobian */
5420:     if (transform) DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, cell, PETSC_TRUE, totDim, &elemMat[cind*totDim*totDim]);
5421:     if (hasPrec) {
5422:       if (hasJac) {
5423:         if (mesh->printFEM > 1) DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);
5424:         DMPlexMatSetClosure(dm, section, globalSection, Jac, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);
5425:       }
5426:       if (mesh->printFEM > 1) DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatP[cind*totDim*totDim]);
5427:       DMPlexMatSetClosure(dm, section, globalSection, JacP, cell, &elemMatP[cind*totDim*totDim], ADD_VALUES);
5428:     } else {
5429:       if (hasJac) {
5430:         if (mesh->printFEM > 1) DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);
5431:         DMPlexMatSetClosure(dm, section, globalSection, JacP, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);
5432:       }
5433:     }
5434:   }
5435:   ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);
5436:   if (hasFV) MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE);
5437:   PetscFree5(u,u_t,elemMat,elemMatP,elemMatD);
5438:   if (dmAux) {
5439:     PetscFree(a);
5440:     DMDestroy(&plex);
5441:   }
5442:   /* Compute boundary integrals */
5443:   DMPlexComputeBdJacobian_Internal(dm, X, X_t, t, X_tShift, Jac, JacP, user);
5444:   /* Assemble matrix */
5445:   if (hasJac && hasPrec) {
5446:     MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY);
5447:     MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY);
5448:   }
5449:   MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);
5450:   MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);
5451:   PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);
5452:   return 0;
5453: }

5455: PetscErrorCode DMPlexComputeJacobian_Hybrid_Internal(DM dm, PetscFormKey key[], IS cellIS, PetscReal t, PetscReal X_tShift, Vec locX, Vec locX_t, Mat Jac, Mat JacP, void *user)
5456: {
5457:   DM_Plex         *mesh          = (DM_Plex *) dm->data;
5458:   const char      *name          = "Hybrid Jacobian";
5459:   DM               dmAux[3]      = {NULL, NULL, NULL};
5460:   DMLabel          ghostLabel    = NULL;
5461:   DM               plex          = NULL;
5462:   DM               plexA         = NULL;
5463:   PetscDS          ds            = NULL;
5464:   PetscDS          dsAux[3]      = {NULL, NULL, NULL};
5465:   Vec              locA[3]       = {NULL, NULL, NULL};
5466:   PetscSection     section       = NULL;
5467:   PetscSection     sectionAux[3] = {NULL, NULL, NULL};
5468:   DMField          coordField    = NULL;
5469:   PetscScalar     *u = NULL, *u_t, *a[3];
5470:   PetscScalar     *elemMat, *elemMatP;
5471:   PetscSection     globalSection;
5472:   IS               chunkIS;
5473:   const PetscInt  *cells;
5474:   PetscInt        *faces;
5475:   PetscInt         cStart, cEnd, numCells;
5476:   PetscInt         Nf, fieldI, fieldJ, totDim, totDimAux[3], numChunks, cellChunkSize, chunk;
5477:   PetscInt         maxDegree = PETSC_MAX_INT;
5478:   PetscQuadrature  affineQuad = NULL, *quads = NULL;
5479:   PetscFEGeom     *affineGeom = NULL, **geoms = NULL;
5480:   PetscBool        hasBdJac, hasBdPrec;

5482:   if ((key[0].label == key[1].label) && (key[0].value == key[1].value) && (key[0].part == key[1].part)) {
5483:     const char *name;
5484:     PetscObjectGetName((PetscObject) key[0].label, &name);
5485:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Form keys for each side of a cohesive surface must be different (%s, %D, %D)", name, key[0].value, key[0].part);
5486:   }
5487:   PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);
5488:   ISGetLocalSize(cellIS, &numCells);
5489:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
5490:   DMConvert(dm, DMPLEX, &plex);
5491:   DMGetSection(dm, &section);
5492:   DMGetGlobalSection(dm, &globalSection);
5493:   DMGetLabel(dm, "ghost", &ghostLabel);
5494:   DMGetCellDS(dm, cStart, &ds);
5495:   PetscDSGetNumFields(ds, &Nf);
5496:   PetscDSGetTotalDimension(ds, &totDim);
5497:   PetscDSHasBdJacobian(ds, &hasBdJac);
5498:   PetscDSHasBdJacobianPreconditioner(ds, &hasBdPrec);
5499:   DMGetAuxiliaryVec(dm, key[2].label, key[2].value, key[2].part, &locA[2]);
5500:   if (locA[2]) {
5501:     VecGetDM(locA[2], &dmAux[2]);
5502:     DMConvert(dmAux[2], DMPLEX, &plexA);
5503:     DMGetSection(dmAux[2], &sectionAux[2]);
5504:     DMGetCellDS(dmAux[2], cStart, &dsAux[2]);
5505:     PetscDSGetTotalDimension(dsAux[2], &totDimAux[2]);
5506:     {
5507:       const PetscInt *cone;
5508:       PetscInt        c;

5510:       DMPlexGetCone(dm, cStart, &cone);
5511:       for (c = 0; c < 2; ++c) {
5512:         const PetscInt *support;
5513:         PetscInt ssize, s;

5515:         DMPlexGetSupport(dm, cone[c], &support);
5516:         DMPlexGetSupportSize(dm, cone[c], &ssize);
5518:         if      (support[0] == cStart) s = 1;
5519:         else if (support[1] == cStart) s = 0;
5520:         else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %D does not have cell %D in its support", cone[c], cStart);
5521:         DMGetAuxiliaryVec(dm, key[c].label, key[c].value, key[c].part, &locA[c]);
5522:         if (locA[c]) VecGetDM(locA[c], &dmAux[c]);
5523:         else         {dmAux[c] = dmAux[2];}
5524:         DMGetCellDS(dmAux[c], support[s], &dsAux[c]);
5525:         PetscDSGetTotalDimension(dsAux[c], &totDimAux[c]);
5526:       }
5527:     }
5528:   }
5529:   DMGetCoordinateField(dm, &coordField);
5530:   DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree);
5531:   if (maxDegree > 1) {
5532:     PetscInt f;
5533:     PetscCalloc2(Nf, &quads, Nf, &geoms);
5534:     for (f = 0; f < Nf; ++f) {
5535:       PetscFE fe;

5537:       PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);
5538:       if (fe) {
5539:         PetscFEGetQuadrature(fe, &quads[f]);
5540:         PetscObjectReference((PetscObject) quads[f]);
5541:       }
5542:     }
5543:   }
5544:   cellChunkSize = numCells;
5545:   numChunks     = !numCells ? 0 : PetscCeilReal(((PetscReal) numCells)/cellChunkSize);
5546:   PetscCalloc1(2*cellChunkSize, &faces);
5547:   ISCreateGeneral(PETSC_COMM_SELF, cellChunkSize, faces, PETSC_USE_POINTER, &chunkIS);
5548:   DMPlexGetCellFields(dm, cellIS, locX, locX_t, locA[2], &u, &u_t, &a[2]);
5549:   DMPlexGetHybridAuxFields(dm, dmAux, dsAux, cellIS, locA, a);
5550:   DMGetWorkArray(dm, hasBdJac  ? cellChunkSize*totDim*totDim : 0, MPIU_SCALAR, &elemMat);
5551:   DMGetWorkArray(dm, hasBdPrec ? cellChunkSize*totDim*totDim : 0, MPIU_SCALAR, &elemMatP);
5552:   for (chunk = 0; chunk < numChunks; ++chunk) {
5553:     PetscInt cS = cStart+chunk*cellChunkSize, cE = PetscMin(cS+cellChunkSize, cEnd), numCells = cE - cS, c;

5555:     if (hasBdJac)  PetscMemzero(elemMat,  numCells*totDim*totDim * sizeof(PetscScalar));
5556:     if (hasBdPrec) PetscMemzero(elemMatP, numCells*totDim*totDim * sizeof(PetscScalar));
5557:     /* Get faces */
5558:     for (c = cS; c < cE; ++c) {
5559:       const PetscInt  cell = cells ? cells[c] : c;
5560:       const PetscInt *cone;
5561:       DMPlexGetCone(plex, cell, &cone);
5562:       faces[(c-cS)*2+0] = cone[0];
5563:       faces[(c-cS)*2+1] = cone[1];
5564:     }
5565:     ISGeneralSetIndices(chunkIS, cellChunkSize, faces, PETSC_USE_POINTER);
5566:     if (maxDegree <= 1) {
5567:       if (!affineQuad) DMFieldCreateDefaultQuadrature(coordField, chunkIS, &affineQuad);
5568:       if (affineQuad)  DMSNESGetFEGeom(coordField, chunkIS, affineQuad, PETSC_TRUE, &affineGeom);
5569:     } else {
5570:       PetscInt f;
5571:       for (f = 0; f < Nf; ++f) {
5572:         if (quads[f]) DMSNESGetFEGeom(coordField, chunkIS, quads[f], PETSC_TRUE, &geoms[f]);
5573:       }
5574:     }

5576:     for (fieldI = 0; fieldI < Nf; ++fieldI) {
5577:       PetscFE         feI;
5578:       PetscFEGeom    *geom = affineGeom ? affineGeom : geoms[fieldI];
5579:       PetscFEGeom    *chunkGeom = NULL, *remGeom = NULL;
5580:       PetscQuadrature quad = affineQuad ? affineQuad : quads[fieldI];
5581:       PetscInt        numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset, Nq, Nb;
5582:       PetscBool       isCohesiveField;

5584:       PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);
5585:       if (!feI) continue;
5586:       PetscFEGetTileSizes(feI, NULL, &numBlocks, NULL, &numBatches);
5587:       PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);
5588:       PetscFEGetDimension(feI, &Nb);
5589:       blockSize = Nb;
5590:       batchSize = numBlocks * blockSize;
5591:       PetscFESetTileSizes(feI, blockSize, numBlocks, batchSize, numBatches);
5592:       numChunks = numCells / (numBatches*batchSize);
5593:       Ne        = numChunks*numBatches*batchSize;
5594:       Nr        = numCells % (numBatches*batchSize);
5595:       offset    = numCells - Nr;
5596:       PetscFEGeomGetChunk(geom,0,offset,&chunkGeom);
5597:       PetscFEGeomGetChunk(geom,offset,numCells,&remGeom);
5598:       PetscDSGetCohesive(ds, fieldI, &isCohesiveField);
5599:       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
5600:         PetscFE feJ;

5602:         PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);
5603:         if (!feJ) continue;
5604:         key[0].field = fieldI*Nf+fieldJ;
5605:         key[1].field = fieldI*Nf+fieldJ;
5606:         key[2].field = fieldI*Nf+fieldJ;
5607:         if (hasBdJac) {
5608:           PetscFEIntegrateHybridJacobian(ds, PETSCFE_JACOBIAN, key[0], 0, Ne, chunkGeom, u, u_t, dsAux[0], a[0], t, X_tShift, elemMat);
5609:           PetscFEIntegrateHybridJacobian(ds, PETSCFE_JACOBIAN, key[0], 0, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, dsAux[0], &a[0][offset*totDimAux[0]], t, X_tShift, &elemMat[offset*totDim*totDim]);
5610:           PetscFEIntegrateHybridJacobian(ds, PETSCFE_JACOBIAN, key[1], 1, Ne, chunkGeom, u, u_t, dsAux[1], a[1], t, X_tShift, elemMat);
5611:           PetscFEIntegrateHybridJacobian(ds, PETSCFE_JACOBIAN, key[1], 1, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, dsAux[1], &a[1][offset*totDimAux[1]], t, X_tShift, &elemMat[offset*totDim*totDim]);
5612:         }
5613:         if (hasBdPrec) {
5614:           PetscFEIntegrateHybridJacobian(ds, PETSCFE_JACOBIAN_PRE, key[0], 0, Ne, chunkGeom, u, u_t, dsAux[0], a[0], t, X_tShift, elemMatP);
5615:           PetscFEIntegrateHybridJacobian(ds, PETSCFE_JACOBIAN_PRE, key[0], 0, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, dsAux[0], &a[0][offset*totDimAux[0]], t, X_tShift, &elemMatP[offset*totDim*totDim]);
5616:           PetscFEIntegrateHybridJacobian(ds, PETSCFE_JACOBIAN_PRE, key[1], 1, Ne, chunkGeom, u, u_t, dsAux[1], a[1], t, X_tShift, elemMatP);
5617:           PetscFEIntegrateHybridJacobian(ds, PETSCFE_JACOBIAN_PRE, key[1], 1, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, dsAux[1], &a[1][offset*totDimAux[1]], t, X_tShift, &elemMatP[offset*totDim*totDim]);
5618:         }
5619:         if (hasBdJac) {
5620:           PetscFEIntegrateHybridJacobian(ds, PETSCFE_JACOBIAN, key[2], 2, Ne, chunkGeom, u, u_t, dsAux[2], a[2], t, X_tShift, elemMat);
5621:           PetscFEIntegrateHybridJacobian(ds, PETSCFE_JACOBIAN, key[2], 2, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, dsAux[2], &a[2][offset*totDimAux[2]], t, X_tShift, &elemMat[offset*totDim*totDim]);
5622:         }
5623:         if (hasBdPrec) {
5624:           PetscFEIntegrateHybridJacobian(ds, PETSCFE_JACOBIAN_PRE, key[2], 2, Ne, chunkGeom, u, u_t, dsAux[2], a[2], t, X_tShift, elemMatP);
5625:           PetscFEIntegrateHybridJacobian(ds, PETSCFE_JACOBIAN_PRE, key[2], 2, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, dsAux[2], &a[2][offset*totDimAux[2]], t, X_tShift, &elemMatP[offset*totDim*totDim]);
5626:         }
5627:       }
5628:       PetscFEGeomRestoreChunk(geom,offset,numCells,&remGeom);
5629:       PetscFEGeomRestoreChunk(geom,0,offset,&chunkGeom);
5630:     }
5631:     /* Insert values into matrix */
5632:     for (c = cS; c < cE; ++c) {
5633:       const PetscInt cell = cells ? cells[c] : c;
5634:       const PetscInt cind = c - cS;

5636:       if (hasBdPrec) {
5637:         if (hasBdJac) {
5638:           if (mesh->printFEM > 1) DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);
5639:           DMPlexMatSetClosure(plex, section, globalSection, Jac, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);
5640:         }
5641:         if (mesh->printFEM > 1) DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatP[cind*totDim*totDim]);
5642:         DMPlexMatSetClosure(plex, section, globalSection, JacP, cell, &elemMatP[cind*totDim*totDim], ADD_VALUES);
5643:       } else if (hasBdJac) {
5644:         if (mesh->printFEM > 1) DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);
5645:         DMPlexMatSetClosure(plex, section, globalSection, JacP, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);
5646:       }
5647:     }
5648:   }
5649:   DMPlexRestoreCellFields(dm, cellIS, locX, locX_t, locA[2], &u, &u_t, &a[2]);
5650:   DMPlexRestoreHybridAuxFields(dmAux, dsAux, cellIS, locA, a);
5651:   DMRestoreWorkArray(dm, hasBdJac  ? cellChunkSize*totDim*totDim : 0, MPIU_SCALAR, &elemMat);
5652:   DMRestoreWorkArray(dm, hasBdPrec ? cellChunkSize*totDim*totDim : 0, MPIU_SCALAR, &elemMatP);
5653:   PetscFree(faces);
5654:   ISDestroy(&chunkIS);
5655:   ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);
5656:   if (maxDegree <= 1) {
5657:     DMSNESRestoreFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);
5658:     PetscQuadratureDestroy(&affineQuad);
5659:   } else {
5660:     PetscInt f;
5661:     for (f = 0; f < Nf; ++f) {
5662:       if (geoms) DMSNESRestoreFEGeom(coordField,cellIS,quads[f],PETSC_FALSE, &geoms[f]);
5663:       if (quads) PetscQuadratureDestroy(&quads[f]);
5664:     }
5665:     PetscFree2(quads,geoms);
5666:   }
5667:   if (dmAux[2]) DMDestroy(&plexA);
5668:   DMDestroy(&plex);
5669:   /* Assemble matrix */
5670:   if (hasBdJac && hasBdPrec) {
5671:     MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY);
5672:     MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY);
5673:   }
5674:   MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);
5675:   MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);
5676:   PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);
5677:   return 0;
5678: }

5680: /*
5681:   DMPlexComputeJacobian_Action_Internal - Form the local portion of the Jacobian action Z = J(X) Y at the local solution X using pointwise functions specified by the user.

5683:   Input Parameters:
5684: + dm     - The mesh
5685: . key    - The PetscWeakFormKey indcating where integration should happen
5686: . cellIS - The cells to integrate over
5687: . t      - The time
5688: . X_tShift - The multiplier for the Jacobian with repsect to X_t
5689: . X      - Local solution vector
5690: . X_t    - Time-derivative of the local solution vector
5691: . Y      - Local input vector
5692: - user   - the user context

5694:   Output Parameter:
5695: . Z - Local output vector

5697:   Note:
5698:   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
5699:   like a GPU, or vectorize on a multicore machine.
5700: */
5701: PetscErrorCode DMPlexComputeJacobian_Action_Internal(DM dm, PetscFormKey key, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Vec Y, Vec Z, void *user)
5702: {
5703:   DM_Plex        *mesh  = (DM_Plex *) dm->data;
5704:   const char     *name  = "Jacobian";
5705:   DM              dmAux = NULL, plex, plexAux = NULL;
5706:   DMEnclosureType encAux;
5707:   Vec             A;
5708:   DMField         coordField;
5709:   PetscDS         prob, probAux = NULL;
5710:   PetscQuadrature quad;
5711:   PetscSection    section, globalSection, sectionAux;
5712:   PetscScalar    *elemMat, *elemMatD, *u, *u_t, *a = NULL, *y, *z;
5713:   const PetscInt *cells;
5714:   PetscInt        Nf, fieldI, fieldJ;
5715:   PetscInt        totDim, totDimAux = 0, cStart, cEnd, numCells, c;
5716:   PetscBool       hasDyn;

5718:   PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);
5719:   DMConvert(dm, DMPLEX, &plex);
5720:   if (!cellIS) {
5721:     PetscInt depth;

5723:     DMPlexGetDepth(plex, &depth);
5724:     DMGetStratumIS(plex, "dim", depth, &cellIS);
5725:     if (!cellIS) DMGetStratumIS(plex, "depth", depth, &cellIS);
5726:   } else {
5727:     PetscObjectReference((PetscObject) cellIS);
5728:   }
5729:   ISGetLocalSize(cellIS, &numCells);
5730:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
5731:   DMGetLocalSection(dm, &section);
5732:   DMGetGlobalSection(dm, &globalSection);
5733:   DMGetCellDS(dm, cells ? cells[cStart] : cStart, &prob);
5734:   PetscDSGetNumFields(prob, &Nf);
5735:   PetscDSGetTotalDimension(prob, &totDim);
5736:   PetscDSHasDynamicJacobian(prob, &hasDyn);
5737:   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
5738:   DMGetAuxiliaryVec(dm, key.label, key.value, key.part, &A);
5739:   if (A) {
5740:     VecGetDM(A, &dmAux);
5741:     DMGetEnclosureRelation(dmAux, dm, &encAux);
5742:     DMConvert(dmAux, DMPLEX, &plexAux);
5743:     DMGetLocalSection(plexAux, &sectionAux);
5744:     DMGetDS(dmAux, &probAux);
5745:     PetscDSGetTotalDimension(probAux, &totDimAux);
5746:   }
5747:   VecSet(Z, 0.0);
5748:   PetscMalloc6(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*totDim*totDim,&elemMat,hasDyn ? numCells*totDim*totDim : 0, &elemMatD,numCells*totDim,&y,totDim,&z);
5749:   if (dmAux) PetscMalloc1(numCells*totDimAux, &a);
5750:   DMGetCoordinateField(dm, &coordField);
5751:   for (c = cStart; c < cEnd; ++c) {
5752:     const PetscInt cell = cells ? cells[c] : c;
5753:     const PetscInt cind = c - cStart;
5754:     PetscScalar   *x = NULL,  *x_t = NULL;
5755:     PetscInt       i;

5757:     DMPlexVecGetClosure(plex, section, X, cell, NULL, &x);
5758:     for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i];
5759:     DMPlexVecRestoreClosure(plex, section, X, cell, NULL, &x);
5760:     if (X_t) {
5761:       DMPlexVecGetClosure(plex, section, X_t, cell, NULL, &x_t);
5762:       for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i];
5763:       DMPlexVecRestoreClosure(plex, section, X_t, cell, NULL, &x_t);
5764:     }
5765:     if (dmAux) {
5766:       PetscInt subcell;
5767:       DMGetEnclosurePoint(dmAux, dm, encAux, cell, &subcell);
5768:       DMPlexVecGetClosure(plexAux, sectionAux, A, subcell, NULL, &x);
5769:       for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i];
5770:       DMPlexVecRestoreClosure(plexAux, sectionAux, A, subcell, NULL, &x);
5771:     }
5772:     DMPlexVecGetClosure(plex, section, Y, cell, NULL, &x);
5773:     for (i = 0; i < totDim; ++i) y[cind*totDim+i] = x[i];
5774:     DMPlexVecRestoreClosure(plex, section, Y, cell, NULL, &x);
5775:   }
5776:   PetscArrayzero(elemMat, numCells*totDim*totDim);
5777:   if (hasDyn)  PetscArrayzero(elemMatD, numCells*totDim*totDim);
5778:   for (fieldI = 0; fieldI < Nf; ++fieldI) {
5779:     PetscFE  fe;
5780:     PetscInt Nb;
5781:     /* Conforming batches */
5782:     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
5783:     /* Remainder */
5784:     PetscInt Nr, offset, Nq;
5785:     PetscQuadrature qGeom = NULL;
5786:     PetscInt    maxDegree;
5787:     PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;

5789:     PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
5790:     PetscFEGetQuadrature(fe, &quad);
5791:     PetscFEGetDimension(fe, &Nb);
5792:     PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
5793:     DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);
5794:     if (maxDegree <= 1) DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);
5795:     if (!qGeom) {
5796:       PetscFEGetQuadrature(fe,&qGeom);
5797:       PetscObjectReference((PetscObject)qGeom);
5798:     }
5799:     PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
5800:     DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
5801:     blockSize = Nb;
5802:     batchSize = numBlocks * blockSize;
5803:     PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
5804:     numChunks = numCells / (numBatches*batchSize);
5805:     Ne        = numChunks*numBatches*batchSize;
5806:     Nr        = numCells % (numBatches*batchSize);
5807:     offset    = numCells - Nr;
5808:     PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);
5809:     PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&remGeom);
5810:     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
5811:       key.field = fieldI*Nf + fieldJ;
5812:       PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);
5813:       PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, key, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMat[offset*totDim*totDim]);
5814:       if (hasDyn) {
5815:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD);
5816:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, key, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatD[offset*totDim*totDim]);
5817:       }
5818:     }
5819:     PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&remGeom);
5820:     PetscFEGeomRestoreChunk(cgeomFEM,0,offset,&chunkGeom);
5821:     DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
5822:     PetscQuadratureDestroy(&qGeom);
5823:   }
5824:   if (hasDyn) {
5825:     for (c = 0; c < numCells*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];
5826:   }
5827:   for (c = cStart; c < cEnd; ++c) {
5828:     const PetscInt     cell = cells ? cells[c] : c;
5829:     const PetscInt     cind = c - cStart;
5830:     const PetscBLASInt M = totDim, one = 1;
5831:     const PetscScalar  a = 1.0, b = 0.0;

5833:     PetscStackCallBLAS("BLASgemv", BLASgemv_("N", &M, &M, &a, &elemMat[cind*totDim*totDim], &M, &y[cind*totDim], &one, &b, z, &one));
5834:     if (mesh->printFEM > 1) {
5835:       DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[cind*totDim*totDim]);
5836:       DMPrintCellVector(c, "Y",  totDim, &y[cind*totDim]);
5837:       DMPrintCellVector(c, "Z",  totDim, z);
5838:     }
5839:     DMPlexVecSetClosure(dm, section, Z, cell, z, ADD_VALUES);
5840:   }
5841:   PetscFree6(u,u_t,elemMat,elemMatD,y,z);
5842:   if (mesh->printFEM) {
5843:     PetscPrintf(PetscObjectComm((PetscObject)Z), "Z:\n");
5844:     VecView(Z, NULL);
5845:   }
5846:   PetscFree(a);
5847:   ISDestroy(&cellIS);
5848:   DMDestroy(&plexAux);
5849:   DMDestroy(&plex);
5850:   PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);
5851:   return 0;
5852: }