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: PetscFunctionBegin;
24: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
25: if (isPlex) {
26: *plex = dm;
27: PetscCall(PetscObjectReference((PetscObject)dm));
28: } else {
29: PetscCall(PetscObjectQuery((PetscObject)dm, "dm_plex", (PetscObject *)plex));
30: if (!*plex) {
31: PetscCall(DMConvert(dm, DMPLEX, plex));
32: PetscCall(PetscObjectCompose((PetscObject)dm, "dm_plex", (PetscObject)*plex));
33: } else {
34: PetscCall(PetscObjectReference((PetscObject)*plex));
35: }
36: if (copy) {
37: DMSubDomainHookLink link;
39: PetscCall(DMCopyDS(dm, PETSC_DETERMINE, PETSC_DETERMINE, *plex));
40: PetscCall(DMCopyAuxiliaryVec(dm, *plex));
41: /* Run the subdomain hook (this will copy the DMSNES/DMTS) */
42: for (link = dm->subdomainhook; link; link = link->next) {
43: if (link->ddhook) PetscCall((*link->ddhook)(dm, *plex, link->ctx));
44: }
45: }
46: }
47: PetscFunctionReturn(PETSC_SUCCESS);
48: }
50: static PetscErrorCode PetscContainerCtxDestroy_PetscFEGeom(void **ctx)
51: {
52: PetscFEGeom *geom = (PetscFEGeom *)*ctx;
54: PetscFunctionBegin;
55: PetscCall(PetscFEGeomDestroy(&geom));
56: PetscFunctionReturn(PETSC_SUCCESS);
57: }
59: static PetscErrorCode DMPlexGetFEGeom(DMField coordField, IS pointIS, PetscQuadrature quad, PetscFEGeomMode mode, PetscFEGeom **geom)
60: {
61: char composeStr[33] = {0};
62: PetscObjectId id;
63: PetscContainer container;
65: PetscFunctionBegin;
66: PetscCall(PetscObjectGetId((PetscObject)quad, &id));
67: PetscCall(PetscSNPrintf(composeStr, 32, "DMPlexGetFEGeom_%" PetscInt64_FMT "\n", id));
68: PetscCall(PetscObjectQuery((PetscObject)pointIS, composeStr, (PetscObject *)&container));
69: if (container) {
70: PetscCall(PetscContainerGetPointer(container, (void **)geom));
71: } else {
72: PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, mode, geom));
73: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
74: PetscCall(PetscContainerSetPointer(container, (void *)*geom));
75: PetscCall(PetscContainerSetCtxDestroy(container, PetscContainerCtxDestroy_PetscFEGeom));
76: PetscCall(PetscObjectCompose((PetscObject)pointIS, composeStr, (PetscObject)container));
77: PetscCall(PetscContainerDestroy(&container));
78: }
79: PetscFunctionReturn(PETSC_SUCCESS);
80: }
82: static PetscErrorCode DMPlexRestoreFEGeom(DMField coordField, IS pointIS, PetscQuadrature quad, PetscFEGeomMode mode, PetscFEGeom **geom)
83: {
84: PetscFunctionBegin;
85: *geom = NULL;
86: PetscFunctionReturn(PETSC_SUCCESS);
87: }
89: /*@
90: DMPlexGetScale - Get the scale for the specified fundamental unit
92: Not Collective
94: Input Parameters:
95: + dm - the `DM`
96: - unit - The SI unit
98: Output Parameter:
99: . scale - The value used to scale all quantities with this unit
101: Level: advanced
103: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetScale()`, `PetscUnit`
104: @*/
105: PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
106: {
107: DM_Plex *mesh = (DM_Plex *)dm->data;
109: PetscFunctionBegin;
111: PetscAssertPointer(scale, 3);
112: *scale = mesh->scale[unit];
113: PetscFunctionReturn(PETSC_SUCCESS);
114: }
116: /*@
117: DMPlexSetScale - Set the scale for the specified fundamental unit
119: Not Collective
121: Input Parameters:
122: + dm - the `DM`
123: . unit - The SI unit
124: - scale - The value used to scale all quantities with this unit
126: Level: advanced
128: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetScale()`, `PetscUnit`
129: @*/
130: PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
131: {
132: DM_Plex *mesh = (DM_Plex *)dm->data;
134: PetscFunctionBegin;
136: mesh->scale[unit] = scale;
137: PetscFunctionReturn(PETSC_SUCCESS);
138: }
140: PetscErrorCode DMPlexGetUseCeed_Plex(DM dm, PetscBool *useCeed)
141: {
142: DM_Plex *mesh = (DM_Plex *)dm->data;
144: PetscFunctionBegin;
145: *useCeed = mesh->useCeed;
146: PetscFunctionReturn(PETSC_SUCCESS);
147: }
148: PetscErrorCode DMPlexSetUseCeed_Plex(DM dm, PetscBool useCeed)
149: {
150: DM_Plex *mesh = (DM_Plex *)dm->data;
152: PetscFunctionBegin;
153: mesh->useCeed = useCeed;
154: PetscFunctionReturn(PETSC_SUCCESS);
155: }
157: /*@
158: DMPlexGetUseCeed - Get flag for using the LibCEED backend
160: Not collective
162: Input Parameter:
163: . dm - The `DM`
165: Output Parameter:
166: . useCeed - The flag
168: Level: intermediate
170: .seealso: `DMPlexSetUseCeed()`
171: @*/
172: PetscErrorCode DMPlexGetUseCeed(DM dm, PetscBool *useCeed)
173: {
174: PetscFunctionBegin;
176: PetscAssertPointer(useCeed, 2);
177: *useCeed = PETSC_FALSE;
178: PetscTryMethod(dm, "DMPlexGetUseCeed_C", (DM, PetscBool *), (dm, useCeed));
179: PetscFunctionReturn(PETSC_SUCCESS);
180: }
182: /*@
183: DMPlexSetUseCeed - Set flag for using the LibCEED backend
185: Not collective
187: Input Parameters:
188: + dm - The `DM`
189: - useCeed - The flag
191: Level: intermediate
193: .seealso: `DMPlexGetUseCeed()`
194: @*/
195: PetscErrorCode DMPlexSetUseCeed(DM dm, PetscBool useCeed)
196: {
197: PetscFunctionBegin;
200: PetscUseMethod(dm, "DMPlexSetUseCeed_C", (DM, PetscBool), (dm, useCeed));
201: PetscFunctionReturn(PETSC_SUCCESS);
202: }
204: /*@
205: DMPlexGetUseMatClosurePermutation - Get flag for using a closure permutation for matrix insertion
207: Not collective
209: Input Parameter:
210: . dm - The `DM`
212: Output Parameter:
213: . useClPerm - The flag
215: Level: intermediate
217: .seealso: `DMPlexSetUseMatClosurePermutation()`
218: @*/
219: PetscErrorCode DMPlexGetUseMatClosurePermutation(DM dm, PetscBool *useClPerm)
220: {
221: DM_Plex *mesh = (DM_Plex *)dm->data;
223: PetscFunctionBegin;
225: PetscAssertPointer(useClPerm, 2);
226: *useClPerm = mesh->useMatClPerm;
227: PetscFunctionReturn(PETSC_SUCCESS);
228: }
230: /*@
231: DMPlexSetUseMatClosurePermutation - Set flag for using a closure permutation for matrix insertion
233: Not collective
235: Input Parameters:
236: + dm - The `DM`
237: - useClPerm - The flag
239: Level: intermediate
241: .seealso: `DMPlexGetUseMatClosurePermutation()`
242: @*/
243: PetscErrorCode DMPlexSetUseMatClosurePermutation(DM dm, PetscBool useClPerm)
244: {
245: DM_Plex *mesh = (DM_Plex *)dm->data;
247: PetscFunctionBegin;
250: mesh->useMatClPerm = useClPerm;
251: PetscFunctionReturn(PETSC_SUCCESS);
252: }
254: static PetscErrorCode DMPlexProjectRigidBody_Private(PetscInt dim, PetscReal t, const PetscReal X[], PetscInt Nc, PetscScalar *mode, void *ctx)
255: {
256: const PetscInt eps[3][3][3] = {
257: {{0, 0, 0}, {0, 0, 1}, {0, -1, 0}},
258: {{0, 0, -1}, {0, 0, 0}, {1, 0, 0} },
259: {{0, 1, 0}, {-1, 0, 0}, {0, 0, 0} }
260: };
261: PetscInt *ctxInt = (PetscInt *)ctx;
262: PetscInt dim2 = ctxInt[0];
263: PetscInt d = ctxInt[1];
264: PetscInt i, j, k = dim > 2 ? d - dim : d;
266: PetscFunctionBegin;
267: PetscCheck(dim == dim2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input dimension %" PetscInt_FMT " does not match context dimension %" PetscInt_FMT, dim, dim2);
268: for (i = 0; i < dim; i++) mode[i] = 0.;
269: if (d < dim) {
270: mode[d] = 1.; /* Translation along axis d */
271: } else {
272: for (i = 0; i < dim; i++) {
273: for (j = 0; j < dim; j++) mode[j] += eps[i][j][k] * X[i]; /* Rotation about axis d */
274: }
275: }
276: PetscFunctionReturn(PETSC_SUCCESS);
277: }
279: /*@
280: DMPlexCreateRigidBody - For the default global section, create rigid body modes by function space interpolation
282: Collective
284: Input Parameters:
285: + dm - the `DM`
286: - field - The field number for the rigid body space, or 0 for the default
288: Output Parameter:
289: . sp - the null space
291: Level: advanced
293: Note:
294: This is necessary to provide a suitable coarse space for algebraic multigrid
296: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `MatNullSpaceCreate()`, `PCGAMG`
297: @*/
298: PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscInt field, MatNullSpace *sp)
299: {
300: PetscErrorCode (**func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *);
301: MPI_Comm comm;
302: Vec mode[6];
303: PetscSection section, globalSection;
304: PetscInt dim, dimEmbed, Nf, n, m, mmin, d, i, j;
305: void **ctxs;
307: PetscFunctionBegin;
308: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
309: PetscCall(DMGetDimension(dm, &dim));
310: PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
311: PetscCall(DMGetNumFields(dm, &Nf));
312: PetscCheck(!Nf || !(field < 0 || field >= Nf), comm, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " is not in [0, %" PetscInt_FMT ")", field, Nf);
313: if (dim == 1 && Nf < 2) {
314: PetscCall(MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp));
315: PetscFunctionReturn(PETSC_SUCCESS);
316: }
317: PetscCall(DMGetLocalSection(dm, §ion));
318: PetscCall(DMGetGlobalSection(dm, &globalSection));
319: PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &n));
320: PetscCall(PetscCalloc2(Nf, &func, Nf, &ctxs));
321: m = (dim * (dim + 1)) / 2;
322: PetscCall(VecCreate(comm, &mode[0]));
323: PetscCall(VecSetType(mode[0], dm->vectype));
324: PetscCall(VecSetSizes(mode[0], n, PETSC_DETERMINE));
325: PetscCall(VecSetUp(mode[0]));
326: PetscCall(VecGetSize(mode[0], &n));
327: mmin = PetscMin(m, n);
328: func[field] = DMPlexProjectRigidBody_Private;
329: for (i = 1; i < m; ++i) PetscCall(VecDuplicate(mode[0], &mode[i]));
330: for (d = 0; d < m; d++) {
331: PetscInt ctx[2];
333: ctxs[field] = (void *)(&ctx[0]);
334: ctx[0] = dimEmbed;
335: ctx[1] = d;
336: PetscCall(DMProjectFunction(dm, 0.0, func, ctxs, INSERT_VALUES, mode[d]));
337: }
338: /* Orthonormalize system */
339: for (i = 0; i < mmin; ++i) {
340: PetscScalar dots[6];
342: PetscCall(VecNormalize(mode[i], NULL));
343: PetscCall(VecMDot(mode[i], mmin - i - 1, mode + i + 1, dots + i + 1));
344: for (j = i + 1; j < mmin; ++j) {
345: dots[j] *= -1.0;
346: PetscCall(VecAXPY(mode[j], dots[j], mode[i]));
347: }
348: }
349: PetscCall(MatNullSpaceCreate(comm, PETSC_FALSE, mmin, mode, sp));
350: for (i = 0; i < m; ++i) PetscCall(VecDestroy(&mode[i]));
351: PetscCall(PetscFree2(func, ctxs));
352: PetscFunctionReturn(PETSC_SUCCESS);
353: }
355: /*@
356: DMPlexCreateRigidBodies - For the default global section, create rigid body modes by function space interpolation
358: Collective
360: Input Parameters:
361: + dm - the `DM`
362: . nb - The number of bodies
363: . label - The `DMLabel` marking each domain
364: . nids - The number of ids per body
365: - ids - An array of the label ids in sequence for each domain
367: Output Parameter:
368: . sp - the null space
370: Level: advanced
372: Note:
373: This is necessary to provide a suitable coarse space for algebraic multigrid
375: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `MatNullSpaceCreate()`
376: @*/
377: PetscErrorCode DMPlexCreateRigidBodies(DM dm, PetscInt nb, DMLabel label, const PetscInt nids[], const PetscInt ids[], MatNullSpace *sp)
378: {
379: MPI_Comm comm;
380: PetscSection section, globalSection;
381: Vec *mode;
382: PetscScalar *dots;
383: PetscInt dim, dimEmbed, n, m, b, d, i, j, off;
385: PetscFunctionBegin;
386: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
387: PetscCall(DMGetDimension(dm, &dim));
388: PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
389: PetscCall(DMGetLocalSection(dm, §ion));
390: PetscCall(DMGetGlobalSection(dm, &globalSection));
391: PetscCall(PetscSectionGetConstrainedStorageSize(globalSection, &n));
392: m = nb * (dim * (dim + 1)) / 2;
393: PetscCall(PetscMalloc2(m, &mode, m, &dots));
394: PetscCall(VecCreate(comm, &mode[0]));
395: PetscCall(VecSetSizes(mode[0], n, PETSC_DETERMINE));
396: PetscCall(VecSetUp(mode[0]));
397: for (i = 1; i < m; ++i) PetscCall(VecDuplicate(mode[0], &mode[i]));
398: for (b = 0, off = 0; b < nb; ++b) {
399: for (d = 0; d < m / nb; ++d) {
400: PetscInt ctx[2];
401: PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *) = DMPlexProjectRigidBody_Private;
402: void *voidctx = (void *)(&ctx[0]);
404: ctx[0] = dimEmbed;
405: ctx[1] = d;
406: PetscCall(DMProjectFunctionLabel(dm, 0.0, label, nids[b], &ids[off], 0, NULL, &func, &voidctx, INSERT_VALUES, mode[d]));
407: off += nids[b];
408: }
409: }
410: /* Orthonormalize system */
411: for (i = 0; i < m; ++i) {
412: PetscScalar dots[6];
414: PetscCall(VecNormalize(mode[i], NULL));
415: PetscCall(VecMDot(mode[i], m - i - 1, mode + i + 1, dots + i + 1));
416: for (j = i + 1; j < m; ++j) {
417: dots[j] *= -1.0;
418: PetscCall(VecAXPY(mode[j], dots[j], mode[i]));
419: }
420: }
421: PetscCall(MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp));
422: for (i = 0; i < m; ++i) PetscCall(VecDestroy(&mode[i]));
423: PetscCall(PetscFree2(mode, dots));
424: PetscFunctionReturn(PETSC_SUCCESS);
425: }
427: /*@
428: DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs
429: are computed by associating the basis function with one of the mesh points in its transitively-closed support, and
430: evaluating the dual space basis of that point.
432: Input Parameters:
433: + dm - the `DMPLEX` object
434: - height - the maximum projection height >= 0
436: Level: advanced
438: Notes:
439: A basis function is associated with the point in its transitively-closed support whose mesh
440: height is highest (w.r.t. DAG height), but not greater than the maximum projection height,
441: which is set with this function. By default, the maximum projection height is zero, which
442: means that only mesh cells are used to project basis functions. A height of one, for
443: example, evaluates a cell-interior basis functions using its cells dual space basis, but all
444: other basis functions with the dual space basis of a face.
446: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetMaxProjectionHeight()`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabelLocal()`
447: @*/
448: PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height)
449: {
450: DM_Plex *plex = (DM_Plex *)dm->data;
452: PetscFunctionBegin;
454: plex->maxProjectionHeight = height;
455: PetscFunctionReturn(PETSC_SUCCESS);
456: }
458: /*@
459: DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in
460: DMPlexProjectXXXLocal() functions.
462: Input Parameter:
463: . dm - the `DMPLEX` object
465: Output Parameter:
466: . height - the maximum projection height
468: Level: intermediate
470: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetMaxProjectionHeight()`, `DMProjectFunctionLocal()`, `DMProjectFunctionLabelLocal()`
471: @*/
472: PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height)
473: {
474: DM_Plex *plex = (DM_Plex *)dm->data;
476: PetscFunctionBegin;
478: *height = plex->maxProjectionHeight;
479: PetscFunctionReturn(PETSC_SUCCESS);
480: }
482: typedef struct {
483: PetscReal alpha; /* The first Euler angle, and in 2D the only one */
484: PetscReal beta; /* The second Euler angle */
485: PetscReal gamma; /* The third Euler angle */
486: PetscInt dim; /* The dimension of R */
487: PetscScalar *R; /* The rotation matrix, transforming a vector in the local basis to the global basis */
488: PetscScalar *RT; /* The transposed rotation matrix, transforming a vector in the global basis to the local basis */
489: } RotCtx;
491: /*
492: Note: Following https://en.wikipedia.org/wiki/Euler_angles, we will specify Euler angles by extrinsic rotations, meaning that
493: 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:
494: $ The XYZ system rotates about the z axis by alpha. The X axis is now at angle alpha with respect to the x axis.
495: $ 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.
496: $ The XYZ system rotates a third time about the z axis by gamma.
497: */
498: static PetscErrorCode DMPlexBasisTransformSetUp_Rotation_Internal(DM dm, void *ctx)
499: {
500: RotCtx *rc = (RotCtx *)ctx;
501: PetscInt dim = rc->dim;
502: PetscReal c1, s1, c2, s2, c3, s3;
504: PetscFunctionBegin;
505: PetscCall(PetscMalloc2(PetscSqr(dim), &rc->R, PetscSqr(dim), &rc->RT));
506: switch (dim) {
507: case 2:
508: c1 = PetscCosReal(rc->alpha);
509: s1 = PetscSinReal(rc->alpha);
510: rc->R[0] = c1;
511: rc->R[1] = s1;
512: rc->R[2] = -s1;
513: rc->R[3] = c1;
514: PetscCall(PetscArraycpy(rc->RT, rc->R, PetscSqr(dim)));
515: DMPlex_Transpose2D_Internal(rc->RT);
516: break;
517: case 3:
518: c1 = PetscCosReal(rc->alpha);
519: s1 = PetscSinReal(rc->alpha);
520: c2 = PetscCosReal(rc->beta);
521: s2 = PetscSinReal(rc->beta);
522: c3 = PetscCosReal(rc->gamma);
523: s3 = PetscSinReal(rc->gamma);
524: rc->R[0] = c1 * c3 - c2 * s1 * s3;
525: rc->R[1] = c3 * s1 + c1 * c2 * s3;
526: rc->R[2] = s2 * s3;
527: rc->R[3] = -c1 * s3 - c2 * c3 * s1;
528: rc->R[4] = c1 * c2 * c3 - s1 * s3;
529: rc->R[5] = c3 * s2;
530: rc->R[6] = s1 * s2;
531: rc->R[7] = -c1 * s2;
532: rc->R[8] = c2;
533: PetscCall(PetscArraycpy(rc->RT, rc->R, PetscSqr(dim)));
534: DMPlex_Transpose3D_Internal(rc->RT);
535: break;
536: default:
537: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported", dim);
538: }
539: PetscFunctionReturn(PETSC_SUCCESS);
540: }
542: static PetscErrorCode DMPlexBasisTransformDestroy_Rotation_Internal(DM dm, void *ctx)
543: {
544: RotCtx *rc = (RotCtx *)ctx;
546: PetscFunctionBegin;
547: PetscCall(PetscFree2(rc->R, rc->RT));
548: PetscCall(PetscFree(rc));
549: PetscFunctionReturn(PETSC_SUCCESS);
550: }
552: static PetscErrorCode DMPlexBasisTransformGetMatrix_Rotation_Internal(DM dm, const PetscReal x[], PetscBool l2g, const PetscScalar **A, void *ctx)
553: {
554: RotCtx *rc = (RotCtx *)ctx;
556: PetscFunctionBeginHot;
557: PetscAssertPointer(ctx, 5);
558: if (l2g) {
559: *A = rc->R;
560: } else {
561: *A = rc->RT;
562: }
563: PetscFunctionReturn(PETSC_SUCCESS);
564: }
566: PetscErrorCode DMPlexBasisTransformApplyReal_Internal(DM dm, const PetscReal x[], PetscBool l2g, PetscInt dim, const PetscReal *y, PetscReal *z, void *ctx)
567: {
568: PetscFunctionBegin;
569: #if defined(PETSC_USE_COMPLEX)
570: switch (dim) {
571: case 2: {
572: PetscScalar yt[2] = {y[0], y[1]}, zt[2] = {0.0, 0.0};
574: PetscCall(DMPlexBasisTransformApply_Internal(dm, x, l2g, dim, yt, zt, ctx));
575: z[0] = PetscRealPart(zt[0]);
576: z[1] = PetscRealPart(zt[1]);
577: } break;
578: case 3: {
579: PetscScalar yt[3] = {y[0], y[1], y[2]}, zt[3] = {0.0, 0.0, 0.0};
581: PetscCall(DMPlexBasisTransformApply_Internal(dm, x, l2g, dim, yt, zt, ctx));
582: z[0] = PetscRealPart(zt[0]);
583: z[1] = PetscRealPart(zt[1]);
584: z[2] = PetscRealPart(zt[2]);
585: } break;
586: }
587: #else
588: PetscCall(DMPlexBasisTransformApply_Internal(dm, x, l2g, dim, y, z, ctx));
589: #endif
590: PetscFunctionReturn(PETSC_SUCCESS);
591: }
593: PetscErrorCode DMPlexBasisTransformApply_Internal(DM dm, const PetscReal x[], PetscBool l2g, PetscInt dim, const PetscScalar *y, PetscScalar *z, void *ctx)
594: {
595: const PetscScalar *A;
597: PetscFunctionBeginHot;
598: PetscCall((*dm->transformGetMatrix)(dm, x, l2g, &A, ctx));
599: switch (dim) {
600: case 2:
601: DMPlex_Mult2D_Internal(A, 1, y, z);
602: break;
603: case 3:
604: DMPlex_Mult3D_Internal(A, 1, y, z);
605: break;
606: }
607: PetscFunctionReturn(PETSC_SUCCESS);
608: }
610: static PetscErrorCode DMPlexBasisTransformField_Internal(DM dm, DM tdm, Vec tv, PetscInt p, PetscInt f, PetscBool l2g, PetscScalar *a)
611: {
612: PetscSection ts;
613: const PetscScalar *ta, *tva;
614: PetscInt dof;
616: PetscFunctionBeginHot;
617: PetscCall(DMGetLocalSection(tdm, &ts));
618: PetscCall(PetscSectionGetFieldDof(ts, p, f, &dof));
619: PetscCall(VecGetArrayRead(tv, &ta));
620: PetscCall(DMPlexPointLocalFieldRead(tdm, p, f, ta, &tva));
621: if (l2g) {
622: switch (dof) {
623: case 4:
624: DMPlex_Mult2D_Internal(tva, 1, a, a);
625: break;
626: case 9:
627: DMPlex_Mult3D_Internal(tva, 1, a, a);
628: break;
629: }
630: } else {
631: switch (dof) {
632: case 4:
633: DMPlex_MultTranspose2D_Internal(tva, 1, a, a);
634: break;
635: case 9:
636: DMPlex_MultTranspose3D_Internal(tva, 1, a, a);
637: break;
638: }
639: }
640: PetscCall(VecRestoreArrayRead(tv, &ta));
641: PetscFunctionReturn(PETSC_SUCCESS);
642: }
644: static PetscErrorCode DMPlexBasisTransformFieldTensor_Internal(DM dm, DM tdm, Vec tv, PetscInt pf, PetscInt f, PetscInt pg, PetscInt g, PetscBool l2g, PetscInt lda, PetscScalar *a)
645: {
646: PetscSection s, ts;
647: const PetscScalar *ta, *tvaf, *tvag;
648: PetscInt fdof, gdof, fpdof, gpdof;
650: PetscFunctionBeginHot;
651: PetscCall(DMGetLocalSection(dm, &s));
652: PetscCall(DMGetLocalSection(tdm, &ts));
653: PetscCall(PetscSectionGetFieldDof(s, pf, f, &fpdof));
654: PetscCall(PetscSectionGetFieldDof(s, pg, g, &gpdof));
655: PetscCall(PetscSectionGetFieldDof(ts, pf, f, &fdof));
656: PetscCall(PetscSectionGetFieldDof(ts, pg, g, &gdof));
657: PetscCall(VecGetArrayRead(tv, &ta));
658: PetscCall(DMPlexPointLocalFieldRead(tdm, pf, f, ta, &tvaf));
659: PetscCall(DMPlexPointLocalFieldRead(tdm, pg, g, ta, &tvag));
660: if (l2g) {
661: switch (fdof) {
662: case 4:
663: DMPlex_MatMult2D_Internal(tvaf, gpdof, lda, a, a);
664: break;
665: case 9:
666: DMPlex_MatMult3D_Internal(tvaf, gpdof, lda, a, a);
667: break;
668: }
669: switch (gdof) {
670: case 4:
671: DMPlex_MatMultTransposeLeft2D_Internal(tvag, fpdof, lda, a, a);
672: break;
673: case 9:
674: DMPlex_MatMultTransposeLeft3D_Internal(tvag, fpdof, lda, a, a);
675: break;
676: }
677: } else {
678: switch (fdof) {
679: case 4:
680: DMPlex_MatMultTranspose2D_Internal(tvaf, gpdof, lda, a, a);
681: break;
682: case 9:
683: DMPlex_MatMultTranspose3D_Internal(tvaf, gpdof, lda, a, a);
684: break;
685: }
686: switch (gdof) {
687: case 4:
688: DMPlex_MatMultLeft2D_Internal(tvag, fpdof, lda, a, a);
689: break;
690: case 9:
691: DMPlex_MatMultLeft3D_Internal(tvag, fpdof, lda, a, a);
692: break;
693: }
694: }
695: PetscCall(VecRestoreArrayRead(tv, &ta));
696: PetscFunctionReturn(PETSC_SUCCESS);
697: }
699: PetscErrorCode DMPlexBasisTransformPoint_Internal(DM dm, DM tdm, Vec tv, PetscInt p, PetscBool fieldActive[], PetscBool l2g, PetscScalar *a)
700: {
701: PetscSection s;
702: PetscSection clSection;
703: IS clPoints;
704: const PetscInt *clp;
705: PetscInt *points = NULL;
706: PetscInt Nf, f, Np, cp, dof, d = 0;
708: PetscFunctionBegin;
709: PetscCall(DMGetLocalSection(dm, &s));
710: PetscCall(PetscSectionGetNumFields(s, &Nf));
711: PetscCall(DMPlexGetCompressedClosure(dm, s, p, 0, &Np, &points, &clSection, &clPoints, &clp));
712: for (f = 0; f < Nf; ++f) {
713: for (cp = 0; cp < Np * 2; cp += 2) {
714: PetscCall(PetscSectionGetFieldDof(s, points[cp], f, &dof));
715: if (!dof) continue;
716: if (fieldActive[f]) PetscCall(DMPlexBasisTransformField_Internal(dm, tdm, tv, points[cp], f, l2g, &a[d]));
717: d += dof;
718: }
719: }
720: PetscCall(DMPlexRestoreCompressedClosure(dm, s, p, &Np, &points, &clSection, &clPoints, &clp));
721: PetscFunctionReturn(PETSC_SUCCESS);
722: }
724: PetscErrorCode DMPlexBasisTransformPointTensor_Internal(DM dm, DM tdm, Vec tv, PetscInt p, PetscBool l2g, PetscInt lda, PetscScalar *a)
725: {
726: PetscSection s;
727: PetscSection clSection;
728: IS clPoints;
729: const PetscInt *clp;
730: PetscInt *points = NULL;
731: PetscInt Nf, f, g, Np, cpf, cpg, fdof, gdof, r, c = 0;
733: PetscFunctionBegin;
734: PetscCall(DMGetLocalSection(dm, &s));
735: PetscCall(PetscSectionGetNumFields(s, &Nf));
736: PetscCall(DMPlexGetCompressedClosure(dm, s, p, 0, &Np, &points, &clSection, &clPoints, &clp));
737: for (f = 0, r = 0; f < Nf; ++f) {
738: for (cpf = 0; cpf < Np * 2; cpf += 2) {
739: PetscCall(PetscSectionGetFieldDof(s, points[cpf], f, &fdof));
740: for (g = 0, c = 0; g < Nf; ++g) {
741: for (cpg = 0; cpg < Np * 2; cpg += 2) {
742: PetscCall(PetscSectionGetFieldDof(s, points[cpg], g, &gdof));
743: PetscCall(DMPlexBasisTransformFieldTensor_Internal(dm, tdm, tv, points[cpf], f, points[cpg], g, l2g, lda, &a[r * lda + c]));
744: c += gdof;
745: }
746: }
747: PetscCheck(c == lda, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of columns %" PetscInt_FMT " should be %" PetscInt_FMT, c, lda);
748: r += fdof;
749: }
750: }
751: PetscCheck(r == lda, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of rows %" PetscInt_FMT " should be %" PetscInt_FMT, c, lda);
752: PetscCall(DMPlexRestoreCompressedClosure(dm, s, p, &Np, &points, &clSection, &clPoints, &clp));
753: PetscFunctionReturn(PETSC_SUCCESS);
754: }
756: static PetscErrorCode DMPlexBasisTransform_Internal(DM dm, Vec lv, PetscBool l2g)
757: {
758: DM tdm;
759: Vec tv;
760: PetscSection ts, s;
761: const PetscScalar *ta;
762: PetscScalar *a, *va;
763: PetscInt pStart, pEnd, p, Nf, f;
765: PetscFunctionBegin;
766: PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
767: PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
768: PetscCall(DMGetLocalSection(tdm, &ts));
769: PetscCall(DMGetLocalSection(dm, &s));
770: PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
771: PetscCall(PetscSectionGetNumFields(s, &Nf));
772: PetscCall(VecGetArray(lv, &a));
773: PetscCall(VecGetArrayRead(tv, &ta));
774: for (p = pStart; p < pEnd; ++p) {
775: for (f = 0; f < Nf; ++f) {
776: PetscCall(DMPlexPointLocalFieldRef(dm, p, f, a, &va));
777: PetscCall(DMPlexBasisTransformField_Internal(dm, tdm, tv, p, f, l2g, va));
778: }
779: }
780: PetscCall(VecRestoreArray(lv, &a));
781: PetscCall(VecRestoreArrayRead(tv, &ta));
782: PetscFunctionReturn(PETSC_SUCCESS);
783: }
785: /*@
786: DMPlexGlobalToLocalBasis - Transform the values in the given local vector from the global basis to the local basis
788: Input Parameters:
789: + dm - The `DM`
790: - lv - A local vector with values in the global basis
792: Output Parameter:
793: . lv - A local vector with values in the local basis
795: Level: developer
797: Note:
798: 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.
800: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexLocalToGlobalBasis()`, `DMGetLocalSection()`, `DMPlexCreateBasisRotation()`
801: @*/
802: PetscErrorCode DMPlexGlobalToLocalBasis(DM dm, Vec lv)
803: {
804: PetscFunctionBegin;
807: PetscCall(DMPlexBasisTransform_Internal(dm, lv, PETSC_FALSE));
808: PetscFunctionReturn(PETSC_SUCCESS);
809: }
811: /*@
812: DMPlexLocalToGlobalBasis - Transform the values in the given local vector from the local basis to the global basis
814: Input Parameters:
815: + dm - The `DM`
816: - lv - A local vector with values in the local basis
818: Output Parameter:
819: . lv - A local vector with values in the global basis
821: Level: developer
823: Note:
824: 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.
826: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGlobalToLocalBasis()`, `DMGetLocalSection()`, `DMPlexCreateBasisRotation()`
827: @*/
828: PetscErrorCode DMPlexLocalToGlobalBasis(DM dm, Vec lv)
829: {
830: PetscFunctionBegin;
833: PetscCall(DMPlexBasisTransform_Internal(dm, lv, PETSC_TRUE));
834: PetscFunctionReturn(PETSC_SUCCESS);
835: }
837: /*@
838: DMPlexCreateBasisRotation - Create an internal transformation from the global basis, used to specify boundary conditions
839: and global solutions, to a local basis, appropriate for discretization integrals and assembly.
841: Input Parameters:
842: + dm - The `DM`
843: . alpha - The first Euler angle, and in 2D the only one
844: . beta - The second Euler angle
845: - gamma - The third Euler angle
847: Level: developer
849: Note:
850: Following https://en.wikipedia.org/wiki/Euler_angles, we will specify Euler angles by extrinsic rotations, meaning that
851: 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
852: .vb
853: The XYZ system rotates about the z axis by alpha. The X axis is now at angle alpha with respect to the x axis.
854: 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.
855: The XYZ system rotates a third time about the z axis by gamma.
856: .ve
858: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGlobalToLocalBasis()`, `DMPlexLocalToGlobalBasis()`
859: @*/
860: PetscErrorCode DMPlexCreateBasisRotation(DM dm, PetscReal alpha, PetscReal beta, PetscReal gamma)
861: {
862: RotCtx *rc;
863: PetscInt cdim;
865: PetscFunctionBegin;
866: PetscCall(DMGetCoordinateDim(dm, &cdim));
867: PetscCall(PetscMalloc1(1, &rc));
868: dm->transformCtx = rc;
869: dm->transformSetUp = DMPlexBasisTransformSetUp_Rotation_Internal;
870: dm->transformDestroy = DMPlexBasisTransformDestroy_Rotation_Internal;
871: dm->transformGetMatrix = DMPlexBasisTransformGetMatrix_Rotation_Internal;
872: rc->dim = cdim;
873: rc->alpha = alpha;
874: rc->beta = beta;
875: rc->gamma = gamma;
876: PetscCall((*dm->transformSetUp)(dm, dm->transformCtx));
877: PetscCall(DMConstructBasisTransform_Internal(dm));
878: PetscFunctionReturn(PETSC_SUCCESS);
879: }
881: /*@C
882: DMPlexInsertBoundaryValuesEssential - Insert boundary values into a local vector using a function of the coordinates
884: Input Parameters:
885: + dm - The `DM`, with a `PetscDS` that matches the problem being constrained
886: . time - The time
887: . field - The field to constrain
888: . Nc - The number of constrained field components, or 0 for all components
889: . comps - An array of constrained component numbers, or `NULL` for all components
890: . label - The `DMLabel` defining constrained points
891: . numids - The number of `DMLabel` ids for constrained points
892: . ids - An array of ids for constrained points
893: . func - A pointwise function giving boundary values
894: - ctx - An optional user context for bcFunc
896: Output Parameter:
897: . locX - A local vector to receives the boundary values
899: Level: developer
901: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMLabel`, `DMPlexInsertBoundaryValuesEssentialField()`, `DMPlexInsertBoundaryValuesEssentialBdField()`, `DMAddBoundary()`
902: @*/
903: 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)
904: {
905: PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
906: void **ctxs;
907: PetscInt numFields;
909: PetscFunctionBegin;
910: PetscCall(DMGetNumFields(dm, &numFields));
911: PetscCall(PetscCalloc2(numFields, &funcs, numFields, &ctxs));
912: funcs[field] = func;
913: ctxs[field] = ctx;
914: PetscCall(DMProjectFunctionLabelLocal(dm, time, label, numids, ids, Nc, comps, funcs, ctxs, INSERT_BC_VALUES, locX));
915: PetscCall(PetscFree2(funcs, ctxs));
916: PetscFunctionReturn(PETSC_SUCCESS);
917: }
919: /*@C
920: DMPlexInsertBoundaryValuesEssentialField - Insert boundary values into a local vector using a function of the coordinates and field data
922: Input Parameters:
923: + dm - The `DM`, with a `PetscDS` that matches the problem being constrained
924: . time - The time
925: . locU - A local vector with the input solution values
926: . field - The field to constrain
927: . Nc - The number of constrained field components, or 0 for all components
928: . comps - An array of constrained component numbers, or `NULL` for all components
929: . label - The `DMLabel` defining constrained points
930: . numids - The number of `DMLabel` ids for constrained points
931: . ids - An array of ids for constrained points
932: . func - A pointwise function giving boundary values
933: - ctx - An optional user context for bcFunc
935: Output Parameter:
936: . locX - A local vector to receives the boundary values
938: Level: developer
940: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexInsertBoundaryValuesEssential()`, `DMPlexInsertBoundaryValuesEssentialBdField()`, `DMAddBoundary()`
941: @*/
942: PetscErrorCode DMPlexInsertBoundaryValuesEssentialField(DM dm, PetscReal time, Vec locU, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[], void (*func)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void *ctx, Vec locX)
943: {
944: void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
945: void **ctxs;
946: PetscInt numFields;
948: PetscFunctionBegin;
949: PetscCall(DMGetNumFields(dm, &numFields));
950: PetscCall(PetscCalloc2(numFields, &funcs, numFields, &ctxs));
951: funcs[field] = func;
952: ctxs[field] = ctx;
953: PetscCall(DMProjectFieldLabelLocal(dm, time, label, numids, ids, Nc, comps, locU, funcs, INSERT_BC_VALUES, locX));
954: PetscCall(PetscFree2(funcs, ctxs));
955: PetscFunctionReturn(PETSC_SUCCESS);
956: }
958: /*@C
959: DMPlexInsertBoundaryValuesEssentialBdField - Insert boundary values into a local vector using a function of the coordinates and boundary field data
961: Collective
963: Input Parameters:
964: + dm - The `DM`, with a `PetscDS` that matches the problem being constrained
965: . time - The time
966: . locU - A local vector with the input solution values
967: . field - The field to constrain
968: . Nc - The number of constrained field components, or 0 for all components
969: . comps - An array of constrained component numbers, or `NULL` for all components
970: . label - The `DMLabel` defining constrained points
971: . numids - The number of `DMLabel` ids for constrained points
972: . ids - An array of ids for constrained points
973: . func - A pointwise function giving boundary values, the calling sequence is given in `DMProjectBdFieldLabelLocal()`
974: - ctx - An optional user context for `func`
976: Output Parameter:
977: . locX - A local vector to receive the boundary values
979: Level: developer
981: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectBdFieldLabelLocal()`, `DMPlexInsertBoundaryValuesEssential()`, `DMPlexInsertBoundaryValuesEssentialField()`, `DMAddBoundary()`
982: @*/
983: PetscErrorCode DMPlexInsertBoundaryValuesEssentialBdField(DM dm, PetscReal time, Vec locU, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[], void (*func)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void *ctx, Vec locX)
984: {
985: void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
986: void **ctxs;
987: PetscInt numFields;
989: PetscFunctionBegin;
990: PetscCall(DMGetNumFields(dm, &numFields));
991: PetscCall(PetscCalloc2(numFields, &funcs, numFields, &ctxs));
992: funcs[field] = func;
993: ctxs[field] = ctx;
994: PetscCall(DMProjectBdFieldLabelLocal(dm, time, label, numids, ids, Nc, comps, locU, funcs, INSERT_BC_VALUES, locX));
995: PetscCall(PetscFree2(funcs, ctxs));
996: PetscFunctionReturn(PETSC_SUCCESS);
997: }
999: /*@C
1000: DMPlexInsertBoundaryValuesRiemann - Insert boundary values into a local vector
1002: Input Parameters:
1003: + dm - The `DM`, with a `PetscDS` that matches the problem being constrained
1004: . time - The time
1005: . faceGeometry - A vector with the FVM face geometry information
1006: . cellGeometry - A vector with the FVM cell geometry information
1007: . Grad - A vector with the FVM cell gradient information
1008: . field - The field to constrain
1009: . Nc - The number of constrained field components, or 0 for all components
1010: . comps - An array of constrained component numbers, or `NULL` for all components
1011: . label - The `DMLabel` defining constrained points
1012: . numids - The number of `DMLabel` ids for constrained points
1013: . ids - An array of ids for constrained points
1014: . func - A pointwise function giving boundary values
1015: - ctx - An optional user context for bcFunc
1017: Output Parameter:
1018: . locX - A local vector to receives the boundary values
1020: Level: developer
1022: Note:
1023: This implementation currently ignores the numcomps/comps argument from `DMAddBoundary()`
1025: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexInsertBoundaryValuesEssential()`, `DMPlexInsertBoundaryValuesEssentialField()`, `DMAddBoundary()`
1026: @*/
1027: 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[], PetscErrorCode (*func)(PetscReal, const PetscReal *, const PetscReal *, const PetscScalar *, PetscScalar *, void *), void *ctx, Vec locX)
1028: {
1029: PetscDS prob;
1030: PetscSF sf;
1031: DM dmFace, dmCell, dmGrad;
1032: const PetscScalar *facegeom, *cellgeom = NULL, *grad;
1033: const PetscInt *leaves;
1034: PetscScalar *x, *fx;
1035: PetscInt dim, nleaves, loc, fStart, fEnd, pdim, i;
1036: PetscErrorCode ierru = PETSC_SUCCESS;
1038: PetscFunctionBegin;
1039: PetscCall(DMGetPointSF(dm, &sf));
1040: PetscCall(PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL));
1041: nleaves = PetscMax(0, nleaves);
1042: PetscCall(DMGetDimension(dm, &dim));
1043: PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
1044: PetscCall(DMGetDS(dm, &prob));
1045: PetscCall(VecGetDM(faceGeometry, &dmFace));
1046: PetscCall(VecGetArrayRead(faceGeometry, &facegeom));
1047: if (cellGeometry) {
1048: PetscCall(VecGetDM(cellGeometry, &dmCell));
1049: PetscCall(VecGetArrayRead(cellGeometry, &cellgeom));
1050: }
1051: if (Grad) {
1052: PetscFV fv;
1054: PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fv));
1055: PetscCall(VecGetDM(Grad, &dmGrad));
1056: PetscCall(VecGetArrayRead(Grad, &grad));
1057: PetscCall(PetscFVGetNumComponents(fv, &pdim));
1058: PetscCall(DMGetWorkArray(dm, pdim, MPIU_SCALAR, &fx));
1059: }
1060: PetscCall(VecGetArray(locX, &x));
1061: for (i = 0; i < numids; ++i) {
1062: IS faceIS;
1063: const PetscInt *faces;
1064: PetscInt numFaces, f;
1066: PetscCall(DMLabelGetStratumIS(label, ids[i], &faceIS));
1067: if (!faceIS) continue; /* No points with that id on this process */
1068: PetscCall(ISGetLocalSize(faceIS, &numFaces));
1069: PetscCall(ISGetIndices(faceIS, &faces));
1070: for (f = 0; f < numFaces; ++f) {
1071: const PetscInt face = faces[f], *cells;
1072: PetscFVFaceGeom *fg;
1074: if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
1075: PetscCall(PetscFindInt(face, nleaves, (PetscInt *)leaves, &loc));
1076: if (loc >= 0) continue;
1077: PetscCall(DMPlexPointLocalRead(dmFace, face, facegeom, &fg));
1078: PetscCall(DMPlexGetSupport(dm, face, &cells));
1079: if (Grad) {
1080: PetscFVCellGeom *cg;
1081: PetscScalar *cx, *cgrad;
1082: PetscScalar *xG;
1083: PetscReal dx[3];
1084: PetscInt d;
1086: PetscCall(DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg));
1087: PetscCall(DMPlexPointLocalRead(dm, cells[0], x, &cx));
1088: PetscCall(DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad));
1089: PetscCall(DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG));
1090: DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
1091: for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d * dim], dx);
1092: PetscCall((*func)(time, fg->centroid, fg->normal, fx, xG, ctx));
1093: } else {
1094: PetscScalar *xI;
1095: PetscScalar *xG;
1097: PetscCall(DMPlexPointLocalRead(dm, cells[0], x, &xI));
1098: PetscCall(DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG));
1099: ierru = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);
1100: if (ierru) {
1101: PetscCall(ISRestoreIndices(faceIS, &faces));
1102: PetscCall(ISDestroy(&faceIS));
1103: goto cleanup;
1104: }
1105: }
1106: }
1107: PetscCall(ISRestoreIndices(faceIS, &faces));
1108: PetscCall(ISDestroy(&faceIS));
1109: }
1110: cleanup:
1111: PetscCall(VecRestoreArray(locX, &x));
1112: if (Grad) {
1113: PetscCall(DMRestoreWorkArray(dm, pdim, MPIU_SCALAR, &fx));
1114: PetscCall(VecRestoreArrayRead(Grad, &grad));
1115: }
1116: if (cellGeometry) PetscCall(VecRestoreArrayRead(cellGeometry, &cellgeom));
1117: PetscCall(VecRestoreArrayRead(faceGeometry, &facegeom));
1118: PetscCall(ierru);
1119: PetscFunctionReturn(PETSC_SUCCESS);
1120: }
1122: static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
1123: {
1124: PetscInt c;
1125: for (c = 0; c < Nc; ++c) u[c] = 0.0;
1126: return PETSC_SUCCESS;
1127: }
1129: PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
1130: {
1131: PetscObject isZero;
1132: PetscDS prob;
1133: PetscInt numBd, b;
1135: PetscFunctionBegin;
1136: PetscCall(DMGetDS(dm, &prob));
1137: PetscCall(PetscDSGetNumBoundary(prob, &numBd));
1138: PetscCall(PetscObjectQuery((PetscObject)locX, "__Vec_bc_zero__", &isZero));
1139: PetscCall(PetscDSUpdateBoundaryLabels(prob, dm));
1140: for (b = 0; b < numBd; ++b) {
1141: PetscWeakForm wf;
1142: DMBoundaryConditionType type;
1143: const char *name;
1144: DMLabel label;
1145: PetscInt field, Nc;
1146: const PetscInt *comps;
1147: PetscObject obj;
1148: PetscClassId id;
1149: void (*bvfunc)(void);
1150: PetscInt numids;
1151: const PetscInt *ids;
1152: void *ctx;
1154: PetscCall(PetscDSGetBoundary(prob, b, &wf, &type, &name, &label, &numids, &ids, &field, &Nc, &comps, &bvfunc, NULL, &ctx));
1155: if (insertEssential != (type & DM_BC_ESSENTIAL)) continue;
1156: PetscCall(DMGetField(dm, field, NULL, &obj));
1157: PetscCall(PetscObjectGetClassId(obj, &id));
1158: if (id == PETSCFE_CLASSID) {
1159: switch (type) {
1160: /* for FEM, there is no insertion to be done for non-essential boundary conditions */
1161: case DM_BC_ESSENTIAL: {
1162: PetscSimplePointFn *func = (PetscSimplePointFn *)bvfunc;
1164: if (isZero) func = zero;
1165: PetscCall(DMPlexLabelAddCells(dm, label));
1166: PetscCall(DMPlexInsertBoundaryValuesEssential(dm, time, field, Nc, comps, label, numids, ids, func, ctx, locX));
1167: PetscCall(DMPlexLabelClearCells(dm, label));
1168: } break;
1169: case DM_BC_ESSENTIAL_FIELD: {
1170: PetscPointFn *func = (PetscPointFn *)bvfunc;
1172: PetscCall(DMPlexLabelAddCells(dm, label));
1173: PetscCall(DMPlexInsertBoundaryValuesEssentialField(dm, time, locX, field, Nc, comps, label, numids, ids, func, ctx, locX));
1174: PetscCall(DMPlexLabelClearCells(dm, label));
1175: } break;
1176: default:
1177: break;
1178: }
1179: } else if (id == PETSCFV_CLASSID) {
1180: {
1181: PetscErrorCode (*func)(PetscReal, const PetscReal *, const PetscReal *, const PetscScalar *, PetscScalar *, void *) = (PetscErrorCode (*)(PetscReal, const PetscReal *, const PetscReal *, const PetscScalar *, PetscScalar *, void *))bvfunc;
1183: if (!faceGeomFVM) continue;
1184: PetscCall(DMPlexInsertBoundaryValuesRiemann(dm, time, faceGeomFVM, cellGeomFVM, gradFVM, field, Nc, comps, label, numids, ids, func, ctx, locX));
1185: }
1186: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1187: }
1188: PetscFunctionReturn(PETSC_SUCCESS);
1189: }
1191: PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues_Plex(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
1192: {
1193: PetscObject isZero;
1194: PetscDS prob;
1195: PetscInt numBd, b;
1197: PetscFunctionBegin;
1198: if (!locX) PetscFunctionReturn(PETSC_SUCCESS);
1199: PetscCall(DMGetDS(dm, &prob));
1200: PetscCall(PetscDSGetNumBoundary(prob, &numBd));
1201: PetscCall(PetscObjectQuery((PetscObject)locX, "__Vec_bc_zero__", &isZero));
1202: for (b = 0; b < numBd; ++b) {
1203: PetscWeakForm wf;
1204: DMBoundaryConditionType type;
1205: const char *name;
1206: DMLabel label;
1207: PetscInt field, Nc;
1208: const PetscInt *comps;
1209: PetscObject obj;
1210: PetscClassId id;
1211: PetscInt numids;
1212: const PetscInt *ids;
1213: void (*bvfunc)(void);
1214: void *ctx;
1216: PetscCall(PetscDSGetBoundary(prob, b, &wf, &type, &name, &label, &numids, &ids, &field, &Nc, &comps, NULL, &bvfunc, &ctx));
1217: if (insertEssential != (type & DM_BC_ESSENTIAL)) continue;
1218: PetscCall(DMGetField(dm, field, NULL, &obj));
1219: PetscCall(PetscObjectGetClassId(obj, &id));
1220: if (id == PETSCFE_CLASSID) {
1221: switch (type) {
1222: /* for FEM, there is no insertion to be done for non-essential boundary conditions */
1223: case DM_BC_ESSENTIAL: {
1224: PetscSimplePointFn *func_t = (PetscSimplePointFn *)bvfunc;
1226: if (isZero) func_t = zero;
1227: PetscCall(DMPlexLabelAddCells(dm, label));
1228: PetscCall(DMPlexInsertBoundaryValuesEssential(dm, time, field, Nc, comps, label, numids, ids, func_t, ctx, locX));
1229: PetscCall(DMPlexLabelClearCells(dm, label));
1230: } break;
1231: case DM_BC_ESSENTIAL_FIELD: {
1232: PetscPointFn *func_t = (PetscPointFn *)bvfunc;
1234: PetscCall(DMPlexLabelAddCells(dm, label));
1235: PetscCall(DMPlexInsertBoundaryValuesEssentialField(dm, time, locX, field, Nc, comps, label, numids, ids, func_t, ctx, locX));
1236: PetscCall(DMPlexLabelClearCells(dm, label));
1237: } break;
1238: default:
1239: break;
1240: }
1241: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1242: }
1243: PetscFunctionReturn(PETSC_SUCCESS);
1244: }
1246: PetscErrorCode DMPlexInsertBounds_Plex(DM dm, PetscBool lower, PetscReal time, Vec locB)
1247: {
1248: PetscDS ds;
1249: PetscInt numBd;
1251: PetscFunctionBegin;
1252: PetscCall(DMGetDS(dm, &ds));
1253: PetscCall(PetscDSGetNumBoundary(ds, &numBd));
1254: PetscCall(PetscDSUpdateBoundaryLabels(ds, dm));
1255: for (PetscInt b = 0; b < numBd; ++b) {
1256: PetscWeakForm wf;
1257: DMBoundaryConditionType type;
1258: const char *name;
1259: DMLabel label;
1260: PetscInt numids;
1261: const PetscInt *ids;
1262: PetscInt field, Nc;
1263: const PetscInt *comps;
1264: void (*bvfunc)(void);
1265: void *ctx;
1267: PetscCall(PetscDSGetBoundary(ds, b, &wf, &type, &name, &label, &numids, &ids, &field, &Nc, &comps, &bvfunc, NULL, &ctx));
1268: if (lower && type != DM_BC_LOWER_BOUND) continue;
1269: if (!lower && type != DM_BC_UPPER_BOUND) continue;
1270: PetscCall(DMPlexLabelAddCells(dm, label));
1271: {
1272: PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
1273: void **ctxs;
1274: PetscInt Nf;
1276: PetscCall(DMGetNumFields(dm, &Nf));
1277: PetscCall(PetscCalloc2(Nf, &funcs, Nf, &ctxs));
1278: funcs[field] = (PetscSimplePointFn *)bvfunc;
1279: ctxs[field] = ctx;
1280: PetscCall(DMProjectFunctionLabelLocal(dm, time, label, numids, ids, Nc, comps, funcs, ctxs, INSERT_ALL_VALUES, locB));
1281: PetscCall(PetscFree2(funcs, ctxs));
1282: }
1283: PetscCall(DMPlexLabelClearCells(dm, label));
1284: }
1285: PetscFunctionReturn(PETSC_SUCCESS);
1286: }
1288: /*@
1289: DMPlexInsertBoundaryValues - Puts coefficients which represent boundary values into the local solution vector
1291: Not Collective
1293: Input Parameters:
1294: + dm - The `DM`
1295: . insertEssential - Should I insert essential (e.g. Dirichlet) or inessential (e.g. Neumann) boundary conditions
1296: . time - The time
1297: . faceGeomFVM - Face geometry data for FV discretizations
1298: . cellGeomFVM - Cell geometry data for FV discretizations
1299: - gradFVM - Gradient reconstruction data for FV discretizations
1301: Output Parameter:
1302: . locX - Solution updated with boundary values
1304: Level: intermediate
1306: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectFunctionLabelLocal()`, `DMAddBoundary()`
1307: @*/
1308: PetscErrorCode DMPlexInsertBoundaryValues(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
1309: {
1310: PetscFunctionBegin;
1316: PetscTryMethod(dm, "DMPlexInsertBoundaryValues_C", (DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec), (dm, insertEssential, locX, time, faceGeomFVM, cellGeomFVM, gradFVM));
1317: PetscFunctionReturn(PETSC_SUCCESS);
1318: }
1320: /*@
1321: DMPlexInsertTimeDerivativeBoundaryValues - Puts coefficients which represent boundary values of the time derivative into the local solution vector
1323: Input Parameters:
1324: + dm - The `DM`
1325: . insertEssential - Should I insert essential (e.g. Dirichlet) or inessential (e.g. Neumann) boundary conditions
1326: . time - The time
1327: . faceGeomFVM - Face geometry data for FV discretizations
1328: . cellGeomFVM - Cell geometry data for FV discretizations
1329: - gradFVM - Gradient reconstruction data for FV discretizations
1331: Output Parameter:
1332: . locX_t - Solution updated with boundary values
1334: Level: developer
1336: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectFunctionLabelLocal()`
1337: @*/
1338: PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues(DM dm, PetscBool insertEssential, Vec locX_t, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
1339: {
1340: PetscFunctionBegin;
1346: PetscTryMethod(dm, "DMPlexInsertTimeDerivativeBoundaryValues_C", (DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec), (dm, insertEssential, locX_t, time, faceGeomFVM, cellGeomFVM, gradFVM));
1347: PetscFunctionReturn(PETSC_SUCCESS);
1348: }
1350: /*@
1351: DMPlexInsertBounds - Puts coefficients which represent solution bounds into the local bounds vector
1353: Not Collective
1355: Input Parameters:
1356: + dm - The `DM`
1357: . lower - If `PETSC_TRUE` use `DM_BC_LOWER_BOUND` conditions, otherwise use `DM_BC_UPPER_BOUND`
1358: - time - The time
1360: Output Parameter:
1361: . locB - Bounds vector updated with new bounds
1363: Level: intermediate
1365: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectFunctionLabelLocal()`, `PetscDSAddBoundary()`
1366: @*/
1367: PetscErrorCode DMPlexInsertBounds(DM dm, PetscBool lower, PetscReal time, Vec locB)
1368: {
1369: PetscFunctionBegin;
1372: PetscTryMethod(dm, "DMPlexInsertBounds_C", (DM, PetscBool, PetscReal, Vec), (dm, lower, time, locB));
1373: PetscFunctionReturn(PETSC_SUCCESS);
1374: }
1376: // Handle non-essential (e.g. outflow) boundary values
1377: PetscErrorCode DMPlexInsertBoundaryValuesFVM(DM dm, PetscFV fv, Vec locX, PetscReal time, Vec *locGradient)
1378: {
1379: DM dmGrad;
1380: Vec cellGeometryFVM, faceGeometryFVM, locGrad = NULL;
1382: PetscFunctionBegin;
1386: if (locGradient) {
1387: PetscAssertPointer(locGradient, 5);
1388: *locGradient = NULL;
1389: }
1390: PetscCall(DMPlexGetGeometryFVM(dm, &faceGeometryFVM, &cellGeometryFVM, NULL));
1391: /* Reconstruct and limit cell gradients */
1392: PetscCall(DMPlexGetGradientDM(dm, fv, &dmGrad));
1393: if (dmGrad) {
1394: Vec grad;
1395: PetscInt fStart, fEnd;
1397: PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
1398: PetscCall(DMGetGlobalVector(dmGrad, &grad));
1399: PetscCall(DMPlexReconstructGradients_Internal(dm, fv, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad));
1400: /* Communicate gradient values */
1401: PetscCall(DMGetLocalVector(dmGrad, &locGrad));
1402: PetscCall(DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad));
1403: PetscCall(DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad));
1404: PetscCall(DMRestoreGlobalVector(dmGrad, &grad));
1405: }
1406: PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_FALSE, locX, time, faceGeometryFVM, cellGeometryFVM, locGrad));
1407: if (locGradient) *locGradient = locGrad;
1408: else if (locGrad) PetscCall(DMRestoreLocalVector(dmGrad, &locGrad));
1409: PetscFunctionReturn(PETSC_SUCCESS);
1410: }
1412: PetscErrorCode DMComputeL2Diff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
1413: {
1414: Vec localX;
1416: PetscFunctionBegin;
1417: PetscCall(DMGetLocalVector(dm, &localX));
1418: PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localX, time, NULL, NULL, NULL));
1419: PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX));
1420: PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX));
1421: PetscCall(DMPlexComputeL2DiffLocal(dm, time, funcs, ctxs, localX, diff));
1422: PetscCall(DMRestoreLocalVector(dm, &localX));
1423: PetscFunctionReturn(PETSC_SUCCESS);
1424: }
1426: /*@C
1427: DMPlexComputeL2DiffLocal - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
1429: Collective
1431: Input Parameters:
1432: + dm - The `DM`
1433: . time - The time
1434: . funcs - The functions to evaluate for each field component
1435: . ctxs - Optional array of contexts to pass to each function, or `NULL`.
1436: - localX - The coefficient vector u_h, a local vector
1438: Output Parameter:
1439: . diff - The diff ||u - u_h||_2
1441: Level: developer
1443: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectFunction()`, `DMComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
1444: @*/
1445: PetscErrorCode DMPlexComputeL2DiffLocal(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec localX, PetscReal *diff)
1446: {
1447: const PetscInt debug = ((DM_Plex *)dm->data)->printL2;
1448: DM tdm;
1449: Vec tv;
1450: PetscSection section;
1451: PetscQuadrature quad;
1452: PetscFEGeom fegeom;
1453: PetscScalar *funcVal, *interpolant;
1454: PetscReal *coords, *gcoords;
1455: PetscReal localDiff = 0.0;
1456: const PetscReal *quadWeights;
1457: PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cellHeight, cStart, cEnd, c, field, fieldOffset;
1458: PetscBool transform;
1460: PetscFunctionBegin;
1461: PetscCall(DMGetDimension(dm, &dim));
1462: PetscCall(DMGetCoordinateDim(dm, &coordDim));
1463: fegeom.dimEmbed = coordDim;
1464: PetscCall(DMGetLocalSection(dm, §ion));
1465: PetscCall(PetscSectionGetNumFields(section, &numFields));
1466: PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
1467: PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
1468: PetscCall(DMHasBasisTransform(dm, &transform));
1469: PetscCheck(numFields, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fields is zero!");
1470: for (field = 0; field < numFields; ++field) {
1471: PetscObject obj;
1472: PetscClassId id;
1473: PetscInt Nc;
1475: PetscCall(DMGetField(dm, field, NULL, &obj));
1476: PetscCall(PetscObjectGetClassId(obj, &id));
1477: if (id == PETSCFE_CLASSID) {
1478: PetscFE fe = (PetscFE)obj;
1480: PetscCall(PetscFEGetQuadrature(fe, &quad));
1481: PetscCall(PetscFEGetNumComponents(fe, &Nc));
1482: } else if (id == PETSCFV_CLASSID) {
1483: PetscFV fv = (PetscFV)obj;
1485: PetscCall(PetscFVGetQuadrature(fv, &quad));
1486: PetscCall(PetscFVGetNumComponents(fv, &Nc));
1487: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1488: numComponents += Nc;
1489: }
1490: PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights));
1491: PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents);
1492: PetscCall(PetscMalloc6(numComponents, &funcVal, numComponents, &interpolant, coordDim * (Nq + 1), &coords, Nq, &fegeom.detJ, coordDim * coordDim * Nq, &fegeom.J, coordDim * coordDim * Nq, &fegeom.invJ));
1493: PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
1494: PetscCall(DMPlexGetSimplexOrBoxCells(dm, cellHeight, &cStart, &cEnd));
1495: for (c = cStart; c < cEnd; ++c) {
1496: PetscScalar *x = NULL;
1497: PetscReal elemDiff = 0.0;
1498: PetscInt qc = 0;
1500: PetscCall(DMPlexComputeCellGeometryFEM(dm, c, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
1501: PetscCall(DMPlexVecGetOrientedClosure(dm, NULL, PETSC_FALSE, localX, c, 0, NULL, &x));
1503: for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1504: PetscObject obj;
1505: PetscClassId id;
1506: void *const ctx = ctxs ? ctxs[field] : NULL;
1507: PetscInt Nb, Nc, q, fc;
1509: PetscCall(DMGetField(dm, field, NULL, &obj));
1510: PetscCall(PetscObjectGetClassId(obj, &id));
1511: if (id == PETSCFE_CLASSID) {
1512: PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
1513: PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb));
1514: } else if (id == PETSCFV_CLASSID) {
1515: PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
1516: Nb = 1;
1517: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1518: if (debug) {
1519: char title[1024];
1520: PetscCall(PetscSNPrintf(title, 1023, "Solution for Field %" PetscInt_FMT, field));
1521: PetscCall(DMPrintCellVector(c, title, Nb, &x[fieldOffset]));
1522: }
1523: for (q = 0; q < Nq; ++q) {
1524: PetscFEGeom qgeom;
1525: PetscErrorCode ierr;
1527: qgeom.dimEmbed = fegeom.dimEmbed;
1528: qgeom.J = &fegeom.J[q * coordDim * coordDim];
1529: qgeom.invJ = &fegeom.invJ[q * coordDim * coordDim];
1530: qgeom.detJ = &fegeom.detJ[q];
1531: PetscCheck(fegeom.detJ[q] > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT ", point %" PetscInt_FMT, (double)fegeom.detJ[q], c, q);
1532: if (transform) {
1533: gcoords = &coords[coordDim * Nq];
1534: PetscCall(DMPlexBasisTransformApplyReal_Internal(dm, &coords[coordDim * q], PETSC_TRUE, coordDim, &coords[coordDim * q], gcoords, dm->transformCtx));
1535: } else {
1536: gcoords = &coords[coordDim * q];
1537: }
1538: PetscCall(PetscArrayzero(funcVal, Nc));
1539: ierr = (*funcs[field])(coordDim, time, gcoords, Nc, funcVal, ctx);
1540: if (ierr) {
1541: PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x));
1542: PetscCall(DMRestoreLocalVector(dm, &localX));
1543: PetscCall(PetscFree6(funcVal, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
1544: }
1545: if (transform) PetscCall(DMPlexBasisTransformApply_Internal(dm, &coords[coordDim * q], PETSC_FALSE, Nc, funcVal, funcVal, dm->transformCtx));
1546: if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolate_Static((PetscFE)obj, &x[fieldOffset], &qgeom, q, interpolant));
1547: else if (id == PETSCFV_CLASSID) PetscCall(PetscFVInterpolate_Static((PetscFV)obj, &x[fieldOffset], q, interpolant));
1548: else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1549: for (fc = 0; fc < Nc; ++fc) {
1550: const PetscReal wt = quadWeights[q * qNc + (qNc == 1 ? 0 : qc + fc)];
1551: if (debug)
1552: PetscCall(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),
1553: (double)(PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc])) * wt * fegeom.detJ[q]), (double)PetscRealPart(interpolant[fc]), (double)PetscRealPart(funcVal[fc])));
1554: elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc])) * wt * fegeom.detJ[q];
1555: }
1556: }
1557: fieldOffset += Nb;
1558: qc += Nc;
1559: }
1560: PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x));
1561: if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " elem %" PetscInt_FMT " diff %g\n", c, (double)elemDiff));
1562: localDiff += elemDiff;
1563: }
1564: PetscCall(PetscFree6(funcVal, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
1565: PetscCallMPI(MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm)));
1566: *diff = PetscSqrtReal(*diff);
1567: PetscFunctionReturn(PETSC_SUCCESS);
1568: }
1570: 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)
1571: {
1572: const PetscInt debug = ((DM_Plex *)dm->data)->printL2;
1573: DM tdm;
1574: PetscSection section;
1575: PetscQuadrature quad;
1576: Vec localX, tv;
1577: PetscScalar *funcVal, *interpolant;
1578: const PetscReal *quadWeights;
1579: PetscFEGeom fegeom;
1580: PetscReal *coords, *gcoords;
1581: PetscReal localDiff = 0.0;
1582: PetscInt dim, coordDim, qNc = 0, Nq = 0, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset;
1583: PetscBool transform;
1585: PetscFunctionBegin;
1586: PetscCall(DMGetDimension(dm, &dim));
1587: PetscCall(DMGetCoordinateDim(dm, &coordDim));
1588: fegeom.dimEmbed = coordDim;
1589: PetscCall(DMGetLocalSection(dm, §ion));
1590: PetscCall(PetscSectionGetNumFields(section, &numFields));
1591: PetscCall(DMGetLocalVector(dm, &localX));
1592: PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX));
1593: PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX));
1594: PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
1595: PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
1596: PetscCall(DMHasBasisTransform(dm, &transform));
1597: for (field = 0; field < numFields; ++field) {
1598: PetscFE fe;
1599: PetscInt Nc;
1601: PetscCall(DMGetField(dm, field, NULL, (PetscObject *)&fe));
1602: PetscCall(PetscFEGetQuadrature(fe, &quad));
1603: PetscCall(PetscFEGetNumComponents(fe, &Nc));
1604: numComponents += Nc;
1605: }
1606: PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights));
1607: PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents);
1608: /* PetscCall(DMProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX)); */
1609: PetscCall(PetscMalloc6(numComponents, &funcVal, coordDim * (Nq + 1), &coords, coordDim * coordDim * Nq, &fegeom.J, coordDim * coordDim * Nq, &fegeom.invJ, numComponents * coordDim, &interpolant, Nq, &fegeom.detJ));
1610: PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
1611: for (c = cStart; c < cEnd; ++c) {
1612: PetscScalar *x = NULL;
1613: PetscReal elemDiff = 0.0;
1614: PetscInt qc = 0;
1616: PetscCall(DMPlexComputeCellGeometryFEM(dm, c, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
1617: PetscCall(DMPlexVecGetOrientedClosure(dm, NULL, PETSC_FALSE, localX, c, 0, NULL, &x));
1619: for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1620: PetscFE fe;
1621: void *const ctx = ctxs ? ctxs[field] : NULL;
1622: PetscInt Nb, Nc, q, fc;
1624: PetscCall(DMGetField(dm, field, NULL, (PetscObject *)&fe));
1625: PetscCall(PetscFEGetDimension(fe, &Nb));
1626: PetscCall(PetscFEGetNumComponents(fe, &Nc));
1627: if (debug) {
1628: char title[1024];
1629: PetscCall(PetscSNPrintf(title, 1023, "Solution for Field %" PetscInt_FMT, field));
1630: PetscCall(DMPrintCellVector(c, title, Nb, &x[fieldOffset]));
1631: }
1632: for (q = 0; q < Nq; ++q) {
1633: PetscFEGeom qgeom;
1634: PetscErrorCode ierr;
1636: qgeom.dimEmbed = fegeom.dimEmbed;
1637: qgeom.J = &fegeom.J[q * coordDim * coordDim];
1638: qgeom.invJ = &fegeom.invJ[q * coordDim * coordDim];
1639: qgeom.detJ = &fegeom.detJ[q];
1640: PetscCheck(fegeom.detJ[q] > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT ", quadrature points %" PetscInt_FMT, (double)fegeom.detJ[q], c, q);
1641: if (transform) {
1642: gcoords = &coords[coordDim * Nq];
1643: PetscCall(DMPlexBasisTransformApplyReal_Internal(dm, &coords[coordDim * q], PETSC_TRUE, coordDim, &coords[coordDim * q], gcoords, dm->transformCtx));
1644: } else {
1645: gcoords = &coords[coordDim * q];
1646: }
1647: PetscCall(PetscArrayzero(funcVal, Nc));
1648: ierr = (*funcs[field])(coordDim, time, gcoords, n, Nc, funcVal, ctx);
1649: if (ierr) {
1650: PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x));
1651: PetscCall(DMRestoreLocalVector(dm, &localX));
1652: PetscCall(PetscFree6(funcVal, coords, fegeom.J, fegeom.invJ, interpolant, fegeom.detJ));
1653: }
1654: if (transform) PetscCall(DMPlexBasisTransformApply_Internal(dm, &coords[coordDim * q], PETSC_FALSE, Nc, funcVal, funcVal, dm->transformCtx));
1655: PetscCall(PetscFEInterpolateGradient_Static(fe, 1, &x[fieldOffset], &qgeom, q, interpolant));
1656: /* Overwrite with the dot product if the normal is given */
1657: if (n) {
1658: for (fc = 0; fc < Nc; ++fc) {
1659: PetscScalar sum = 0.0;
1660: PetscInt d;
1661: for (d = 0; d < dim; ++d) sum += interpolant[fc * dim + d] * n[d];
1662: interpolant[fc] = sum;
1663: }
1664: }
1665: for (fc = 0; fc < Nc; ++fc) {
1666: const PetscReal wt = quadWeights[q * qNc + (qNc == 1 ? 0 : qc + fc)];
1667: if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " elem %" PetscInt_FMT " fieldDer %" PetscInt_FMT ",%" PetscInt_FMT " diff %g\n", c, field, fc, (double)(PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc])) * wt * fegeom.detJ[q])));
1668: elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc])) * wt * fegeom.detJ[q];
1669: }
1670: }
1671: fieldOffset += Nb;
1672: qc += Nc;
1673: }
1674: PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x));
1675: if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " elem %" PetscInt_FMT " diff %g\n", c, (double)elemDiff));
1676: localDiff += elemDiff;
1677: }
1678: PetscCall(PetscFree6(funcVal, coords, fegeom.J, fegeom.invJ, interpolant, fegeom.detJ));
1679: PetscCall(DMRestoreLocalVector(dm, &localX));
1680: PetscCallMPI(MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm)));
1681: *diff = PetscSqrtReal(*diff);
1682: PetscFunctionReturn(PETSC_SUCCESS);
1683: }
1685: PetscErrorCode DMComputeL2FieldDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
1686: {
1687: const PetscInt debug = ((DM_Plex *)dm->data)->printL2;
1688: DM tdm;
1689: DMLabel depthLabel;
1690: PetscSection section;
1691: Vec localX, tv;
1692: PetscReal *localDiff;
1693: PetscInt dim, depth, dE, Nf, f, Nds, s;
1694: PetscBool transform;
1696: PetscFunctionBegin;
1697: PetscCall(DMGetDimension(dm, &dim));
1698: PetscCall(DMGetCoordinateDim(dm, &dE));
1699: PetscCall(DMGetLocalSection(dm, §ion));
1700: PetscCall(DMGetLocalVector(dm, &localX));
1701: PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
1702: PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
1703: PetscCall(DMHasBasisTransform(dm, &transform));
1704: PetscCall(DMGetNumFields(dm, &Nf));
1705: PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
1706: PetscCall(DMLabelGetNumValues(depthLabel, &depth));
1708: PetscCall(VecSet(localX, 0.0));
1709: PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX));
1710: PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX));
1711: PetscCall(DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX));
1712: PetscCall(DMGetNumDS(dm, &Nds));
1713: PetscCall(PetscCalloc1(Nf, &localDiff));
1714: for (s = 0; s < Nds; ++s) {
1715: PetscDS ds;
1716: DMLabel label;
1717: IS fieldIS, pointIS;
1718: const PetscInt *fields, *points = NULL;
1719: PetscQuadrature quad;
1720: const PetscReal *quadPoints, *quadWeights;
1721: PetscFEGeom fegeom;
1722: PetscReal *coords, *gcoords;
1723: PetscScalar *funcVal, *interpolant;
1724: PetscBool isCohesive;
1725: PetscInt qNc, Nq, totNc, cStart = 0, cEnd, c, dsNf;
1727: PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL));
1728: PetscCall(ISGetIndices(fieldIS, &fields));
1729: PetscCall(PetscDSIsCohesive(ds, &isCohesive));
1730: PetscCall(PetscDSGetNumFields(ds, &dsNf));
1731: PetscCall(PetscDSGetTotalComponents(ds, &totNc));
1732: PetscCall(PetscDSGetQuadrature(ds, &quad));
1733: PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
1734: PetscCheck(!(qNc != 1) || !(qNc != totNc), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, totNc);
1735: PetscCall(PetscCalloc6(totNc, &funcVal, totNc, &interpolant, dE * (Nq + 1), &coords, Nq, &fegeom.detJ, dE * dE * Nq, &fegeom.J, dE * dE * Nq, &fegeom.invJ));
1736: if (!label) {
1737: PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
1738: } else {
1739: PetscCall(DMLabelGetStratumIS(label, 1, &pointIS));
1740: PetscCall(ISGetLocalSize(pointIS, &cEnd));
1741: PetscCall(ISGetIndices(pointIS, &points));
1742: }
1743: for (c = cStart; c < cEnd; ++c) {
1744: const PetscInt cell = points ? points[c] : c;
1745: PetscScalar *x = NULL;
1746: const PetscInt *cone;
1747: PetscInt qc = 0, fOff = 0, dep;
1749: PetscCall(DMLabelGetValue(depthLabel, cell, &dep));
1750: if (dep != depth - 1) continue;
1751: if (isCohesive) {
1752: PetscCall(DMPlexGetCone(dm, cell, &cone));
1753: PetscCall(DMPlexComputeCellGeometryFEM(dm, cone[0], quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
1754: } else {
1755: PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
1756: }
1757: PetscCall(DMPlexVecGetOrientedClosure(dm, NULL, PETSC_FALSE, localX, cell, 0, NULL, &x));
1758: for (f = 0; f < dsNf; ++f) {
1759: PetscObject obj;
1760: PetscClassId id;
1761: void *const ctx = ctxs ? ctxs[fields[f]] : NULL;
1762: PetscInt Nb, Nc, q, fc;
1763: PetscReal elemDiff = 0.0;
1764: PetscBool cohesive;
1766: PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
1767: PetscCall(PetscDSGetDiscretization(ds, f, &obj));
1768: PetscCall(PetscObjectGetClassId(obj, &id));
1769: if (id == PETSCFE_CLASSID) {
1770: PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
1771: PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb));
1772: } else if (id == PETSCFV_CLASSID) {
1773: PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
1774: Nb = 1;
1775: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, fields[f]);
1776: if (isCohesive && !cohesive) {
1777: fOff += Nb * 2;
1778: qc += Nc;
1779: continue;
1780: }
1781: if (debug) {
1782: char title[1024];
1783: PetscCall(PetscSNPrintf(title, 1023, "Solution for Field %" PetscInt_FMT, fields[f]));
1784: PetscCall(DMPrintCellVector(cell, title, Nb, &x[fOff]));
1785: }
1786: for (q = 0; q < Nq; ++q) {
1787: PetscFEGeom qgeom;
1788: PetscErrorCode ierr;
1790: qgeom.dimEmbed = fegeom.dimEmbed;
1791: qgeom.J = &fegeom.J[q * dE * dE];
1792: qgeom.invJ = &fegeom.invJ[q * dE * dE];
1793: qgeom.detJ = &fegeom.detJ[q];
1794: PetscCheck(fegeom.detJ[q] > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for cell %" PetscInt_FMT ", quadrature point %" PetscInt_FMT, (double)fegeom.detJ[q], cell, q);
1795: if (transform) {
1796: gcoords = &coords[dE * Nq];
1797: PetscCall(DMPlexBasisTransformApplyReal_Internal(dm, &coords[dE * q], PETSC_TRUE, dE, &coords[dE * q], gcoords, dm->transformCtx));
1798: } else {
1799: gcoords = &coords[dE * q];
1800: }
1801: for (fc = 0; fc < Nc; ++fc) funcVal[fc] = 0.;
1802: ierr = (*funcs[fields[f]])(dE, time, gcoords, Nc, funcVal, ctx);
1803: if (ierr) {
1804: PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, cell, NULL, &x));
1805: PetscCall(DMRestoreLocalVector(dm, &localX));
1806: PetscCall(PetscFree6(funcVal, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
1807: }
1808: if (transform) PetscCall(DMPlexBasisTransformApply_Internal(dm, &coords[dE * q], PETSC_FALSE, Nc, funcVal, funcVal, dm->transformCtx));
1809: /* Call once for each face, except for lagrange field */
1810: if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolate_Static((PetscFE)obj, &x[fOff], &qgeom, q, interpolant));
1811: else if (id == PETSCFV_CLASSID) PetscCall(PetscFVInterpolate_Static((PetscFV)obj, &x[fOff], q, interpolant));
1812: else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, fields[f]);
1813: for (fc = 0; fc < Nc; ++fc) {
1814: const PetscReal wt = quadWeights[q * qNc + (qNc == 1 ? 0 : qc + fc)];
1815: if (debug)
1816: PetscCall(PetscPrintf(PETSC_COMM_SELF, " cell %" PetscInt_FMT " field %" PetscInt_FMT ",%" PetscInt_FMT " 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),
1817: (double)(PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc])) * wt * fegeom.detJ[q])));
1818: elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc])) * wt * fegeom.detJ[q];
1819: }
1820: }
1821: fOff += Nb;
1822: qc += Nc;
1823: localDiff[fields[f]] += elemDiff;
1824: if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " cell %" PetscInt_FMT " field %" PetscInt_FMT " cum diff %g\n", cell, fields[f], (double)localDiff[fields[f]]));
1825: }
1826: PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, cell, NULL, &x));
1827: }
1828: if (label) {
1829: PetscCall(ISRestoreIndices(pointIS, &points));
1830: PetscCall(ISDestroy(&pointIS));
1831: }
1832: PetscCall(ISRestoreIndices(fieldIS, &fields));
1833: PetscCall(PetscFree6(funcVal, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
1834: }
1835: PetscCall(DMRestoreLocalVector(dm, &localX));
1836: PetscCallMPI(MPIU_Allreduce(localDiff, diff, Nf, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm)));
1837: PetscCall(PetscFree(localDiff));
1838: for (f = 0; f < Nf; ++f) diff[f] = PetscSqrtReal(diff[f]);
1839: PetscFunctionReturn(PETSC_SUCCESS);
1840: }
1842: /*@C
1843: 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.
1845: Collective
1847: Input Parameters:
1848: + dm - The `DM`
1849: . time - The time
1850: . funcs - The functions to evaluate for each field component: `NULL` means that component does not contribute to error calculation
1851: . ctxs - Optional array of contexts to pass to each function, or `NULL`.
1852: - X - The coefficient vector u_h
1854: Output Parameter:
1855: . D - A `Vec` which holds the difference ||u - u_h||_2 for each cell
1857: Level: developer
1859: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectFunction()`, `DMComputeL2Diff()`, `DMPlexComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
1860: @*/
1861: PetscErrorCode DMPlexComputeL2DiffVec(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, Vec D)
1862: {
1863: PetscSection section;
1864: PetscQuadrature quad;
1865: Vec localX;
1866: PetscFEGeom fegeom;
1867: PetscScalar *funcVal, *interpolant;
1868: PetscReal *coords;
1869: const PetscReal *quadPoints, *quadWeights;
1870: PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, c, field, fieldOffset;
1872: PetscFunctionBegin;
1873: PetscCall(VecSet(D, 0.0));
1874: PetscCall(DMGetDimension(dm, &dim));
1875: PetscCall(DMGetCoordinateDim(dm, &coordDim));
1876: PetscCall(DMGetLocalSection(dm, §ion));
1877: PetscCall(PetscSectionGetNumFields(section, &numFields));
1878: PetscCall(DMGetLocalVector(dm, &localX));
1879: PetscCall(DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX));
1880: PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX));
1881: PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX));
1882: for (field = 0; field < numFields; ++field) {
1883: PetscObject obj;
1884: PetscClassId id;
1885: PetscInt Nc;
1887: PetscCall(DMGetField(dm, field, NULL, &obj));
1888: PetscCall(PetscObjectGetClassId(obj, &id));
1889: if (id == PETSCFE_CLASSID) {
1890: PetscFE fe = (PetscFE)obj;
1892: PetscCall(PetscFEGetQuadrature(fe, &quad));
1893: PetscCall(PetscFEGetNumComponents(fe, &Nc));
1894: } else if (id == PETSCFV_CLASSID) {
1895: PetscFV fv = (PetscFV)obj;
1897: PetscCall(PetscFVGetQuadrature(fv, &quad));
1898: PetscCall(PetscFVGetNumComponents(fv, &Nc));
1899: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1900: numComponents += Nc;
1901: }
1902: PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
1903: PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents);
1904: PetscCall(PetscMalloc6(numComponents, &funcVal, numComponents, &interpolant, coordDim * Nq, &coords, Nq, &fegeom.detJ, coordDim * coordDim * Nq, &fegeom.J, coordDim * coordDim * Nq, &fegeom.invJ));
1905: PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
1906: for (c = cStart; c < cEnd; ++c) {
1907: PetscScalar *x = NULL;
1908: PetscScalar elemDiff = 0.0;
1909: PetscInt qc = 0;
1911: PetscCall(DMPlexComputeCellGeometryFEM(dm, c, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
1912: PetscCall(DMPlexVecGetOrientedClosure(dm, NULL, PETSC_FALSE, localX, c, 0, NULL, &x));
1914: for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1915: PetscObject obj;
1916: PetscClassId id;
1917: void *const ctx = ctxs ? ctxs[field] : NULL;
1918: PetscInt Nb, Nc, q, fc;
1920: PetscCall(DMGetField(dm, field, NULL, &obj));
1921: PetscCall(PetscObjectGetClassId(obj, &id));
1922: if (id == PETSCFE_CLASSID) {
1923: PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
1924: PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb));
1925: } else if (id == PETSCFV_CLASSID) {
1926: PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
1927: Nb = 1;
1928: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1929: if (funcs[field]) {
1930: for (q = 0; q < Nq; ++q) {
1931: PetscFEGeom qgeom;
1933: qgeom.dimEmbed = fegeom.dimEmbed;
1934: qgeom.J = &fegeom.J[q * coordDim * coordDim];
1935: qgeom.invJ = &fegeom.invJ[q * coordDim * coordDim];
1936: qgeom.detJ = &fegeom.detJ[q];
1937: PetscCheck(fegeom.detJ[q] > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT ", quadrature points %" PetscInt_FMT, (double)fegeom.detJ[q], c, q);
1938: PetscCall((*funcs[field])(coordDim, time, &coords[q * coordDim], Nc, funcVal, ctx));
1939: #if defined(needs_fix_with_return_code_argument)
1940: if (ierr) {
1941: PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x));
1942: PetscCall(PetscFree6(funcVal, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
1943: PetscCall(DMRestoreLocalVector(dm, &localX));
1944: }
1945: #endif
1946: if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolate_Static((PetscFE)obj, &x[fieldOffset], &qgeom, q, interpolant));
1947: else if (id == PETSCFV_CLASSID) PetscCall(PetscFVInterpolate_Static((PetscFV)obj, &x[fieldOffset], q, interpolant));
1948: else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
1949: for (fc = 0; fc < Nc; ++fc) {
1950: const PetscReal wt = quadWeights[q * qNc + (qNc == 1 ? 0 : qc + fc)];
1951: elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc])) * wt * fegeom.detJ[q];
1952: }
1953: }
1954: }
1955: fieldOffset += Nb;
1956: qc += Nc;
1957: }
1958: PetscCall(DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x));
1959: PetscCall(VecSetValue(D, c - cStart, elemDiff, INSERT_VALUES));
1960: }
1961: PetscCall(PetscFree6(funcVal, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
1962: PetscCall(DMRestoreLocalVector(dm, &localX));
1963: PetscCall(VecSqrtAbs(D));
1964: PetscFunctionReturn(PETSC_SUCCESS);
1965: }
1967: /*@
1968: DMPlexComputeL2FluxDiffVecLocal - This function computes the integral of the difference between the gradient of field `f`in `u` and field `mf` in `mu`
1970: Collective
1972: Input Parameters:
1973: + lu - The local `Vec` containing the primal solution
1974: . f - The field number for the potential
1975: . lmu - The local `Vec` containing the mixed solution
1976: - mf - The field number for the flux
1978: Output Parameter:
1979: . eFlux - A global `Vec` which holds $||\nabla u_f - \mu_{mf}||$
1981: Level: advanced
1983: Notes:
1984: We assume that the `DM` for each solution has the same topology, geometry, and quadrature.
1986: This is usually used to get an error estimate for the primal solution, using the flux from a mixed solution.
1988: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexComputeL2FluxDiffVec()`, `DMProjectFunction()`, `DMComputeL2Diff()`, `DMPlexComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
1989: @*/
1990: PetscErrorCode DMPlexComputeL2FluxDiffVecLocal(Vec lu, PetscInt f, Vec lmu, PetscInt mf, Vec eFlux)
1991: {
1992: DM dm, mdm, edm;
1993: PetscFE fe, mfe;
1994: PetscFEGeom fegeom;
1995: PetscQuadrature quad;
1996: const PetscReal *quadWeights;
1997: PetscReal *coords;
1998: PetscScalar *interpolant, *minterpolant, *earray;
1999: PetscInt cdim, mcdim, cStart, cEnd, Nc, mNc, qNc, Nq;
2000: MPI_Comm comm;
2002: PetscFunctionBegin;
2003: PetscCall(VecGetDM(lu, &dm));
2004: PetscCall(VecGetDM(lmu, &mdm));
2005: PetscCall(VecGetDM(eFlux, &edm));
2006: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2007: PetscCall(VecSet(eFlux, 0.0));
2009: // Check if the both problems are on the same mesh
2010: PetscCall(DMGetCoordinateDim(dm, &cdim));
2011: PetscCall(DMGetCoordinateDim(mdm, &mcdim));
2012: PetscCheck(cdim == mcdim, comm, PETSC_ERR_ARG_SIZ, "primal coordinate Dim %" PetscInt_FMT " != %" PetscInt_FMT " mixed coordinate Dim", cdim, mcdim);
2013: fegeom.dimEmbed = cdim;
2015: PetscCall(DMGetField(dm, f, NULL, (PetscObject *)&fe));
2016: PetscCall(DMGetField(mdm, mf, NULL, (PetscObject *)&mfe));
2017: PetscCall(PetscFEGetNumComponents(fe, &Nc));
2018: PetscCall(PetscFEGetNumComponents(mfe, &mNc));
2019: PetscCall(PetscFEGetQuadrature(fe, &quad));
2020: PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights));
2021: PetscCheck(qNc == 1 || qNc == mNc, comm, PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, mNc);
2023: PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
2024: PetscCall(VecGetArrayWrite(eFlux, &earray));
2025: PetscCall(PetscMalloc6(Nc * cdim, &interpolant, mNc * cdim, &minterpolant, cdim * (Nq + 1), &coords, cdim * cdim * Nq, &fegeom.J, cdim * cdim * Nq, &fegeom.invJ, Nq, &fegeom.detJ));
2026: for (PetscInt c = cStart; c < cEnd; ++c) {
2027: PetscScalar *x = NULL;
2028: PetscScalar *mx = NULL;
2029: PetscScalar *eval = NULL;
2030: PetscReal fluxElemDiff = 0.0;
2032: PetscCall(DMPlexComputeCellGeometryFEM(dm, c, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
2033: PetscCall(DMPlexVecGetClosure(dm, NULL, lu, c, NULL, &x));
2034: PetscCall(DMPlexVecGetClosure(mdm, NULL, lmu, c, NULL, &mx));
2036: for (PetscInt q = 0; q < Nq; ++q) {
2037: PetscFEGeom qgeom;
2039: qgeom.dimEmbed = fegeom.dimEmbed;
2040: qgeom.J = &fegeom.J[q * cdim * cdim];
2041: qgeom.invJ = &fegeom.invJ[q * cdim * cdim];
2042: qgeom.detJ = &fegeom.detJ[q];
2044: PetscCheck(fegeom.detJ[q] > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT ", quadrature points %" PetscInt_FMT, (double)fegeom.detJ[q], c, q);
2046: PetscCall(PetscFEInterpolate_Static(mfe, &mx[0], &qgeom, q, minterpolant));
2047: PetscCall(PetscFEInterpolateGradient_Static(fe, 1, &x[0], &qgeom, q, interpolant));
2049: /* Now take the elementwise difference and store that in a vector. */
2050: for (PetscInt fc = 0; fc < mNc; ++fc) {
2051: const PetscReal wt = quadWeights[q * qNc + (qNc == 1 ? 0 : fc)];
2052: fluxElemDiff += PetscSqr(PetscRealPart(interpolant[fc] - minterpolant[fc])) * wt * fegeom.detJ[q];
2053: }
2054: }
2055: PetscCall(DMPlexVecRestoreClosure(dm, NULL, lu, c, NULL, &x));
2056: PetscCall(DMPlexVecRestoreClosure(mdm, NULL, lmu, c, NULL, &mx));
2057: PetscCall(DMPlexPointGlobalRef(edm, c, earray, (void *)&eval));
2058: if (eval) eval[0] = fluxElemDiff;
2059: }
2060: PetscCall(PetscFree6(interpolant, minterpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
2061: PetscCall(VecRestoreArrayWrite(eFlux, &earray));
2063: PetscCall(VecAssemblyBegin(eFlux));
2064: PetscCall(VecAssemblyEnd(eFlux));
2065: PetscCall(VecSqrtAbs(eFlux));
2066: PetscFunctionReturn(PETSC_SUCCESS);
2067: }
2069: /*@
2070: DMPlexComputeL2FluxDiffVec - This function computes the integral of the difference between the gradient of field `f`in `u` and field `mf` in `mu`
2072: Collective
2074: Input Parameters:
2075: + u - The global `Vec` containing the primal solution
2076: . f - The field number for the potential
2077: . mu - The global `Vec` containing the mixed solution
2078: - mf - The field number for the flux
2080: Output Parameter:
2081: . eFlux - A global `Vec` which holds $||\nabla u_f - \mu_{mf}||$
2083: Level: advanced
2085: Notes:
2086: We assume that the `DM` for each solution has the same topology, geometry, and quadrature.
2088: This is usually used to get an error estimate for the primal solution, using the flux from a mixed solution.
2090: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexComputeL2FluxDiffVecLocal()`, `DMProjectFunction()`, `DMComputeL2Diff()`, `DMPlexComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
2091: @*/
2092: PetscErrorCode DMPlexComputeL2FluxDiffVec(Vec u, PetscInt f, Vec mu, PetscInt mf, Vec eFlux)
2093: {
2094: DM dm, mdm;
2095: Vec lu, lmu;
2097: PetscFunctionBegin;
2098: PetscCall(VecGetDM(u, &dm));
2099: PetscCall(DMGetLocalVector(dm, &lu));
2100: PetscCall(DMGlobalToLocal(dm, u, INSERT_VALUES, lu));
2101: PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, lu, 0.0, NULL, NULL, NULL));
2103: PetscCall(VecGetDM(mu, &mdm));
2104: PetscCall(DMGetLocalVector(mdm, &lmu));
2105: PetscCall(DMGlobalToLocal(mdm, mu, INSERT_VALUES, lmu));
2106: PetscCall(DMPlexInsertBoundaryValues(mdm, PETSC_TRUE, lmu, 0.0, NULL, NULL, NULL));
2108: PetscCall(DMPlexComputeL2FluxDiffVecLocal(lu, f, lmu, mf, eFlux));
2110: PetscCall(DMRestoreLocalVector(dm, &lu));
2111: PetscCall(DMRestoreLocalVector(mdm, &lmu));
2112: PetscFunctionReturn(PETSC_SUCCESS);
2113: }
2115: /*@
2116: DMPlexComputeClementInterpolant - This function computes the L2 projection of the cellwise values of a function u onto P1
2118: Collective
2120: Input Parameters:
2121: + dm - The `DM`
2122: - locX - The coefficient vector u_h
2124: Output Parameter:
2125: . locC - A `Vec` which holds the Clement interpolant of the function
2127: Level: developer
2129: Note:
2130: $ 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
2132: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectFunction()`, `DMComputeL2Diff()`, `DMPlexComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
2133: @*/
2134: PetscErrorCode DMPlexComputeClementInterpolant(DM dm, Vec locX, Vec locC)
2135: {
2136: PetscInt debug = ((DM_Plex *)dm->data)->printFEM;
2137: DM dmc;
2138: PetscQuadrature quad;
2139: PetscScalar *interpolant, *valsum;
2140: PetscFEGeom fegeom;
2141: PetscReal *coords;
2142: const PetscReal *quadPoints, *quadWeights;
2143: PetscInt dim, cdim, Nf, f, Nc = 0, Nq, qNc, cStart, cEnd, vStart, vEnd, v;
2145: PetscFunctionBegin;
2146: PetscCall(PetscCitationsRegister(ClementCitation, &Clementcite));
2147: PetscCall(VecGetDM(locC, &dmc));
2148: PetscCall(VecSet(locC, 0.0));
2149: PetscCall(DMGetDimension(dm, &dim));
2150: PetscCall(DMGetCoordinateDim(dm, &cdim));
2151: fegeom.dimEmbed = cdim;
2152: PetscCall(DMGetNumFields(dm, &Nf));
2153: PetscCheck(Nf > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fields is zero!");
2154: for (f = 0; f < Nf; ++f) {
2155: PetscObject obj;
2156: PetscClassId id;
2157: PetscInt fNc;
2159: PetscCall(DMGetField(dm, f, NULL, &obj));
2160: PetscCall(PetscObjectGetClassId(obj, &id));
2161: if (id == PETSCFE_CLASSID) {
2162: PetscFE fe = (PetscFE)obj;
2164: PetscCall(PetscFEGetQuadrature(fe, &quad));
2165: PetscCall(PetscFEGetNumComponents(fe, &fNc));
2166: } else if (id == PETSCFV_CLASSID) {
2167: PetscFV fv = (PetscFV)obj;
2169: PetscCall(PetscFVGetQuadrature(fv, &quad));
2170: PetscCall(PetscFVGetNumComponents(fv, &fNc));
2171: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
2172: Nc += fNc;
2173: }
2174: PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
2175: PetscCheck(qNc == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " > 1", qNc);
2176: PetscCall(PetscMalloc6(Nc * 2, &valsum, Nc, &interpolant, cdim * Nq, &coords, Nq, &fegeom.detJ, cdim * cdim * Nq, &fegeom.J, cdim * cdim * Nq, &fegeom.invJ));
2177: PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2178: PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
2179: for (v = vStart; v < vEnd; ++v) {
2180: PetscScalar volsum = 0.0;
2181: PetscInt *star = NULL;
2182: PetscInt starSize, st, fc;
2184: PetscCall(PetscArrayzero(valsum, Nc));
2185: PetscCall(DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
2186: for (st = 0; st < starSize * 2; st += 2) {
2187: const PetscInt cell = star[st];
2188: PetscScalar *val = &valsum[Nc];
2189: PetscScalar *x = NULL;
2190: PetscReal vol = 0.0;
2191: PetscInt foff = 0;
2193: if ((cell < cStart) || (cell >= cEnd)) continue;
2194: PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
2195: PetscCall(DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x));
2196: for (f = 0; f < Nf; ++f) {
2197: PetscObject obj;
2198: PetscClassId id;
2199: PetscInt Nb, fNc, q;
2201: PetscCall(PetscArrayzero(val, Nc));
2202: PetscCall(DMGetField(dm, f, NULL, &obj));
2203: PetscCall(PetscObjectGetClassId(obj, &id));
2204: if (id == PETSCFE_CLASSID) {
2205: PetscCall(PetscFEGetNumComponents((PetscFE)obj, &fNc));
2206: PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb));
2207: } else if (id == PETSCFV_CLASSID) {
2208: PetscCall(PetscFVGetNumComponents((PetscFV)obj, &fNc));
2209: Nb = 1;
2210: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
2211: for (q = 0; q < Nq; ++q) {
2212: const PetscReal wt = quadWeights[q] * fegeom.detJ[q];
2213: PetscFEGeom qgeom;
2215: qgeom.dimEmbed = fegeom.dimEmbed;
2216: qgeom.J = &fegeom.J[q * cdim * cdim];
2217: qgeom.invJ = &fegeom.invJ[q * cdim * cdim];
2218: qgeom.detJ = &fegeom.detJ[q];
2219: PetscCheck(fegeom.detJ[q] > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT ", quadrature points %" PetscInt_FMT, (double)fegeom.detJ[q], cell, q);
2220: PetscCheck(id == PETSCFE_CLASSID, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
2221: PetscCall(PetscFEInterpolate_Static((PetscFE)obj, &x[foff], &qgeom, q, interpolant));
2222: for (fc = 0; fc < fNc; ++fc) val[foff + fc] += interpolant[fc] * wt;
2223: vol += wt;
2224: }
2225: foff += Nb;
2226: }
2227: PetscCall(DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x));
2228: for (fc = 0; fc < Nc; ++fc) valsum[fc] += val[fc];
2229: volsum += vol;
2230: if (debug) {
2231: PetscCall(PetscPrintf(PETSC_COMM_SELF, "Vertex %" PetscInt_FMT " Cell %" PetscInt_FMT " value: [", v, cell));
2232: for (fc = 0; fc < Nc; ++fc) {
2233: if (fc) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
2234: PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(val[fc])));
2235: }
2236: PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
2237: }
2238: }
2239: for (fc = 0; fc < Nc; ++fc) valsum[fc] /= volsum;
2240: PetscCall(DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
2241: PetscCall(DMPlexVecSetClosure(dmc, NULL, locC, v, valsum, INSERT_VALUES));
2242: }
2243: PetscCall(PetscFree6(valsum, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
2244: PetscFunctionReturn(PETSC_SUCCESS);
2245: }
2247: /*@
2248: DMPlexComputeGradientClementInterpolant - This function computes the L2 projection of the cellwise gradient of a function u onto P1
2250: Collective
2252: Input Parameters:
2253: + dm - The `DM`
2254: - locX - The coefficient vector u_h
2256: Output Parameter:
2257: . locC - A `Vec` which holds the Clement interpolant of the gradient
2259: Level: developer
2261: Note:
2262: $\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
2264: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMProjectFunction()`, `DMComputeL2Diff()`, `DMPlexComputeL2FieldDiff()`, `DMComputeL2GradientDiff()`
2265: @*/
2266: PetscErrorCode DMPlexComputeGradientClementInterpolant(DM dm, Vec locX, Vec locC)
2267: {
2268: DM_Plex *mesh = (DM_Plex *)dm->data;
2269: PetscInt debug = mesh->printFEM;
2270: DM dmC;
2271: PetscQuadrature quad;
2272: PetscScalar *interpolant, *gradsum;
2273: PetscFEGeom fegeom;
2274: PetscReal *coords;
2275: const PetscReal *quadPoints, *quadWeights;
2276: PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, vStart, vEnd, v, field, fieldOffset;
2278: PetscFunctionBegin;
2279: PetscCall(PetscCitationsRegister(ClementCitation, &Clementcite));
2280: PetscCall(VecGetDM(locC, &dmC));
2281: PetscCall(VecSet(locC, 0.0));
2282: PetscCall(DMGetDimension(dm, &dim));
2283: PetscCall(DMGetCoordinateDim(dm, &coordDim));
2284: fegeom.dimEmbed = coordDim;
2285: PetscCall(DMGetNumFields(dm, &numFields));
2286: PetscCheck(numFields, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fields is zero!");
2287: for (field = 0; field < numFields; ++field) {
2288: PetscObject obj;
2289: PetscClassId id;
2290: PetscInt Nc;
2292: PetscCall(DMGetField(dm, field, NULL, &obj));
2293: PetscCall(PetscObjectGetClassId(obj, &id));
2294: if (id == PETSCFE_CLASSID) {
2295: PetscFE fe = (PetscFE)obj;
2297: PetscCall(PetscFEGetQuadrature(fe, &quad));
2298: PetscCall(PetscFEGetNumComponents(fe, &Nc));
2299: } else if (id == PETSCFV_CLASSID) {
2300: PetscFV fv = (PetscFV)obj;
2302: PetscCall(PetscFVGetQuadrature(fv, &quad));
2303: PetscCall(PetscFVGetNumComponents(fv, &Nc));
2304: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
2305: numComponents += Nc;
2306: }
2307: PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
2308: PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents);
2309: PetscCall(PetscMalloc6(coordDim * numComponents * 2, &gradsum, coordDim * numComponents, &interpolant, coordDim * Nq, &coords, Nq, &fegeom.detJ, coordDim * coordDim * Nq, &fegeom.J, coordDim * coordDim * Nq, &fegeom.invJ));
2310: PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2311: PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
2312: for (v = vStart; v < vEnd; ++v) {
2313: PetscScalar volsum = 0.0;
2314: PetscInt *star = NULL;
2315: PetscInt starSize, st, d, fc;
2317: PetscCall(PetscArrayzero(gradsum, coordDim * numComponents));
2318: PetscCall(DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
2319: for (st = 0; st < starSize * 2; st += 2) {
2320: const PetscInt cell = star[st];
2321: PetscScalar *grad = &gradsum[coordDim * numComponents];
2322: PetscScalar *x = NULL;
2323: PetscReal vol = 0.0;
2325: if ((cell < cStart) || (cell >= cEnd)) continue;
2326: PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
2327: PetscCall(DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x));
2328: for (field = 0, fieldOffset = 0; field < numFields; ++field) {
2329: PetscObject obj;
2330: PetscClassId id;
2331: PetscInt Nb, Nc, q, qc = 0;
2333: PetscCall(PetscArrayzero(grad, coordDim * numComponents));
2334: PetscCall(DMGetField(dm, field, NULL, &obj));
2335: PetscCall(PetscObjectGetClassId(obj, &id));
2336: if (id == PETSCFE_CLASSID) {
2337: PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
2338: PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb));
2339: } else if (id == PETSCFV_CLASSID) {
2340: PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
2341: Nb = 1;
2342: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
2343: for (q = 0; q < Nq; ++q) {
2344: PetscFEGeom qgeom;
2346: qgeom.dimEmbed = fegeom.dimEmbed;
2347: qgeom.J = &fegeom.J[q * coordDim * coordDim];
2348: qgeom.invJ = &fegeom.invJ[q * coordDim * coordDim];
2349: qgeom.detJ = &fegeom.detJ[q];
2350: PetscCheck(fegeom.detJ[q] > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT ", quadrature points %" PetscInt_FMT, (double)fegeom.detJ[q], cell, q);
2351: PetscCheck(id == PETSCFE_CLASSID, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
2352: PetscCall(PetscFEInterpolateGradient_Static((PetscFE)obj, 1, &x[fieldOffset], &qgeom, q, interpolant));
2353: for (fc = 0; fc < Nc; ++fc) {
2354: const PetscReal wt = quadWeights[q * qNc + qc];
2356: for (d = 0; d < coordDim; ++d) grad[fc * coordDim + d] += interpolant[fc * dim + d] * wt * fegeom.detJ[q];
2357: }
2358: vol += quadWeights[q * qNc] * fegeom.detJ[q];
2359: }
2360: fieldOffset += Nb;
2361: qc += Nc;
2362: }
2363: PetscCall(DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x));
2364: for (fc = 0; fc < numComponents; ++fc) {
2365: for (d = 0; d < coordDim; ++d) gradsum[fc * coordDim + d] += grad[fc * coordDim + d];
2366: }
2367: volsum += vol;
2368: if (debug) {
2369: PetscCall(PetscPrintf(PETSC_COMM_SELF, "Vertex %" PetscInt_FMT " Cell %" PetscInt_FMT " gradient: [", v, cell));
2370: for (fc = 0; fc < numComponents; ++fc) {
2371: for (d = 0; d < coordDim; ++d) {
2372: if (fc || d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
2373: PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc * coordDim + d])));
2374: }
2375: }
2376: PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
2377: }
2378: }
2379: for (fc = 0; fc < numComponents; ++fc) {
2380: for (d = 0; d < coordDim; ++d) gradsum[fc * coordDim + d] /= volsum;
2381: }
2382: PetscCall(DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
2383: PetscCall(DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES));
2384: }
2385: PetscCall(PetscFree6(gradsum, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
2386: PetscFunctionReturn(PETSC_SUCCESS);
2387: }
2389: PetscErrorCode DMPlexComputeIntegral_Internal(DM dm, Vec locX, PetscInt cStart, PetscInt cEnd, PetscScalar *cintegral, void *user)
2390: {
2391: DM dmAux = NULL, plexA = NULL;
2392: PetscDS prob, probAux = NULL;
2393: PetscSection section, sectionAux;
2394: Vec locA;
2395: PetscInt dim, numCells = cEnd - cStart, c, f;
2396: PetscBool useFVM = PETSC_FALSE;
2397: /* DS */
2398: PetscInt Nf, totDim, *uOff, *uOff_x, numConstants;
2399: PetscInt NfAux, totDimAux, *aOff;
2400: PetscScalar *u, *a = NULL;
2401: const PetscScalar *constants;
2402: /* Geometry */
2403: PetscFEGeom *cgeomFEM;
2404: DM dmGrad;
2405: PetscQuadrature affineQuad = NULL;
2406: Vec cellGeometryFVM = NULL, faceGeometryFVM = NULL, locGrad = NULL;
2407: PetscFVCellGeom *cgeomFVM;
2408: const PetscScalar *lgrad;
2409: PetscInt maxDegree;
2410: DMField coordField;
2411: IS cellIS;
2413: PetscFunctionBegin;
2414: PetscCall(DMGetDS(dm, &prob));
2415: PetscCall(DMGetDimension(dm, &dim));
2416: PetscCall(DMGetLocalSection(dm, §ion));
2417: PetscCall(DMGetNumFields(dm, &Nf));
2418: /* Determine which discretizations we have */
2419: for (f = 0; f < Nf; ++f) {
2420: PetscObject obj;
2421: PetscClassId id;
2423: PetscCall(PetscDSGetDiscretization(prob, f, &obj));
2424: PetscCall(PetscObjectGetClassId(obj, &id));
2425: if (id == PETSCFV_CLASSID) useFVM = PETSC_TRUE;
2426: }
2427: /* Read DS information */
2428: PetscCall(PetscDSGetTotalDimension(prob, &totDim));
2429: PetscCall(PetscDSGetComponentOffsets(prob, &uOff));
2430: PetscCall(PetscDSGetComponentDerivativeOffsets(prob, &uOff_x));
2431: PetscCall(ISCreateStride(PETSC_COMM_SELF, numCells, cStart, 1, &cellIS));
2432: PetscCall(PetscDSGetConstants(prob, &numConstants, &constants));
2433: /* Read Auxiliary DS information */
2434: PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &locA));
2435: if (locA) {
2436: PetscCall(VecGetDM(locA, &dmAux));
2437: PetscCall(DMConvert(dmAux, DMPLEX, &plexA));
2438: PetscCall(DMGetDS(dmAux, &probAux));
2439: PetscCall(PetscDSGetNumFields(probAux, &NfAux));
2440: PetscCall(DMGetLocalSection(dmAux, §ionAux));
2441: PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
2442: PetscCall(PetscDSGetComponentOffsets(probAux, &aOff));
2443: }
2444: /* Allocate data arrays */
2445: PetscCall(PetscCalloc1(numCells * totDim, &u));
2446: if (dmAux) PetscCall(PetscMalloc1(numCells * totDimAux, &a));
2447: /* Read out geometry */
2448: PetscCall(DMGetCoordinateField(dm, &coordField));
2449: PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
2450: if (maxDegree <= 1) {
2451: PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, &affineQuad));
2452: if (affineQuad) PetscCall(DMFieldCreateFEGeom(coordField, cellIS, affineQuad, PETSC_FEGEOM_BASIC, &cgeomFEM));
2453: }
2454: if (useFVM) {
2455: PetscFV fv = NULL;
2456: Vec grad;
2457: PetscInt fStart, fEnd;
2458: PetscBool compGrad;
2460: for (f = 0; f < Nf; ++f) {
2461: PetscObject obj;
2462: PetscClassId id;
2464: PetscCall(PetscDSGetDiscretization(prob, f, &obj));
2465: PetscCall(PetscObjectGetClassId(obj, &id));
2466: if (id == PETSCFV_CLASSID) {
2467: fv = (PetscFV)obj;
2468: break;
2469: }
2470: }
2471: PetscCall(PetscFVGetComputeGradients(fv, &compGrad));
2472: PetscCall(PetscFVSetComputeGradients(fv, PETSC_TRUE));
2473: PetscCall(DMPlexComputeGeometryFVM(dm, &cellGeometryFVM, &faceGeometryFVM));
2474: PetscCall(DMPlexComputeGradientFVM(dm, fv, faceGeometryFVM, cellGeometryFVM, &dmGrad));
2475: PetscCall(PetscFVSetComputeGradients(fv, compGrad));
2476: PetscCall(VecGetArrayRead(cellGeometryFVM, (const PetscScalar **)&cgeomFVM));
2477: /* Reconstruct and limit cell gradients */
2478: PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
2479: PetscCall(DMGetGlobalVector(dmGrad, &grad));
2480: PetscCall(DMPlexReconstructGradients_Internal(dm, fv, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad));
2481: /* Communicate gradient values */
2482: PetscCall(DMGetLocalVector(dmGrad, &locGrad));
2483: PetscCall(DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad));
2484: PetscCall(DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad));
2485: PetscCall(DMRestoreGlobalVector(dmGrad, &grad));
2486: /* Handle non-essential (e.g. outflow) boundary values */
2487: PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_FALSE, locX, 0.0, faceGeometryFVM, cellGeometryFVM, locGrad));
2488: PetscCall(VecGetArrayRead(locGrad, &lgrad));
2489: }
2490: /* Read out data from inputs */
2491: for (c = cStart; c < cEnd; ++c) {
2492: PetscScalar *x = NULL;
2493: PetscInt i;
2495: PetscCall(DMPlexVecGetClosure(dm, section, locX, c, NULL, &x));
2496: for (i = 0; i < totDim; ++i) u[c * totDim + i] = x[i];
2497: PetscCall(DMPlexVecRestoreClosure(dm, section, locX, c, NULL, &x));
2498: if (dmAux) {
2499: PetscCall(DMPlexVecGetClosure(plexA, sectionAux, locA, c, NULL, &x));
2500: for (i = 0; i < totDimAux; ++i) a[c * totDimAux + i] = x[i];
2501: PetscCall(DMPlexVecRestoreClosure(plexA, sectionAux, locA, c, NULL, &x));
2502: }
2503: }
2504: /* Do integration for each field */
2505: for (f = 0; f < Nf; ++f) {
2506: PetscObject obj;
2507: PetscClassId id;
2508: PetscInt numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;
2510: PetscCall(PetscDSGetDiscretization(prob, f, &obj));
2511: PetscCall(PetscObjectGetClassId(obj, &id));
2512: if (id == PETSCFE_CLASSID) {
2513: PetscFE fe = (PetscFE)obj;
2514: PetscQuadrature q;
2515: PetscFEGeom *chunkGeom = NULL;
2516: PetscInt Nq, Nb;
2518: PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
2519: PetscCall(PetscFEGetQuadrature(fe, &q));
2520: PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL));
2521: PetscCall(PetscFEGetDimension(fe, &Nb));
2522: blockSize = Nb * Nq;
2523: batchSize = numBlocks * blockSize;
2524: PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
2525: numChunks = numCells / (numBatches * batchSize);
2526: Ne = numChunks * numBatches * batchSize;
2527: Nr = numCells % (numBatches * batchSize);
2528: offset = numCells - Nr;
2529: if (!affineQuad) PetscCall(DMFieldCreateFEGeom(coordField, cellIS, q, PETSC_FEGEOM_BASIC, &cgeomFEM));
2530: PetscCall(PetscFEGeomGetChunk(cgeomFEM, 0, offset, &chunkGeom));
2531: PetscCall(PetscFEIntegrate(prob, f, Ne, chunkGeom, u, probAux, a, cintegral));
2532: PetscCall(PetscFEGeomGetChunk(cgeomFEM, offset, numCells, &chunkGeom));
2533: PetscCall(PetscFEIntegrate(prob, f, Nr, chunkGeom, &u[offset * totDim], probAux, PetscSafePointerPlusOffset(a, offset * totDimAux), &cintegral[offset * Nf]));
2534: PetscCall(PetscFEGeomRestoreChunk(cgeomFEM, offset, numCells, &chunkGeom));
2535: if (!affineQuad) PetscCall(PetscFEGeomDestroy(&cgeomFEM));
2536: } else if (id == PETSCFV_CLASSID) {
2537: PetscInt foff;
2538: PetscPointFn *obj_func;
2540: PetscCall(PetscDSGetObjective(prob, f, &obj_func));
2541: PetscCall(PetscDSGetFieldOffset(prob, f, &foff));
2542: if (obj_func) {
2543: for (c = 0; c < numCells; ++c) {
2544: PetscScalar *u_x;
2545: PetscScalar lint = 0.;
2547: PetscCall(DMPlexPointLocalRead(dmGrad, c, lgrad, &u_x));
2548: obj_func(dim, Nf, NfAux, uOff, uOff_x, &u[totDim * c + foff], NULL, u_x, aOff, NULL, PetscSafePointerPlusOffset(a, totDimAux * c), NULL, NULL, 0.0, cgeomFVM[c].centroid, numConstants, constants, &lint);
2549: cintegral[c * Nf + f] += PetscRealPart(lint) * cgeomFVM[c].volume;
2550: }
2551: }
2552: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
2553: }
2554: /* Cleanup data arrays */
2555: if (useFVM) {
2556: PetscCall(VecRestoreArrayRead(locGrad, &lgrad));
2557: PetscCall(VecRestoreArrayRead(cellGeometryFVM, (const PetscScalar **)&cgeomFVM));
2558: PetscCall(DMRestoreLocalVector(dmGrad, &locGrad));
2559: PetscCall(VecDestroy(&faceGeometryFVM));
2560: PetscCall(VecDestroy(&cellGeometryFVM));
2561: PetscCall(DMDestroy(&dmGrad));
2562: }
2563: if (dmAux) PetscCall(PetscFree(a));
2564: PetscCall(DMDestroy(&plexA));
2565: PetscCall(PetscFree(u));
2566: /* Cleanup */
2567: if (affineQuad) PetscCall(PetscFEGeomDestroy(&cgeomFEM));
2568: PetscCall(PetscQuadratureDestroy(&affineQuad));
2569: PetscCall(ISDestroy(&cellIS));
2570: PetscFunctionReturn(PETSC_SUCCESS);
2571: }
2573: /*@
2574: DMPlexComputeIntegralFEM - Form the integral over the domain from the global input X using pointwise functions specified by the user
2576: Input Parameters:
2577: + dm - The mesh
2578: . X - Global input vector
2579: - user - The user context
2581: Output Parameter:
2582: . integral - Integral for each field
2584: Level: developer
2586: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSNESComputeResidualFEM()`
2587: @*/
2588: PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscScalar *integral, void *user)
2589: {
2590: PetscInt printFEM;
2591: PetscScalar *cintegral, *lintegral;
2592: PetscInt Nf, f, cellHeight, cStart, cEnd, cell;
2593: Vec locX;
2595: PetscFunctionBegin;
2598: PetscAssertPointer(integral, 3);
2599: PetscCall(PetscLogEventBegin(DMPLEX_IntegralFEM, dm, 0, 0, 0));
2600: PetscCall(DMPlexConvertPlex(dm, &dm, PETSC_TRUE));
2601: PetscCall(DMGetNumFields(dm, &Nf));
2602: PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
2603: PetscCall(DMPlexGetSimplexOrBoxCells(dm, cellHeight, &cStart, &cEnd));
2604: /* TODO Introduce a loop over large chunks (right now this is a single chunk) */
2605: PetscCall(PetscCalloc2(Nf, &lintegral, (cEnd - cStart) * Nf, &cintegral));
2606: /* Get local solution with boundary values */
2607: PetscCall(DMGetLocalVector(dm, &locX));
2608: PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locX, 0.0, NULL, NULL, NULL));
2609: PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
2610: PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
2611: PetscCall(DMPlexComputeIntegral_Internal(dm, locX, cStart, cEnd, cintegral, user));
2612: PetscCall(DMRestoreLocalVector(dm, &locX));
2613: printFEM = ((DM_Plex *)dm->data)->printFEM;
2614: /* Sum up values */
2615: for (cell = cStart; cell < cEnd; ++cell) {
2616: const PetscInt c = cell - cStart;
2618: if (printFEM > 1) PetscCall(DMPrintCellVector(cell, "Cell Integral", Nf, &cintegral[c * Nf]));
2619: for (f = 0; f < Nf; ++f) lintegral[f] += cintegral[c * Nf + f];
2620: }
2621: PetscCallMPI(MPIU_Allreduce(lintegral, integral, Nf, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)dm)));
2622: if (printFEM) {
2623: PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Integral:"));
2624: for (f = 0; f < Nf; ++f) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), " %g", (double)PetscRealPart(integral[f])));
2625: PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "\n"));
2626: }
2627: PetscCall(PetscFree2(lintegral, cintegral));
2628: PetscCall(PetscLogEventEnd(DMPLEX_IntegralFEM, dm, 0, 0, 0));
2629: PetscCall(DMDestroy(&dm));
2630: PetscFunctionReturn(PETSC_SUCCESS);
2631: }
2633: /*@
2634: DMPlexComputeCellwiseIntegralFEM - Form the vector of cellwise integrals F from the global input X using pointwise functions specified by the user
2636: Input Parameters:
2637: + dm - The mesh
2638: . X - Global input vector
2639: - user - The user context
2641: Output Parameter:
2642: . F - Cellwise integrals for each field
2644: Level: developer
2646: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSNESComputeResidualFEM()`
2647: @*/
2648: PetscErrorCode DMPlexComputeCellwiseIntegralFEM(DM dm, Vec X, Vec F, void *user)
2649: {
2650: PetscInt printFEM;
2651: DM dmF;
2652: PetscSection sectionF = NULL;
2653: PetscScalar *cintegral, *af;
2654: PetscInt Nf, f, cellHeight, cStart, cEnd, cell, n;
2655: Vec locX;
2657: PetscFunctionBegin;
2661: PetscCall(PetscLogEventBegin(DMPLEX_IntegralFEM, dm, 0, 0, 0));
2662: PetscCall(DMPlexConvertPlex(dm, &dm, PETSC_TRUE));
2663: PetscCall(DMGetNumFields(dm, &Nf));
2664: PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight));
2665: PetscCall(DMPlexGetSimplexOrBoxCells(dm, cellHeight, &cStart, &cEnd));
2666: /* TODO Introduce a loop over large chunks (right now this is a single chunk) */
2667: PetscCall(PetscCalloc1((cEnd - cStart) * Nf, &cintegral));
2668: /* Get local solution with boundary values */
2669: PetscCall(DMGetLocalVector(dm, &locX));
2670: PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locX, 0.0, NULL, NULL, NULL));
2671: PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
2672: PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
2673: PetscCall(DMPlexComputeIntegral_Internal(dm, locX, cStart, cEnd, cintegral, user));
2674: PetscCall(DMRestoreLocalVector(dm, &locX));
2675: /* Put values in F */
2676: PetscCall(VecGetArray(F, &af));
2677: PetscCall(VecGetDM(F, &dmF));
2678: if (dmF) PetscCall(DMGetLocalSection(dmF, §ionF));
2679: PetscCall(VecGetLocalSize(F, &n));
2680: PetscCheck(n >= (cEnd - cStart) * Nf, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Vector size %" PetscInt_FMT " < %" PetscInt_FMT, n, (cEnd - cStart) * Nf);
2681: printFEM = ((DM_Plex *)dm->data)->printFEM;
2682: for (cell = cStart; cell < cEnd; ++cell) {
2683: const PetscInt c = cell - cStart;
2684: PetscInt dof = Nf, off = c * Nf;
2686: if (printFEM > 1) PetscCall(DMPrintCellVector(cell, "Cell Integral", Nf, &cintegral[c * Nf]));
2687: if (sectionF) {
2688: PetscCall(PetscSectionGetDof(sectionF, cell, &dof));
2689: PetscCall(PetscSectionGetOffset(sectionF, cell, &off));
2690: }
2691: PetscCheck(dof == Nf, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of cell dofs %" PetscInt_FMT " != %" PetscInt_FMT, dof, Nf);
2692: for (f = 0; f < Nf; ++f) af[off + f] = cintegral[c * Nf + f];
2693: }
2694: PetscCall(VecRestoreArray(F, &af));
2695: PetscCall(PetscFree(cintegral));
2696: PetscCall(PetscLogEventEnd(DMPLEX_IntegralFEM, dm, 0, 0, 0));
2697: PetscCall(DMDestroy(&dm));
2698: PetscFunctionReturn(PETSC_SUCCESS);
2699: }
2701: static PetscErrorCode DMPlexComputeBdIntegral_Internal(DM dm, Vec locX, IS pointIS, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscScalar *fintegral, void *user)
2702: {
2703: DM plex = NULL, plexA = NULL;
2704: DMEnclosureType encAux;
2705: PetscDS prob, probAux = NULL;
2706: PetscSection section, sectionAux = NULL;
2707: Vec locA = NULL;
2708: DMField coordField;
2709: PetscInt Nf, totDim, *uOff, *uOff_x;
2710: PetscInt NfAux = 0, totDimAux = 0, *aOff = NULL;
2711: PetscScalar *u, *a = NULL;
2712: const PetscScalar *constants;
2713: PetscInt numConstants, f;
2715: PetscFunctionBegin;
2716: PetscCall(DMGetCoordinateField(dm, &coordField));
2717: PetscCall(DMConvert(dm, DMPLEX, &plex));
2718: PetscCall(DMGetDS(dm, &prob));
2719: PetscCall(DMGetLocalSection(dm, §ion));
2720: PetscCall(PetscSectionGetNumFields(section, &Nf));
2721: /* Determine which discretizations we have */
2722: for (f = 0; f < Nf; ++f) {
2723: PetscObject obj;
2724: PetscClassId id;
2726: PetscCall(PetscDSGetDiscretization(prob, f, &obj));
2727: PetscCall(PetscObjectGetClassId(obj, &id));
2728: PetscCheck(id != PETSCFV_CLASSID, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not supported for FVM (field %" PetscInt_FMT ")", f);
2729: }
2730: /* Read DS information */
2731: PetscCall(PetscDSGetTotalDimension(prob, &totDim));
2732: PetscCall(PetscDSGetComponentOffsets(prob, &uOff));
2733: PetscCall(PetscDSGetComponentDerivativeOffsets(prob, &uOff_x));
2734: PetscCall(PetscDSGetConstants(prob, &numConstants, &constants));
2735: /* Read Auxiliary DS information */
2736: PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &locA));
2737: if (locA) {
2738: DM dmAux;
2740: PetscCall(VecGetDM(locA, &dmAux));
2741: PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
2742: PetscCall(DMConvert(dmAux, DMPLEX, &plexA));
2743: PetscCall(DMGetDS(dmAux, &probAux));
2744: PetscCall(PetscDSGetNumFields(probAux, &NfAux));
2745: PetscCall(DMGetLocalSection(dmAux, §ionAux));
2746: PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
2747: PetscCall(PetscDSGetComponentOffsets(probAux, &aOff));
2748: }
2749: /* Integrate over points */
2750: {
2751: PetscFEGeom *fgeom, *chunkGeom = NULL;
2752: PetscInt maxDegree;
2753: PetscQuadrature qGeom = NULL;
2754: const PetscInt *points;
2755: PetscInt numFaces, face, Nq, field;
2756: PetscInt numChunks, chunkSize, chunk, Nr, offset;
2758: PetscCall(ISGetLocalSize(pointIS, &numFaces));
2759: PetscCall(ISGetIndices(pointIS, &points));
2760: PetscCall(PetscCalloc2(numFaces * totDim, &u, (locA ? (size_t)numFaces * totDimAux : 0), &a));
2761: PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
2762: for (face = 0; face < numFaces; ++face) {
2763: const PetscInt point = points[face], *support;
2764: PetscScalar *x = NULL;
2766: PetscCall(DMPlexGetSupport(dm, point, &support));
2767: PetscCall(DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x));
2768: for (PetscInt i = 0; i < totDim; ++i) u[face * totDim + i] = x[i];
2769: PetscCall(DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x));
2770: if (locA) {
2771: PetscInt subp;
2772: PetscCall(DMGetEnclosurePoint(plexA, dm, encAux, support[0], &subp));
2773: PetscCall(DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x));
2774: for (PetscInt i = 0; i < totDimAux; ++i) a[f * totDimAux + i] = x[i];
2775: PetscCall(DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x));
2776: }
2777: }
2778: for (field = 0; field < Nf; ++field) {
2779: PetscFE fe;
2781: PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe));
2782: if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &qGeom));
2783: if (!qGeom) {
2784: PetscCall(PetscFEGetFaceQuadrature(fe, &qGeom));
2785: PetscCall(PetscObjectReference((PetscObject)qGeom));
2786: }
2787: PetscCall(PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL));
2788: PetscCall(DMPlexGetFEGeom(coordField, pointIS, qGeom, PETSC_FEGEOM_BOUNDARY, &fgeom));
2789: /* Get blocking */
2790: {
2791: PetscQuadrature q;
2792: PetscInt numBatches, batchSize, numBlocks, blockSize;
2793: PetscInt Nq, Nb;
2795: PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
2796: PetscCall(PetscFEGetQuadrature(fe, &q));
2797: PetscCall(PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL));
2798: PetscCall(PetscFEGetDimension(fe, &Nb));
2799: blockSize = Nb * Nq;
2800: batchSize = numBlocks * blockSize;
2801: chunkSize = numBatches * batchSize;
2802: PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
2803: numChunks = numFaces / chunkSize;
2804: Nr = numFaces % chunkSize;
2805: offset = numFaces - Nr;
2806: }
2807: /* Do integration for each field */
2808: for (chunk = 0; chunk < numChunks; ++chunk) {
2809: PetscCall(PetscFEGeomGetChunk(fgeom, chunk * chunkSize, (chunk + 1) * chunkSize, &chunkGeom));
2810: PetscCall(PetscFEIntegrateBd(prob, field, funcs[field], chunkSize, chunkGeom, &u[chunk * chunkSize * totDim], probAux, PetscSafePointerPlusOffset(a, chunk * chunkSize * totDimAux), &fintegral[chunk * chunkSize * Nf]));
2811: PetscCall(PetscFEGeomRestoreChunk(fgeom, 0, offset, &chunkGeom));
2812: }
2813: PetscCall(PetscFEGeomGetChunk(fgeom, offset, numFaces, &chunkGeom));
2814: PetscCall(PetscFEIntegrateBd(prob, field, funcs[field], Nr, chunkGeom, &u[offset * totDim], probAux, PetscSafePointerPlusOffset(a, offset * totDimAux), &fintegral[offset * Nf]));
2815: PetscCall(PetscFEGeomRestoreChunk(fgeom, offset, numFaces, &chunkGeom));
2816: /* Cleanup data arrays */
2817: PetscCall(DMPlexRestoreFEGeom(coordField, pointIS, qGeom, PETSC_FEGEOM_BOUNDARY, &fgeom));
2818: PetscCall(PetscQuadratureDestroy(&qGeom));
2819: }
2820: PetscCall(PetscFree2(u, a));
2821: PetscCall(ISRestoreIndices(pointIS, &points));
2822: }
2823: if (plex) PetscCall(DMDestroy(&plex));
2824: if (plexA) PetscCall(DMDestroy(&plexA));
2825: PetscFunctionReturn(PETSC_SUCCESS);
2826: }
2828: /*@C
2829: DMPlexComputeBdIntegral - Form the integral over the specified boundary from the global input X using pointwise functions specified by the user
2831: Input Parameters:
2832: + dm - The mesh
2833: . X - Global input vector
2834: . label - The boundary `DMLabel`
2835: . numVals - The number of label values to use, or `PETSC_DETERMINE` for all values
2836: . vals - The label values to use, or NULL for all values
2837: . funcs - The functions to integrate along the boundary for each field
2838: - user - The user context
2840: Output Parameter:
2841: . integral - Integral for each field
2843: Level: developer
2845: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexComputeIntegralFEM()`, `DMPlexComputeBdResidualFEM()`
2846: @*/
2847: PetscErrorCode DMPlexComputeBdIntegral(DM dm, Vec X, DMLabel label, PetscInt numVals, const PetscInt vals[], void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscScalar *integral, void *user)
2848: {
2849: Vec locX;
2850: PetscSection section;
2851: DMLabel depthLabel;
2852: IS facetIS;
2853: PetscInt dim, Nf, f, v;
2855: PetscFunctionBegin;
2859: if (vals) PetscAssertPointer(vals, 5);
2860: PetscAssertPointer(integral, 7);
2861: PetscCall(PetscLogEventBegin(DMPLEX_IntegralFEM, dm, 0, 0, 0));
2862: PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
2863: PetscCall(DMGetDimension(dm, &dim));
2864: PetscCall(DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS));
2865: PetscCall(DMGetLocalSection(dm, §ion));
2866: PetscCall(PetscSectionGetNumFields(section, &Nf));
2867: /* Get local solution with boundary values */
2868: PetscCall(DMGetLocalVector(dm, &locX));
2869: PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locX, 0.0, NULL, NULL, NULL));
2870: PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX));
2871: PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX));
2872: /* Loop over label values */
2873: PetscCall(PetscArrayzero(integral, Nf));
2874: for (v = 0; v < numVals; ++v) {
2875: IS pointIS;
2876: PetscInt numFaces, face;
2877: PetscScalar *fintegral;
2879: PetscCall(DMLabelGetStratumIS(label, vals[v], &pointIS));
2880: if (!pointIS) continue; /* No points with that id on this process */
2881: {
2882: IS isectIS;
2884: /* TODO: Special cases of ISIntersect where it is quick to check a priori if one is a superset of the other */
2885: PetscCall(ISIntersect_Caching_Internal(facetIS, pointIS, &isectIS));
2886: PetscCall(ISDestroy(&pointIS));
2887: pointIS = isectIS;
2888: }
2889: PetscCall(ISGetLocalSize(pointIS, &numFaces));
2890: PetscCall(PetscCalloc1(numFaces * Nf, &fintegral));
2891: PetscCall(DMPlexComputeBdIntegral_Internal(dm, locX, pointIS, funcs, fintegral, user));
2892: /* Sum point contributions into integral */
2893: for (f = 0; f < Nf; ++f)
2894: for (face = 0; face < numFaces; ++face) integral[f] += fintegral[face * Nf + f];
2895: PetscCall(PetscFree(fintegral));
2896: PetscCall(ISDestroy(&pointIS));
2897: }
2898: PetscCall(DMRestoreLocalVector(dm, &locX));
2899: PetscCall(ISDestroy(&facetIS));
2900: PetscCall(PetscLogEventEnd(DMPLEX_IntegralFEM, dm, 0, 0, 0));
2901: PetscFunctionReturn(PETSC_SUCCESS);
2902: }
2904: /*@
2905: DMPlexComputeInterpolatorNested - Form the local portion of the interpolation matrix from the coarse `DM` to a uniformly refined `DM`.
2907: Input Parameters:
2908: + dmc - The coarse mesh
2909: . dmf - The fine mesh
2910: . isRefined - Flag indicating regular refinement, rather than the same topology
2911: - user - The user context
2913: Output Parameter:
2914: . In - The interpolation matrix
2916: Level: developer
2918: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexComputeInterpolatorGeneral()`
2919: @*/
2920: PetscErrorCode DMPlexComputeInterpolatorNested(DM dmc, DM dmf, PetscBool isRefined, Mat In, void *user)
2921: {
2922: DM_Plex *mesh = (DM_Plex *)dmc->data;
2923: const char *name = "Interpolator";
2924: PetscFE *feRef;
2925: PetscFV *fvRef;
2926: PetscSection fsection, fglobalSection;
2927: PetscSection csection, cglobalSection;
2928: PetscScalar *elemMat;
2929: PetscInt dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
2930: PetscInt cTotDim = 0, rTotDim = 0;
2932: PetscFunctionBegin;
2933: PetscCall(PetscLogEventBegin(DMPLEX_InterpolatorFEM, dmc, dmf, 0, 0));
2934: PetscCall(DMGetDimension(dmf, &dim));
2935: PetscCall(DMGetLocalSection(dmf, &fsection));
2936: PetscCall(DMGetGlobalSection(dmf, &fglobalSection));
2937: PetscCall(DMGetLocalSection(dmc, &csection));
2938: PetscCall(DMGetGlobalSection(dmc, &cglobalSection));
2939: PetscCall(PetscSectionGetNumFields(fsection, &Nf));
2940: PetscCall(DMPlexGetSimplexOrBoxCells(dmc, 0, &cStart, &cEnd));
2941: PetscCall(PetscCalloc2(Nf, &feRef, Nf, &fvRef));
2942: for (f = 0; f < Nf; ++f) {
2943: PetscObject obj, objc;
2944: PetscClassId id, idc;
2945: PetscInt rNb = 0, Nc = 0, cNb = 0;
2947: PetscCall(DMGetField(dmf, f, NULL, &obj));
2948: PetscCall(PetscObjectGetClassId(obj, &id));
2949: if (id == PETSCFE_CLASSID) {
2950: PetscFE fe = (PetscFE)obj;
2952: if (isRefined) {
2953: PetscCall(PetscFERefine(fe, &feRef[f]));
2954: } else {
2955: PetscCall(PetscObjectReference((PetscObject)fe));
2956: feRef[f] = fe;
2957: }
2958: PetscCall(PetscFEGetDimension(feRef[f], &rNb));
2959: PetscCall(PetscFEGetNumComponents(fe, &Nc));
2960: } else if (id == PETSCFV_CLASSID) {
2961: PetscFV fv = (PetscFV)obj;
2962: PetscDualSpace Q;
2964: if (isRefined) {
2965: PetscCall(PetscFVRefine(fv, &fvRef[f]));
2966: } else {
2967: PetscCall(PetscObjectReference((PetscObject)fv));
2968: fvRef[f] = fv;
2969: }
2970: PetscCall(PetscFVGetDualSpace(fvRef[f], &Q));
2971: PetscCall(PetscDualSpaceGetDimension(Q, &rNb));
2972: PetscCall(PetscFVGetDualSpace(fv, &Q));
2973: PetscCall(PetscFVGetNumComponents(fv, &Nc));
2974: }
2975: PetscCall(DMGetField(dmc, f, NULL, &objc));
2976: PetscCall(PetscObjectGetClassId(objc, &idc));
2977: if (idc == PETSCFE_CLASSID) {
2978: PetscFE fe = (PetscFE)objc;
2980: PetscCall(PetscFEGetDimension(fe, &cNb));
2981: } else if (id == PETSCFV_CLASSID) {
2982: PetscFV fv = (PetscFV)obj;
2983: PetscDualSpace Q;
2985: PetscCall(PetscFVGetDualSpace(fv, &Q));
2986: PetscCall(PetscDualSpaceGetDimension(Q, &cNb));
2987: }
2988: rTotDim += rNb;
2989: cTotDim += cNb;
2990: }
2991: PetscCall(PetscMalloc1(rTotDim * cTotDim, &elemMat));
2992: PetscCall(PetscArrayzero(elemMat, rTotDim * cTotDim));
2993: for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
2994: PetscDualSpace Qref;
2995: PetscQuadrature f;
2996: const PetscReal *qpoints, *qweights;
2997: PetscReal *points;
2998: PetscInt npoints = 0, Nc, Np, fpdim, i, k, p, d;
3000: /* Compose points from all dual basis functionals */
3001: if (feRef[fieldI]) {
3002: PetscCall(PetscFEGetDualSpace(feRef[fieldI], &Qref));
3003: PetscCall(PetscFEGetNumComponents(feRef[fieldI], &Nc));
3004: } else {
3005: PetscCall(PetscFVGetDualSpace(fvRef[fieldI], &Qref));
3006: PetscCall(PetscFVGetNumComponents(fvRef[fieldI], &Nc));
3007: }
3008: PetscCall(PetscDualSpaceGetDimension(Qref, &fpdim));
3009: for (i = 0; i < fpdim; ++i) {
3010: PetscCall(PetscDualSpaceGetFunctional(Qref, i, &f));
3011: PetscCall(PetscQuadratureGetData(f, NULL, NULL, &Np, NULL, NULL));
3012: npoints += Np;
3013: }
3014: PetscCall(PetscMalloc1(npoints * dim, &points));
3015: for (i = 0, k = 0; i < fpdim; ++i) {
3016: PetscCall(PetscDualSpaceGetFunctional(Qref, i, &f));
3017: PetscCall(PetscQuadratureGetData(f, NULL, NULL, &Np, &qpoints, NULL));
3018: for (p = 0; p < Np; ++p, ++k)
3019: for (d = 0; d < dim; ++d) points[k * dim + d] = qpoints[p * dim + d];
3020: }
3022: for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
3023: PetscObject obj;
3024: PetscClassId id;
3025: PetscInt NcJ = 0, cpdim = 0, j, qNc;
3027: PetscCall(DMGetField(dmc, fieldJ, NULL, &obj));
3028: PetscCall(PetscObjectGetClassId(obj, &id));
3029: if (id == PETSCFE_CLASSID) {
3030: PetscFE fe = (PetscFE)obj;
3031: PetscTabulation T = NULL;
3033: /* Evaluate basis at points */
3034: PetscCall(PetscFEGetNumComponents(fe, &NcJ));
3035: PetscCall(PetscFEGetDimension(fe, &cpdim));
3036: /* For now, fields only interpolate themselves */
3037: if (fieldI == fieldJ) {
3038: PetscCheck(Nc == NcJ, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %" PetscInt_FMT " does not match coarse field %" PetscInt_FMT, Nc, NcJ);
3039: PetscCall(PetscFECreateTabulation(fe, 1, npoints, points, 0, &T));
3040: for (i = 0, k = 0; i < fpdim; ++i) {
3041: PetscCall(PetscDualSpaceGetFunctional(Qref, i, &f));
3042: PetscCall(PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights));
3043: PetscCheck(qNc == NcJ, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %" PetscInt_FMT " does not match coarse field %" PetscInt_FMT, qNc, NcJ);
3044: for (p = 0; p < Np; ++p, ++k) {
3045: for (j = 0; j < cpdim; ++j) {
3046: /*
3047: cTotDim: Total columns in element interpolation matrix, sum of number of dual basis functionals in each field
3048: offsetI, offsetJ: Offsets into the larger element interpolation matrix for different fields
3049: fpdim, i, cpdim, j: Dofs for fine and coarse grids, correspond to dual space basis functionals
3050: qNC, Nc, Ncj, c: Number of components in this field
3051: Np, p: Number of quad points in the fine grid functional i
3052: k: i*Np + p, overall point number for the interpolation
3053: */
3054: 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];
3055: }
3056: }
3057: }
3058: PetscCall(PetscTabulationDestroy(&T));
3059: }
3060: } else if (id == PETSCFV_CLASSID) {
3061: PetscFV fv = (PetscFV)obj;
3063: /* Evaluate constant function at points */
3064: PetscCall(PetscFVGetNumComponents(fv, &NcJ));
3065: cpdim = 1;
3066: /* For now, fields only interpolate themselves */
3067: if (fieldI == fieldJ) {
3068: PetscCheck(Nc == NcJ, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %" PetscInt_FMT " does not match coarse field %" PetscInt_FMT, Nc, NcJ);
3069: for (i = 0, k = 0; i < fpdim; ++i) {
3070: PetscCall(PetscDualSpaceGetFunctional(Qref, i, &f));
3071: PetscCall(PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights));
3072: PetscCheck(qNc == NcJ, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %" PetscInt_FMT " does not match coarse field %" PetscInt_FMT, qNc, NcJ);
3073: for (p = 0; p < Np; ++p, ++k) {
3074: for (j = 0; j < cpdim; ++j) {
3075: for (c = 0; c < Nc; ++c) elemMat[(offsetI + i) * cTotDim + offsetJ + j] += 1.0 * qweights[p * qNc + c];
3076: }
3077: }
3078: }
3079: }
3080: }
3081: offsetJ += cpdim;
3082: }
3083: offsetI += fpdim;
3084: PetscCall(PetscFree(points));
3085: }
3086: if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat));
3087: /* Preallocate matrix */
3088: {
3089: Mat preallocator;
3090: PetscScalar *vals;
3091: PetscInt *cellCIndices, *cellFIndices;
3092: PetscInt locRows, locCols, cell;
3094: PetscCall(MatGetLocalSize(In, &locRows, &locCols));
3095: PetscCall(MatCreate(PetscObjectComm((PetscObject)In), &preallocator));
3096: PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
3097: PetscCall(MatSetSizes(preallocator, locRows, locCols, PETSC_DETERMINE, PETSC_DETERMINE));
3098: PetscCall(MatSetUp(preallocator));
3099: PetscCall(PetscCalloc3(rTotDim * cTotDim, &vals, cTotDim, &cellCIndices, rTotDim, &cellFIndices));
3100: if (locRows || locCols) {
3101: for (cell = cStart; cell < cEnd; ++cell) {
3102: if (isRefined) {
3103: PetscCall(DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices));
3104: PetscCall(MatSetValues(preallocator, rTotDim, cellFIndices, cTotDim, cellCIndices, vals, INSERT_VALUES));
3105: } else {
3106: PetscCall(DMPlexMatSetClosureGeneral(dmf, fsection, fglobalSection, PETSC_FALSE, dmc, csection, cglobalSection, PETSC_FALSE, preallocator, cell, vals, INSERT_VALUES));
3107: }
3108: }
3109: }
3110: PetscCall(PetscFree3(vals, cellCIndices, cellFIndices));
3111: PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
3112: PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
3113: PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, In));
3114: PetscCall(MatDestroy(&preallocator));
3115: }
3116: /* Fill matrix */
3117: PetscCall(MatZeroEntries(In));
3118: for (c = cStart; c < cEnd; ++c) {
3119: if (isRefined) {
3120: PetscCall(DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES));
3121: } else {
3122: PetscCall(DMPlexMatSetClosureGeneral(dmf, fsection, fglobalSection, PETSC_FALSE, dmc, csection, cglobalSection, PETSC_FALSE, In, c, elemMat, INSERT_VALUES));
3123: }
3124: }
3125: for (f = 0; f < Nf; ++f) PetscCall(PetscFEDestroy(&feRef[f]));
3126: PetscCall(PetscFree2(feRef, fvRef));
3127: PetscCall(PetscFree(elemMat));
3128: PetscCall(MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY));
3129: PetscCall(MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY));
3130: if (mesh->printFEM > 1) {
3131: PetscCall(PetscPrintf(PetscObjectComm((PetscObject)In), "%s:\n", name));
3132: PetscCall(MatFilter(In, 1.0e-10, PETSC_FALSE, PETSC_FALSE));
3133: PetscCall(MatView(In, NULL));
3134: }
3135: PetscCall(PetscLogEventEnd(DMPLEX_InterpolatorFEM, dmc, dmf, 0, 0));
3136: PetscFunctionReturn(PETSC_SUCCESS);
3137: }
3139: PetscErrorCode DMPlexComputeMassMatrixNested(DM dmc, DM dmf, Mat mass, void *user)
3140: {
3141: SETERRQ(PetscObjectComm((PetscObject)dmc), PETSC_ERR_SUP, "Laziness");
3142: }
3144: /*@
3145: DMPlexComputeInterpolatorGeneral - Form the local portion of the interpolation matrix from the coarse `DM` to a non-nested fine `DM`.
3147: Input Parameters:
3148: + dmf - The fine mesh
3149: . dmc - The coarse mesh
3150: - user - The user context
3152: Output Parameter:
3153: . In - The interpolation matrix
3155: Level: developer
3157: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexComputeInterpolatorNested()`
3158: @*/
3159: PetscErrorCode DMPlexComputeInterpolatorGeneral(DM dmc, DM dmf, Mat In, void *user)
3160: {
3161: DM_Plex *mesh = (DM_Plex *)dmf->data;
3162: const char *name = "Interpolator";
3163: PetscDS prob;
3164: Mat interp;
3165: PetscSection fsection, globalFSection;
3166: PetscSection csection, globalCSection;
3167: PetscInt locRows, locCols;
3168: PetscReal *x, *v0, *J, *invJ, detJ;
3169: PetscReal *v0c, *Jc, *invJc, detJc;
3170: PetscScalar *elemMat;
3171: PetscInt dim, Nf, field, totDim, cStart, cEnd, cell, ccell, s;
3173: PetscFunctionBegin;
3174: PetscCall(PetscLogEventBegin(DMPLEX_InterpolatorFEM, dmc, dmf, 0, 0));
3175: PetscCall(DMGetCoordinateDim(dmc, &dim));
3176: PetscCall(DMGetDS(dmc, &prob));
3177: PetscCall(PetscDSGetWorkspace(prob, &x, NULL, NULL, NULL, NULL));
3178: PetscCall(PetscDSGetNumFields(prob, &Nf));
3179: PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ));
3180: PetscCall(PetscMalloc3(dim, &v0c, dim * dim, &Jc, dim * dim, &invJc));
3181: PetscCall(DMGetLocalSection(dmf, &fsection));
3182: PetscCall(DMGetGlobalSection(dmf, &globalFSection));
3183: PetscCall(DMGetLocalSection(dmc, &csection));
3184: PetscCall(DMGetGlobalSection(dmc, &globalCSection));
3185: PetscCall(DMPlexGetSimplexOrBoxCells(dmf, 0, &cStart, &cEnd));
3186: PetscCall(PetscDSGetTotalDimension(prob, &totDim));
3187: PetscCall(PetscMalloc1(totDim, &elemMat));
3189: PetscCall(MatGetLocalSize(In, &locRows, &locCols));
3190: PetscCall(MatCreate(PetscObjectComm((PetscObject)In), &interp));
3191: PetscCall(MatSetType(interp, MATPREALLOCATOR));
3192: PetscCall(MatSetSizes(interp, locRows, locCols, PETSC_DETERMINE, PETSC_DETERMINE));
3193: PetscCall(MatSetUp(interp));
3194: for (s = 0; s < 2; ++s) {
3195: for (field = 0; field < Nf; ++field) {
3196: PetscObject obj;
3197: PetscClassId id;
3198: PetscDualSpace Q = NULL;
3199: PetscTabulation T = NULL;
3200: PetscQuadrature f;
3201: const PetscReal *qpoints, *qweights;
3202: PetscInt Nc, qNc, Np, fpdim, off, i, d;
3204: PetscCall(PetscDSGetFieldOffset(prob, field, &off));
3205: PetscCall(PetscDSGetDiscretization(prob, field, &obj));
3206: PetscCall(PetscObjectGetClassId(obj, &id));
3207: if (id == PETSCFE_CLASSID) {
3208: PetscFE fe = (PetscFE)obj;
3210: PetscCall(PetscFEGetDualSpace(fe, &Q));
3211: PetscCall(PetscFEGetNumComponents(fe, &Nc));
3212: if (s) PetscCall(PetscFECreateTabulation(fe, 1, 1, x, 0, &T));
3213: } else if (id == PETSCFV_CLASSID) {
3214: PetscFV fv = (PetscFV)obj;
3216: PetscCall(PetscFVGetDualSpace(fv, &Q));
3217: Nc = 1;
3218: } else SETERRQ(PetscObjectComm((PetscObject)dmc), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
3219: PetscCall(PetscDualSpaceGetDimension(Q, &fpdim));
3220: /* For each fine grid cell */
3221: for (cell = cStart; cell < cEnd; ++cell) {
3222: PetscInt *findices, *cindices;
3223: PetscInt numFIndices, numCIndices;
3225: PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3226: PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
3227: PetscCheck(numFIndices == totDim, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %" PetscInt_FMT " != %" PetscInt_FMT " dual basis vecs", numFIndices, totDim);
3228: for (i = 0; i < fpdim; ++i) {
3229: Vec pointVec;
3230: PetscScalar *pV;
3231: PetscSF coarseCellSF = NULL;
3232: const PetscSFNode *coarseCells;
3233: PetscInt numCoarseCells, cpdim, row = findices[i + off], q, c, j;
3235: /* Get points from the dual basis functional quadrature */
3236: PetscCall(PetscDualSpaceGetFunctional(Q, i, &f));
3237: PetscCall(PetscQuadratureGetData(f, NULL, &qNc, &Np, &qpoints, &qweights));
3238: PetscCheck(qNc == Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %" PetscInt_FMT " does not match coarse field %" PetscInt_FMT, qNc, Nc);
3239: PetscCall(VecCreateSeq(PETSC_COMM_SELF, Np * dim, &pointVec));
3240: PetscCall(VecSetBlockSize(pointVec, dim));
3241: PetscCall(VecGetArray(pointVec, &pV));
3242: for (q = 0; q < Np; ++q) {
3243: const PetscReal xi0[3] = {-1., -1., -1.};
3245: /* Transform point to real space */
3246: CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q * dim], x);
3247: for (d = 0; d < dim; ++d) pV[q * dim + d] = x[d];
3248: }
3249: PetscCall(VecRestoreArray(pointVec, &pV));
3250: /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
3251: /* OPT: Read this out from preallocation information */
3252: PetscCall(DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF));
3253: /* Update preallocation info */
3254: PetscCall(PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells));
3255: PetscCheck(numCoarseCells == Np, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not all closure points located");
3256: PetscCall(VecGetArray(pointVec, &pV));
3257: for (ccell = 0; ccell < numCoarseCells; ++ccell) {
3258: PetscReal pVReal[3];
3259: const PetscReal xi0[3] = {-1., -1., -1.};
3261: PetscCall(DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL));
3262: if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetDimension((PetscFE)obj, &cpdim));
3263: else cpdim = 1;
3265: if (s) {
3266: /* Transform points from real space to coarse reference space */
3267: PetscCall(DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc));
3268: for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell * dim + d]);
3269: CoordinatesRealToRef(dim, dim, xi0, v0c, invJc, pVReal, x);
3271: if (id == PETSCFE_CLASSID) {
3272: /* Evaluate coarse basis on contained point */
3273: PetscCall(PetscFEComputeTabulation((PetscFE)obj, 1, x, 0, T));
3274: PetscCall(PetscArrayzero(elemMat, cpdim));
3275: /* Get elemMat entries by multiplying by weight */
3276: for (j = 0; j < cpdim; ++j) {
3277: for (c = 0; c < Nc; ++c) elemMat[j] += T->T[0][j * Nc + c] * qweights[ccell * qNc + c];
3278: }
3279: } else {
3280: for (j = 0; j < cpdim; ++j) {
3281: for (c = 0; c < Nc; ++c) elemMat[j] += 1.0 * qweights[ccell * qNc + c];
3282: }
3283: }
3284: if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
3285: }
3286: /* Update interpolator */
3287: PetscCheck(numCIndices == totDim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %" PetscInt_FMT " != %" PetscInt_FMT, numCIndices, totDim);
3288: PetscCall(MatSetValues(interp, 1, &row, cpdim, &cindices[off], elemMat, INSERT_VALUES));
3289: PetscCall(DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL));
3290: }
3291: PetscCall(VecRestoreArray(pointVec, &pV));
3292: PetscCall(PetscSFDestroy(&coarseCellSF));
3293: PetscCall(VecDestroy(&pointVec));
3294: }
3295: PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3296: }
3297: if (s && id == PETSCFE_CLASSID) PetscCall(PetscTabulationDestroy(&T));
3298: }
3299: if (!s) {
3300: PetscCall(MatAssemblyBegin(interp, MAT_FINAL_ASSEMBLY));
3301: PetscCall(MatAssemblyEnd(interp, MAT_FINAL_ASSEMBLY));
3302: PetscCall(MatPreallocatorPreallocate(interp, PETSC_TRUE, In));
3303: PetscCall(MatDestroy(&interp));
3304: interp = In;
3305: }
3306: }
3307: PetscCall(PetscFree3(v0, J, invJ));
3308: PetscCall(PetscFree3(v0c, Jc, invJc));
3309: PetscCall(PetscFree(elemMat));
3310: PetscCall(MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY));
3311: PetscCall(MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY));
3312: PetscCall(PetscLogEventEnd(DMPLEX_InterpolatorFEM, dmc, dmf, 0, 0));
3313: PetscFunctionReturn(PETSC_SUCCESS);
3314: }
3316: /*@
3317: DMPlexComputeMassMatrixGeneral - Form the local portion of the mass matrix from the coarse `DM` to a non-nested fine `DM`.
3319: Input Parameters:
3320: + dmf - The fine mesh
3321: . dmc - The coarse mesh
3322: - user - The user context
3324: Output Parameter:
3325: . mass - The mass matrix
3327: Level: developer
3329: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexComputeMassMatrixNested()`, `DMPlexComputeInterpolatorNested()`, `DMPlexComputeInterpolatorGeneral()`
3330: @*/
3331: PetscErrorCode DMPlexComputeMassMatrixGeneral(DM dmc, DM dmf, Mat mass, void *user)
3332: {
3333: DM_Plex *mesh = (DM_Plex *)dmf->data;
3334: const char *name = "Mass Matrix";
3335: PetscDS prob;
3336: PetscSection fsection, csection, globalFSection, globalCSection;
3337: PetscHSetIJ ht;
3338: PetscLayout rLayout;
3339: PetscInt *dnz, *onz;
3340: PetscInt locRows, rStart, rEnd;
3341: PetscReal *x, *v0, *J, *invJ, detJ;
3342: PetscReal *v0c, *Jc, *invJc, detJc;
3343: PetscScalar *elemMat;
3344: PetscInt dim, Nf, field, totDim, cStart, cEnd, cell, ccell;
3346: PetscFunctionBegin;
3347: PetscCall(DMGetCoordinateDim(dmc, &dim));
3348: PetscCall(DMGetDS(dmc, &prob));
3349: PetscCall(PetscDSGetWorkspace(prob, &x, NULL, NULL, NULL, NULL));
3350: PetscCall(PetscDSGetNumFields(prob, &Nf));
3351: PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ));
3352: PetscCall(PetscMalloc3(dim, &v0c, dim * dim, &Jc, dim * dim, &invJc));
3353: PetscCall(DMGetLocalSection(dmf, &fsection));
3354: PetscCall(DMGetGlobalSection(dmf, &globalFSection));
3355: PetscCall(DMGetLocalSection(dmc, &csection));
3356: PetscCall(DMGetGlobalSection(dmc, &globalCSection));
3357: PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd));
3358: PetscCall(PetscDSGetTotalDimension(prob, &totDim));
3359: PetscCall(PetscMalloc1(totDim, &elemMat));
3361: PetscCall(MatGetLocalSize(mass, &locRows, NULL));
3362: PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)mass), &rLayout));
3363: PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
3364: PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
3365: PetscCall(PetscLayoutSetUp(rLayout));
3366: PetscCall(PetscLayoutGetRange(rLayout, &rStart, &rEnd));
3367: PetscCall(PetscLayoutDestroy(&rLayout));
3368: PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz));
3369: PetscCall(PetscHSetIJCreate(&ht));
3370: for (field = 0; field < Nf; ++field) {
3371: PetscObject obj;
3372: PetscClassId id;
3373: PetscQuadrature quad;
3374: const PetscReal *qpoints;
3375: PetscInt Nq, Nc, i, d;
3377: PetscCall(PetscDSGetDiscretization(prob, field, &obj));
3378: PetscCall(PetscObjectGetClassId(obj, &id));
3379: if (id == PETSCFE_CLASSID) PetscCall(PetscFEGetQuadrature((PetscFE)obj, &quad));
3380: else PetscCall(PetscFVGetQuadrature((PetscFV)obj, &quad));
3381: PetscCall(PetscQuadratureGetData(quad, NULL, &Nc, &Nq, &qpoints, NULL));
3382: /* For each fine grid cell */
3383: for (cell = cStart; cell < cEnd; ++cell) {
3384: Vec pointVec;
3385: PetscScalar *pV;
3386: PetscSF coarseCellSF = NULL;
3387: const PetscSFNode *coarseCells;
3388: PetscInt numCoarseCells, q, c;
3389: PetscInt *findices, *cindices;
3390: PetscInt numFIndices, numCIndices;
3392: PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3393: PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
3394: /* Get points from the quadrature */
3395: PetscCall(VecCreateSeq(PETSC_COMM_SELF, Nq * dim, &pointVec));
3396: PetscCall(VecSetBlockSize(pointVec, dim));
3397: PetscCall(VecGetArray(pointVec, &pV));
3398: for (q = 0; q < Nq; ++q) {
3399: const PetscReal xi0[3] = {-1., -1., -1.};
3401: /* Transform point to real space */
3402: CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q * dim], x);
3403: for (d = 0; d < dim; ++d) pV[q * dim + d] = x[d];
3404: }
3405: PetscCall(VecRestoreArray(pointVec, &pV));
3406: /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
3407: PetscCall(DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF));
3408: PetscCall(PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view"));
3409: /* Update preallocation info */
3410: PetscCall(PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells));
3411: PetscCheck(numCoarseCells == Nq, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not all closure points located");
3412: {
3413: PetscHashIJKey key;
3414: PetscBool missing;
3416: for (i = 0; i < numFIndices; ++i) {
3417: key.i = findices[i];
3418: if (key.i >= 0) {
3419: /* Get indices for coarse elements */
3420: for (ccell = 0; ccell < numCoarseCells; ++ccell) {
3421: PetscCall(DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL));
3422: for (c = 0; c < numCIndices; ++c) {
3423: key.j = cindices[c];
3424: if (key.j < 0) continue;
3425: PetscCall(PetscHSetIJQueryAdd(ht, key, &missing));
3426: if (missing) {
3427: if ((key.j >= rStart) && (key.j < rEnd)) ++dnz[key.i - rStart];
3428: else ++onz[key.i - rStart];
3429: }
3430: }
3431: PetscCall(DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL));
3432: }
3433: }
3434: }
3435: }
3436: PetscCall(PetscSFDestroy(&coarseCellSF));
3437: PetscCall(VecDestroy(&pointVec));
3438: PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3439: }
3440: }
3441: PetscCall(PetscHSetIJDestroy(&ht));
3442: PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL));
3443: PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
3444: PetscCall(PetscFree2(dnz, onz));
3445: for (field = 0; field < Nf; ++field) {
3446: PetscObject obj;
3447: PetscClassId id;
3448: PetscTabulation T, Tfine;
3449: PetscQuadrature quad;
3450: const PetscReal *qpoints, *qweights;
3451: PetscInt Nq, Nc, i, d;
3453: PetscCall(PetscDSGetDiscretization(prob, field, &obj));
3454: PetscCall(PetscObjectGetClassId(obj, &id));
3455: if (id == PETSCFE_CLASSID) {
3456: PetscCall(PetscFEGetQuadrature((PetscFE)obj, &quad));
3457: PetscCall(PetscFEGetCellTabulation((PetscFE)obj, 1, &Tfine));
3458: PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, 1, x, 0, &T));
3459: } else {
3460: PetscCall(PetscFVGetQuadrature((PetscFV)obj, &quad));
3461: }
3462: PetscCall(PetscQuadratureGetData(quad, NULL, &Nc, &Nq, &qpoints, &qweights));
3463: /* For each fine grid cell */
3464: for (cell = cStart; cell < cEnd; ++cell) {
3465: Vec pointVec;
3466: PetscScalar *pV;
3467: PetscSF coarseCellSF = NULL;
3468: const PetscSFNode *coarseCells;
3469: PetscInt numCoarseCells, cpdim, q, c, j;
3470: PetscInt *findices, *cindices;
3471: PetscInt numFIndices, numCIndices;
3473: PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3474: PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ));
3475: /* Get points from the quadrature */
3476: PetscCall(VecCreateSeq(PETSC_COMM_SELF, Nq * dim, &pointVec));
3477: PetscCall(VecSetBlockSize(pointVec, dim));
3478: PetscCall(VecGetArray(pointVec, &pV));
3479: for (q = 0; q < Nq; ++q) {
3480: const PetscReal xi0[3] = {-1., -1., -1.};
3482: /* Transform point to real space */
3483: CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q * dim], x);
3484: for (d = 0; d < dim; ++d) pV[q * dim + d] = x[d];
3485: }
3486: PetscCall(VecRestoreArray(pointVec, &pV));
3487: /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
3488: PetscCall(DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF));
3489: /* Update matrix */
3490: PetscCall(PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells));
3491: PetscCheck(numCoarseCells == Nq, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not all closure points located");
3492: PetscCall(VecGetArray(pointVec, &pV));
3493: for (ccell = 0; ccell < numCoarseCells; ++ccell) {
3494: PetscReal pVReal[3];
3495: const PetscReal xi0[3] = {-1., -1., -1.};
3497: PetscCall(DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL));
3498: /* Transform points from real space to coarse reference space */
3499: PetscCall(DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc));
3500: for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell * dim + d]);
3501: CoordinatesRealToRef(dim, dim, xi0, v0c, invJc, pVReal, x);
3503: if (id == PETSCFE_CLASSID) {
3504: PetscFE fe = (PetscFE)obj;
3506: /* Evaluate coarse basis on contained point */
3507: PetscCall(PetscFEGetDimension(fe, &cpdim));
3508: PetscCall(PetscFEComputeTabulation(fe, 1, x, 0, T));
3509: /* Get elemMat entries by multiplying by weight */
3510: for (i = 0; i < numFIndices; ++i) {
3511: PetscCall(PetscArrayzero(elemMat, cpdim));
3512: for (j = 0; j < cpdim; ++j) {
3513: 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;
3514: }
3515: /* Update interpolator */
3516: if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
3517: PetscCheck(numCIndices == cpdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %" PetscInt_FMT " != %" PetscInt_FMT, numCIndices, cpdim);
3518: PetscCall(MatSetValues(mass, 1, &findices[i], numCIndices, cindices, elemMat, ADD_VALUES));
3519: }
3520: } else {
3521: cpdim = 1;
3522: for (i = 0; i < numFIndices; ++i) {
3523: PetscCall(PetscArrayzero(elemMat, cpdim));
3524: for (j = 0; j < cpdim; ++j) {
3525: for (c = 0; c < Nc; ++c) elemMat[j] += 1.0 * 1.0 * qweights[ccell * Nc + c] * detJ;
3526: }
3527: /* Update interpolator */
3528: if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat));
3529: PetscCall(PetscPrintf(PETSC_COMM_SELF, "Nq: %" PetscInt_FMT " %" PetscInt_FMT " Nf: %" PetscInt_FMT " %" PetscInt_FMT " Nc: %" PetscInt_FMT " %" PetscInt_FMT "\n", ccell, Nq, i, numFIndices, j, numCIndices));
3530: PetscCheck(numCIndices == cpdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %" PetscInt_FMT " != %" PetscInt_FMT, numCIndices, cpdim);
3531: PetscCall(MatSetValues(mass, 1, &findices[i], numCIndices, cindices, elemMat, ADD_VALUES));
3532: }
3533: }
3534: PetscCall(DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, PETSC_FALSE, &numCIndices, &cindices, NULL, NULL));
3535: }
3536: PetscCall(VecRestoreArray(pointVec, &pV));
3537: PetscCall(PetscSFDestroy(&coarseCellSF));
3538: PetscCall(VecDestroy(&pointVec));
3539: PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL));
3540: }
3541: if (id == PETSCFE_CLASSID) PetscCall(PetscTabulationDestroy(&T));
3542: }
3543: PetscCall(PetscFree3(v0, J, invJ));
3544: PetscCall(PetscFree3(v0c, Jc, invJc));
3545: PetscCall(PetscFree(elemMat));
3546: PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY));
3547: PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY));
3548: PetscFunctionReturn(PETSC_SUCCESS);
3549: }
3551: /*@
3552: DMPlexComputeInjectorFEM - Compute a mapping from coarse unknowns to fine unknowns
3554: Input Parameters:
3555: + dmc - The coarse mesh
3556: . dmf - The fine mesh
3557: - user - The user context
3559: Output Parameter:
3560: . sc - The mapping
3562: Level: developer
3564: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexComputeInterpolatorNested()`
3565: @*/
3566: PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
3567: {
3568: PetscDS prob;
3569: PetscFE *feRef;
3570: PetscFV *fvRef;
3571: Vec fv, cv;
3572: IS fis, cis;
3573: PetscSection fsection, fglobalSection, csection, cglobalSection;
3574: PetscInt *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
3575: PetscInt cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, c, dim, d, startC, endC, offsetC, offsetF, m;
3576: PetscBool *needAvg;
3578: PetscFunctionBegin;
3579: PetscCall(PetscLogEventBegin(DMPLEX_InjectorFEM, dmc, dmf, 0, 0));
3580: PetscCall(DMGetDimension(dmf, &dim));
3581: PetscCall(DMGetLocalSection(dmf, &fsection));
3582: PetscCall(DMGetGlobalSection(dmf, &fglobalSection));
3583: PetscCall(DMGetLocalSection(dmc, &csection));
3584: PetscCall(DMGetGlobalSection(dmc, &cglobalSection));
3585: PetscCall(PetscSectionGetNumFields(fsection, &Nf));
3586: PetscCall(DMPlexGetSimplexOrBoxCells(dmc, 0, &cStart, &cEnd));
3587: PetscCall(DMGetDS(dmc, &prob));
3588: PetscCall(PetscCalloc3(Nf, &feRef, Nf, &fvRef, Nf, &needAvg));
3589: for (f = 0; f < Nf; ++f) {
3590: PetscObject obj;
3591: PetscClassId id;
3592: PetscInt fNb = 0, Nc = 0;
3594: PetscCall(PetscDSGetDiscretization(prob, f, &obj));
3595: PetscCall(PetscObjectGetClassId(obj, &id));
3596: if (id == PETSCFE_CLASSID) {
3597: PetscFE fe = (PetscFE)obj;
3598: PetscSpace sp;
3599: PetscInt maxDegree;
3601: PetscCall(PetscFERefine(fe, &feRef[f]));
3602: PetscCall(PetscFEGetDimension(feRef[f], &fNb));
3603: PetscCall(PetscFEGetNumComponents(fe, &Nc));
3604: PetscCall(PetscFEGetBasisSpace(fe, &sp));
3605: PetscCall(PetscSpaceGetDegree(sp, NULL, &maxDegree));
3606: if (!maxDegree) needAvg[f] = PETSC_TRUE;
3607: } else if (id == PETSCFV_CLASSID) {
3608: PetscFV fv = (PetscFV)obj;
3609: PetscDualSpace Q;
3611: PetscCall(PetscFVRefine(fv, &fvRef[f]));
3612: PetscCall(PetscFVGetDualSpace(fvRef[f], &Q));
3613: PetscCall(PetscDualSpaceGetDimension(Q, &fNb));
3614: PetscCall(PetscFVGetNumComponents(fv, &Nc));
3615: needAvg[f] = PETSC_TRUE;
3616: }
3617: fTotDim += fNb;
3618: }
3619: PetscCall(PetscDSGetTotalDimension(prob, &cTotDim));
3620: PetscCall(PetscMalloc1(cTotDim, &cmap));
3621: for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
3622: PetscFE feC;
3623: PetscFV fvC;
3624: PetscDualSpace QF, QC;
3625: PetscInt order = -1, NcF, NcC, fpdim, cpdim;
3627: if (feRef[field]) {
3628: PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&feC));
3629: PetscCall(PetscFEGetNumComponents(feC, &NcC));
3630: PetscCall(PetscFEGetNumComponents(feRef[field], &NcF));
3631: PetscCall(PetscFEGetDualSpace(feRef[field], &QF));
3632: PetscCall(PetscDualSpaceGetOrder(QF, &order));
3633: PetscCall(PetscDualSpaceGetDimension(QF, &fpdim));
3634: PetscCall(PetscFEGetDualSpace(feC, &QC));
3635: PetscCall(PetscDualSpaceGetDimension(QC, &cpdim));
3636: } else {
3637: PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fvC));
3638: PetscCall(PetscFVGetNumComponents(fvC, &NcC));
3639: PetscCall(PetscFVGetNumComponents(fvRef[field], &NcF));
3640: PetscCall(PetscFVGetDualSpace(fvRef[field], &QF));
3641: PetscCall(PetscDualSpaceGetDimension(QF, &fpdim));
3642: PetscCall(PetscFVGetDualSpace(fvC, &QC));
3643: PetscCall(PetscDualSpaceGetDimension(QC, &cpdim));
3644: }
3645: PetscCheck(NcF == NcC, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %" PetscInt_FMT " does not match coarse field %" PetscInt_FMT, NcF, NcC);
3646: for (c = 0; c < cpdim; ++c) {
3647: PetscQuadrature cfunc;
3648: const PetscReal *cqpoints, *cqweights;
3649: PetscInt NqcC, NpC;
3650: PetscBool found = PETSC_FALSE;
3652: PetscCall(PetscDualSpaceGetFunctional(QC, c, &cfunc));
3653: PetscCall(PetscQuadratureGetData(cfunc, NULL, &NqcC, &NpC, &cqpoints, &cqweights));
3654: PetscCheck(NqcC == NcC, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of quadrature components %" PetscInt_FMT " must match number of field components %" PetscInt_FMT, NqcC, NcC);
3655: PetscCheck(NpC == 1 || !feRef[field], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
3656: for (f = 0; f < fpdim; ++f) {
3657: PetscQuadrature ffunc;
3658: const PetscReal *fqpoints, *fqweights;
3659: PetscReal sum = 0.0;
3660: PetscInt NqcF, NpF;
3662: PetscCall(PetscDualSpaceGetFunctional(QF, f, &ffunc));
3663: PetscCall(PetscQuadratureGetData(ffunc, NULL, &NqcF, &NpF, &fqpoints, &fqweights));
3664: PetscCheck(NqcF == NcF, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of quadrature components %" PetscInt_FMT " must match number of field components %" PetscInt_FMT, NqcF, NcF);
3665: if (NpC != NpF) continue;
3666: for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
3667: if (sum > 1.0e-9) continue;
3668: for (d = 0; d < NcC; ++d) sum += PetscAbsReal(cqweights[d] * fqweights[d]);
3669: if (sum < 1.0e-9) continue;
3670: cmap[offsetC + c] = offsetF + f;
3671: found = PETSC_TRUE;
3672: break;
3673: }
3674: if (!found) {
3675: /* TODO We really want the average here, but some asshole put VecScatter in the interface */
3676: PetscCheck(fvRef[field] || (feRef[field] && order == 0), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate matching functional for injection");
3677: cmap[offsetC + c] = offsetF + 0;
3678: }
3679: }
3680: offsetC += cpdim;
3681: offsetF += fpdim;
3682: }
3683: for (f = 0; f < Nf; ++f) {
3684: PetscCall(PetscFEDestroy(&feRef[f]));
3685: PetscCall(PetscFVDestroy(&fvRef[f]));
3686: }
3687: PetscCall(PetscFree3(feRef, fvRef, needAvg));
3689: PetscCall(DMGetGlobalVector(dmf, &fv));
3690: PetscCall(DMGetGlobalVector(dmc, &cv));
3691: PetscCall(VecGetOwnershipRange(cv, &startC, &endC));
3692: PetscCall(PetscSectionGetConstrainedStorageSize(cglobalSection, &m));
3693: PetscCall(PetscMalloc2(cTotDim, &cellCIndices, fTotDim, &cellFIndices));
3694: PetscCall(PetscMalloc1(m, &cindices));
3695: PetscCall(PetscMalloc1(m, &findices));
3696: for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
3697: for (c = cStart; c < cEnd; ++c) {
3698: PetscCall(DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices));
3699: for (d = 0; d < cTotDim; ++d) {
3700: if ((cellCIndices[d] < startC) || (cellCIndices[d] >= endC)) continue;
3701: PetscCheck(!(findices[cellCIndices[d] - startC] >= 0) || !(findices[cellCIndices[d] - startC] != cellFIndices[cmap[d]]), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %" PetscInt_FMT " Coarse dof %" PetscInt_FMT " maps to both %" PetscInt_FMT " and %" PetscInt_FMT, c, cindices[cellCIndices[d] - startC], findices[cellCIndices[d] - startC], cellFIndices[cmap[d]]);
3702: cindices[cellCIndices[d] - startC] = cellCIndices[d];
3703: findices[cellCIndices[d] - startC] = cellFIndices[cmap[d]];
3704: }
3705: }
3706: PetscCall(PetscFree(cmap));
3707: PetscCall(PetscFree2(cellCIndices, cellFIndices));
3709: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis));
3710: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis));
3711: PetscCall(VecScatterCreate(cv, cis, fv, fis, sc));
3712: PetscCall(ISDestroy(&cis));
3713: PetscCall(ISDestroy(&fis));
3714: PetscCall(DMRestoreGlobalVector(dmf, &fv));
3715: PetscCall(DMRestoreGlobalVector(dmc, &cv));
3716: PetscCall(PetscLogEventEnd(DMPLEX_InjectorFEM, dmc, dmf, 0, 0));
3717: PetscFunctionReturn(PETSC_SUCCESS);
3718: }
3720: /*@C
3721: DMPlexGetCellFields - Retrieve the field values values for a chunk of cells
3723: Input Parameters:
3724: + dm - The `DM`
3725: . cellIS - The cells to include
3726: . locX - A local vector with the solution fields
3727: . locX_t - A local vector with solution field time derivatives, or `NULL`
3728: - locA - A local vector with auxiliary fields, or `NULL`
3730: Output Parameters:
3731: + u - The field coefficients
3732: . u_t - The fields derivative coefficients
3733: - a - The auxiliary field coefficients
3735: Level: developer
3737: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetFaceFields()`
3738: @*/
3739: PetscErrorCode DMPlexGetCellFields(DM dm, IS cellIS, Vec locX, PeOp Vec locX_t, PeOp Vec locA, PetscScalar *u[], PetscScalar *u_t[], PetscScalar *a[])
3740: {
3741: DM plex, plexA = NULL;
3742: DMEnclosureType encAux;
3743: PetscSection section, sectionAux;
3744: PetscDS prob;
3745: const PetscInt *cells;
3746: PetscInt cStart, cEnd, numCells, totDim, totDimAux, c;
3748: PetscFunctionBegin;
3753: PetscAssertPointer(u, 6);
3754: PetscAssertPointer(u_t, 7);
3755: PetscAssertPointer(a, 8);
3756: PetscCall(DMPlexConvertPlex(dm, &plex, PETSC_FALSE));
3757: PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
3758: PetscCall(DMGetLocalSection(dm, §ion));
3759: PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &prob, NULL));
3760: PetscCall(PetscDSGetTotalDimension(prob, &totDim));
3761: if (locA) {
3762: DM dmAux;
3763: PetscDS probAux;
3765: PetscCall(VecGetDM(locA, &dmAux));
3766: PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
3767: PetscCall(DMPlexConvertPlex(dmAux, &plexA, PETSC_FALSE));
3768: PetscCall(DMGetLocalSection(dmAux, §ionAux));
3769: PetscCall(DMGetDS(dmAux, &probAux));
3770: PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
3771: }
3772: numCells = cEnd - cStart;
3773: PetscCall(DMGetWorkArray(dm, numCells * totDim, MPIU_SCALAR, u));
3774: if (locX_t) PetscCall(DMGetWorkArray(dm, numCells * totDim, MPIU_SCALAR, u_t));
3775: else *u_t = NULL;
3776: if (locA) PetscCall(DMGetWorkArray(dm, numCells * totDimAux, MPIU_SCALAR, a));
3777: else *a = NULL;
3778: for (c = cStart; c < cEnd; ++c) {
3779: const PetscInt cell = cells ? cells[c] : c;
3780: const PetscInt cind = c - cStart;
3781: PetscScalar *x = NULL, *x_t = NULL, *ul = *u, *ul_t = *u_t, *al = *a;
3782: PetscInt i;
3784: PetscCall(DMPlexVecGetClosure(plex, section, locX, cell, NULL, &x));
3785: for (i = 0; i < totDim; ++i) ul[cind * totDim + i] = x[i];
3786: PetscCall(DMPlexVecRestoreClosure(plex, section, locX, cell, NULL, &x));
3787: if (locX_t) {
3788: PetscCall(DMPlexVecGetClosure(plex, section, locX_t, cell, NULL, &x_t));
3789: for (i = 0; i < totDim; ++i) ul_t[cind * totDim + i] = x_t[i];
3790: PetscCall(DMPlexVecRestoreClosure(plex, section, locX_t, cell, NULL, &x_t));
3791: }
3792: if (locA) {
3793: PetscInt subcell;
3794: PetscCall(DMGetEnclosurePoint(plexA, dm, encAux, cell, &subcell));
3795: PetscCall(DMPlexVecGetClosure(plexA, sectionAux, locA, subcell, NULL, &x));
3796: for (i = 0; i < totDimAux; ++i) al[cind * totDimAux + i] = x[i];
3797: PetscCall(DMPlexVecRestoreClosure(plexA, sectionAux, locA, subcell, NULL, &x));
3798: }
3799: }
3800: PetscCall(DMDestroy(&plex));
3801: if (locA) PetscCall(DMDestroy(&plexA));
3802: PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
3803: PetscFunctionReturn(PETSC_SUCCESS);
3804: }
3806: /*@C
3807: DMPlexRestoreCellFields - Restore the field values values for a chunk of cells
3809: Input Parameters:
3810: + dm - The `DM`
3811: . cellIS - The cells to include
3812: . locX - A local vector with the solution fields
3813: . locX_t - A local vector with solution field time derivatives, or `NULL`
3814: - locA - A local vector with auxiliary fields, or `NULL`
3816: Output Parameters:
3817: + u - The field coefficients
3818: . u_t - The fields derivative coefficients
3819: - a - The auxiliary field coefficients
3821: Level: developer
3823: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetFaceFields()`
3824: @*/
3825: PetscErrorCode DMPlexRestoreCellFields(DM dm, IS cellIS, Vec locX, PeOp Vec locX_t, PeOp Vec locA, PetscScalar *u[], PetscScalar *u_t[], PetscScalar *a[])
3826: {
3827: PetscFunctionBegin;
3828: PetscCall(DMRestoreWorkArray(dm, 0, MPIU_SCALAR, u));
3829: if (locX_t) PetscCall(DMRestoreWorkArray(dm, 0, MPIU_SCALAR, u_t));
3830: if (locA) PetscCall(DMRestoreWorkArray(dm, 0, MPIU_SCALAR, a));
3831: PetscFunctionReturn(PETSC_SUCCESS);
3832: }
3834: static PetscErrorCode DMPlexGetHybridCellFields(DM dm, IS cellIS, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a)
3835: {
3836: DM plex, plexA = NULL;
3837: DMEnclosureType encAux;
3838: PetscSection section, sectionAux;
3839: PetscDS ds, dsIn;
3840: const PetscInt *cells;
3841: PetscInt cStart, cEnd, numCells, c, totDim, totDimAux, Nf, f;
3843: PetscFunctionBegin;
3849: PetscAssertPointer(u, 6);
3850: PetscAssertPointer(u_t, 7);
3851: PetscAssertPointer(a, 8);
3852: PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
3853: numCells = cEnd - cStart;
3854: PetscCall(DMPlexConvertPlex(dm, &plex, PETSC_FALSE));
3855: PetscCall(DMGetLocalSection(dm, §ion));
3856: PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &ds, &dsIn));
3857: PetscCall(PetscDSGetNumFields(dsIn, &Nf));
3858: PetscCall(PetscDSGetTotalDimension(dsIn, &totDim));
3859: if (locA) {
3860: DM dmAux;
3861: PetscDS probAux;
3863: PetscCall(VecGetDM(locA, &dmAux));
3864: PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
3865: PetscCall(DMPlexConvertPlex(dmAux, &plexA, PETSC_FALSE));
3866: PetscCall(DMGetLocalSection(dmAux, §ionAux));
3867: PetscCall(DMGetDS(dmAux, &probAux));
3868: PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
3869: }
3870: PetscCall(DMGetWorkArray(dm, numCells * totDim, MPIU_SCALAR, u));
3871: if (locX_t) PetscCall(DMGetWorkArray(dm, numCells * totDim, MPIU_SCALAR, u_t));
3872: else {
3873: *u_t = NULL;
3874: }
3875: if (locA) PetscCall(DMGetWorkArray(dm, numCells * totDimAux, MPIU_SCALAR, a));
3876: else {
3877: *a = NULL;
3878: }
3879: // Loop over cohesive cells
3880: for (c = cStart; c < cEnd; ++c) {
3881: const PetscInt cell = cells ? cells[c] : c;
3882: const PetscInt cind = c - cStart;
3883: PetscScalar *xf = NULL, *xc = NULL, *x = NULL, *xf_t = NULL, *xc_t = NULL;
3884: PetscScalar *ul = &(*u)[cind * totDim], *ul_t = PetscSafePointerPlusOffset(*u_t, cind * totDim);
3885: const PetscInt *cone, *ornt;
3886: PetscInt Nx = 0, Nxf, s;
3888: PetscCall(DMPlexGetCone(dm, cell, &cone));
3889: PetscCall(DMPlexGetConeOrientation(dm, cell, &ornt));
3890: // Put in cohesive unknowns
3891: PetscCall(DMPlexVecGetClosure(plex, section, locX, cell, &Nxf, &xf));
3892: if (locX_t) PetscCall(DMPlexVecGetClosure(plex, section, locX_t, cell, NULL, &xf_t));
3893: for (f = 0; f < Nf; ++f) {
3894: PetscInt fdofIn, foff, foffIn;
3895: PetscBool cohesive;
3897: PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive));
3898: if (!cohesive) continue;
3899: PetscCall(PetscDSGetFieldSize(dsIn, f, &fdofIn));
3900: PetscCall(PetscDSGetFieldOffsetCohesive(ds, f, &foff));
3901: PetscCall(PetscDSGetFieldOffsetCohesive(dsIn, f, &foffIn));
3902: for (PetscInt i = 0; i < fdofIn; ++i) ul[foffIn + i] = xf[foff + i];
3903: if (locX_t)
3904: for (PetscInt i = 0; i < fdofIn; ++i) ul_t[foffIn + i] = xf_t[foff + i];
3905: Nx += fdofIn;
3906: }
3907: PetscCall(DMPlexVecRestoreClosure(plex, section, locX, cell, &Nxf, &xf));
3908: if (locX_t) PetscCall(DMPlexVecRestoreClosure(plex, section, locX_t, cell, NULL, &xf_t));
3909: // Loop over sides of surface
3910: for (s = 0; s < 2; ++s) {
3911: const PetscInt *support;
3912: const PetscInt face = cone[s];
3913: PetscDS dsC;
3914: PetscInt ssize, ncell, Nxc;
3916: // I don't think I need the face to have 0 orientation in the hybrid cell
3917: //PetscCheck(!ornt[s], PETSC_COMM_SELF, PETSC_ERR_SUP, "Face %" PetscInt_FMT " in hybrid cell %" PetscInt_FMT " has orientation %" PetscInt_FMT " != 0", face, cell, ornt[s]);
3918: PetscCall(DMPlexGetSupport(dm, face, &support));
3919: PetscCall(DMPlexGetSupportSize(dm, face, &ssize));
3920: if (support[0] == cell) ncell = support[1];
3921: else if (support[1] == cell) ncell = support[0];
3922: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " does not have cell %" PetscInt_FMT " in its support", face, cell);
3923: // Get closure of both face and cell, stick in cell for normal fields and face for cohesive fields
3924: PetscCall(DMGetCellDS(dm, ncell, &dsC, NULL));
3925: PetscCall(DMPlexVecGetClosure(plex, section, locX, ncell, &Nxc, &xc));
3926: if (locX_t) PetscCall(DMPlexVecGetClosure(plex, section, locX_t, ncell, NULL, &xc_t));
3927: for (f = 0; f < Nf; ++f) {
3928: PetscInt fdofIn, foffIn, foff;
3929: PetscBool cohesive;
3931: PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive));
3932: if (cohesive) continue;
3933: PetscCall(PetscDSGetFieldSize(dsIn, f, &fdofIn));
3934: PetscCall(PetscDSGetFieldOffset(dsC, f, &foff));
3935: PetscCall(PetscDSGetFieldOffsetCohesive(dsIn, f, &foffIn));
3936: for (PetscInt i = 0; i < fdofIn; ++i) ul[foffIn + s * fdofIn + i] = xc[foff + i];
3937: if (locX_t)
3938: for (PetscInt i = 0; i < fdofIn; ++i) ul_t[foffIn + s * fdofIn + i] = xc_t[foff + i];
3939: Nx += fdofIn;
3940: }
3941: PetscCall(DMPlexVecRestoreClosure(plex, section, locX, ncell, &Nxc, &xc));
3942: if (locX_t) PetscCall(DMPlexVecRestoreClosure(plex, section, locX_t, ncell, NULL, &xc_t));
3943: }
3944: PetscCheck(Nx == totDim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Closure size %" PetscInt_FMT " for cell %" PetscInt_FMT " does not match DS size %" PetscInt_FMT, Nx, cell, totDim);
3946: if (locA) {
3947: PetscScalar *al = &(*a)[cind * totDimAux];
3948: PetscInt subcell;
3950: PetscCall(DMGetEnclosurePoint(plexA, dm, encAux, cell, &subcell));
3951: PetscCall(DMPlexVecGetClosure(plexA, sectionAux, locA, subcell, &Nx, &x));
3952: PetscCheck(Nx == totDimAux, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Closure size %" PetscInt_FMT " for subcell %" PetscInt_FMT "does not match DS size %" PetscInt_FMT, Nx, subcell, totDimAux);
3953: for (PetscInt i = 0; i < totDimAux; ++i) al[i] = x[i];
3954: PetscCall(DMPlexVecRestoreClosure(plexA, sectionAux, locA, subcell, &Nx, &x));
3955: }
3956: }
3957: PetscCall(DMDestroy(&plex));
3958: PetscCall(DMDestroy(&plexA));
3959: PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
3960: PetscFunctionReturn(PETSC_SUCCESS);
3961: }
3963: /*
3964: DMPlexGetHybridFields - Get the field values for the negative side (s = 0) and positive side (s = 1) of the interface
3966: Input Parameters:
3967: + dm - The full domain DM
3968: . dmX - An array of DM for the field, say an auxiliary DM, indexed by s
3969: . dsX - An array of PetscDS for the field, indexed by s
3970: . cellIS - The interface cells for which we want values
3971: . locX - An array of local vectors with the field values, indexed by s
3972: - useCell - Flag to have values come from neighboring cell rather than endcap face
3974: Output Parameter:
3975: . x - An array of field values, indexed by s
3977: Note:
3978: The arrays in `x` will be allocated using `DMGetWorkArray()`, and must be returned using `DMPlexRestoreHybridFields()`.
3980: Level: advanced
3982: .seealso: `DMPlexRestoreHybridFields()`, `DMGetWorkArray()`
3983: */
3984: static PetscErrorCode DMPlexGetHybridFields(DM dm, DM dmX[], PetscDS dsX[], IS cellIS, Vec locX[], PetscBool useCell, PetscScalar *x[])
3985: {
3986: DM plexX[2];
3987: DMEnclosureType encX[2];
3988: PetscSection sectionX[2];
3989: const PetscInt *cells;
3990: PetscInt cStart, cEnd, numCells, c, s, totDimX[2];
3992: PetscFunctionBegin;
3993: PetscAssertPointer(locX, 5);
3994: if (!locX[0] || !locX[1]) PetscFunctionReturn(PETSC_SUCCESS);
3995: PetscAssertPointer(dmX, 2);
3996: PetscAssertPointer(dsX, 3);
3998: PetscAssertPointer(x, 7);
3999: PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
4000: numCells = cEnd - cStart;
4001: for (s = 0; s < 2; ++s) {
4005: PetscCall(DMPlexConvertPlex(dmX[s], &plexX[s], PETSC_FALSE));
4006: PetscCall(DMGetEnclosureRelation(dmX[s], dm, &encX[s]));
4007: PetscCall(DMGetLocalSection(dmX[s], §ionX[s]));
4008: PetscCall(PetscDSGetTotalDimension(dsX[s], &totDimX[s]));
4009: PetscCall(DMGetWorkArray(dmX[s], numCells * totDimX[s], MPIU_SCALAR, &x[s]));
4010: }
4011: for (c = cStart; c < cEnd; ++c) {
4012: const PetscInt cell = cells ? cells[c] : c;
4013: const PetscInt cind = c - cStart;
4014: const PetscInt *cone, *ornt;
4016: PetscCall(DMPlexGetCone(dm, cell, &cone));
4017: PetscCall(DMPlexGetConeOrientation(dm, cell, &ornt));
4018: //PetscCheck(!ornt[0], PETSC_COMM_SELF, PETSC_ERR_SUP, "Face %" PetscInt_FMT " in hybrid cell %" PetscInt_FMT " has orientation %" PetscInt_FMT " != 0", cone[0], cell, ornt[0]);
4019: for (s = 0; s < 2; ++s) {
4020: const PetscInt tdX = totDimX[s];
4021: PetscScalar *closure = NULL, *xl = &x[s][cind * tdX];
4022: PetscInt face = cone[s], point = face, subpoint, Nx, i;
4024: if (useCell) {
4025: const PetscInt *support;
4026: PetscInt ssize;
4028: PetscCall(DMPlexGetSupport(dm, face, &support));
4029: PetscCall(DMPlexGetSupportSize(dm, face, &ssize));
4030: PetscCheck(ssize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " from cell %" PetscInt_FMT " has support size %" PetscInt_FMT " != 2", face, cell, ssize);
4031: if (support[0] == cell) point = support[1];
4032: else if (support[1] == cell) point = support[0];
4033: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " does not have cell %" PetscInt_FMT " in its support", face, cell);
4034: }
4035: PetscCall(DMGetEnclosurePoint(plexX[s], dm, encX[s], point, &subpoint));
4036: PetscCall(DMPlexVecGetOrientedClosure(plexX[s], sectionX[s], PETSC_FALSE, locX[s], subpoint, ornt[s], &Nx, &closure));
4037: PetscCheck(Nx == tdX, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Closure size %" PetscInt_FMT " for subpoint %" PetscInt_FMT " does not match DS size %" PetscInt_FMT, Nx, subpoint, tdX);
4038: for (i = 0; i < Nx; ++i) xl[i] = closure[i];
4039: PetscCall(DMPlexVecRestoreClosure(plexX[s], sectionX[s], locX[s], subpoint, &Nx, &closure));
4040: }
4041: }
4042: for (s = 0; s < 2; ++s) PetscCall(DMDestroy(&plexX[s]));
4043: PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
4044: PetscFunctionReturn(PETSC_SUCCESS);
4045: }
4047: static PetscErrorCode DMPlexRestoreHybridFields(DM dm, DM dmX[], PetscDS dsX[], IS cellIS, Vec locX[], PetscBool useCell, PetscScalar *x[])
4048: {
4049: PetscFunctionBegin;
4050: if (!locX[0] || !locX[1]) PetscFunctionReturn(PETSC_SUCCESS);
4051: PetscCall(DMRestoreWorkArray(dmX[0], 0, MPIU_SCALAR, &x[0]));
4052: PetscCall(DMRestoreWorkArray(dmX[1], 0, MPIU_SCALAR, &x[1]));
4053: PetscFunctionReturn(PETSC_SUCCESS);
4054: }
4056: /*@C
4057: DMPlexGetFaceFields - Retrieve the field values values for a chunk of faces
4059: Input Parameters:
4060: + dm - The `DM`
4061: . fStart - The first face to include
4062: . fEnd - The first face to exclude
4063: . locX - A local vector with the solution fields
4064: . locX_t - A local vector with solution field time derivatives, or `NULL`
4065: . faceGeometry - A local vector with face geometry
4066: . cellGeometry - A local vector with cell geometry
4067: - locGrad - A local vector with field gradients, or `NULL`
4069: Output Parameters:
4070: + Nface - The number of faces with field values
4071: . uL - The field values at the left side of the face
4072: - uR - The field values at the right side of the face
4074: Level: developer
4076: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCellFields()`
4077: @*/
4078: PetscErrorCode DMPlexGetFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, PeOp Vec locX_t, Vec faceGeometry, Vec cellGeometry, PeOp Vec locGrad, PetscInt *Nface, PetscScalar *uL[], PetscScalar *uR[])
4079: {
4080: DM dmFace, dmCell, dmGrad = NULL;
4081: PetscSection section;
4082: PetscDS prob;
4083: DMLabel ghostLabel;
4084: const PetscScalar *facegeom, *cellgeom, *x, *lgrad;
4085: PetscBool *isFE;
4086: PetscInt dim, Nf, f, Nc, numFaces = fEnd - fStart, iface, face;
4088: PetscFunctionBegin;
4095: PetscAssertPointer(uL, 10);
4096: PetscAssertPointer(uR, 11);
4097: PetscCall(DMGetDimension(dm, &dim));
4098: PetscCall(DMGetDS(dm, &prob));
4099: PetscCall(DMGetLocalSection(dm, §ion));
4100: PetscCall(PetscDSGetNumFields(prob, &Nf));
4101: PetscCall(PetscDSGetTotalComponents(prob, &Nc));
4102: PetscCall(PetscMalloc1(Nf, &isFE));
4103: for (f = 0; f < Nf; ++f) {
4104: PetscObject obj;
4105: PetscClassId id;
4107: PetscCall(PetscDSGetDiscretization(prob, f, &obj));
4108: PetscCall(PetscObjectGetClassId(obj, &id));
4109: if (id == PETSCFE_CLASSID) {
4110: isFE[f] = PETSC_TRUE;
4111: } else if (id == PETSCFV_CLASSID) {
4112: isFE[f] = PETSC_FALSE;
4113: } else {
4114: isFE[f] = PETSC_FALSE;
4115: }
4116: }
4117: PetscCall(DMGetLabel(dm, "ghost", &ghostLabel));
4118: PetscCall(VecGetArrayRead(locX, &x));
4119: PetscCall(VecGetDM(faceGeometry, &dmFace));
4120: PetscCall(VecGetArrayRead(faceGeometry, &facegeom));
4121: PetscCall(VecGetDM(cellGeometry, &dmCell));
4122: PetscCall(VecGetArrayRead(cellGeometry, &cellgeom));
4123: if (locGrad) {
4124: PetscCall(VecGetDM(locGrad, &dmGrad));
4125: PetscCall(VecGetArrayRead(locGrad, &lgrad));
4126: }
4127: PetscCall(DMGetWorkArray(dm, numFaces * Nc, MPIU_SCALAR, uL));
4128: PetscCall(DMGetWorkArray(dm, numFaces * Nc, MPIU_SCALAR, uR));
4129: /* Right now just eat the extra work for FE (could make a cell loop) */
4130: for (face = fStart, iface = 0; face < fEnd; ++face) {
4131: const PetscInt *cells;
4132: PetscFVFaceGeom *fg;
4133: PetscFVCellGeom *cgL, *cgR;
4134: PetscScalar *xL, *xR, *gL, *gR;
4135: PetscScalar *uLl = *uL, *uRl = *uR;
4136: PetscInt ghost, nsupp, nchild;
4138: PetscCall(DMLabelGetValue(ghostLabel, face, &ghost));
4139: PetscCall(DMPlexGetSupportSize(dm, face, &nsupp));
4140: PetscCall(DMPlexGetTreeChildren(dm, face, &nchild, NULL));
4141: if (ghost >= 0 || nsupp > 2 || nchild > 0) continue;
4142: PetscCall(DMPlexPointLocalRead(dmFace, face, facegeom, &fg));
4143: PetscCall(DMPlexGetSupport(dm, face, &cells));
4144: PetscCall(DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL));
4145: PetscCall(DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR));
4146: for (f = 0; f < Nf; ++f) {
4147: PetscInt off;
4149: PetscCall(PetscDSGetComponentOffset(prob, f, &off));
4150: if (isFE[f]) {
4151: const PetscInt *cone;
4152: PetscInt comp, coneSizeL, coneSizeR, faceLocL, faceLocR, ldof, rdof, d;
4154: xL = xR = NULL;
4155: PetscCall(PetscSectionGetFieldComponents(section, f, &comp));
4156: PetscCall(DMPlexVecGetClosure(dm, section, locX, cells[0], &ldof, &xL));
4157: PetscCall(DMPlexVecGetClosure(dm, section, locX, cells[1], &rdof, &xR));
4158: PetscCall(DMPlexGetCone(dm, cells[0], &cone));
4159: PetscCall(DMPlexGetConeSize(dm, cells[0], &coneSizeL));
4160: for (faceLocL = 0; faceLocL < coneSizeL; ++faceLocL)
4161: if (cone[faceLocL] == face) break;
4162: PetscCall(DMPlexGetCone(dm, cells[1], &cone));
4163: PetscCall(DMPlexGetConeSize(dm, cells[1], &coneSizeR));
4164: for (faceLocR = 0; faceLocR < coneSizeR; ++faceLocR)
4165: if (cone[faceLocR] == face) break;
4166: PetscCheck(faceLocL != coneSizeL || faceLocR != coneSizeR, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %" PetscInt_FMT " in cone of cell %" PetscInt_FMT " or cell %" PetscInt_FMT, face, cells[0], cells[1]);
4167: /* Check that FEM field has values in the right cell (sometimes its an FV ghost cell) */
4168: /* TODO: this is a hack that might not be right for nonconforming */
4169: if (faceLocL < coneSizeL) {
4170: PetscCall(PetscFEEvaluateFaceFields_Internal(prob, f, faceLocL, xL, &uLl[iface * Nc + off]));
4171: if (rdof == ldof && faceLocR < coneSizeR) PetscCall(PetscFEEvaluateFaceFields_Internal(prob, f, faceLocR, xR, &uRl[iface * Nc + off]));
4172: else {
4173: for (d = 0; d < comp; ++d) uRl[iface * Nc + off + d] = uLl[iface * Nc + off + d];
4174: }
4175: } else {
4176: PetscCall(PetscFEEvaluateFaceFields_Internal(prob, f, faceLocR, xR, &uRl[iface * Nc + off]));
4177: PetscCall(PetscSectionGetFieldComponents(section, f, &comp));
4178: for (d = 0; d < comp; ++d) uLl[iface * Nc + off + d] = uRl[iface * Nc + off + d];
4179: }
4180: PetscCall(DMPlexVecRestoreClosure(dm, section, locX, cells[0], &ldof, &xL));
4181: PetscCall(DMPlexVecRestoreClosure(dm, section, locX, cells[1], &rdof, &xR));
4182: } else {
4183: PetscFV fv;
4184: PetscInt numComp, c;
4186: PetscCall(PetscDSGetDiscretization(prob, f, (PetscObject *)&fv));
4187: PetscCall(PetscFVGetNumComponents(fv, &numComp));
4188: PetscCall(DMPlexPointLocalFieldRead(dm, cells[0], f, x, &xL));
4189: PetscCall(DMPlexPointLocalFieldRead(dm, cells[1], f, x, &xR));
4190: if (dmGrad) {
4191: PetscReal dxL[3], dxR[3];
4193: PetscCall(DMPlexPointLocalRead(dmGrad, cells[0], lgrad, &gL));
4194: PetscCall(DMPlexPointLocalRead(dmGrad, cells[1], lgrad, &gR));
4195: DMPlex_WaxpyD_Internal(dim, -1, cgL->centroid, fg->centroid, dxL);
4196: DMPlex_WaxpyD_Internal(dim, -1, cgR->centroid, fg->centroid, dxR);
4197: for (c = 0; c < numComp; ++c) {
4198: uLl[iface * Nc + off + c] = xL[c] + DMPlex_DotD_Internal(dim, &gL[c * dim], dxL);
4199: uRl[iface * Nc + off + c] = xR[c] + DMPlex_DotD_Internal(dim, &gR[c * dim], dxR);
4200: }
4201: } else {
4202: for (c = 0; c < numComp; ++c) {
4203: uLl[iface * Nc + off + c] = xL[c];
4204: uRl[iface * Nc + off + c] = xR[c];
4205: }
4206: }
4207: }
4208: }
4209: ++iface;
4210: }
4211: *Nface = iface;
4212: PetscCall(VecRestoreArrayRead(locX, &x));
4213: PetscCall(VecRestoreArrayRead(faceGeometry, &facegeom));
4214: PetscCall(VecRestoreArrayRead(cellGeometry, &cellgeom));
4215: if (locGrad) PetscCall(VecRestoreArrayRead(locGrad, &lgrad));
4216: PetscCall(PetscFree(isFE));
4217: PetscFunctionReturn(PETSC_SUCCESS);
4218: }
4220: /*@C
4221: DMPlexRestoreFaceFields - Restore the field values values for a chunk of faces
4223: Input Parameters:
4224: + dm - The `DM`
4225: . fStart - The first face to include
4226: . fEnd - The first face to exclude
4227: . locX - A local vector with the solution fields
4228: . locX_t - A local vector with solution field time derivatives, or `NULL`
4229: . faceGeometry - A local vector with face geometry
4230: . cellGeometry - A local vector with cell geometry
4231: - locGrad - A local vector with field gradients, or `NULL`
4233: Output Parameters:
4234: + Nface - The number of faces with field values
4235: . uL - The field values at the left side of the face
4236: - uR - The field values at the right side of the face
4238: Level: developer
4240: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetFaceFields()`
4241: @*/
4242: PetscErrorCode DMPlexRestoreFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, PeOp Vec locX_t, Vec faceGeometry, Vec cellGeometry, PeOp Vec locGrad, PetscInt *Nface, PetscScalar *uL[], PetscScalar *uR[])
4243: {
4244: PetscFunctionBegin;
4245: PetscCall(DMRestoreWorkArray(dm, 0, MPIU_SCALAR, uL));
4246: PetscCall(DMRestoreWorkArray(dm, 0, MPIU_SCALAR, uR));
4247: PetscFunctionReturn(PETSC_SUCCESS);
4248: }
4250: /*@C
4251: DMPlexGetFaceGeometry - Retrieve the geometric values for a chunk of faces
4253: Input Parameters:
4254: + dm - The `DM`
4255: . fStart - The first face to include
4256: . fEnd - The first face to exclude
4257: . faceGeometry - A local vector with face geometry
4258: - cellGeometry - A local vector with cell geometry
4260: Output Parameters:
4261: + Nface - The number of faces with field values
4262: . fgeom - The face centroid and normals
4263: - vol - The cell volumes
4265: Level: developer
4267: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetCellFields()`
4268: @*/
4269: PetscErrorCode DMPlexGetFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscInt *Nface, PetscFVFaceGeom *fgeom[], PetscReal *vol[])
4270: {
4271: DM dmFace, dmCell;
4272: DMLabel ghostLabel;
4273: const PetscScalar *facegeom, *cellgeom;
4274: PetscInt dim, numFaces = fEnd - fStart, iface, face;
4276: PetscFunctionBegin;
4280: PetscAssertPointer(fgeom, 7);
4281: PetscAssertPointer(vol, 8);
4282: PetscCall(DMGetDimension(dm, &dim));
4283: PetscCall(DMGetLabel(dm, "ghost", &ghostLabel));
4284: PetscCall(VecGetDM(faceGeometry, &dmFace));
4285: PetscCall(VecGetArrayRead(faceGeometry, &facegeom));
4286: PetscCall(VecGetDM(cellGeometry, &dmCell));
4287: PetscCall(VecGetArrayRead(cellGeometry, &cellgeom));
4288: PetscCall(PetscMalloc1(numFaces, fgeom));
4289: PetscCall(DMGetWorkArray(dm, numFaces * 2, MPIU_SCALAR, vol));
4290: for (face = fStart, iface = 0; face < fEnd; ++face) {
4291: const PetscInt *cells;
4292: PetscFVFaceGeom *fg;
4293: PetscFVCellGeom *cgL, *cgR;
4294: PetscFVFaceGeom *fgeoml = *fgeom;
4295: PetscReal *voll = *vol;
4296: PetscInt ghost, d, nchild, nsupp;
4298: PetscCall(DMLabelGetValue(ghostLabel, face, &ghost));
4299: PetscCall(DMPlexGetSupportSize(dm, face, &nsupp));
4300: PetscCall(DMPlexGetTreeChildren(dm, face, &nchild, NULL));
4301: if (ghost >= 0 || nsupp > 2 || nchild > 0) continue;
4302: PetscCall(DMPlexPointLocalRead(dmFace, face, facegeom, &fg));
4303: PetscCall(DMPlexGetSupport(dm, face, &cells));
4304: PetscCall(DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL));
4305: PetscCall(DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR));
4306: for (d = 0; d < dim; ++d) {
4307: fgeoml[iface].centroid[d] = fg->centroid[d];
4308: fgeoml[iface].normal[d] = fg->normal[d];
4309: }
4310: voll[iface * 2 + 0] = cgL->volume;
4311: voll[iface * 2 + 1] = cgR->volume;
4312: ++iface;
4313: }
4314: *Nface = iface;
4315: PetscCall(VecRestoreArrayRead(faceGeometry, &facegeom));
4316: PetscCall(VecRestoreArrayRead(cellGeometry, &cellgeom));
4317: PetscFunctionReturn(PETSC_SUCCESS);
4318: }
4320: /*@C
4321: DMPlexRestoreFaceGeometry - Restore the field values values for a chunk of faces
4323: Input Parameters:
4324: + dm - The `DM`
4325: . fStart - The first face to include
4326: . fEnd - The first face to exclude
4327: . faceGeometry - A local vector with face geometry
4328: - cellGeometry - A local vector with cell geometry
4330: Output Parameters:
4331: + Nface - The number of faces with field values
4332: . fgeom - The face centroid and normals
4333: - vol - The cell volumes
4335: Level: developer
4337: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetFaceFields()`
4338: @*/
4339: PetscErrorCode DMPlexRestoreFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscInt *Nface, PetscFVFaceGeom *fgeom[], PetscReal *vol[])
4340: {
4341: PetscFunctionBegin;
4342: PetscCall(PetscFree(*fgeom));
4343: PetscCall(DMRestoreWorkArray(dm, 0, MPIU_REAL, vol));
4344: PetscFunctionReturn(PETSC_SUCCESS);
4345: }
4347: PetscErrorCode DMSNESGetFEGeom(DMField coordField, IS pointIS, PetscQuadrature quad, PetscFEGeomMode mode, PetscFEGeom **geom)
4348: {
4349: char composeStr[33] = {0};
4350: PetscObjectId id;
4351: PetscContainer container;
4353: PetscFunctionBegin;
4354: PetscCall(PetscObjectGetId((PetscObject)quad, &id));
4355: PetscCall(PetscSNPrintf(composeStr, 32, "DMSNESGetFEGeom_%" PetscInt64_FMT "\n", id));
4356: PetscCall(PetscObjectQuery((PetscObject)pointIS, composeStr, (PetscObject *)&container));
4357: if (container) {
4358: PetscCall(PetscContainerGetPointer(container, (void **)geom));
4359: } else {
4360: PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, mode, geom));
4361: PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
4362: PetscCall(PetscContainerSetPointer(container, (void *)*geom));
4363: PetscCall(PetscContainerSetCtxDestroy(container, PetscContainerCtxDestroy_PetscFEGeom));
4364: PetscCall(PetscObjectCompose((PetscObject)pointIS, composeStr, (PetscObject)container));
4365: PetscCall(PetscContainerDestroy(&container));
4366: }
4367: PetscFunctionReturn(PETSC_SUCCESS);
4368: }
4370: PetscErrorCode DMSNESRestoreFEGeom(DMField coordField, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
4371: {
4372: PetscFunctionBegin;
4373: *geom = NULL;
4374: PetscFunctionReturn(PETSC_SUCCESS);
4375: }
4377: PetscErrorCode DMPlexComputeResidual_Patch_Internal(DM dm, PetscSection section, IS cellIS, PetscReal t, Vec locX, Vec locX_t, Vec locF, void *user)
4378: {
4379: DM_Plex *mesh = (DM_Plex *)dm->data;
4380: const char *name = "Residual";
4381: DM dmAux = NULL;
4382: DMLabel ghostLabel = NULL;
4383: PetscDS prob = NULL;
4384: PetscDS probAux = NULL;
4385: PetscBool useFEM = PETSC_FALSE;
4386: PetscBool isImplicit = (locX_t || t == PETSC_MIN_REAL) ? PETSC_TRUE : PETSC_FALSE;
4387: DMField coordField = NULL;
4388: Vec locA;
4389: PetscScalar *u = NULL, *u_t, *a, *uL = NULL, *uR = NULL;
4390: IS chunkIS;
4391: const PetscInt *cells;
4392: PetscInt cStart, cEnd, numCells;
4393: PetscInt Nf, f, totDim, totDimAux, numChunks, cellChunkSize, chunk, fStart, fEnd;
4394: PetscInt maxDegree = PETSC_INT_MAX;
4395: PetscFormKey key;
4396: PetscQuadrature affineQuad = NULL, *quads = NULL;
4397: PetscFEGeom *affineGeom = NULL, **geoms = NULL;
4399: PetscFunctionBegin;
4400: PetscCall(PetscLogEventBegin(DMPLEX_ResidualFEM, dm, 0, 0, 0));
4401: /* FEM+FVM */
4402: /* 1: Get sizes from dm and dmAux */
4403: PetscCall(DMGetLabel(dm, "ghost", &ghostLabel));
4404: PetscCall(DMGetDS(dm, &prob));
4405: PetscCall(PetscDSGetNumFields(prob, &Nf));
4406: PetscCall(PetscDSGetTotalDimension(prob, &totDim));
4407: PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &locA));
4408: if (locA) {
4409: PetscCall(VecGetDM(locA, &dmAux));
4410: PetscCall(DMGetDS(dmAux, &probAux));
4411: PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
4412: }
4413: /* 2: Get geometric data */
4414: for (f = 0; f < Nf; ++f) {
4415: PetscObject obj;
4416: PetscClassId id;
4417: PetscBool fimp;
4419: PetscCall(PetscDSGetImplicit(prob, f, &fimp));
4420: if (isImplicit != fimp) continue;
4421: PetscCall(PetscDSGetDiscretization(prob, f, &obj));
4422: PetscCall(PetscObjectGetClassId(obj, &id));
4423: if (id == PETSCFE_CLASSID) useFEM = PETSC_TRUE;
4424: PetscCheck(id != PETSCFV_CLASSID, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Use of FVM with PCPATCH not yet implemented");
4425: }
4426: if (useFEM) {
4427: PetscCall(DMGetCoordinateField(dm, &coordField));
4428: PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
4429: if (maxDegree <= 1) {
4430: PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, &affineQuad));
4431: if (affineQuad) PetscCall(DMSNESGetFEGeom(coordField, cellIS, affineQuad, PETSC_FEGEOM_BASIC, &affineGeom));
4432: } else {
4433: PetscCall(PetscCalloc2(Nf, &quads, Nf, &geoms));
4434: for (f = 0; f < Nf; ++f) {
4435: PetscObject obj;
4436: PetscClassId id;
4437: PetscBool fimp;
4439: PetscCall(PetscDSGetImplicit(prob, f, &fimp));
4440: if (isImplicit != fimp) continue;
4441: PetscCall(PetscDSGetDiscretization(prob, f, &obj));
4442: PetscCall(PetscObjectGetClassId(obj, &id));
4443: if (id == PETSCFE_CLASSID) {
4444: PetscFE fe = (PetscFE)obj;
4446: PetscCall(PetscFEGetQuadrature(fe, &quads[f]));
4447: PetscCall(PetscObjectReference((PetscObject)quads[f]));
4448: PetscCall(DMSNESGetFEGeom(coordField, cellIS, quads[f], PETSC_FEGEOM_BASIC, &geoms[f]));
4449: }
4450: }
4451: }
4452: }
4453: /* Loop over chunks */
4454: PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
4455: PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
4456: if (useFEM) PetscCall(ISCreate(PETSC_COMM_SELF, &chunkIS));
4457: numCells = cEnd - cStart;
4458: numChunks = 1;
4459: cellChunkSize = numCells / numChunks;
4460: numChunks = PetscMin(1, numCells);
4461: key.label = NULL;
4462: key.value = 0;
4463: key.part = 0;
4464: for (chunk = 0; chunk < numChunks; ++chunk) {
4465: PetscScalar *elemVec, *fluxL = NULL, *fluxR = NULL;
4466: PetscReal *vol = NULL;
4467: PetscFVFaceGeom *fgeom = NULL;
4468: PetscInt cS = cStart + chunk * cellChunkSize, cE = PetscMin(cS + cellChunkSize, cEnd), numCells = cE - cS, c;
4469: PetscInt numFaces = 0;
4471: /* Extract field coefficients */
4472: if (useFEM) {
4473: PetscCall(ISGetPointSubrange(chunkIS, cS, cE, cells));
4474: PetscCall(DMPlexGetCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a));
4475: PetscCall(DMGetWorkArray(dm, numCells * totDim, MPIU_SCALAR, &elemVec));
4476: PetscCall(PetscArrayzero(elemVec, numCells * totDim));
4477: }
4478: /* TODO We will interlace both our field coefficients (u, u_t, uL, uR, etc.) and our output (elemVec, fL, fR). I think this works */
4479: /* Loop over fields */
4480: for (f = 0; f < Nf; ++f) {
4481: PetscObject obj;
4482: PetscClassId id;
4483: PetscBool fimp;
4484: PetscInt numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;
4486: key.field = f;
4487: PetscCall(PetscDSGetImplicit(prob, f, &fimp));
4488: if (isImplicit != fimp) continue;
4489: PetscCall(PetscDSGetDiscretization(prob, f, &obj));
4490: PetscCall(PetscObjectGetClassId(obj, &id));
4491: if (id == PETSCFE_CLASSID) {
4492: PetscFE fe = (PetscFE)obj;
4493: PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f];
4494: PetscFEGeom *chunkGeom = NULL;
4495: PetscQuadrature quad = affineQuad ? affineQuad : quads[f];
4496: PetscInt Nq, Nb;
4498: PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
4499: PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
4500: PetscCall(PetscFEGetDimension(fe, &Nb));
4501: blockSize = Nb;
4502: batchSize = numBlocks * blockSize;
4503: PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
4504: numChunks = numCells / (numBatches * batchSize);
4505: Ne = numChunks * numBatches * batchSize;
4506: Nr = numCells % (numBatches * batchSize);
4507: offset = numCells - Nr;
4508: /* Integrate FE residual to get elemVec (need fields at quadrature points) */
4509: /* 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) */
4510: PetscCall(PetscFEGeomGetChunk(geom, 0, offset, &chunkGeom));
4511: PetscCall(PetscFEIntegrateResidual(prob, key, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec));
4512: PetscCall(PetscFEGeomGetChunk(geom, offset, numCells, &chunkGeom));
4513: PetscCall(PetscFEIntegrateResidual(prob, key, Nr, chunkGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), probAux, &a[offset * totDimAux], t, &elemVec[offset * totDim]));
4514: PetscCall(PetscFEGeomRestoreChunk(geom, offset, numCells, &chunkGeom));
4515: } else if (id == PETSCFV_CLASSID) {
4516: PetscFV fv = (PetscFV)obj;
4518: Ne = numFaces;
4519: /* Riemann solve over faces (need fields at face centroids) */
4520: /* We need to evaluate FE fields at those coordinates */
4521: PetscCall(PetscFVIntegrateRHSFunction(fv, prob, f, Ne, fgeom, vol, uL, uR, fluxL, fluxR));
4522: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
4523: }
4524: /* Loop over domain */
4525: if (useFEM) {
4526: /* Add elemVec to locX */
4527: for (c = cS; c < cE; ++c) {
4528: const PetscInt cell = cells ? cells[c] : c;
4529: const PetscInt cind = c - cStart;
4531: if (mesh->printFEM > 1) PetscCall(DMPrintCellVector(cell, name, totDim, &elemVec[cind * totDim]));
4532: if (ghostLabel) {
4533: PetscInt ghostVal;
4535: PetscCall(DMLabelGetValue(ghostLabel, cell, &ghostVal));
4536: if (ghostVal > 0) continue;
4537: }
4538: PetscCall(DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cind * totDim], ADD_ALL_VALUES));
4539: }
4540: }
4541: /* Handle time derivative */
4542: if (locX_t) {
4543: PetscScalar *x_t, *fa;
4545: PetscCall(VecGetArray(locF, &fa));
4546: PetscCall(VecGetArray(locX_t, &x_t));
4547: for (f = 0; f < Nf; ++f) {
4548: PetscFV fv;
4549: PetscObject obj;
4550: PetscClassId id;
4551: PetscInt pdim, d;
4553: PetscCall(PetscDSGetDiscretization(prob, f, &obj));
4554: PetscCall(PetscObjectGetClassId(obj, &id));
4555: if (id != PETSCFV_CLASSID) continue;
4556: fv = (PetscFV)obj;
4557: PetscCall(PetscFVGetNumComponents(fv, &pdim));
4558: for (c = cS; c < cE; ++c) {
4559: const PetscInt cell = cells ? cells[c] : c;
4560: PetscScalar *u_t, *r;
4562: if (ghostLabel) {
4563: PetscInt ghostVal;
4565: PetscCall(DMLabelGetValue(ghostLabel, cell, &ghostVal));
4566: if (ghostVal > 0) continue;
4567: }
4568: PetscCall(DMPlexPointLocalFieldRead(dm, cell, f, x_t, &u_t));
4569: PetscCall(DMPlexPointLocalFieldRef(dm, cell, f, fa, &r));
4570: for (d = 0; d < pdim; ++d) r[d] += u_t[d];
4571: }
4572: }
4573: PetscCall(VecRestoreArray(locX_t, &x_t));
4574: PetscCall(VecRestoreArray(locF, &fa));
4575: }
4576: if (useFEM) {
4577: PetscCall(DMPlexRestoreCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a));
4578: PetscCall(DMRestoreWorkArray(dm, numCells * totDim, MPIU_SCALAR, &elemVec));
4579: }
4580: }
4581: if (useFEM) PetscCall(ISDestroy(&chunkIS));
4582: PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
4583: /* TODO Could include boundary residual here (see DMPlexComputeResidualByKey) */
4584: if (useFEM) {
4585: if (maxDegree <= 1) {
4586: PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, affineQuad, PETSC_FALSE, &affineGeom));
4587: PetscCall(PetscQuadratureDestroy(&affineQuad));
4588: } else {
4589: for (f = 0; f < Nf; ++f) {
4590: PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, quads[f], PETSC_FALSE, &geoms[f]));
4591: PetscCall(PetscQuadratureDestroy(&quads[f]));
4592: }
4593: PetscCall(PetscFree2(quads, geoms));
4594: }
4595: }
4596: PetscCall(PetscLogEventEnd(DMPLEX_ResidualFEM, dm, 0, 0, 0));
4597: PetscFunctionReturn(PETSC_SUCCESS);
4598: }
4600: /*
4601: 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
4603: X - The local solution vector
4604: X_t - The local solution time derivative vector, or NULL
4605: */
4606: PetscErrorCode DMPlexComputeJacobian_Patch_Internal(DM dm, PetscSection section, PetscSection globalSection, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Mat Jac, Mat JacP, void *ctx)
4607: {
4608: DM_Plex *mesh = (DM_Plex *)dm->data;
4609: const char *name = "Jacobian", *nameP = "JacobianPre";
4610: DM dmAux = NULL;
4611: PetscDS prob, probAux = NULL;
4612: PetscSection sectionAux = NULL;
4613: Vec A;
4614: DMField coordField;
4615: PetscFEGeom *cgeomFEM;
4616: PetscQuadrature qGeom = NULL;
4617: Mat J = Jac, JP = JacP;
4618: PetscScalar *work, *u = NULL, *u_t = NULL, *a = NULL, *elemMat = NULL, *elemMatP = NULL, *elemMatD = NULL;
4619: PetscBool hasJac, hasPrec, hasDyn, assembleJac, *isFE, hasFV = PETSC_FALSE;
4620: const PetscInt *cells;
4621: PetscFormKey key;
4622: PetscInt Nf, fieldI, fieldJ, maxDegree, numCells, cStart, cEnd, numChunks, chunkSize, chunk, totDim, totDimAux = 0, sz, wsz, off = 0, offCell = 0;
4624: PetscFunctionBegin;
4625: PetscCall(ISGetLocalSize(cellIS, &numCells));
4626: PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
4627: PetscCall(PetscLogEventBegin(DMPLEX_JacobianFEM, dm, 0, 0, 0));
4628: PetscCall(DMGetDS(dm, &prob));
4629: PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &A));
4630: if (A) {
4631: PetscCall(VecGetDM(A, &dmAux));
4632: PetscCall(DMGetLocalSection(dmAux, §ionAux));
4633: PetscCall(DMGetDS(dmAux, &probAux));
4634: }
4635: /* Get flags */
4636: PetscCall(PetscDSGetNumFields(prob, &Nf));
4637: PetscCall(DMGetWorkArray(dm, Nf, MPIU_BOOL, &isFE));
4638: for (fieldI = 0; fieldI < Nf; ++fieldI) {
4639: PetscObject disc;
4640: PetscClassId id;
4641: PetscCall(PetscDSGetDiscretization(prob, fieldI, &disc));
4642: PetscCall(PetscObjectGetClassId(disc, &id));
4643: if (id == PETSCFE_CLASSID) {
4644: isFE[fieldI] = PETSC_TRUE;
4645: } else if (id == PETSCFV_CLASSID) {
4646: hasFV = PETSC_TRUE;
4647: isFE[fieldI] = PETSC_FALSE;
4648: }
4649: }
4650: PetscCall(PetscDSHasJacobian(prob, &hasJac));
4651: PetscCall(PetscDSHasJacobianPreconditioner(prob, &hasPrec));
4652: PetscCall(PetscDSHasDynamicJacobian(prob, &hasDyn));
4653: assembleJac = hasJac && hasPrec && (Jac != JacP) ? PETSC_TRUE : PETSC_FALSE;
4654: hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
4655: if (hasFV) PetscCall(MatSetOption(JP, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE)); /* No allocated space for FV stuff, so ignore the zero entries */
4656: PetscCall(PetscDSGetTotalDimension(prob, &totDim));
4657: if (probAux) PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
4658: /* Compute batch sizes */
4659: if (isFE[0]) {
4660: PetscFE fe;
4661: PetscQuadrature q;
4662: PetscInt numQuadPoints, numBatches, batchSize, numBlocks, blockSize, Nb;
4664: PetscCall(PetscDSGetDiscretization(prob, 0, (PetscObject *)&fe));
4665: PetscCall(PetscFEGetQuadrature(fe, &q));
4666: PetscCall(PetscQuadratureGetData(q, NULL, NULL, &numQuadPoints, NULL, NULL));
4667: PetscCall(PetscFEGetDimension(fe, &Nb));
4668: PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
4669: blockSize = Nb * numQuadPoints;
4670: batchSize = numBlocks * blockSize;
4671: chunkSize = numBatches * batchSize;
4672: numChunks = numCells / chunkSize + numCells % chunkSize;
4673: PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
4674: } else {
4675: chunkSize = numCells;
4676: numChunks = 1;
4677: }
4678: /* Get work space */
4679: wsz = (((X ? 1 : 0) + (X_t ? 1 : 0) + (dmAux ? 1 : 0)) * totDim + ((hasJac ? 1 : 0) + (hasPrec ? 1 : 0) + (hasDyn ? 1 : 0)) * totDim * totDim) * chunkSize;
4680: PetscCall(DMGetWorkArray(dm, wsz, MPIU_SCALAR, &work));
4681: PetscCall(PetscArrayzero(work, wsz));
4682: off = 0;
4683: u = X ? (sz = chunkSize * totDim, off += sz, work + off - sz) : NULL;
4684: u_t = X_t ? (sz = chunkSize * totDim, off += sz, work + off - sz) : NULL;
4685: a = dmAux ? (sz = chunkSize * totDimAux, off += sz, work + off - sz) : NULL;
4686: elemMat = hasJac ? (sz = chunkSize * totDim * totDim, off += sz, work + off - sz) : NULL;
4687: elemMatP = hasPrec ? (sz = chunkSize * totDim * totDim, off += sz, work + off - sz) : NULL;
4688: elemMatD = hasDyn ? (sz = chunkSize * totDim * totDim, off += sz, work + off - sz) : NULL;
4689: PetscCheck(off == wsz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error is workspace size %" PetscInt_FMT " should be %" PetscInt_FMT, off, wsz);
4690: /* Setup geometry */
4691: PetscCall(DMGetCoordinateField(dm, &coordField));
4692: PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
4693: if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, &qGeom));
4694: if (!qGeom) {
4695: PetscFE fe;
4697: PetscCall(PetscDSGetDiscretization(prob, 0, (PetscObject *)&fe));
4698: PetscCall(PetscFEGetQuadrature(fe, &qGeom));
4699: PetscCall(PetscObjectReference((PetscObject)qGeom));
4700: }
4701: PetscCall(DMSNESGetFEGeom(coordField, cellIS, qGeom, PETSC_FEGEOM_BASIC, &cgeomFEM));
4702: /* Compute volume integrals */
4703: if (assembleJac) PetscCall(MatZeroEntries(J));
4704: PetscCall(MatZeroEntries(JP));
4705: key.label = NULL;
4706: key.value = 0;
4707: key.part = 0;
4708: for (chunk = 0; chunk < numChunks; ++chunk, offCell += chunkSize) {
4709: const PetscInt Ncell = PetscMin(chunkSize, numCells - offCell);
4710: PetscInt c;
4712: /* Extract values */
4713: for (c = 0; c < Ncell; ++c) {
4714: const PetscInt cell = cells ? cells[c + offCell] : c + offCell;
4715: PetscScalar *x = NULL, *x_t = NULL;
4716: PetscInt i;
4718: if (X) {
4719: PetscCall(DMPlexVecGetClosure(dm, section, X, cell, NULL, &x));
4720: for (i = 0; i < totDim; ++i) u[c * totDim + i] = x[i];
4721: PetscCall(DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x));
4722: }
4723: if (X_t) {
4724: PetscCall(DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t));
4725: for (i = 0; i < totDim; ++i) u_t[c * totDim + i] = x_t[i];
4726: PetscCall(DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t));
4727: }
4728: if (dmAux) {
4729: PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, A, cell, NULL, &x));
4730: for (i = 0; i < totDimAux; ++i) a[c * totDimAux + i] = x[i];
4731: PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, A, cell, NULL, &x));
4732: }
4733: }
4734: for (fieldI = 0; fieldI < Nf; ++fieldI) {
4735: PetscFE fe;
4736: PetscCall(PetscDSGetDiscretization(prob, fieldI, (PetscObject *)&fe));
4737: for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
4738: key.field = fieldI * Nf + fieldJ;
4739: if (hasJac) PetscCall(PetscFEIntegrateJacobian(prob, prob, PETSCFE_JACOBIAN, key, Ncell, cgeomFEM, u, u_t, probAux, a, t, X_tShift, elemMat));
4740: if (hasPrec) PetscCall(PetscFEIntegrateJacobian(prob, prob, PETSCFE_JACOBIAN_PRE, key, Ncell, cgeomFEM, u, u_t, probAux, a, t, X_tShift, elemMatP));
4741: if (hasDyn) PetscCall(PetscFEIntegrateJacobian(prob, prob, PETSCFE_JACOBIAN_DYN, key, Ncell, cgeomFEM, u, u_t, probAux, a, t, X_tShift, elemMatD));
4742: }
4743: /* For finite volume, add the identity */
4744: if (!isFE[fieldI]) {
4745: PetscFV fv;
4746: PetscInt eOffset = 0, Nc, fc, foff;
4748: PetscCall(PetscDSGetFieldOffset(prob, fieldI, &foff));
4749: PetscCall(PetscDSGetDiscretization(prob, fieldI, (PetscObject *)&fv));
4750: PetscCall(PetscFVGetNumComponents(fv, &Nc));
4751: for (c = 0; c < chunkSize; ++c, eOffset += totDim * totDim) {
4752: for (fc = 0; fc < Nc; ++fc) {
4753: const PetscInt i = foff + fc;
4754: if (hasJac) elemMat[eOffset + i * totDim + i] = 1.0;
4755: if (hasPrec) elemMatP[eOffset + i * totDim + i] = 1.0;
4756: }
4757: }
4758: }
4759: }
4760: /* Add contribution from X_t */
4761: if (hasDyn) {
4762: for (c = 0; c < chunkSize * totDim * totDim; ++c) elemMat[c] += X_tShift * elemMatD[c];
4763: }
4764: /* Insert values into matrix */
4765: for (c = 0; c < Ncell; ++c) {
4766: const PetscInt cell = cells ? cells[c + offCell] : c + offCell;
4767: if (mesh->printFEM > 1) {
4768: if (hasJac) PetscCall(DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[(c - cStart) * totDim * totDim]));
4769: if (hasPrec) PetscCall(DMPrintCellMatrix(cell, nameP, totDim, totDim, &elemMatP[(c - cStart) * totDim * totDim]));
4770: }
4771: if (assembleJac) PetscCall(DMPlexMatSetClosure_Internal(dm, section, globalSection, mesh->useMatClPerm, Jac, cell, &elemMat[(c - cStart) * totDim * totDim], ADD_VALUES));
4772: PetscCall(DMPlexMatSetClosure_Internal(dm, section, globalSection, mesh->useMatClPerm, JP, cell, &elemMat[(c - cStart) * totDim * totDim], ADD_VALUES));
4773: }
4774: }
4775: /* Cleanup */
4776: PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, qGeom, PETSC_FALSE, &cgeomFEM));
4777: PetscCall(PetscQuadratureDestroy(&qGeom));
4778: if (hasFV) PetscCall(MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE));
4779: PetscCall(DMRestoreWorkArray(dm, Nf, MPIU_BOOL, &isFE));
4780: PetscCall(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));
4781: /* Compute boundary integrals */
4782: /* PetscCall(DMPlexComputeBdJacobian_Internal(dm, X, X_t, t, X_tShift, Jac, JacP, ctx)); */
4783: /* Assemble matrix */
4784: if (assembleJac) {
4785: PetscCall(MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY));
4786: PetscCall(MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY));
4787: }
4788: PetscCall(MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY));
4789: PetscCall(MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY));
4790: PetscCall(PetscLogEventEnd(DMPLEX_JacobianFEM, dm, 0, 0, 0));
4791: PetscFunctionReturn(PETSC_SUCCESS);
4792: }
4794: /* FEM Assembly Function */
4796: static PetscErrorCode DMConvertPlex_Internal(DM dm, DM *plex, PetscBool copy)
4797: {
4798: PetscBool isPlex;
4800: PetscFunctionBegin;
4801: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
4802: if (isPlex) {
4803: *plex = dm;
4804: PetscCall(PetscObjectReference((PetscObject)dm));
4805: } else {
4806: PetscCall(PetscObjectQuery((PetscObject)dm, "dm_plex", (PetscObject *)plex));
4807: if (!*plex) {
4808: PetscCall(DMConvert(dm, DMPLEX, plex));
4809: PetscCall(PetscObjectCompose((PetscObject)dm, "dm_plex", (PetscObject)*plex));
4810: } else {
4811: PetscCall(PetscObjectReference((PetscObject)*plex));
4812: }
4813: if (copy) PetscCall(DMCopyAuxiliaryVec(dm, *plex));
4814: }
4815: PetscFunctionReturn(PETSC_SUCCESS);
4816: }
4818: /*@
4819: DMPlexGetGeometryFVM - Return precomputed geometric data
4821: Collective
4823: Input Parameter:
4824: . dm - The `DM`
4826: Output Parameters:
4827: + facegeom - The values precomputed from face geometry
4828: . cellgeom - The values precomputed from cell geometry
4829: - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell, or `NULL` if not needed
4831: Level: developer
4833: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMTSSetRHSFunctionLocal()`
4834: @*/
4835: PetscErrorCode DMPlexGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PeOp PetscReal *minRadius)
4836: {
4837: DM plex;
4839: PetscFunctionBegin;
4841: PetscCall(DMConvertPlex_Internal(dm, &plex, PETSC_TRUE));
4842: PetscCall(DMPlexGetDataFVM(plex, NULL, cellgeom, facegeom, NULL));
4843: if (minRadius) PetscCall(DMPlexGetMinRadius(plex, minRadius));
4844: PetscCall(DMDestroy(&plex));
4845: PetscFunctionReturn(PETSC_SUCCESS);
4846: }
4848: /*@
4849: DMPlexGetGradientDM - Return gradient data layout
4851: Collective
4853: Input Parameters:
4854: + dm - The `DM`
4855: - fv - The `PetscFV`
4857: Output Parameter:
4858: . dmGrad - The layout for gradient values
4860: Level: developer
4862: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetGeometryFVM()`
4863: @*/
4864: PetscErrorCode DMPlexGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
4865: {
4866: DM plex;
4867: PetscBool computeGradients;
4869: PetscFunctionBegin;
4872: PetscAssertPointer(dmGrad, 3);
4873: PetscCall(PetscFVGetComputeGradients(fv, &computeGradients));
4874: if (!computeGradients) {
4875: *dmGrad = NULL;
4876: PetscFunctionReturn(PETSC_SUCCESS);
4877: }
4878: PetscCall(DMConvertPlex_Internal(dm, &plex, PETSC_TRUE));
4879: PetscCall(DMPlexGetDataFVM(plex, fv, NULL, NULL, dmGrad));
4880: PetscCall(DMDestroy(&plex));
4881: PetscFunctionReturn(PETSC_SUCCESS);
4882: }
4884: /*@
4885: DMPlexComputeBdResidualSingleByKey - Compute the local boundary residual for terms matching the input key
4887: Not collective
4889: Input Parameters:
4890: + dm - The output `DM`
4891: . wf - The `PetscWeakForm` holding forms on this boundary
4892: . key - The `PetscFormKey` indicating what should be integrated
4893: . facetIS - The `IS` giving a set of faces to integrate over
4894: . locX - The local solution
4895: . locX_t - The time derivative of the local solution, or `NULL` for time-independent problems
4896: . t - The time
4897: - coordField - The `DMField` object with coordinates for these faces
4899: Output Parameter:
4900: . locF - The local residual
4902: Level: developer
4904: .seealso: `DMPlexComputeBdResidualSingle()`, `DMPlexComputeJacobianByKey()`, `DMPlexComputeResidualHybridByKey()`, `DMPlexComputeJacobianHybridByKey()`, `PetscFormKey`
4905: @*/
4906: PetscErrorCode DMPlexComputeBdResidualSingleByKey(DM dm, PetscWeakForm wf, PetscFormKey key, IS facetIS, Vec locX, Vec locX_t, PetscReal t, DMField coordField, Vec locF)
4907: {
4908: DM_Plex *mesh = (DM_Plex *)dm->data;
4909: DM plex = NULL, plexA = NULL;
4910: const char *name = "BdResidual";
4911: DMEnclosureType encAux;
4912: PetscDS prob, probAux = NULL;
4913: PetscSection section, sectionAux = NULL;
4914: Vec locA = NULL;
4915: PetscScalar *u = NULL, *u_t = NULL, *a = NULL, *elemVec = NULL;
4916: PetscInt totDim, totDimAux = 0;
4918: PetscFunctionBegin;
4919: PetscCall(DMConvert(dm, DMPLEX, &plex));
4920: PetscCall(DMGetLocalSection(dm, §ion));
4921: PetscCall(DMGetDS(dm, &prob));
4922: PetscCall(PetscDSGetTotalDimension(prob, &totDim));
4923: PetscCall(DMGetAuxiliaryVec(dm, key.label, key.value, key.part, &locA));
4924: if (locA) {
4925: DM dmAux;
4927: PetscCall(VecGetDM(locA, &dmAux));
4928: PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
4929: PetscCall(DMConvert(dmAux, DMPLEX, &plexA));
4930: PetscCall(DMGetDS(plexA, &probAux));
4931: PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
4932: PetscCall(DMGetLocalSection(plexA, §ionAux));
4933: }
4934: {
4935: PetscFEGeom *fgeom;
4936: PetscInt maxDegree;
4937: PetscQuadrature qGeom = NULL;
4938: IS pointIS;
4939: const PetscInt *points;
4940: PetscInt numFaces, face, Nq;
4942: PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
4943: if (!pointIS) goto end; /* No points with that id on this process */
4944: {
4945: IS isectIS;
4947: /* TODO: Special cases of ISIntersect where it is quick to check a priori if one is a superset of the other */
4948: PetscCall(ISIntersect_Caching_Internal(facetIS, pointIS, &isectIS));
4949: PetscCall(ISDestroy(&pointIS));
4950: pointIS = isectIS;
4951: }
4952: PetscCall(ISGetLocalSize(pointIS, &numFaces));
4953: PetscCall(ISGetIndices(pointIS, &points));
4954: PetscCall(PetscMalloc4(numFaces * totDim, &u, (locX_t ? (size_t)numFaces * totDim : 0), &u_t, numFaces * totDim, &elemVec, (locA ? (size_t)numFaces * totDimAux : 0), &a));
4955: PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
4956: if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &qGeom));
4957: if (!qGeom) {
4958: PetscFE fe;
4960: PetscCall(PetscDSGetDiscretization(prob, key.field, (PetscObject *)&fe));
4961: PetscCall(PetscFEGetFaceQuadrature(fe, &qGeom));
4962: PetscCall(PetscObjectReference((PetscObject)qGeom));
4963: }
4964: PetscCall(PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL));
4965: PetscCall(DMSNESGetFEGeom(coordField, pointIS, qGeom, PETSC_FEGEOM_BOUNDARY, &fgeom));
4966: for (face = 0; face < numFaces; ++face) {
4967: const PetscInt point = points[face], *support;
4968: PetscScalar *x = NULL;
4969: PetscInt i;
4971: PetscCall(DMPlexGetSupport(dm, point, &support));
4972: PetscCall(DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x));
4973: for (i = 0; i < totDim; ++i) u[face * totDim + i] = x[i];
4974: PetscCall(DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x));
4975: if (locX_t) {
4976: PetscCall(DMPlexVecGetClosure(plex, section, locX_t, support[0], NULL, &x));
4977: for (i = 0; i < totDim; ++i) u_t[face * totDim + i] = x[i];
4978: PetscCall(DMPlexVecRestoreClosure(plex, section, locX_t, support[0], NULL, &x));
4979: }
4980: if (locA) {
4981: PetscInt subp;
4983: PetscCall(DMGetEnclosurePoint(plexA, dm, encAux, support[0], &subp));
4984: PetscCall(DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x));
4985: for (i = 0; i < totDimAux; ++i) a[face * totDimAux + i] = x[i];
4986: PetscCall(DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x));
4987: }
4988: }
4989: PetscCall(PetscArrayzero(elemVec, numFaces * totDim));
4990: {
4991: PetscFE fe;
4992: PetscInt Nb;
4993: PetscFEGeom *chunkGeom = NULL;
4994: /* Conforming batches */
4995: PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
4996: /* Remainder */
4997: PetscInt Nr, offset;
4999: PetscCall(PetscDSGetDiscretization(prob, key.field, (PetscObject *)&fe));
5000: PetscCall(PetscFEGetDimension(fe, &Nb));
5001: PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
5002: /* TODO: documentation is unclear about what is going on with these numbers: how should Nb / Nq factor in ? */
5003: blockSize = Nb;
5004: batchSize = numBlocks * blockSize;
5005: PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
5006: numChunks = numFaces / (numBatches * batchSize);
5007: Ne = numChunks * numBatches * batchSize;
5008: Nr = numFaces % (numBatches * batchSize);
5009: offset = numFaces - Nr;
5010: PetscCall(PetscFEGeomGetChunk(fgeom, 0, offset, &chunkGeom));
5011: PetscCall(PetscFEIntegrateBdResidual(prob, wf, key, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec));
5012: PetscCall(PetscFEGeomRestoreChunk(fgeom, 0, offset, &chunkGeom));
5013: PetscCall(PetscFEGeomGetChunk(fgeom, offset, numFaces, &chunkGeom));
5014: PetscCall(PetscFEIntegrateBdResidual(prob, wf, key, Nr, chunkGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), probAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, &elemVec[offset * totDim]));
5015: PetscCall(PetscFEGeomRestoreChunk(fgeom, offset, numFaces, &chunkGeom));
5016: }
5017: for (face = 0; face < numFaces; ++face) {
5018: const PetscInt point = points[face], *support;
5020: if (mesh->printFEM > 1) PetscCall(DMPrintCellVector(point, name, totDim, &elemVec[face * totDim]));
5021: PetscCall(DMPlexGetSupport(plex, point, &support));
5022: PetscCall(DMPlexVecSetClosure(plex, NULL, locF, support[0], &elemVec[face * totDim], ADD_ALL_VALUES));
5023: }
5024: PetscCall(DMSNESRestoreFEGeom(coordField, pointIS, qGeom, PETSC_TRUE, &fgeom));
5025: PetscCall(PetscQuadratureDestroy(&qGeom));
5026: PetscCall(ISRestoreIndices(pointIS, &points));
5027: PetscCall(ISDestroy(&pointIS));
5028: PetscCall(PetscFree4(u, u_t, elemVec, a));
5029: }
5030: end:
5031: if (mesh->printFEM) {
5032: PetscSection s;
5033: Vec locFbc;
5034: PetscInt pStart, pEnd, maxDof;
5035: PetscScalar *zeroes;
5037: PetscCall(DMGetLocalSection(dm, &s));
5038: PetscCall(VecDuplicate(locF, &locFbc));
5039: PetscCall(VecCopy(locF, locFbc));
5040: PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
5041: PetscCall(PetscSectionGetMaxDof(s, &maxDof));
5042: PetscCall(PetscCalloc1(maxDof, &zeroes));
5043: for (PetscInt p = pStart; p < pEnd; p++) PetscCall(VecSetValuesSection(locFbc, s, p, zeroes, INSERT_BC_VALUES));
5044: PetscCall(PetscFree(zeroes));
5045: PetscCall(DMPrintLocalVec(dm, name, mesh->printTol, locFbc));
5046: PetscCall(VecDestroy(&locFbc));
5047: }
5048: PetscCall(DMDestroy(&plex));
5049: PetscCall(DMDestroy(&plexA));
5050: PetscFunctionReturn(PETSC_SUCCESS);
5051: }
5053: /*@
5054: DMPlexComputeBdResidualSingle - Compute the local boundary residual
5056: Not collective
5058: Input Parameters:
5059: + dm - The output `DM`
5060: . wf - The `PetscWeakForm` holding forms on this boundary
5061: . key - The `PetscFormKey` indicating what should be integrated
5062: . locX - The local solution
5063: . locX_t - The time derivative of the local solution, or `NULL` for time-independent problems
5064: - t - The time
5066: Output Parameter:
5067: . locF - The local residual
5069: Level: developer
5071: .seealso: `DMPlexComputeBdResidualSingleByKey()`, `DMPlexComputeJacobianByKey()`, `DMPlexComputeResidualHybridByKey()`, `DMPlexComputeJacobianHybridByKey()`, `PetscFormKey`
5072: @*/
5073: PetscErrorCode DMPlexComputeBdResidualSingle(DM dm, PetscWeakForm wf, PetscFormKey key, Vec locX, Vec locX_t, PetscReal t, Vec locF)
5074: {
5075: DMField coordField;
5076: DMLabel depthLabel;
5077: IS facetIS;
5078: PetscInt dim;
5080: PetscFunctionBegin;
5081: PetscCall(DMGetDimension(dm, &dim));
5082: PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
5083: PetscCall(DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS));
5084: PetscCall(DMGetCoordinateField(dm, &coordField));
5085: PetscCall(DMPlexComputeBdResidualSingleByKey(dm, wf, key, facetIS, locX, locX_t, t, coordField, locF));
5086: PetscCall(ISDestroy(&facetIS));
5087: PetscFunctionReturn(PETSC_SUCCESS);
5088: }
5090: static PetscErrorCode DMPlexComputeBdResidual_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
5091: {
5092: PetscDS prob;
5093: PetscInt numBd, bd;
5094: DMField coordField = NULL;
5095: IS facetIS = NULL;
5096: DMLabel depthLabel;
5097: PetscInt dim;
5099: PetscFunctionBegin;
5100: PetscCall(DMGetDS(dm, &prob));
5101: PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
5102: PetscCall(DMGetDimension(dm, &dim));
5103: PetscCall(DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS));
5104: PetscCall(PetscDSGetNumBoundary(prob, &numBd));
5105: for (bd = 0; bd < numBd; ++bd) {
5106: PetscWeakForm wf;
5107: DMBoundaryConditionType type;
5108: DMLabel label;
5109: const PetscInt *values;
5110: PetscInt field, numValues, v;
5111: PetscObject obj;
5112: PetscClassId id;
5113: PetscFormKey key;
5115: PetscCall(PetscDSGetBoundary(prob, bd, &wf, &type, NULL, &label, &numValues, &values, &field, NULL, NULL, NULL, NULL, NULL));
5116: if (type & DM_BC_ESSENTIAL) continue;
5117: PetscCall(PetscDSGetDiscretization(prob, field, &obj));
5118: PetscCall(PetscObjectGetClassId(obj, &id));
5119: if (id != PETSCFE_CLASSID) continue;
5120: if (!facetIS) {
5121: DMLabel depthLabel;
5122: PetscInt dim;
5124: PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
5125: PetscCall(DMGetDimension(dm, &dim));
5126: PetscCall(DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS));
5127: }
5128: PetscCall(DMGetCoordinateField(dm, &coordField));
5129: for (v = 0; v < numValues; ++v) {
5130: key.label = label;
5131: key.value = values[v];
5132: key.field = field;
5133: key.part = 0;
5134: PetscCall(DMPlexComputeBdResidualSingleByKey(dm, wf, key, facetIS, locX, locX_t, t, coordField, locF));
5135: }
5136: }
5137: PetscCall(ISDestroy(&facetIS));
5138: PetscFunctionReturn(PETSC_SUCCESS);
5139: }
5141: /*@
5142: DMPlexComputeResidualByKey - Compute the local residual for terms matching the input key
5144: Collective
5146: Input Parameters:
5147: + dm - The output `DM`
5148: . key - The `PetscFormKey` indicating what should be integrated
5149: . cellIS - The `IS` giving a set of cells to integrate over
5150: . time - The time, or `PETSC_MIN_REAL` to include implicit terms in a time-independent problems
5151: . locX - The local solution
5152: . locX_t - The time derivative of the local solution, or `NULL` for time-independent problems
5153: . t - The time
5154: - user - An optional user context, passed to the pointwise functions
5156: Output Parameter:
5157: . locF - The local residual
5159: Level: developer
5161: .seealso: `DMPlexComputeJacobianByKey()`, `DMPlexComputeResidualHybridByKey()`, `DMPlexComputeJacobianHybridByKey()`, `PetscFormKey`
5162: @*/
5163: PetscErrorCode DMPlexComputeResidualByKey(DM dm, PetscFormKey key, IS cellIS, PetscReal time, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
5164: {
5165: DM_Plex *mesh = (DM_Plex *)dm->data;
5166: const char *name = "Residual";
5167: DM dmAux = NULL;
5168: DM dmGrad = NULL;
5169: DMLabel ghostLabel = NULL;
5170: PetscDS ds = NULL;
5171: PetscDS dsAux = NULL;
5172: PetscSection section = NULL;
5173: PetscBool useFEM = PETSC_FALSE;
5174: PetscBool useFVM = PETSC_FALSE;
5175: PetscBool isImplicit = (locX_t || time == PETSC_MIN_REAL) ? PETSC_TRUE : PETSC_FALSE;
5176: PetscFV fvm = NULL;
5177: DMField coordField = NULL;
5178: Vec locA, cellGeometryFVM = NULL, faceGeometryFVM = NULL, locGrad = NULL;
5179: PetscScalar *u = NULL, *u_t, *a, *uL, *uR;
5180: IS chunkIS;
5181: const PetscInt *cells;
5182: PetscInt cStart, cEnd, numCells;
5183: PetscInt Nf, f, totDim, totDimAux, numChunks, cellChunkSize, faceChunkSize, chunk, fStart, fEnd;
5184: PetscInt maxDegree = PETSC_INT_MAX;
5185: PetscQuadrature affineQuad = NULL, *quads = NULL;
5186: PetscFEGeom *affineGeom = NULL, **geoms = NULL;
5188: PetscFunctionBegin;
5189: PetscCall(PetscLogEventBegin(DMPLEX_ResidualFEM, dm, 0, 0, 0));
5190: if (!cellIS) goto end;
5191: PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
5192: if (cStart >= cEnd) goto end;
5193: /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */
5194: /* TODO The FVM geometry is over-manipulated. Make the precalc functions return exactly what we need */
5195: /* FEM+FVM */
5196: PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
5197: /* 1: Get sizes from dm and dmAux */
5198: PetscCall(DMGetLocalSection(dm, §ion));
5199: PetscCall(DMGetLabel(dm, "ghost", &ghostLabel));
5200: PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &ds, NULL));
5201: PetscCall(PetscDSGetNumFields(ds, &Nf));
5202: PetscCall(PetscDSGetTotalDimension(ds, &totDim));
5203: PetscCall(DMGetAuxiliaryVec(dm, key.label, key.value, key.part, &locA));
5204: if (locA) {
5205: PetscInt subcell;
5206: PetscCall(VecGetDM(locA, &dmAux));
5207: PetscCall(DMGetEnclosurePoint(dmAux, dm, DM_ENC_UNKNOWN, cells ? cells[cStart] : cStart, &subcell));
5208: PetscCall(DMGetCellDS(dmAux, subcell, &dsAux, NULL));
5209: PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
5210: }
5211: /* 2: Get geometric data */
5212: for (f = 0; f < Nf; ++f) {
5213: PetscObject obj;
5214: PetscClassId id;
5215: PetscBool fimp;
5217: PetscCall(PetscDSGetImplicit(ds, f, &fimp));
5218: if (isImplicit != fimp) continue;
5219: PetscCall(PetscDSGetDiscretization(ds, f, &obj));
5220: PetscCall(PetscObjectGetClassId(obj, &id));
5221: if (id == PETSCFE_CLASSID) useFEM = PETSC_TRUE;
5222: if (id == PETSCFV_CLASSID) {
5223: useFVM = PETSC_TRUE;
5224: fvm = (PetscFV)obj;
5225: }
5226: }
5227: if (useFEM) {
5228: PetscCall(DMGetCoordinateField(dm, &coordField));
5229: PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
5230: if (maxDegree <= 1) {
5231: PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, &affineQuad));
5232: if (affineQuad) PetscCall(DMSNESGetFEGeom(coordField, cellIS, affineQuad, PETSC_FEGEOM_BASIC, &affineGeom));
5233: } else {
5234: PetscCall(PetscCalloc2(Nf, &quads, Nf, &geoms));
5235: for (f = 0; f < Nf; ++f) {
5236: PetscObject obj;
5237: PetscClassId id;
5238: PetscBool fimp;
5240: PetscCall(PetscDSGetImplicit(ds, f, &fimp));
5241: if (isImplicit != fimp) continue;
5242: PetscCall(PetscDSGetDiscretization(ds, f, &obj));
5243: PetscCall(PetscObjectGetClassId(obj, &id));
5244: if (id == PETSCFE_CLASSID) {
5245: PetscFE fe = (PetscFE)obj;
5247: PetscCall(PetscFEGetQuadrature(fe, &quads[f]));
5248: PetscCall(PetscObjectReference((PetscObject)quads[f]));
5249: PetscCall(DMSNESGetFEGeom(coordField, cellIS, quads[f], PETSC_FEGEOM_BASIC, &geoms[f]));
5250: }
5251: }
5252: }
5253: }
5254: // Handle non-essential (e.g. outflow) boundary values
5255: if (useFVM) {
5256: PetscCall(DMPlexInsertBoundaryValuesFVM(dm, fvm, locX, time, &locGrad));
5257: PetscCall(DMPlexGetGeometryFVM(dm, &faceGeometryFVM, &cellGeometryFVM, NULL));
5258: PetscCall(DMPlexGetGradientDM(dm, fvm, &dmGrad));
5259: }
5260: /* Loop over chunks */
5261: if (useFEM) PetscCall(ISCreate(PETSC_COMM_SELF, &chunkIS));
5262: numCells = cEnd - cStart;
5263: numChunks = 1;
5264: cellChunkSize = numCells / numChunks;
5265: faceChunkSize = (fEnd - fStart) / numChunks;
5266: numChunks = PetscMin(1, numCells);
5267: for (chunk = 0; chunk < numChunks; ++chunk) {
5268: PetscScalar *elemVec, *fluxL, *fluxR;
5269: PetscReal *vol;
5270: PetscFVFaceGeom *fgeom;
5271: PetscInt cS = cStart + chunk * cellChunkSize, cE = PetscMin(cS + cellChunkSize, cEnd), numCells = cE - cS, c;
5272: PetscInt fS = fStart + chunk * faceChunkSize, fE = PetscMin(fS + faceChunkSize, fEnd), numFaces = 0, face;
5274: /* Extract field coefficients */
5275: if (useFEM) {
5276: PetscCall(ISGetPointSubrange(chunkIS, cS, cE, cells));
5277: PetscCall(DMPlexGetCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a));
5278: PetscCall(DMGetWorkArray(dm, numCells * totDim, MPIU_SCALAR, &elemVec));
5279: PetscCall(PetscArrayzero(elemVec, numCells * totDim));
5280: }
5281: if (useFVM) {
5282: PetscCall(DMPlexGetFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &numFaces, &uL, &uR));
5283: PetscCall(DMPlexGetFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &numFaces, &fgeom, &vol));
5284: PetscCall(DMGetWorkArray(dm, numFaces * totDim, MPIU_SCALAR, &fluxL));
5285: PetscCall(DMGetWorkArray(dm, numFaces * totDim, MPIU_SCALAR, &fluxR));
5286: PetscCall(PetscArrayzero(fluxL, numFaces * totDim));
5287: PetscCall(PetscArrayzero(fluxR, numFaces * totDim));
5288: }
5289: /* TODO We will interlace both our field coefficients (u, u_t, uL, uR, etc.) and our output (elemVec, fL, fR). I think this works */
5290: /* Loop over fields */
5291: for (f = 0; f < Nf; ++f) {
5292: PetscObject obj;
5293: PetscClassId id;
5294: PetscBool fimp;
5295: PetscInt numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;
5297: key.field = f;
5298: PetscCall(PetscDSGetImplicit(ds, f, &fimp));
5299: if (isImplicit != fimp) continue;
5300: PetscCall(PetscDSGetDiscretization(ds, f, &obj));
5301: PetscCall(PetscObjectGetClassId(obj, &id));
5302: if (id == PETSCFE_CLASSID) {
5303: PetscFE fe = (PetscFE)obj;
5304: PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f];
5305: PetscFEGeom *chunkGeom = NULL;
5306: PetscQuadrature quad = affineQuad ? affineQuad : quads[f];
5307: PetscInt Nq, Nb;
5309: PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
5310: PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
5311: PetscCall(PetscFEGetDimension(fe, &Nb));
5312: blockSize = Nb;
5313: batchSize = numBlocks * blockSize;
5314: PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
5315: numChunks = numCells / (numBatches * batchSize);
5316: Ne = numChunks * numBatches * batchSize;
5317: Nr = numCells % (numBatches * batchSize);
5318: offset = numCells - Nr;
5319: /* Integrate FE residual to get elemVec (need fields at quadrature points) */
5320: /* 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) */
5321: PetscCall(PetscFEGeomGetChunk(geom, 0, offset, &chunkGeom));
5322: PetscCall(PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, u_t, dsAux, a, t, elemVec));
5323: PetscCall(PetscFEGeomGetChunk(geom, offset, numCells, &chunkGeom));
5324: PetscCall(PetscFEIntegrateResidual(ds, key, Nr, chunkGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), dsAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, &elemVec[offset * totDim]));
5325: PetscCall(PetscFEGeomRestoreChunk(geom, offset, numCells, &chunkGeom));
5326: } else if (id == PETSCFV_CLASSID) {
5327: PetscFV fv = (PetscFV)obj;
5329: Ne = numFaces;
5330: /* Riemann solve over faces (need fields at face centroids) */
5331: /* We need to evaluate FE fields at those coordinates */
5332: PetscCall(PetscFVIntegrateRHSFunction(fv, ds, f, Ne, fgeom, vol, uL, uR, fluxL, fluxR));
5333: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f);
5334: }
5335: /* Loop over domain */
5336: if (useFEM) {
5337: /* Add elemVec to locX */
5338: for (c = cS; c < cE; ++c) {
5339: const PetscInt cell = cells ? cells[c] : c;
5340: const PetscInt cind = c - cStart;
5342: if (mesh->printFEM > 1) PetscCall(DMPrintCellVector(cell, name, totDim, &elemVec[cind * totDim]));
5343: if (ghostLabel) {
5344: PetscInt ghostVal;
5346: PetscCall(DMLabelGetValue(ghostLabel, cell, &ghostVal));
5347: if (ghostVal > 0) continue;
5348: }
5349: PetscCall(DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cind * totDim], ADD_ALL_VALUES));
5350: }
5351: }
5352: if (useFVM) {
5353: PetscScalar *fa;
5354: PetscInt iface;
5356: PetscCall(VecGetArray(locF, &fa));
5357: for (f = 0; f < Nf; ++f) {
5358: PetscFV fv;
5359: PetscObject obj;
5360: PetscClassId id;
5361: PetscInt cdim, foff, pdim;
5363: PetscCall(DMGetCoordinateDim(dm, &cdim));
5364: PetscCall(PetscDSGetDiscretization(ds, f, &obj));
5365: PetscCall(PetscDSGetFieldOffset(ds, f, &foff));
5366: PetscCall(PetscObjectGetClassId(obj, &id));
5367: if (id != PETSCFV_CLASSID) continue;
5368: fv = (PetscFV)obj;
5369: PetscCall(PetscFVGetNumComponents(fv, &pdim));
5370: /* Accumulate fluxes to cells */
5371: for (face = fS, iface = 0; face < fE; ++face) {
5372: const PetscInt *scells;
5373: PetscScalar *fL = NULL, *fR = NULL;
5374: PetscInt ghost, d, nsupp, nchild;
5376: PetscCall(DMLabelGetValue(ghostLabel, face, &ghost));
5377: PetscCall(DMPlexGetSupportSize(dm, face, &nsupp));
5378: PetscCall(DMPlexGetTreeChildren(dm, face, &nchild, NULL));
5379: if (ghost >= 0 || nsupp > 2 || nchild > 0) continue;
5380: PetscCall(DMPlexGetSupport(dm, face, &scells));
5381: PetscCall(DMLabelGetValue(ghostLabel, scells[0], &ghost));
5382: if (ghost <= 0) PetscCall(DMPlexPointLocalFieldRef(dm, scells[0], f, fa, &fL));
5383: PetscCall(DMLabelGetValue(ghostLabel, scells[1], &ghost));
5384: if (ghost <= 0) PetscCall(DMPlexPointLocalFieldRef(dm, scells[1], f, fa, &fR));
5385: if (mesh->printFVM > 1) {
5386: PetscCall(DMPrintCellVectorReal(face, "Residual: normal", cdim, fgeom[iface].normal));
5387: PetscCall(DMPrintCellVector(face, "Residual: left state", pdim, &uL[iface * totDim + foff]));
5388: PetscCall(DMPrintCellVector(face, "Residual: right state", pdim, &uR[iface * totDim + foff]));
5389: PetscCall(DMPrintCellVector(face, "Residual: left flux", pdim, &fluxL[iface * totDim + foff]));
5390: PetscCall(DMPrintCellVector(face, "Residual: right flux", pdim, &fluxR[iface * totDim + foff]));
5391: }
5392: for (d = 0; d < pdim; ++d) {
5393: if (fL) fL[d] -= fluxL[iface * totDim + foff + d];
5394: if (fR) fR[d] += fluxR[iface * totDim + foff + d];
5395: }
5396: ++iface;
5397: }
5398: }
5399: PetscCall(VecRestoreArray(locF, &fa));
5400: }
5401: /* Handle time derivative */
5402: if (locX_t) {
5403: PetscScalar *x_t, *fa;
5405: PetscCall(VecGetArray(locF, &fa));
5406: PetscCall(VecGetArray(locX_t, &x_t));
5407: for (f = 0; f < Nf; ++f) {
5408: PetscFV fv;
5409: PetscObject obj;
5410: PetscClassId id;
5411: PetscInt pdim, d;
5413: PetscCall(PetscDSGetDiscretization(ds, f, &obj));
5414: PetscCall(PetscObjectGetClassId(obj, &id));
5415: if (id != PETSCFV_CLASSID) continue;
5416: fv = (PetscFV)obj;
5417: PetscCall(PetscFVGetNumComponents(fv, &pdim));
5418: for (c = cS; c < cE; ++c) {
5419: const PetscInt cell = cells ? cells[c] : c;
5420: PetscScalar *u_t, *r;
5422: if (ghostLabel) {
5423: PetscInt ghostVal;
5425: PetscCall(DMLabelGetValue(ghostLabel, cell, &ghostVal));
5426: if (ghostVal > 0) continue;
5427: }
5428: PetscCall(DMPlexPointLocalFieldRead(dm, cell, f, x_t, &u_t));
5429: PetscCall(DMPlexPointLocalFieldRef(dm, cell, f, fa, &r));
5430: for (d = 0; d < pdim; ++d) r[d] += u_t[d];
5431: }
5432: }
5433: PetscCall(VecRestoreArray(locX_t, &x_t));
5434: PetscCall(VecRestoreArray(locF, &fa));
5435: }
5436: if (useFEM) {
5437: PetscCall(DMPlexRestoreCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a));
5438: PetscCall(DMRestoreWorkArray(dm, numCells * totDim, MPIU_SCALAR, &elemVec));
5439: }
5440: if (useFVM) {
5441: PetscCall(DMPlexRestoreFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &numFaces, &uL, &uR));
5442: PetscCall(DMPlexRestoreFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &numFaces, &fgeom, &vol));
5443: PetscCall(DMRestoreWorkArray(dm, numFaces * totDim, MPIU_SCALAR, &fluxL));
5444: PetscCall(DMRestoreWorkArray(dm, numFaces * totDim, MPIU_SCALAR, &fluxR));
5445: if (dmGrad) PetscCall(DMRestoreLocalVector(dmGrad, &locGrad));
5446: }
5447: }
5448: if (useFEM) PetscCall(ISDestroy(&chunkIS));
5449: PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
5451: if (useFEM) {
5452: PetscCall(DMPlexComputeBdResidual_Internal(dm, locX, locX_t, t, locF, user));
5454: if (maxDegree <= 1) {
5455: PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, affineQuad, PETSC_FALSE, &affineGeom));
5456: PetscCall(PetscQuadratureDestroy(&affineQuad));
5457: } else {
5458: for (f = 0; f < Nf; ++f) {
5459: PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, quads[f], PETSC_FALSE, &geoms[f]));
5460: PetscCall(PetscQuadratureDestroy(&quads[f]));
5461: }
5462: PetscCall(PetscFree2(quads, geoms));
5463: }
5464: }
5466: /* FEM */
5467: /* 1: Get sizes from dm and dmAux */
5468: /* 2: Get geometric data */
5469: /* 3: Handle boundary values */
5470: /* 4: Loop over domain */
5471: /* Extract coefficients */
5472: /* Loop over fields */
5473: /* Set tiling for FE*/
5474: /* Integrate FE residual to get elemVec */
5475: /* Loop over subdomain */
5476: /* Loop over quad points */
5477: /* Transform coords to real space */
5478: /* Evaluate field and aux fields at point */
5479: /* Evaluate residual at point */
5480: /* Transform residual to real space */
5481: /* Add residual to elemVec */
5482: /* Loop over domain */
5483: /* Add elemVec to locX */
5485: /* FVM */
5486: /* Get geometric data */
5487: /* If using gradients */
5488: /* Compute gradient data */
5489: /* Loop over domain faces */
5490: /* Count computational faces */
5491: /* Reconstruct cell gradient */
5492: /* Loop over domain cells */
5493: /* Limit cell gradients */
5494: /* Handle boundary values */
5495: /* Loop over domain faces */
5496: /* Read out field, centroid, normal, volume for each side of face */
5497: /* Riemann solve over faces */
5498: /* Loop over domain faces */
5499: /* Accumulate fluxes to cells */
5500: /* TODO Change printFEM to printDisc here */
5501: if (mesh->printFEM) {
5502: Vec locFbc;
5503: PetscInt pStart, pEnd, p, maxDof;
5504: PetscScalar *zeroes;
5506: PetscCall(VecDuplicate(locF, &locFbc));
5507: PetscCall(VecCopy(locF, locFbc));
5508: PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
5509: PetscCall(PetscSectionGetMaxDof(section, &maxDof));
5510: PetscCall(PetscCalloc1(maxDof, &zeroes));
5511: for (p = pStart; p < pEnd; p++) PetscCall(VecSetValuesSection(locFbc, section, p, zeroes, INSERT_BC_VALUES));
5512: PetscCall(PetscFree(zeroes));
5513: PetscCall(DMPrintLocalVec(dm, name, mesh->printTol, locFbc));
5514: PetscCall(VecDestroy(&locFbc));
5515: }
5516: end:
5517: PetscCall(PetscLogEventEnd(DMPLEX_ResidualFEM, dm, 0, 0, 0));
5518: PetscFunctionReturn(PETSC_SUCCESS);
5519: }
5521: /*@
5522: DMPlexComputeResidualHybridByKey - Compute the local residual over hybrid cells for terms matching the input key
5524: Collective
5526: Input Parameters:
5527: + dm - The output `DM`
5528: . key - The `PetscFormKey` array (left cell, right cell, cohesive cell) indicating what should be integrated
5529: . cellIS - The `IS` give a set of cells to integrate over
5530: . time - The time, or `PETSC_MIN_REAL` to include implicit terms in a time-independent problems
5531: . locX - The local solution
5532: . locX_t - The time derivative of the local solution, or `NULL` for time-independent problems
5533: . t - The time
5534: - user - An optional user context, passed to the pointwise functions
5536: Output Parameter:
5537: . locF - The local residual
5539: Level: developer
5541: .seealso: `DMPlexComputeResidualByKey()`, `DMPlexComputeJacobianByKey()`, `DMPlexComputeJacobianHybridByKey()`, `PetscFormKey`
5542: @*/
5543: PetscErrorCode DMPlexComputeResidualHybridByKey(DM dm, PetscFormKey key[], IS cellIS, PetscReal time, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
5544: {
5545: DM_Plex *mesh = (DM_Plex *)dm->data;
5546: const char *name = "Hybrid Residual";
5547: DM dmAux[3] = {NULL, NULL, NULL};
5548: DMLabel ghostLabel = NULL;
5549: PetscDS ds = NULL;
5550: PetscDS dsIn = NULL;
5551: PetscDS dsAux[3] = {NULL, NULL, NULL};
5552: Vec locA[3] = {NULL, NULL, NULL};
5553: DM dmScale[3] = {NULL, NULL, NULL};
5554: PetscDS dsScale[3] = {NULL, NULL, NULL};
5555: Vec locS[3] = {NULL, NULL, NULL};
5556: PetscSection section = NULL;
5557: DMField coordField = NULL;
5558: PetscScalar *a[3] = {NULL, NULL, NULL};
5559: PetscScalar *s[3] = {NULL, NULL, NULL};
5560: PetscScalar *u = NULL, *u_t;
5561: PetscScalar *elemVecNeg, *elemVecPos, *elemVecCoh;
5562: IS chunkISF, chunkISN;
5563: const PetscInt *cells;
5564: PetscInt *faces, *neighbors;
5565: PetscInt cStart, cEnd, numCells;
5566: PetscInt Nf, f, totDim, totDimIn, totDimAux[3], totDimScale[3], numChunks, cellChunkSize, chunk;
5567: PetscInt maxDegree = PETSC_INT_MAX;
5568: PetscQuadrature affineQuadF = NULL, *quadsF = NULL;
5569: PetscFEGeom *affineGeomF = NULL, **geomsF = NULL;
5570: PetscQuadrature affineQuadN = NULL, *quadsN = NULL;
5571: PetscFEGeom *affineGeomN = NULL, **geomsN = NULL;
5573: PetscFunctionBegin;
5574: PetscCall(PetscLogEventBegin(DMPLEX_ResidualFEM, dm, 0, 0, 0));
5575: if (!cellIS) goto end;
5576: PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
5577: PetscCall(ISGetLocalSize(cellIS, &numCells));
5578: if (cStart >= cEnd) goto end;
5579: if ((key[0].label == key[1].label) && (key[0].value == key[1].value) && (key[0].part == key[1].part)) {
5580: const char *name;
5581: PetscCall(PetscObjectGetName((PetscObject)key[0].label, &name));
5582: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Form keys for each side of a cohesive surface must be different (%s, %" PetscInt_FMT ", %" PetscInt_FMT ")", name, key[0].value, key[0].part);
5583: }
5584: /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */
5585: /* FEM */
5586: /* 1: Get sizes from dm and dmAux */
5587: PetscCall(DMGetLocalSection(dm, §ion));
5588: PetscCall(DMGetLabel(dm, "ghost", &ghostLabel));
5589: PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &ds, &dsIn));
5590: PetscCall(PetscDSGetNumFields(ds, &Nf));
5591: PetscCall(PetscDSGetTotalDimension(ds, &totDim));
5592: PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn));
5593: PetscCall(DMGetAuxiliaryVec(dm, key[2].label, key[2].value, key[2].part, &locA[2]));
5594: if (locA[2]) {
5595: const PetscInt cellStart = cells ? cells[cStart] : cStart;
5597: PetscCall(VecGetDM(locA[2], &dmAux[2]));
5598: PetscCall(DMGetCellDS(dmAux[2], cellStart, &dsAux[2], NULL));
5599: PetscCall(PetscDSGetTotalDimension(dsAux[2], &totDimAux[2]));
5600: {
5601: const PetscInt *cone;
5602: PetscInt c;
5604: PetscCall(DMPlexGetCone(dm, cellStart, &cone));
5605: for (c = 0; c < 2; ++c) {
5606: const PetscInt *support;
5607: PetscInt ssize, s;
5609: PetscCall(DMPlexGetSupport(dm, cone[c], &support));
5610: PetscCall(DMPlexGetSupportSize(dm, cone[c], &ssize));
5611: PetscCheck(ssize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " from cell %" PetscInt_FMT " has support size %" PetscInt_FMT " != 2", cone[c], cellStart, ssize);
5612: if (support[0] == cellStart) s = 1;
5613: else if (support[1] == cellStart) s = 0;
5614: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " does not have cell %" PetscInt_FMT " in its support", cone[c], cellStart);
5615: PetscCall(DMGetAuxiliaryVec(dm, key[c].label, key[c].value, key[c].part, &locA[c]));
5616: PetscCheck(locA[c], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must have auxiliary vector for (%p, %" PetscInt_FMT ", %" PetscInt_FMT ")", (void *)key[c].label, key[c].value, key[c].part);
5617: if (locA[c]) PetscCall(VecGetDM(locA[c], &dmAux[c]));
5618: else dmAux[c] = dmAux[2];
5619: PetscCall(DMGetCellDS(dmAux[c], support[s], &dsAux[c], NULL));
5620: PetscCall(PetscDSGetTotalDimension(dsAux[c], &totDimAux[c]));
5621: }
5622: }
5623: }
5624: /* Handle mass matrix scaling
5625: The field in key[2] is the field to be scaled, and the scaling field is the first in the dsScale */
5626: PetscCall(DMGetAuxiliaryVec(dm, key[2].label, -key[2].value, key[2].part, &locS[2]));
5627: if (locS[2]) {
5628: const PetscInt cellStart = cells ? cells[cStart] : cStart;
5629: PetscInt Nb, Nbs;
5631: PetscCall(VecGetDM(locS[2], &dmScale[2]));
5632: PetscCall(DMGetCellDS(dmScale[2], cellStart, &dsScale[2], NULL));
5633: PetscCall(PetscDSGetTotalDimension(dsScale[2], &totDimScale[2]));
5634: // BRAD: This is not set correctly
5635: key[2].field = 2;
5636: PetscCall(PetscDSGetFieldSize(ds, key[2].field, &Nb));
5637: PetscCall(PetscDSGetFieldSize(dsScale[2], 0, &Nbs));
5638: PetscCheck(Nb == Nbs, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field %" PetscInt_FMT " of size %" PetscInt_FMT " cannot be scaled by field of size %" PetscInt_FMT, key[2].field, Nb, Nbs);
5639: {
5640: const PetscInt *cone;
5641: PetscInt c;
5643: locS[1] = locS[0] = locS[2];
5644: dmScale[1] = dmScale[0] = dmScale[2];
5645: PetscCall(DMPlexGetCone(dm, cellStart, &cone));
5646: for (c = 0; c < 2; ++c) {
5647: const PetscInt *support;
5648: PetscInt ssize, s;
5650: PetscCall(DMPlexGetSupport(dm, cone[c], &support));
5651: PetscCall(DMPlexGetSupportSize(dm, cone[c], &ssize));
5652: PetscCheck(ssize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " from cell %" PetscInt_FMT " has support size %" PetscInt_FMT " != 2", cone[c], cellStart, ssize);
5653: if (support[0] == cellStart) s = 1;
5654: else if (support[1] == cellStart) s = 0;
5655: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " does not have cell %" PetscInt_FMT " in its support", cone[c], cellStart);
5656: PetscCall(DMGetCellDS(dmScale[c], support[s], &dsScale[c], NULL));
5657: PetscCall(PetscDSGetTotalDimension(dsScale[c], &totDimScale[c]));
5658: }
5659: }
5660: }
5661: /* 2: Setup geometric data */
5662: PetscCall(DMGetCoordinateField(dm, &coordField));
5663: PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
5664: if (maxDegree > 1) {
5665: PetscCall(PetscCalloc4(Nf, &quadsF, Nf, &geomsF, Nf, &quadsN, Nf, &geomsN));
5666: for (f = 0; f < Nf; ++f) {
5667: PetscFE fe;
5668: PetscBool isCohesiveField;
5670: PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
5671: if (fe) {
5672: PetscCall(PetscFEGetQuadrature(fe, &quadsF[f]));
5673: PetscCall(PetscObjectReference((PetscObject)quadsF[f]));
5674: }
5675: PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fe));
5676: PetscCall(PetscDSGetCohesive(dsIn, f, &isCohesiveField));
5677: if (fe) {
5678: if (isCohesiveField) {
5679: for (PetscInt g = 0; g < Nf; ++g) {
5680: PetscCall(PetscDSGetDiscretization(dsIn, g, (PetscObject *)&fe));
5681: PetscCall(PetscDSGetCohesive(dsIn, g, &isCohesiveField));
5682: if (!isCohesiveField) break;
5683: }
5684: }
5685: PetscCall(PetscFEGetQuadrature(fe, &quadsN[f]));
5686: PetscCall(PetscObjectReference((PetscObject)quadsN[f]));
5687: }
5688: }
5689: }
5690: /* Loop over chunks */
5691: cellChunkSize = numCells;
5692: numChunks = !numCells ? 0 : PetscCeilReal(((PetscReal)numCells) / cellChunkSize);
5693: PetscCall(PetscCalloc2(2 * cellChunkSize, &faces, 2 * cellChunkSize, &neighbors));
5694: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 2 * cellChunkSize, faces, PETSC_USE_POINTER, &chunkISF));
5695: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 2 * cellChunkSize, neighbors, PETSC_USE_POINTER, &chunkISN));
5696: /* Extract field coefficients */
5697: /* NOTE This needs the end cap faces to have identical orientations */
5698: PetscCall(DMPlexGetHybridCellFields(dm, cellIS, locX, locX_t, locA[2], &u, &u_t, &a[2]));
5699: PetscCall(DMPlexGetHybridFields(dm, dmAux, dsAux, cellIS, locA, PETSC_TRUE, a));
5700: PetscCall(DMPlexGetHybridFields(dm, dmScale, dsScale, cellIS, locS, PETSC_TRUE, s));
5701: PetscCall(DMGetWorkArray(dm, cellChunkSize * totDim, MPIU_SCALAR, &elemVecNeg));
5702: PetscCall(DMGetWorkArray(dm, cellChunkSize * totDim, MPIU_SCALAR, &elemVecPos));
5703: PetscCall(DMGetWorkArray(dm, cellChunkSize * totDim, MPIU_SCALAR, &elemVecCoh));
5704: for (chunk = 0; chunk < numChunks; ++chunk) {
5705: PetscInt cS = cStart + chunk * cellChunkSize, cE = PetscMin(cS + cellChunkSize, cEnd), numCells = cE - cS, c;
5707: PetscCall(PetscArrayzero(elemVecNeg, cellChunkSize * totDim));
5708: PetscCall(PetscArrayzero(elemVecPos, cellChunkSize * totDim));
5709: PetscCall(PetscArrayzero(elemVecCoh, cellChunkSize * totDim));
5710: /* Get faces and neighbors */
5711: for (c = cS; c < cE; ++c) {
5712: const PetscInt cell = cells ? cells[c] : c;
5713: const PetscInt *cone, *support;
5714: PetscCall(DMPlexGetCone(dm, cell, &cone));
5715: faces[(c - cS) * 2 + 0] = cone[0];
5716: faces[(c - cS) * 2 + 1] = cone[1];
5717: PetscCall(DMPlexGetSupport(dm, cone[0], &support));
5718: neighbors[(c - cS) * 2 + 0] = support[0] == cell ? support[1] : support[0];
5719: PetscCall(DMPlexGetSupport(dm, cone[1], &support));
5720: neighbors[(c - cS) * 2 + 1] = support[0] == cell ? support[1] : support[0];
5721: }
5722: PetscCall(ISGeneralSetIndices(chunkISF, 2 * cellChunkSize, faces, PETSC_USE_POINTER));
5723: PetscCall(ISGeneralSetIndices(chunkISN, 2 * cellChunkSize, neighbors, PETSC_USE_POINTER));
5724: /* Get geometric data */
5725: if (maxDegree <= 1) {
5726: if (!affineQuadF) PetscCall(DMFieldCreateDefaultQuadrature(coordField, chunkISF, &affineQuadF));
5727: if (affineQuadF) PetscCall(DMSNESGetFEGeom(coordField, chunkISF, affineQuadF, PETSC_FEGEOM_COHESIVE, &affineGeomF));
5728: if (!affineQuadN) {
5729: PetscInt dim;
5730: PetscCall(PetscQuadratureGetData(affineQuadF, &dim, NULL, NULL, NULL, NULL));
5731: PetscCall(DMFieldCreateDefaultFaceQuadrature(coordField, chunkISN, &affineQuadN));
5732: PetscCall(PetscQuadratureSetData(affineQuadN, dim + 1, PETSC_DECIDE, PETSC_DECIDE, NULL, NULL));
5733: }
5734: if (affineQuadN) PetscCall(DMSNESGetFEGeom(coordField, chunkISN, affineQuadN, PETSC_FEGEOM_BASIC, &affineGeomN));
5735: } else {
5736: for (f = 0; f < Nf; ++f) {
5737: if (quadsF[f]) PetscCall(DMSNESGetFEGeom(coordField, chunkISF, quadsF[f], PETSC_FEGEOM_COHESIVE, &geomsF[f]));
5738: if (quadsN[f]) PetscCall(DMSNESGetFEGeom(coordField, chunkISN, quadsN[f], PETSC_FEGEOM_BASIC, &geomsN[f]));
5739: }
5740: }
5741: /* Loop over fields */
5742: for (f = 0; f < Nf; ++f) {
5743: PetscFE fe;
5744: PetscFEGeom *geomF = affineGeomF ? affineGeomF : geomsF[f];
5745: PetscFEGeom *chunkGeomF = NULL, *remGeomF = NULL;
5746: PetscFEGeom *geomN = affineGeomN ? affineGeomN : geomsN[f];
5747: PetscFEGeom *chunkGeomN = NULL, *remGeomN = NULL;
5748: PetscQuadrature quadF = affineQuadF ? affineQuadF : quadsF[f];
5749: PetscInt numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset, Nq, Nb;
5750: PetscBool isCohesiveField;
5752: PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
5753: if (!fe) continue;
5754: PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
5755: PetscCall(PetscQuadratureGetData(quadF, NULL, NULL, &Nq, NULL, NULL));
5756: PetscCall(PetscFEGetDimension(fe, &Nb));
5757: blockSize = Nb;
5758: batchSize = numBlocks * blockSize;
5759: PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
5760: numChunks = numCells / (numBatches * batchSize);
5761: Ne = numChunks * numBatches * batchSize;
5762: Nr = numCells % (numBatches * batchSize);
5763: offset = numCells - Nr;
5764: PetscCall(PetscFEGeomGetChunk(geomF, 0, offset * 2, &chunkGeomF));
5765: PetscCall(PetscFEGeomGetChunk(geomF, offset * 2, numCells * 2, &remGeomF));
5766: PetscCall(PetscFEGeomGetChunk(geomN, 0, offset * 2, &chunkGeomN));
5767: PetscCall(PetscFEGeomGetChunk(geomN, offset * 2, numCells * 2, &remGeomN));
5768: PetscCall(PetscDSGetCohesive(ds, f, &isCohesiveField));
5769: // TODO Do I need to set isCohesive on the chunks?
5770: key[0].field = f;
5771: key[1].field = f;
5772: key[2].field = f;
5773: PetscCall(PetscFEIntegrateHybridResidual(ds, dsIn, key[0], 0, Ne, chunkGeomF, chunkGeomN, u, u_t, dsAux[0], a[0], t, elemVecNeg));
5774: PetscCall(PetscFEIntegrateHybridResidual(ds, dsIn, key[0], 0, Nr, remGeomF, remGeomN, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[0], PetscSafePointerPlusOffset(a[0], offset * totDimAux[0]), t, &elemVecNeg[offset * totDim]));
5775: PetscCall(PetscFEIntegrateHybridResidual(ds, dsIn, key[1], 1, Ne, chunkGeomF, chunkGeomN, u, u_t, dsAux[1], a[1], t, elemVecPos));
5776: PetscCall(PetscFEIntegrateHybridResidual(ds, dsIn, key[1], 1, Nr, remGeomF, remGeomN, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[1], PetscSafePointerPlusOffset(a[1], offset * totDimAux[1]), t, &elemVecPos[offset * totDim]));
5777: PetscCall(PetscFEIntegrateHybridResidual(ds, dsIn, key[2], 2, Ne, chunkGeomF, chunkGeomN, u, u_t, dsAux[2], a[2], t, elemVecCoh));
5778: PetscCall(PetscFEIntegrateHybridResidual(ds, dsIn, key[2], 2, Nr, remGeomF, remGeomN, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[2], PetscSafePointerPlusOffset(a[2], offset * totDimAux[2]), t, &elemVecCoh[offset * totDim]));
5779: PetscCall(PetscFEGeomRestoreChunk(geomF, offset, numCells, &remGeomF));
5780: PetscCall(PetscFEGeomRestoreChunk(geomF, 0, offset, &chunkGeomF));
5781: PetscCall(PetscFEGeomRestoreChunk(geomN, offset, numCells, &remGeomN));
5782: PetscCall(PetscFEGeomRestoreChunk(geomN, 0, offset, &chunkGeomN));
5783: }
5784: /* Add elemVec to locX */
5785: for (c = cS; c < cE; ++c) {
5786: const PetscInt cell = cells ? cells[c] : c;
5787: const PetscInt cind = c - cStart;
5788: PetscInt i;
5790: /* Scale element values */
5791: if (locS[0]) {
5792: PetscInt Nb, off = cind * totDim, soff = cind * totDimScale[0];
5793: PetscBool cohesive;
5795: for (f = 0; f < Nf; ++f) {
5796: PetscCall(PetscDSGetFieldSize(ds, f, &Nb));
5797: PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
5798: if (f == key[2].field) {
5799: PetscCheck(cohesive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Scaling should not happen for face fields");
5800: // No cohesive scaling field is currently input
5801: for (i = 0; i < Nb; ++i) elemVecCoh[off + i] += s[0][soff + i] * elemVecNeg[off + i] + s[1][soff + i] * elemVecPos[off + i];
5802: off += Nb;
5803: } else {
5804: const PetscInt N = cohesive ? Nb : Nb * 2;
5806: for (i = 0; i < N; ++i) elemVecCoh[off + i] += elemVecNeg[off + i] + elemVecPos[off + i];
5807: off += N;
5808: }
5809: }
5810: } else {
5811: for (i = cind * totDim; i < (cind + 1) * totDim; ++i) elemVecCoh[i] += elemVecNeg[i] + elemVecPos[i];
5812: }
5813: if (mesh->printFEM > 1) PetscCall(DMPrintCellVector(cell, name, totDim, &elemVecCoh[cind * totDim]));
5814: if (ghostLabel) {
5815: PetscInt ghostVal;
5817: PetscCall(DMLabelGetValue(ghostLabel, cell, &ghostVal));
5818: if (ghostVal > 0) continue;
5819: }
5820: PetscCall(DMPlexVecSetClosure(dm, section, locF, cell, &elemVecCoh[cind * totDim], ADD_ALL_VALUES));
5821: }
5822: }
5823: PetscCall(DMPlexRestoreCellFields(dm, cellIS, locX, locX_t, locA[2], &u, &u_t, &a[2]));
5824: PetscCall(DMPlexRestoreHybridFields(dm, dmAux, dsAux, cellIS, locA, PETSC_TRUE, a));
5825: PetscCall(DMPlexRestoreHybridFields(dm, dmScale, dsScale, cellIS, locS, PETSC_TRUE, s));
5826: PetscCall(DMRestoreWorkArray(dm, numCells * totDim, MPIU_SCALAR, &elemVecNeg));
5827: PetscCall(DMRestoreWorkArray(dm, numCells * totDim, MPIU_SCALAR, &elemVecPos));
5828: PetscCall(DMRestoreWorkArray(dm, numCells * totDim, MPIU_SCALAR, &elemVecCoh));
5829: PetscCall(PetscFree2(faces, neighbors));
5830: PetscCall(ISDestroy(&chunkISF));
5831: PetscCall(ISDestroy(&chunkISN));
5832: PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
5833: if (maxDegree <= 1) {
5834: PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, affineQuadF, PETSC_FALSE, &affineGeomF));
5835: PetscCall(PetscQuadratureDestroy(&affineQuadF));
5836: PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, affineQuadN, PETSC_FALSE, &affineGeomN));
5837: PetscCall(PetscQuadratureDestroy(&affineQuadN));
5838: } else {
5839: for (f = 0; f < Nf; ++f) {
5840: if (geomsF) PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, quadsF[f], PETSC_FALSE, &geomsF[f]));
5841: if (quadsF) PetscCall(PetscQuadratureDestroy(&quadsF[f]));
5842: if (geomsN) PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, quadsN[f], PETSC_FALSE, &geomsN[f]));
5843: if (quadsN) PetscCall(PetscQuadratureDestroy(&quadsN[f]));
5844: }
5845: PetscCall(PetscFree4(quadsF, geomsF, quadsN, geomsN));
5846: }
5847: if (mesh->printFEM) {
5848: Vec locFbc;
5849: PetscInt pStart, pEnd, p, maxDof;
5850: PetscScalar *zeroes;
5852: PetscCall(VecDuplicate(locF, &locFbc));
5853: PetscCall(VecCopy(locF, locFbc));
5854: PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
5855: PetscCall(PetscSectionGetMaxDof(section, &maxDof));
5856: PetscCall(PetscCalloc1(maxDof, &zeroes));
5857: for (p = pStart; p < pEnd; p++) PetscCall(VecSetValuesSection(locFbc, section, p, zeroes, INSERT_BC_VALUES));
5858: PetscCall(PetscFree(zeroes));
5859: PetscCall(DMPrintLocalVec(dm, name, mesh->printTol, locFbc));
5860: PetscCall(VecDestroy(&locFbc));
5861: }
5862: end:
5863: PetscCall(PetscLogEventEnd(DMPLEX_ResidualFEM, dm, 0, 0, 0));
5864: PetscFunctionReturn(PETSC_SUCCESS);
5865: }
5867: /*@
5868: DMPlexComputeBdJacobianSingleByLabel - Compute the local boundary Jacobian for terms matching the input label
5870: Not collective
5872: Input Parameters:
5873: + dm - The output `DM`
5874: . wf - The `PetscWeakForm` holding forms on this boundary
5875: . label - The `DMLabel` indicating what faces should be integrated over
5876: . numValues - The number of label values
5877: . values - The array of label values
5878: . fieldI - The test field for these integrals
5879: . facetIS - The `IS` giving the set of possible faces to integrate over (intersected with the label)
5880: . locX - The local solution
5881: . locX_t - The time derivative of the local solution, or `NULL` for time-independent problems
5882: . t - The time
5883: . coordField - The `DMField` object with coordinates for these faces
5884: - X_tShift - The multiplier for dF/dxdot
5886: Output Parameters:
5887: + Jac - The local Jacobian
5888: - JacP - The local Jacobian preconditioner
5890: Level: developer
5892: .seealso: `DMPlexComputeBdJacobianSingle()`, `DMPlexComputeJacobianByKey()`, `DMPlexComputeResidualHybridByKey()`, `DMPlexComputeJacobianHybridByKey()`, `PetscFormKey`
5893: @*/
5894: PetscErrorCode DMPlexComputeBdJacobianSingleByLabel(DM dm, PetscWeakForm wf, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt fieldI, IS facetIS, Vec locX, Vec locX_t, PetscReal t, DMField coordField, PetscReal X_tShift, Mat Jac, Mat JacP)
5895: {
5896: DM_Plex *mesh = (DM_Plex *)dm->data;
5897: DM plex = NULL, plexA = NULL, tdm;
5898: DMEnclosureType encAux;
5899: PetscDS ds, dsAux = NULL;
5900: PetscSection section, sectionAux = NULL;
5901: PetscSection globalSection;
5902: Vec locA = NULL, tv;
5903: PetscScalar *u = NULL, *u_t = NULL, *a = NULL, *elemMat = NULL, *elemMatP = NULL;
5904: PetscInt v;
5905: PetscInt Nf, totDim, totDimAux = 0;
5906: PetscBool hasJac = PETSC_FALSE, hasPrec = PETSC_FALSE, transform;
5908: PetscFunctionBegin;
5909: PetscCall(DMHasBasisTransform(dm, &transform));
5910: PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
5911: PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
5912: PetscCall(DMGetLocalSection(dm, §ion));
5913: PetscCall(DMGetDS(dm, &ds));
5914: PetscCall(PetscDSGetNumFields(ds, &Nf));
5915: PetscCall(PetscDSGetTotalDimension(ds, &totDim));
5916: PetscCall(PetscWeakFormHasBdJacobian(wf, &hasJac));
5917: PetscCall(PetscWeakFormHasBdJacobianPreconditioner(wf, &hasPrec));
5918: if (!hasJac && !hasPrec) PetscFunctionReturn(PETSC_SUCCESS);
5919: PetscCall(DMConvert(dm, DMPLEX, &plex));
5920: PetscCall(DMGetAuxiliaryVec(dm, label, values[0], 0, &locA));
5921: if (locA) {
5922: DM dmAux;
5924: PetscCall(VecGetDM(locA, &dmAux));
5925: PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
5926: PetscCall(DMConvert(dmAux, DMPLEX, &plexA));
5927: PetscCall(DMGetDS(plexA, &dsAux));
5928: PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
5929: PetscCall(DMGetLocalSection(plexA, §ionAux));
5930: }
5932: PetscCall(DMGetGlobalSection(dm, &globalSection));
5933: for (v = 0; v < numValues; ++v) {
5934: PetscFEGeom *fgeom;
5935: PetscInt maxDegree;
5936: PetscQuadrature qGeom = NULL;
5937: IS pointIS;
5938: const PetscInt *points;
5939: PetscFormKey key;
5940: PetscInt numFaces, face, Nq;
5942: key.label = label;
5943: key.value = values[v];
5944: key.part = 0;
5945: PetscCall(DMLabelGetStratumIS(label, values[v], &pointIS));
5946: if (!pointIS) continue; /* No points with that id on this process */
5947: {
5948: IS isectIS;
5950: /* TODO: Special cases of ISIntersect where it is quick to check a prior if one is a superset of the other */
5951: PetscCall(ISIntersect_Caching_Internal(facetIS, pointIS, &isectIS));
5952: PetscCall(ISDestroy(&pointIS));
5953: pointIS = isectIS;
5954: }
5955: PetscCall(ISGetLocalSize(pointIS, &numFaces));
5956: PetscCall(ISGetIndices(pointIS, &points));
5957: PetscCall(PetscMalloc5(numFaces * totDim, &u, (locX_t ? (size_t)numFaces * totDim : 0), &u_t, (hasJac ? (size_t)numFaces * totDim * totDim : 0), &elemMat, (hasPrec ? (size_t)numFaces * totDim * totDim : 0), &elemMatP, (locA ? (size_t)numFaces * totDimAux : 0), &a));
5958: PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
5959: if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &qGeom));
5960: if (!qGeom) {
5961: PetscFE fe;
5963: PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&fe));
5964: PetscCall(PetscFEGetFaceQuadrature(fe, &qGeom));
5965: PetscCall(PetscObjectReference((PetscObject)qGeom));
5966: }
5967: PetscCall(PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL));
5968: PetscCall(DMSNESGetFEGeom(coordField, pointIS, qGeom, PETSC_FEGEOM_BOUNDARY, &fgeom));
5969: for (face = 0; face < numFaces; ++face) {
5970: const PetscInt point = points[face], *support;
5971: PetscScalar *x = NULL;
5972: PetscInt i;
5974: PetscCall(DMPlexGetSupport(dm, point, &support));
5975: PetscCall(DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x));
5976: for (i = 0; i < totDim; ++i) u[face * totDim + i] = x[i];
5977: PetscCall(DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x));
5978: if (locX_t) {
5979: PetscCall(DMPlexVecGetClosure(plex, section, locX_t, support[0], NULL, &x));
5980: for (i = 0; i < totDim; ++i) u_t[face * totDim + i] = x[i];
5981: PetscCall(DMPlexVecRestoreClosure(plex, section, locX_t, support[0], NULL, &x));
5982: }
5983: if (locA) {
5984: PetscInt subp;
5985: PetscCall(DMGetEnclosurePoint(plexA, dm, encAux, support[0], &subp));
5986: PetscCall(DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x));
5987: for (i = 0; i < totDimAux; ++i) a[face * totDimAux + i] = x[i];
5988: PetscCall(DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x));
5989: }
5990: }
5991: if (elemMat) PetscCall(PetscArrayzero(elemMat, numFaces * totDim * totDim));
5992: if (elemMatP) PetscCall(PetscArrayzero(elemMatP, numFaces * totDim * totDim));
5993: {
5994: PetscFE fe;
5995: PetscInt Nb;
5996: /* Conforming batches */
5997: PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
5998: /* Remainder */
5999: PetscFEGeom *chunkGeom = NULL;
6000: PetscInt fieldJ, Nr, offset;
6002: PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&fe));
6003: PetscCall(PetscFEGetDimension(fe, &Nb));
6004: PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
6005: blockSize = Nb;
6006: batchSize = numBlocks * blockSize;
6007: PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
6008: numChunks = numFaces / (numBatches * batchSize);
6009: Ne = numChunks * numBatches * batchSize;
6010: Nr = numFaces % (numBatches * batchSize);
6011: offset = numFaces - Nr;
6012: PetscCall(PetscFEGeomGetChunk(fgeom, 0, offset, &chunkGeom));
6013: for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
6014: key.field = fieldI * Nf + fieldJ;
6015: if (hasJac) PetscCall(PetscFEIntegrateBdJacobian(ds, wf, PETSCFE_JACOBIAN, key, Ne, chunkGeom, u, u_t, dsAux, a, t, X_tShift, elemMat));
6016: if (hasPrec) PetscCall(PetscFEIntegrateBdJacobian(ds, wf, PETSCFE_JACOBIAN_PRE, key, Ne, chunkGeom, u, u_t, dsAux, a, t, X_tShift, elemMatP));
6017: }
6018: PetscCall(PetscFEGeomGetChunk(fgeom, offset, numFaces, &chunkGeom));
6019: for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
6020: key.field = fieldI * Nf + fieldJ;
6021: if (hasJac)
6022: PetscCall(PetscFEIntegrateBdJacobian(ds, wf, PETSCFE_JACOBIAN, key, Nr, chunkGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), dsAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMat[offset * totDim * totDim]));
6023: if (hasPrec)
6024: PetscCall(PetscFEIntegrateBdJacobian(ds, wf, PETSCFE_JACOBIAN_PRE, key, Nr, chunkGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), dsAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMatP[offset * totDim * totDim]));
6025: }
6026: PetscCall(PetscFEGeomRestoreChunk(fgeom, offset, numFaces, &chunkGeom));
6027: }
6028: for (face = 0; face < numFaces; ++face) {
6029: const PetscInt point = points[face], *support;
6031: /* Transform to global basis before insertion in Jacobian */
6032: PetscCall(DMPlexGetSupport(plex, point, &support));
6033: if (hasJac && transform) PetscCall(DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, support[0], PETSC_TRUE, totDim, &elemMat[face * totDim * totDim]));
6034: if (hasPrec && transform) PetscCall(DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, support[0], PETSC_TRUE, totDim, &elemMatP[face * totDim * totDim]));
6035: if (hasPrec) {
6036: if (hasJac) {
6037: if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(point, "BdJacobian", totDim, totDim, &elemMat[face * totDim * totDim]));
6038: PetscCall(DMPlexMatSetClosure_Internal(plex, section, globalSection, mesh->useMatClPerm, Jac, support[0], &elemMat[face * totDim * totDim], ADD_VALUES));
6039: }
6040: if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(point, "BdJacobian", totDim, totDim, &elemMatP[face * totDim * totDim]));
6041: PetscCall(DMPlexMatSetClosure_Internal(plex, section, globalSection, mesh->useMatClPerm, JacP, support[0], &elemMatP[face * totDim * totDim], ADD_VALUES));
6042: } else {
6043: if (hasJac) {
6044: if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(point, "BdJacobian", totDim, totDim, &elemMat[face * totDim * totDim]));
6045: PetscCall(DMPlexMatSetClosure_Internal(plex, section, globalSection, mesh->useMatClPerm, Jac, support[0], &elemMat[face * totDim * totDim], ADD_VALUES));
6046: }
6047: }
6048: }
6049: PetscCall(DMSNESRestoreFEGeom(coordField, pointIS, qGeom, PETSC_TRUE, &fgeom));
6050: PetscCall(PetscQuadratureDestroy(&qGeom));
6051: PetscCall(ISRestoreIndices(pointIS, &points));
6052: PetscCall(ISDestroy(&pointIS));
6053: PetscCall(PetscFree5(u, u_t, elemMat, elemMatP, a));
6054: }
6055: if (plex) PetscCall(DMDestroy(&plex));
6056: if (plexA) PetscCall(DMDestroy(&plexA));
6057: PetscFunctionReturn(PETSC_SUCCESS);
6058: }
6060: /*@
6061: DMPlexComputeBdJacobianSingle - Compute the local boundary Jacobian
6063: Not collective
6065: Input Parameters:
6066: + dm - The output `DM`
6067: . wf - The `PetscWeakForm` holding forms on this boundary
6068: . label - The `DMLabel` indicating what faces should be integrated over
6069: . numValues - The number of label values
6070: . values - The array of label values
6071: . fieldI - The test field for these integrals
6072: . locX - The local solution
6073: . locX_t - The time derivative of the local solution, or `NULL` for time-independent problems
6074: . t - The time
6075: - X_tShift - The multiplier for dF/dxdot
6077: Output Parameters:
6078: + Jac - The local Jacobian
6079: - JacP - The local Jacobian preconditioner
6081: Level: developer
6083: .seealso: `DMPlexComputeBdJacobianSingleByLabel()`, `DMPlexComputeJacobianByKey()`, `DMPlexComputeResidualHybridByKey()`, `DMPlexComputeJacobianHybridByKey()`, `PetscFormKey`
6084: @*/
6085: PetscErrorCode DMPlexComputeBdJacobianSingle(DM dm, PetscWeakForm wf, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt fieldI, Vec locX, Vec locX_t, PetscReal t, PetscReal X_tShift, Mat Jac, Mat JacP)
6086: {
6087: DMField coordField;
6088: DMLabel depthLabel;
6089: IS facetIS;
6090: PetscInt dim;
6092: PetscFunctionBegin;
6093: PetscCall(DMGetDimension(dm, &dim));
6094: PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
6095: PetscCall(DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS));
6096: PetscCall(DMGetCoordinateField(dm, &coordField));
6097: PetscCall(DMPlexComputeBdJacobianSingleByLabel(dm, wf, label, numValues, values, fieldI, facetIS, locX, locX_t, t, coordField, X_tShift, Jac, JacP));
6098: PetscCall(ISDestroy(&facetIS));
6099: PetscFunctionReturn(PETSC_SUCCESS);
6100: }
6102: static PetscErrorCode DMPlexComputeBdJacobian_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
6103: {
6104: PetscDS prob;
6105: PetscInt dim, numBd, bd;
6106: DMLabel depthLabel;
6107: DMField coordField = NULL;
6108: IS facetIS;
6110: PetscFunctionBegin;
6111: PetscCall(DMGetDS(dm, &prob));
6112: PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
6113: PetscCall(DMGetDimension(dm, &dim));
6114: PetscCall(DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS));
6115: PetscCall(PetscDSGetNumBoundary(prob, &numBd));
6116: PetscCall(DMGetCoordinateField(dm, &coordField));
6117: for (bd = 0; bd < numBd; ++bd) {
6118: PetscWeakForm wf;
6119: DMBoundaryConditionType type;
6120: DMLabel label;
6121: const PetscInt *values;
6122: PetscInt fieldI, numValues;
6123: PetscObject obj;
6124: PetscClassId id;
6126: PetscCall(PetscDSGetBoundary(prob, bd, &wf, &type, NULL, &label, &numValues, &values, &fieldI, NULL, NULL, NULL, NULL, NULL));
6127: if (type & DM_BC_ESSENTIAL) continue;
6128: PetscCall(PetscDSGetDiscretization(prob, fieldI, &obj));
6129: PetscCall(PetscObjectGetClassId(obj, &id));
6130: if (id != PETSCFE_CLASSID) continue;
6131: PetscCall(DMPlexComputeBdJacobianSingleByLabel(dm, wf, label, numValues, values, fieldI, facetIS, locX, locX_t, t, coordField, X_tShift, Jac, JacP));
6132: }
6133: PetscCall(ISDestroy(&facetIS));
6134: PetscFunctionReturn(PETSC_SUCCESS);
6135: }
6137: /*@
6138: DMPlexComputeJacobianByKey - Compute the local Jacobian for terms matching the input key
6140: Collective
6142: Input Parameters:
6143: + dm - The output `DM`
6144: . key - The `PetscFormKey` indicating what should be integrated
6145: . cellIS - The `IS` give a set of cells to integrate over
6146: . t - The time
6147: . X_tShift - The multiplier for the Jacobian with respect to $X_t$
6148: . locX - The local solution
6149: . locX_t - The time derivative of the local solution, or `NULL` for time-independent problems
6150: - user - An optional user context, passed to the pointwise functions
6152: Output Parameters:
6153: + Jac - The local Jacobian
6154: - JacP - The local Jacobian preconditioner
6156: Level: developer
6158: .seealso: `DMPlexComputeResidualByKey()`, `DMPlexComputeResidualHybridByKey()`, `DMPlexComputeJacobianHybridByKey()`, `PetscFormKey`
6159: @*/
6160: PetscErrorCode DMPlexComputeJacobianByKey(DM dm, PetscFormKey key, IS cellIS, PetscReal t, PetscReal X_tShift, Vec locX, Vec locX_t, Mat Jac, Mat JacP, void *user)
6161: {
6162: DM_Plex *mesh = (DM_Plex *)dm->data;
6163: const char *name = "Jacobian";
6164: DM dmAux = NULL, plex, tdm;
6165: DMEnclosureType encAux;
6166: Vec A, tv;
6167: DMField coordField;
6168: PetscDS prob, probAux = NULL;
6169: PetscSection section, globalSection, sectionAux;
6170: PetscScalar *elemMat, *elemMatP, *elemMatD, *u, *u_t, *a = NULL;
6171: const PetscInt *cells;
6172: PetscInt Nf, fieldI, fieldJ;
6173: PetscInt totDim, totDimAux = 0, cStart, cEnd, numCells, c;
6174: PetscBool hasJac = PETSC_FALSE, hasPrec = PETSC_FALSE, hasDyn, hasFV = PETSC_FALSE, transform;
6176: PetscFunctionBegin;
6177: PetscCall(PetscLogEventBegin(DMPLEX_JacobianFEM, dm, 0, 0, 0));
6178: PetscCall(DMGetLocalSection(dm, §ion));
6179: PetscCall(DMGetGlobalSection(dm, &globalSection));
6180: PetscCall(DMGetAuxiliaryVec(dm, key.label, key.value, key.part, &A));
6181: if (A) {
6182: PetscCall(VecGetDM(A, &dmAux));
6183: PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
6184: PetscCall(DMConvert(dmAux, DMPLEX, &plex));
6185: PetscCall(DMGetLocalSection(plex, §ionAux));
6186: PetscCall(DMGetDS(dmAux, &probAux));
6187: PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
6188: }
6189: PetscCall(DMGetCoordinateField(dm, &coordField));
6190: if (!cellIS) goto end;
6191: PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
6192: PetscCall(ISGetLocalSize(cellIS, &numCells));
6193: if (cStart >= cEnd) goto end;
6194: PetscCall(DMHasBasisTransform(dm, &transform));
6195: PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
6196: PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
6197: PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &prob, NULL));
6198: PetscCall(PetscDSGetNumFields(prob, &Nf));
6199: PetscCall(PetscDSGetTotalDimension(prob, &totDim));
6200: PetscCall(PetscDSHasJacobian(prob, &hasJac));
6201: PetscCall(PetscDSHasJacobianPreconditioner(prob, &hasPrec));
6202: /* user passed in the same matrix, avoid double contributions and
6203: only assemble the Jacobian */
6204: if (hasJac && Jac == JacP) hasPrec = PETSC_FALSE;
6205: PetscCall(PetscDSHasDynamicJacobian(prob, &hasDyn));
6206: hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
6207: PetscCall(PetscMalloc5(numCells * totDim, &u, (locX_t ? (size_t)numCells * totDim : 0), &u_t, (hasJac ? (size_t)numCells * totDim * totDim : 0), &elemMat, (hasPrec ? (size_t)numCells * totDim * totDim : 0), &elemMatP, (hasDyn ? (size_t)numCells * totDim * totDim : 0), &elemMatD));
6208: if (dmAux) PetscCall(PetscMalloc1(numCells * totDimAux, &a));
6209: for (c = cStart; c < cEnd; ++c) {
6210: const PetscInt cell = cells ? cells[c] : c;
6211: const PetscInt cind = c - cStart;
6212: PetscScalar *x = NULL, *x_t = NULL;
6213: PetscInt i;
6215: PetscCall(DMPlexVecGetClosure(dm, section, locX, cell, NULL, &x));
6216: for (i = 0; i < totDim; ++i) u[cind * totDim + i] = x[i];
6217: PetscCall(DMPlexVecRestoreClosure(dm, section, locX, cell, NULL, &x));
6218: if (locX_t) {
6219: PetscCall(DMPlexVecGetClosure(dm, section, locX_t, cell, NULL, &x_t));
6220: for (i = 0; i < totDim; ++i) u_t[cind * totDim + i] = x_t[i];
6221: PetscCall(DMPlexVecRestoreClosure(dm, section, locX_t, cell, NULL, &x_t));
6222: }
6223: if (dmAux) {
6224: PetscInt subcell;
6225: PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, cell, &subcell));
6226: PetscCall(DMPlexVecGetClosure(plex, sectionAux, A, subcell, NULL, &x));
6227: for (i = 0; i < totDimAux; ++i) a[cind * totDimAux + i] = x[i];
6228: PetscCall(DMPlexVecRestoreClosure(plex, sectionAux, A, subcell, NULL, &x));
6229: }
6230: }
6231: if (hasJac) PetscCall(PetscArrayzero(elemMat, numCells * totDim * totDim));
6232: if (hasPrec) PetscCall(PetscArrayzero(elemMatP, numCells * totDim * totDim));
6233: if (hasDyn) PetscCall(PetscArrayzero(elemMatD, numCells * totDim * totDim));
6234: for (fieldI = 0; fieldI < Nf; ++fieldI) {
6235: PetscClassId id;
6236: PetscFE fe;
6237: PetscQuadrature qGeom = NULL;
6238: PetscInt Nb;
6239: /* Conforming batches */
6240: PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
6241: /* Remainder */
6242: PetscInt Nr, offset, Nq;
6243: PetscInt maxDegree;
6244: PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;
6246: PetscCall(PetscDSGetDiscretization(prob, fieldI, (PetscObject *)&fe));
6247: PetscCall(PetscObjectGetClassId((PetscObject)fe, &id));
6248: if (id == PETSCFV_CLASSID) {
6249: hasFV = PETSC_TRUE;
6250: continue;
6251: }
6252: PetscCall(PetscFEGetDimension(fe, &Nb));
6253: PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
6254: PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
6255: if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, &qGeom));
6256: if (!qGeom) {
6257: PetscCall(PetscFEGetQuadrature(fe, &qGeom));
6258: PetscCall(PetscObjectReference((PetscObject)qGeom));
6259: }
6260: PetscCall(PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL));
6261: PetscCall(DMSNESGetFEGeom(coordField, cellIS, qGeom, PETSC_FEGEOM_BASIC, &cgeomFEM));
6262: blockSize = Nb;
6263: batchSize = numBlocks * blockSize;
6264: PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
6265: numChunks = numCells / (numBatches * batchSize);
6266: Ne = numChunks * numBatches * batchSize;
6267: Nr = numCells % (numBatches * batchSize);
6268: offset = numCells - Nr;
6269: PetscCall(PetscFEGeomGetChunk(cgeomFEM, 0, offset, &chunkGeom));
6270: PetscCall(PetscFEGeomGetChunk(cgeomFEM, offset, numCells, &remGeom));
6271: for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
6272: key.field = fieldI * Nf + fieldJ;
6273: if (hasJac) {
6274: PetscCall(PetscFEIntegrateJacobian(prob, prob, PETSCFE_JACOBIAN, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat));
6275: PetscCall(PetscFEIntegrateJacobian(prob, prob, PETSCFE_JACOBIAN, key, Nr, remGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), probAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMat[offset * totDim * totDim]));
6276: }
6277: if (hasPrec) {
6278: PetscCall(PetscFEIntegrateJacobian(prob, prob, PETSCFE_JACOBIAN_PRE, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatP));
6279: PetscCall(PetscFEIntegrateJacobian(prob, prob, PETSCFE_JACOBIAN_PRE, key, Nr, remGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), probAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMatP[offset * totDim * totDim]));
6280: }
6281: if (hasDyn) {
6282: PetscCall(PetscFEIntegrateJacobian(prob, prob, PETSCFE_JACOBIAN_DYN, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD));
6283: PetscCall(PetscFEIntegrateJacobian(prob, prob, PETSCFE_JACOBIAN_DYN, key, Nr, remGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), probAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMatD[offset * totDim * totDim]));
6284: }
6285: }
6286: PetscCall(PetscFEGeomRestoreChunk(cgeomFEM, offset, numCells, &remGeom));
6287: PetscCall(PetscFEGeomRestoreChunk(cgeomFEM, 0, offset, &chunkGeom));
6288: PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, qGeom, PETSC_FALSE, &cgeomFEM));
6289: PetscCall(PetscQuadratureDestroy(&qGeom));
6290: }
6291: /* Add contribution from X_t */
6292: if (hasDyn) {
6293: for (c = 0; c < numCells * totDim * totDim; ++c) elemMat[c] += X_tShift * elemMatD[c];
6294: }
6295: if (hasFV) {
6296: PetscClassId id;
6297: PetscFV fv;
6298: PetscInt offsetI, NcI, NbI = 1, fc, f;
6300: for (fieldI = 0; fieldI < Nf; ++fieldI) {
6301: PetscCall(PetscDSGetDiscretization(prob, fieldI, (PetscObject *)&fv));
6302: PetscCall(PetscDSGetFieldOffset(prob, fieldI, &offsetI));
6303: PetscCall(PetscObjectGetClassId((PetscObject)fv, &id));
6304: if (id != PETSCFV_CLASSID) continue;
6305: /* Put in the weighted identity */
6306: PetscCall(PetscFVGetNumComponents(fv, &NcI));
6307: for (c = cStart; c < cEnd; ++c) {
6308: const PetscInt cind = c - cStart;
6309: const PetscInt eOffset = cind * totDim * totDim;
6310: PetscReal vol;
6312: PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL));
6313: for (fc = 0; fc < NcI; ++fc) {
6314: for (f = 0; f < NbI; ++f) {
6315: const PetscInt i = offsetI + f * NcI + fc;
6316: if (hasPrec) {
6317: if (hasJac) elemMat[eOffset + i * totDim + i] = vol;
6318: elemMatP[eOffset + i * totDim + i] = vol;
6319: } else {
6320: elemMat[eOffset + i * totDim + i] = vol;
6321: }
6322: }
6323: }
6324: }
6325: }
6326: /* No allocated space for FV stuff, so ignore the zero entries */
6327: PetscCall(MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
6328: }
6329: /* Insert values into matrix */
6330: for (c = cStart; c < cEnd; ++c) {
6331: const PetscInt cell = cells ? cells[c] : c;
6332: const PetscInt cind = c - cStart;
6334: /* Transform to global basis before insertion in Jacobian */
6335: if (transform) PetscCall(DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, cell, PETSC_TRUE, totDim, &elemMat[cind * totDim * totDim]));
6336: if (hasPrec) {
6337: if (hasJac) {
6338: if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind * totDim * totDim]));
6339: PetscCall(DMPlexMatSetClosure_Internal(dm, section, globalSection, mesh->useMatClPerm, Jac, cell, &elemMat[cind * totDim * totDim], ADD_VALUES));
6340: }
6341: if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatP[cind * totDim * totDim]));
6342: PetscCall(DMPlexMatSetClosure_Internal(dm, section, globalSection, mesh->useMatClPerm, JacP, cell, &elemMatP[cind * totDim * totDim], ADD_VALUES));
6343: } else {
6344: if (hasJac) {
6345: if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind * totDim * totDim]));
6346: PetscCall(DMPlexMatSetClosure_Internal(dm, section, globalSection, mesh->useMatClPerm, JacP, cell, &elemMat[cind * totDim * totDim], ADD_VALUES));
6347: }
6348: }
6349: }
6350: PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
6351: if (hasFV) PetscCall(MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE));
6352: PetscCall(PetscFree5(u, u_t, elemMat, elemMatP, elemMatD));
6353: if (dmAux) PetscCall(PetscFree(a));
6354: /* Compute boundary integrals */
6355: PetscCall(DMPlexComputeBdJacobian_Internal(dm, locX, locX_t, t, X_tShift, Jac, JacP, user));
6356: /* Assemble matrix */
6357: end: {
6358: PetscBool assOp = hasJac && hasPrec ? PETSC_TRUE : PETSC_FALSE, gassOp;
6360: if (dmAux) PetscCall(DMDestroy(&plex));
6361: PetscCallMPI(MPIU_Allreduce(&assOp, &gassOp, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)dm)));
6362: if (hasJac && hasPrec) {
6363: PetscCall(MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY));
6364: PetscCall(MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY));
6365: }
6366: }
6367: PetscCall(MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY));
6368: PetscCall(MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY));
6369: PetscCall(PetscLogEventEnd(DMPLEX_JacobianFEM, dm, 0, 0, 0));
6370: PetscFunctionReturn(PETSC_SUCCESS);
6371: }
6373: PetscErrorCode DMPlexComputeJacobianByKeyGeneral(DM dmr, DM dmc, PetscFormKey key, IS cellIS, PetscReal t, PetscReal X_tShift, Vec locX, Vec locX_t, Mat Jac, Mat JacP, void *user)
6374: {
6375: DM_Plex *mesh = (DM_Plex *)dmr->data;
6376: const char *name = "Jacobian";
6377: DM dmAux = NULL, plex, tdm;
6378: PetscInt printFEM = mesh->printFEM;
6379: PetscBool clPerm = mesh->useMatClPerm;
6380: DMEnclosureType encAux;
6381: Vec A, tv;
6382: DMField coordField;
6383: PetscDS rds, cds, dsAux = NULL;
6384: PetscSection rsection, rglobalSection, csection, cglobalSection, sectionAux;
6385: PetscScalar *elemMat, *elemMatP, *elemMatD, *u, *u_t, *a = NULL;
6386: const PetscInt *cells;
6387: PetscInt Nf, cNf;
6388: PetscInt totDim, ctotDim, totDimAux = 0, cStart, cEnd, numCells;
6389: PetscBool hasJac = PETSC_FALSE, hasPrec = PETSC_FALSE, hasDyn, hasFV = PETSC_FALSE, transform;
6390: MPI_Comm comm;
6392: PetscFunctionBegin;
6393: PetscCall(PetscObjectGetComm((PetscObject)dmr, &comm));
6394: PetscCall(PetscLogEventBegin(DMPLEX_JacobianFEM, dmr, 0, 0, 0));
6395: PetscCall(DMGetLocalSection(dmr, &rsection));
6396: PetscCall(DMGetGlobalSection(dmr, &rglobalSection));
6397: PetscCall(DMGetLocalSection(dmc, &csection));
6398: PetscCall(DMGetGlobalSection(dmc, &cglobalSection));
6399: PetscCall(DMGetAuxiliaryVec(dmr, key.label, key.value, key.part, &A));
6400: if (A) {
6401: PetscCall(VecGetDM(A, &dmAux));
6402: PetscCall(DMGetEnclosureRelation(dmAux, dmr, &encAux));
6403: PetscCall(DMConvert(dmAux, DMPLEX, &plex));
6404: PetscCall(DMGetLocalSection(plex, §ionAux));
6405: PetscCall(DMGetDS(dmAux, &dsAux));
6406: PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
6407: }
6408: PetscCall(DMGetCoordinateField(dmr, &coordField));
6409: if (!cellIS) goto end;
6410: PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
6411: PetscCall(ISGetLocalSize(cellIS, &numCells));
6412: if (cStart >= cEnd) goto end;
6413: PetscCall(DMHasBasisTransform(dmr, &transform));
6414: PetscCall(DMGetBasisTransformDM_Internal(dmr, &tdm));
6415: PetscCall(DMGetBasisTransformVec_Internal(dmr, &tv));
6416: PetscCall(DMGetCellDS(dmr, cells ? cells[cStart] : cStart, &rds, NULL));
6417: PetscCall(DMGetCellDS(dmc, cells ? cells[cStart] : cStart, &cds, NULL));
6418: PetscCall(PetscDSGetNumFields(rds, &Nf));
6419: PetscCall(PetscDSGetNumFields(cds, &cNf));
6420: PetscCheck(Nf == cNf, comm, PETSC_ERR_ARG_WRONG, "Number of row fields %" PetscInt_FMT " != %" PetscInt_FMT " number of columns field", Nf, cNf);
6421: PetscCall(PetscDSGetTotalDimension(rds, &totDim));
6422: PetscCall(PetscDSGetTotalDimension(cds, &ctotDim));
6423: PetscCall(PetscDSHasJacobian(rds, &hasJac));
6424: PetscCall(PetscDSHasJacobianPreconditioner(rds, &hasPrec));
6425: /* user passed in the same matrix, avoid double contributions and
6426: only assemble the Jacobian */
6427: if (hasJac && Jac == JacP) hasPrec = PETSC_FALSE;
6428: PetscCall(PetscDSHasDynamicJacobian(rds, &hasDyn));
6429: hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
6430: PetscCall(PetscMalloc5(numCells * totDim, &u, (locX_t ? (size_t)numCells * totDim : 0), &u_t, (hasJac ? (size_t)numCells * totDim * ctotDim : 0), &elemMat, (hasPrec ? (size_t)numCells * totDim * ctotDim : 0), &elemMatP, (hasDyn ? (size_t)numCells * totDim * ctotDim : 0), &elemMatD));
6431: if (dmAux) PetscCall(PetscMalloc1(numCells * totDimAux, &a));
6432: for (PetscInt c = cStart; c < cEnd; ++c) {
6433: const PetscInt cell = cells ? cells[c] : c;
6434: const PetscInt cind = c - cStart;
6435: PetscScalar *x = NULL, *x_t = NULL;
6436: PetscInt i;
6438: PetscCall(DMPlexVecGetClosure(dmr, rsection, locX, cell, NULL, &x));
6439: for (i = 0; i < totDim; ++i) u[cind * totDim + i] = x[i];
6440: PetscCall(DMPlexVecRestoreClosure(dmr, rsection, locX, cell, NULL, &x));
6441: if (locX_t) {
6442: PetscCall(DMPlexVecGetClosure(dmr, rsection, locX_t, cell, NULL, &x_t));
6443: for (i = 0; i < totDim; ++i) u_t[cind * totDim + i] = x_t[i];
6444: PetscCall(DMPlexVecRestoreClosure(dmr, rsection, locX_t, cell, NULL, &x_t));
6445: }
6446: if (dmAux) {
6447: PetscInt subcell;
6448: PetscCall(DMGetEnclosurePoint(dmAux, dmr, encAux, cell, &subcell));
6449: PetscCall(DMPlexVecGetClosure(plex, sectionAux, A, subcell, NULL, &x));
6450: for (i = 0; i < totDimAux; ++i) a[cind * totDimAux + i] = x[i];
6451: PetscCall(DMPlexVecRestoreClosure(plex, sectionAux, A, subcell, NULL, &x));
6452: }
6453: }
6454: if (hasJac) PetscCall(PetscArrayzero(elemMat, numCells * totDim * ctotDim));
6455: if (hasPrec) PetscCall(PetscArrayzero(elemMatP, numCells * totDim * ctotDim));
6456: if (hasDyn) PetscCall(PetscArrayzero(elemMatD, numCells * totDim * ctotDim));
6457: for (PetscInt fieldI = 0; fieldI < Nf; ++fieldI) {
6458: PetscClassId id;
6459: PetscFE fe;
6460: PetscQuadrature qGeom = NULL;
6461: PetscInt Nb;
6462: /* Conforming batches */
6463: PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
6464: /* Remainder */
6465: PetscInt Nr, offset, Nq;
6466: PetscInt maxDegree;
6467: PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;
6469: PetscCall(PetscDSGetDiscretization(rds, fieldI, (PetscObject *)&fe));
6470: PetscCall(PetscObjectGetClassId((PetscObject)fe, &id));
6471: if (id == PETSCFV_CLASSID) {
6472: hasFV = PETSC_TRUE;
6473: continue;
6474: }
6475: PetscCall(PetscFEGetDimension(fe, &Nb));
6476: PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
6477: PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
6478: if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, &qGeom));
6479: if (!qGeom) {
6480: PetscCall(PetscFEGetQuadrature(fe, &qGeom));
6481: PetscCall(PetscObjectReference((PetscObject)qGeom));
6482: }
6483: PetscCall(PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL));
6484: PetscCall(DMSNESGetFEGeom(coordField, cellIS, qGeom, PETSC_FEGEOM_BASIC, &cgeomFEM));
6485: blockSize = Nb;
6486: batchSize = numBlocks * blockSize;
6487: PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
6488: numChunks = numCells / (numBatches * batchSize);
6489: Ne = numChunks * numBatches * batchSize;
6490: Nr = numCells % (numBatches * batchSize);
6491: offset = numCells - Nr;
6492: PetscCall(PetscFEGeomGetChunk(cgeomFEM, 0, offset, &chunkGeom));
6493: PetscCall(PetscFEGeomGetChunk(cgeomFEM, offset, numCells, &remGeom));
6494: for (PetscInt fieldJ = 0; fieldJ < Nf; ++fieldJ) {
6495: key.field = fieldI * Nf + fieldJ;
6496: if (hasJac) {
6497: PetscCall(PetscFEIntegrateJacobian(rds, cds, PETSCFE_JACOBIAN, key, Ne, chunkGeom, u, u_t, dsAux, a, t, X_tShift, elemMat));
6498: PetscCall(PetscFEIntegrateJacobian(rds, cds, PETSCFE_JACOBIAN, key, Nr, remGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), dsAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMat[offset * totDim * ctotDim]));
6499: }
6500: if (hasPrec) {
6501: PetscCall(PetscFEIntegrateJacobian(rds, cds, PETSCFE_JACOBIAN_PRE, key, Ne, chunkGeom, u, u_t, dsAux, a, t, X_tShift, elemMatP));
6502: PetscCall(PetscFEIntegrateJacobian(rds, cds, PETSCFE_JACOBIAN_PRE, key, Nr, remGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), dsAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMatP[offset * totDim * ctotDim]));
6503: }
6504: if (hasDyn) {
6505: PetscCall(PetscFEIntegrateJacobian(rds, cds, PETSCFE_JACOBIAN_DYN, key, Ne, chunkGeom, u, u_t, dsAux, a, t, X_tShift, elemMatD));
6506: PetscCall(PetscFEIntegrateJacobian(rds, cds, PETSCFE_JACOBIAN_DYN, key, Nr, remGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), dsAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMatD[offset * totDim * ctotDim]));
6507: }
6508: }
6509: PetscCall(PetscFEGeomRestoreChunk(cgeomFEM, offset, numCells, &remGeom));
6510: PetscCall(PetscFEGeomRestoreChunk(cgeomFEM, 0, offset, &chunkGeom));
6511: PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, qGeom, PETSC_FALSE, &cgeomFEM));
6512: PetscCall(PetscQuadratureDestroy(&qGeom));
6513: }
6514: /* Add contribution from X_t */
6515: if (hasDyn) {
6516: for (PetscInt c = 0; c < numCells * totDim * ctotDim; ++c) elemMat[c] += X_tShift * elemMatD[c];
6517: }
6518: if (hasFV) {
6519: PetscClassId id;
6520: PetscFV fv;
6521: PetscInt offsetI, NcI, NbI = 1;
6523: for (PetscInt fieldI = 0; fieldI < Nf; ++fieldI) {
6524: PetscCall(PetscDSGetDiscretization(rds, fieldI, (PetscObject *)&fv));
6525: PetscCall(PetscDSGetFieldOffset(rds, fieldI, &offsetI));
6526: PetscCall(PetscObjectGetClassId((PetscObject)fv, &id));
6527: if (id != PETSCFV_CLASSID) continue;
6528: /* Put in the weighted identity */
6529: PetscCall(PetscFVGetNumComponents(fv, &NcI));
6530: for (PetscInt c = cStart; c < cEnd; ++c) {
6531: const PetscInt cind = c - cStart;
6532: const PetscInt eOffset = cind * totDim * ctotDim;
6533: PetscReal vol;
6535: PetscCall(DMPlexComputeCellGeometryFVM(dmr, c, &vol, NULL, NULL));
6536: for (PetscInt fc = 0; fc < NcI; ++fc) {
6537: for (PetscInt f = 0; f < NbI; ++f) {
6538: const PetscInt i = offsetI + f * NcI + fc;
6539: if (hasPrec) {
6540: if (hasJac) elemMat[eOffset + i * ctotDim + i] = vol;
6541: elemMatP[eOffset + i * ctotDim + i] = vol;
6542: } else {
6543: elemMat[eOffset + i * ctotDim + i] = vol;
6544: }
6545: }
6546: }
6547: }
6548: }
6549: /* No allocated space for FV stuff, so ignore the zero entries */
6550: PetscCall(MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
6551: }
6552: /* Insert values into matrix */
6553: for (PetscInt c = cStart; c < cEnd; ++c) {
6554: const PetscInt cell = cells ? cells[c] : c;
6555: const PetscInt cind = c - cStart;
6557: /* Transform to global basis before insertion in Jacobian */
6558: if (transform) PetscCall(DMPlexBasisTransformPointTensor_Internal(dmr, tdm, tv, cell, PETSC_TRUE, totDim, &elemMat[cind * totDim * ctotDim]));
6559: if (hasPrec) {
6560: if (hasJac) {
6561: if (printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, ctotDim, &elemMat[cind * totDim * ctotDim]));
6562: PetscCall(DMPlexMatSetClosureGeneral(dmr, rsection, rglobalSection, clPerm, dmc, csection, cglobalSection, clPerm, Jac, cell, &elemMat[cind * totDim * ctotDim], ADD_VALUES));
6563: }
6564: if (printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, ctotDim, &elemMatP[cind * totDim * ctotDim]));
6565: PetscCall(DMPlexMatSetClosureGeneral(dmr, rsection, rglobalSection, clPerm, dmc, csection, cglobalSection, clPerm, JacP, cell, &elemMatP[cind * totDim * ctotDim], ADD_VALUES));
6566: } else {
6567: if (hasJac) {
6568: if (printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, ctotDim, &elemMat[cind * totDim * ctotDim]));
6569: PetscCall(DMPlexMatSetClosureGeneral(dmr, rsection, rglobalSection, clPerm, dmc, csection, cglobalSection, clPerm, JacP, cell, &elemMat[cind * totDim * ctotDim], ADD_VALUES));
6570: }
6571: }
6572: }
6573: PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
6574: if (hasFV) PetscCall(MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE));
6575: PetscCall(PetscFree5(u, u_t, elemMat, elemMatP, elemMatD));
6576: if (dmAux) PetscCall(PetscFree(a));
6577: /* Compute boundary integrals */
6578: PetscCall(DMPlexComputeBdJacobian_Internal(dmr, locX, locX_t, t, X_tShift, Jac, JacP, user));
6579: /* Assemble matrix */
6580: end: {
6581: PetscBool assOp = hasJac && hasPrec ? PETSC_TRUE : PETSC_FALSE, gassOp;
6583: if (dmAux) PetscCall(DMDestroy(&plex));
6584: PetscCallMPI(MPIU_Allreduce(&assOp, &gassOp, 1, MPIU_BOOL, MPI_LOR, comm));
6585: if (hasJac && hasPrec) {
6586: PetscCall(MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY));
6587: PetscCall(MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY));
6588: }
6589: }
6590: PetscCall(MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY));
6591: PetscCall(MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY));
6592: PetscCall(PetscLogEventEnd(DMPLEX_JacobianFEM, dmr, 0, 0, 0));
6593: PetscFunctionReturn(PETSC_SUCCESS);
6594: }
6596: /*@
6597: DMPlexComputeJacobianHybridByKey - Compute the local Jacobian over hybrid cells for terms matching the input key
6599: Collective
6601: Input Parameters:
6602: + dm - The output `DM`
6603: . key - The `PetscFormKey` array (left cell, right cell, cohesive cell) indicating what should be integrated
6604: . cellIS - The `IS` give a set of cells to integrate over
6605: . t - The time
6606: . X_tShift - The multiplier for the Jacobian with respect to $X_t$
6607: . locX - The local solution
6608: . locX_t - The time derivative of the local solution, or `NULL` for time-independent problems
6609: - user - An optional user context, passed to the pointwise functions
6611: Output Parameters:
6612: + Jac - The local Jacobian
6613: - JacP - The local Jacobian preconditioner
6615: Level: developer
6617: .seealso: `DMPlexComputeResidualByKey()`, `DMPlexComputeJacobianByKey()`, `DMPlexComputeResidualHybridByKey()`, `PetscFormKey`
6618: @*/
6619: PetscErrorCode DMPlexComputeJacobianHybridByKey(DM dm, PetscFormKey key[], IS cellIS, PetscReal t, PetscReal X_tShift, Vec locX, Vec locX_t, Mat Jac, Mat JacP, void *user)
6620: {
6621: DM_Plex *mesh = (DM_Plex *)dm->data;
6622: const char *name = "Hybrid Jacobian";
6623: DM dmAux[3] = {NULL, NULL, NULL};
6624: DMLabel ghostLabel = NULL;
6625: DM plex = NULL;
6626: DM plexA = NULL;
6627: PetscDS ds = NULL;
6628: PetscDS dsIn = NULL;
6629: PetscDS dsAux[3] = {NULL, NULL, NULL};
6630: Vec locA[3] = {NULL, NULL, NULL};
6631: DM dmScale[3] = {NULL, NULL, NULL};
6632: PetscDS dsScale[3] = {NULL, NULL, NULL};
6633: Vec locS[3] = {NULL, NULL, NULL};
6634: PetscSection section = NULL;
6635: PetscSection sectionAux[3] = {NULL, NULL, NULL};
6636: DMField coordField = NULL;
6637: PetscScalar *a[3] = {NULL, NULL, NULL};
6638: PetscScalar *s[3] = {NULL, NULL, NULL};
6639: PetscScalar *u = NULL, *u_t;
6640: PetscScalar *elemMatNeg, *elemMatPos, *elemMatCoh;
6641: PetscScalar *elemMatNegP, *elemMatPosP, *elemMatCohP;
6642: PetscSection globalSection;
6643: IS chunkISF, chunkISN;
6644: const PetscInt *cells;
6645: PetscInt *faces, *neighbors;
6646: PetscInt cStart, cEnd, numCells;
6647: PetscInt Nf, fieldI, fieldJ, totDim, totDimIn, totDimAux[3], totDimScale[3], numChunks, cellChunkSize, chunk;
6648: PetscInt maxDegree = PETSC_INT_MAX;
6649: PetscQuadrature affineQuadF = NULL, *quadsF = NULL;
6650: PetscFEGeom *affineGeomF = NULL, **geomsF = NULL;
6651: PetscQuadrature affineQuadN = NULL;
6652: PetscFEGeom *affineGeomN = NULL;
6653: PetscBool hasBdJac, hasBdPrec;
6655: PetscFunctionBegin;
6656: PetscCall(PetscLogEventBegin(DMPLEX_JacobianFEM, dm, 0, 0, 0));
6657: if (!cellIS) goto end;
6658: PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
6659: PetscCall(ISGetLocalSize(cellIS, &numCells));
6660: if (cStart >= cEnd) goto end;
6661: if ((key[0].label == key[1].label) && (key[0].value == key[1].value) && (key[0].part == key[1].part)) {
6662: const char *name;
6663: PetscCall(PetscObjectGetName((PetscObject)key[0].label, &name));
6664: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Form keys for each side of a cohesive surface must be different (%s, %" PetscInt_FMT ", %" PetscInt_FMT ")", name, key[0].value, key[0].part);
6665: }
6666: PetscCall(DMConvert(dm, DMPLEX, &plex));
6667: PetscCall(DMGetLocalSection(dm, §ion));
6668: PetscCall(DMGetGlobalSection(dm, &globalSection));
6669: PetscCall(DMGetLabel(dm, "ghost", &ghostLabel));
6670: PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &ds, &dsIn));
6671: PetscCall(PetscDSGetNumFields(ds, &Nf));
6672: PetscCall(PetscDSGetTotalDimension(ds, &totDim));
6673: PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn));
6674: PetscCall(PetscDSHasBdJacobian(ds, &hasBdJac));
6675: PetscCall(PetscDSHasBdJacobianPreconditioner(ds, &hasBdPrec));
6676: PetscCall(DMGetAuxiliaryVec(dm, key[2].label, key[2].value, key[2].part, &locA[2]));
6677: if (locA[2]) {
6678: const PetscInt cellStart = cells ? cells[cStart] : cStart;
6680: PetscCall(VecGetDM(locA[2], &dmAux[2]));
6681: PetscCall(DMConvert(dmAux[2], DMPLEX, &plexA));
6682: PetscCall(DMGetLocalSection(dmAux[2], §ionAux[2]));
6683: PetscCall(DMGetCellDS(dmAux[2], cellStart, &dsAux[2], NULL));
6684: PetscCall(PetscDSGetTotalDimension(dsAux[2], &totDimAux[2]));
6685: {
6686: const PetscInt *cone;
6687: PetscInt c;
6689: PetscCall(DMPlexGetCone(dm, cellStart, &cone));
6690: for (c = 0; c < 2; ++c) {
6691: const PetscInt *support;
6692: PetscInt ssize, s;
6694: PetscCall(DMPlexGetSupport(dm, cone[c], &support));
6695: PetscCall(DMPlexGetSupportSize(dm, cone[c], &ssize));
6696: PetscCheck(ssize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " from cell %" PetscInt_FMT " has support size %" PetscInt_FMT " != 2", cone[c], cellStart, ssize);
6697: if (support[0] == cellStart) s = 1;
6698: else if (support[1] == cellStart) s = 0;
6699: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " does not have cell %" PetscInt_FMT " in its support", cone[c], cellStart);
6700: PetscCall(DMGetAuxiliaryVec(dm, key[c].label, key[c].value, key[c].part, &locA[c]));
6701: if (locA[c]) PetscCall(VecGetDM(locA[c], &dmAux[c]));
6702: else dmAux[c] = dmAux[2];
6703: PetscCall(DMGetCellDS(dmAux[c], support[s], &dsAux[c], NULL));
6704: PetscCall(PetscDSGetTotalDimension(dsAux[c], &totDimAux[c]));
6705: }
6706: }
6707: }
6708: /* Handle mass matrix scaling
6709: The field in key[2] is the field to be scaled, and the scaling field is the first in the dsScale */
6710: PetscCall(DMGetAuxiliaryVec(dm, key[2].label, -key[2].value, key[2].part, &locS[2]));
6711: if (locS[2]) {
6712: const PetscInt cellStart = cells ? cells[cStart] : cStart;
6713: PetscInt Nb, Nbs;
6715: PetscCall(VecGetDM(locS[2], &dmScale[2]));
6716: PetscCall(DMGetCellDS(dmScale[2], cells ? cells[cStart] : cStart, &dsScale[2], NULL));
6717: PetscCall(PetscDSGetTotalDimension(dsScale[2], &totDimScale[2]));
6718: // BRAD: This is not set correctly
6719: key[2].field = 2;
6720: PetscCall(PetscDSGetFieldSize(ds, key[2].field, &Nb));
6721: PetscCall(PetscDSGetFieldSize(dsScale[2], 0, &Nbs));
6722: PetscCheck(Nb == Nbs, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field %" PetscInt_FMT " of size %" PetscInt_FMT " cannot be scaled by field of size %" PetscInt_FMT, key[2].field, Nb, Nbs);
6723: {
6724: const PetscInt *cone;
6725: PetscInt c;
6727: locS[1] = locS[0] = locS[2];
6728: dmScale[1] = dmScale[0] = dmScale[2];
6729: PetscCall(DMPlexGetCone(dm, cellStart, &cone));
6730: for (c = 0; c < 2; ++c) {
6731: const PetscInt *support;
6732: PetscInt ssize, s;
6734: PetscCall(DMPlexGetSupport(dm, cone[c], &support));
6735: PetscCall(DMPlexGetSupportSize(dm, cone[c], &ssize));
6736: PetscCheck(ssize == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " from cell %" PetscInt_FMT " has support size %" PetscInt_FMT " != 2", cone[c], cellStart, ssize);
6737: if (support[0] == cellStart) s = 1;
6738: else if (support[1] == cellStart) s = 0;
6739: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " does not have cell %" PetscInt_FMT " in its support", cone[c], cellStart);
6740: PetscCall(DMGetCellDS(dmScale[c], support[s], &dsScale[c], NULL));
6741: PetscCall(PetscDSGetTotalDimension(dsScale[c], &totDimScale[c]));
6742: }
6743: }
6744: }
6745: /* 2: Setup geometric data */
6746: PetscCall(DMGetCoordinateField(dm, &coordField));
6747: PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
6748: if (maxDegree > 1) {
6749: PetscInt f;
6750: PetscCall(PetscCalloc2(Nf, &quadsF, Nf, &geomsF));
6751: for (f = 0; f < Nf; ++f) {
6752: PetscFE fe;
6754: PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
6755: if (fe) {
6756: PetscCall(PetscFEGetQuadrature(fe, &quadsF[f]));
6757: PetscCall(PetscObjectReference((PetscObject)quadsF[f]));
6758: }
6759: }
6760: }
6761: /* Loop over chunks */
6762: cellChunkSize = numCells;
6763: numChunks = !numCells ? 0 : PetscCeilReal(((PetscReal)numCells) / cellChunkSize);
6764: PetscCall(PetscCalloc2(2 * cellChunkSize, &faces, 2 * cellChunkSize, &neighbors));
6765: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 2 * cellChunkSize, faces, PETSC_USE_POINTER, &chunkISF));
6766: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 2 * cellChunkSize, neighbors, PETSC_USE_POINTER, &chunkISN));
6767: /* Extract field coefficients */
6768: /* NOTE This needs the end cap faces to have identical orientations */
6769: PetscCall(DMPlexGetHybridCellFields(dm, cellIS, locX, locX_t, locA[2], &u, &u_t, &a[2]));
6770: PetscCall(DMPlexGetHybridFields(dm, dmAux, dsAux, cellIS, locA, PETSC_TRUE, a));
6771: PetscCall(DMPlexGetHybridFields(dm, dmScale, dsScale, cellIS, locS, PETSC_TRUE, s));
6772: PetscCall(DMGetWorkArray(dm, hasBdJac ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatNeg));
6773: PetscCall(DMGetWorkArray(dm, hasBdJac ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatPos));
6774: PetscCall(DMGetWorkArray(dm, hasBdJac ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatCoh));
6775: PetscCall(DMGetWorkArray(dm, hasBdPrec ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatNegP));
6776: PetscCall(DMGetWorkArray(dm, hasBdPrec ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatPosP));
6777: PetscCall(DMGetWorkArray(dm, hasBdPrec ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatCohP));
6778: for (chunk = 0; chunk < numChunks; ++chunk) {
6779: PetscInt cS = cStart + chunk * cellChunkSize, cE = PetscMin(cS + cellChunkSize, cEnd), numCells = cE - cS, c;
6781: if (hasBdJac) {
6782: PetscCall(PetscArrayzero(elemMatNeg, cellChunkSize * totDim * totDim));
6783: PetscCall(PetscArrayzero(elemMatPos, cellChunkSize * totDim * totDim));
6784: PetscCall(PetscArrayzero(elemMatCoh, cellChunkSize * totDim * totDim));
6785: }
6786: if (hasBdPrec) {
6787: PetscCall(PetscArrayzero(elemMatNegP, cellChunkSize * totDim * totDim));
6788: PetscCall(PetscArrayzero(elemMatPosP, cellChunkSize * totDim * totDim));
6789: PetscCall(PetscArrayzero(elemMatCohP, cellChunkSize * totDim * totDim));
6790: }
6791: /* Get faces */
6792: for (c = cS; c < cE; ++c) {
6793: const PetscInt cell = cells ? cells[c] : c;
6794: const PetscInt *cone, *support;
6795: PetscCall(DMPlexGetCone(plex, cell, &cone));
6796: faces[(c - cS) * 2 + 0] = cone[0];
6797: faces[(c - cS) * 2 + 1] = cone[1];
6798: PetscCall(DMPlexGetSupport(dm, cone[0], &support));
6799: neighbors[(c - cS) * 2 + 0] = support[0] == cell ? support[1] : support[0];
6800: PetscCall(DMPlexGetSupport(dm, cone[1], &support));
6801: neighbors[(c - cS) * 2 + 1] = support[0] == cell ? support[1] : support[0];
6802: }
6803: PetscCall(ISGeneralSetIndices(chunkISF, 2 * cellChunkSize, faces, PETSC_USE_POINTER));
6804: PetscCall(ISGeneralSetIndices(chunkISN, 2 * cellChunkSize, neighbors, PETSC_USE_POINTER));
6805: if (maxDegree <= 1) {
6806: if (!affineQuadF) PetscCall(DMFieldCreateDefaultQuadrature(coordField, chunkISF, &affineQuadF));
6807: if (affineQuadF) PetscCall(DMSNESGetFEGeom(coordField, chunkISF, affineQuadF, PETSC_FEGEOM_COHESIVE, &affineGeomF));
6808: if (!affineQuadN) {
6809: PetscInt dim;
6810: PetscCall(PetscQuadratureGetData(affineQuadF, &dim, NULL, NULL, NULL, NULL));
6811: PetscCall(DMFieldCreateDefaultFaceQuadrature(coordField, chunkISN, &affineQuadN));
6812: PetscCall(PetscQuadratureSetData(affineQuadN, dim + 1, PETSC_DECIDE, PETSC_DECIDE, NULL, NULL));
6813: }
6814: if (affineQuadN) PetscCall(DMSNESGetFEGeom(coordField, chunkISN, affineQuadN, PETSC_FEGEOM_BASIC, &affineGeomN));
6815: } else {
6816: PetscInt f;
6817: for (f = 0; f < Nf; ++f) {
6818: if (quadsF[f]) PetscCall(DMSNESGetFEGeom(coordField, chunkISF, quadsF[f], PETSC_FEGEOM_COHESIVE, &geomsF[f]));
6819: }
6820: }
6822: for (fieldI = 0; fieldI < Nf; ++fieldI) {
6823: PetscFE feI;
6824: PetscFEGeom *geomF = affineGeomF ? affineGeomF : geomsF[fieldI];
6825: PetscFEGeom *chunkGeomF = NULL, *remGeomF = NULL;
6826: PetscFEGeom *geomN = affineGeomN ? affineGeomN : geomsF[fieldI];
6827: PetscFEGeom *chunkGeomN = NULL, *remGeomN = NULL;
6828: PetscQuadrature quadF = affineQuadF ? affineQuadF : quadsF[fieldI];
6829: PetscInt numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset, Nq, Nb;
6830: PetscBool isCohesiveField;
6832: PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
6833: if (!feI) continue;
6834: PetscCall(PetscFEGetTileSizes(feI, NULL, &numBlocks, NULL, &numBatches));
6835: PetscCall(PetscQuadratureGetData(quadF, NULL, NULL, &Nq, NULL, NULL));
6836: PetscCall(PetscFEGetDimension(feI, &Nb));
6837: blockSize = Nb;
6838: batchSize = numBlocks * blockSize;
6839: PetscCall(PetscFESetTileSizes(feI, blockSize, numBlocks, batchSize, numBatches));
6840: numChunks = numCells / (numBatches * batchSize);
6841: Ne = numChunks * numBatches * batchSize;
6842: Nr = numCells % (numBatches * batchSize);
6843: offset = numCells - Nr;
6844: PetscCall(PetscFEGeomGetChunk(geomF, 0, offset * 2, &chunkGeomF));
6845: PetscCall(PetscFEGeomGetChunk(geomF, offset * 2, numCells * 2, &remGeomF));
6846: PetscCall(PetscFEGeomGetChunk(geomN, 0, offset * 2, &chunkGeomN));
6847: PetscCall(PetscFEGeomGetChunk(geomN, offset * 2, numCells * 2, &remGeomN));
6848: PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveField));
6849: for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
6850: PetscFE feJ;
6852: PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
6853: if (!feJ) continue;
6854: key[0].field = fieldI * Nf + fieldJ;
6855: key[1].field = fieldI * Nf + fieldJ;
6856: key[2].field = fieldI * Nf + fieldJ;
6857: if (hasBdJac) {
6858: PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN, key[0], 0, Ne, chunkGeomF, chunkGeomN, u, u_t, dsAux[0], a[0], t, X_tShift, elemMatNeg));
6859: PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN, key[0], 0, Nr, remGeomF, remGeomN, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[0], PetscSafePointerPlusOffset(a[0], offset * totDimAux[0]), t, X_tShift, &elemMatNeg[offset * totDim * totDim]));
6860: PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN, key[1], 1, Ne, chunkGeomF, chunkGeomN, u, u_t, dsAux[1], a[1], t, X_tShift, elemMatPos));
6861: PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN, key[1], 1, Nr, remGeomF, remGeomN, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[1], PetscSafePointerPlusOffset(a[1], offset * totDimAux[1]), t, X_tShift, &elemMatPos[offset * totDim * totDim]));
6862: }
6863: if (hasBdPrec) {
6864: PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN_PRE, key[0], 0, Ne, chunkGeomF, chunkGeomN, u, u_t, dsAux[0], a[0], t, X_tShift, elemMatNegP));
6865: PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN_PRE, key[0], 0, Nr, remGeomF, remGeomN, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[0], &a[0][offset * totDimAux[0]], t, X_tShift, &elemMatNegP[offset * totDim * totDim]));
6866: PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN_PRE, key[1], 1, Ne, chunkGeomF, chunkGeomN, u, u_t, dsAux[1], a[1], t, X_tShift, elemMatPosP));
6867: PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN_PRE, key[1], 1, Nr, remGeomF, remGeomN, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[1], &a[1][offset * totDimAux[1]], t, X_tShift, &elemMatPosP[offset * totDim * totDim]));
6868: }
6869: if (hasBdJac) {
6870: PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN, key[2], 2, Ne, chunkGeomF, chunkGeomN, u, u_t, dsAux[2], a[2], t, X_tShift, elemMatCoh));
6871: PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN, key[2], 2, Nr, remGeomF, remGeomN, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[2], PetscSafePointerPlusOffset(a[2], offset * totDimAux[2]), t, X_tShift, &elemMatCoh[offset * totDim * totDim]));
6872: }
6873: if (hasBdPrec) {
6874: PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN_PRE, key[2], 2, Ne, chunkGeomF, chunkGeomN, u, u_t, dsAux[2], a[2], t, X_tShift, elemMatCohP));
6875: PetscCall(PetscFEIntegrateHybridJacobian(ds, dsIn, PETSCFE_JACOBIAN_PRE, key[2], 2, Nr, remGeomF, remGeomN, &u[offset * totDimIn], PetscSafePointerPlusOffset(u_t, offset * totDimIn), dsAux[2], &a[2][offset * totDimAux[2]], t, X_tShift, &elemMatCohP[offset * totDim * totDim]));
6876: }
6877: }
6878: PetscCall(PetscFEGeomRestoreChunk(geomF, offset, numCells, &remGeomF));
6879: PetscCall(PetscFEGeomRestoreChunk(geomF, 0, offset, &chunkGeomF));
6880: PetscCall(PetscFEGeomRestoreChunk(geomN, offset, numCells, &remGeomN));
6881: PetscCall(PetscFEGeomRestoreChunk(geomN, 0, offset, &chunkGeomN));
6882: }
6883: /* Insert values into matrix */
6884: for (c = cS; c < cE; ++c) {
6885: const PetscInt cell = cells ? cells[c] : c;
6886: const PetscInt cind = c - cS, coff = cind * totDim * totDim;
6887: PetscInt i, j;
6889: /* Scale element values */
6890: if (locS[0]) {
6891: PetscInt Nb, soff = cind * totDimScale[0], off = 0;
6892: PetscBool cohesive;
6894: for (fieldI = 0; fieldI < Nf; ++fieldI) {
6895: PetscCall(PetscDSGetFieldSize(ds, fieldI, &Nb));
6896: PetscCall(PetscDSGetCohesive(ds, fieldI, &cohesive));
6898: if (fieldI == key[2].field) {
6899: PetscCheck(cohesive, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Scaling should not happen for face fields");
6900: for (i = 0; i < Nb; ++i) {
6901: for (j = 0; j < totDim; ++j) elemMatCoh[coff + (off + i) * totDim + j] += s[0][soff + i] * elemMatNeg[coff + (off + i) * totDim + j] + s[1][soff + i] * elemMatPos[coff + (off + i) * totDim + j];
6902: if (hasBdPrec)
6903: for (j = 0; j < totDim; ++j) elemMatCohP[coff + (off + i) * totDim + j] += s[0][soff + i] * elemMatNegP[coff + (off + i) * totDim + j] + s[1][soff + i] * elemMatPosP[coff + (off + i) * totDim + j];
6904: }
6905: off += Nb;
6906: } else {
6907: const PetscInt N = cohesive ? Nb : Nb * 2;
6909: for (i = 0; i < N; ++i) {
6910: for (j = 0; j < totDim; ++j) elemMatCoh[coff + (off + i) * totDim + j] += elemMatNeg[coff + (off + i) * totDim + j] + elemMatPos[coff + (off + i) * totDim + j];
6911: if (hasBdPrec)
6912: for (j = 0; j < totDim; ++j) elemMatCohP[coff + (off + i) * totDim + j] += elemMatNegP[coff + (off + i) * totDim + j] + elemMatPosP[coff + (off + i) * totDim + j];
6913: }
6914: off += N;
6915: }
6916: }
6917: } else {
6918: for (i = 0; i < totDim * totDim; ++i) elemMatCoh[coff + i] += elemMatNeg[coff + i] + elemMatPos[coff + i];
6919: if (hasBdPrec)
6920: for (i = 0; i < totDim * totDim; ++i) elemMatCohP[coff + i] += elemMatNegP[coff + i] + elemMatPosP[coff + i];
6921: }
6922: if (hasBdPrec) {
6923: if (hasBdJac) {
6924: if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatCoh[cind * totDim * totDim]));
6925: PetscCall(DMPlexMatSetClosure_Internal(plex, section, globalSection, mesh->useMatClPerm, Jac, cell, &elemMatCoh[cind * totDim * totDim], ADD_VALUES));
6926: }
6927: if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatCohP[cind * totDim * totDim]));
6928: PetscCall(DMPlexMatSetClosure(plex, section, globalSection, JacP, cell, &elemMatCohP[cind * totDim * totDim], ADD_VALUES));
6929: } else if (hasBdJac) {
6930: if (mesh->printFEM > 1) PetscCall(DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatCoh[cind * totDim * totDim]));
6931: PetscCall(DMPlexMatSetClosure_Internal(plex, section, globalSection, mesh->useMatClPerm, JacP, cell, &elemMatCoh[cind * totDim * totDim], ADD_VALUES));
6932: }
6933: }
6934: }
6935: PetscCall(DMPlexRestoreCellFields(dm, cellIS, locX, locX_t, locA[2], &u, &u_t, &a[2]));
6936: PetscCall(DMPlexRestoreHybridFields(dm, dmAux, dsAux, cellIS, locA, PETSC_TRUE, a));
6937: PetscCall(DMRestoreWorkArray(dm, hasBdJac ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatNeg));
6938: PetscCall(DMRestoreWorkArray(dm, hasBdJac ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatPos));
6939: PetscCall(DMRestoreWorkArray(dm, hasBdJac ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatCoh));
6940: PetscCall(DMRestoreWorkArray(dm, hasBdPrec ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatNegP));
6941: PetscCall(DMRestoreWorkArray(dm, hasBdPrec ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatPosP));
6942: PetscCall(DMRestoreWorkArray(dm, hasBdPrec ? cellChunkSize * totDim * totDim : 0, MPIU_SCALAR, &elemMatCohP));
6943: PetscCall(PetscFree2(faces, neighbors));
6944: PetscCall(ISDestroy(&chunkISF));
6945: PetscCall(ISDestroy(&chunkISN));
6946: PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
6947: if (maxDegree <= 1) {
6948: PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, affineQuadF, PETSC_FALSE, &affineGeomF));
6949: PetscCall(PetscQuadratureDestroy(&affineQuadF));
6950: PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, affineQuadN, PETSC_FALSE, &affineGeomN));
6951: PetscCall(PetscQuadratureDestroy(&affineQuadN));
6952: } else {
6953: PetscInt f;
6954: for (f = 0; f < Nf; ++f) {
6955: if (geomsF) PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, quadsF[f], PETSC_FALSE, &geomsF[f]));
6956: if (quadsF) PetscCall(PetscQuadratureDestroy(&quadsF[f]));
6957: }
6958: PetscCall(PetscFree2(quadsF, geomsF));
6959: }
6960: if (dmAux[2]) PetscCall(DMDestroy(&plexA));
6961: PetscCall(DMDestroy(&plex));
6962: end:
6963: PetscCall(PetscLogEventEnd(DMPLEX_JacobianFEM, dm, 0, 0, 0));
6964: PetscFunctionReturn(PETSC_SUCCESS);
6965: }
6967: /*@
6968: DMPlexComputeJacobianActionByKey - Compute the local Jacobian for terms matching the input key
6970: Collective
6972: Input Parameters:
6973: + dm - The output `DM`
6974: . key - The `PetscFormKey` indicating what should be integrated
6975: . cellIS - The `IS` give a set of cells to integrate over
6976: . t - The time
6977: . X_tShift - The multiplier for the Jacobian with respect to $X_t$
6978: . locX - The local solution
6979: . locX_t - The time derivative of the local solution, or `NULL` for time-independent problems
6980: . locY - The local vector acted on by J
6981: - user - An optional user context, passed to the pointwise functions
6983: Output Parameter:
6984: . locF - The local residual F = J(X) Y
6986: Level: developer
6988: .seealso: `DMPlexComputeResidualByKey()`, `DMPlexComputeJacobianByKey()`, `DMPlexComputeResidualHybridByKey()`, `DMPlexComputeJacobianHybridByKey()`, `PetscFormKey`
6989: @*/
6990: PetscErrorCode DMPlexComputeJacobianActionByKey(DM dm, PetscFormKey key, IS cellIS, PetscReal t, PetscReal X_tShift, Vec locX, Vec locX_t, Vec locY, Vec locF, void *user)
6991: {
6992: DM_Plex *mesh = (DM_Plex *)dm->data;
6993: const char *name = "Jacobian";
6994: DM dmAux = NULL, plex, plexAux = NULL;
6995: DMEnclosureType encAux;
6996: Vec A;
6997: DMField coordField;
6998: PetscDS prob, probAux = NULL;
6999: PetscQuadrature quad;
7000: PetscSection section, globalSection, sectionAux;
7001: PetscScalar *elemMat, *elemMatD, *u, *u_t, *a = NULL, *y, *z;
7002: const PetscInt *cells;
7003: PetscInt Nf, fieldI, fieldJ;
7004: PetscInt totDim, totDimAux = 0, cStart, cEnd, numCells, c;
7005: PetscBool hasDyn;
7007: PetscFunctionBegin;
7008: PetscCall(PetscLogEventBegin(DMPLEX_JacobianFEM, dm, 0, 0, 0));
7009: PetscCall(DMConvert(dm, DMPLEX, &plex));
7010: PetscCall(ISGetLocalSize(cellIS, &numCells));
7011: PetscCall(ISGetPointRange(cellIS, &cStart, &cEnd, &cells));
7012: PetscCall(DMGetLocalSection(dm, §ion));
7013: PetscCall(DMGetGlobalSection(dm, &globalSection));
7014: PetscCall(DMGetCellDS(dm, cells ? cells[cStart] : cStart, &prob, NULL));
7015: PetscCall(PetscDSGetNumFields(prob, &Nf));
7016: PetscCall(PetscDSGetTotalDimension(prob, &totDim));
7017: PetscCall(PetscDSHasDynamicJacobian(prob, &hasDyn));
7018: hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
7019: PetscCall(DMGetAuxiliaryVec(dm, key.label, key.value, key.part, &A));
7020: if (A) {
7021: PetscCall(VecGetDM(A, &dmAux));
7022: PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
7023: PetscCall(DMConvert(dmAux, DMPLEX, &plexAux));
7024: PetscCall(DMGetLocalSection(plexAux, §ionAux));
7025: PetscCall(DMGetDS(dmAux, &probAux));
7026: PetscCall(PetscDSGetTotalDimension(probAux, &totDimAux));
7027: }
7028: PetscCall(VecSet(locF, 0.0));
7029: PetscCall(PetscMalloc6(numCells * totDim, &u, (locX_t ? (size_t)numCells * totDim : 0), &u_t, numCells * totDim * totDim, &elemMat, (hasDyn ? (size_t)numCells * totDim * totDim : 0), &elemMatD, numCells * totDim, &y, totDim, &z));
7030: if (dmAux) PetscCall(PetscMalloc1(numCells * totDimAux, &a));
7031: PetscCall(DMGetCoordinateField(dm, &coordField));
7032: for (c = cStart; c < cEnd; ++c) {
7033: const PetscInt cell = cells ? cells[c] : c;
7034: const PetscInt cind = c - cStart;
7035: PetscScalar *x = NULL, *x_t = NULL;
7036: PetscInt i;
7038: PetscCall(DMPlexVecGetClosure(plex, section, locX, cell, NULL, &x));
7039: for (i = 0; i < totDim; ++i) u[cind * totDim + i] = x[i];
7040: PetscCall(DMPlexVecRestoreClosure(plex, section, locX, cell, NULL, &x));
7041: if (locX_t) {
7042: PetscCall(DMPlexVecGetClosure(plex, section, locX_t, cell, NULL, &x_t));
7043: for (i = 0; i < totDim; ++i) u_t[cind * totDim + i] = x_t[i];
7044: PetscCall(DMPlexVecRestoreClosure(plex, section, locX_t, cell, NULL, &x_t));
7045: }
7046: if (dmAux) {
7047: PetscInt subcell;
7048: PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, cell, &subcell));
7049: PetscCall(DMPlexVecGetClosure(plexAux, sectionAux, A, subcell, NULL, &x));
7050: for (i = 0; i < totDimAux; ++i) a[cind * totDimAux + i] = x[i];
7051: PetscCall(DMPlexVecRestoreClosure(plexAux, sectionAux, A, subcell, NULL, &x));
7052: }
7053: PetscCall(DMPlexVecGetClosure(plex, section, locY, cell, NULL, &x));
7054: for (i = 0; i < totDim; ++i) y[cind * totDim + i] = x[i];
7055: PetscCall(DMPlexVecRestoreClosure(plex, section, locY, cell, NULL, &x));
7056: }
7057: PetscCall(PetscArrayzero(elemMat, numCells * totDim * totDim));
7058: if (hasDyn) PetscCall(PetscArrayzero(elemMatD, numCells * totDim * totDim));
7059: for (fieldI = 0; fieldI < Nf; ++fieldI) {
7060: PetscFE fe;
7061: PetscInt Nb;
7062: /* Conforming batches */
7063: PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
7064: /* Remainder */
7065: PetscInt Nr, offset, Nq;
7066: PetscQuadrature qGeom = NULL;
7067: PetscInt maxDegree;
7068: PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;
7070: PetscCall(PetscDSGetDiscretization(prob, fieldI, (PetscObject *)&fe));
7071: PetscCall(PetscFEGetQuadrature(fe, &quad));
7072: PetscCall(PetscFEGetDimension(fe, &Nb));
7073: PetscCall(PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches));
7074: PetscCall(DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree));
7075: if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, cellIS, &qGeom));
7076: if (!qGeom) {
7077: PetscCall(PetscFEGetQuadrature(fe, &qGeom));
7078: PetscCall(PetscObjectReference((PetscObject)qGeom));
7079: }
7080: PetscCall(PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL));
7081: PetscCall(DMSNESGetFEGeom(coordField, cellIS, qGeom, PETSC_FEGEOM_BASIC, &cgeomFEM));
7082: blockSize = Nb;
7083: batchSize = numBlocks * blockSize;
7084: PetscCall(PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches));
7085: numChunks = numCells / (numBatches * batchSize);
7086: Ne = numChunks * numBatches * batchSize;
7087: Nr = numCells % (numBatches * batchSize);
7088: offset = numCells - Nr;
7089: PetscCall(PetscFEGeomGetChunk(cgeomFEM, 0, offset, &chunkGeom));
7090: PetscCall(PetscFEGeomGetChunk(cgeomFEM, offset, numCells, &remGeom));
7091: for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
7092: key.field = fieldI * Nf + fieldJ;
7093: PetscCall(PetscFEIntegrateJacobian(prob, prob, PETSCFE_JACOBIAN, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat));
7094: PetscCall(PetscFEIntegrateJacobian(prob, prob, PETSCFE_JACOBIAN, key, Nr, remGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), probAux, PetscSafePointerPlusOffset(a, offset * totDimAux), t, X_tShift, &elemMat[offset * totDim * totDim]));
7095: if (hasDyn) {
7096: PetscCall(PetscFEIntegrateJacobian(prob, prob, PETSCFE_JACOBIAN_DYN, key, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD));
7097: PetscCall(PetscFEIntegrateJacobian(prob, prob, PETSCFE_JACOBIAN_DYN, key, Nr, remGeom, &u[offset * totDim], PetscSafePointerPlusOffset(u_t, offset * totDim), probAux, &a[offset * totDimAux], t, X_tShift, &elemMatD[offset * totDim * totDim]));
7098: }
7099: }
7100: PetscCall(PetscFEGeomRestoreChunk(cgeomFEM, offset, numCells, &remGeom));
7101: PetscCall(PetscFEGeomRestoreChunk(cgeomFEM, 0, offset, &chunkGeom));
7102: PetscCall(DMSNESRestoreFEGeom(coordField, cellIS, qGeom, PETSC_FALSE, &cgeomFEM));
7103: PetscCall(PetscQuadratureDestroy(&qGeom));
7104: }
7105: if (hasDyn) {
7106: for (c = 0; c < numCells * totDim * totDim; ++c) elemMat[c] += X_tShift * elemMatD[c];
7107: }
7108: for (c = cStart; c < cEnd; ++c) {
7109: const PetscInt cell = cells ? cells[c] : c;
7110: const PetscInt cind = c - cStart;
7111: const PetscBLASInt one = 1;
7112: PetscBLASInt M;
7113: const PetscScalar a = 1.0, b = 0.0;
7115: PetscCall(PetscBLASIntCast(totDim, &M));
7116: PetscCallBLAS("BLASgemv", BLASgemv_("N", &M, &M, &a, &elemMat[cind * totDim * totDim], &M, &y[cind * totDim], &one, &b, z, &one));
7117: if (mesh->printFEM > 1) {
7118: PetscCall(DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[cind * totDim * totDim]));
7119: PetscCall(DMPrintCellVector(c, "Y", totDim, &y[cind * totDim]));
7120: PetscCall(DMPrintCellVector(c, "Z", totDim, z));
7121: }
7122: PetscCall(DMPlexVecSetClosure(dm, section, locF, cell, z, ADD_VALUES));
7123: }
7124: PetscCall(PetscFree6(u, u_t, elemMat, elemMatD, y, z));
7125: if (mesh->printFEM) {
7126: PetscCall(PetscPrintf(PetscObjectComm((PetscObject)locF), "Z:\n"));
7127: PetscCall(VecView(locF, NULL));
7128: }
7129: PetscCall(ISRestorePointRange(cellIS, &cStart, &cEnd, &cells));
7130: PetscCall(PetscFree(a));
7131: PetscCall(DMDestroy(&plexAux));
7132: PetscCall(DMDestroy(&plex));
7133: PetscCall(PetscLogEventEnd(DMPLEX_JacobianFEM, dm, 0, 0, 0));
7134: PetscFunctionReturn(PETSC_SUCCESS);
7135: }
7137: static void f0_1(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
7138: {
7139: f0[0] = u[0];
7140: }
7142: static void f0_x(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
7143: {
7144: f0[0] = x[(int)PetscRealPart(constants[0])] * u[0];
7145: }
7147: static void f0_x2(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
7148: {
7149: PetscInt d;
7151: f0[0] = 0.0;
7152: for (d = 0; d < dim; ++d) f0[0] += PetscSqr(x[d]) * u[0];
7153: }
7155: /*@
7156: DMPlexComputeMoments - Compute the first three moments for a field
7158: Noncollective
7160: Input Parameters:
7161: + dm - the `DMPLEX`
7162: - u - the field
7164: Output Parameter:
7165: . moments - the field moments
7167: Level: intermediate
7169: Note:
7170: The `moments` array should be of length cdim + 2, where cdim is the number of components for the coordinate field.
7172: .seealso: `DM`, `DMPLEX`, `DMSwarmComputeMoments()`
7173: @*/
7174: PetscErrorCode DMPlexComputeMoments(DM dm, Vec u, PetscReal moments[])
7175: {
7176: PetscDS ds;
7177: PetscScalar mom, constants[1];
7178: const PetscScalar *oldConstants;
7179: PetscInt cdim, Nf, field = 0, Ncon;
7180: MPI_Comm comm;
7181: void *user;
7183: PetscFunctionBeginUser;
7184: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
7185: PetscCall(DMGetCoordinateDim(dm, &cdim));
7186: PetscCall(DMGetApplicationContext(dm, &user));
7187: PetscCall(DMGetDS(dm, &ds));
7188: PetscCall(PetscDSGetNumFields(ds, &Nf));
7189: PetscCall(PetscDSGetConstants(ds, &Ncon, &oldConstants));
7190: PetscCheck(Nf == 1, comm, PETSC_ERR_ARG_WRONG, "We currently only support 1 field, not %" PetscInt_FMT, Nf);
7191: PetscCall(PetscDSSetObjective(ds, field, &f0_1));
7192: PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
7193: moments[0] = PetscRealPart(mom);
7194: for (PetscInt c = 0; c < cdim; ++c) {
7195: constants[0] = c;
7196: PetscCall(PetscDSSetConstants(ds, 1, constants));
7197: PetscCall(PetscDSSetObjective(ds, field, &f0_x));
7198: PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
7199: moments[c + 1] = PetscRealPart(mom);
7200: }
7201: PetscCall(PetscDSSetObjective(ds, field, &f0_x2));
7202: PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
7203: moments[cdim + 1] = PetscRealPart(mom);
7204: PetscCall(PetscDSSetConstants(ds, Ncon, (PetscScalar *)oldConstants));
7205: PetscFunctionReturn(PETSC_SUCCESS);
7206: }