The Open FUSION Toolkit 1.0.0-beta6
Modeling tools for plasma and fusion research and engineering
Loading...
Searching...
No Matches
oft_trimesh_type Module Reference

Detailed Description

Triangular boundary mesh definitions.

  • MPI seam
  • MPI I/O index
  • MPI global context
  • Base mesh linkage
  • MPI preallocated Send/Recv
  • High order tet representation
  • Mesh container

Global Tet variables

  • Tetrahedra edge and face lists
  • Grounding point position
Author
Chris Hansen
Date
June 2010

Data Types

type  oft_trimesh
 Triangular surface mesh type. More...

Functions/Subroutines

subroutine tri_2d_grid (order, xnodes, inodesf)
 Needs docs.
real(r8) function, dimension(3) trimesh_glogphys (self, cell, 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 spatial hessian matrices for a given cell at a given logical position.
integer(i4) function trimesh_in_cell (self, f, tol)
 Test if logical position lies within the base cell.
subroutine trimesh_invert_cell (self, cell)
 Turn cell "inside out", used to ensure consistent orientations.
subroutine trimesh_jacinv (a, c, j)
 Invert a 2x2 matrix.
subroutine trimesh_jacl (self, cell, f, gop, j)
 Linear implementation of @trimesh_jacobian.
subroutine trimesh_jacobian (self, cell, f, gop, j)
 Compute the spatial jacobian matrix and its determinant for a given cell at a given logical position.
subroutine trimesh_load (self, filename, ascii)
 Load trimesh from transfer file.
real(r8) function, dimension(3) trimesh_log2phys (self, cell, f)
 Map from logical to physical coordinates in a given cell.
subroutine trimesh_m3inv (a, c, j)
 Invert a 3x3 matrix.
subroutine trimesh_norm (self, cell, f, n)
 Get unit normal for surface at a given point in a given cell.
subroutine trimesh_phys2log (self, cell, pt, f)
 Map from physical to logical coordinates in a given cell.
subroutine trimesh_quad_rule (self, order, quad_rule)
 Retrieve suitable quadrature rule for triangular mesh with given order.
subroutine trimesh_save (self, filename, ascii)
 Save trimesh from transfer file.
subroutine trimesh_set_order (self, order)
 Set maximum order of spatial mapping.
subroutine trimesh_setup (self, cad_type, has_parent)
 Setup mesh with implementation specifics (cell_np, cell_ne, etc.).
subroutine trimesh_tang (self, cell, f, t)
 Get tangent basis set for surface at a given point in a given cell.
subroutine trimesh_tessellate (self, rtmp, lctmp, order)
 Tessellate mesh onto lagrange FE nodes of specified order (usually for plotting).
integer(i4) function, dimension(2) trimesh_tessellated_sizes (self)
 Get sizes of arrays returned by trimesh_tessellate.
subroutine trimesh_vlog (self, i, f)
 Get position in logical space of vertex i.

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 )

Needs docs.

◆ trimesh_glogphys()

real(r8) function, dimension(3) trimesh_glogphys ( class(oft_trimesh), intent(in) self,
integer, intent(in) cell,
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 object
[in]cellIndex of cell for evaulation
[in]jNeeds docs
[in]fLogical coordinate in cell [4]
Returns
\( \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 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]

◆ 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 )

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

◆ trimesh_invert_cell()

subroutine trimesh_invert_cell ( class(oft_trimesh), intent(inout) self,
integer(i4), intent(in) cell )

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

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

◆ 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) cell,
real(r8), dimension(:), intent(in) f,
real(r8), dimension(:,:), intent(out) gop,
real(r8), intent(out) j )

Linear implementation of @trimesh_jacobian.

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

◆ 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 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 [3]
[out]gopJacobian matrix \( (\frac{\partial x_i}{\partial \lambda_j})^{-1} \) [3,4]
[out]jJacobian of transformation from logical to physical coordinates

◆ trimesh_load()

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

Load trimesh from transfer file.

Parameters
[in,out]selfMesh object
[in]filenameFile to load mesh from
[in]asciiMesh stored in ASCII format?

◆ 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 in a given cell.

Parameters
[in]selfMesh object
[in]cellIndex of cell for evaulation
[in]fLogical coordinate in cell [4]
Returns
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) cell,
real(r8), dimension(:), intent(in) f,
real(r8), dimension(3), intent(out) n )

Get unit normal for surface at a given point in a given cell.

Parameters
[in]selfMesh object
[in]cellCell containing point
[in]fLogical coordinates in cell
[out]nUnit normal [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 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]

◆ 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 )

Retrieve suitable quadrature rule for triangular mesh with given order.

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

◆ trimesh_save()

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

Save trimesh from transfer file.

Parameters
[in]selfMesh object
[in]filenameFile to save mesh to
[in]asciiSave file in ASCII format?

◆ trimesh_set_order()

subroutine trimesh_set_order ( class(oft_trimesh), 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

◆ trimesh_setup()

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

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

Parameters
[in,out]selfMesh object
[in]cad_typeCAD/mesh interface ID number
[in]has_parentIs this mesh the/a surface of a volume mesh?

◆ trimesh_tang()

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

Get tangent basis set for surface at a given point in a given cell.

Parameters
[in]selfMesh object
[in]cellCell containing point
[in]fLogical coordinates in cell
[out]tUnit tangent basis set [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 )

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 [selfncp,:]
[in]orderTessellation order

◆ trimesh_tessellated_sizes()

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

Get sizes of arrays returned by trimesh_tessellate.

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

◆ trimesh_vlog()

subroutine trimesh_vlog ( class(oft_trimesh), 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

◆ tri_ed

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

Triangle edge list.