OGS
|
Definition at line 60 of file SmallDeformationFEM.h.
#include <SmallDeformationFEM.h>
Public Types | |
using | ShapeMatricesType |
using | NodalMatrixType = typename ShapeMatricesType::NodalMatrixType |
using | NodalVectorType = typename ShapeMatricesType::NodalVectorType |
using | ShapeMatrices = typename ShapeMatricesType::ShapeMatrices |
using | BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim> |
using | BMatrixType = typename BMatricesType::BMatrixType |
using | BBarMatrixType = typename BMatricesType::BBarMatrixType |
using | StiffnessMatrixType = typename BMatricesType::StiffnessMatrixType |
using | NodalForceVectorType = typename BMatricesType::NodalForceVectorType |
using | NodalDisplacementVectorType |
using | GMatricesType = GMatrixPolicyType<ShapeFunction, DisplacementDim> |
using | GradientVectorType = typename GMatricesType::GradientVectorType |
using | GradientMatrixType = typename GMatricesType::GradientMatrixType |
using | IpData |
Public Member Functions | |
SmallDeformationLocalAssembler (SmallDeformationLocalAssembler const &)=delete | |
SmallDeformationLocalAssembler (SmallDeformationLocalAssembler &&)=delete | |
SmallDeformationLocalAssembler (MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, SmallDeformationProcessData< DisplacementDim > &process_data) | |
void | initializeConcrete () override |
ConstitutiveRelations::ConstitutiveData< DisplacementDim > | updateConstitutiveRelations (BMatrixType const &B, Eigen::Ref< Eigen::VectorXd const > const &u, Eigen::Ref< Eigen::VectorXd const > const &u_prev, ParameterLib::SpatialPosition const &x_position, double const t, double const dt, typename ConstitutiveRelations::ConstitutiveSetting< DisplacementDim > &CS, MaterialPropertyLib::Medium const &medium, typename ConstitutiveRelations::StatefulData< DisplacementDim > ¤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< BBarMatrixType > | getDilatationalBBarMatrix () const |
void | assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override |
void | postTimestepConcrete (Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, double const t, double const dt, int const) override |
std::vector< double > const & | getMaterialForces (std::vector< double > const &local_x, std::vector< double > &nodal_values) override |
Eigen::Map< const Eigen::RowVectorXd > | getShapeMatrix (const unsigned integration_point) const override |
Provides the shape matrix at the given integration point. | |
Public Member Functions inherited from ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim > | |
SmallDeformationLocalAssemblerInterface (MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, SmallDeformationProcessData< DisplacementDim > &process_data) | |
std::size_t | setIPDataInitialConditions (std::string_view name, double const *values, int const integration_order) |
Returns number of read integration points. | |
void | computeSecondaryVariableConcrete (double const, double const, Eigen::VectorXd const &, Eigen::VectorXd const &) override |
unsigned | getNumberOfIntegrationPoints () const |
int | getMaterialID () const |
std::vector< double > | getMaterialStateVariableInternalState (std::function< std::span< double >(typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables &)> const &get_values_span, int const &n_components) const |
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & | getMaterialStateVariablesAt (unsigned integration_point) const |
Public Member Functions inherited from ProcessLib::LocalAssemblerInterface | |
virtual | ~LocalAssemblerInterface ()=default |
virtual void | setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, int const process_id) |
virtual void | initialize (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table) |
virtual void | preAssemble (double const, double const, std::vector< double > const &) |
virtual void | assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) |
virtual void | assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) |
virtual void | computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id) |
virtual void | preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t) |
virtual void | postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id) |
void | postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id) |
virtual Eigen::Vector3d | getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const |
virtual Eigen::Vector3d | getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double > > const &) const |
Fits to staggered scheme. | |
Public Member Functions inherited from ProcessLib::SmallDeformation::MaterialForcesInterface | |
virtual | ~MaterialForcesInterface ()=default |
Public Member Functions inherited from NumLib::ExtrapolatableElement | |
virtual | ~ExtrapolatableElement ()=default |
Static Public Attributes | |
static constexpr auto & | N_u_op |
Static Private Member Functions | |
static constexpr auto | localDOF (std::vector< double > const &x) |
Private Attributes | |
std::vector< IpData, Eigen::aligned_allocator< IpData > > | _ip_data |
SecondaryData< typename ShapeMatrices::ShapeType > | _secondary_data |
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::BBarMatrixType = typename BMatricesType::BBarMatrixType |
Definition at line 72 of file SmallDeformationFEM.h.
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim> |
Definition at line 69 of file SmallDeformationFEM.h.
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::BMatrixType = typename BMatricesType::BMatrixType |
Definition at line 71 of file SmallDeformationFEM.h.
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GMatricesType = GMatrixPolicyType<ShapeFunction, DisplacementDim> |
Definition at line 78 of file SmallDeformationFEM.h.
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GradientMatrixType = typename GMatricesType::GradientMatrixType |
Definition at line 80 of file SmallDeformationFEM.h.
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GradientVectorType = typename GMatricesType::GradientVectorType |
Definition at line 79 of file SmallDeformationFEM.h.
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::IpData |
Definition at line 81 of file SmallDeformationFEM.h.
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalDisplacementVectorType |
Definition at line 75 of file SmallDeformationFEM.h.
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalForceVectorType = typename BMatricesType::NodalForceVectorType |
Definition at line 74 of file SmallDeformationFEM.h.
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType |
Definition at line 66 of file SmallDeformationFEM.h.
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType |
Definition at line 67 of file SmallDeformationFEM.h.
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices |
Definition at line 68 of file SmallDeformationFEM.h.
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatricesType |
Definition at line 64 of file SmallDeformationFEM.h.
using ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::StiffnessMatrixType = typename BMatricesType::StiffnessMatrixType |
Definition at line 73 of file SmallDeformationFEM.h.
|
delete |
|
delete |
|
inline |
Definition at line 91 of file SmallDeformationFEM.h.
References ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data, NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, and ProcessLib::SmallDeformation::SecondaryData< ShapeMatrixType >::N.
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 203 of file SmallDeformationFEM.h.
References OGS_FATAL.
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 229 of file SmallDeformationFEM.h.
References ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::LinearBMatrix::computeBMatrixPossiblyWithBbar(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::constitutive_setting, MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::getDilatationalBBarMatrix(), MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateXCoordinate(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::localDOF(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::N_u_op, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::output_data_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), and ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::updateConstitutiveRelations().
|
inline |
Definition at line 215 of file SmallDeformationFEM.h.
References ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::LinearBMatrix::computeDilatationalBbar(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, and ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::process_data_.
Referenced by ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian(), and ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete().
|
inlineoverridevirtual |
Implements ProcessLib::SmallDeformation::MaterialForcesInterface.
Definition at line 339 of file SmallDeformationFEM.h.
References ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, ProcessLib::SmallDeformation::getMaterialForces(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, and ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::output_data_.
|
inlineoverridevirtual |
Provides the shape matrix at the given integration point.
Implements NumLib::ExtrapolatableElement.
Definition at line 353 of file SmallDeformationFEM.h.
References ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data, and ProcessLib::SmallDeformation::SecondaryData< ShapeMatrixType >::N.
|
inlineoverridevirtual |
Set initial stress from parameter.
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 126 of file SmallDeformationFEM.h.
References ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::solid_material_, and MathLib::KelvinVector::symmetricTensorToKelvinVector().
|
inlinestaticconstexprprivate |
Definition at line 363 of file SmallDeformationFEM.h.
References NumLib::localDOF().
Referenced by ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian().
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 290 of file SmallDeformationFEM.h.
References ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::LinearBMatrix::computeBMatrixPossiblyWithBbar(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::constitutive_setting, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::element_, ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::getDilatationalBBarMatrix(), MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateXCoordinate(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::output_data_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), and ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::updateConstitutiveRelations().
|
inline |
Definition at line 167 of file SmallDeformationFEM.h.
References ProcessLib::SmallDeformation::ConstitutiveRelations::ConstitutiveSetting< DisplacementDim >::eval(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, and ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::solid_material_.
Referenced by ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian(), and ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete().
|
private |
Definition at line 370 of file SmallDeformationFEM.h.
Referenced by ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::SmallDeformationLocalAssembler(), ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian(), ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::getDilatationalBBarMatrix(), ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::getMaterialForces(), ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete(), and ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete().
|
private |
Definition at line 372 of file SmallDeformationFEM.h.
Referenced by ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::SmallDeformationLocalAssembler(), and ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::getShapeMatrix().
|
staticconstexpr |
Definition at line 84 of file SmallDeformationFEM.h.
Referenced by ProcessLib::SmallDeformation::SmallDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian().