The Open FUSION Toolkit 1.0.0-8905cc5
Modeling tools for plasma and fusion research and engineering
Loading...
Searching...
No Matches
Data Types | Functions/Subroutines | Variables
oft_tetmesh_type Module Reference

Detailed Description

Tetrahedral mesh structure definitions.

Global Tet variables

Author
Chris Hansen
Date
June 2010

Data Types

type  oft_tetmesh
 Tetrahedral volume mesh type. More...
 

Functions/Subroutines

subroutine tet_3d_grid (order, xnodes, inodesf, inodesc)
 Needs docs.
 
subroutine tetmesh_ctang (self, cell, ind, f, tang)
 Compute the curve tangent vector for a given edge on a cell.
 
subroutine tetmesh_g2inv (jfull, g2op)
 Needs docs.
 
subroutine tetmesh_get_surf_map (self, face, cell, lmap)
 Get mapping between boundary and volume logical coordinates.
 
subroutine tetmesh_hessian (self, cell, f, g2op, k)
 Compute the spatial hessian matrices for a given cell at a given logical position.
 
integer(i4) function tetmesh_in_cell (self, f, tol)
 Test if logical position lies within the base cell.
 
subroutine tetmesh_invert_cell (self, cell)
 Turn cell "inside out", used to ensure consistent orientations.
 
subroutine tetmesh_jacinv (jfull, gop, jac)
 Invert a 3x3 matrix.
 
subroutine tetmesh_jacl (self, cell, gop, j)
 Linear implementation of @tetmesh_jacobian.
 
subroutine tetmesh_jacobian (self, cell, f, gop, j)
 Compute the spatial jacobian matrix and its determinant for a given cell at a given logical position.
 
real(r8) function, dimension(3) tetmesh_log2phys (self, cell, f)
 Map from logical to physical coordinates in a given cell.
 
subroutine tetmesh_phys2log (self, cell, pt, f)
 Map from physical to logical coordinates in a given cell.
 
subroutine tetmesh_quad_rule (self, order, quad_rule)
 Retrieve suitable quadrature rule for mesh with given order.
 
subroutine tetmesh_set_order (self, order)
 Set maximum order of spatial mapping.
 
subroutine tetmesh_setup (self, cad_type)
 Setup mesh with implementation specifics (cell_np, cell_ne, etc.)
 
subroutine tetmesh_snormal (self, cell, ind, f, norm)
 Compute the surface normal vector for a given face on a cell.
 
subroutine tetmesh_surf_to_vol (self, fsurf, lmap, fvol)
 Map between surface and volume logical coordinates.
 
subroutine tetmesh_tessellate (self, rtmp, lctmp, order)
 Tessellate mesh onto lagrange FE nodes of specified order (usually for plotting)
 
integer(i4) function, dimension(2) tetmesh_tessellated_sizes (self)
 Get sizes of arrays returned by tetmesh_tessellate.
 
subroutine tetmesh_vlog (self, i, f)
 Get position in logical space of vertex i
 

Variables

integer(i4), private active_cell = 0
 Active cell for high order find_cell.
 
class(oft_tetmesh), pointer, private active_mesh => NULL()
 Active mesh for high order find_cell.
 
real(r8), dimension(3), private active_pt = 0.d0
 Active point for high order find_cell.
 
real(r8), parameter, private ho_find_du =1.d-6
 Step size used for jacobian eval during high order find_cell.
 
integer(i4), parameter, private ho_find_nsteps =100
 Maximum number of steps during high order find_cell.
 
real(r8), parameter, private ho_find_tol =1.d-6
 Convergence tolerance for high order find_cell.
 
integer(i4), dimension(2, 6), parameter tet_ed =RESHAPE((/1,4, 2,4, 3,4, 2,3, 3,1, 1,2/), (/2,6/))
 Tetrahedron edge list.
 
