OGS
ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int DisplacementDim>
class ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >

Definition at line 60 of file SmallDeformationFEM.h.

#include <SmallDeformationFEM.h>

Inheritance diagram for ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >:
[legend]

Public Types

using ShapeMatricesType
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
using BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>
using BMatrixType = typename BMatricesType::BMatrixType
using BBarMatrixType = typename BMatricesType::BBarMatrixType
using StiffnessMatrixType = typename BMatricesType::StiffnessMatrixType
using NodalForceVectorType = typename BMatricesType::NodalForceVectorType
using NodalDisplacementVectorType
using GMatricesType = GMatrixPolicyType<ShapeFunction, DisplacementDim>
using GradientVectorType = typename GMatricesType::GradientVectorType
using GradientMatrixType = typename GMatricesType::GradientMatrixType
using IpData

Public Member Functions

 SmallDeformationLocalAssembler (SmallDeformationLocalAssembler const &)=delete
 SmallDeformationLocalAssembler (SmallDeformationLocalAssembler &&)=delete
 SmallDeformationLocalAssembler (MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, SmallDeformationProcessData< DisplacementDim > &process_data)
void initializeConcrete () override
void setInitialConditionsConcrete (Eigen::VectorXd const local_x, double const, int const) override
ConstitutiveRelations::ConstitutiveData< DisplacementDim > updateConstitutiveRelations (BMatrixType const &B, Eigen::Ref< Eigen::VectorXd const > const &u, Eigen::Ref< Eigen::VectorXd const > const &u_prev, ParameterLib::SpatialPosition const &x_position, double const t, double const dt, typename ConstitutiveRelations::ConstitutiveSetting< DisplacementDim > &CS, MaterialPropertyLib::Medium const &medium, typename ConstitutiveRelations::StatefulData< DisplacementDim > &current_state, typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > const &prev_state, MaterialStateData< DisplacementDim > &material_state, typename ConstitutiveRelations::OutputData< DisplacementDim > &output_data) const
void assemble (double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
std::optional< BBarMatrixTypegetDilatationalBBarMatrix () const
void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
void postTimestepConcrete (Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, double const t, double const dt, int const) override
std::vector< double > const & getMaterialForces (std::vector< double > const &local_x, std::vector< double > &nodal_values) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
Public Member Functions inherited from ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >
 SmallDeformationLocalAssemblerInterface (MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, SmallDeformationProcessData< DisplacementDim > &process_data)
std::size_t setIPDataInitialConditions (std::string_view name, double const *values, int const integration_order)
 Returns number of read integration points.
void computeSecondaryVariableConcrete (double const, double const, Eigen::VectorXd const &, Eigen::VectorXd const &) override
unsigned getNumberOfIntegrationPoints () const
int getMaterialID () const
std::vector< double > getMaterialStateVariableInternalState (std::function< std::span< double >(typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables &)> const &get_values_span, int const &n_components) const
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt (unsigned integration_point) const
Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
virtual void setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, int const process_id)
virtual void initialize (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
virtual void preAssemble (double const, double const, std::vector< double > const &)
virtual void assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
virtual void assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
virtual void computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
virtual void preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
virtual void postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
void postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double > > const &) const
 Fits to staggered scheme.
virtual std::optional< VectorSegmentgetVectorDeformationSegment () const
Public Member Functions inherited from ProcessLib::SmallDeformation::MaterialForcesInterface
virtual ~MaterialForcesInterface ()=default
Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default

Static Public Attributes

static constexpr auto & N_u_op

Static Private Member Functions

static constexpr auto localDOF (std::vector< double > const &x)

Private Attributes

std::vector< IpData, Eigen::aligned_allocator< IpData > > ip_data_
SecondaryData< typename ShapeMatrices::ShapeType > secondary_data_

Additional Inherited Members

Static Public Member Functions inherited from ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >
static auto getReflectionDataForOutput ()
Protected Attributes inherited from ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >
SmallDeformationProcessData< DisplacementDim > & process_data_
std::vector< MaterialStateData< DisplacementDim > > material_states_
std::vector< typename ConstitutiveRelations::StatefulData< DisplacementDim > > current_states_
std::vector< typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > > prev_states_
std::vector< typename ConstitutiveRelations::OutputData< DisplacementDim > > output_data_
NumLib::GenericIntegrationMethod const & integration_method_
MeshLib::Element const & element_
bool const is_axially_symmetric_
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material_
ConstitutiveRelations::ConstitutiveSetting< DisplacementDim > constitutive_setting

Member Typedef Documentation

◆ BBarMatrixType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::BBarMatrixType = typename BMatricesType::BBarMatrixType

Definition at line 72 of file SmallDeformationFEM.h.

◆ BMatricesType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 69 of file SmallDeformationFEM.h.

◆ BMatrixType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::BMatrixType = typename BMatricesType::BMatrixType

Definition at line 71 of file SmallDeformationFEM.h.

◆ GMatricesType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GMatricesType = GMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 78 of file SmallDeformationFEM.h.

◆ GradientMatrixType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GradientMatrixType = typename GMatricesType::GradientMatrixType

Definition at line 80 of file SmallDeformationFEM.h.

◆ GradientVectorType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GradientVectorType = typename GMatricesType::GradientVectorType

Definition at line 79 of file SmallDeformationFEM.h.

◆ IpData

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::IpData

◆ NodalDisplacementVectorType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalDisplacementVectorType
Initial value:
VectorType< _number_of_dof > NodalForceVectorType
Rhs residual.

Definition at line 75 of file SmallDeformationFEM.h.

◆ NodalForceVectorType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalForceVectorType = typename BMatricesType::NodalForceVectorType

Definition at line 74 of file SmallDeformationFEM.h.

◆ NodalMatrixType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType

Definition at line 66 of file SmallDeformationFEM.h.

◆ NodalVectorType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType

Definition at line 67 of file SmallDeformationFEM.h.

◆ ShapeMatrices

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices

Definition at line 68 of file SmallDeformationFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatricesType
Initial value:
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType

Definition at line 64 of file SmallDeformationFEM.h.

◆ StiffnessMatrixType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::StiffnessMatrixType = typename BMatricesType::StiffnessMatrixType

Definition at line 73 of file SmallDeformationFEM.h.

Constructor & Destructor Documentation

◆ SmallDeformationLocalAssembler() [1/3]

template<typename ShapeFunction, int DisplacementDim>
ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::SmallDeformationLocalAssembler ( SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim > const & )
delete

◆ SmallDeformationLocalAssembler() [2/3]

template<typename ShapeFunction, int DisplacementDim>
ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::SmallDeformationLocalAssembler ( SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim > && )
delete

◆ SmallDeformationLocalAssembler() [3/3]

template<typename ShapeFunction, int DisplacementDim>
ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::SmallDeformationLocalAssembler ( MeshLib::Element const & e,
std::size_t const ,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
SmallDeformationProcessData< DisplacementDim > & process_data )
inline

Definition at line 91 of file SmallDeformationFEM.h.

99 {
100 unsigned const n_integration_points =
101 this->integration_method_.getNumberOfPoints();
102
105
106 auto const shape_matrices =
110
111 for (unsigned ip = 0; ip < n_integration_points; ip++)
112 {
113 auto& ip_data = ip_data_[ip];
114 auto const& sm = shape_matrices[ip];
115 ip_data_[ip].integration_weight =
116 this->integration_method_.getWeightedPoint(ip).getWeight() *
117 sm.integralMeasure * sm.detJ;
118
119 ip_data.N_u = sm.N;
120 ip_data.dNdx_u = sm.dNdx;
121
123 }
124 }
SecondaryData< typename ShapeMatrices::ShapeType > secondary_data_
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
std::vector< IpData, Eigen::aligned_allocator< IpData > > ip_data_
SmallDeformationLocalAssemblerInterface(MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, SmallDeformationProcessData< DisplacementDim > &process_data)

References ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::SmallDeformationLocalAssemblerInterface(), NumLib::initShapeMatrices(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, ip_data_, and secondary_data_.

Member Function Documentation

◆ assemble()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assemble ( double const ,
double const ,
std::vector< double > const & ,
std::vector< double > const & ,
std::vector< double > & ,
std::vector< double > & ,
std::vector< double > &  )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 235 of file SmallDeformationFEM.h.

241 {
242 OGS_FATAL(
243 "SmallDeformationLocalAssembler: assembly without jacobian is not "
244 "implemented.");
245 }
#define OGS_FATAL(...)
Definition Error.h:26

References OGS_FATAL.

◆ assembleWithJacobian()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian ( double const t,
double const dt,
std::vector< double > const & local_x,
std::vector< double > const & local_x_prev,
std::vector< double > & local_b_data,
std::vector< double > & local_Jac_data )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 261 of file SmallDeformationFEM.h.

266 {
267 auto const local_matrix_size = local_x.size();
268
271
274
275 auto [u] = localDOF(local_x);
276 auto [u_prev] = localDOF(local_x_prev);
277
278 unsigned const n_integration_points =
279 this->integration_method_.getNumberOfPoints();
280
283 auto const& medium =
284 *this->process_data_.media_map.getMedium(this->element_.getID());
285
287
288 for (unsigned ip = 0; ip < n_integration_points; ip++)
289 {
290 auto const& w = ip_data_[ip].integration_weight;
291 auto const& N = ip_data_[ip].N_u;
292 auto const& dNdx = ip_data_[ip].dNdx_u;
293
295 std::nullopt, this->element_.getID(),
299 this->element_, N))};
300
301 auto const x_coord =
304 this->element_, N);
309
310 auto const CD = updateConstitutiveRelations(
312 this->current_states_[ip], this->prev_states_[ip],
313 this->material_states_[ip], this->output_data_[ip]);
314
315 auto const& sigma = this->current_states_[ip].stress_data.sigma;
316 auto const& b = *CD.volumetric_body_force;
317 auto const& C = CD.s_mech_data.stiffness_tensor;
318
319 local_b.noalias() -=
320 (B.transpose() * sigma - N_u_op(N).transpose() * b) * w;
321 local_Jac.noalias() += B.transpose() * C * B * w;
322 }
323 }
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
ConstitutiveRelations::ConstitutiveData< DisplacementDim > updateConstitutiveRelations(BMatrixType const &B, Eigen::Ref< Eigen::VectorXd const > const &u, Eigen::Ref< Eigen::VectorXd const > const &u_prev, ParameterLib::SpatialPosition const &x_position, double const t, double const dt, typename ConstitutiveRelations::ConstitutiveSetting< DisplacementDim > &CS, MaterialPropertyLib::Medium const &medium, typename ConstitutiveRelations::StatefulData< DisplacementDim > &current_state, typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > const &prev_state, MaterialStateData< DisplacementDim > &material_state, typename ConstitutiveRelations::OutputData< DisplacementDim > &output_data) const
std::optional< BBarMatrixType > getDilatationalBBarMatrix() const
static constexpr auto localDOF(std::vector< double > const &x)
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
ConstitutiveRelations::ConstitutiveSetting< DisplacementDim > constitutive_setting
std::vector< typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > > prev_states_
std::vector< typename ConstitutiveRelations::StatefulData< DisplacementDim > > current_states_
std::vector< typename ConstitutiveRelations::OutputData< DisplacementDim > > output_data_
std::vector< MaterialStateData< DisplacementDim > > material_states_

References ProcessLib::LinearBMatrix::computeBMatrixPossiblyWithBbar(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::constitutive_setting, MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, getDilatationalBBarMatrix(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), NumLib::interpolateXCoordinate(), ip_data_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, localDOF(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, N_u_op, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::output_data_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, and updateConstitutiveRelations().

◆ getDilatationalBBarMatrix()

◆ getMaterialForces()

template<typename ShapeFunction, int DisplacementDim>
std::vector< double > const & ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::getMaterialForces ( std::vector< double > const & local_x,
std::vector< double > & nodal_values )
inlineoverridevirtual

Implements ProcessLib::SmallDeformation::MaterialForcesInterface.

Definition at line 378 of file SmallDeformationFEM.h.

References ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, ProcessLib::SmallDeformation::getMaterialForces(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, ip_data_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, and ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::output_data_.

◆ getShapeMatrix()

template<typename ShapeFunction, int DisplacementDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 392 of file SmallDeformationFEM.h.

394 {
395 auto const& N = secondary_data_.N[integration_point];
396
397 // assumes N is stored contiguously in memory
398 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
399 }

References secondary_data_.

◆ initializeConcrete()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete ( )
inlineoverridevirtual

Set initial stress from parameter.

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 126 of file SmallDeformationFEM.h.

127 {
128 unsigned const n_integration_points =
129 this->integration_method_.getNumberOfPoints();
130 for (unsigned ip = 0; ip < n_integration_points; ip++)
131 {
132 auto const& ip_data = ip_data_[ip];
133
135 std::nullopt, this->element_.getID(),
139 this->element_, ip_data.N_u))};
140
142 if (this->process_data_.initial_stress != nullptr)
143 {
144 this->current_states_[ip].stress_data.sigma.noalias() =
146 DisplacementDim>((*this->process_data_.initial_stress)(
148 double>::quiet_NaN() /* time independent */,
149 x_position));
150 }
151
152 double const t = 0; // TODO (naumov) pass t from top
153 auto& material_state = this->material_states_[ip];
154 this->solid_material_.initializeInternalStateVariables(
155 t, x_position, *material_state.material_state_variables);
156
157 material_state.pushBackState();
158 }
159
160 for (unsigned ip = 0; ip < n_integration_points; ip++)
161 {
162 this->prev_states_[ip] = this->current_states_[ip];
163 }
164 }
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material_

References ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), ip_data_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::solid_material_, and MathLib::KelvinVector::symmetricTensorToKelvinVector().

