OGS
TESLocalAssembler.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <memory>
14#include <vector>
15
20#include "TESAssemblyParams.h"
22
23namespace ProcessLib
24{
25namespace TES
26{
29{
30public:
31 ~TESLocalAssemblerInterface() override = default;
32
33 virtual bool checkBounds(std::vector<double> const& local_x,
34 std::vector<double> const& local_x_prev_ts) = 0;
35
36 virtual std::vector<double> const& getIntPtSolidDensity(
37 const double t,
38 std::vector<GlobalVector*> const& x,
39 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
40 std::vector<double>& cache) const = 0;
41
42 virtual std::vector<double> const& getIntPtLoading(
43 const double t,
44 std::vector<GlobalVector*> const& x,
45 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
46 std::vector<double>& cache) const = 0;
47
48 virtual std::vector<double> const& getIntPtReactionDampingFactor(
49 const double t,
50 std::vector<GlobalVector*> const& x,
51 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
52 std::vector<double>& cache) const = 0;
53
54 virtual std::vector<double> const& getIntPtReactionRate(
55 const double t,
56 std::vector<GlobalVector*> const& x,
57 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
58 std::vector<double>& cache) const = 0;
59
60 virtual std::vector<double> const& getIntPtDarcyVelocity(
61 const double t,
62 std::vector<GlobalVector*> const& x,
63 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
64 std::vector<double>& cache) const = 0;
65};
66
67template <typename ShapeFunction_, int GlobalDim>
69{
70public:
71 using ShapeFunction = ShapeFunction_;
74
76 MeshLib::Element const& e,
77 std::size_t const local_matrix_size,
78 NumLib::GenericIntegrationMethod const& integration_method,
79 bool const is_axially_symmetric,
80 AssemblyParams const& asm_params);
81
82 void assemble(double const t, double const dt,
83 std::vector<double> const& local_x,
84 std::vector<double> const& local_x_prev,
85 std::vector<double>& local_M_data,
86 std::vector<double>& local_K_data,
87 std::vector<double>& local_b_data) override;
88
89 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
90 const unsigned integration_point) const override
91 {
92 auto const& N = _shape_matrices[integration_point].N;
93
94 // assumes N is stored contiguously in memory
95 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
96 }
97
98 bool checkBounds(std::vector<double> const& local_x,
99 std::vector<double> const& local_x_prev_ts) override;
100
101 std::vector<double> const& getIntPtSolidDensity(
102 const double t,
103 std::vector<GlobalVector*> const& x,
104 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
105 std::vector<double>& cache) const override;
106
107 std::vector<double> const& getIntPtLoading(
108 const double t,
109 std::vector<GlobalVector*> const& x,
110 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
111 std::vector<double>& cache) const override;
112
113 std::vector<double> const& getIntPtReactionDampingFactor(
114 const double t,
115 std::vector<GlobalVector*> const& x,
116 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
117 std::vector<double>& cache) const override;
118
119 std::vector<double> const& getIntPtReactionRate(
120 const double t,
121 std::vector<GlobalVector*> const& x,
122 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
123 std::vector<double>& cache) const override;
124
125 std::vector<double> const& getIntPtDarcyVelocity(
126 const double t,
127 std::vector<GlobalVector*> const& x,
128 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
129 std::vector<double>& cache) const override;
130
131private:
133
135
136 std::vector<ShapeMatrices, Eigen::aligned_allocator<ShapeMatrices>>
138
139 using LAT = LocalAssemblerTraits<ShapeMatricesType, ShapeFunction::NPOINTS,
140 NODAL_DOF, GlobalDim>;
141
143
146
147 static_assert(std::is_same_v<NodalMatrixType, typename LAT::LocalMatrix>,
148 "local matrix and data traits matrix do not coincide");
149 static_assert(std::is_same_v<NodalVectorType, typename LAT::LocalVector>,
150 "local vector and data traits vector do not coincide");
151};
152
153} // namespace TES
154} // namespace ProcessLib
155
Definition of the Element class.
virtual std::vector< double > const & getIntPtSolidDensity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtLoading(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtReactionDampingFactor(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual bool checkBounds(std::vector< double > const &local_x, std::vector< double > const &local_x_prev_ts)=0
virtual std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtReactionRate(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
NumLib::GenericIntegrationMethod const & _integration_method
std::vector< double > const & getIntPtReactionDampingFactor(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
std::vector< double > const & getIntPtSolidDensity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
typename LAT::LocalMatrix NodalMatrixType
typename LAT::LocalVector NodalVectorType
TESLocalAssembler(MeshLib::Element const &e, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, AssemblyParams const &asm_params)
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
TESLocalAssemblerInner< LAT > _d
std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::vector< double > const & getIntPtLoading(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::vector< double > const & getIntPtReactionRate(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
bool checkBounds(std::vector< double > const &local_x, std::vector< double > const &local_x_prev_ts) override
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
const unsigned NODAL_DOF
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix