1c4762a1bSJed Brown static char help[] = "Tests dual space symmetry.\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscfe.h>
4c4762a1bSJed Brown #include <petscdmplex.h>
5c4762a1bSJed Brown
CheckSymmetry(PetscInt dim,PetscInt order,PetscBool tensor)6d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckSymmetry(PetscInt dim, PetscInt order, PetscBool tensor)
7d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown DM dm;
9c4762a1bSJed Brown PetscDualSpace sp;
10c4762a1bSJed Brown PetscInt nFunc, *ids, *idsCopy, *idsCopy2, i, closureSize, *closure = NULL, offset, depth;
11c4762a1bSJed Brown DMLabel depthLabel;
12c4762a1bSJed Brown PetscBool printed = PETSC_FALSE;
13c4762a1bSJed Brown PetscScalar *vals, *valsCopy, *valsCopy2;
14c4762a1bSJed Brown const PetscInt *numDofs;
15c4762a1bSJed Brown const PetscInt ***perms = NULL;
16c4762a1bSJed Brown const PetscScalar ***flips = NULL;
17c4762a1bSJed Brown
18c4762a1bSJed Brown PetscFunctionBegin;
199566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreate(PETSC_COMM_SELF, &sp));
209566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, tensor ? PETSC_FALSE : PETSC_TRUE), &dm));
219566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetType(sp, PETSCDUALSPACELAGRANGE));
229566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetDM(sp, dm));
239566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetOrder(sp, order));
249566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetContinuity(sp, PETSC_TRUE));
259566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetTensor(sp, tensor));
269566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetFromOptions(sp));
279566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetUp(sp));
289566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(sp, &nFunc));
299566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetSymmetries(sp, &perms, &flips));
30c4762a1bSJed Brown if (!perms && !flips) {
319566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&sp));
329566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
34c4762a1bSJed Brown }
359566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(nFunc, &ids, nFunc, &idsCopy, nFunc, &idsCopy2, nFunc * dim, &vals, nFunc * dim, &valsCopy, nFunc * dim, &valsCopy2));
36c4762a1bSJed Brown for (i = 0; i < nFunc; i++) ids[i] = idsCopy2[i] = i;
37c4762a1bSJed Brown for (i = 0; i < nFunc; i++) {
38c4762a1bSJed Brown PetscQuadrature q;
39c4762a1bSJed Brown PetscInt numPoints, Nc, j;
40c4762a1bSJed Brown const PetscReal *points;
41c4762a1bSJed Brown const PetscReal *weights;
42c4762a1bSJed Brown
439566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetFunctional(sp, i, &q));
449566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(q, NULL, &Nc, &numPoints, &points, &weights));
4563a3b9bcSJacob Faibussowitsch PetscCheck(Nc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support scalar quadrature, not %" PetscInt_FMT " components", Nc);
46c4762a1bSJed Brown for (j = 0; j < dim; j++) vals[dim * i + j] = valsCopy2[dim * i + j] = (PetscScalar)points[j];
47c4762a1bSJed Brown }
489566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetNumDof(sp, &numDofs));
499566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth));
509566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, 0, PETSC_TRUE, &closureSize, &closure));
519566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
52c4762a1bSJed Brown for (i = 0, offset = 0; i < closureSize; i++, offset += numDofs[depth]) {
53b5a892a1SMatthew G. Knepley PetscInt point = closure[2 * i], numFaces, j;
54c4762a1bSJed Brown const PetscInt **pointPerms = perms ? perms[i] : NULL;
55c4762a1bSJed Brown const PetscScalar **pointFlips = flips ? flips[i] : NULL;
56c4762a1bSJed Brown PetscBool anyPrinted = PETSC_FALSE;
57c4762a1bSJed Brown
58c4762a1bSJed Brown if (!pointPerms && !pointFlips) continue;
599566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(depthLabel, point, &depth));
60b5a892a1SMatthew G. Knepley {
61b5a892a1SMatthew G. Knepley DMPolytopeType ct;
62b5a892a1SMatthew G. Knepley /* The number of arrangements is no longer based on the number of faces */
639566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, point, &ct));
64*85036b15SMatthew G. Knepley numFaces = DMPolytopeTypeGetNumArrangements(ct) / 2;
65b5a892a1SMatthew G. Knepley }
66b5a892a1SMatthew G. Knepley for (j = -numFaces; j < numFaces; j++) {
67c4762a1bSJed Brown PetscInt k, l;
68c4762a1bSJed Brown const PetscInt *perm = pointPerms ? pointPerms[j] : NULL;
69c4762a1bSJed Brown const PetscScalar *flip = pointFlips ? pointFlips[j] : NULL;
70c4762a1bSJed Brown
71c4762a1bSJed Brown for (k = 0; k < numDofs[depth]; k++) {
72c4762a1bSJed Brown PetscInt kLocal = perm ? perm[k] : k;
73c4762a1bSJed Brown
74c4762a1bSJed Brown idsCopy[kLocal] = ids[offset + k];
75ad540459SPierre Jolivet for (l = 0; l < dim; l++) valsCopy[kLocal * dim + l] = vals[(offset + k) * dim + l] * (flip ? flip[kLocal] : 1.);
76c4762a1bSJed Brown }
77c4762a1bSJed Brown if (!printed && numDofs[depth] > 1) {
78c4762a1bSJed Brown IS is;
79c4762a1bSJed Brown Vec vec;
80c4762a1bSJed Brown char name[256];
81c4762a1bSJed Brown
82c4762a1bSJed Brown anyPrinted = PETSC_TRUE;
8363a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, 256, "%" PetscInt_FMT "D, %s, Order %" PetscInt_FMT ", Point %" PetscInt_FMT " Symmetry %" PetscInt_FMT, dim, tensor ? "Tensor" : "Simplex", order, point, j));
849566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numDofs[depth], idsCopy, PETSC_USE_POINTER, &is));
859566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)is, name));
869566063dSJacob Faibussowitsch PetscCall(ISView(is, PETSC_VIEWER_STDOUT_SELF));
879566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is));
889566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, dim, numDofs[depth] * dim, valsCopy, &vec));
899566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)vec, name));
909566063dSJacob Faibussowitsch PetscCall(VecView(vec, PETSC_VIEWER_STDOUT_SELF));
919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec));
92c4762a1bSJed Brown }
93c4762a1bSJed Brown for (k = 0; k < numDofs[depth]; k++) {
94c4762a1bSJed Brown PetscInt kLocal = perm ? perm[k] : k;
95c4762a1bSJed Brown
96c4762a1bSJed Brown idsCopy2[offset + k] = idsCopy[kLocal];
97ad540459SPierre Jolivet for (l = 0; l < dim; l++) valsCopy2[(offset + k) * dim + l] = valsCopy[kLocal * dim + l] * (flip ? PetscConj(flip[kLocal]) : 1.);
98c4762a1bSJed Brown }
99c4762a1bSJed Brown for (k = 0; k < nFunc; k++) {
10063a3b9bcSJacob Faibussowitsch PetscCheck(idsCopy2[k] == ids[k], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Symmetry failure: %" PetscInt_FMT "D, %s, point %" PetscInt_FMT ", symmetry %" PetscInt_FMT ", order %" PetscInt_FMT ", functional %" PetscInt_FMT ": (%" PetscInt_FMT " != %" PetscInt_FMT ")", dim, tensor ? "Tensor" : "Simplex", point, j, order, k, ids[k], k);
101c4762a1bSJed Brown for (l = 0; l < dim; l++) {
10263a3b9bcSJacob Faibussowitsch PetscCheck(valsCopy2[dim * k + l] == vals[dim * k + l], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Symmetry failure: %" PetscInt_FMT "D, %s, point %" PetscInt_FMT ", symmetry %" PetscInt_FMT ", order %" PetscInt_FMT ", functional %" PetscInt_FMT ", component %" PetscInt_FMT ": (%g != %g)", dim, tensor ? "Tensor" : "Simplex", point, j, order, k, l, (double)PetscAbsScalar(valsCopy2[dim * k + l]), (double)PetscAbsScalar(vals[dim * k + l]));
103c4762a1bSJed Brown }
104c4762a1bSJed Brown }
105c4762a1bSJed Brown }
106c4762a1bSJed Brown if (anyPrinted && !printed) printed = PETSC_TRUE;
107c4762a1bSJed Brown }
1089566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, &closureSize, &closure));
1099566063dSJacob Faibussowitsch PetscCall(PetscFree6(ids, idsCopy, idsCopy2, vals, valsCopy, valsCopy2));
1109566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&sp));
1119566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
1123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
113c4762a1bSJed Brown }
114c4762a1bSJed Brown
main(int argc,char ** argv)115d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
116d71ae5a4SJacob Faibussowitsch {
117c4762a1bSJed Brown PetscInt dim, order, tensor;
118c4762a1bSJed Brown
119327415f7SBarry Smith PetscFunctionBeginUser;
1209566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
121c4762a1bSJed Brown for (tensor = 0; tensor < 2; tensor++) {
122c4762a1bSJed Brown for (dim = 1; dim <= 3; dim++) {
123c4762a1bSJed Brown if (dim == 1 && tensor) continue;
12448a46eb9SPierre Jolivet for (order = 0; order <= (tensor ? 5 : 6); order++) PetscCall(CheckSymmetry(dim, order, tensor ? PETSC_TRUE : PETSC_FALSE));
125c4762a1bSJed Brown }
126c4762a1bSJed Brown }
1279566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
128b122ec5aSJacob Faibussowitsch return 0;
129c4762a1bSJed Brown }
130c4762a1bSJed Brown
131c4762a1bSJed Brown /*TEST
132c4762a1bSJed Brown test:
133c4762a1bSJed Brown suffix: 0
134c4762a1bSJed Brown TEST*/
135