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 74 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 GlobalDimNodalMatrixType
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 VectorTypeForFbar = typename BMatricesType::NodalForceVectorType
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
double computeOutputStrainData (double const detF0, BMatrixType const &B, GradientVectorType const &grad_u, Eigen::Ref< Eigen::VectorXd const > const &u, typename ConstitutiveRelations::OutputData< DisplacementDim > &output_data) const
ConstitutiveRelations::ConstitutiveData< DisplacementDim > updateConstitutiveRelations (double const alpha, GradientMatrixType const &G, 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< std::tuple< double, VectorTypeForFbar > > computeFBarVariables (bool const compute_detF0_only, Eigen::VectorXd const &local_x) 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
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_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 int getNumberOfVectorElementsForDeformation () const
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 85 of file LargeDeformationFEM.h.

◆ BMatrixType

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

Definition at line 87 of file LargeDeformationFEM.h.

◆ GlobalDimNodalMatrixType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GlobalDimNodalMatrixType
Initial value:
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType

Definition at line 81 of file LargeDeformationFEM.h.

◆ GMatricesType

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

Definition at line 93 of file LargeDeformationFEM.h.

◆ GradientMatrixType

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

Definition at line 95 of file LargeDeformationFEM.h.

◆ GradientVectorType

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

Definition at line 94 of file LargeDeformationFEM.h.

◆ IpData

◆ 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 90 of file LargeDeformationFEM.h.

◆ NodalForceVectorType

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

Definition at line 89 of file LargeDeformationFEM.h.

◆ NodalMatrixType

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

Definition at line 80 of file LargeDeformationFEM.h.

◆ NodalVectorType

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

Definition at line 83 of file LargeDeformationFEM.h.

◆ ShapeMatrices

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

Definition at line 84 of file LargeDeformationFEM.h.

◆ ShapeMatricesType

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

Definition at line 78 of file LargeDeformationFEM.h.

◆ StiffnessMatrixType

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

Definition at line 88 of file LargeDeformationFEM.h.

◆ VectorTypeForFbar

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

Definition at line 97 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.

117 {
118 unsigned const n_integration_points =
119 this->integration_method_.getNumberOfPoints();
120
123
124 auto const shape_matrices =
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 =
134 this->integration_method_.getWeightedPoint(ip).getWeight() *
135 sm.integralMeasure * sm.detJ;
136
137 ip_data.N_u = sm.N;
138 ip_data.dNdx_u = sm.dNdx;
139
141 }
142 }
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
LargeDeformationLocalAssemblerInterface(MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, LargeDeformationProcessData< DisplacementDim > &process_data)

References ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::LargeDeformationLocalAssemblerInterface(), _ip_data, _secondary_data, NumLib::initShapeMatrices(), and ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_.

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 280 of file LargeDeformationFEM.h.

286 {
287 OGS_FATAL(
288 "LargeDeformationLocalAssembler: assembly without jacobian is not "
289 "implemented.");
290 }
#define OGS_FATAL(...)
Definition Error.h:19

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 > & local_b_data,
std::vector< double > & local_Jac_data )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 319 of file LargeDeformationFEM.h.

