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