OGS
ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int DisplacementDim>
class ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >

Definition at line 79 of file LargeDeformationFEM.h.

#include <LargeDeformationFEM.h>

Inheritance diagram for ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< 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

 LargeDeformationLocalAssembler (LargeDeformationLocalAssembler const &)=delete
 
 LargeDeformationLocalAssembler (LargeDeformationLocalAssembler &&)=delete
 
 LargeDeformationLocalAssembler (MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, LargeDeformationProcessData< 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
 
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::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >
 LargeDeformationLocalAssemblerInterface (MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, LargeDeformationProcessData< DisplacementDim > &process_data)
 
std::size_t setIPDataInitialConditions (std::string_view name, double const *values, int const integration_order)
 Returns number of read integration points.
 
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 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
 

Static Private Attributes

static const int displacement_size
 

Additional Inherited Members

- Static Public Member Functions inherited from ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >
static auto getReflectionDataForOutput ()
 
- Protected Attributes inherited from ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >
LargeDeformationProcessData< 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_
 
ConstitutiveRelations::SolidConstitutiveRelation< DisplacementDim > const & solid_material_
 

Member Typedef Documentation

◆ BMatricesType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 88 of file LargeDeformationFEM.h.

◆ BMatrixType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::BMatrixType = typename BMatricesType::BMatrixType

Definition at line 90 of file LargeDeformationFEM.h.

◆ GMatricesType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GMatricesType = GMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 96 of file LargeDeformationFEM.h.

◆ GradientMatrixType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GradientMatrixType = typename GMatricesType::GradientMatrixType

Definition at line 98 of file LargeDeformationFEM.h.

◆ GradientVectorType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GradientVectorType = typename GMatricesType::GradientVectorType

Definition at line 97 of file LargeDeformationFEM.h.

◆ IpData

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

Definition at line 99 of file LargeDeformationFEM.h.

◆ NodalDisplacementVectorType

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

Definition at line 93 of file LargeDeformationFEM.h.

◆ NodalForceVectorType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalForceVectorType = typename BMatricesType::NodalForceVectorType

Definition at line 92 of file LargeDeformationFEM.h.

◆ NodalMatrixType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType

Definition at line 85 of file LargeDeformationFEM.h.

◆ NodalVectorType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType

Definition at line 86 of file LargeDeformationFEM.h.

◆ ShapeMatrices

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices

Definition at line 87 of file LargeDeformationFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatricesType

◆ StiffnessMatrixType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::StiffnessMatrixType = typename BMatricesType::StiffnessMatrixType

Definition at line 91 of file LargeDeformationFEM.h.

Constructor & Destructor Documentation

◆ LargeDeformationLocalAssembler() [1/3]

template<typename ShapeFunction , int DisplacementDim>
ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::LargeDeformationLocalAssembler ( LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim > const & )
delete

◆ LargeDeformationLocalAssembler() [2/3]

template<typename ShapeFunction , int DisplacementDim>
ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::LargeDeformationLocalAssembler ( LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim > && )
delete

◆ LargeDeformationLocalAssembler() [3/3]

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

Definition at line 109 of file LargeDeformationFEM.h.

115 : LargeDeformationLocalAssemblerInterface<DisplacementDim>(
116 e, integration_method, is_axially_symmetric, process_data)
117 {
118 unsigned const n_integration_points =
120
121 _ip_data.resize(n_integration_points);
122 _secondary_data.N.resize(n_integration_points);
123
124 auto const shape_matrices =
126 DisplacementDim>(
127 e, is_axially_symmetric, this->integration_method_);
128
129 for (unsigned ip = 0; ip < n_integration_points; ip++)
130 {
131 auto& ip_data = _ip_data[ip];
132 auto const& sm = shape_matrices[ip];
133 _ip_data[ip].integration_weight =
135 sm.integralMeasure * sm.detJ;
136
137 ip_data.N = sm.N;
138 ip_data.dNdx = sm.dNdx;
139
140 _secondary_data.N[ip] = shape_matrices[ip].N;
141 }
142 }
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::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data, NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, and ProcessLib::HeatTransportBHE::SecondaryData< ShapeMatrixType >::N.

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< 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 257 of file LargeDeformationFEM.h.

263 {
264 OGS_FATAL(
265 "LargeDeformationLocalAssembler: assembly without jacobian is not "
266 "implemented.");
267 }
#define OGS_FATAL(...)
Definition Error.h:26

References OGS_FATAL.

◆ assembleWithJacobian()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< 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 269 of file LargeDeformationFEM.h.

276 {
277 auto const local_matrix_size = local_x.size();
278
279 auto local_Jac = MathLib::createZeroedMatrix<StiffnessMatrixType>(
280 local_Jac_data, local_matrix_size, local_matrix_size);
281
282 auto local_b = MathLib::createZeroedVector<NodalDisplacementVectorType>(
283 local_b_data, local_matrix_size);
284
285 auto [u] = localDOF(local_x);
286 auto [u_prev] = localDOF(local_x_prev);
287
288 unsigned const n_integration_points =
290
292 x_position.setElementID(this->element_.getID());
293
294 typename ConstitutiveRelations::ConstitutiveSetting<DisplacementDim>
295 constitutive_setting;
296 auto const& medium =
297 *this->process_data_.media_map.getMedium(this->element_.getID());
298
299 for (unsigned ip = 0; ip < n_integration_points; ip++)
300 {
301 x_position.setIntegrationPoint(ip);
302 auto const& w = _ip_data[ip].integration_weight;
303 auto const& N = _ip_data[ip].N;
304 auto const& dNdx = _ip_data[ip].dNdx;
305
306 auto const x_coord =
307 NumLib::interpolateXCoordinate<ShapeFunction,
309 this->element_, N);
310
311 // For the 2D case the 33-component is needed (and the four entries
312 // of the non-symmetric matrix); In 3d there are nine entries.
313 GradientMatrixType G(DisplacementDim * DisplacementDim +
314 (DisplacementDim == 2 ? 1 : 0),
315 ShapeFunction::NPOINTS * DisplacementDim);
316 Deformation::computeGMatrix<DisplacementDim,
317 ShapeFunction::NPOINTS>(
318 dNdx, G, this->is_axially_symmetric_, N, x_coord);
319
320 GradientVectorType const grad_u = G * u;
321
322 auto const B_0 = LinearBMatrix::computeBMatrix<
323 DisplacementDim, ShapeFunction::NPOINTS,
325 dNdx, N, x_coord, this->is_axially_symmetric_);
326
327 auto const B_N = NonLinearBMatrix::computeBMatrix<
328 DisplacementDim, ShapeFunction::NPOINTS,
330 dNdx, N, grad_u, x_coord, this->is_axially_symmetric_);
331
332 auto const B = B_0 + B_N;
333
334 auto const CD = updateConstitutiveRelations(
335 u, u_prev, x_position, t, dt, _ip_data[ip],
336 constitutive_setting, medium, this->current_states_[ip],
337 this->prev_states_[ip], this->material_states_[ip],
338 this->output_data_[ip]);
339
340 auto const& sigma = this->current_states_[ip].stress_data.sigma;
341 auto const& b = *CD.volumetric_body_force;
342 auto const& C = CD.s_mech_data.stiffness_tensor;
343
344 local_b.noalias() -=
345 (B.transpose() * sigma - N_u_op(N).transpose() * b) * w;
346
347 auto const sigma_geom =
348 computeSigmaGeom<DisplacementDim, ShapeMatricesType>(
350
351 local_Jac.noalias() += G.transpose() * sigma_geom * G * w;
352
353 local_Jac.noalias() += B.transpose() * C * B * w;
354 }
355 }
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
typename GMatricesType::GradientVectorType GradientVectorType
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)
typename GMatricesType::GradientMatrixType GradientMatrixType
Eigen::Matrix< double, 3, 3 > kelvinVectorToTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
double interpolateXCoordinate(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
void computeGMatrix(DNDX_Type const &dNdx, GMatrixType &g_matrix, const bool is_axially_symmetric, N_Type const &N, const double radius)
Fills a G-matrix based on given shape function dN/dx values.
Definition GMatrix.h:25
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.
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, GradientVectorType const &grad_u, const double radius, const bool is_axially_symmetric)
std::vector< MaterialStateData< DisplacementDim > > material_states_
std::vector< typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > > prev_states_
std::vector< typename ConstitutiveRelations::StatefulData< DisplacementDim > > current_states_
std::vector< typename ConstitutiveRelations::OutputData< DisplacementDim > > output_data_

References ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::NonLinearBMatrix::computeBMatrix(), ProcessLib::Deformation::computeGMatrix(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::element_, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateXCoordinate(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, MathLib::KelvinVector::kelvinVectorToTensor(), ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::localDOF(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::N_u_op, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::output_data_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), and ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::updateConstitutiveRelations().

◆ getShapeMatrix()

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

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 393 of file LargeDeformationFEM.h.

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

References ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data, and ProcessLib::HeatTransportBHE::SecondaryData< ShapeMatrixType >::N.

◆ initializeConcrete()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete ( )
inlineoverridevirtual

Set initial stress from parameter.

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 144 of file LargeDeformationFEM.h.

145 {
146 unsigned const n_integration_points =
148 for (unsigned ip = 0; ip < n_integration_points; ip++)
149 {
150 auto const& ip_data = _ip_data[ip];
151
152 ParameterLib::SpatialPosition const x_position{
153 std::nullopt, this->element_.getID(), ip,
155 NumLib::interpolateCoordinates<ShapeFunction,
157 this->element_, ip_data.N))};
158
160 if (this->process_data_.initial_stress != nullptr)
161 {
162 this->current_states_[ip].stress_data.sigma.noalias() =
164 DisplacementDim>((*this->process_data_.initial_stress)(
165 std::numeric_limits<
166 double>::quiet_NaN() /* time independent */,
167 x_position));
168 }
169
170 double const t = 0; // TODO (naumov) pass t from top
171 auto& material_state = this->material_states_[ip];
172 this->solid_material_.initializeInternalStateVariables(
173 t, x_position, *material_state.material_state_variables);
174
175 material_state.pushBackState();
176 this->prev_states_[ip] = this->current_states_[ip];
177 }
178 }
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)
ConstitutiveRelations::SolidConstitutiveRelation< DisplacementDim > const & solid_material_

References ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::element_, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::solid_material_, and MathLib::KelvinVector::symmetricTensorToKelvinVector().

◆ localDOF()

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

◆ postTimestepConcrete()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< 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 357 of file LargeDeformationFEM.h.

361 {
362 unsigned const n_integration_points =
364
366 x_position.setElementID(this->element_.getID());
367
368 typename ConstitutiveRelations::ConstitutiveSetting<DisplacementDim>
369 constitutive_setting;
370
371 auto& medium =
372 *this->process_data_.media_map.getMedium(this->element_.getID());
373
374 for (unsigned ip = 0; ip < n_integration_points; ip++)
375 {
376 x_position.setIntegrationPoint(ip);
377
379 local_x, local_x_prev, x_position, t, dt, _ip_data[ip],
380 constitutive_setting, medium, this->current_states_[ip],
381 this->prev_states_[ip], this->material_states_[ip],
382 this->output_data_[ip]);
383
384 this->material_states_[ip].pushBackState();
385 }
386
387 for (unsigned ip = 0; ip < n_integration_points; ip++)
388 {
389 this->prev_states_[ip] = this->current_states_[ip];
390 }
391 }

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

◆ updateConstitutiveRelations()

template<typename ShapeFunction , int DisplacementDim>
ConstitutiveRelations::ConstitutiveData< DisplacementDim > ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< 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 181 of file LargeDeformationFEM.h.

