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 59 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 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
 
ConstitutiveRelations::ConstitutiveData< DisplacementDim > updateConstitutiveRelations (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, IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim > &ip_data, 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
 
void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &, std::vector< double > &, 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_M_data, std::vector< double > &local_K_data, 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.
 
- 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

◆ BMatricesType

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

Definition at line 68 of file SmallDeformationFEM.h.

◆ BMatrixType

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

Definition at line 70 of file SmallDeformationFEM.h.

◆ GMatricesType

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

Definition at line 76 of file SmallDeformationFEM.h.

◆ GradientMatrixType

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

Definition at line 78 of file SmallDeformationFEM.h.

◆ GradientVectorType

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

Definition at line 77 of file SmallDeformationFEM.h.

◆ IpData

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::IpData
Initial value:
IntegrationPointData<BMatricesType, ShapeMatricesType, DisplacementDim>

Definition at line 79 of file SmallDeformationFEM.h.

◆ 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 73 of file SmallDeformationFEM.h.

◆ NodalForceVectorType

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

Definition at line 72 of file SmallDeformationFEM.h.

◆ NodalMatrixType

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

Definition at line 65 of file SmallDeformationFEM.h.

◆ NodalVectorType

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

Definition at line 66 of file SmallDeformationFEM.h.

◆ ShapeMatrices

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

Definition at line 67 of file SmallDeformationFEM.h.

◆ ShapeMatricesType

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

◆ StiffnessMatrixType

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

Definition at line 71 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 89 of file SmallDeformationFEM.h.

95 : SmallDeformationLocalAssemblerInterface<DisplacementDim>(
96 e, integration_method, is_axially_symmetric, process_data)
97 {
98 unsigned const n_integration_points =
100
101 _ip_data.resize(n_integration_points);
102 _secondary_data.N.resize(n_integration_points);
103
104 auto const shape_matrices =
106 DisplacementDim>(
107 e, is_axially_symmetric, this->integration_method_);
108
109 for (unsigned ip = 0; ip < n_integration_points; ip++)
110 {
111 auto& ip_data = _ip_data[ip];
112 auto const& sm = shape_matrices[ip];
113 _ip_data[ip].integration_weight =
115 sm.integralMeasure * sm.detJ;
116
117 ip_data.N = sm.N;
118 ip_data.dNdx = sm.dNdx;
119
120 _secondary_data.N[ip] = shape_matrices[ip].N;
121 }
122 }
double getWeight() const
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N

References ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data, NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, and ProcessLib::HeatTransportBHE::SecondaryData< ShapeMatrixType >::N.

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 214 of file SmallDeformationFEM.h.

220 {
221 OGS_FATAL(
222 "SmallDeformationLocalAssembler: assembly without jacobian is not "
223 "implemented.");
224 }
#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 > & ,
std::vector< double > & ,
std::vector< double > & local_b_data,
std::vector< double > & local_Jac_data )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 226 of file SmallDeformationFEM.h.

233 {
234 auto const local_matrix_size = local_x.size();
235
236 auto local_Jac = MathLib::createZeroedMatrix<StiffnessMatrixType>(
237 local_Jac_data, local_matrix_size, local_matrix_size);
238
239 auto local_b = MathLib::createZeroedVector<NodalDisplacementVectorType>(
240 local_b_data, local_matrix_size);
241
242 auto [u] = localDOF(local_x);
243 auto [u_prev] = localDOF(local_x_prev);
244
245 unsigned const n_integration_points =
247
249 x_position.setElementID(this->element_.getID());
250
251 typename ConstitutiveRelations::ConstitutiveSetting<DisplacementDim>
253 auto const& medium =
254 *this->process_data_.media_map.getMedium(this->element_.getID());
255
256 for (unsigned ip = 0; ip < n_integration_points; ip++)
257 {
258 x_position.setIntegrationPoint(ip);
259 auto const& w = _ip_data[ip].integration_weight;
260 auto const& N = _ip_data[ip].N;
261 auto const& dNdx = _ip_data[ip].dNdx;
262
263 auto const x_coord =
264 NumLib::interpolateXCoordinate<ShapeFunction,
266 this->element_, N);
267 auto const B = LinearBMatrix::computeBMatrix<
268 DisplacementDim, ShapeFunction::NPOINTS,
270 dNdx, N, x_coord, this->is_axially_symmetric_);
271
272 auto const CD = updateConstitutiveRelations(
273 u, u_prev, x_position, t, dt, _ip_data[ip],
274 constitutive_setting, medium, this->current_states_[ip],
275 this->prev_states_[ip], this->material_states_[ip],
276 this->output_data_[ip]);
277
278 auto const& sigma = this->current_states_[ip].stress_data.sigma;
279 auto const& b = *CD.volumetric_body_force;
280 auto const& C = CD.s_mech_data.stiffness_tensor;
281
282 local_b.noalias() -=
283 (B.transpose() * sigma - N_u_op(N).transpose() * b) * w;
284 local_Jac.noalias() += B.transpose() * C * B * w;
285 }
286 }
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
ConstitutiveRelations::ConstitutiveData< DisplacementDim > updateConstitutiveRelations(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, IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim > &ip_data, 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
static constexpr auto localDOF(std::vector< double > const &x)
double interpolateXCoordinate(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, const double radius, const bool is_axially_symmetric)
Fills a B-matrix based on given shape function dN/dx values.
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::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::constitutive_setting, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateXCoordinate(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::localDOF(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::N_u_op, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::output_data_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), and ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::updateConstitutiveRelations().

◆ 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 324 of file SmallDeformationFEM.h.

327 {
329 DisplacementDim, ShapeFunction, ShapeMatricesType,
332 GradientMatrixType>(local_x, nodal_values,
334 this->current_states_, this->output_data_,
335 this->element_, this->is_axially_symmetric_);
336 }
typename BMatricesType::NodalForceVectorType NodalDisplacementVectorType
typename GMatricesType::GradientMatrixType GradientMatrixType
typename GMatricesType::GradientVectorType GradientVectorType
std::vector< double > const & getMaterialForces(std::vector< double > const &local_x, std::vector< double > &nodal_values, IntegrationMethod const &_integration_method, IPData const &_ip_data, StressData const &stress_data, OutputData const &output_data, MeshLib::Element const &element, bool const is_axially_symmetric)

References ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, ProcessLib::SmallDeformation::getMaterialForces(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, 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 338 of file SmallDeformationFEM.h.

340 {
341 auto const& N = _secondary_data.N[integration_point];
342
343 // assumes N is stored contiguously in memory
344 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
345 }

References ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data, and ProcessLib::HeatTransportBHE::SecondaryData< ShapeMatrixType >::N.

◆ 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 124 of file SmallDeformationFEM.h.

125 {
126 unsigned const n_integration_points =
128 for (unsigned ip = 0; ip < n_integration_points; ip++)
129 {
130 auto const& ip_data = _ip_data[ip];
131
132 ParameterLib::SpatialPosition const x_position{
133 std::nullopt, this->element_.getID(), ip,
135 NumLib::interpolateCoordinates<ShapeFunction,
137 this->element_, ip_data.N))};
138
140 if (this->process_data_.initial_stress != nullptr)
141 {
142 this->current_states_[ip].stress_data.sigma.noalias() =
144 DisplacementDim>((*this->process_data_.initial_stress)(
145 std::numeric_limits<
146 double>::quiet_NaN() /* time independent */,
147 x_position));
148 }
149
150 double const t = 0; // TODO (naumov) pass t from top
151 auto& material_state = this->material_states_[ip];
152 this->solid_material_.initializeInternalStateVariables(
153 t, x_position, *material_state.material_state_variables);
154
155 material_state.pushBackState();
156 }
157
158 for (unsigned ip = 0; ip < n_integration_points; ip++)
159 {
160 this->prev_states_[ip] = this->current_states_[ip];
161 }
162 }
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
MaterialLib::Solids::MechanicsBase< DisplacementDim > const & solid_material_

References ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), 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>
static constexpr auto ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::localDOF ( std::vector< double > const & x)
inlinestaticconstexprprivate

