xref: /petsc/src/dm/dt/tests/ex16.c (revision 2dce792e531186164765a9583d36d03ffc15e9ea)
1*2dce792eSToby Isaac const char help[] = "Test PETSCFEVECTOR";
2*2dce792eSToby Isaac 
3*2dce792eSToby Isaac #include <petscfe.h>
4*2dce792eSToby Isaac 
5*2dce792eSToby Isaac static PetscErrorCode PetscFEVectorTest(PetscFE orig_fe, PetscInt n_copies, PetscBool interleave_basis, PetscBool interleave_components)
6*2dce792eSToby Isaac {
7*2dce792eSToby Isaac   PetscFE          vec_fe, dup_fe;
8*2dce792eSToby Isaac   PetscQuadrature  quad;
9*2dce792eSToby Isaac   PetscInt         num_points;
10*2dce792eSToby Isaac   const PetscReal *points;
11*2dce792eSToby Isaac   PetscViewer      viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)orig_fe));
12*2dce792eSToby Isaac   PetscTabulation  orig_T, vec_T, dup_T;
13*2dce792eSToby Isaac   PetscSpace       space;
14*2dce792eSToby Isaac   PetscInt         Nb, vNb, vNb_s, vNb_d, Nc, vNc, cdim;
15*2dce792eSToby Isaac   PetscDualSpace   dual_space, dup_dual_space;
16*2dce792eSToby Isaac   PetscBool        ib_s, ic_s, ib_d, ic_d;
17*2dce792eSToby Isaac 
18*2dce792eSToby Isaac   PetscFunctionBegin;
19*2dce792eSToby Isaac   PetscCall(PetscFEGetQuadrature(orig_fe, &quad));
20*2dce792eSToby Isaac   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &num_points, &points, NULL));
21*2dce792eSToby Isaac   PetscCall(PetscFECreateVector(orig_fe, n_copies, interleave_basis, interleave_components, &vec_fe));
22*2dce792eSToby Isaac   PetscCall(PetscFEGetBasisSpace(vec_fe, &space));
23*2dce792eSToby Isaac   PetscCall(PetscFEGetDualSpace(vec_fe, &dual_space));
24*2dce792eSToby Isaac   PetscCall(PetscObjectSetName((PetscObject)vec_fe, "vector fe"));
25*2dce792eSToby Isaac   PetscCall(PetscObjectSetName((PetscObject)space, "vector basis space"));
26*2dce792eSToby Isaac   PetscCall(PetscObjectSetName((PetscObject)dual_space, "vector dual space"));
27*2dce792eSToby Isaac   PetscCall(PetscFEView(vec_fe, viewer));
28*2dce792eSToby Isaac   PetscCall(PetscFECreateTabulation(orig_fe, 1, num_points, points, 1, &orig_T));
29*2dce792eSToby Isaac   PetscCall(PetscFECreateTabulation(vec_fe, 1, num_points, points, 1, &vec_T));
30*2dce792eSToby Isaac   PetscCall(PetscFEGetDimension(orig_fe, &Nb));
31*2dce792eSToby Isaac   PetscCall(PetscFEGetDimension(vec_fe, &vNb));
32*2dce792eSToby Isaac   PetscCall(PetscFEGetNumComponents(orig_fe, &Nc));
33*2dce792eSToby Isaac   PetscCall(PetscFEGetNumComponents(vec_fe, &vNc));
34*2dce792eSToby Isaac   PetscCall(PetscFEGetSpatialDimension(orig_fe, &cdim));
35*2dce792eSToby Isaac   {
36*2dce792eSToby Isaac     PetscInt *pre_image;
37*2dce792eSToby Isaac     PetscInt  c_stride = interleave_components ? n_copies : 1;
38*2dce792eSToby Isaac     PetscInt  c_incr   = interleave_components ? 1 : Nc;
39*2dce792eSToby Isaac 
40*2dce792eSToby Isaac     PetscCall(PetscMalloc1(vNb, &pre_image));
41*2dce792eSToby Isaac     for (PetscInt e = 0; e < vNb; e++) pre_image[e] = -1;
42*2dce792eSToby Isaac     for (PetscInt copy = 0, coffset = 0; copy < n_copies; copy++, coffset += c_incr) {
43*2dce792eSToby Isaac       for (PetscInt b = 0; b < Nb; b++) {
44*2dce792eSToby Isaac         for (PetscInt e = 0; e < vNb; e++) {
45*2dce792eSToby Isaac           PetscReal err = 0.0;
46*2dce792eSToby Isaac 
47*2dce792eSToby Isaac           for (PetscInt k = 0; k <= orig_T->K; k++) {
48*2dce792eSToby Isaac             const PetscReal *s_Tk   = orig_T->T[k];
49*2dce792eSToby Isaac             const PetscReal *v_Tk   = vec_T->T[k];
50*2dce792eSToby Isaac             PetscInt         dblock = PetscPowInt(cdim, k);
51*2dce792eSToby Isaac 
52*2dce792eSToby Isaac             for (PetscInt p = 0; p < num_points; p++) {
53*2dce792eSToby Isaac               const PetscReal *s_Tp = &s_Tk[(p * Nb + b) * Nc * dblock];
54*2dce792eSToby Isaac               const PetscReal *v_Tp = &v_Tk[(p * vNb + e) * vNc * dblock];
55*2dce792eSToby Isaac               for (PetscInt c = 0; c < Nc; c++) {
56*2dce792eSToby Isaac                 PetscInt         vc   = coffset + c * c_stride;
57*2dce792eSToby Isaac                 const PetscReal *s_Tc = &s_Tp[c * dblock];
58*2dce792eSToby Isaac                 const PetscReal *v_Tc = &v_Tp[vc * dblock];
59*2dce792eSToby Isaac                 for (PetscInt d = 0; d < PetscPowInt(cdim, k); d++) err = PetscMax(err, PetscAbsReal(s_Tc[d] - v_Tc[d]));
60*2dce792eSToby Isaac               }
61*2dce792eSToby Isaac             }
62*2dce792eSToby Isaac           }
63*2dce792eSToby Isaac           if (err < PETSC_SMALL) {
64*2dce792eSToby Isaac             PetscCheck(pre_image[e] == -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Original basis %d and %d both match to vector basis %d\n", (int)pre_image[e], (int)b, (int)e);
65*2dce792eSToby Isaac             pre_image[e] = b;
66*2dce792eSToby Isaac           }
67*2dce792eSToby Isaac         }
68*2dce792eSToby Isaac       }
69*2dce792eSToby Isaac     }
70*2dce792eSToby Isaac     for (PetscInt e = 0; e < vNb; e++) PetscCheck(pre_image[e] >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "No original basis matched to %d\n", (int)e);
71*2dce792eSToby Isaac     PetscCall(PetscViewerASCIIPrintf(viewer, "Vector basis to original basis:"));
72*2dce792eSToby Isaac     for (PetscInt e = 0; e < vNb; e++) {
73*2dce792eSToby Isaac       if (!(e % 16)) PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
74*2dce792eSToby Isaac       PetscCall(PetscViewerASCIIPrintf(viewer, " %3d", (int)pre_image[e]));
75*2dce792eSToby Isaac     }
76*2dce792eSToby Isaac     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
77*2dce792eSToby Isaac     PetscCall(PetscFree(pre_image));
78*2dce792eSToby Isaac   }
79*2dce792eSToby Isaac   PetscCall(PetscSpaceSumGetInterleave(space, &ib_s, &ic_s));
80*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetInterleave(dual_space, &ib_d, &ic_d));
81*2dce792eSToby Isaac   PetscCheck(ib_s == interleave_basis, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interleave basis of space does not match");
82*2dce792eSToby Isaac   PetscCheck(ic_s == interleave_components, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interleave components of space does not match");
83*2dce792eSToby Isaac   PetscCheck(ib_d == interleave_basis, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interleave basis of dual space does not match");
84*2dce792eSToby Isaac   PetscCheck(ic_d == interleave_components, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interleave components of dual space does not match");
85*2dce792eSToby Isaac   PetscCall(PetscSpaceGetDimension(space, &vNb_s));
86*2dce792eSToby Isaac   PetscCall(PetscDualSpaceGetDimension(dual_space, &vNb_d));
87*2dce792eSToby Isaac   PetscCheck(vNb_s == vNb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dimension of space does not match");
88*2dce792eSToby Isaac   PetscCheck(vNb_s == vNb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dimension of dual space does not match");
89*2dce792eSToby Isaac   PetscCall(PetscObjectReference((PetscObject)space));
90*2dce792eSToby Isaac   PetscCall(PetscDualSpaceDuplicate(dual_space, &dup_dual_space)); // not necessary just testing interface
91*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSetUp(dup_dual_space));
92*2dce792eSToby Isaac   PetscCall(PetscFECreateFromSpaces(space, dup_dual_space, NULL, NULL, &dup_fe));
93*2dce792eSToby Isaac   PetscCall(PetscFECreateTabulation(dup_fe, 1, num_points, points, 1, &dup_T));
94*2dce792eSToby Isaac   {
95*2dce792eSToby Isaac     PetscReal err = 0.0;
96*2dce792eSToby Isaac 
97*2dce792eSToby Isaac     for (PetscInt k = 0; k <= vec_T->K; k++) {
98*2dce792eSToby Isaac       PetscInt dblock = PetscPowInt(cdim, k);
99*2dce792eSToby Isaac       PetscInt size   = num_points * vNb * vNc * dblock;
100*2dce792eSToby Isaac       for (PetscInt i = 0; i < size; i++) err = PetscMax(err, PetscAbsReal(vec_T->T[k][i] - dup_T->T[k][i]));
101*2dce792eSToby Isaac     }
102*2dce792eSToby Isaac     PetscCheck(err < PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error between direct tabulation and indirect tabulation: %g\n", (double)err);
103*2dce792eSToby Isaac   }
104*2dce792eSToby Isaac   PetscCall(PetscTabulationDestroy(&dup_T));
105*2dce792eSToby Isaac   PetscCall(PetscTabulationDestroy(&vec_T));
106*2dce792eSToby Isaac   PetscCall(PetscTabulationDestroy(&orig_T));
107*2dce792eSToby Isaac   PetscCall(PetscFEDestroy(&dup_fe));
108*2dce792eSToby Isaac   PetscCall(PetscFEDestroy(&vec_fe));
109*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
110*2dce792eSToby Isaac }
111*2dce792eSToby Isaac 
112*2dce792eSToby Isaac int main(int argc, char **argv)
113*2dce792eSToby Isaac {
114*2dce792eSToby Isaac   PetscFE     scalar, vector;
115*2dce792eSToby Isaac   PetscViewer viewer;
116*2dce792eSToby Isaac 
117*2dce792eSToby Isaac   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
118*2dce792eSToby Isaac   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, 3, 1, PETSC_TRUE, 3, PETSC_DETERMINE, &scalar));
119*2dce792eSToby Isaac   viewer = PETSC_VIEWER_STDOUT_SELF;
120*2dce792eSToby Isaac   PetscCall(PetscObjectSetName((PetscObject)scalar, "base FE (scalar)"));
121*2dce792eSToby Isaac   PetscCall(PetscFEView(scalar, viewer));
122*2dce792eSToby Isaac   PetscCall(PetscViewerASCIIPushTab(viewer));
123*2dce792eSToby Isaac   for (PetscInt n_copies = 1; n_copies <= 3; n_copies++) {
124*2dce792eSToby Isaac     PetscCall(PetscFEVectorTest(scalar, n_copies, PETSC_FALSE, PETSC_FALSE));
125*2dce792eSToby Isaac     PetscCall(PetscFEVectorTest(scalar, n_copies, PETSC_FALSE, PETSC_TRUE));
126*2dce792eSToby Isaac     PetscCall(PetscFEVectorTest(scalar, n_copies, PETSC_TRUE, PETSC_FALSE));
127*2dce792eSToby Isaac     PetscCall(PetscFEVectorTest(scalar, n_copies, PETSC_TRUE, PETSC_TRUE));
128*2dce792eSToby Isaac   }
129*2dce792eSToby Isaac   PetscCall(PetscViewerASCIIPopTab(viewer));
130*2dce792eSToby Isaac   PetscCall(PetscFEDestroy(&scalar));
131*2dce792eSToby Isaac   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, 3, 3, PETSC_TRUE, 3, PETSC_DETERMINE, &vector));
132*2dce792eSToby Isaac   PetscCall(PetscObjectSetName((PetscObject)vector, "base FE (vector)"));
133*2dce792eSToby Isaac   PetscCall(PetscFEView(vector, viewer));
134*2dce792eSToby Isaac   PetscCall(PetscViewerASCIIPushTab(viewer));
135*2dce792eSToby Isaac   for (PetscInt n_copies = 1; n_copies <= 3; n_copies++) {
136*2dce792eSToby Isaac     PetscCall(PetscFEVectorTest(vector, n_copies, PETSC_FALSE, PETSC_FALSE));
137*2dce792eSToby Isaac     PetscCall(PetscFEVectorTest(vector, n_copies, PETSC_FALSE, PETSC_TRUE));
138*2dce792eSToby Isaac     PetscCall(PetscFEVectorTest(vector, n_copies, PETSC_TRUE, PETSC_FALSE));
139*2dce792eSToby Isaac     PetscCall(PetscFEVectorTest(vector, n_copies, PETSC_TRUE, PETSC_TRUE));
140*2dce792eSToby Isaac   }
141*2dce792eSToby Isaac   PetscCall(PetscViewerASCIIPopTab(viewer));
142*2dce792eSToby Isaac   PetscCall(PetscFEDestroy(&vector));
143*2dce792eSToby Isaac   PetscCall(PetscFinalize());
144*2dce792eSToby Isaac   return 0;
145*2dce792eSToby Isaac }
146*2dce792eSToby Isaac 
147*2dce792eSToby Isaac /*TEST
148*2dce792eSToby Isaac 
149*2dce792eSToby Isaac   test:
150*2dce792eSToby Isaac     suffix: 0
151*2dce792eSToby Isaac 
152*2dce792eSToby Isaac TEST*/
153