OGS
MaterialLib/MPL/Property.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 <Eigen/Core>
7#include <array>
8#include <string>
9#include <typeinfo>
10#include <variant>
11
13#include "BaseLib/Error.h"
15#include "PropertyType.h"
16#include "VariableType.h"
17
18namespace MaterialPropertyLib
19{
20class Medium;
21class Phase;
22class Component;
23
25 std::variant<double, Eigen::Matrix<double, 2, 1>,
26 Eigen::Matrix<double, 3, 1>, Eigen::Matrix<double, 2, 2>,
27 Eigen::Matrix<double, 3, 3>, Eigen::Matrix<double, 4, 1>,
28 Eigen::Matrix<double, 6, 1>, Eigen::MatrixXd>;
29
38PropertyDataType fromVector(std::vector<double> const& values);
39
44{
45public:
46#ifndef NDEBUG
47 virtual ~Property()
48 {
49 if (property_used)
50 {
51 DBUG("Property is used: '{:s}'", description());
52 }
53 else
54 {
55 WARN("Property is not used: '{:s}'", description());
56 }
57 }
58#else
59 virtual ~Property() = default;
60#endif
61
65 ParameterLib::SpatialPosition const& pos, double const t) const;
66
69 virtual PropertyDataType value() const;
73 virtual PropertyDataType value(VariableArray const& variable_array,
74 VariableArray const& variable_array_prev,
76 double const t, double const dt) const;
80 virtual PropertyDataType value(VariableArray const& variable_array,
82 double const t, double const dt) const;
86 virtual PropertyDataType dValue(VariableArray const& variable_array,
87 VariableArray const& variable_array_prev,
88 Variable const variable,
90 double const t, double const dt) const;
94 virtual PropertyDataType dValue(VariableArray const& variable_array,
95 Variable const variable,
97 double const t, double const dt) const;
100 virtual PropertyDataType d2Value(VariableArray const& variable_array,
101 Variable const variable1,
102 Variable const variable2,
104 double const t, double const dt) const;
105
108 virtual void setProperties(
109 std::vector<std::unique_ptr<Phase>> const& phases);
110
111 void setScale(std::variant<Medium*, Phase*, Component*> scale)
112 {
113 scale_ = scale;
114 checkScale();
115 };
116
117 template <typename T>
119 double const t) const
120 {
121 try
122 {
123 return std::get<T>(initialValue(pos, t));
124 }
125 catch (std::bad_variant_access const&)
126 {
127 OGS_FATAL(
128 "The initial value of {:s} does not hold requested type '{:s}' "
129 "but a {:s}.",
130 description(),
132 property_data_type_names_[initialValue(pos, t).index()]);
133 }
134 }
135
136 template <typename T>
137 T value() const
138 {
139 try
140 {
141#ifndef NDEBUG
142 property_used = true;
143#endif
144 return std::get<T>(value());
145 }
146 catch (std::bad_variant_access const&)
147 {
148 OGS_FATAL(
149 "The value of {:s} does not hold requested type '{:s}' but a "
150 "{:s}.",
151 description(),
154 }
155 }
156
157 template <typename T>
158 T value(VariableArray const& variable_array,
159 VariableArray const& variable_array_prev,
160 ParameterLib::SpatialPosition const& pos, double const t,
161 double const dt) const
162 {
163 try
164 {
165#ifndef NDEBUG
166 property_used = true;
167#endif
168 return std::get<T>(
169 value(variable_array, variable_array_prev, pos, t, dt));
170 }
171 catch (std::bad_variant_access const&)
172 {
173 OGS_FATAL(
174 "The value of {:s} is not of the requested type '{:s}' but a "
175 "{:s}.",
176 description(),
178 property_data_type_names_[value(variable_array,
179 variable_array_prev, pos, t, dt)
180 .index()]);
181 }
182 }
183 template <typename T>
184 T value(VariableArray const& variable_array,
185 ParameterLib::SpatialPosition const& pos, double const t,
186 double const dt) const
187 {
188 try
189 {
190#ifndef NDEBUG
191 property_used = true;
192#endif
193 return std::get<T>(value(variable_array, pos, t, dt));
194 }
195 catch (std::bad_variant_access const&)
196 {
197 OGS_FATAL(
198 "The value of {:s} is not of the requested type '{:s}' but a "
199 "{:s}.",
200 description(),
202 property_data_type_names_[value(variable_array, pos, t, dt)
203 .index()]);
204 }
205 }
206
207 template <typename T>
208 T dValue(VariableArray const& variable_array,
209 VariableArray const& variable_array_prev, Variable const variable,
210 ParameterLib::SpatialPosition const& pos, double const t,
211 double const dt) const
212 {
213 try
214 {
215#ifndef NDEBUG
216 property_used = true;
217#endif
218 return std::get<T>(dValue(variable_array, variable_array_prev,
219 variable, pos, t, dt));
220 }
221 catch (std::bad_variant_access const&)
222 {
223 OGS_FATAL(
224 "The first derivative value of {:s} is not of the requested "
225 "type '{:s}' but a {:s}.",
226 description(),
229 [dValue(variable_array, variable, pos, t, dt).index()]);
230 }
231 }
232 template <typename T>
233 T dValue(VariableArray const& variable_array, Variable const variable,
234 ParameterLib::SpatialPosition const& pos, double const t,
235 double const dt) const
236 {
237 try
238 {
239#ifndef NDEBUG
240 property_used = true;
241#endif
242 return std::get<T>(dValue(variable_array, variable, pos, t, dt));
243 }
244 catch (std::bad_variant_access const&)
245 {
246 OGS_FATAL(
247 "The first derivative value of {:s} is not of the requested "
248 "type '{:s}' but a {:s}.",
249 description(),
252 [dValue(variable_array, variable, pos, t, dt).index()]);
253 }
254 }
255 template <typename T>
256 T d2Value(VariableArray const& variable_array, Variable const& variable1,
257 Variable const& variable2,
258 ParameterLib::SpatialPosition const& pos, double const t,
259 double const dt) const
260 {
261 try
262 {
263#ifndef NDEBUG
264 property_used = true;
265#endif
266 return std::get<T>(
267 d2Value(variable_array, variable1, variable2, pos, t, dt));
268 }
269 catch (std::bad_variant_access const&)
270 {
271 OGS_FATAL(
272 "The second derivative value of {:s} is not of the requested "
273 "type '{:s}' but a {:s}.",
274 description(),
276 property_data_type_names_[d2Value(variable_array, variable1,
277 variable2, pos, t, dt)
278 .index()]);
279 }
280 }
281
282protected:
283 std::string name_;
290 std::variant<Medium*, Phase*, Component*> scale_;
291
292private:
293 virtual void checkScale() const
294 {
295 // Empty check for properties which can be defined on every scale,
296 // medium, phase or component
297 }
298 std::string description() const;
299#ifndef NDEBUG
300 mutable bool property_used = false;
301#endif
302
303private:
305 static constexpr std::array property_data_type_names_ = {
306 "scalar", "2-vector", "3-vector",
307 "2x2-matrix", "3x3-matrix", "2D-Kelvin vector",
308 "3D-Kelvin vector", "dynamic matrix type"};
309 static_assert(property_data_type_names_.size() ==
310 std::variant_size_v<PropertyDataType>,
311 "The array of property data type names has different size "
312 "than the PropertyDataType variant type.");
313};
314
316 PropertyArray& properties,
317 PropertyArray& new_properties,
318 std::variant<Medium*, Phase*, Component*>
319 scale_pointer)
320{
321 for (std::size_t i = 0; i < properties.size(); ++i)
322 {
323 if (new_properties[i] != nullptr)
324 {
325 properties[i] = std::move(new_properties[i]);
326 properties[i]->setScale(scale_pointer);
327 }
328 }
329}
330
332 PropertyArray& properties,
333 std::vector<std::unique_ptr<Phase>> const& phases)
334{
335 for (auto& p : properties)
336 {
337 if (p != nullptr)
338 {
339 p->setProperties(phases);
340 }
341 }
342}
343
344} // namespace MaterialPropertyLib
#define OGS_FATAL(...)
Definition Error.h:19
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
This class defines components (substances).
Definition Component.h:18
virtual void setProperties(std::vector< std::unique_ptr< Phase > > const &phases)
Default implementation:
virtual PropertyDataType d2Value(VariableArray const &variable_array, Variable const variable1, Variable const variable2, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
Default implementation: 2nd derivative of any constant property is zero.
static constexpr std::array property_data_type_names_
Corresponds to the PropertyDataType.
T dValue(VariableArray const &variable_array, VariableArray const &variable_array_prev, Variable const variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
T value(VariableArray const &variable_array, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
T d2Value(VariableArray const &variable_array, Variable const &variable1, Variable const &variable2, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
PropertyDataType value_
The single value of a property.
T dValue(VariableArray const &variable_array, Variable const variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
T initialValue(ParameterLib::SpatialPosition const &pos, double const t) const
virtual PropertyDataType value() const
void setScale(std::variant< Medium *, Phase *, Component * > scale)
T value(VariableArray const &variable_array, VariableArray const &variable_array_prev, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
std::variant< Medium *, Phase *, Component * > scale_
virtual PropertyDataType initialValue(ParameterLib::SpatialPosition const &pos, double const t) const
virtual PropertyDataType dValue(VariableArray const &variable_array, VariableArray const &variable_array_prev, Variable const variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const
std::string typeToString()
void overwriteExistingProperties(PropertyArray &properties, PropertyArray &new_properties, std::variant< Medium *, Phase *, Component * > scale_pointer)
PropertyDataType fromVector(std::vector< double > const &values)
std::array< std::unique_ptr< Property >, PropertyType::number_of_properties > PropertyArray
void updatePropertiesForAllPhases(PropertyArray &properties, std::vector< std::unique_ptr< Phase > > const &phases)
std::variant< double, Eigen::Matrix< double, 2, 1 >, Eigen::Matrix< double, 3, 1 >, Eigen::Matrix< double, 2, 2 >, Eigen::Matrix< double, 3, 3 >, Eigen::Matrix< double, 4, 1 >, Eigen::Matrix< double, 6, 1 >, Eigen::MatrixXd > PropertyDataType