xref: /petsc/src/dm/dt/tests/ex5.c (revision 061e922f3926be00487707c73b78dd3d40309129)
1c4762a1bSJed Brown static char help[] = "Tests affine subspaces.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscfe.h>
4c4762a1bSJed Brown #include <petscdmplex.h>
5c4762a1bSJed Brown #include <petscdmshell.h>
6c4762a1bSJed Brown 
main(int argc,char ** argv)7*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
8*d71ae5a4SJacob Faibussowitsch {
9c4762a1bSJed Brown   DM             dm;
10c4762a1bSJed Brown   PetscFE        fe;
11c4762a1bSJed Brown   PetscSpace     space;
12c4762a1bSJed Brown   PetscDualSpace dualspace, dualsubspace;
13c4762a1bSJed Brown   PetscInt       dim = 2, Nc = 3, cStart, cEnd;
14c4762a1bSJed Brown   PetscBool      simplex = PETSC_TRUE;
15c4762a1bSJed Brown   MPI_Comm       comm;
16c4762a1bSJed Brown 
17327415f7SBarry Smith   PetscFunctionBeginUser;
189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
19c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
20d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Options for subspace test", "none");
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-dim", "The spatial dimension", "ex5.c", dim, &dim, NULL, 1, 3));
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-simplex", "Test simplex element", "ex5.c", simplex, &simplex, NULL));
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_comp", "Number of components in space", "ex5.c", Nc, &Nc, NULL, 1));
24d0609cedSBarry Smith   PetscOptionsEnd();
259566063dSJacob Faibussowitsch   PetscCall(DMShellCreate(comm, &dm));
269566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, Nc, simplex, NULL, PETSC_DEFAULT, &fe));
279566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
289566063dSJacob Faibussowitsch   PetscCall(PetscFESetName(fe, "solution"));
299566063dSJacob Faibussowitsch   PetscCall(PetscFEGetBasisSpace(fe, &space));
309566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(space, &Nc));
319566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(fe, &dualspace));
329566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetHeightSubspace(dualspace, 1, &dualsubspace));
339566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(dualspace, &dm));
349566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
35c4762a1bSJed Brown   if (cEnd > cStart) {
36c4762a1bSJed Brown     PetscInt coneSize;
37c4762a1bSJed Brown 
389566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(dm, cStart, &coneSize));
39c4762a1bSJed Brown     if (coneSize) {
40c4762a1bSJed Brown       PetscFE         traceFE;
41c4762a1bSJed Brown       const PetscInt *cone;
42c4762a1bSJed Brown       PetscInt        point, nSub, nFull;
43c4762a1bSJed Brown       PetscReal       xi0[3] = {-1., -1., -1.};
44c4762a1bSJed Brown       PetscScalar    *outSub, *outFull;
45c4762a1bSJed Brown       PetscReal      *testSub, *testFull;
46c4762a1bSJed Brown       PetscTabulation Tsub, Tfull;
47c4762a1bSJed Brown       PetscReal       J[9], detJ;
48c4762a1bSJed Brown       PetscInt        i, j;
49c4762a1bSJed Brown       PetscSection    sectionFull;
50c4762a1bSJed Brown       Vec             vecFull;
51c4762a1bSJed Brown       PetscScalar    *arrayFull, *arraySub;
52c4762a1bSJed Brown       PetscReal       err;
53c4762a1bSJed Brown       PetscRandom     rand;
54c4762a1bSJed Brown 
559566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, cStart, &cone));
56c4762a1bSJed Brown       point = cone[0];
579566063dSJacob Faibussowitsch       PetscCall(PetscFECreatePointTrace(fe, point, &traceFE));
589566063dSJacob Faibussowitsch       PetscCall(PetscFESetUp(traceFE));
599566063dSJacob Faibussowitsch       PetscCall(PetscFEViewFromOptions(traceFE, NULL, "-trace_fe_view"));
609566063dSJacob Faibussowitsch       PetscCall(PetscMalloc4(dim - 1, &testSub, dim, &testFull, Nc, &outSub, Nc, &outFull));
619566063dSJacob Faibussowitsch       PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
629566063dSJacob Faibussowitsch       PetscCall(PetscRandomSetFromOptions(rand));
639566063dSJacob Faibussowitsch       PetscCall(PetscRandomSetInterval(rand, -1., 1.));
64c4762a1bSJed Brown       /* create a random point in the trace domain */
6548a46eb9SPierre Jolivet       for (i = 0; i < dim - 1; i++) PetscCall(PetscRandomGetValueReal(rand, &testSub[i]));
669566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dm, point, NULL, testFull, J, NULL, &detJ));
67c4762a1bSJed Brown       /* project it into the full domain */
68c4762a1bSJed Brown       for (i = 0; i < dim; i++) {
69c4762a1bSJed Brown         for (j = 0; j < dim - 1; j++) testFull[i] += J[i * dim + j] * (testSub[j] - xi0[j]);
70c4762a1bSJed Brown       }
71c4762a1bSJed Brown       /* create a random vector in the full domain */
729566063dSJacob Faibussowitsch       PetscCall(PetscFEGetDimension(fe, &nFull));
739566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF, nFull, &vecFull));
749566063dSJacob Faibussowitsch       PetscCall(VecGetArray(vecFull, &arrayFull));
7548a46eb9SPierre Jolivet       for (i = 0; i < nFull; i++) PetscCall(PetscRandomGetValue(rand, &arrayFull[i]));
769566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(vecFull, &arrayFull));
77c4762a1bSJed Brown       /* create a vector on the trace domain */
789566063dSJacob Faibussowitsch       PetscCall(PetscFEGetDimension(traceFE, &nSub));
79c4762a1bSJed Brown       /* get the subset of the original finite element space that is supported on the trace space */
809566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetSection(dualspace, &sectionFull));
819566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetUp(sectionFull));
82c4762a1bSJed Brown       /* get the trace degrees of freedom */
839566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nSub, &arraySub));
849566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(dm, sectionFull, vecFull, point, &nSub, &arraySub));
85c4762a1bSJed Brown       /* get the tabulations */
869566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation(traceFE, 1, 1, testSub, 0, &Tsub));
879566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation(fe, 1, 1, testFull, 0, &Tfull));
88c4762a1bSJed Brown       for (i = 0; i < Nc; i++) {
89c4762a1bSJed Brown         outSub[i] = 0.0;
90ad540459SPierre Jolivet         for (j = 0; j < nSub; j++) outSub[i] += Tsub->T[0][j * Nc + i] * arraySub[j];
91c4762a1bSJed Brown       }
929566063dSJacob Faibussowitsch       PetscCall(VecGetArray(vecFull, &arrayFull));
93c4762a1bSJed Brown       err = 0.0;
94c4762a1bSJed Brown       for (i = 0; i < Nc; i++) {
95c4762a1bSJed Brown         PetscScalar diff;
96c4762a1bSJed Brown 
97c4762a1bSJed Brown         outFull[i] = 0.0;
98ad540459SPierre Jolivet         for (j = 0; j < nFull; j++) outFull[i] += Tfull->T[0][j * Nc + i] * arrayFull[j];
99c4762a1bSJed Brown         diff = outFull[i] - outSub[i];
100c4762a1bSJed Brown         err += PetscRealPart(PetscConj(diff) * diff);
101c4762a1bSJed Brown       }
102c4762a1bSJed Brown       err = PetscSqrtReal(err);
1037a46b595SBarry Smith       PetscCheck(err <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Trace FE error %g", (double)err);
1049566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(vecFull, &arrayFull));
1059566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&Tfull));
1069566063dSJacob Faibussowitsch       PetscCall(PetscTabulationDestroy(&Tsub));
107c4762a1bSJed Brown       /* clean up */
1089566063dSJacob Faibussowitsch       PetscCall(PetscFree(arraySub));
1099566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&vecFull));
1109566063dSJacob Faibussowitsch       PetscCall(PetscRandomDestroy(&rand));
1119566063dSJacob Faibussowitsch       PetscCall(PetscFree4(testSub, testFull, outSub, outFull));
1129566063dSJacob Faibussowitsch       PetscCall(PetscFEDestroy(&traceFE));
113c4762a1bSJed Brown     }
114c4762a1bSJed Brown   }
1159566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
1169566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
117b122ec5aSJacob Faibussowitsch   return 0;
118c4762a1bSJed Brown }
119c4762a1bSJed Brown 
120c4762a1bSJed Brown /*TEST
121c4762a1bSJed Brown   test:
122c4762a1bSJed Brown     suffix: 0
123c4762a1bSJed Brown     args: -petscspace_degree 1 -trace_fe_view
124c4762a1bSJed Brown TEST*/
125