◆ 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 288 of file SmallDeformationFEM.h.

292 {
293 unsigned const n_integration_points =
295
297 x_position.setElementID(this->element_.getID());
298
299 typename ConstitutiveRelations::ConstitutiveSetting<DisplacementDim>
301
302 auto const& medium =
303 *this->process_data_.media_map.getMedium(this->element_.getID());
304
305 for (unsigned ip = 0; ip < n_integration_points; ip++)
306 {
307 x_position.setIntegrationPoint(ip);
308
310 local_x, local_x_prev, x_position, t, dt, _ip_data[ip],
311 constitutive_setting, medium, this->current_states_[ip],
312 this->prev_states_[ip], this->material_states_[ip],
313 this->output_data_[ip]);
314
315 this->material_states_[ip].pushBackState();
316 }
317
318 for (unsigned ip = 0; ip < n_integration_points; ip++)
319 {
320 this->prev_states_[ip] = this->current_states_[ip];
321 }
322 }

References ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::constitutive_setting, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::output_data_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), and ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::updateConstitutiveRelations().

◆ updateConstitutiveRelations()

template<typename ShapeFunction , int DisplacementDim>
ConstitutiveRelations::ConstitutiveData< DisplacementDim > ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::updateConstitutiveRelations ( 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,
IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim > & ip_data,
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

Definition at line 165 of file SmallDeformationFEM.h.

182 {
183 auto const& N = ip_data.N;
184 auto const& dNdx = ip_data.dNdx;
185 auto const x_coord =
186 NumLib::interpolateXCoordinate<ShapeFunction, ShapeMatricesType>(
187 this->element_, N);
188 auto const B =
189 LinearBMatrix::computeBMatrix<DisplacementDim,
190 ShapeFunction::NPOINTS,
192 dNdx, N, x_coord, this->is_axially_symmetric_);
193
194 double const T_ref =
195 this->process_data_.reference_temperature
196 ? (*this->process_data_.reference_temperature)(t, x_position)[0]
197 : std::numeric_limits<double>::quiet_NaN();
198
199 typename ConstitutiveRelations::ConstitutiveModels<DisplacementDim>
200 models{this->process_data_, this->solid_material_};
201 typename ConstitutiveRelations::ConstitutiveTempData<DisplacementDim>
202 tmp;
203 typename ConstitutiveRelations::ConstitutiveData<DisplacementDim> CD;
204
205 CS.eval(models, t, dt, x_position, //
206 medium, //
207 T_ref, B * u, B * u_prev, //
208 current_state, prev_state, material_state, tmp, output_data,
209 CD);
210
211 return CD;
212 }

References ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::SmallDeformation::IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim >::dNdx, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, ProcessLib::SmallDeformation::ConstitutiveRelations::ConstitutiveSetting< DisplacementDim >::eval(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, ProcessLib::SmallDeformation::IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim >::N, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, and ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::solid_material_.

Referenced by ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian(), and ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete().

Member Data Documentation

◆ _ip_data

◆ _secondary_data

◆ N_u_op

template<typename ShapeFunction , int DisplacementDim>
constexpr 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 82 of file SmallDeformationFEM.h.

Referenced by ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian().


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