OGS
LocalAssemblerInterface.h
Go to the documentation of this file.
1
11#pragma once
12
23
24namespace ProcessLib
25{
26namespace RichardsMechanics
27{
28template <int DisplacementDim>
31{
33 MeshLib::Element const& e,
34 NumLib::GenericIntegrationMethod const& integration_method,
35 bool const is_axially_symmetric,
37 : process_data_(process_data),
38 integration_method_(integration_method),
39 element_(e),
40 is_axially_symmetric_(is_axially_symmetric),
41 solid_material_(MaterialLib::Solids::selectSolidConstitutiveRelation(
42 process_data_.solid_materials, process_data_.material_ids,
43 e.getID()))
44 {
45 unsigned const n_integration_points =
47
48 current_states_.resize(n_integration_points);
49 prev_states_.resize(n_integration_points);
50 output_data_.resize(n_integration_points);
51
52 material_states_.reserve(n_integration_points);
53
54 for (unsigned ip = 0; ip < n_integration_points; ++ip)
55 {
56 material_states_.emplace_back(
57 solid_material_.createMaterialStateVariables());
58
59 // Set initial strain field to zero.
60 std::get<StrainData<DisplacementDim>>(current_states_[ip]).eps =
62 }
63 }
64
65 std::size_t setIPDataInitialConditions(std::string_view name,
66 double const* values,
67 int const integration_order)
68 {
69 if (integration_order !=
70 static_cast<int>(integration_method_.getIntegrationOrder()))
71 {
73 "Setting integration point initial conditions; The integration "
74 "order of the local assembler for element {:d} is different "
75 "from the integration order in the initial condition.",
76 element_.getID());
77 }
78
79 if (name == "sigma" && process_data_.initial_stress.value != nullptr)
80 {
82 "Setting initial conditions for stress from integration "
83 "point data and from a parameter '{:s}' is not possible "
84 "simultaneously.",
85 process_data_.initial_stress.value->name);
86 }
87
88 if (name.starts_with("material_state_variable_"))
89 {
90 name.remove_prefix(24);
91
92 auto const& internal_variables =
93 solid_material_.getInternalVariables();
94 if (auto const iv = std::find_if(
95 begin(internal_variables), end(internal_variables),
96 [&name](auto const& iv) { return iv.name == name; });
97 iv != end(internal_variables))
98 {
99 DBUG("Setting material state variable '{:s}'", name);
100 return ProcessLib::
102 values, material_states_,
104 DisplacementDim>::material_state_variables,
105 iv->reference);
106 }
107 return 0;
108 }
109
110 // TODO this logic could be pulled out of the local assembler into the
111 // process. That might lead to a slightly better performance due to less
112 // string comparisons.
114 name, values, current_states_);
115 }
116
118 std::function<std::span<double>(
120 MaterialStateVariables&)> const& get_values_span,
121 int const& n_components) const
122 {
126 DisplacementDim>::material_state_variables,
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
151 {
153
155 &Self::current_states_, &Self::output_data_);
156 }
157
158 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
159 Eigen::VectorXd const& /*local_x_prev*/,
160 double const /*t*/, double const /*dt*/,
161 int const /*process_id*/) override final
162 {
163 unsigned const n_integration_points =
165
166 for (auto& s : material_states_)
167 {
168 s.pushBackState();
169 }
170
171 for (unsigned ip = 0; ip < n_integration_points; ip++)
172 {
174 }
175 }
176
177protected:
182
184
185 std::vector<StatefulData<DisplacementDim>> current_states_;
186 std::vector<StatefulDataPrev<DisplacementDim>> prev_states_;
187 std::vector<
190 std::vector<OutputData<DisplacementDim>> output_data_;
191};
192} // namespace RichardsMechanics
193} // namespace ProcessLib
#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:27
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< StatefulData< DisplacementDim > > current_states_
LocalAssemblerInterface(MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, RichardsMechanicsProcessData< DisplacementDim > &process_data)
std::vector< StatefulDataPrev< DisplacementDim > > prev_states_
RichardsMechanicsProcessData< DisplacementDim > & process_data_
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt(unsigned integration_point) const
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material_
std::vector< ProcessLib::ThermoRichardsMechanics::MaterialStateData< DisplacementDim > > material_states_
void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override final
std::vector< OutputData< DisplacementDim > > output_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
NumLib::GenericIntegrationMethod const & integration_method_
std::size_t setIPDataInitialConditions(std::string_view name, double const *values, int const integration_order)