integer(i4), dimension(3, 4), parameter tet_fc =RESHAPE((/2,3,4,3,1,4,1,2,4,1,2,3/), (/3,4/))
 Tetrahedron face list.
 
integer(i4), dimension(3, 4), parameter tet_fe =RESHAPE((/2,3,4, 1,3,5, 1,2,6, 4,5,6/), (/3,4/))
 Tetrahedron face edge list.
 
integer(i4), dimension(2, 3), parameter tri_ed =RESHAPE((/3,2,1,3,2,1/), (/2,3/))
 Triangle edge list.
 

Function/Subroutine Documentation

◆ tet_3d_grid()

subroutine tet_3d_grid ( integer(i4), intent(in)  order,
real(r8), dimension(:), intent(out), pointer  xnodes,
integer(i4), dimension(:,:,:), intent(out), pointer  inodesf,
integer(i4), dimension(:,:), intent(out), pointer  inodesc 
)

Needs docs.

◆ tetmesh_ctang()

subroutine tetmesh_ctang ( class(oft_tetmesh), intent(in)  self,
integer(i4), intent(in)  cell,
integer(i4), intent(in)  ind,
real(r8), dimension(:), intent(in)  f,
real(r8), dimension(3), intent(out)  tang 
)

Compute the curve tangent vector for a given edge on a cell.

If edge is not a global boundary edge the function returns with tang = 0

Note
The logical position in the cell must be on the chosen edge for this subroutine to return a meaningful result
Parameters
[in]selfMesh object
[in]cellIndex of cell
[in]indIndex of edge within cell
[in]fLogical coordinate in cell [4]
[out]tangUnit vector tangent to the edge [3]

◆ tetmesh_g2inv()

subroutine tetmesh_g2inv ( real(r8), dimension(3,4), intent(in)  jfull,
real(r8), dimension(6,10), intent(out)  g2op 
)

Needs docs.

◆ tetmesh_get_surf_map()

subroutine tetmesh_get_surf_map ( class(oft_tetmesh), intent(in)  self,
integer(i4), intent(in)  face,
integer(i4), intent(out)  cell,
integer(i4), dimension(3), intent(out)  lmap 
)

Get mapping between boundary and volume logical coordinates.

Parameters
[in]selfMesh object
[in]faceIndex of face on boundary mesh
[out]cellCell containing face
[out]lmapCoordinate mapping

◆ tetmesh_hessian()

subroutine tetmesh_hessian ( class(oft_tetmesh), intent(in)  self,
integer(i4), intent(in)  cell,
real(r8), dimension(:), intent(in)  f,
real(r8), dimension(:,:), intent(out)  g2op,
real(r8), dimension(:,:), intent(out)  k 
)

Compute the spatial hessian matrices for a given cell at a given logical position.

Parameters
[in]selfMesh object
[in]cellIndex of cell for evaulation
[in]fLogical coordinate in cell [4]
[out]g2opSecond order Jacobian matrix \( (\frac{\partial x_i}{\partial \lambda_l} \frac{\partial x_j}{\partial \lambda_k})^{-1} \)
[out]kGradient correction matrix \( \frac{\partial^2 x_i}{\partial \lambda_k \partial \lambda_l}\) [10,3]

◆ tetmesh_in_cell()

integer(i4) function tetmesh_in_cell ( class(oft_tetmesh), intent(in)  self,
real(r8), dimension(:), intent(in)  f,
real(r8), intent(in)  tol 
)

Test if logical position lies within the base cell.

Returns
Position f is inside the base cell?
Parameters
[in]selfMesh object
[in]fLogical coordinate to evaluate
[in]tolTolerance for test

◆ tetmesh_invert_cell()

subroutine tetmesh_invert_cell ( class(oft_tetmesh), intent(inout)  self,
integer(i4), intent(in)  cell 
)

Turn cell "inside out", used to ensure consistent orientations.

Parameters
[in,out]selfMesh object
[in]cellIndex of cell to invert

◆ tetmesh_jacinv()

subroutine tetmesh_jacinv ( real(r8), dimension(3,4), intent(in)  jfull,
real(r8), dimension(3,4), intent(out)  gop,
real(r8), intent(out)  jac 
)

Invert a 3x3 matrix.

Parameters
[in]jfullMatrix to invert
[out]gop\( A^{-1} \)
[out]jac|A|

◆ tetmesh_jacl()

subroutine tetmesh_jacl ( class(oft_tetmesh), intent(in)  self,
integer(i4), intent(in)  cell,
real(r8), dimension(:,:), intent(out)  gop,
real(r8), intent(out)  j 
)

Linear implementation of @tetmesh_jacobian.

Parameters
[in]selfMesh object
[in]cellIndex of cell for evaulation
[out]gopJacobian matrix \( (\frac{\partial x_i}{\partial \lambda_j})^{-1} \) [3,4]
[out]jJacobian of transformation from logical to physical coordinates

◆ tetmesh_jacobian()

subroutine tetmesh_jacobian ( class(oft_tetmesh), intent(in)  self,
integer(i4), intent(in)  cell,
real(r8), dimension(:), intent(in)  f,
real(r8), dimension(:,:), intent(out)  gop,
real(r8), intent(out)  j 
)

Compute the spatial jacobian matrix and its determinant for a given cell at a given logical position.

Parameters
[in]selfMesh object
[in]cellIndex of cell for evaulation
[in]fLogical coordinate in cell [4]
[out]gopJacobian matrix \( (\frac{\partial x_i}{\partial \lambda_j})^{-1} \) [3,4]
[out]jJacobian of transformation from logical to physical coordinates

◆ tetmesh_log2phys()

real(r8) function, dimension(3) tetmesh_log2phys ( class(oft_tetmesh), intent(in)  self,
integer, intent(in)  cell,
real(r8), dimension(:), intent(in)  f 
)

Map from logical to physical coordinates in a given cell.

Parameters
[in]selfMesh object
[in]cellIndex of cell for evaulation
[in]fLogical coordinate in cell [4]
Returns
Physical position [3]

◆ tetmesh_phys2log()

subroutine tetmesh_phys2log ( class(oft_tetmesh), intent(in)  self,
integer(i4), intent(in)  cell,
real(r8), dimension(3), intent(in)  pt,
real(r8), dimension(:), intent(out)  f 
)

Map from physical to logical coordinates in a given cell.

Parameters
[in]selfMesh object
[in]cellIndex of cell for evaulation
[in]ptPhysical position [3]
[out]fLogical coordinates within the cell [4]

◆ tetmesh_quad_rule()

subroutine tetmesh_quad_rule ( class(oft_tetmesh), intent(in)  self,
integer(i4), intent(in)  order,
type(oft_quad_type), intent(out)  quad_rule 
)

Retrieve suitable quadrature rule for mesh with given order.

Parameters
[in]selfMesh object
[in]orderDesired order of quadrature rule
[out]quad_ruleResulting quadrature rule

◆ tetmesh_set_order()

subroutine tetmesh_set_order ( class(oft_tetmesh), intent(inout)  self,
integer(i4), intent(in)  order 
)

Set maximum order of spatial mapping.

Parameters
[in,out]selfMesh object
[in]orderMaximum order of spatial mapping

◆ tetmesh_setup()

subroutine tetmesh_setup ( class(oft_tetmesh), intent(inout)  self,
integer(i4), intent(in)  cad_type 
)

Setup mesh with implementation specifics (cell_np, cell_ne, etc.)

Parameters
[in,out]selfMesh object
[in]cad_typeCAD/mesh interface ID number

◆ tetmesh_snormal()

subroutine tetmesh_snormal ( class(oft_tetmesh), intent(in)  self,
integer(i4), intent(in)  cell,
integer(i4), intent(in)  ind,
real(r8), dimension(:), intent(in)  f,
real(r8), dimension(3), intent(out)  norm 
)

