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

Detailed Description

Triangular boundary mesh definitions.

Global Tet variables

Author
Chris Hansen
Date
June 2010

Data Types

type  oft_trimesh
 Triangular boundary mesh type. More...
 

Functions/Subroutines

subroutine tri_2d_grid (order, xnodes, inodesf)
 
integer(i4) function, dimension(2) trimesh_get_io_sizes (self)
 Get variable sizes following tessellation.
 
real(r8) function, dimension(3) trimesh_glogphys (self, face, j, f)
 Compute the partial derivative of the physical coordinates with a specific logical coordinate.
 
subroutine trimesh_hessian (self, cell, f, g2op, k)
 Compute the second order jacobians for a grid cell.
 
integer(i4) function trimesh_in_cell (self, f, tol)
 Logical locations of vertices.
 
subroutine trimesh_invert_face (self, i)
 Invert cell to produce positive volume.
 
subroutine trimesh_jacinv (a, c, j)
 Invert a 2x2 matrix.
 
subroutine trimesh_jacl (self, i, f, gop, j)
 Compute the jacobian matrix and its determinant for a linear element.
 
subroutine trimesh_jacobian (self, cell, f, gop, j)
 Compute the jacobian matrix and its determinant for a grid cell.
 
subroutine trimesh_load (self, filename)
 Load trimesh from transfer file.
 
real(r8) function, dimension(3) trimesh_log2phys (self, cell, f)
 Map from logical to physical coordinates.
 
subroutine trimesh_m3inv (a, c, j)
 Invert a 3x3 matrix.
 
subroutine trimesh_norm (self, i, f, n)
 Compute the unit normal vector to a face.
 
subroutine trimesh_phys2log (self, cell, pt, f)
 Map from physical to logical coordinates for a linear element.
 
subroutine trimesh_quad_rule (self, order, quad_rule)
 Create quadrature rule for tetrahedra.
 
subroutine trimesh_save (self, filename)
 Load trimesh from transfer file.
 
subroutine trimesh_set_order (self, order)
 Load trimesh from transfer file.
 
subroutine trimesh_setup (self, cad_type, has_parent)
 Load trimesh from transfer file.
 
subroutine trimesh_tang (self, i, f, t)
 Compute a orthonormal set of axis tangential to a face.
 
subroutine trimesh_tessellate (self, rtmp, lctmp, order)
 Driver for order specific tessellation subroutines.
 
subroutine trimesh_vlog (self, i, f)
 Logical locations of vertices.
 

Variables

integer(i4), dimension(2, 3), parameter tri_ed =RESHAPE((/3,2,1,3,2,1/), (/2,3/))
 Triangle edge list.
 

Function/Subroutine Documentation

◆ tri_2d_grid()

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

◆ trimesh_get_io_sizes()

integer(i4) function, dimension(2) trimesh_get_io_sizes ( class(oft_trimesh), intent(in)  self)

Get variable sizes following tessellation.

◆ trimesh_glogphys()

real(r8) function, dimension(3) trimesh_glogphys ( class(oft_trimesh), intent(in)  self,
integer, intent(in)  face,
integer, intent(in)  j,
real(r8), dimension(3), intent(in)  f 
)

Compute the partial derivative of the physical coordinates with a specific logical coordinate.

Driver function calls mapping specific function depending on mesh order.

Parameters
[in]selfMesh containing cell
[in]iIndex of cell for evaulation
[in]kLogical coordinate for differentiation
[in]fLogical coordinate in cell [4]
Returns
pt \( \frac{\partial r}{\partial f_k} \) [3]

◆ trimesh_hessian()

subroutine trimesh_hessian ( class(oft_trimesh), 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 second order jacobians for a grid cell.

Parameters
[in]selfMesh containing cell
[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} \) [6,6]
[out]KGradient correction matrix \( \frac{\partial^2 x_i}{\partial \lambda_k \partial \lambda_l}\) [6,3]

◆ trimesh_in_cell()

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

Logical locations of vertices.

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

◆ trimesh_invert_face()

subroutine trimesh_invert_face ( class(oft_trimesh), intent(inout)  self,
integer(i4), intent(in)  i 
)

Invert cell to produce positive volume.

◆ trimesh_jacinv()

