1243b716aSBarry Smith! Contributed by Noem T 2243b716aSBarry Smith#include <petsc/finclude/petscdmplex.h> 3243b716aSBarry Smith#include <petsc/finclude/petscdm.h> 4*c5e229c2SMartin Diehlprogram test 5243b716aSBarry Smith use PETScDM 6243b716aSBarry Smith use PETScDMplex 7243b716aSBarry Smith implicit none 8243b716aSBarry Smith DM :: dm 9243b716aSBarry Smith PetscInt :: numCells = 1 10243b716aSBarry Smith PetscInt :: cStart 11243b716aSBarry Smith PetscInt :: numVertices = 3, numCorners = 3 12243b716aSBarry Smith PetscErrorCode :: ierr 13243b716aSBarry Smith PetscInt, parameter :: numDim = 2 14243b716aSBarry Smith PetscReal, parameter :: eps = 5.0*epsilon(1.0) 15243b716aSBarry Smith PetscCallA(PetscInitialize(ierr)) 16243b716aSBarry Smith triangle: block 17243b716aSBarry Smith 18243b716aSBarry Smith! 19243b716aSBarry Smith! Single triangle element 20243b716aSBarry Smith! Corner coordinates (1, 1), (5, 5), (3, 6) 21243b716aSBarry Smith! Using a 3-point quadrature rule 22243b716aSBarry Smith! 23243b716aSBarry Smith PetscInt :: i 24243b716aSBarry Smith PetscInt :: zero = 0 25243b716aSBarry Smith 26243b716aSBarry Smith! Number of quadrature points 27243b716aSBarry Smith PetscInt, parameter :: tria_qpts = 3 28243b716aSBarry Smith! Quadrature order 29243b716aSBarry Smith PetscInt, parameter :: tria_qorder = 2 30243b716aSBarry Smith! Mapped quadrature points 31243b716aSBarry Smith PetscReal, parameter :: tria_v(tria_qpts*numDim) = [2.0, 2.5, 3.0, 5.0, 4.0, 4.5] 32243b716aSBarry Smith! Jacobian (constant, repeated tria_qpts times) 33243b716aSBarry Smith PetscReal, parameter :: tria_J(tria_qpts*numDim**2) = [([2.0, 1.0, 2.0, 2.5], i=1, tria_qpts)] 34243b716aSBarry Smith 35243b716aSBarry Smith PetscQuadrature :: q 36243b716aSBarry Smith PetscReal :: J(tria_qpts*numDim**2), invJ(tria_qpts*numDim**2), v(tria_qpts*numDim), detJ(tria_qpts) 37243b716aSBarry Smith PetscReal :: vertexCoords(6) = 1.0*[1.0, 1.0, 5.0, 5.0, 3.0, 6.0] 38243b716aSBarry Smith PetscInt :: cells(3) = [0, 1, 2] 39243b716aSBarry Smith 40243b716aSBarry Smith PetscCallA(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, numDim, numCells, numVertices, numCorners, PETSC_FALSE, cells, numDim, vertexCoords, dm, ierr)) 41243b716aSBarry Smith PetscCallA(PetscDTSimplexQuadrature(numDim, tria_qorder, PETSCDTSIMPLEXQUAD_DEFAULT, q, ierr)) 42243b716aSBarry Smith PetscCallA(DMPlexGetHeightStratum(dm, zero, cStart, PETSC_NULL_INTEGER, ierr)) 43243b716aSBarry Smith PetscCallA(DMPlexComputeCellGeometryFEM(dm, cStart, q, v, J, invJ, detJ, ierr)) 44243b716aSBarry Smith PetscCheckA(all(abs(v - tria_v) < eps), PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Wrong mapped quadrature points (triangle)') 45243b716aSBarry Smith PetscCheckA(all(abs(J - tria_J) < eps), PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Wrong jacobian (triangle)') 46243b716aSBarry Smith PetscCallA(PetscQuadratureDestroy(q, ierr)) 47243b716aSBarry Smith PetscCallA(DMDestroy(dm, ierr)) 48243b716aSBarry Smith end block triangle 49243b716aSBarry Smith 50243b716aSBarry Smith quadrilateral: block 51243b716aSBarry Smith 52243b716aSBarry Smith! 53243b716aSBarry Smith! Single quadrilateral element 54243b716aSBarry Smith! Corner coordinates (-3, -2), (3, -1), (2, 4), (-2, 3) 55243b716aSBarry Smith! Using a 4-point (2x2) Gauss quadrature rule 56243b716aSBarry Smith! 57243b716aSBarry Smith 58243b716aSBarry Smith! Number of quadrature points 59243b716aSBarry Smith PetscInt, parameter :: quad_qpts = 4 60243b716aSBarry Smith 61243b716aSBarry Smith PetscQuadrature :: q 62243b716aSBarry Smith PetscReal :: vertexCoords(8) = [-3.0, -2.0, 3.0, -1.0, 2.0, 4.0, -2.0, 3.0] 63243b716aSBarry Smith PetscReal :: J(quad_qpts*numDim**2), invJ(quad_qpts*numDim**2), v(quad_qpts*numDim), detJ(quad_qpts) 64243b716aSBarry Smith PetscReal :: a = -1.0, b = 1.0 65243b716aSBarry Smith PetscInt :: cells(4) = [0, 1, 2, 3] 66243b716aSBarry Smith PetscInt :: nc = 1, npoints = 2 67243b716aSBarry Smith PetscInt :: k 68243b716aSBarry Smith PetscInt :: zero = 0 69243b716aSBarry Smith 70243b716aSBarry Smith numCells = 1 71243b716aSBarry Smith numCorners = 4 72243b716aSBarry Smith numVertices = 4 73243b716aSBarry Smith 74243b716aSBarry Smith PetscCallA(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, numDim, numCells, numVertices, numCorners, PETSC_FALSE, cells, numDim, vertexCoords, dm, ierr)) 75243b716aSBarry Smith PetscCallA(PetscDTGaussTensorQuadrature(numDim, nc, npoints, a, b, q, ierr)) 76243b716aSBarry Smith PetscCallA(DMPlexGetHeightStratum(dm, zero, cStart, PETSC_NULL_INTEGER, ierr)) 77243b716aSBarry Smith PetscCallA(DMPlexComputeCellGeometryFEM(dm, cStart, q, v, J, invJ, detJ, ierr)) 78243b716aSBarry Smith do k = 1, quad_qpts 79243b716aSBarry Smith PetscCheckA(all(abs(v((k - 1)*numDim + 1:k*numDim) - quad_v(Gauss_rs(k))) < eps), PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Wrong mapped quadrature points (quadrilateral)') 80243b716aSBarry Smith PetscCheckA(all(abs(J((k - 1)*numDim**2 + 1:k*numDim**2) - quad_J(Gauss_rs(k))) < eps), PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Wrong jacobian (quadrilateral)') 81243b716aSBarry Smith end do 82243b716aSBarry Smith PetscCallA(PetscQuadratureDestroy(q, ierr)) 83243b716aSBarry Smith PetscCallA(DMDestroy(dm, ierr)) 84243b716aSBarry Smith end block quadrilateral 85243b716aSBarry Smith PetscCallA(PetscFinalize(ierr)) 86243b716aSBarry Smith 87243b716aSBarry Smithcontains 88243b716aSBarry Smith! Quadrature points in a quadrilateral in [-1,+1] 89243b716aSBarry Smith function Gauss_rs(n) 90243b716aSBarry Smith PetscInt, intent(in) :: n 91243b716aSBarry Smith PetscReal :: Gauss_rs(2) 92243b716aSBarry Smith 93243b716aSBarry Smith PetscReal, parameter :: p = 1.0/sqrt(3.0) 94243b716aSBarry Smith 95243b716aSBarry Smith select case (n) 96243b716aSBarry Smith case (1) 97243b716aSBarry Smith Gauss_rs = [-p, -p] 98243b716aSBarry Smith case (2) 99243b716aSBarry Smith Gauss_rs = [-p, +p] 100243b716aSBarry Smith case (3) 101243b716aSBarry Smith Gauss_rs = [+p, -p] 102243b716aSBarry Smith case (4) 103243b716aSBarry Smith Gauss_rs = [+p, +p] 104243b716aSBarry Smith end select 105243b716aSBarry Smith end function Gauss_rs 106243b716aSBarry Smith! Mapped quadrature points 107243b716aSBarry Smith function quad_v(rs) 108243b716aSBarry Smith PetscReal, intent(in) :: rs(2) 109243b716aSBarry Smith PetscReal :: quad_v(2) 110243b716aSBarry Smith 111243b716aSBarry Smith PetscReal :: r, s 112243b716aSBarry Smith 113243b716aSBarry Smith r = rs(1) 114243b716aSBarry Smith s = rs(2) 115243b716aSBarry Smith quad_v(1) = -0.5*r*(s - 5) ! sum N_i * x_i 116243b716aSBarry Smith quad_v(2) = 0.5*(r + 5*s + 2) ! sum N_i * y_i 117243b716aSBarry Smith 118243b716aSBarry Smith end function quad_v 119243b716aSBarry Smith! Jacobian 120243b716aSBarry Smith function quad_J(rs) 121243b716aSBarry Smith PetscReal, intent(in) :: rs(2) 122243b716aSBarry Smith PetscReal :: quad_J(4) 123243b716aSBarry Smith 124243b716aSBarry Smith PetscReal :: r, s 125243b716aSBarry Smith PetscReal :: pfive = .5, twopfive = 2.5 126243b716aSBarry Smith 127243b716aSBarry Smith r = rs(1) 128243b716aSBarry Smith s = rs(2) 129243b716aSBarry Smith quad_J = [-0.5*(s - 5), -0.5*r, pfive, twopfive] 130243b716aSBarry Smith end function quad_J 131243b716aSBarry Smithend program test 132243b716aSBarry Smith 133243b716aSBarry Smith! /*TEST 134243b716aSBarry Smith! 135243b716aSBarry Smith! test: 136e0008caeSPierre Jolivet! output_file: output/empty.out 137243b716aSBarry Smith! 138243b716aSBarry Smith! TEST*/ 139