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 79371c9d4SSatish Balay int main(int argc, char **argv) { 8c4762a1bSJed Brown DM dm; 9c4762a1bSJed Brown PetscFE fe; 10c4762a1bSJed Brown PetscSpace space; 11c4762a1bSJed Brown PetscDualSpace dualspace, dualsubspace; 12c4762a1bSJed Brown PetscInt dim = 2, Nc = 3, cStart, cEnd; 13c4762a1bSJed Brown PetscBool simplex = PETSC_TRUE; 14c4762a1bSJed Brown MPI_Comm comm; 15c4762a1bSJed Brown 16327415f7SBarry Smith PetscFunctionBeginUser; 179566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 18c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 19d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Options for subspace test", "none"); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsRangeInt("-dim", "The spatial dimension", "ex5.c", dim, &dim, NULL, 1, 3)); 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-simplex", "Test simplex element", "ex5.c", simplex, &simplex, NULL)); 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-num_comp", "Number of components in space", "ex5.c", Nc, &Nc, NULL, 1)); 23d0609cedSBarry Smith PetscOptionsEnd(); 249566063dSJacob Faibussowitsch PetscCall(DMShellCreate(comm, &dm)); 259566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, Nc, simplex, NULL, PETSC_DEFAULT, &fe)); 269566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 279566063dSJacob Faibussowitsch PetscCall(PetscFESetName(fe, "solution")); 289566063dSJacob Faibussowitsch PetscCall(PetscFEGetBasisSpace(fe, &space)); 299566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumComponents(space, &Nc)); 309566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(fe, &dualspace)); 319566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetHeightSubspace(dualspace, 1, &dualsubspace)); 329566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(dualspace, &dm)); 339566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 34c4762a1bSJed Brown if (cEnd > cStart) { 35c4762a1bSJed Brown PetscInt coneSize; 36c4762a1bSJed Brown 379566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cStart, &coneSize)); 38c4762a1bSJed Brown if (coneSize) { 39c4762a1bSJed Brown PetscFE traceFE; 40c4762a1bSJed Brown const PetscInt *cone; 41c4762a1bSJed Brown PetscInt point, nSub, nFull; 42c4762a1bSJed Brown PetscReal xi0[3] = {-1., -1., -1.}; 43c4762a1bSJed Brown PetscScalar *outSub, *outFull; 44c4762a1bSJed Brown PetscReal *testSub, *testFull; 45c4762a1bSJed Brown PetscTabulation Tsub, Tfull; 46c4762a1bSJed Brown PetscReal J[9], detJ; 47c4762a1bSJed Brown PetscInt i, j; 48c4762a1bSJed Brown PetscSection sectionFull; 49c4762a1bSJed Brown Vec vecFull; 50c4762a1bSJed Brown PetscScalar *arrayFull, *arraySub; 51c4762a1bSJed Brown PetscReal err; 52c4762a1bSJed Brown PetscRandom rand; 53c4762a1bSJed Brown 549566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cStart, &cone)); 55c4762a1bSJed Brown point = cone[0]; 569566063dSJacob Faibussowitsch PetscCall(PetscFECreatePointTrace(fe, point, &traceFE)); 579566063dSJacob Faibussowitsch PetscCall(PetscFESetUp(traceFE)); 589566063dSJacob Faibussowitsch PetscCall(PetscFEViewFromOptions(traceFE, NULL, "-trace_fe_view")); 599566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dim - 1, &testSub, dim, &testFull, Nc, &outSub, Nc, &outFull)); 609566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand)); 619566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 629566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand, -1., 1.)); 63c4762a1bSJed Brown /* create a random point in the trace domain */ 64*48a46eb9SPierre Jolivet for (i = 0; i < dim - 1; i++) PetscCall(PetscRandomGetValueReal(rand, &testSub[i])); 659566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, point, NULL, testFull, J, NULL, &detJ)); 66c4762a1bSJed Brown /* project it into the full domain */ 67c4762a1bSJed Brown for (i = 0; i < dim; i++) { 68c4762a1bSJed Brown for (j = 0; j < dim - 1; j++) testFull[i] += J[i * dim + j] * (testSub[j] - xi0[j]); 69c4762a1bSJed Brown } 70c4762a1bSJed Brown /* create a random vector in the full domain */ 719566063dSJacob Faibussowitsch PetscCall(PetscFEGetDimension(fe, &nFull)); 729566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, nFull, &vecFull)); 739566063dSJacob Faibussowitsch PetscCall(VecGetArray(vecFull, &arrayFull)); 74*48a46eb9SPierre Jolivet for (i = 0; i < nFull; i++) PetscCall(PetscRandomGetValue(rand, &arrayFull[i])); 759566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vecFull, &arrayFull)); 76c4762a1bSJed Brown /* create a vector on the trace domain */ 779566063dSJacob Faibussowitsch PetscCall(PetscFEGetDimension(traceFE, &nSub)); 78c4762a1bSJed Brown /* get the subset of the original finite element space that is supported on the trace space */ 799566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetSection(dualspace, §ionFull)); 809566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(sectionFull)); 81c4762a1bSJed Brown /* get the trace degrees of freedom */ 829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nSub, &arraySub)); 839566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, sectionFull, vecFull, point, &nSub, &arraySub)); 84c4762a1bSJed Brown /* get the tabulations */ 859566063dSJacob Faibussowitsch PetscCall(PetscFECreateTabulation(traceFE, 1, 1, testSub, 0, &Tsub)); 869566063dSJacob Faibussowitsch PetscCall(PetscFECreateTabulation(fe, 1, 1, testFull, 0, &Tfull)); 87c4762a1bSJed Brown for (i = 0; i < Nc; i++) { 88c4762a1bSJed Brown outSub[i] = 0.0; 899371c9d4SSatish Balay for (j = 0; j < nSub; j++) { outSub[i] += Tsub->T[0][j * Nc + i] * arraySub[j]; } 90c4762a1bSJed Brown } 919566063dSJacob Faibussowitsch PetscCall(VecGetArray(vecFull, &arrayFull)); 92c4762a1bSJed Brown err = 0.0; 93c4762a1bSJed Brown for (i = 0; i < Nc; i++) { 94c4762a1bSJed Brown PetscScalar diff; 95c4762a1bSJed Brown 96c4762a1bSJed Brown outFull[i] = 0.0; 979371c9d4SSatish Balay for (j = 0; j < nFull; j++) { outFull[i] += Tfull->T[0][j * Nc + i] * arrayFull[j]; } 98c4762a1bSJed Brown diff = outFull[i] - outSub[i]; 99c4762a1bSJed Brown err += PetscRealPart(PetscConj(diff) * diff); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown err = PetscSqrtReal(err); 1027a46b595SBarry Smith PetscCheck(err <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Trace FE error %g", (double)err); 1039566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vecFull, &arrayFull)); 1049566063dSJacob Faibussowitsch PetscCall(PetscTabulationDestroy(&Tfull)); 1059566063dSJacob Faibussowitsch PetscCall(PetscTabulationDestroy(&Tsub)); 106c4762a1bSJed Brown /* clean up */ 1079566063dSJacob Faibussowitsch PetscCall(PetscFree(arraySub)); 1089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecFull)); 1099566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 1109566063dSJacob Faibussowitsch PetscCall(PetscFree4(testSub, testFull, outSub, outFull)); 1119566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&traceFE)); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown } 1149566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 1159566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 116b122ec5aSJacob Faibussowitsch return 0; 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown /*TEST 120c4762a1bSJed Brown test: 121c4762a1bSJed Brown suffix: 0 122c4762a1bSJed Brown args: -petscspace_degree 1 -trace_fe_view 123c4762a1bSJed Brown TEST*/ 124