324 {
325 auto const local_matrix_size = local_x.size();
326
329
332
333 auto [u] = localDOF(local_x);
334 auto [u_prev] = localDOF(local_x_prev);
335
336 unsigned const n_integration_points =
337 this->integration_method_.getNumberOfPoints();
338
341 auto const& medium =
342 *this->process_data_.media_map.getMedium(this->element_.getID());
343
344 auto const f_bar_variables =
345 computeFBarVariables(/*compute_detF0_only?*/ false, u);
346
347 for (unsigned ip = 0; ip < n_integration_points; ip++)
348 {
349 auto const& w = _ip_data[ip].integration_weight;
350 auto const& N = _ip_data[ip].N_u;
351 auto const& dNdx = _ip_data[ip].dNdx_u;
352
354 std::nullopt, this->element_.getID(),
358 this->element_, N))};
359
360 auto const x_coord = x_position.getCoordinates().value()[0];
361
362 // For the 2D case the 33-component is needed (and the four entries
363 // of the non-symmetric matrix); In 3d there are nine entries.
365 (DisplacementDim == 2 ? 1 : 0),
370
371 GradientVectorType const grad_u = G * u;
372
377
382
383 double const detF0 =
385 double const alpha = computeOutputStrainData(
386 detF0, B_0 + 0.5 * B_N, grad_u, u, this->output_data_[ip]);
387
388 auto const CD = updateConstitutiveRelations(
390 medium, this->current_states_[ip], this->prev_states_[ip],
391 this->material_states_[ip], this->output_data_[ip]);
392
393 BMatrixType B = B_0 + B_N;
394
395 auto const& sigma =
397 .sigma;
399 auto const& C =
403
404 auto const sigma_geom =
407
408 if (!f_bar_variables)
409 {
410 local_b.noalias() -=
411 (B.transpose() * sigma - N_u_op(N).transpose() * b) * w;
412
413 local_Jac.noalias() += G.transpose() * sigma_geom * G * w;
414
415 local_Jac.noalias() += B.transpose() * C * B * w;
416 }
417 else // with F bar
418 {
419 double const alpha_p2 = alpha * alpha;
420 local_Jac.noalias() +=
421 alpha_p2 * G.transpose() * sigma_geom * G * w;
422
424
425 auto const& F =
427 this->output_data_[ip])
434
435 auto const q0_q = q0 - q;
436
437 auto const S_q = sigma * (q0_q).transpose();
438
439 local_Jac.noalias() +=
440 2.0 * alpha_p2 * S_q.transpose() * B * w / DisplacementDim;
441
442 auto const twoEplsI =
444 alpha_p2,
446 this->output_data_[ip])
447 .eps,
449
450 // B bar:
451 BMatrixType const Bbar =
452 alpha_p2 *
453 (B + twoEplsI * (q0_q).transpose() / DisplacementDim);
454
455 local_Jac.noalias() +=
456 2.0 * Bbar.transpose() * S_q * w / DisplacementDim;
457
458 local_Jac.noalias() += Bbar.transpose() * C * Bbar * w;
459
460 local_Jac.noalias() -=
461 alpha_p2 * sigma.dot(twoEplsI) *
462 (q0 * q0.transpose() - q * q.transpose()) * w /
464
465 local_b.noalias() -=
466 (Bbar.transpose() * sigma - N_u_op(N).transpose() * b) * w;
467 }
468 }
469 }
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
std::optional< std::tuple< double, VectorTypeForFbar > > computeFBarVariables(bool const compute_detF0_only, Eigen::VectorXd const &local_x) const
typename GMatricesType::GradientVectorType GradientVectorType
double computeOutputStrainData(double const detF0, BMatrixType const &B, GradientVectorType const &grad_u, Eigen::Ref< Eigen::VectorXd const > const &u, typename ConstitutiveRelations::OutputData< DisplacementDim > &output_data) const
static constexpr auto localDOF(std::vector< double > const &x)
typename BMatricesType::NodalForceVectorType VectorTypeForFbar
ConstitutiveRelations::ConstitutiveData< DisplacementDim > updateConstitutiveRelations(double const alpha, GradientMatrixType const &G, 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
typename GMatricesType::GradientMatrixType GradientMatrixType
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)
Eigen::Matrix< double, MPL::tensorSize(DisplacementDim), MPL::tensorSize(DisplacementDim)> computeSigmaGeom(Eigen::Matrix3d const &sigma_tensor)
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > compute2EPlusI(double const alpha_p2, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps_bar, bool const is_axially_symmetric)
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 _ip_data, ProcessLib::NonLinearFbar::compute2EPlusI(), ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::NonLinearBMatrix::computeBMatrix(), computeFBarVariables(), ProcessLib::Deformation::computeGMatrix(), computeOutputStrainData(), ProcessLib::NonLinearFbar::computeQVector(), ProcessLib::LargeDeformation::computeSigmaGeom(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::element_, ParameterLib::SpatialPosition::getCoordinates(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, MathLib::KelvinVector::kelvinVectorToTensor(), localDOF(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, N_u_op, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::output_data_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, and updateConstitutiveRelations().

◆ computeFBarVariables()

template<typename ShapeFunction, int DisplacementDim>
std::optional< std::tuple< double, VectorTypeForFbar > > ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::computeFBarVariables ( bool const compute_detF0_only,
Eigen::VectorXd const & local_x ) const
inline

Definition at line 292 of file LargeDeformationFEM.h.

References _ip_data, ProcessLib::NonLinearFbar::computeFBarInitialVariablesAverage(), ProcessLib::NonLinearFbar::computeFBarInitialVariablesElementCenter(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::element_, ProcessLib::NonLinearFbar::ELEMENT_AVERAGE, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, ProcessLib::NonLinearFbar::NONE, and ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::process_data_.

Referenced by assembleWithJacobian(), and postTimestepConcrete().

◆ computeOutputStrainData()

template<typename ShapeFunction, int DisplacementDim>
double ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::computeOutputStrainData ( double const detF0,
BMatrixType const & B,
GradientVectorType const & grad_u,
Eigen::Ref< Eigen::VectorXd const > const & u,
typename ConstitutiveRelations::OutputData< DisplacementDim > & output_data ) const
inline

Definition at line 185 of file LargeDeformationFEM.h.

193 {
194 // Note: Here B=B_{linear}+B_{nonlinear}/2 For Green-Lagrange strain.
198
199 eps_data.eps = B * u;
200 deformation_gradient_data.deformation_gradient =
202 deformation_gradient_data.volume_ratio =
204 deformation_gradient_data.deformation_gradient);
205
206 if (this->process_data_.bar_det_f_type ==
208 {
209 return 1.0;
210 }
211
212 double const detF = deformation_gradient_data.volume_ratio;
213 double const J_ratio = detF0 / detF;
214
215 if (J_ratio < .0)
216 {
217 OGS_FATAL(
218 "det(F0)/det(F) is negative with det(F0) = {:g}, and det(F) = "
219 "{:g} ",
220 detF0, detF);
221 }
222
223 double const alpha =
225
226 deformation_gradient_data.deformation_gradient
228 deformation_gradient_data.volume_ratio *=
230
231 double const alpha_p2 = alpha * alpha;
232 eps_data.eps = alpha_p2 * eps_data.eps +
233 0.5 * (alpha_p2 - 1) *
236 return alpha;
237 }
double determinant(Eigen::MatrixBase< Derived > const &tensor)
Computes determinant of a vectorized tensor.

References MathLib::VectorizedTensor::determinant(), MathLib::VectorizedTensor::identity(), ProcessLib::NonLinearFbar::identityForF(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, ProcessLib::NonLinearFbar::NONE, OGS_FATAL, and ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::process_data_.

Referenced by assembleWithJacobian(), and postTimestepConcrete().

◆ 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 544 of file LargeDeformationFEM.h.

546 {
547 auto const& N = _secondary_data.N[integration_point];
548
549 // assumes N is stored contiguously in memory
550 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
551 }

References _secondary_data.

◆ 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 =
147 this->integration_method_.getNumberOfPoints();
148 for (unsigned ip = 0; ip < n_integration_points; ip++)
149 {
150 auto const& ip_data = _ip_data[ip];
151
153 std::nullopt, this->element_.getID(),
157 this->element_, ip_data.N_u))};
158
160 if (this->process_data_.initial_stress != nullptr)
161 {
163 .sigma.noalias() =
165 DisplacementDim>((*this->process_data_.initial_stress)(
167 double>::quiet_NaN() /* time independent */,
168 x_position));
169 }
170
171 double const t = 0; // TODO (naumov) pass t from top
172 auto& material_state = this->material_states_[ip];
173 this->solid_material_.initializeInternalStateVariables(
174 t, x_position, *material_state.material_state_variables);
175
179
180 material_state.pushBackState();
181 this->prev_states_[ip] = this->current_states_[ip];
182 }
183 }
ConstitutiveRelations::SolidConstitutiveRelation< DisplacementDim > const & solid_material_

References _ip_data, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::element_, ProcessLib::LargeDeformation::ConstitutiveRelations::ConstitutiveSetting< DisplacementDim >::init(), 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>
constexpr auto ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::localDOF ( std::vector< double > const & x)
inlinestaticconstexprprivate

Definition at line 554 of file LargeDeformationFEM.h.

References NumLib::localDOF().

Referenced by assembleWithJacobian().

◆ 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 471 of file LargeDeformationFEM.h.

475 {
476 unsigned const n_integration_points =
477 this->integration_method_.getNumberOfPoints();
478
481
482 auto& medium =
483 *this->process_data_.media_map.getMedium(this->element_.getID());
484
485 auto const f_bar_variables =
486 computeFBarVariables(/*compute_detF0_only?*/ true, local_x);
487
488 for (unsigned ip = 0; ip < n_integration_points; ip++)
489 {
490 auto const& N = _ip_data[ip].N_u;
491 auto const& dNdx = _ip_data[ip].dNdx_u;
492
494 std::nullopt, this->element_.getID(),
498 this->element_, N))};
499
500 auto const x_coord = x_position.getCoordinates().value()[0];
501
502 // For the 2D case the 33-component is needed (and the four entries
503 // of the non-symmetric matrix); In 3d there are nine entries.
505 (DisplacementDim == 2 ? 1 : 0),
510
512
513 // B for Green-Lagrange strain.
518
519 B.noalias() += 0.5 * NonLinearBMatrix::computeBMatrix<
522 dNdx, N, grad_u, x_coord,
524
525 double const detF0 =
527 double const alpha = computeOutputStrainData(
528 detF0, B, grad_u, local_x, this->output_data_[ip]);
529
532 medium, this->current_states_[ip], this->prev_states_[ip],
533 this->material_states_[ip], this->output_data_[ip]);
534
535 this->material_states_[ip].pushBackState();
536 }
537
538 for (unsigned ip = 0; ip < n_integration_points; ip++)
539 {
540 this->prev_states_[ip] = this->current_states_[ip];
541 }
542 }

