OGS
TemperatureDependentDiffusion.cpp
Go to the documentation of this file.
1
11
12#include <algorithm>
13#include <cmath>
14#include <iterator>
15
17
18namespace MaterialPropertyLib
19{
21{
22 if (!std::holds_alternative<Component*>(scale_))
23 {
25 "The property 'TemperatureDependentDiffusion' is "
26 "implemented on the 'component' scale only.");
27 }
28}
29
31 VariableArray const& variable_array,
32 ParameterLib::SpatialPosition const& pos, double const t,
33 double const /*dt*/) const
34{
35 auto const T = variable_array.temperature;
36 double const gas_constant = MaterialLib::PhysicalConstant::IdealGasConstant;
37 double const Arrhenius_exponent =
38 std::exp(Ea_ / gas_constant * (1 / T0_ - 1 / T));
39
40 auto const D0_data = D0_(t, pos);
41 std::vector<double> D;
42 std::transform(D0_data.cbegin(), D0_data.cend(), std::back_inserter(D),
43 [&Arrhenius_exponent](double const D0_component)
44 { return D0_component * Arrhenius_exponent; });
45
46 return fromVector(D);
47}
48} // namespace MaterialPropertyLib
#define OGS_FATAL(...)
Definition Error.h:26
virtual PropertyDataType value() const
Definition Property.cpp:76
std::variant< Medium *, Phase *, Component * > scale_
Definition Property.h:297
ParameterLib::Parameter< double > const & D0_
the molecular diffusion at the reference temperature
double const Ea_
the activition energy for diffusion
PropertyDataType fromVector(std::vector< double > const &values)
Definition Property.cpp:23
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
Definition Property.h:31