xref: /petsc/src/dm/dt/tests/ex4.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
69371c9d4SSatish Balay static PetscErrorCode CheckSymmetry(PetscInt dim, PetscInt order, PetscBool tensor) {
7c4762a1bSJed Brown   DM                   dm;
8c4762a1bSJed Brown   PetscDualSpace       sp;
9c4762a1bSJed Brown   PetscInt             nFunc, *ids, *idsCopy, *idsCopy2, i, closureSize, *closure = NULL, offset, depth;
10c4762a1bSJed Brown   DMLabel              depthLabel;
11c4762a1bSJed Brown   PetscBool            printed = PETSC_FALSE;
12c4762a1bSJed Brown   PetscScalar         *vals, *valsCopy, *valsCopy2;
13c4762a1bSJed Brown   const PetscInt      *numDofs;
14c4762a1bSJed Brown   const PetscInt    ***perms = NULL;
15c4762a1bSJed Brown   const PetscScalar ***flips = NULL;
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   PetscFunctionBegin;
189566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceCreate(PETSC_COMM_SELF, &sp));
199566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, tensor ? PETSC_FALSE : PETSC_TRUE), &dm));
209566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetType(sp, PETSCDUALSPACELAGRANGE));
219566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetDM(sp, dm));
229566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetOrder(sp, order));
239566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeSetContinuity(sp, PETSC_TRUE));
249566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeSetTensor(sp, tensor));
259566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetFromOptions(sp));
269566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetUp(sp));
279566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(sp, &nFunc));
289566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetSymmetries(sp, &perms, &flips));
29c4762a1bSJed Brown   if (!perms && !flips) {
309566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceDestroy(&sp));
319566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dm));
32c4762a1bSJed Brown     PetscFunctionReturn(0);
33c4762a1bSJed Brown   }
349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc6(nFunc, &ids, nFunc, &idsCopy, nFunc, &idsCopy2, nFunc * dim, &vals, nFunc * dim, &valsCopy, nFunc * dim, &valsCopy2));
35c4762a1bSJed Brown   for (i = 0; i < nFunc; i++) ids[i] = idsCopy2[i] = i;
36c4762a1bSJed Brown   for (i = 0; i < nFunc; i++) {
37c4762a1bSJed Brown     PetscQuadrature  q;
38c4762a1bSJed Brown     PetscInt         numPoints, Nc, j;
39c4762a1bSJed Brown     const PetscReal *points;
40c4762a1bSJed Brown     const PetscReal *weights;
41c4762a1bSJed Brown 
429566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetFunctional(sp, i, &q));
439566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(q, NULL, &Nc, &numPoints, &points, &weights));
4463a3b9bcSJacob Faibussowitsch     PetscCheck(Nc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support scalar quadrature, not %" PetscInt_FMT " components", Nc);
45c4762a1bSJed Brown     for (j = 0; j < dim; j++) vals[dim * i + j] = valsCopy2[dim * i + j] = (PetscScalar)points[j];
46c4762a1bSJed Brown   }
479566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetNumDof(sp, &numDofs));
489566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
499566063dSJacob Faibussowitsch   PetscCall(DMPlexGetTransitiveClosure(dm, 0, PETSC_TRUE, &closureSize, &closure));
509566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
51c4762a1bSJed Brown   for (i = 0, offset = 0; i < closureSize; i++, offset += numDofs[depth]) {
52b5a892a1SMatthew G. Knepley     PetscInt            point      = closure[2 * i], numFaces, j;
53c4762a1bSJed Brown     const PetscInt    **pointPerms = perms ? perms[i] : NULL;
54c4762a1bSJed Brown     const PetscScalar **pointFlips = flips ? flips[i] : NULL;
55c4762a1bSJed Brown     PetscBool           anyPrinted = PETSC_FALSE;
56c4762a1bSJed Brown 
57c4762a1bSJed Brown     if (!pointPerms && !pointFlips) continue;
589566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(depthLabel, point, &depth));
59b5a892a1SMatthew G. Knepley     {
60b5a892a1SMatthew G. Knepley       DMPolytopeType ct;
61b5a892a1SMatthew G. Knepley       /* The number of arrangements is no longer based on the number of faces */
629566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(dm, point, &ct));
63b5a892a1SMatthew G. Knepley       numFaces = DMPolytopeTypeGetNumArrangments(ct) / 2;
64b5a892a1SMatthew G. Knepley     }
65b5a892a1SMatthew G. Knepley     for (j = -numFaces; j < numFaces; j++) {
66c4762a1bSJed Brown       PetscInt           k, l;
67c4762a1bSJed Brown       const PetscInt    *perm = pointPerms ? pointPerms[j] : NULL;
68c4762a1bSJed Brown       const PetscScalar *flip = pointFlips ? pointFlips[j] : NULL;
69c4762a1bSJed Brown 
70c4762a1bSJed Brown       for (k = 0; k < numDofs[depth]; k++) {
71c4762a1bSJed Brown         PetscInt kLocal = perm ? perm[k] : k;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown         idsCopy[kLocal] = ids[offset + k];
749371c9d4SSatish Balay         for (l = 0; l < dim; l++) { valsCopy[kLocal * dim + l] = vals[(offset + k) * dim + l] * (flip ? flip[kLocal] : 1.); }
75c4762a1bSJed Brown       }
76c4762a1bSJed Brown       if (!printed && numDofs[depth] > 1) {
77c4762a1bSJed Brown         IS   is;
78c4762a1bSJed Brown         Vec  vec;
79c4762a1bSJed Brown         char name[256];
80c4762a1bSJed Brown 
81c4762a1bSJed Brown         anyPrinted = PETSC_TRUE;
8263a3b9bcSJacob 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));
839566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numDofs[depth], idsCopy, PETSC_USE_POINTER, &is));
849566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetName((PetscObject)is, name));
859566063dSJacob Faibussowitsch         PetscCall(ISView(is, PETSC_VIEWER_STDOUT_SELF));
869566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&is));
879566063dSJacob Faibussowitsch         PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, dim, numDofs[depth] * dim, valsCopy, &vec));
889566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetName((PetscObject)vec, name));
899566063dSJacob Faibussowitsch         PetscCall(VecView(vec, PETSC_VIEWER_STDOUT_SELF));
909566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&vec));
91c4762a1bSJed Brown       }
92c4762a1bSJed Brown       for (k = 0; k < numDofs[depth]; k++) {
93c4762a1bSJed Brown         PetscInt kLocal = perm ? perm[k] : k;
94c4762a1bSJed Brown 
95c4762a1bSJed Brown         idsCopy2[offset + k] = idsCopy[kLocal];
969371c9d4SSatish Balay         for (l = 0; l < dim; l++) { valsCopy2[(offset + k) * dim + l] = valsCopy[kLocal * dim + l] * (flip ? PetscConj(flip[kLocal]) : 1.); }
97c4762a1bSJed Brown       }
98c4762a1bSJed Brown       for (k = 0; k < nFunc; k++) {
9963a3b9bcSJacob 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);
100c4762a1bSJed Brown         for (l = 0; l < dim; l++) {
10163a3b9bcSJacob 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]));
102c4762a1bSJed Brown         }
103c4762a1bSJed Brown       }
104c4762a1bSJed Brown     }
105c4762a1bSJed Brown     if (anyPrinted && !printed) printed = PETSC_TRUE;
106c4762a1bSJed Brown   }
1079566063dSJacob Faibussowitsch   PetscCall(DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, &closureSize, &closure));
1089566063dSJacob Faibussowitsch   PetscCall(PetscFree6(ids, idsCopy, idsCopy2, vals, valsCopy, valsCopy2));
1099566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDestroy(&sp));
1109566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
111c4762a1bSJed Brown   PetscFunctionReturn(0);
112c4762a1bSJed Brown }
113c4762a1bSJed Brown 
1149371c9d4SSatish Balay int main(int argc, char **argv) {
115c4762a1bSJed Brown   PetscInt dim, order, tensor;
116c4762a1bSJed Brown 
117327415f7SBarry Smith   PetscFunctionBeginUser;
1189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
119c4762a1bSJed Brown   for (tensor = 0; tensor < 2; tensor++) {
120c4762a1bSJed Brown     for (dim = 1; dim <= 3; dim++) {
121c4762a1bSJed Brown       if (dim == 1 && tensor) continue;
122*48a46eb9SPierre Jolivet       for (order = 0; order <= (tensor ? 5 : 6); order++) PetscCall(CheckSymmetry(dim, order, tensor ? PETSC_TRUE : PETSC_FALSE));
123c4762a1bSJed Brown     }
124c4762a1bSJed Brown   }
1259566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
126b122ec5aSJacob Faibussowitsch   return 0;
127c4762a1bSJed Brown }
128c4762a1bSJed Brown 
129c4762a1bSJed Brown /*TEST
130c4762a1bSJed Brown   test:
131c4762a1bSJed Brown     suffix: 0
132c4762a1bSJed Brown TEST*/
133