OGS
LocalAssemblerInterface.h
Go to the documentation of this file.
1
11#pragma once
12
22#include "TH2MProcessData.h"
23
24namespace ProcessLib
25{
26namespace TH2M
27{
28template <int DisplacementDim>
29struct LocalAssemblerInterface : public ProcessLib::LocalAssemblerInterface,
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 element_.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 for (unsigned ip = 0; ip < n_integration_points; ++ip)
54 {
55 material_states_.emplace_back(
56 solid_material_.createMaterialStateVariables());
57
58 // Set initial stress field to zero. Might be overwritten by
59 // integration point data or initial stress.
60 current_states_[ip].eff_stress_data.sigma =
62 }
63 }
64
65 virtual std::size_t setIPDataInitialConditions(
66 std::string_view name, double const* values,
67 int const integration_order) = 0;
68
69 // TODO move to NumLib::ExtrapolatableElement
74
75 int getMaterialID() const
76 {
77 return process_data_.material_ids == nullptr
78 ? 0
79 : (*process_data_.material_ids)[element_.getID()];
80 }
81
83 std::function<std::span<double>(
85 DisplacementDim>::MaterialStateVariables&)> const&
86 get_values_span,
87 int const& n_components) const
88 {
92 DisplacementDim>::material_state_variables,
93 get_values_span, n_components);
94 }
95
97 DisplacementDim>::MaterialStateVariables const&
98 getMaterialStateVariablesAt(unsigned integration_point) const
99 {
100 return *material_states_[integration_point].material_state_variables;
101 }
102
104 {
106
108 &Self::current_states_, &Self::output_data_);
109 }
110
112
113 std::vector<typename ConstitutiveRelations::StatefulData<DisplacementDim>>
114 current_states_; // TODO maybe do not store but rather re-evaluate for
115 // state update
116 std::vector<
119
120 // Material state is special, because it contains both the current and the
121 // old state.
122 std::vector<ConstitutiveRelations::MaterialStateData<DisplacementDim>>
124 std::vector<ConstitutiveRelations::OutputData<DisplacementDim>>
126
132};
133
134} // namespace TH2M
135} // namespace ProcessLib
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
auto reflectWithoutName(Accessors &&... accessors)
KV::KelvinVectorType< DisplacementDim > KelvinVector
Definition Base.h:26
std::vector< double > getIntegrationPointDataMaterialStateVariables(IntegrationPointDataVector const &ip_data_vector, MemberType member, std::function< std::span< double >(MaterialStateVariables &)> get_values_span, int const n_components)
NumLib::GenericIntegrationMethod const & integration_method_
LocalAssemblerInterface(MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, TH2MProcessData< DisplacementDim > &process_data)
std::vector< double > getMaterialStateVariableInternalState(std::function< std::span< double >(typename ConstitutiveRelations::SolidConstitutiveRelation< DisplacementDim >::MaterialStateVariables &)> const &get_values_span, int const &n_components) const
TH2MProcessData< DisplacementDim > & process_data_
virtual std::size_t setIPDataInitialConditions(std::string_view name, double const *values, int const integration_order)=0
ConstitutiveRelations::SolidConstitutiveRelation< DisplacementDim > const & solid_material_
std::vector< typename ConstitutiveRelations::StatefulData< DisplacementDim > > current_states_
std::vector< typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > > prev_states_
ConstitutiveRelations::SolidConstitutiveRelation< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt(unsigned integration_point) const
std::vector< ConstitutiveRelations::OutputData< DisplacementDim > > output_data_
std::vector< ConstitutiveRelations::MaterialStateData< DisplacementDim > > material_states_