1*c4762a1bSJed Brown static char help[] = "Tests affine subspaces.\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscfe.h> 4*c4762a1bSJed Brown #include <petscdmplex.h> 5*c4762a1bSJed Brown #include <petscdmshell.h> 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown int main(int argc, char **argv) 8*c4762a1bSJed Brown { 9*c4762a1bSJed Brown DM dm; 10*c4762a1bSJed Brown PetscFE fe; 11*c4762a1bSJed Brown PetscSpace space; 12*c4762a1bSJed Brown PetscDualSpace dualspace, dualsubspace; 13*c4762a1bSJed Brown PetscInt dim = 2, Nc = 3, cStart, cEnd; 14*c4762a1bSJed Brown PetscBool simplex = PETSC_TRUE; 15*c4762a1bSJed Brown MPI_Comm comm; 16*c4762a1bSJed Brown PetscErrorCode ierr; 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 19*c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 20*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm,"","Options for subspace test","none");CHKERRQ(ierr); 21*c4762a1bSJed Brown ierr = PetscOptionsRangeInt("-dim", "The spatial dimension","ex5.c",dim,&dim,NULL,1,3);CHKERRQ(ierr); 22*c4762a1bSJed Brown ierr = PetscOptionsBool("-simplex", "Test simplex element","ex5.c",simplex,&simplex,NULL);CHKERRQ(ierr); 23*c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-num_comp", "Number of components in space","ex5.c",Nc,&Nc,NULL,1);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = PetscOptionsEnd(); 25*c4762a1bSJed Brown ierr = DMShellCreate(comm,&dm);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = PetscFECreateDefault(comm,dim,Nc,simplex,NULL,PETSC_DEFAULT,&fe);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = PetscFESetName(fe, "solution");CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = PetscFEGetBasisSpace(fe,&space);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = PetscSpaceGetNumComponents(space,&Nc);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = PetscFEGetDualSpace(fe,&dualspace);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = PetscDualSpaceGetHeightSubspace(dualspace,1,&dualsubspace);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = PetscDualSpaceGetDM(dualspace,&dm);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 35*c4762a1bSJed Brown if (cEnd > cStart) { 36*c4762a1bSJed Brown PetscInt coneSize; 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown ierr = DMPlexGetConeSize(dm,cStart,&coneSize);CHKERRQ(ierr); 39*c4762a1bSJed Brown if (coneSize) { 40*c4762a1bSJed Brown PetscFE traceFE; 41*c4762a1bSJed Brown const PetscInt *cone; 42*c4762a1bSJed Brown PetscInt point, nSub, nFull; 43*c4762a1bSJed Brown PetscReal xi0[3] = {-1., -1., -1.}; 44*c4762a1bSJed Brown PetscScalar *outSub, *outFull; 45*c4762a1bSJed Brown PetscReal *testSub, *testFull; 46*c4762a1bSJed Brown PetscTabulation Tsub, Tfull; 47*c4762a1bSJed Brown PetscReal J[9], detJ; 48*c4762a1bSJed Brown PetscInt i, j; 49*c4762a1bSJed Brown PetscSection sectionFull; 50*c4762a1bSJed Brown Vec vecFull; 51*c4762a1bSJed Brown PetscScalar *arrayFull, *arraySub; 52*c4762a1bSJed Brown PetscReal err; 53*c4762a1bSJed Brown PetscRandom rand; 54*c4762a1bSJed Brown 55*c4762a1bSJed Brown ierr = DMPlexGetCone(dm,cStart,&cone);CHKERRQ(ierr); 56*c4762a1bSJed Brown point = cone[0]; 57*c4762a1bSJed Brown ierr = PetscFECreatePointTrace(fe,point,&traceFE);CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = PetscFESetUp(traceFE);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = PetscFEViewFromOptions(traceFE,NULL,"-trace_fe_view");CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = PetscMalloc4(dim - 1,&testSub,dim,&testFull,Nc,&outSub,Nc,&outFull);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_SELF,&rand);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = PetscRandomSetInterval(rand,-1.,1.);CHKERRQ(ierr); 64*c4762a1bSJed Brown /* create a random point in the trace domain */ 65*c4762a1bSJed Brown for (i = 0; i < dim - 1; i++) { 66*c4762a1bSJed Brown ierr = PetscRandomGetValueReal(rand,&testSub[i]);CHKERRQ(ierr); 67*c4762a1bSJed Brown } 68*c4762a1bSJed Brown ierr = DMPlexComputeCellGeometryFEM(dm,point,NULL,testFull,J,NULL,&detJ);CHKERRQ(ierr); 69*c4762a1bSJed Brown /* project it into the full domain */ 70*c4762a1bSJed Brown for (i = 0; i < dim; i++) { 71*c4762a1bSJed Brown for (j = 0; j < dim - 1; j++) testFull[i] += J[i * dim + j] * (testSub[j] - xi0[j]); 72*c4762a1bSJed Brown } 73*c4762a1bSJed Brown /* create a random vector in the full domain */ 74*c4762a1bSJed Brown ierr = PetscFEGetDimension(fe,&nFull);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,nFull,&vecFull);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = VecGetArray(vecFull,&arrayFull);CHKERRQ(ierr); 77*c4762a1bSJed Brown for (i = 0; i < nFull; i++) { 78*c4762a1bSJed Brown ierr = PetscRandomGetValue(rand,&arrayFull[i]);CHKERRQ(ierr); 79*c4762a1bSJed Brown } 80*c4762a1bSJed Brown ierr = VecRestoreArray(vecFull,&arrayFull);CHKERRQ(ierr); 81*c4762a1bSJed Brown /* create a vector on the trace domain */ 82*c4762a1bSJed Brown ierr = PetscFEGetDimension(traceFE,&nSub);CHKERRQ(ierr); 83*c4762a1bSJed Brown /* get the subset of the original finite element space that is supported on the trace space */ 84*c4762a1bSJed Brown ierr = PetscDualSpaceCreateSection(dualspace,§ionFull);CHKERRQ(ierr); 85*c4762a1bSJed Brown ierr = PetscSectionSetUp(sectionFull);CHKERRQ(ierr); 86*c4762a1bSJed Brown /* get the trace degrees of freedom */ 87*c4762a1bSJed Brown ierr = PetscMalloc1(nSub,&arraySub);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = DMPlexVecGetClosure(dm,sectionFull,vecFull,point,&nSub,&arraySub);CHKERRQ(ierr); 89*c4762a1bSJed Brown /* get the tabulations */ 90*c4762a1bSJed Brown ierr = PetscFECreateTabulation(traceFE,1,1,testSub,0,&Tsub);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = PetscFECreateTabulation(fe,1,1,testFull,0,&Tfull);CHKERRQ(ierr); 92*c4762a1bSJed Brown for (i = 0; i < Nc; i++) { 93*c4762a1bSJed Brown outSub[i] = 0.0; 94*c4762a1bSJed Brown for (j = 0; j < nSub; j++) { 95*c4762a1bSJed Brown outSub[i] += Tsub->T[0][j * Nc + i] * arraySub[j]; 96*c4762a1bSJed Brown } 97*c4762a1bSJed Brown } 98*c4762a1bSJed Brown ierr = VecGetArray(vecFull,&arrayFull);CHKERRQ(ierr); 99*c4762a1bSJed Brown err = 0.0; 100*c4762a1bSJed Brown for (i = 0; i < Nc; i++) { 101*c4762a1bSJed Brown PetscScalar diff; 102*c4762a1bSJed Brown 103*c4762a1bSJed Brown outFull[i] = 0.0; 104*c4762a1bSJed Brown for (j = 0; j < nFull; j++) { 105*c4762a1bSJed Brown outFull[i] += Tfull->T[0][j * Nc + i] * arrayFull[j]; 106*c4762a1bSJed Brown } 107*c4762a1bSJed Brown diff = outFull[i] - outSub[i]; 108*c4762a1bSJed Brown err += PetscRealPart(PetscConj(diff) * diff); 109*c4762a1bSJed Brown } 110*c4762a1bSJed Brown err = PetscSqrtReal(err); 111*c4762a1bSJed Brown if (err > PETSC_SMALL) { 112*c4762a1bSJed Brown SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Trace FE error %g\n",err); 113*c4762a1bSJed Brown } 114*c4762a1bSJed Brown ierr = VecRestoreArray(vecFull,&arrayFull);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = PetscTabulationDestroy(&Tfull);CHKERRQ(ierr); 116*c4762a1bSJed Brown ierr = PetscTabulationDestroy(&Tsub);CHKERRQ(ierr); 117*c4762a1bSJed Brown /* clean up */ 118*c4762a1bSJed Brown ierr = PetscFree(arraySub);CHKERRQ(ierr); 119*c4762a1bSJed Brown ierr = PetscSectionDestroy(§ionFull);CHKERRQ(ierr); 120*c4762a1bSJed Brown ierr = VecDestroy(&vecFull);CHKERRQ(ierr); 121*c4762a1bSJed Brown ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = PetscFree4(testSub,testFull,outSub,outFull);CHKERRQ(ierr); 123*c4762a1bSJed Brown ierr = PetscFEDestroy(&traceFE);CHKERRQ(ierr); 124*c4762a1bSJed Brown } 125*c4762a1bSJed Brown } 126*c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 127*c4762a1bSJed Brown ierr = PetscFinalize(); 128*c4762a1bSJed Brown return ierr; 129*c4762a1bSJed Brown } 130*c4762a1bSJed Brown 131*c4762a1bSJed Brown /*TEST 132*c4762a1bSJed Brown test: 133*c4762a1bSJed Brown suffix: 0 134*c4762a1bSJed Brown args: -petscspace_degree 1 -trace_fe_view 135*c4762a1bSJed Brown TEST*/ 136