|
The Open FUSION Toolkit 1.0.0-beta6
Modeling tools for plasma and fusion research and engineering
|
Triangular boundary mesh definitions.
Global Tet variables
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. | |
| 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.
| 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
| [in] | self | Mesh object |
| [in] | cell | Index of cell for evaulation |
| [in] | j | Needs docs |
| [in] | f | Logical coordinate in cell [4] |
| 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.
| [in] | self | Mesh object |
| [in] | cell | Index of cell for evaulation |
| [in] | f | Logical coordinate in cell [4] |
| [out] | g2op | Second order Jacobian matrix \( (\frac{\partial x_i}{\partial \lambda_l} \frac{\partial x_j}{\partial \lambda_k})^{-1} \) |
| [out] | k | Gradient correction matrix \( \frac{\partial^2 x_i}{\partial \lambda_k \partial \lambda_l}\) [10,3] |
| 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.
f is inside the base cell? | [in] | self | Mesh object |
| [in] | f | Logical coordinate to evaluate |
| [in] | tol | Tolerance for test |
| subroutine trimesh_invert_cell | ( | class(oft_trimesh), intent(inout) | self, |
| integer(i4), intent(in) | cell | ||
| ) |
Turn cell "inside out", used to ensure consistent orientations.
| [in,out] | self | Mesh object |
| [in] | cell | Maximum order of spatial mapping |
| 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.
| [in] | a | Matrix to invert |
| [out] | c | \( A^{-1} \) |
| [out] | j | |A| |
| 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.
| [in] | self | Mesh object |
| [in] | cell | Index of cell for evaulation |
| [in] | f | Logical coordinate in cell [3] |
| [out] | gop | Jacobian matrix \( (\frac{\partial x_i}{\partial \lambda_j})^{-1} \) [3,4] |
| [out] | j | Jacobian of transformation from logical to physical coordinates |
| 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.
| [in] | self | Mesh object |
| [in] | cell | Index of cell for evaulation |
| [in] | f | Logical coordinate in cell [3] |
| [out] | gop | Jacobian matrix \( (\frac{\partial x_i}{\partial \lambda_j})^{-1} \) [3,4] |
| [out] | j | Jacobian of transformation from logical to physical coordinates |
| subroutine trimesh_load | ( | class(oft_trimesh), intent(inout) | self, |
| character(len=*), intent(in) | filename, | ||
| logical, intent(in), optional | ascii | ||
| ) |
Load trimesh from transfer file.
| [in,out] | self | Mesh object |
| [in] | filename | File to load mesh from |
| [in] | ascii | Mesh stored in ASCII format? |
| 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.
| [in] | self | Mesh object |
| [in] | cell | Index of cell for evaulation |
| [in] | f | Logical coordinate in cell [4] |
| 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.
| [in] | a | Matrix to invert |
| [out] | c | \( A^{-1} \) |
| [out] | j | |A| |
| 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.
| [in] | self | Mesh object |
| [in] | cell | Cell containing point |
| [in] | f | Logical coordinates in cell |
| [out] | n | Unit normal [3] |
| 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.
| [in] | self | Mesh object |
| [in] | cell | Index of cell for evaulation |
| [in] | pt | Physical position [3] |
| [out] | f | Logical coordinates within the cell [4] |
| 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.
| [in] | self | Mesh object |
| [in] | order | Desired order of quadrature rule |
| [out] | quad_rule | Resulting quadrature rule |
| subroutine trimesh_save | ( | class(oft_trimesh), intent(in) | self, |
| character(len=*), intent(in) | filename, | ||
| logical, intent(in), optional | ascii | ||
| ) |
Save trimesh from transfer file.
| [in] | self | Mesh object |
| [in] | filename | File to save mesh to |
| [in] | ascii | Save file in ASCII format? |
| subroutine trimesh_set_order | ( | class(oft_trimesh), intent(inout) | self, |
| integer(i4), intent(in) | order | ||
| ) |
Set maximum order of spatial mapping.
| [in,out] | self | Mesh object |
| [in] | order | Maximum order of spatial mapping |
| 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.)
| [in,out] | self | Mesh object |
| [in] | cad_type | CAD/mesh interface ID number |
| [in] | has_parent | Is this mesh the/a surface of a volume mesh? |
| 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.
| [in] | self | Mesh object |
| [in] | cell | Cell containing point |
| [in] | f | Logical coordinates in cell |
| [out] | t | Unit tangent basis set [3,2] |
| 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)
| [in] | self | Mesh object |
| [out] | rtmp | Tessellated point list [3,:] |
| [out] | lctmp | Tessellated cell list [selfncp,:] |
| [in] | order | Tessellation order |
| integer(i4) function, dimension(2) trimesh_tessellated_sizes | ( | class(oft_trimesh), intent(in) | self | ) |
Get sizes of arrays returned by trimesh_tessellate.
| [in] | self | Mesh object |
| 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
| [in] | self | Mesh object |
| [in] | i | Vertex to locate |
| [out] | f | Logical coordinates of vertex i |
| integer(i4), dimension(2,3), parameter tri_ed =RESHAPE((/3,2,1,3,2,1/), (/2,3/)) |
Triangle edge list.