OGS
ExtrapolatableElementCollection.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <functional>
14#include <vector>
15
18
19namespace NumLib
20{
21class LocalToGlobalIndexMap;
22
30{
31public:
34 virtual Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
35 std::size_t const id, unsigned const integration_point) const = 0;
36
63 virtual std::vector<double> const& getIntegrationPointValues(
64 std::size_t const id, const double t,
65 std::vector<GlobalVector*> const& x,
66 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
67 std::vector<double>& cache) const = 0;
68
70 virtual std::size_t size() const = 0;
71
73};
74
75template <typename LocalAssemblerCollection>
78{
79public:
82 typename std::decay_t<decltype(*std::declval<LocalAssemblerCollection>()
83 [static_cast<std::size_t>(0)])>;
84
85 static_assert(std::is_base_of<ExtrapolatableElement, LocalAssembler>::value,
86 "Local assemblers used for extrapolation must be derived "
87 "from ExtrapolatableElement.");
88
98 std::function<std::vector<double> const&(
99 LocalAssembler const& loc_asm, const double t,
100 std::vector<GlobalVector*> const& x,
101 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
102 std::vector<double>& cache)>;
103
112 LocalAssemblerCollection const& local_assemblers,
113 IntegrationPointValuesMethod const& integration_point_values_method)
114 : _local_assemblers(local_assemblers),
115 _integration_point_values_method{integration_point_values_method}
116 {
117 }
118
119 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
120 std::size_t const id, unsigned const integration_point) const override
121 {
122 ExtrapolatableElement const& loc_asm = *_local_assemblers[id];
123 return loc_asm.getShapeMatrix(integration_point);
124 }
125
126 std::vector<double> const& getIntegrationPointValues(
127 std::size_t const id, const double t,
128 std::vector<GlobalVector*> const& x,
129 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
130 std::vector<double>& cache) const override
131 {
132 auto const& loc_asm = *_local_assemblers[id];
133 return _integration_point_values_method(loc_asm, t, x, dof_table,
134 cache);
135 }
136
137 std::size_t size() const override { return _local_assemblers.size(); }
138
139private:
140 LocalAssemblerCollection const& _local_assemblers;
142};
143
146template <typename LocalAssemblerCollection,
147 typename IntegrationPointValuesMethod>
149makeExtrapolatable(LocalAssemblerCollection const& local_assemblers,
150 IntegrationPointValuesMethod integration_point_values_method)
151{
153 local_assemblers, integration_point_values_method};
154}
155} // namespace NumLib
virtual Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(std::size_t const id, unsigned const integration_point) const =0
virtual std::size_t size() const =0
Returns the number of elements whose properties shall be extrapolated.
virtual ~ExtrapolatableElementCollection()=default
virtual std::vector< double > const & getIntegrationPointValues(std::size_t const id, const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const =0
Provides the shape matrix at the given integration point.
std::function< std::vector< double > const &( LocalAssembler const &loc_asm, const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache)> IntegrationPointValuesMethod
std::vector< double > const & getIntegrationPointValues(std::size_t const id, const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
typename std::decay_t< decltype(*std::declval< LocalAssemblerCollection >()[static_cast< std::size_t >(0)])> LocalAssembler
LocalAssemblerCollection contains many LocalAssembler's.
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(std::size_t const id, unsigned const integration_point) const override
std::size_t size() const override
Returns the number of elements whose properties shall be extrapolated.
IntegrationPointValuesMethod const _integration_point_values_method
ExtrapolatableLocalAssemblerCollection(LocalAssemblerCollection const &local_assemblers, IntegrationPointValuesMethod const &integration_point_values_method)
ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection > makeExtrapolatable(LocalAssemblerCollection const &local_assemblers, IntegrationPointValuesMethod integration_point_values_method)