Compute the surface normal vector for a given face on a cell.

If face is not a global boundary face the function returns with norm = 0

Note
The logical position in the cell must be on the chosen face for this subroutine, else an error will be thrown
Parameters
[in]selfMesh object
[in]cellIndex of cell
[in]indIndex of edge within cell
[in]fLogical coordinate in cell [4]
[out]normUnit vector normal to the face [3]

◆ tetmesh_surf_to_vol()

subroutine tetmesh_surf_to_vol ( class(oft_tetmesh), intent(in)  self,
real(r8), dimension(:), intent(in)  fsurf,
integer(i4), dimension(3), intent(in)  lmap,
real(r8), dimension(:), intent(out)  fvol 
)

Map between surface and volume logical coordinates.

Parameters
[in]selfMesh object
[in]fsurfSurface coordinates [3]
[in]lmapCoordinate mapping
[out]fvolVolume coordinates [4]

◆ tetmesh_tessellate()

subroutine tetmesh_tessellate ( class(oft_tetmesh), intent(in)  self,
real(r8), dimension(:,:), intent(out), pointer  rtmp,
integer(i4), dimension(:,:), intent(out), pointer  lctmp,
integer(i4), intent(in)  order 
)

Tessellate mesh onto lagrange FE nodes of specified order (usually for plotting)

Note
The maximum tessellation order currently supported is 4 (may be lower for certain mesh types).
Warning
Cell lists are returned with zero based indexing
Parameters
[in]selfMesh object
[out]rtmpTessellated point list [3,:]
[out]lctmpTessellated cell list [4,:]
[in]orderTessellation order

◆ tetmesh_tessellated_sizes()

integer(i4) function, dimension(2) tetmesh_tessellated_sizes ( class(oft_tetmesh), intent(in)  self)

Get sizes of arrays returned by tetmesh_tessellate.

Parameters
[in]selfMesh object
Returns
Array sizes following tessellation [np_tess,nc_tess]

◆ tetmesh_vlog()

subroutine tetmesh_vlog ( class(oft_tetmesh), intent(in)  self,
integer(i4), intent(in)  i,
real(r8), dimension(:), intent(out)  f 
)

Get position in logical space of vertex i

Parameters
[in]selfMesh object
[in]iVertex to locate
[out]fLogical coordinates of vertex i

Variable Documentation

◆ active_cell

integer(i4), private active_cell = 0
private

Active cell for high order find_cell.

◆ active_mesh

class(oft_tetmesh), pointer, private active_mesh => NULL()
private

Active mesh for high order find_cell.

◆ active_pt

real(r8), dimension(3), private active_pt = 0.d0
private

Active point for high order find_cell.

◆ ho_find_du

real(r8), parameter, private ho_find_du =1.d-6
private

Step size used for jacobian eval during high order find_cell.

◆ ho_find_nsteps

integer(i4), parameter, private ho_find_nsteps =100
private

Maximum number of steps during high order find_cell.

◆ ho_find_tol

real(r8), parameter, private ho_find_tol =1.d-6
private

Convergence tolerance for high order find_cell.

◆ tet_ed

integer(i4), dimension(2,6), parameter tet_ed =RESHAPE((/1,4, 2,4, 3,4, 2,3, 3,1, 1,2/), (/2,6/))

Tetrahedron edge list.

◆ tet_fc

integer(i4), dimension(3,4), parameter tet_fc =RESHAPE((/2,3,4,3,1,4,1,2,4,1,2,3/), (/3,4/))

Tetrahedron face list.

◆ tet_fe

integer(i4), dimension(3,4), parameter tet_fe =RESHAPE((/2,3,4, 1,3,5, 1,2,6, 4,5,6/), (/3,4/))

Tetrahedron face edge list.

◆ tri_ed

integer(i4), dimension(2,3), parameter tri_ed =RESHAPE((/3,2,1,3,2,1/), (/2,3/))

Triangle edge list.