OGS
ThermoRichardsMechanics/LocalAssemblerInterface.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
16
18{
19template <int DisplacementDim, typename ConstitutiveTraits>
22{
24 MeshLib::Element const& e,
25 NumLib::GenericIntegrationMethod const& integration_method,
26 bool const is_axially_symmetric,
28 process_data)
29 : process_data_(process_data),
30 integration_method_(integration_method),
31 element_(e),
32 is_axially_symmetric_(is_axially_symmetric),
33 solid_material_(MaterialLib::Solids::selectSolidConstitutiveRelation(
34 process_data_.solid_materials, process_data_.material_ids,
35 e.getID()))
36 {
37 unsigned const n_integration_points =
38 integration_method_.getNumberOfPoints();
39
40 current_states_.resize(n_integration_points);
41 prev_states_.resize(n_integration_points);
42 output_data_.resize(n_integration_points);
43
44 material_states_.reserve(n_integration_points);
45 for (unsigned ip = 0; ip < n_integration_points; ++ip)
46 {
47 material_states_.emplace_back(
48 solid_material_.createMaterialStateVariables());
49
50 // Set initial strain field to zero.
51 std::get<StrainData<DisplacementDim>>(current_states_[ip]).eps =
53 }
54 }
55
56 std::size_t setIPDataInitialConditions(std::string_view name,
57 double const* values,
58 int const integration_order)
59 {
60 if (integration_order !=
61 static_cast<int>(integration_method_.getIntegrationOrder()))
62 {
64 "Setting integration point initial conditions; The integration "
65 "order of the local assembler for element {:d} is different "
66 "from the integration order in the initial condition.",
67 element_.getID());
68 }
69
70 if (name == "sigma" && process_data_.initial_stress.value)
71 {
73 "Setting initial conditions for stress from integration "
74 "point data and from a parameter '{:s}' is not possible "
75 "simultaneously.",
76 process_data_.initial_stress.value->name);
77 }
78
79 // TODO (naumov) this information is runtime information and I'm not
80 // sure how to put it into the reflected data structure. The
81 // reflectWithName function also supports only a single return value.
82 if (name.starts_with("material_state_variable_"))
83 {
84 name.remove_prefix(24);
85
86 auto const& internal_variables =
87 solid_material_.getInternalVariables();
88 if (auto const iv = std::find_if(
89 begin(internal_variables), end(internal_variables),
90 [&name](auto const& iv) { return iv.name == name; });
91 iv != end(internal_variables))
92 {
93 DBUG("Setting material state variable '{:s}'", name);
94 return ProcessLib::
96 values, material_states_,
98 DisplacementDim>::material_state_variables,
99 iv->reference);
100 }
101 return 0;
102 }
103
104 // TODO this logic could be pulled out of the local assembler into the
105 // process. That might lead to a slightly better performance due to less
106 // string comparisons.
108 name, values, current_states_);
109 }
110
112 std::function<std::span<double>(
114 MaterialStateVariables&)> const& get_values_span,
115 int const& n_components) const
116 {
120 get_values_span, n_components);
121 }
122
123 // TODO move to NumLib::ExtrapolatableElement
125 {
126 return integration_method_.getNumberOfPoints();
127 }
128
129 int getMaterialID() const
130 {
131 return process_data_.material_ids == nullptr
132 ? 0
133 : (*process_data_.material_ids)[element_.getID()];
134 }
135
137 DisplacementDim>::MaterialStateVariables const&
138 getMaterialStateVariablesAt(unsigned integration_point) const
139 {
140 return *material_states_[integration_point].material_state_variables;
141 }
142
143 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
144 Eigen::VectorXd const& /*local_x_prev*/,
145 double const /*t*/, double const /*dt*/,
146 int const /*process_id*/) override
147 {
148 unsigned const n_integration_points =
149 integration_method_.getNumberOfPoints();
150
151 for (unsigned ip = 0; ip < n_integration_points; ip++)
152 {
153 // TODO re-evaluate part of the assembly in order to be consistent?
154 material_states_[ip].pushBackState();
155 }
156
157 for (unsigned ip = 0; ip < n_integration_points; ip++)
158 {
160 }
161 }
162
164 {
165 using Self =
167
169 &Self::current_states_, &Self::output_data_);
170 }
171
172protected:
175
176 std::vector<typename ConstitutiveTraits::StatefulData>
177 current_states_; // TODO maybe do not store but rather re-evaluate for
178 // state update
179 std::vector<typename ConstitutiveTraits::StatefulDataPrev> prev_states_;
180
181 // Material state is special, because it contains both the current and the
182 // old state.
183 std::vector<MaterialStateData<DisplacementDim>> material_states_;
184
188
189 typename ConstitutiveTraits::SolidConstitutiveRelation const&
191
192 std::vector<typename ConstitutiveTraits::OutputData> output_data_;
193};
194
195} // namespace ProcessLib::ThermoRichardsMechanics
#define OGS_FATAL(...)
Definition Error.h:19
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
KV::KelvinVectorType< DisplacementDim > KelvinVector
std::size_t reflectSetIPData(std::string_view const name, double const *values, std::vector< IPData > &ip_data_vector)
auto reflectWithoutName(Accessors &&... accessors)
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)
ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > & process_data_
LocalAssemblerInterface(MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, ThermoRichardsMechanicsProcessData< DisplacementDim, ConstitutiveTraits > &process_data)
std::vector< double > getMaterialStateVariableInternalState(std::function< std::span< double >(typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables &)> const &get_values_span, int const &n_components) const
void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
std::size_t setIPDataInitialConditions(std::string_view name, double const *values, int const integration_order)
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt(unsigned integration_point) const