References _ip_data, ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::NonLinearBMatrix::computeBMatrix(), computeFBarVariables(), ProcessLib::Deformation::computeGMatrix(), computeOutputStrainData(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::element_, ParameterLib::SpatialPosition::getCoordinates(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::output_data_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, and updateConstitutiveRelations().

◆ updateConstitutiveRelations()

template<typename ShapeFunction, int DisplacementDim>
ConstitutiveRelations::ConstitutiveData< DisplacementDim > ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::updateConstitutiveRelations ( double const alpha,
GradientMatrixType const & G,
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

Definition at line 240 of file LargeDeformationFEM.h.

255 {
256 double const T_ref =
257 this->process_data_.reference_temperature
258 ? (*this->process_data_.reference_temperature)(t, x_position)[0]
260
261 auto models =
263 this->process_data_, this->solid_material_);
265 tmp;
267
268 CS.eval(
269 models, t, dt, x_position, //
270 medium, //
271 T_ref, //
273 alpha * (G * u_prev +
276
277 return CD;
278 }
ConstitutiveModels< DisplacementDim > createConstitutiveModels(LDProcessData const &process_data, SolidConstitutiveRelation< DisplacementDim > const &solid_material)

References ProcessLib::LargeDeformation::ConstitutiveRelations::createConstitutiveModels(), ProcessLib::LargeDeformation::ConstitutiveRelations::ConstitutiveSetting< DisplacementDim >::eval(), MathLib::VectorizedTensor::identity(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, and ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::solid_material_.

Referenced by assembleWithJacobian(), and postTimestepConcrete().

Member Data Documentation

◆ _ip_data

template<typename ShapeFunction, int DisplacementDim>
std::vector<IpData, Eigen::aligned_allocator<IpData> > ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data
private

◆ _secondary_data

template<typename ShapeFunction, int DisplacementDim>
SecondaryData<typename ShapeMatrices::ShapeType> ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data
private

Definition at line 563 of file LargeDeformationFEM.h.

Referenced by LargeDeformationLocalAssembler(), and getShapeMatrix().

◆ 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 565 of file LargeDeformationFEM.h.

◆ N_u_op

template<typename ShapeFunction, int DisplacementDim>
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 assembleWithJacobian().


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