12dce792eSToby Isaac const char help[] = "Test PETSCFEVECTOR";
22dce792eSToby Isaac
32dce792eSToby Isaac #include <petscfe.h>
42dce792eSToby Isaac
PetscFEVectorTest(PetscFE orig_fe,PetscInt n_copies,PetscBool interleave_basis,PetscBool interleave_components)52dce792eSToby Isaac static PetscErrorCode PetscFEVectorTest(PetscFE orig_fe, PetscInt n_copies, PetscBool interleave_basis, PetscBool interleave_components)
62dce792eSToby Isaac {
72dce792eSToby Isaac PetscFE vec_fe, dup_fe;
82dce792eSToby Isaac PetscQuadrature quad;
92dce792eSToby Isaac PetscInt num_points;
102dce792eSToby Isaac const PetscReal *points;
112dce792eSToby Isaac PetscViewer viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)orig_fe));
122dce792eSToby Isaac PetscTabulation orig_T, vec_T, dup_T;
132dce792eSToby Isaac PetscSpace space;
142dce792eSToby Isaac PetscInt Nb, vNb, vNb_s, vNb_d, Nc, vNc, cdim;
152dce792eSToby Isaac PetscDualSpace dual_space, dup_dual_space;
162dce792eSToby Isaac PetscBool ib_s, ic_s, ib_d, ic_d;
172dce792eSToby Isaac
182dce792eSToby Isaac PetscFunctionBegin;
192dce792eSToby Isaac PetscCall(PetscFEGetQuadrature(orig_fe, &quad));
202dce792eSToby Isaac PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &num_points, &points, NULL));
212dce792eSToby Isaac PetscCall(PetscFECreateVector(orig_fe, n_copies, interleave_basis, interleave_components, &vec_fe));
222dce792eSToby Isaac PetscCall(PetscFEGetBasisSpace(vec_fe, &space));
232dce792eSToby Isaac PetscCall(PetscFEGetDualSpace(vec_fe, &dual_space));
242dce792eSToby Isaac PetscCall(PetscObjectSetName((PetscObject)vec_fe, "vector fe"));
252dce792eSToby Isaac PetscCall(PetscObjectSetName((PetscObject)space, "vector basis space"));
262dce792eSToby Isaac PetscCall(PetscObjectSetName((PetscObject)dual_space, "vector dual space"));
272dce792eSToby Isaac PetscCall(PetscFEView(vec_fe, viewer));
282dce792eSToby Isaac PetscCall(PetscFECreateTabulation(orig_fe, 1, num_points, points, 1, &orig_T));
292dce792eSToby Isaac PetscCall(PetscFECreateTabulation(vec_fe, 1, num_points, points, 1, &vec_T));
302dce792eSToby Isaac PetscCall(PetscFEGetDimension(orig_fe, &Nb));
312dce792eSToby Isaac PetscCall(PetscFEGetDimension(vec_fe, &vNb));
322dce792eSToby Isaac PetscCall(PetscFEGetNumComponents(orig_fe, &Nc));
332dce792eSToby Isaac PetscCall(PetscFEGetNumComponents(vec_fe, &vNc));
342dce792eSToby Isaac PetscCall(PetscFEGetSpatialDimension(orig_fe, &cdim));
352dce792eSToby Isaac {
362dce792eSToby Isaac PetscInt *pre_image;
372dce792eSToby Isaac PetscInt c_stride = interleave_components ? n_copies : 1;
382dce792eSToby Isaac PetscInt c_incr = interleave_components ? 1 : Nc;
392dce792eSToby Isaac
402dce792eSToby Isaac PetscCall(PetscMalloc1(vNb, &pre_image));
412dce792eSToby Isaac for (PetscInt e = 0; e < vNb; e++) pre_image[e] = -1;
422dce792eSToby Isaac for (PetscInt copy = 0, coffset = 0; copy < n_copies; copy++, coffset += c_incr) {
432dce792eSToby Isaac for (PetscInt b = 0; b < Nb; b++) {
442dce792eSToby Isaac for (PetscInt e = 0; e < vNb; e++) {
452dce792eSToby Isaac PetscReal err = 0.0;
462dce792eSToby Isaac
472dce792eSToby Isaac for (PetscInt k = 0; k <= orig_T->K; k++) {
482dce792eSToby Isaac const PetscReal *s_Tk = orig_T->T[k];
492dce792eSToby Isaac const PetscReal *v_Tk = vec_T->T[k];
502dce792eSToby Isaac PetscInt dblock = PetscPowInt(cdim, k);
512dce792eSToby Isaac
522dce792eSToby Isaac for (PetscInt p = 0; p < num_points; p++) {
532dce792eSToby Isaac const PetscReal *s_Tp = &s_Tk[(p * Nb + b) * Nc * dblock];
542dce792eSToby Isaac const PetscReal *v_Tp = &v_Tk[(p * vNb + e) * vNc * dblock];
552dce792eSToby Isaac for (PetscInt c = 0; c < Nc; c++) {
562dce792eSToby Isaac PetscInt vc = coffset + c * c_stride;
572dce792eSToby Isaac const PetscReal *s_Tc = &s_Tp[c * dblock];
582dce792eSToby Isaac const PetscReal *v_Tc = &v_Tp[vc * dblock];
592dce792eSToby Isaac for (PetscInt d = 0; d < PetscPowInt(cdim, k); d++) err = PetscMax(err, PetscAbsReal(s_Tc[d] - v_Tc[d]));
602dce792eSToby Isaac }
612dce792eSToby Isaac }
622dce792eSToby Isaac }
632dce792eSToby Isaac if (err < PETSC_SMALL) {
64*300f1712SStefano Zampini PetscCheck(pre_image[e] == -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Original basis %" PetscInt_FMT " and %" PetscInt_FMT " both match to vector basis %" PetscInt_FMT, pre_image[e], b, e);
652dce792eSToby Isaac pre_image[e] = b;
662dce792eSToby Isaac }
672dce792eSToby Isaac }
682dce792eSToby Isaac }
692dce792eSToby Isaac }
70*300f1712SStefano Zampini for (PetscInt e = 0; e < vNb; e++) PetscCheck(pre_image[e] >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "No original basis matched to %" PetscInt_FMT, e);
712dce792eSToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "Vector basis to original basis:"));
722dce792eSToby Isaac for (PetscInt e = 0; e < vNb; e++) {
732dce792eSToby Isaac if (!(e % 16)) PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
74*300f1712SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " %3" PetscInt_FMT, pre_image[e]));
752dce792eSToby Isaac }
762dce792eSToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
772dce792eSToby Isaac PetscCall(PetscFree(pre_image));
782dce792eSToby Isaac }
792dce792eSToby Isaac PetscCall(PetscSpaceSumGetInterleave(space, &ib_s, &ic_s));
802dce792eSToby Isaac PetscCall(PetscDualSpaceSumGetInterleave(dual_space, &ib_d, &ic_d));
812dce792eSToby Isaac PetscCheck(ib_s == interleave_basis, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interleave basis of space does not match");
822dce792eSToby Isaac PetscCheck(ic_s == interleave_components, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interleave components of space does not match");
832dce792eSToby Isaac PetscCheck(ib_d == interleave_basis, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interleave basis of dual space does not match");
842dce792eSToby Isaac PetscCheck(ic_d == interleave_components, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interleave components of dual space does not match");
852dce792eSToby Isaac PetscCall(PetscSpaceGetDimension(space, &vNb_s));
862dce792eSToby Isaac PetscCall(PetscDualSpaceGetDimension(dual_space, &vNb_d));
872dce792eSToby Isaac PetscCheck(vNb_s == vNb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dimension of space does not match");
882dce792eSToby Isaac PetscCheck(vNb_s == vNb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dimension of dual space does not match");
892dce792eSToby Isaac PetscCall(PetscObjectReference((PetscObject)space));
902dce792eSToby Isaac PetscCall(PetscDualSpaceDuplicate(dual_space, &dup_dual_space)); // not necessary just testing interface
912dce792eSToby Isaac PetscCall(PetscDualSpaceSetUp(dup_dual_space));
922dce792eSToby Isaac PetscCall(PetscFECreateFromSpaces(space, dup_dual_space, NULL, NULL, &dup_fe));
932dce792eSToby Isaac PetscCall(PetscFECreateTabulation(dup_fe, 1, num_points, points, 1, &dup_T));
942dce792eSToby Isaac {
952dce792eSToby Isaac PetscReal err = 0.0;
962dce792eSToby Isaac
972dce792eSToby Isaac for (PetscInt k = 0; k <= vec_T->K; k++) {
982dce792eSToby Isaac PetscInt dblock = PetscPowInt(cdim, k);
992dce792eSToby Isaac PetscInt size = num_points * vNb * vNc * dblock;
1002dce792eSToby Isaac for (PetscInt i = 0; i < size; i++) err = PetscMax(err, PetscAbsReal(vec_T->T[k][i] - dup_T->T[k][i]));
1012dce792eSToby Isaac }
10200045ab3SPierre Jolivet PetscCheck(err < PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error between direct tabulation and indirect tabulation: %g", (double)err);
1032dce792eSToby Isaac }
1042dce792eSToby Isaac PetscCall(PetscTabulationDestroy(&dup_T));
1052dce792eSToby Isaac PetscCall(PetscTabulationDestroy(&vec_T));
1062dce792eSToby Isaac PetscCall(PetscTabulationDestroy(&orig_T));
1072dce792eSToby Isaac PetscCall(PetscFEDestroy(&dup_fe));
1082dce792eSToby Isaac PetscCall(PetscFEDestroy(&vec_fe));
1092dce792eSToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
1102dce792eSToby Isaac }
1112dce792eSToby Isaac
main(int argc,char ** argv)1122dce792eSToby Isaac int main(int argc, char **argv)
1132dce792eSToby Isaac {
1142dce792eSToby Isaac PetscFE scalar, vector;
1152dce792eSToby Isaac PetscViewer viewer;
1162dce792eSToby Isaac
1172dce792eSToby Isaac PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1182dce792eSToby Isaac PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, 3, 1, PETSC_TRUE, 3, PETSC_DETERMINE, &scalar));
1192dce792eSToby Isaac viewer = PETSC_VIEWER_STDOUT_SELF;
1202dce792eSToby Isaac PetscCall(PetscObjectSetName((PetscObject)scalar, "base FE (scalar)"));
1212dce792eSToby Isaac PetscCall(PetscFEView(scalar, viewer));
1222dce792eSToby Isaac PetscCall(PetscViewerASCIIPushTab(viewer));
1232dce792eSToby Isaac for (PetscInt n_copies = 1; n_copies <= 3; n_copies++) {
1242dce792eSToby Isaac PetscCall(PetscFEVectorTest(scalar, n_copies, PETSC_FALSE, PETSC_FALSE));
1252dce792eSToby Isaac PetscCall(PetscFEVectorTest(scalar, n_copies, PETSC_FALSE, PETSC_TRUE));
1262dce792eSToby Isaac PetscCall(PetscFEVectorTest(scalar, n_copies, PETSC_TRUE, PETSC_FALSE));
1272dce792eSToby Isaac PetscCall(PetscFEVectorTest(scalar, n_copies, PETSC_TRUE, PETSC_TRUE));
1282dce792eSToby Isaac }
1292dce792eSToby Isaac PetscCall(PetscViewerASCIIPopTab(viewer));
1302dce792eSToby Isaac PetscCall(PetscFEDestroy(&scalar));
1312dce792eSToby Isaac PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, 3, 3, PETSC_TRUE, 3, PETSC_DETERMINE, &vector));
1322dce792eSToby Isaac PetscCall(PetscObjectSetName((PetscObject)vector, "base FE (vector)"));
1332dce792eSToby Isaac PetscCall(PetscFEView(vector, viewer));
1342dce792eSToby Isaac PetscCall(PetscViewerASCIIPushTab(viewer));
1352dce792eSToby Isaac for (PetscInt n_copies = 1; n_copies <= 3; n_copies++) {
1362dce792eSToby Isaac PetscCall(PetscFEVectorTest(vector, n_copies, PETSC_FALSE, PETSC_FALSE));
1372dce792eSToby Isaac PetscCall(PetscFEVectorTest(vector, n_copies, PETSC_FALSE, PETSC_TRUE));
1382dce792eSToby Isaac PetscCall(PetscFEVectorTest(vector, n_copies, PETSC_TRUE, PETSC_FALSE));
1392dce792eSToby Isaac PetscCall(PetscFEVectorTest(vector, n_copies, PETSC_TRUE, PETSC_TRUE));
1402dce792eSToby Isaac }
1412dce792eSToby Isaac PetscCall(PetscViewerASCIIPopTab(viewer));
1422dce792eSToby Isaac PetscCall(PetscFEDestroy(&vector));
1432dce792eSToby Isaac PetscCall(PetscFinalize());
1442dce792eSToby Isaac return 0;
1452dce792eSToby Isaac }
1462dce792eSToby Isaac
1472dce792eSToby Isaac /*TEST
1482dce792eSToby Isaac
1492dce792eSToby Isaac test:
1502dce792eSToby Isaac suffix: 0
1512dce792eSToby Isaac
1522dce792eSToby Isaac TEST*/
153