OGS
SetOrGetIntegrationPointData.h
Go to the documentation of this file.
1 
12 #pragma once
13 
14 #include <Eigen/Eigen>
15 #include <vector>
16 
17 #include "MathLib/KelvinVector.h"
20 
21 namespace ProcessLib
22 {
23 template <int DisplacementDim, typename IntegrationPointData,
24  typename MemberType>
25 std::vector<double> const& getIntegrationPointKelvinVectorData(
26  std::vector<IntegrationPointData,
27  Eigen::aligned_allocator<IntegrationPointData>> const& ip_data,
28  MemberType member, std::vector<double>& cache)
29 {
30  constexpr int kelvin_vector_size =
32  auto const n_integration_points = ip_data.size();
33 
34  cache.clear();
35  auto cache_mat = MathLib::createZeroedMatrix<Eigen::Matrix<
36  double, kelvin_vector_size, Eigen::Dynamic, Eigen::RowMajor>>(
37  cache, kelvin_vector_size, n_integration_points);
38 
39  for (unsigned ip = 0; ip < n_integration_points; ++ip)
40  {
41  auto const& kelvin_vector = ip_data[ip].*member;
42  cache_mat.col(ip) =
44  }
45 
46  return cache;
47 }
48 
49 template <int DisplacementDim, typename IntegrationPointData,
50  typename MemberType>
52  std::vector<IntegrationPointData,
53  Eigen::aligned_allocator<IntegrationPointData>> const& ip_data,
54  MemberType member)
55 {
56  constexpr int kelvin_vector_size =
58  auto const n_integration_points = ip_data.size();
59 
60  std::vector<double> ip_kelvin_vector_values;
61  auto cache_mat = MathLib::createZeroedMatrix<Eigen::Matrix<
62  double, Eigen::Dynamic, kelvin_vector_size, Eigen::RowMajor>>(
63  ip_kelvin_vector_values, n_integration_points, kelvin_vector_size);
64 
65  for (unsigned ip = 0; ip < n_integration_points; ++ip)
66  {
67  auto const& ip_member = ip_data[ip].*member;
68  cache_mat.row(ip) =
70  }
71 
72  return ip_kelvin_vector_values;
73 }
74 
75 template <int DisplacementDim, typename IntegrationPointData,
76  typename MemberType>
78  double const* values,
79  std::vector<IntegrationPointData,
80  Eigen::aligned_allocator<IntegrationPointData>>& ip_data,
81  MemberType member)
82 {
83  constexpr int kelvin_vector_size =
85  auto const n_integration_points = ip_data.size();
86 
87  auto kelvin_vector_values =
88  Eigen::Map<Eigen::Matrix<double, kelvin_vector_size, Eigen::Dynamic,
89  Eigen::ColMajor> const>(
90  values, kelvin_vector_size, n_integration_points);
91 
92  for (unsigned ip = 0; ip < n_integration_points; ++ip)
93  {
94  ip_data[ip].*member =
96  kelvin_vector_values.col(ip));
97  }
98 
99  return n_integration_points;
100 }
101 
102 template <typename IntegrationPointData, typename MemberType>
103 std::vector<double> const& getIntegrationPointScalarData(
104  std::vector<IntegrationPointData,
105  Eigen::aligned_allocator<IntegrationPointData>> const& ip_data,
106  MemberType member, std::vector<double>& cache)
107 {
108  auto const n_integration_points = ip_data.size();
109 
110  cache.clear();
111  auto cache_mat = MathLib::createZeroedMatrix<
112  Eigen::Matrix<double, 1, Eigen::Dynamic, Eigen::RowMajor>>(
113  cache, 1, n_integration_points);
114 
115  for (unsigned ip = 0; ip < n_integration_points; ++ip)
116  {
117  cache_mat[ip] = ip_data[ip].*member;
118  }
119 
120  return cache;
121 }
122 
123 template <typename IntegrationPointData, typename MemberType>
125  double const* values,
126  std::vector<IntegrationPointData,
127  Eigen::aligned_allocator<IntegrationPointData>>& ip_data,
128  MemberType member)
129 {
130  auto const n_integration_points = ip_data.size();
131 
132  for (unsigned ip = 0; ip < n_integration_points; ++ip)
133  {
134  ip_data[ip].*member = values[ip];
135  }
136  return n_integration_points;
137 }
138 
139 template <typename IntegrationPointData, typename MemberType,
140  typename MaterialStateVariables>
142  std::vector<IntegrationPointData,
143  Eigen::aligned_allocator<IntegrationPointData>> const&
144  ip_data_vector,
145  MemberType member,
146  std::function<BaseLib::DynamicSpan<double>(MaterialStateVariables&)>
147  get_values_span,
148  int const n_components)
149 {
150  std::vector<double> result;
151  result.reserve(ip_data_vector.size() * n_components);
152 
153  for (auto& ip_data : ip_data_vector)
154  {
155  auto const values_span = get_values_span(*(ip_data.*member));
156  assert(values_span.size == static_cast<std::size_t>(n_components));
157 
158  result.insert(end(result), values_span.data,
159  values_span.data + values_span.size);
160  }
161 
162  return result;
163 }
164 
165 template <typename IntegrationPointData, typename MemberType,
166  typename MaterialStateVariables>
168  double const* values,
169  std::vector<IntegrationPointData,
170  Eigen::aligned_allocator<IntegrationPointData>>& ip_data_vector,
171  MemberType member,
172  std::function<BaseLib::DynamicSpan<double>(MaterialStateVariables&)>
173  get_values_span)
174 {
175  auto const n_integration_points = ip_data_vector.size();
176 
177  std::size_t position = 0;
178  for (auto& ip_data : ip_data_vector)
179  {
180  auto const values_span = get_values_span(*(ip_data.*member));
181  std::copy_n(values + position, values_span.size, values_span.data);
182  position += values_span.size;
183  }
184  return n_integration_points;
185 }
186 } // 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.
Definition: KelvinVector.h:23
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
Definition: KelvinVector.h:170
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:32
std::vector< double > const & getIntegrationPointKelvinVectorData(std::vector< IntegrationPointData, Eigen::aligned_allocator< IntegrationPointData >> const &ip_data, MemberType member, std::vector< double > &cache)
std::vector< double > const & getIntegrationPointScalarData(std::vector< IntegrationPointData, Eigen::aligned_allocator< IntegrationPointData >> const &ip_data, MemberType member, std::vector< double > &cache)
std::vector< double > getIntegrationPointDataMaterialStateVariables(std::vector< IntegrationPointData, Eigen::aligned_allocator< IntegrationPointData >> const &ip_data_vector, MemberType member, std::function< BaseLib::DynamicSpan< double >(MaterialStateVariables &)> get_values_span, int const n_components)
std::size_t setIntegrationPointKelvinVectorData(double const *values, std::vector< IntegrationPointData, Eigen::aligned_allocator< IntegrationPointData >> &ip_data, MemberType member)
std::size_t setIntegrationPointScalarData(double const *values, std::vector< IntegrationPointData, Eigen::aligned_allocator< IntegrationPointData >> &ip_data, MemberType member)
std::size_t setIntegrationPointDataMaterialStateVariables(double const *values, std::vector< IntegrationPointData, Eigen::aligned_allocator< IntegrationPointData >> &ip_data_vector, MemberType member, std::function< BaseLib::DynamicSpan< double >(MaterialStateVariables &)> get_values_span)