xref: /petsc/src/dm/dt/tests/ex5.c (revision 98921bda46e76d7aaed9e0138c5ff9d0ce93f355)
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 
7c4762a1bSJed Brown int main(int argc, char **argv)
8c4762a1bSJed Brown {
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   PetscErrorCode ierr;
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
19c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
20c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm,"","Options for subspace test","none");CHKERRQ(ierr);
21c4762a1bSJed Brown   ierr = PetscOptionsRangeInt("-dim", "The spatial dimension","ex5.c",dim,&dim,NULL,1,3);CHKERRQ(ierr);
22c4762a1bSJed Brown   ierr = PetscOptionsBool("-simplex", "Test simplex element","ex5.c",simplex,&simplex,NULL);CHKERRQ(ierr);
23c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-num_comp", "Number of components in space","ex5.c",Nc,&Nc,NULL,1);CHKERRQ(ierr);
241e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
25c4762a1bSJed Brown   ierr = DMShellCreate(comm,&dm);CHKERRQ(ierr);
26c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm,dim,Nc,simplex,NULL,PETSC_DEFAULT,&fe);CHKERRQ(ierr);
27c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
28c4762a1bSJed Brown   ierr = PetscFESetName(fe, "solution");CHKERRQ(ierr);
29c4762a1bSJed Brown   ierr = PetscFEGetBasisSpace(fe,&space);CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = PetscSpaceGetNumComponents(space,&Nc);CHKERRQ(ierr);
31c4762a1bSJed Brown   ierr = PetscFEGetDualSpace(fe,&dualspace);CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr = PetscDualSpaceGetHeightSubspace(dualspace,1,&dualsubspace);CHKERRQ(ierr);
33c4762a1bSJed Brown   ierr = PetscDualSpaceGetDM(dualspace,&dm);CHKERRQ(ierr);
34c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
35c4762a1bSJed Brown   if (cEnd > cStart) {
36c4762a1bSJed Brown     PetscInt coneSize;
37c4762a1bSJed Brown 
38c4762a1bSJed Brown     ierr = DMPlexGetConeSize(dm,cStart,&coneSize);CHKERRQ(ierr);
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 
55c4762a1bSJed Brown       ierr = DMPlexGetCone(dm,cStart,&cone);CHKERRQ(ierr);
56c4762a1bSJed Brown       point = cone[0];
57c4762a1bSJed Brown       ierr = PetscFECreatePointTrace(fe,point,&traceFE);CHKERRQ(ierr);
58c4762a1bSJed Brown       ierr = PetscFESetUp(traceFE);CHKERRQ(ierr);
59c4762a1bSJed Brown       ierr = PetscFEViewFromOptions(traceFE,NULL,"-trace_fe_view");CHKERRQ(ierr);
60c4762a1bSJed Brown       ierr = PetscMalloc4(dim - 1,&testSub,dim,&testFull,Nc,&outSub,Nc,&outFull);CHKERRQ(ierr);
61c4762a1bSJed Brown       ierr = PetscRandomCreate(PETSC_COMM_SELF,&rand);CHKERRQ(ierr);
62c4762a1bSJed Brown       ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr);
63c4762a1bSJed Brown       ierr = PetscRandomSetInterval(rand,-1.,1.);CHKERRQ(ierr);
64c4762a1bSJed Brown       /* create a random point in the trace domain */
65c4762a1bSJed Brown       for (i = 0; i < dim - 1; i++) {
66c4762a1bSJed Brown         ierr = PetscRandomGetValueReal(rand,&testSub[i]);CHKERRQ(ierr);
67c4762a1bSJed Brown       }
68c4762a1bSJed Brown       ierr = DMPlexComputeCellGeometryFEM(dm,point,NULL,testFull,J,NULL,&detJ);CHKERRQ(ierr);
69c4762a1bSJed Brown       /* project it into the full domain */
70c4762a1bSJed Brown       for (i = 0; i < dim; i++) {
71c4762a1bSJed Brown         for (j = 0; j < dim - 1; j++) testFull[i] += J[i * dim + j] * (testSub[j] - xi0[j]);
72c4762a1bSJed Brown       }
73c4762a1bSJed Brown       /* create a random vector in the full domain */
74c4762a1bSJed Brown       ierr = PetscFEGetDimension(fe,&nFull);CHKERRQ(ierr);
75c4762a1bSJed Brown       ierr = VecCreateSeq(PETSC_COMM_SELF,nFull,&vecFull);CHKERRQ(ierr);
76c4762a1bSJed Brown       ierr = VecGetArray(vecFull,&arrayFull);CHKERRQ(ierr);
77c4762a1bSJed Brown       for (i = 0; i < nFull; i++) {
78c4762a1bSJed Brown         ierr = PetscRandomGetValue(rand,&arrayFull[i]);CHKERRQ(ierr);
79c4762a1bSJed Brown       }
80c4762a1bSJed Brown       ierr = VecRestoreArray(vecFull,&arrayFull);CHKERRQ(ierr);
81c4762a1bSJed Brown       /* create a vector on the trace domain */
82c4762a1bSJed Brown       ierr = PetscFEGetDimension(traceFE,&nSub);CHKERRQ(ierr);
83c4762a1bSJed Brown       /* get the subset of the original finite element space that is supported on the trace space */
84d4afb720SToby Isaac       ierr = PetscDualSpaceGetSection(dualspace,&sectionFull);CHKERRQ(ierr);
85c4762a1bSJed Brown       ierr = PetscSectionSetUp(sectionFull);CHKERRQ(ierr);
86c4762a1bSJed Brown       /* get the trace degrees of freedom */
87c4762a1bSJed Brown       ierr = PetscMalloc1(nSub,&arraySub);CHKERRQ(ierr);
88c4762a1bSJed Brown       ierr = DMPlexVecGetClosure(dm,sectionFull,vecFull,point,&nSub,&arraySub);CHKERRQ(ierr);
89c4762a1bSJed Brown       /* get the tabulations */
90c4762a1bSJed Brown       ierr = PetscFECreateTabulation(traceFE,1,1,testSub,0,&Tsub);CHKERRQ(ierr);
91c4762a1bSJed Brown       ierr = PetscFECreateTabulation(fe,1,1,testFull,0,&Tfull);CHKERRQ(ierr);
92c4762a1bSJed Brown       for (i = 0; i < Nc; i++) {
93c4762a1bSJed Brown         outSub[i] = 0.0;
94c4762a1bSJed Brown         for (j = 0; j < nSub; j++) {
95c4762a1bSJed Brown           outSub[i] += Tsub->T[0][j * Nc + i] * arraySub[j];
96c4762a1bSJed Brown         }
97c4762a1bSJed Brown       }
98c4762a1bSJed Brown       ierr = VecGetArray(vecFull,&arrayFull);CHKERRQ(ierr);
99c4762a1bSJed Brown       err = 0.0;
100c4762a1bSJed Brown       for (i = 0; i < Nc; i++) {
101c4762a1bSJed Brown         PetscScalar diff;
102c4762a1bSJed Brown 
103c4762a1bSJed Brown         outFull[i] = 0.0;
104c4762a1bSJed Brown         for (j = 0; j < nFull; j++) {
105c4762a1bSJed Brown           outFull[i] += Tfull->T[0][j * Nc + i] * arrayFull[j];
106c4762a1bSJed Brown         }
107c4762a1bSJed Brown         diff = outFull[i] - outSub[i];
108c4762a1bSJed Brown         err += PetscRealPart(PetscConj(diff) * diff);
109c4762a1bSJed Brown       }
110c4762a1bSJed Brown       err = PetscSqrtReal(err);
111c4762a1bSJed Brown       if (err > PETSC_SMALL) {
112*98921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Trace FE error %g",err);
113c4762a1bSJed Brown       }
114c4762a1bSJed Brown       ierr = VecRestoreArray(vecFull,&arrayFull);CHKERRQ(ierr);
115c4762a1bSJed Brown       ierr = PetscTabulationDestroy(&Tfull);CHKERRQ(ierr);
116c4762a1bSJed Brown       ierr = PetscTabulationDestroy(&Tsub);CHKERRQ(ierr);
117c4762a1bSJed Brown       /* clean up */
118c4762a1bSJed Brown       ierr = PetscFree(arraySub);CHKERRQ(ierr);
119c4762a1bSJed Brown       ierr = VecDestroy(&vecFull);CHKERRQ(ierr);
120c4762a1bSJed Brown       ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
121c4762a1bSJed Brown       ierr = PetscFree4(testSub,testFull,outSub,outFull);CHKERRQ(ierr);
122c4762a1bSJed Brown       ierr = PetscFEDestroy(&traceFE);CHKERRQ(ierr);
123c4762a1bSJed Brown     }
124c4762a1bSJed Brown   }
125c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
126c4762a1bSJed Brown   ierr = PetscFinalize();
127c4762a1bSJed Brown   return ierr;
128c4762a1bSJed Brown }
129c4762a1bSJed Brown 
130c4762a1bSJed Brown /*TEST
131c4762a1bSJed Brown   test:
132c4762a1bSJed Brown     suffix: 0
133c4762a1bSJed Brown     args: -petscspace_degree 1 -trace_fe_view
134c4762a1bSJed Brown TEST*/
135