OGS
SetOrGetIntegrationPointData.h
Go to the documentation of this file.
1
12#pragma once
13
14#include <Eigen/Core>
15#include <span>
16#include <vector>
17
19#include "TransposeInPlace.h"
20
21namespace ProcessLib
22{
23
24template <int DisplacementDim, typename IntegrationPointDataVector,
25 typename IpData, typename MemberType>
26std::vector<double> const& getIntegrationPointDimMatrixData(
27 IntegrationPointDataVector const& ip_data_vector,
28 MemberType IpData::*const member, std::vector<double>& cache)
29{
30 return getIntegrationPointDimMatrixData<DisplacementDim>(
31 ip_data_vector,
32 [member](IpData const& ip_data) -> auto const&
33 { return ip_data.*member; },
34 cache);
35}
36
37template <int DisplacementDim, typename IntegrationPointDataVector,
38 typename Accessor>
39std::vector<double> const& getIntegrationPointDimMatrixData(
40 IntegrationPointDataVector const& ip_data_vector, Accessor&& accessor,
41 std::vector<double>& cache)
42{
43 using AccessorResult = decltype(std::declval<Accessor>()(
44 std::declval<IntegrationPointDataVector>()[0]));
45 static_assert(std::is_lvalue_reference_v<AccessorResult>,
46 "The ip data accessor should return a reference. This check "
47 "prevents accidental copies.");
48
49 auto const n_integration_points = ip_data_vector.size();
50 constexpr int size = DisplacementDim * DisplacementDim;
51
52 cache.clear();
53 auto cache_mat = MathLib::createZeroedMatrix<
54 Eigen::Matrix<double, size, Eigen::Dynamic, Eigen::RowMajor>>(
55 cache, size, n_integration_points);
56
57 for (unsigned ip = 0; ip < n_integration_points; ++ip)
58 {
59 auto const& dim_matix = accessor(ip_data_vector[ip]);
60 cache_mat.col(ip) = dim_matix.template reshaped<Eigen::RowMajor>();
61 }
62
63 return cache;
64}
65
66template <int DisplacementDim, typename IntegrationPointDataVector,
67 typename IpData, typename MemberType>
68std::vector<double> const& getIntegrationPointKelvinVectorData(
69 IntegrationPointDataVector const& ip_data_vector,
70 MemberType IpData::*const member, std::vector<double>& cache)
71{
72 return getIntegrationPointKelvinVectorData<DisplacementDim>(
73 ip_data_vector,
74 [member](IpData const& ip_data) -> auto const&
75 { return ip_data.*member; },
76 cache);
77}
78
79template <int DisplacementDim, typename IntegrationPointDataVector,
80 typename Accessor>
81std::vector<double> const& getIntegrationPointKelvinVectorData(
82 IntegrationPointDataVector const& ip_data_vector, Accessor&& accessor,
83 std::vector<double>& cache)
84{
85 using AccessorResult = decltype(std::declval<Accessor>()(
86 std::declval<IntegrationPointDataVector>()[0]));
87 static_assert(std::is_lvalue_reference_v<AccessorResult>,
88 "The ip data accessor should return a reference. This check "
89 "prevents accidental copies.");
90
91 constexpr int kelvin_vector_size =
93 auto const n_integration_points = ip_data_vector.size();
94
95 cache.clear();
96 auto cache_mat = MathLib::createZeroedMatrix<Eigen::Matrix<
97 double, kelvin_vector_size, Eigen::Dynamic, Eigen::RowMajor>>(
98 cache, kelvin_vector_size, n_integration_points);
99
100 for (unsigned ip = 0; ip < n_integration_points; ++ip)
101 {
102 auto const& kelvin_vector = accessor(ip_data_vector[ip]);
103 cache_mat.col(ip) =
105 }
106
107 return cache;
108}
109
114template <int DisplacementDim, typename IntegrationPointDataVector,
115 typename MemberType>
117 IntegrationPointDataVector const& ip_data_vector, MemberType member)
118{
119 constexpr int kelvin_vector_size =
121
122 return transposeInPlace<kelvin_vector_size>(
123 [&](std::vector<double>& values)
124 {
125 return getIntegrationPointKelvinVectorData<DisplacementDim>(
126 ip_data_vector, member, values);
127 ;
128 });
129}
130
131template <int DisplacementDim, typename IntegrationPointDataVector,
132 typename IpData, typename MemberType>
134 double const* values,
135 IntegrationPointDataVector& ip_data_vector,
136 MemberType IpData::*const member)
137{
138 return setIntegrationPointKelvinVectorData<DisplacementDim>(
139 values, ip_data_vector,
140 [member](IpData& ip_data) -> auto& { return ip_data.*member; });
141}
142
143template <int DisplacementDim, typename IntegrationPointDataVector,
144 typename Accessor>
146 double const* values,
147 IntegrationPointDataVector& ip_data_vector,
148 Accessor&& accessor)
149{
150 using AccessorResult = decltype(std::declval<Accessor>()(
151 std::declval<IntegrationPointDataVector>()[0]));
152 static_assert(std::is_lvalue_reference_v<AccessorResult>,
153 "The ip data accessor must return a reference.");
154
155 constexpr int kelvin_vector_size =
157 auto const n_integration_points = ip_data_vector.size();
158
159 auto kelvin_vector_values =
160 Eigen::Map<Eigen::Matrix<double, kelvin_vector_size, Eigen::Dynamic,
161 Eigen::ColMajor> const>(
162 values, kelvin_vector_size, n_integration_points);
163
164 for (unsigned ip = 0; ip < n_integration_points; ++ip)
165 {
166 accessor(ip_data_vector[ip]) =
168 kelvin_vector_values.col(ip));
169 }
170
171 return n_integration_points;
172}
173
174template <typename IntegrationPointDataVector, typename IpData,
175 typename MemberType>
176std::vector<double> const& getIntegrationPointScalarData(
177 IntegrationPointDataVector const& ip_data_vector,
178 MemberType IpData::*const member, std::vector<double>& cache)
179{
181 ip_data_vector,
182 [member](IpData const& ip_data) -> auto const&
183 { return ip_data.*member; },
184 cache);
185}
186
187template <typename IntegrationPointDataVector, typename Accessor>
188std::vector<double> const& getIntegrationPointScalarData(
189 IntegrationPointDataVector const& ip_data_vector, Accessor&& accessor,
190 std::vector<double>& cache)
191{
192 using AccessorResult = decltype(std::declval<Accessor>()(
193 std::declval<IntegrationPointDataVector>()[0]));
194 static_assert(std::is_lvalue_reference_v<AccessorResult>,
195 "The ip data accessor should return a reference. This check "
196 "prevents accidental copies.");
197
198 auto const n_integration_points = ip_data_vector.size();
199
200 cache.clear();
201 auto cache_mat = MathLib::createZeroedMatrix<
202 Eigen::Matrix<double, 1, Eigen::Dynamic, Eigen::RowMajor>>(
203 cache, 1, n_integration_points);
204
205 for (unsigned ip = 0; ip < n_integration_points; ++ip)
206 {
207 cache_mat[ip] = accessor(ip_data_vector[ip]);
208 }
209
210 return cache;
211}
212
213template <typename IntegrationPointDataVector, typename IpData,
214 typename MemberType>
216 double const* values,
217 IntegrationPointDataVector& ip_data_vector,
218 MemberType IpData::*const member)
219{
220 return setIntegrationPointScalarData(values, ip_data_vector,
221 [member](IpData& ip_data) -> auto&
222 { return ip_data.*member; });
223}
224
225template <typename IntegrationPointDataVector, typename Accessor>
227 double const* values,
228 IntegrationPointDataVector& ip_data_vector,
229 Accessor&& accessor)
230{
231 using AccessorResult = decltype(std::declval<Accessor>()(
232 std::declval<IntegrationPointDataVector>()[0]));
233 static_assert(std::is_lvalue_reference_v<AccessorResult>,
234 "The ip data accessor must return a reference.");
235
236 auto const n_integration_points = ip_data_vector.size();
237
238 for (unsigned ip = 0; ip < n_integration_points; ++ip)
239 {
240 accessor(ip_data_vector[ip]) = values[ip];
241 }
242 return n_integration_points;
243}
244
245template <typename IntegrationPointDataVector, typename MemberType,
246 typename MaterialStateVariables>
248 IntegrationPointDataVector const& ip_data_vector,
249 MemberType member,
250 std::function<std::span<double>(MaterialStateVariables&)>
251 get_values_span,
252 int const n_components)
253{
254 std::vector<double> result;
255 result.reserve(ip_data_vector.size() * n_components);
256
257 for (auto& ip_data : ip_data_vector)
258 {
259 auto const values_span = get_values_span(*(ip_data.*member));
260 assert(values_span.size() == static_cast<std::size_t>(n_components));
261
262 result.insert(end(result), values_span.begin(), values_span.end());
263 }
264
265 return result;
266}
267
268template <typename IntegrationPointDataVector, typename MemberType,
269 typename MaterialStateVariables>
271 double const* values,
272 IntegrationPointDataVector& ip_data_vector,
273 MemberType member,
274 std::function<std::span<double>(MaterialStateVariables&)>
275 get_values_span)
276{
277 auto const n_integration_points = ip_data_vector.size();
278
279 std::size_t position = 0;
280 for (auto const& ip_data : ip_data_vector)
281 {
282 auto const values_span = get_values_span(*(ip_data.*member));
283 std::copy_n(values + position, values_span.size(), values_span.begin());
284 position += values_span.size();
285 }
286 return n_integration_points;
287}
288} // namespace ProcessLib
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
std::vector< double > const & getIntegrationPointScalarData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
std::size_t setIntegrationPointDataMaterialStateVariables(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType member, std::function< std::span< double >(MaterialStateVariables &)> get_values_span)
std::vector< double > getIntegrationPointDataMaterialStateVariables(IntegrationPointDataVector const &ip_data_vector, MemberType member, std::function< std::span< double >(MaterialStateVariables &)> get_values_span, int const n_components)
std::vector< double > const & getIntegrationPointDimMatrixData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
std::vector< double > const & getIntegrationPointKelvinVectorData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
std::size_t setIntegrationPointKelvinVectorData(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType IpData::*const member)
std::size_t setIntegrationPointScalarData(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType IpData::*const member)