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