OGS
ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface Class Referenceabstract

Detailed Description

Definition at line 27 of file SmallDeformationLocalAssemblerInterface.h.

#include <SmallDeformationLocalAssemblerInterface.h>

Inheritance diagram for ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface:
[legend]
Collaboration diagram for ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface:
[legend]

Public Member Functions

 SmallDeformationLocalAssemblerInterface ()
 
 SmallDeformationLocalAssemblerInterface (std::size_t n_local_size, std::vector< unsigned > dofIndex_to_localIndex)
 
void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x_, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
 
virtual void assembleWithJacobian (double const, double const, Eigen::VectorXd const &, Eigen::VectorXd &, Eigen::MatrixXd &)
 
void computeSecondaryVariableConcrete (double const t, double const, Eigen::VectorXd const &local_x, Eigen::VectorXd const &) override
 
virtual std::vector< double > const & getIntPtSigmaXX (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtSigmaYY (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtSigmaZZ (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtSigmaXY (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtSigmaXZ (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtSigmaYZ (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtEpsilonXX (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtEpsilonYY (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtEpsilonZZ (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtEpsilonXY (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtEpsilonXZ (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtEpsilonYZ (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
 
- Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
 
void setInitialConditions (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, bool const use_monolithic_scheme, 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 assemble (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
 
virtual void assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
 
virtual void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, 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 assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, 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_dot, 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, double const t, double const dt)
 
void postNonLinearSolver (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, 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. More...
 
- Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const =0
 Provides the shape matrix at the given integration point. More...
 
virtual ~ExtrapolatableElement ()=default
 

Protected Member Functions

virtual void computeSecondaryVariableConcreteWithVector (double const t, Eigen::VectorXd const &local_u)=0
 

Private Attributes

Eigen::VectorXd _local_u
 
Eigen::VectorXd _local_b
 
Eigen::MatrixXd _local_J
 
std::vector< unsigned > const _dofIndex_to_localIndex
 

Constructor & Destructor Documentation

◆ SmallDeformationLocalAssemblerInterface() [1/2]

ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::SmallDeformationLocalAssemblerInterface ( )
inline

◆ SmallDeformationLocalAssemblerInterface() [2/2]

ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::SmallDeformationLocalAssemblerInterface ( std::size_t  n_local_size,
std::vector< unsigned >  dofIndex_to_localIndex 
)
inline

Member Function Documentation

◆ assembleWithJacobian() [1/2]

void ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::assembleWithJacobian ( double const  t,
double const  dt,
std::vector< double > const &  local_x_,
std::vector< double > const &  ,
std::vector< double > &  ,
std::vector< double > &  ,
std::vector< double > &  local_b_data,
std::vector< double > &  local_Jac_data 
)
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Reimplemented in ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >.

Definition at line 42 of file SmallDeformationLocalAssemblerInterface.h.

49 {
50 auto const local_dof_size = local_x_.size();
51
52 _local_u.setZero();
53 for (unsigned i = 0; i < local_dof_size; i++)
54 {
55 _local_u[_dofIndex_to_localIndex[i]] = local_x_[i];
56 }
57 _local_b.setZero();
58 _local_J.setZero();
59
61
62 local_b_data.resize(local_dof_size);
63 for (unsigned i = 0; i < local_dof_size; i++)
64 {
65 local_b_data[i] = _local_b[_dofIndex_to_localIndex[i]];
66 }
67
68 local_Jac_data.resize(local_dof_size * local_dof_size);
69 for (unsigned i = 0; i < local_dof_size; i++)
70 {
71 for (unsigned j = 0; j < local_dof_size; j++)
72 {
73 local_Jac_data[i * local_dof_size + j] = _local_J(
75 }
76 }
77 }
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x_, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
static const double t

References _dofIndex_to_localIndex, _local_b, _local_J, _local_u, assembleWithJacobian(), and MathLib::t.

Referenced by assembleWithJacobian().

◆ assembleWithJacobian() [2/2]

virtual void ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::assembleWithJacobian ( double const  ,
double const  ,
Eigen::VectorXd const &  ,
Eigen::VectorXd &  ,
Eigen::MatrixXd &   
)
inlinevirtual

◆ computeSecondaryVariableConcrete()

void ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::computeSecondaryVariableConcrete ( double const  t,
double const  ,
Eigen::VectorXd const &  local_x,
Eigen::VectorXd const &   
)
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 89 of file SmallDeformationLocalAssemblerInterface.h.

92 {
93 if (!_dofIndex_to_localIndex.empty())
94 {
95 _local_u.setZero();
96 for (auto i = 0; i < local_x.rows(); i++)
97 {
98 _local_u[_dofIndex_to_localIndex[i]] = local_x[i];
99 }
100 }
101
103 }
virtual void computeSecondaryVariableConcreteWithVector(double const t, Eigen::VectorXd const &local_u)=0

References _dofIndex_to_localIndex, _local_u, computeSecondaryVariableConcreteWithVector(), and MathLib::t.

◆ computeSecondaryVariableConcreteWithVector()

◆ getIntPtEpsilonXX()

◆ getIntPtEpsilonXY()

◆ getIntPtEpsilonXZ()

◆ getIntPtEpsilonYY()

◆ getIntPtEpsilonYZ()

◆ getIntPtEpsilonZZ()

◆ getIntPtSigmaXX()

◆ getIntPtSigmaXY()

◆ getIntPtSigmaXZ()

◆ getIntPtSigmaYY()

◆ getIntPtSigmaYZ()

◆ getIntPtSigmaZZ()

Member Data Documentation

◆ _dofIndex_to_localIndex

std::vector<unsigned> const ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::_dofIndex_to_localIndex
private

◆ _local_b

Eigen::VectorXd ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::_local_b
private

◆ _local_J

Eigen::MatrixXd ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::_local_J
private

◆ _local_u

Eigen::VectorXd ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::_local_u
private

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