198 {
199 auto const& N = ip_data.N;
200 auto const& dNdx = ip_data.dNdx;
201 auto const x_coord =
202 NumLib::interpolateXCoordinate<ShapeFunction, ShapeMatricesType>(
203 this->element_, N);
204
205 // For the 2D case the 33-component is needed (and the four entries
206 // of the non-symmetric matrix); In 3d there are nine entries.
208 DisplacementDim * DisplacementDim + (DisplacementDim == 2 ? 1 : 0),
209 ShapeFunction::NPOINTS * DisplacementDim);
210 Deformation::computeGMatrix<DisplacementDim, ShapeFunction::NPOINTS>(
211 dNdx, G, this->is_axially_symmetric_, N, x_coord);
212
213 GradientVectorType const grad_u = G * u;
214
215 auto const B_0 =
216 LinearBMatrix::computeBMatrix<DisplacementDim,
217 ShapeFunction::NPOINTS,
219 dNdx, N, x_coord, this->is_axially_symmetric_);
220
221 auto const B_N = NonLinearBMatrix::computeBMatrix<
222 DisplacementDim, ShapeFunction::NPOINTS,
223 typename BMatricesType::BMatrixType>(dNdx, N, grad_u, x_coord,
225
226 auto const B = B_0 + 0.5 * B_N; // For Green-Lagrange strain.
227
228 double const T_ref =
229 this->process_data_.reference_temperature
230 ? (*this->process_data_.reference_temperature)(t, x_position)[0]
231 : std::numeric_limits<double>::quiet_NaN();
232
233 typename ConstitutiveRelations::ConstitutiveModels<DisplacementDim>
234 models(this->process_data_, this->solid_material_);
235 typename ConstitutiveRelations::ConstitutiveTempData<DisplacementDim>
236 tmp;
237 typename ConstitutiveRelations::ConstitutiveData<DisplacementDim> CD;
238
239 output_data.eps_data.eps = B * u;
240 output_data.deformation_gradient_data.deformation_gradient =
241 grad_u + MathLib::VectorizedTensor::identity<DisplacementDim>();
242 output_data.deformation_gradient_data.volume_ratio =
244 output_data.deformation_gradient_data.deformation_gradient);
245
246 CS.eval(
247 models, t, dt, x_position, //
248 medium, //
249 T_ref, //
250 output_data.deformation_gradient_data, //
251 G * u_prev + MathLib::VectorizedTensor::identity<DisplacementDim>(),
252 current_state, prev_state, material_state, tmp, CD);
253
254 return CD;
255 }
double determinant(Eigen::MatrixBase< Derived > const &tensor)
Computes determinant of a vectorized tensor.

References ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::NonLinearBMatrix::computeBMatrix(), ProcessLib::LargeDeformation::ConstitutiveRelations::OutputData< DisplacementDim >::deformation_gradient_data, MathLib::VectorizedTensor::determinant(), ProcessLib::LargeDeformation::IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim >::dNdx, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::element_, ProcessLib::LargeDeformation::ConstitutiveRelations::OutputData< DisplacementDim >::eps_data, ProcessLib::LargeDeformation::ConstitutiveRelations::ConstitutiveSetting< DisplacementDim >::eval(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, ProcessLib::LargeDeformation::IntegrationPointData< BMatricesType, ShapeMatricesType, DisplacementDim >::N, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, and ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::solid_material_.

Referenced by ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian(), and ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete().

Member Data Documentation

◆ _ip_data

◆ _secondary_data

◆ displacement_size

template<typename ShapeFunction , int DisplacementDim>
const int ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::displacement_size
staticprivate
Initial value:
=
ShapeFunction::NPOINTS * DisplacementDim

Definition at line 414 of file LargeDeformationFEM.h.

◆ N_u_op

template<typename ShapeFunction , int DisplacementDim>
constexpr auto& ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< 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 102 of file LargeDeformationFEM.h.

Referenced by ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian().


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