subroutine trimesh_jacinv ( real(r8), dimension(2,2), intent(in)  a,
real(r8), dimension(2,2), intent(out)  c,
real(r8), intent(out)  j 
)

Invert a 2x2 matrix.

Parameters
[in]AMatrix to invert
[out]C\( A^{-1} \)
[out]j|A|

◆ trimesh_jacl()

subroutine trimesh_jacl ( class(oft_trimesh), intent(in)  self,
integer(i4), intent(in)  i,
real(r8), dimension(3), intent(in)  f,
real(r8), dimension(3,3), intent(out)  gop,
real(r8), intent(out)  j 
)

Compute the jacobian matrix and its determinant for a linear element.

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

◆ trimesh_jacobian()

subroutine trimesh_jacobian ( class(oft_trimesh), 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 jacobian matrix and its determinant for a grid cell.

Driver function calls mapping specific function depending on mesh order.

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

◆ trimesh_load()

subroutine trimesh_load ( class(oft_trimesh), intent(inout)  self,
character(len=*), intent(in)  filename 
)

Load trimesh from transfer file.

◆ trimesh_log2phys()

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

Map from logical to physical coordinates.

Driver function calls mapping specific function depending on mesh order.

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

◆ trimesh_m3inv()

subroutine trimesh_m3inv ( real(8), dimension(3,3), intent(in)  a,
real(8), dimension(3,3), intent(out)  c,
real(8), intent(out)  j 
)

Invert a 3x3 matrix.

Parameters
[in]AMatrix to invert
[out]C\( A^{-1} \)
[out]j|A|

◆ trimesh_norm()

subroutine trimesh_norm ( class(oft_trimesh), intent(in), target  self,
integer(i4), intent(in)  i,
real(r8), dimension(:), intent(in)  f,
real(r8), dimension(3), intent(out)  n 
)

Compute the unit normal vector to a face.

Parameters
[in]selfMesh containing face
[in]iIndex of face for evaulation
[in]fLogical coordinate on face [3]
[out]nNormal vector [3]

◆ trimesh_phys2log()

subroutine trimesh_phys2log ( class(oft_trimesh), 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 for a linear element.

Parameters
[in]selfMesh containing cell
[in]cellIndex of cell for evaulation
[in]ptPhysical position [3]
Returns
f Logical coordinates within the cell [4]

◆ trimesh_quad_rule()

subroutine trimesh_quad_rule ( class(oft_trimesh), intent(in)  self,
integer(i4), intent(in)  order,
type(oft_quad_type), intent(out)  quad_rule 
)

Create quadrature rule for tetrahedra.

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

◆ trimesh_save()

subroutine trimesh_save ( class(oft_trimesh), intent(in)  self,
character(len=*), intent(in)  filename 
)

Load trimesh from transfer file.

◆ trimesh_set_order()

subroutine trimesh_set_order ( class(oft_trimesh), intent(inout)  self,
integer(i4), intent(in)  order 
)

Load trimesh from transfer file.

◆ trimesh_setup()

subroutine trimesh_setup ( class(oft_trimesh), intent(inout)  self,
integer(i4), intent(in)  cad_type,
logical, intent(in)  has_parent 
)

Load trimesh from transfer file.

◆ trimesh_tang()

subroutine trimesh_tang ( class(oft_trimesh), intent(in), target  self,
integer(i4), intent(in)  i,
real(r8), dimension(:), intent(in)  f,
real(r8), dimension(3,2), intent(out)  t 
)

Compute a orthonormal set of axis tangential to a face.

Parameters
[in]selfMesh containing face
[in]iIndex of face for evaulation
[in]fLogical coordinate on face [3]
[out]tOrthonormal vectors [3,2]

◆ trimesh_tessellate()

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

Driver for order specific tessellation subroutines.

Note
Cell lists are returned with zero based indexing.
Warning
The maximum tessellation order currently supported is 4.
Parameters
[in]selfMesh to tessellate
[out]rtmpTessellated mesh points [3,np_tess]
[out]lctmpTessellated face list [3,nf_tess]
[in]orderDesired tessellation order

◆ trimesh_vlog()

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

Logical locations of vertices.

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

Variable Documentation

◆ tri_ed

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

Triangle edge list.