xref: /petsc/src/dm/impls/plex/tests/dmplexcomputecellgeometryfem.F90 (revision e0008cae9d8b220d3a191ff57f84cb339f1acef2)
1243b716aSBarry Smith!     Contributed by Noem T
2243b716aSBarry Smithprogram test
3243b716aSBarry Smith#include <petsc/finclude/petscdmplex.h>
4243b716aSBarry Smith#include <petsc/finclude/petscdm.h>
5243b716aSBarry Smithuse PETScDM
6243b716aSBarry Smithuse PETScDMplex
7243b716aSBarry Smithimplicit none
8243b716aSBarry SmithDM       :: dm
9243b716aSBarry SmithPetscInt :: numCells = 1
10243b716aSBarry SmithPetscInt :: cStart
11243b716aSBarry SmithPetscInt :: numVertices = 3, numCorners = 3
12243b716aSBarry SmithPetscErrorCode :: ierr
13243b716aSBarry SmithPetscInt, parameter :: numDim = 2
14243b716aSBarry SmithPetscReal,    parameter :: eps = 5.0 * epsilon(1.0)
15243b716aSBarry SmithPetscCallA(PetscInitialize(ierr))
16243b716aSBarry Smithtriangle: 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 SmithPetscInt :: i
24243b716aSBarry SmithPetscInt :: zero = 0
25243b716aSBarry Smith
26243b716aSBarry Smith! Number of quadrature points
27243b716aSBarry SmithPetscInt, parameter :: tria_qpts = 3
28243b716aSBarry Smith! Quadrature order
29243b716aSBarry SmithPetscInt, parameter :: tria_qorder = 2
30243b716aSBarry Smith! Mapped quadrature points
31243b716aSBarry SmithPetscReal, 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 SmithPetscReal, parameter :: tria_J(tria_qpts*numDim**2) = [([2.0, 1.0, 2.0, 2.5], i = 1, tria_qpts)]
34243b716aSBarry Smith
35243b716aSBarry SmithPetscQuadrature :: q
36243b716aSBarry SmithPetscReal :: J(tria_qpts*numDim**2), invJ(tria_qpts*numDim**2), v(tria_qpts*numDim), detJ(tria_qpts)
37243b716aSBarry SmithPetscReal :: vertexCoords(6) = 1.0 * [1.0, 1.0, 5.0, 5.0, 3.0, 6.0]
38243b716aSBarry SmithPetscInt  :: cells(3) = [0, 1, 2]
39243b716aSBarry Smith
40243b716aSBarry SmithPetscCallA(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD,numDim,numCells,numVertices,numCorners,PETSC_FALSE,cells,numDim,vertexCoords,dm,ierr))
41243b716aSBarry SmithPetscCallA(PetscDTSimplexQuadrature(numDim, tria_qorder, PETSCDTSIMPLEXQUAD_DEFAULT, q, ierr))
42243b716aSBarry SmithPetscCallA(DMPlexGetHeightStratum(dm, zero, cStart, PETSC_NULL_INTEGER, ierr))
43243b716aSBarry SmithPetscCallA(DMPlexComputeCellGeometryFEM(dm, cStart, q, v, J, invJ, detJ, ierr))
44243b716aSBarry SmithPetscCheckA(all(abs(v - tria_v) < eps), PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Wrong mapped quadrature points (triangle)')
45243b716aSBarry SmithPetscCheckA(all(abs(J - tria_J) < eps), PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Wrong jacobian (triangle)')
46243b716aSBarry SmithPetscCallA(PetscQuadratureDestroy(q, ierr))
47243b716aSBarry SmithPetscCallA(DMDestroy(dm,ierr))
48243b716aSBarry Smithend block triangle
49243b716aSBarry Smith
50243b716aSBarry Smithquadrilateral: 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 SmithPetscInt, parameter :: quad_qpts = 4
60243b716aSBarry Smith
61243b716aSBarry SmithPetscQuadrature :: q
62243b716aSBarry SmithPetscReal :: vertexCoords(8) = [-3.0, -2.0, 3.0, -1.0, 2.0, 4.0, -2.0, 3.0]
63243b716aSBarry SmithPetscReal :: J(quad_qpts*numDim**2), invJ(quad_qpts*numDim**2), v(quad_qpts*numDim), detJ(quad_qpts)
64243b716aSBarry SmithPetscReal :: a = -1.0, b = 1.0
65243b716aSBarry SmithPetscInt  :: cells(4) = [0, 1, 2, 3]
66243b716aSBarry SmithPetscInt  :: nc = 1, npoints = 2
67243b716aSBarry SmithPetscInt  :: k
68243b716aSBarry SmithPetscInt :: zero = 0
69243b716aSBarry Smith
70243b716aSBarry SmithnumCells = 1
71243b716aSBarry SmithnumCorners = 4
72243b716aSBarry SmithnumVertices = 4
73243b716aSBarry Smith
74243b716aSBarry SmithPetscCallA(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD,numDim,numCells,numVertices, numCorners,PETSC_FALSE,cells,numDim,vertexCoords,dm,ierr))
75243b716aSBarry SmithPetscCallA(PetscDTGaussTensorQuadrature(numDim, nc, npoints, a, b, q, ierr))
76243b716aSBarry SmithPetscCallA(DMPlexGetHeightStratum(dm, zero, cStart, PETSC_NULL_INTEGER, ierr))
77243b716aSBarry SmithPetscCallA(DMPlexComputeCellGeometryFEM(dm, cStart, q, v, J, invJ, detJ, ierr))
78243b716aSBarry Smithdo 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 Smithend do
82243b716aSBarry SmithPetscCallA(PetscQuadratureDestroy(q, ierr))
83243b716aSBarry SmithPetscCallA(DMDestroy(dm,ierr))
84243b716aSBarry Smithend block quadrilateral
85243b716aSBarry SmithPetscCallA(PetscFinalize(ierr))
86243b716aSBarry Smith
87243b716aSBarry Smithcontains
88243b716aSBarry Smith! Quadrature points in a quadrilateral in [-1,+1]
89243b716aSBarry Smithfunction Gauss_rs(n)
90243b716aSBarry SmithPetscInt, intent(in) :: n
91243b716aSBarry SmithPetscReal :: 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 Smithend select
105243b716aSBarry Smithend function Gauss_rs
106243b716aSBarry Smith! Mapped quadrature points
107243b716aSBarry Smithfunction quad_v(rs)
108243b716aSBarry SmithPetscReal, intent(in) :: rs(2)
109243b716aSBarry SmithPetscReal :: 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 Smithend function quad_v
119243b716aSBarry Smith! Jacobian
120243b716aSBarry Smithfunction quad_J(rs)
121243b716aSBarry SmithPetscReal, intent(in) :: rs(2)
122243b716aSBarry SmithPetscReal :: 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 Smithend function quad_J
131243b716aSBarry Smithend program test
132243b716aSBarry Smith
133243b716aSBarry Smith! /*TEST
134243b716aSBarry Smith!
135243b716aSBarry Smith! test:
136*e0008caeSPierre Jolivet!   output_file: output/empty.out
137243b716aSBarry Smith!
138243b716aSBarry Smith! TEST*/
139