|
The Open FUSION Toolkit 1.0.0-beta5
Modeling tools for plasma and fusion research and engineering
|
Hexhedral mesh definitions.
Data Types | |
| type | oft_hexmesh |
| Tetrahedral Mesh type. More... | |
Functions/Subroutines | |
| subroutine | hex_3d_grid (order, xnodes, inodesp, inodese, inodesf, inodesc) |
| pure real(r8) function, dimension(6) | hex_get_bary (flog) |
| pure real(r8) function, dimension(3) | hex_get_bary_cgop (cglog, i, j) |
| pure real(r8) function, dimension(3, 6) | hex_get_bary_gop (glog) |
| subroutine | hex_grid_forient (oflag, order, finds) |
| subroutine | hex_grid_forient_inv (oflag, order, finds) |
| subroutine | hexmesh_ctang (self, i, ind, f, tang) |
| Compute the curve tangent vector for a given edge on a cell. | |
| subroutine | hexmesh_g2inv (jac, a) |
| Invert a 3x3 matrix. | |
| integer(i4) function, dimension(2) | hexmesh_get_io_sizes (self) |
| Get variable sizes following tessellation. | |
| subroutine | hexmesh_get_surf_map (self, face, cell, lmap) |
| Get mapping between boundary and volume logical coordinates. | |
| subroutine | hexmesh_hessian (self, i, f, g2op, k) |
| Compute the second order jacobians for a linear element. | |
| integer(i4) function | hexmesh_in_cell (self, f, tol) |
| Logical locations of vertices. | |
| subroutine | hexmesh_invert_cell (self, i) |
| Invert cell to produce positive volume. | |
| subroutine | hexmesh_jacinv (a, c, j) |
| Invert a 3x3 matrix. | |
| subroutine | hexmesh_jacobian (self, cell, f, gop, j) |
| Compute the jacobian matrix and its determinant for a quadratic element. | |
| real(r8) function, dimension(3) | hexmesh_log2phys (self, cell, f) |
| Map from logical to physical coordinates for a quadratic element. | |
| subroutine | hexmesh_phys2log (self, i, pt, f) |
| Map from physical to logical coordinates. | |
| real(r8) function, dimension(4) | hexmesh_phys2logho (self, i, pt) |
| Map from physical to logical coordinates for general high order elements. | |
| subroutine | hexmesh_quad_rule (self, order, quad_rule) |
| Create quadrature rule for tetrahedra. | |
| subroutine | hexmesh_set_order (self, order) |
| Load trimesh from transfer file. | |
| subroutine | hexmesh_setup (self, cad_type) |
| Load trimesh from transfer file. | |
| subroutine | hexmesh_snormal (self, i, ind, f, norm) |
| Compute the surface normal vector for a given face on a cell. | |
| subroutine | hexmesh_surf_to_vol (self, fsurf, lmap, fvol) |
| Map between surface and volume logical coordinates. | |
| subroutine | hexmesh_tessellate (self, rtmp, lctmp, order) |
| Driver for order specific tessellation subroutines. | |
| subroutine | hexmesh_vlog (self, i, f) |
| Logical locations of vertices. | |
| subroutine | tm_findcell_error (m, n, uv, err, iflag) |
| Evalute the error between a logical point and the current active point. | |
Variables | |
| integer(i4), private | active_cell = 0 |
| Active cell for high order find_cell. | |
| class(oft_hexmesh), 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. | |
| integer(i4), dimension(2, 12), parameter | hex_bary_ecoords = RESHAPE((/ 5,3, 2,4, 3,5, 4,2, 1,6, 1,6, 1,6, 1,6, 5,3, 2,4, 3,5, 4,2/), (/2,12/)) |
| integer(i4), dimension(2, 12), parameter | hex_bary_efcoords = RESHAPE((/ 1,2, 1,3, 1,4, 1,5, 2,5, 2,3, 3,4, 4,5, 2,6, 3,6, 4,6, 5,6/), (/2,12/)) |
| integer(i4), dimension(4, 6), parameter | hex_bary_fcoords = RESHAPE((/ 5,2,3,4, 1,5,6,3, 1,2,6,4, 1,3,6,5, 2,1,4,6, 2,5,4,3/), (/4,6/)) |
| integer(i4), dimension(6), parameter | hex_bary_map = (/-3,-2,1, 2,-1,3/) |
| integer(i4), dimension(3, 8), parameter | hex_bary_pfcoords = RESHAPE((/ 1,2,5, 1,2,3, 1,3,4, 1,4,5, 2,5,6, 2,3,6, 3,4,6, 4,5,6/), (/3,8/)) |
| integer(i4), dimension(2, 12), parameter | hex_ed =RESHAPE((/ 1,2, 2,3, 3,4, 4,1, 1,5, 2,6, 3,7, 4,8, 5,6, 6,7, 7,8, 8,5/), (/2,12/)) |
| integer(i4), dimension(4, 6), parameter | hex_fc =RESHAPE((/ 1,2,3,4, 1,5,6,2, 2,6,7,3, 3,7,8,4, 1,4,8,5, 5,8,7,6/), (/4,6/)) |
| integer(i4), dimension(4, 6), parameter | hex_fe =RESHAPE((/ 1,2,3,4, 5,9,-6,-1, 6,10,-7,-2, 7,11,-8,-3, -4,8,12,-5, -12,-11,-10,-9/), (/4,6/)) |
| 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(3, 8), parameter | inodes1p = RESHAPE((/ 1,1,1, 2,1,1, 2,2,1, 1,2,1, 1,1,2, 2,1,2, 2,2,2, 1,2,2/), (/3,8/)) |
| integer(i4), dimension(3, 12), parameter | inodes2e = RESHAPE((/ 2,1,1, 3,2,1, 2,3,1, 1,2,1, 1,1,2, 3,1,2, 3,3,2, 1,3,2, 2,1,3, 3,2,3, 2,3,3, 1,2,3/), (/3,12/)) |
| integer(i4), dimension(3, 6), parameter | inodes2f = RESHAPE((/ 2,2,1, 2,1,2, 3,2,2, 2,3,2, 1,2,2, 2,2,3/), (/3,6/)) |
| integer(i4), dimension(3, 8), parameter | inodes2p = RESHAPE((/ 1,1,1, 3,1,1, 3,3,1, 1,3,1, 1,1,3, 3,1,3, 3,3,3, 1,3,3/), (/3,8/)) |
| integer(i4), dimension(3, 8), parameter | inodesp_base = RESHAPE((/ 0,0,0, 1,0,0, 1,1,0, 0,1,0, 0,0,1, 1,0,1, 1,1,1, 0,1,1/), (/3,8/)) |
| integer(i4), dimension(2, 4), parameter | quad_ed =RESHAPE((/1,2, 2,3, 3,4, 4,1/), (/2,4/)) |
| Quad edge list. | |
| subroutine hex_3d_grid | ( | integer(i4), intent(in) | order, |
| real(r8), dimension(:), intent(out), pointer | xnodes, | ||
| integer(i4), dimension(3,8), intent(out) | inodesp, | ||
| integer(i4), dimension(:,:,:), intent(out), pointer | inodese, | ||
| integer(i4), dimension(:,:,:), intent(out), pointer | inodesf, | ||
| integer(i4), dimension(:,:), intent(out), pointer | inodesc | ||
| ) |
| pure real(r8) function, dimension(3) hex_get_bary_cgop | ( | real(r8), dimension(3,3), intent(in) | cglog, |
| integer(i4), intent(in) | i, | ||
| integer(i4), intent(in) | j | ||
| ) |
| pure real(r8) function, dimension(3,6) hex_get_bary_gop | ( | real(r8), dimension(:,:), intent(in) | glog | ) |
| subroutine hex_grid_forient | ( | integer(i4), intent(in) | oflag, |
| integer(i4), intent(in) | order, | ||
| integer(i4), dimension(:), intent(inout) | finds | ||
| ) |
| subroutine hex_grid_forient_inv | ( | integer(i4), intent(in) | oflag, |
| integer(i4), intent(in) | order, | ||
| integer(i4), dimension(:), intent(inout) | finds | ||
| ) |
| subroutine hexmesh_ctang | ( | class(oft_hexmesh), intent(in) | self, |
| integer(i4), intent(in) | i, | ||
| 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
| [in] | self | Mesh containing face |
| [in] | i | Index of cell |
| [in] | ind | Index of edge within cell |
| [in] | f | Logical coordinate in cell [4] |
| [out] | tang | Unit vector tangent to the edge [3] |
| subroutine hexmesh_g2inv | ( | real(r8), dimension(3,3), intent(in) | jac, |
| real(r8), dimension(6,6), intent(out) | a | ||
| ) |
Invert a 3x3 matrix.
| [in] | A | Matrix to invert |
| [out] | C | \( A^{-1} \) |
| [out] | j | |A| |
| integer(i4) function, dimension(2) hexmesh_get_io_sizes | ( | class(oft_hexmesh), intent(in) | self | ) |
Get variable sizes following tessellation.
| subroutine hexmesh_get_surf_map | ( | class(oft_hexmesh), 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.
| [in] | face | Index of face on boundary mesh |
| [out] | cell | Cell containing face |
| [out] | lmap | Coordinate mapping |
| subroutine hexmesh_hessian | ( | class(oft_hexmesh), intent(in) | self, |
| integer(i4), intent(in) | i, | ||
| 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 linear element.
| [in] | self | Mesh containing cell |
| [in] | i | Index of cell for evaulation |
| [out] | g2op | Second order Jacobian matrix \( (\frac{\partial x_i}{\partial \lambda_l} \frac{\partial x_j}{\partial \lambda_k})^{-1} \) [6,10] |
| [out] | K | Gradient correction matrix \( \frac{\partial^2 x_i}{\partial \lambda_k \partial \lambda_l}\) [10,3] |
| integer(i4) function hexmesh_in_cell | ( | class(oft_hexmesh), intent(in) | self, |
| real(r8), dimension(:), intent(in) | f, | ||
| real(r8), intent(in) | tol | ||
| ) |
Logical locations of vertices.
| [in] | self | Mesh containing cell |
| [in] | i | Index of cell for evaulation |
| [in] | f | Logical coordinate in cell [4] |
| subroutine hexmesh_invert_cell | ( | class(oft_hexmesh), intent(inout) | self, |
| integer(i4), intent(in) | i | ||
| ) |
Invert cell to produce positive volume.
| subroutine hexmesh_jacinv | ( | real(r8), dimension(3,3), intent(in) | a, |
| real(r8), dimension(3,3), intent(out) | c, | ||
| real(r8), intent(out) | j | ||
| ) |
Invert a 3x3 matrix.
| [in] | A | Matrix to invert |
| [out] | C | \( A^{-1} \) |
| [out] | j | |A| |
| subroutine hexmesh_jacobian | ( | class(oft_hexmesh), 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 quadratic element.
| [in] | self | Mesh containing cell |
| [in] | i | Index of cell for evaulation |
| [in] | f | Logical coordinate in cell [4] |
| [out] | gop | Jacobian matrix \( (\frac{\partial x_i}{\partial \lambda_j})^{-1} \) [3,4] |
| [out] | j | Jacobian of transformation from logical to physical coordinates |
| real(r8) function, dimension(3) hexmesh_log2phys | ( | class(oft_hexmesh), intent(in) | self, |
| integer(i4), intent(in) | cell, | ||
| real(r8), dimension(:), intent(in) | f | ||
| ) |
Map from logical to physical coordinates for a quadratic element.
| [in] | self | Mesh containing cell |
| [in] | i | Index of cell for evaulation |
| [in] | f | Logical coordinate in cell [4] |
| subroutine hexmesh_phys2log | ( | class(oft_hexmesh), intent(in) | self, |
| integer(i4), intent(in) | i, | ||
| real(r8), dimension(3), intent(in) | pt, | ||
| real(r8), dimension(:), intent(out) | f | ||
| ) |
Map from physical to logical coordinates.
Driver function calls mapping specific function depending on mesh order.
| [in] | self | Mesh containing cell |
| [in] | i | Index of cell for evaulation |
| [in] | pt | Physical position [3] |
| real(r8) function, dimension(4) hexmesh_phys2logho | ( | class(oft_hexmesh), intent(in), target | self, |
| integer(i4), intent(in) | i, | ||
| real(r8), dimension(3), intent(in) | pt | ||
| ) |
Map from physical to logical coordinates for general high order elements.
The MINPACK package is used with step size given by ho_find_du. The convergence tolerance is set by the variable ho_find_tol.
| [in] | self | Mesh containing cell |
| [in] | i | Index of cell for evaulation |
| [in] | pt | Physical position [3] |
| subroutine hexmesh_quad_rule | ( | class(oft_hexmesh), intent(in) | self, |
| integer(i4), intent(in) | order, | ||
| type(oft_quad_type), intent(out) | quad_rule | ||
| ) |
Create quadrature rule for tetrahedra.
| [in] | self | Mesh containing cell |
| [in] | i | Index of cell for evaulation |
| [in] | f | Logical coordinate in cell [4] |
| subroutine hexmesh_set_order | ( | class(oft_hexmesh), intent(inout) | self, |
| integer(i4), intent(in) | order | ||
| ) |
Load trimesh from transfer file.
| subroutine hexmesh_setup | ( | class(oft_hexmesh), intent(inout) | self, |
| integer(i4), intent(in) | cad_type | ||
| ) |
Load trimesh from transfer file.
| subroutine hexmesh_snormal | ( | class(oft_hexmesh), intent(in) | self, |
| integer(i4), intent(in) | i, | ||
| 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
| [in] | self | Mesh containing face |
| [in] | i | Index of cell |
| [in] | ind | Index of face within cell |
| [in] | f | Logical coordinate in cell [4] |
| [out] | norm | Unit vector normal to the face [3] |
| subroutine hexmesh_surf_to_vol | ( | class(oft_hexmesh), 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.
| [in] | fsurf | Surface coordinates [2] |
| [in] | lmap | Coordinate mapping |
| [out] | fvol | Volume coordinates [3] |
| subroutine hexmesh_tessellate | ( | class(oft_hexmesh), 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.
| [in] | self | Mesh to tessellate |
| [out] | rtmp | Tessellated mesh points [3,np_tess] |
| [out] | lctmp | Tessellated cell list [4,nc_tess] |
| [in] | order | Desired tessellation order |
| subroutine hexmesh_vlog | ( | class(oft_hexmesh), intent(in) | self, |
| integer(i4), intent(in) | i, | ||
| real(r8), dimension(:), intent(out) | f | ||
| ) |
Logical locations of vertices.
| [in] | self | Mesh containing cell |
| [in] | i | Index of cell for evaulation |
| [in] | f | Logical coordinate in cell [4] |
| subroutine tm_findcell_error | ( | integer(i4), intent(in) | m, |
| integer(i4), intent(in) | n, | ||
| real(r8), dimension(n), intent(in) | uv, | ||
| real(r8), dimension(m), intent(out) | err, | ||
| integer(i4), intent(inout) | iflag | ||
| ) |
Evalute the error between a logical point and the current active point.
| [in] | m | Number of spatial dimensions (3) |
| [in] | n | Number of parametric dimensions (3) |
| [in] | uv | Parametric possition [n] |
| [out] | err | Error vector between current and desired point [3] |
| [in,out] | iflag | Unused flag |
|
private |
Active cell for high order find_cell.
|
private |
Active mesh for high order find_cell.
|
private |
Active point for high order find_cell.
| integer(i4), dimension(2,12), parameter hex_bary_ecoords = RESHAPE((/ 5,3, 2,4, 3,5, 4,2, 1,6, 1,6, 1,6, 1,6, 5,3, 2,4, 3,5, 4,2/), (/2,12/)) |
| integer(i4), dimension(2,12), parameter hex_bary_efcoords = RESHAPE((/ 1,2, 1,3, 1,4, 1,5, 2,5, 2,3, 3,4, 4,5, 2,6, 3,6, 4,6, 5,6/), (/2,12/)) |
| integer(i4), dimension(4,6), parameter hex_bary_fcoords = RESHAPE((/ 5,2,3,4, 1,5,6,3, 1,2,6,4, 1,3,6,5, 2,1,4,6, 2,5,4,3/), (/4,6/)) |
| integer(i4), dimension(6), parameter hex_bary_map = (/-3,-2,1, 2,-1,3/) |
| integer(i4), dimension(3,8), parameter hex_bary_pfcoords = RESHAPE((/ 1,2,5, 1,2,3, 1,3,4, 1,4,5, 2,5,6, 2,3,6, 3,4,6, 4,5,6/), (/3,8/)) |
| integer(i4), dimension(2,12), parameter hex_ed =RESHAPE((/ 1,2, 2,3, 3,4, 4,1, 1,5, 2,6, 3,7, 4,8, 5,6, 6,7, 7,8, 8,5/), (/2,12/)) |
| integer(i4), dimension(4,6), parameter hex_fc =RESHAPE((/ 1,2,3,4, 1,5,6,2, 2,6,7,3, 3,7,8,4, 1,4,8,5, 5,8,7,6/), (/4,6/)) |
| integer(i4), dimension(4,6), parameter hex_fe =RESHAPE((/ 1,2,3,4, 5,9,-6,-1, 6,10,-7,-2, 7,11,-8,-3, -4,8,12,-5, -12,-11,-10,-9/), (/4,6/)) |
|
private |
Step size used for jacobian eval during high order find_cell.
|
private |
Maximum number of steps during high order find_cell.
|
private |
Convergence tolerance for high order find_cell.
| integer(i4), dimension(3,8), parameter inodes1p = RESHAPE((/ 1,1,1, 2,1,1, 2,2,1, 1,2,1, 1,1,2, 2,1,2, 2,2,2, 1,2,2/), (/3,8/)) |
| integer(i4), dimension(3,12), parameter inodes2e = RESHAPE((/ 2,1,1, 3,2,1, 2,3,1, 1,2,1, 1,1,2, 3,1,2, 3,3,2, 1,3,2, 2,1,3, 3,2,3, 2,3,3, 1,2,3/), (/3,12/)) |
| integer(i4), dimension(3,6), parameter inodes2f = RESHAPE((/ 2,2,1, 2,1,2, 3,2,2, 2,3,2, 1,2,2, 2,2,3/), (/3,6/)) |
| integer(i4), dimension(3,8), parameter inodes2p = RESHAPE((/ 1,1,1, 3,1,1, 3,3,1, 1,3,1, 1,1,3, 3,1,3, 3,3,3, 1,3,3/), (/3,8/)) |
| integer(i4), dimension(3,8), parameter inodesp_base = RESHAPE((/ 0,0,0, 1,0,0, 1,1,0, 0,1,0, 0,0,1, 1,0,1, 1,1,1, 0,1,1/), (/3,8/)) |
| integer(i4), dimension(2,4), parameter quad_ed =RESHAPE((/1,2, 2,3, 3,4, 4,1/), (/2,4/)) |
Quad edge list.