◆ localDOF()

template<typename ShapeFunction, int DisplacementDim>
constexpr auto ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::localDOF ( std::vector< double > const & x)
inlinestaticconstexprprivate

Definition at line 402 of file SmallDeformationFEM.h.

References NumLib::localDOF().

Referenced by assembleWithJacobian().

◆ postTimestepConcrete()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete ( Eigen::VectorXd const & local_x,
Eigen::VectorXd const & local_x_prev,
double const t,
double const dt,
int const  )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 325 of file SmallDeformationFEM.h.

329 {
330 unsigned const n_integration_points =
331 this->integration_method_.getNumberOfPoints();
332
335
336 auto const& medium =
337 *this->process_data_.media_map.getMedium(this->element_.getID());
338
340
341 for (unsigned ip = 0; ip < n_integration_points; ip++)
342 {
343 auto const& N = ip_data_[ip].N_u;
344 auto const& dNdx = ip_data_[ip].dNdx_u;
345
346 auto const x_coord =
349 this->element_, N);
350
352 std::nullopt, this->element_.getID(),
356 this->element_, N))};
357
362
366 this->prev_states_[ip], this->material_states_[ip],
367 this->output_data_[ip]);
368
369 this->material_states_[ip].pushBackState();
370 }
371
372 for (unsigned ip = 0; ip < n_integration_points; ip++)
373 {
374 this->prev_states_[ip] = this->current_states_[ip];
375 }
376 }

References ProcessLib::LinearBMatrix::computeBMatrixPossiblyWithBbar(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::constitutive_setting, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, getDilatationalBBarMatrix(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), NumLib::interpolateXCoordinate(), ip_data_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::output_data_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, and updateConstitutiveRelations().

◆ setInitialConditionsConcrete()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::setInitialConditionsConcrete ( Eigen::VectorXd const local_x,
double const ,
int const  )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 166 of file SmallDeformationFEM.h.

169 {
170 unsigned const n_integration_points =
171 this->integration_method_.getNumberOfPoints();
173
174 for (unsigned ip = 0; ip < n_integration_points; ip++)
175 {
176 auto const& N = ip_data_[ip].N_u;
177 auto const& dNdx = ip_data_[ip].dNdx_u;
179 std::nullopt, this->element_.getID(),
183 this->element_, N))};
184
185 auto const x_coord =
188 this->element_, N);
193
194 this->output_data_[ip].eps_data.eps.noalias() = B * local_x;
195 }
196 }

References ProcessLib::LinearBMatrix::computeBMatrixPossiblyWithBbar(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, getDilatationalBBarMatrix(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), NumLib::interpolateXCoordinate(), ip_data_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, and ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::output_data_.

◆ updateConstitutiveRelations()

template<typename ShapeFunction, int DisplacementDim>
ConstitutiveRelations::ConstitutiveData< DisplacementDim > ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::updateConstitutiveRelations ( BMatrixType const & B,
Eigen::Ref< Eigen::VectorXd const > const & u,
Eigen::Ref< Eigen::VectorXd const > const & u_prev,
ParameterLib::SpatialPosition const & x_position,
double const t,
double const dt,
typename ConstitutiveRelations::ConstitutiveSetting< DisplacementDim > & CS,
MaterialPropertyLib::Medium const & medium,
typename ConstitutiveRelations::StatefulData< DisplacementDim > & current_state,
typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > const & prev_state,
MaterialStateData< DisplacementDim > & material_state,
typename ConstitutiveRelations::OutputData< DisplacementDim > & output_data ) const
inline

Member Data Documentation

◆ ip_data_

template<typename ShapeFunction, int DisplacementDim>
std::vector<IpData, Eigen::aligned_allocator<IpData> > ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::ip_data_
private

◆ N_u_op

template<typename ShapeFunction, int DisplacementDim>
auto& ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::N_u_op
staticconstexpr
Initial value:
DisplacementDim, typename ShapeMatricesType::NodalRowVectorType>
constexpr Eigen::CwiseNullaryOp< EigenBlockMatrixViewFunctor< D, M >, typename EigenBlockMatrixViewFunctor< D, M >::Matrix > eigenBlockMatrixView(const Eigen::MatrixBase< M > &matrix)
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType

Definition at line 84 of file SmallDeformationFEM.h.

Referenced by assembleWithJacobian().

◆ secondary_data_

template<typename ShapeFunction, int DisplacementDim>
SecondaryData<typename ShapeMatrices::ShapeType> ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::secondary_data_
private

Definition at line 411 of file SmallDeformationFEM.h.

Referenced by SmallDeformationLocalAssembler(), and getShapeMatrix().


The documentation for this class was generated from the following file: