![]() |
OGS
|
Definition at line 81 of file LargeDeformationFEM.h.
#include <LargeDeformationFEM.h>
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 > ¤t_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. | |
![]() | |
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 |
![]() | |
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 | ~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 auto | getReflectionDataForOutput () |
![]() | |
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_ |
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim> |
Definition at line 92 of file LargeDeformationFEM.h.
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::BMatrixType = typename BMatricesType::BMatrixType |
Definition at line 94 of file LargeDeformationFEM.h.
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GlobalDimNodalMatrixType |
Definition at line 88 of file LargeDeformationFEM.h.
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GMatricesType = GMatrixPolicyType<ShapeFunction, DisplacementDim> |
Definition at line 100 of file LargeDeformationFEM.h.
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GradientMatrixType = typename GMatricesType::GradientMatrixType |
Definition at line 102 of file LargeDeformationFEM.h.
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GradientVectorType = typename GMatricesType::GradientVectorType |
Definition at line 101 of file LargeDeformationFEM.h.
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::IpData |
Definition at line 106 of file LargeDeformationFEM.h.
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalDisplacementVectorType |
Definition at line 97 of file LargeDeformationFEM.h.
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalForceVectorType = typename BMatricesType::NodalForceVectorType |
Definition at line 96 of file LargeDeformationFEM.h.
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType |
Definition at line 87 of file LargeDeformationFEM.h.
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType |
Definition at line 90 of file LargeDeformationFEM.h.
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices |
Definition at line 91 of file LargeDeformationFEM.h.
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatricesType |
Definition at line 85 of file LargeDeformationFEM.h.
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::StiffnessMatrixType = typename BMatricesType::StiffnessMatrixType |
Definition at line 95 of file LargeDeformationFEM.h.
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::VectorTypeForFbar = typename BMatricesType::NodalForceVectorType |
Definition at line 104 of file LargeDeformationFEM.h.
|
delete |
|
delete |
|
inline |
Definition at line 116 of file LargeDeformationFEM.h.
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::LargeDeformation::SecondaryData< ShapeMatrixType >::N.
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 278 of file LargeDeformationFEM.h.
References OGS_FATAL.
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 317 of file LargeDeformationFEM.h.
References ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::NonLinearFbar::compute2EPlusI(), ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::NonLinearBMatrix::computeBMatrix(), ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::computeFBarVariables(), ProcessLib::Deformation::computeGMatrix(), ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::computeOutputStrainData(), ProcessLib::NonLinearFbar::computeQVector(), ProcessLib::LargeDeformation::computeSigmaGeom(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), 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(), and ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::updateConstitutiveRelations().
|
inline |
Definition at line 290 of file LargeDeformationFEM.h.
References ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_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 ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian(), and ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete().
|
inline |
Definition at line 187 of file LargeDeformationFEM.h.
References ProcessLib::LargeDeformation::ConstitutiveRelations::OutputData< DisplacementDim >::deformation_gradient_data, MathLib::VectorizedTensor::determinant(), ProcessLib::LargeDeformation::ConstitutiveRelations::OutputData< DisplacementDim >::eps_data, 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 ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian(), and ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete().
|
inlineoverridevirtual |
Provides the shape matrix at the given integration point.
Implements NumLib::ExtrapolatableElement.
Definition at line 530 of file LargeDeformationFEM.h.
References ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data, and ProcessLib::LargeDeformation::SecondaryData< ShapeMatrixType >::N.
|
inlineoverridevirtual |
Set initial stress from parameter.
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 151 of file LargeDeformationFEM.h.
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().
|
inlinestaticconstexprprivate |
Definition at line 540 of file LargeDeformationFEM.h.
References NumLib::localDOF().
Referenced by ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian().
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 459 of file LargeDeformationFEM.h.
References ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::NonLinearBMatrix::computeBMatrix(), ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::computeFBarVariables(), ProcessLib::Deformation::computeGMatrix(), ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::computeOutputStrainData(), 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_, 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(), and ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::updateConstitutiveRelations().
|
inline |
Definition at line 239 of file LargeDeformationFEM.h.
References ProcessLib::LargeDeformation::ConstitutiveRelations::OutputData< DisplacementDim >::deformation_gradient_data, ProcessLib::LargeDeformation::ConstitutiveRelations::ConstitutiveSetting< DisplacementDim >::eval(), MathLib::VectorizedTensor::identity(), 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().
|
private |
Definition at line 547 of file LargeDeformationFEM.h.
Referenced by ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::LargeDeformationLocalAssembler(), ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian(), ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::computeFBarVariables(), ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete(), and ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete().
|
private |
Definition at line 549 of file LargeDeformationFEM.h.
Referenced by ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::LargeDeformationLocalAssembler(), and ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::getShapeMatrix().
|
staticprivate |
Definition at line 551 of file LargeDeformationFEM.h.
|
staticconstexpr |
Definition at line 109 of file LargeDeformationFEM.h.
Referenced by ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian().