OGS
LocalAssemblerInterface.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <vector>
14
26
27namespace ProcessLib
28{
29namespace LargeDeformation
30{
31template <int DisplacementDim>
35{
37 MeshLib::Element const& e,
38 NumLib::GenericIntegrationMethod const& integration_method,
39 bool const is_axially_symmetric,
41 : process_data_(process_data),
42 integration_method_(integration_method),
43 element_(e),
44 is_axially_symmetric_(is_axially_symmetric),
45 solid_material_(MaterialLib::Solids::selectSolidConstitutiveRelation(
46 process_data_.solid_materials, process_data_.material_ids,
47 element_.getID()))
48 {
49 unsigned const n_integration_points =
51
52 material_states_.reserve(n_integration_points);
53 for (unsigned ip = 0; ip < n_integration_points; ++ip)
54 {
55 material_states_.emplace_back(
56 solid_material_.createMaterialStateVariables());
57 }
58
59 current_states_.resize(n_integration_points);
60 prev_states_.resize(n_integration_points);
61 output_data_.resize(n_integration_points);
62 }
64 std::size_t setIPDataInitialConditions(std::string_view name,
65 double const* values,
66 int const integration_order)
67 {
68 if (integration_order !=
69 static_cast<int>(integration_method_.getIntegrationOrder()))
70 {
72 "Setting integration point initial conditions; The integration "
73 "order of the local assembler for element {:d} is different "
74 "from the integration order in the initial condition.",
75 element_.getID());
76 }
77
78 if (name.starts_with("material_state_variable_"))
79 {
80 name.remove_prefix(24);
81
82 auto const& internal_variables =
83 solid_material_.getInternalVariables();
84 if (auto const iv = std::find_if(
85 begin(internal_variables), end(internal_variables),
86 [&name](auto const& iv) { return iv.name == name; });
87 iv != end(internal_variables))
88 {
89 DBUG("Setting material state variable '{:s}'", name);
90 return ProcessLib::
92 values, material_states_,
94 DisplacementDim>::material_state_variables,
95 iv->reference);
96 }
97
98 WARN(
99 "Could not find variable {:s} in solid material model's "
100 "internal variables.",
101 name);
102 return 0;
103 }
104
105 // TODO this logic could be pulled out of the local assembler into the
106 // process. That might lead to a slightly better performance due to less
107 // string comparisons.
109 name, values, current_states_);
110 }
111
112 // TODO move to NumLib::ExtrapolatableElement
114 {
116 }
117
118 int getMaterialID() const
119 {
120 return process_data_.material_ids == nullptr
121 ? 0
122 : (*process_data_.material_ids)[element_.getID()];
123 }
124
126 std::function<std::span<double>(
128 MaterialStateVariables&)> const& get_values_span,
129 int const& n_components) const
130 {
134 get_values_span, n_components);
135 }
136
138 DisplacementDim>::MaterialStateVariables const&
139 getMaterialStateVariablesAt(unsigned integration_point) const
140 {
141 return *material_states_[integration_point].material_state_variables;
142 }
143
145 {
147
149 &Self::current_states_, &Self::output_data_);
150 }
151
152protected:
154
155 // Material state is special, because it contains both the current and the
156 // old state.
157 std::vector<MaterialStateData<DisplacementDim>> material_states_;
158 std::vector<typename ConstitutiveRelations::StatefulData<DisplacementDim>>
159 current_states_; // TODO maybe do not store but rather re-evaluate for
160 // state update
161 std::vector<
164 std::vector<typename ConstitutiveRelations::OutputData<DisplacementDim>>
166
171 DisplacementDim> const& solid_material_;
172};
173
174} // namespace LargeDeformation
175} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
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)
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< MaterialStateData< DisplacementDim > > material_states_
std::vector< typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > > prev_states_
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt(unsigned integration_point) const
std::size_t setIPDataInitialConditions(std::string_view name, double const *values, int const integration_order)
Returns number of read integration points.
std::vector< double > getMaterialStateVariableInternalState(std::function< std::span< double >(typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables &)> const &get_values_span, int const &n_components) const
LargeDeformationLocalAssemblerInterface(MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, LargeDeformationProcessData< DisplacementDim > &process_data)
std::vector< typename ConstitutiveRelations::StatefulData< DisplacementDim > > current_states_
std::vector< typename ConstitutiveRelations::OutputData< DisplacementDim > > output_data_
ConstitutiveRelations::SolidConstitutiveRelation< DisplacementDim > const & solid_material_