|
The Open FUSION Toolkit 1.0.0-beta5
Modeling tools for plasma and fusion research and engineering
|
Quadralateral boundary mesh definitions.
Data Types | |
| type | oft_quadmesh |
| Quadralateral boundary mesh type. More... | |
Functions/Subroutines | |
| subroutine | quad_2d_grid (order, xnodes, inodesp, inodese, inodesf) |
| pure real(r8) function, dimension(4) | quad_get_bary (flog) |
| subroutine | quad_grid_orient (oflag, order, finds) |
| subroutine | quad_grid_orient_inv (oflag, order, finds) |
| integer(i4) function, dimension(2) | quadmesh_get_io_sizes (self) |
| Get variable sizes following tessellation. | |
| real(r8) function, dimension(3) | quadmesh_glogphys (self, face, j, f) |
| Compute the partial derivative of the physical coordinates with a specific logical coordinate. | |
| subroutine | quadmesh_hessian (self, cell, f, g2op, k) |
| Compute the second order jacobians for a grid cell. | |
| integer(i4) function | quadmesh_in_cell (self, f, tol) |
| Logical locations of vertices. | |
| subroutine | quadmesh_invert_face (self, i) |
| Invert cell to produce positive volume. | |
| subroutine | quadmesh_jacinv (a, c, j) |
| Invert a 2x2 matrix. | |
| subroutine | quadmesh_jacobian (self, cell, f, gop, j) |
| Compute the jacobian matrix and its determinant for a grid cell. | |
| subroutine | quadmesh_load (self, filename) |
| Load quadmesh from transfer file. | |
| real(r8) function, dimension(3) | quadmesh_log2phys (self, cell, f) |
| Map from logical to physical coordinates. | |
| subroutine | quadmesh_norm (self, i, f, n) |
| Compute the unit normal vector to a face. | |
| subroutine | quadmesh_phys2log (self, cell, pt, f) |
| Map from physical to logical coordinates. | |
| subroutine | quadmesh_quad_rule (self, order, quad_rule) |
| Create quadrature rule for tetrahedra. | |
| subroutine | quadmesh_save (self, filename) |
| Load trimesh from transfer file. | |
| subroutine | quadmesh_set_order (self, order) |
| Load trimesh from transfer file. | |
| subroutine | quadmesh_setup (self, cad_type, has_parent) |
| Load quadmesh from transfer file. | |
| subroutine | quadmesh_tang (self, i, f, t) |
| Compute a orthonormal set of axis tangential to a face. | |
| subroutine | quadmesh_tessellate (self, rtmp, lctmp, order) |
| Driver for order specific tessellation subroutines. | |
| subroutine | quadmesh_vlog (self, i, f) |
| Logical locations of vertices. | |
Variables | |
| integer(i4), dimension(2, 4), parameter | inodes1p = RESHAPE((/ 1,1, 2,1, 2,2, 1,2/), (/2,4/)) |
| integer(i4), dimension(2, 4), parameter | inodes2e = RESHAPE((/ 2,1, 3,2, 2,3, 1,2/), (/2,4/)) |
| integer(i4), dimension(2), parameter | inodes2f = (/2,2/) |
| integer(i4), dimension(2, 4), parameter | inodes2p = RESHAPE((/ 1,1, 3,1, 3,3, 1,3/), (/2,4/)) |
| integer(i4), dimension(2, 4), parameter | inodesp_base = RESHAPE((/ 0,0, 1,0, 1,1, 0,1/), (/2,4/)) |
| integer(i4), dimension(4), parameter | quad_bary_map = [-2,1,2,-1] |
| integer(i4), dimension(2, 4), parameter | quad_ed =RESHAPE((/1,2, 2,3, 3,4, 4,1/), (/2,4/)) |
| Quad edge list. | |
| subroutine quad_2d_grid | ( | integer(i4), intent(in) | order, |
| real(r8), dimension(:), intent(out), pointer | xnodes, | ||
| integer(i4), dimension(2,4), intent(out) | inodesp, | ||
| integer(i4), dimension(:,:,:), intent(out), pointer | inodese, | ||
| integer(i4), dimension(:,:), intent(out), pointer | inodesf | ||
| ) |
| subroutine quad_grid_orient | ( | integer(i4), intent(in) | oflag, |
| integer(i4), intent(in) | order, | ||
| integer(i4), dimension(:), intent(inout) | finds | ||
| ) |
| subroutine quad_grid_orient_inv | ( | integer(i4), intent(in) | oflag, |
| integer(i4), intent(in) | order, | ||
| integer(i4), dimension(:), intent(inout) | finds | ||
| ) |
| integer(i4) function, dimension(2) quadmesh_get_io_sizes | ( | class(oft_quadmesh), intent(in) | self | ) |
Get variable sizes following tessellation.
| real(r8) function, dimension(3) quadmesh_glogphys | ( | class(oft_quadmesh), intent(in) | self, |
| integer, intent(in) | face, | ||
| integer, intent(in) | j, | ||
| real(r8), dimension(2), 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 containing cell |
| [in] | i | Index of cell for evaulation |
| [in] | k | Logical coordinate for differentiation |
| [in] | f | Logical coordinate in cell [4] |
| subroutine quadmesh_hessian | ( | class(oft_quadmesh), 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.
| [in] | self | Mesh containing cell |
| [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} \) [6,6] |
| [out] | K | Gradient correction matrix \( \frac{\partial^2 x_i}{\partial \lambda_k \partial \lambda_l}\) [6,3] |
| integer(i4) function quadmesh_in_cell | ( | class(oft_quadmesh), 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 quadmesh_invert_face | ( | class(oft_quadmesh), intent(inout) | self, |
| integer(i4), intent(in) | i | ||
| ) |
Invert cell to produce positive volume.
| subroutine quadmesh_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 quadmesh_jacobian | ( | class(oft_quadmesh), 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.
| [in] | self | Mesh containing cell |
| [in] | cell | Index of cell for evaulation |
| [in] | f | Logical coordinate in cell [4] |
| [out] | gop | Jacobian matrix ( \( \frac{\partial x_i}{\partial \lambda_j} \) ) [3,4] (optional) |
| [out] | j | Jacobian of transformation from logical to physical coordinates (optional) |
| subroutine quadmesh_load | ( | class(oft_quadmesh), intent(inout) | self, |
| character(len=*), intent(in) | filename | ||
| ) |
Load quadmesh from transfer file.
| real(r8) function, dimension(3) quadmesh_log2phys | ( | class(oft_quadmesh), 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.
| [in] | self | Mesh containing cell |
| [in] | cell | Index of cell for evaulation |
| [in] | f | Logical coordinate in cell [4] |
| subroutine quadmesh_norm | ( | class(oft_quadmesh), 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.
| [in] | self | Mesh containing face |
| [in] | i | Index of face for evaulation |
| [in] | f | Logical coordinate on face [3] |
| [out] | n | Normal vector [3] |
| subroutine quadmesh_phys2log | ( | class(oft_quadmesh), 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] | self | Mesh containing cell |
| [in] | cell | Index of cell for evaulation |
| [in] | pt | Physical position [3] |
| subroutine quadmesh_quad_rule | ( | class(oft_quadmesh), 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 quadmesh_save | ( | class(oft_quadmesh), intent(in) | self, |
| character(len=*), intent(in) | filename | ||
| ) |
Load trimesh from transfer file.
| subroutine quadmesh_set_order | ( | class(oft_quadmesh), intent(inout) | self, |
| integer(i4), intent(in) | order | ||
| ) |
Load trimesh from transfer file.
| subroutine quadmesh_setup | ( | class(oft_quadmesh), intent(inout) | self, |
| integer(i4), intent(in) | cad_type, | ||
| logical, intent(in) | has_parent | ||
| ) |
Load quadmesh from transfer file.
| subroutine quadmesh_tang | ( | class(oft_quadmesh), 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.
| [in] | self | Mesh containing face |
| [in] | i | Index of face for evaulation |
| [in] | f | Logical coordinate on face [3] |
| [out] | t | Orthonormal vectors [3,2] |
| subroutine quadmesh_tessellate | ( | class(oft_quadmesh), 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 quadmesh_vlog | ( | class(oft_quadmesh), 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] |
| integer(i4), dimension(2,4), parameter inodes1p = RESHAPE((/ 1,1, 2,1, 2,2, 1,2/), (/2,4/)) |
| integer(i4), dimension(2,4), parameter inodes2e = RESHAPE((/ 2,1, 3,2, 2,3, 1,2/), (/2,4/)) |
| integer(i4), dimension(2), parameter inodes2f = (/2,2/) |
| integer(i4), dimension(2,4), parameter inodes2p = RESHAPE((/ 1,1, 3,1, 3,3, 1,3/), (/2,4/)) |
| integer(i4), dimension(2,4), parameter inodesp_base = RESHAPE((/ 0,0, 1,0, 1,1, 0,1/), (/2,4/)) |
| integer(i4), dimension(4), parameter quad_bary_map = [-2,1,2,-1] |
| integer(i4), dimension(2,4), parameter quad_ed =RESHAPE((/1,2, 2,3, 3,4, 4,1/), (/2,4/)) |
Quad edge list.