 |
waves
Basic FE playground
|
Go to the documentation of this file. 1 #ifndef KATOPTRON_ELEMENTMATRICES_H
2 #define KATOPTRON_ELEMENTMATRICES_H
4 #include "Stokhos_Sacado_Kokkos_MP_Vector.hpp"
13 #include <Kokkos_Core.hpp>
14 #include "Kokkos_ViewFactory.hpp"
15 #include "wCacheHex8.h"
16 #include "wCacheQuad4.h"
17 #include "wCacheTetra4.h"
18 #include "wCacheTri3.h"
45 template <
typename scalar,
int element_type,
int element_size>
49 typedef Kokkos::View<scalar ***, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace>
view_type_3D;
50 typedef Kokkos::View<scalar **, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace>
view_type_2D;
56 Teuchos::RCP<katoptron::Domain<scalar>>
domain;
63 void compute(
int e,
int i_thread);
69 #endif //KATOPTRON_ELEMENTMATRICES_H
Class which is used to store all the information related to the discretized domain:
Definition: Domain.h:26
view_type_2D material
Definition: ElementMatrices.h:54
view_type_3D S
Definition: ElementMatrices.h:52
Class used to compute the element matrices:
Definition: ElementMatrices.h:46
~ElementMatrices()
Definition: ElementMatrices.h:61
view_type_3D L
Definition: ElementMatrices.h:52
size_t numPrimalDPN
Definition: ElementMatrices.h:58
view_type_2D getS(size_t i_thread) const
Return the coupling element matrix computed by the thread i_thread.
Definition: ElementMatrices.hpp:296
view_type_2D getK(size_t i_thread) const
Return the element stifness matrix computed by the thread i_thread.
Definition: ElementMatrices.hpp:264
ElementMatrices(Teuchos::RCP< katoptron::Domain< scalar >> domain, size_t numPrimalDPN, size_t pool_size)
ElementMatrices constructor.
Definition: ElementMatrices.hpp:17
Teuchos::RCP< katoptron::Domain< scalar > > domain
Definition: ElementMatrices.h:56
Kokkos::View< scalar ***, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace > view_type_3D
Definition: ElementMatrices.h:49
view_type_2D getL(size_t i_thread) const
Return the thermal element matrix computed by the thread i_thread.
Definition: ElementMatrices.hpp:280
view_type_3D K
Definition: ElementMatrices.h:52
Kokkos::View< scalar **, Kokkos::LayoutRight, Kokkos::DefaultExecutionSpace > view_type_2D
Definition: ElementMatrices.h:50
void compute(int e, int i_thread)
Compute all the matrices and store them as member data. This function must be called before getK,...
Definition: ElementMatrices.hpp:43
Base class for the element computations.
Definition: ElementComputation.h:30