Class used to compute the element matrices:
More...
#include <ElementMatrices.h>
|
typedef Kokkos::View< scalar ***, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace > | view_type_3D |
|
typedef Kokkos::View< scalar **, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace > | view_type_2D |
|
template<typename scalar, int element_type, int element_size>
class ElementMatrices< scalar, element_type, element_size >
Class used to compute the element matrices:
-
The element stifness matrix \(\boldsymbol{K}\),
-
The thermal element matrix \(\boldsymbol{L}\),
-
The coupling element matrix \(\boldsymbol{S}\).
◆ view_type_2D
template<typename scalar , int element_type, int element_size>
typedef Kokkos::View<scalar **, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace> ElementMatrices< scalar, element_type, element_size >::view_type_2D |
◆ view_type_3D
template<typename scalar , int element_type, int element_size>
typedef Kokkos::View<scalar ***, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace> ElementMatrices< scalar, element_type, element_size >::view_type_3D |
◆ ElementMatrices()
template<typename scalar , int element_type, int element_size>
ElementMatrices constructor.
Arguments:
-
domain: an RCP to a Domain object,
-
numPrimalDPN: the number of degrees of freedom per node (without taking into account the Lagrange multipliers),
-
pool_size: the number of threads that will be used to compute the element matrices.
The number of used threads is specify to allocate distinct memory for each thread to ensure that the threads do not use the same memory addresses.
◆ ~ElementMatrices()
template<typename scalar , int element_type, int element_size>
◆ compute()
template<typename scalar , int element_type, int element_size>
void ElementMatrices< scalar, element_type, element_size >::compute |
( |
int |
e, |
|
|
int |
i_thread |
|
) |
| |
Compute all the matrices and store them as member data. This function must be called before getK, getL, or getS.
The computations have been gathered into one function to reuse reusable data such as Jacobian computations.
Arguments:
-
e: the local ID of the element for which the matrices have to be computed,
-
i_thread: the ID of the thread used to compute the element matrices.
◆ getK()
template<typename scalar , int element_type, int element_size>
Kokkos::View< scalar **, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace > ElementMatrices< scalar, element_type, element_size >::getK |
( |
size_t |
i_thread | ) |
const |
Return the element stifness matrix \(\boldsymbol{K}\) computed by the thread i_thread.
Argument:
-
i_thread: the ID of the thread used to compute the element matrices.
◆ getL()
template<typename scalar , int element_type, int element_size>
Kokkos::View< scalar **, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace > ElementMatrices< scalar, element_type, element_size >::getL |
( |
size_t |
i_thread | ) |
const |
Return the thermal element matrix \(\boldsymbol{L}\) computed by the thread i_thread.
Argument:
-
i_thread: the ID of the thread used to compute the element matrices.
◆ getS()
template<typename scalar , int element_type, int element_size>
Kokkos::View< scalar **, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace > ElementMatrices< scalar, element_type, element_size >::getS |
( |
size_t |
i_thread | ) |
const |
Return the coupling element matrix \(\boldsymbol{S}\) computed by the thread i_thread.
Argument:
-
i_thread: the ID of the thread used to compute the element matrices.
◆ domain
template<typename scalar , int element_type, int element_size>
template<typename scalar , int element_type, int element_size>
template<typename scalar , int element_type, int element_size>
◆ material
template<typename scalar , int element_type, int element_size>
◆ numPrimalDPN
template<typename scalar , int element_type, int element_size>
template<typename scalar , int element_type, int element_size>
The documentation for this class was generated from the following files: