OGS
ExtrapolatableElementCollection.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <functional>
7#include <vector>
8
11
12namespace NumLib
13{
15
23{
24public:
27 virtual Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
28 std::size_t const id, unsigned const integration_point) const = 0;
29
56 virtual std::vector<double> const& getIntegrationPointValues(
57 std::size_t const id, const double t,
58 std::vector<GlobalVector*> const& x,
59 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
60 std::vector<double>& cache) const = 0;
61
63 virtual std::size_t size() const = 0;
64
66};
67
68template <typename LocalAssemblerCollection>
71{
72public:
75 typename std::decay_t<decltype(*std::declval<LocalAssemblerCollection>()
76 [static_cast<std::size_t>(0)])>;
77
78 static_assert(std::is_base_of<ExtrapolatableElement, LocalAssembler>::value,
79 "Local assemblers used for extrapolation must be derived "
80 "from ExtrapolatableElement.");
81
91 std::function<std::vector<double> const&(
92 LocalAssembler const& loc_asm, const double t,
93 std::vector<GlobalVector*> const& x,
94 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
95 std::vector<double>& cache)>;
96
105 LocalAssemblerCollection const& local_assemblers,
106 IntegrationPointValuesMethod const& integration_point_values_method)
107 : _local_assemblers(local_assemblers),
108 _integration_point_values_method{integration_point_values_method}
109 {
110 }
111
112 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
113 std::size_t const id, unsigned const integration_point) const override
114 {
115 ExtrapolatableElement const& loc_asm = *_local_assemblers[id];
116 return loc_asm.getShapeMatrix(integration_point);
117 }
118
119 std::vector<double> const& getIntegrationPointValues(
120 std::size_t const id, 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 auto const& loc_asm = *_local_assemblers[id];
126 return _integration_point_values_method(loc_asm, t, x, dof_table,
127 cache);
128 }
129
130 std::size_t size() const override { return _local_assemblers.size(); }
131
132private:
133 LocalAssemblerCollection const& _local_assemblers;
135};
136
139template <typename LocalAssemblerCollection,
140 typename IntegrationPointValuesMethod>
142makeExtrapolatable(LocalAssemblerCollection const& local_assemblers,
143 IntegrationPointValuesMethod integration_point_values_method)
144{
146 local_assemblers, integration_point_values_method};
147}
148} // 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)