OGS
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface Class Referenceabstract

Detailed Description

Definition at line 35 of file HydroMechanicsLocalAssemblerInterface.h.

#include <HydroMechanicsLocalAssemblerInterface.h>

Inheritance diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface:
[legend]
Collaboration diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface:
[legend]

Public Member Functions

 HydroMechanicsLocalAssemblerInterface (MeshLib::Element const &element, bool const is_axially_symmetric, std::size_t n_local_size, std::vector< unsigned > dofIndex_to_localIndex)
 
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 &, const double t, double const dt, int const) override
 
virtual std::vector< double > const & getIntPtSigma (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtEpsilon (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtDarcyVelocity (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtFractureVelocity (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtFractureStress (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtFractureAperture (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
 
virtual std::vector< double > const & getIntPtFracturePermeability (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
 
- 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 Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const =0
 Provides the shape matrix at the given integration point.
 
virtual ~ExtrapolatableElement ()=default
 

Protected Member Functions

virtual void assembleWithJacobianConcrete (double const t, double const dt, Eigen::VectorXd const &local_u, Eigen::VectorXd const &local_u_prev, Eigen::VectorXd &local_b, Eigen::MatrixXd &local_J)=0
 
virtual void postTimestepConcreteWithVector (double const t, double const dt, Eigen::VectorXd const &local_u)=0
 

Protected Attributes

MeshLib::Element const & _element
 
bool const _is_axially_symmetric
 

Private Attributes

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

Constructor & Destructor Documentation

◆ HydroMechanicsLocalAssemblerInterface()

ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::HydroMechanicsLocalAssemblerInterface ( MeshLib::Element const & element,
bool const is_axially_symmetric,
std::size_t n_local_size,
std::vector< unsigned > dofIndex_to_localIndex )
inline

Definition at line 40 of file HydroMechanicsLocalAssemblerInterface.h.

45 : _element(element),
46 _is_axially_symmetric(is_axially_symmetric),
47 _dofIndex_to_localIndex(std::move(dofIndex_to_localIndex))
48 {
49 _local_u.resize(n_local_size);
50 _local_u_prev.resize(n_local_size);
51 _local_b.resize(_local_u.size());
52 _local_J.resize(_local_u.size(), _local_u.size());
53 }

References _local_b, _local_J, _local_u, and _local_u_prev.

Member Function Documentation

◆ assemble()

void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::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 55 of file HydroMechanicsLocalAssemblerInterface.h.

61 {
63 "HydroMechanicsLocalAssembler: assembly without jacobian is not "
64 "implemented.");
65 }
#define OGS_FATAL(...)
Definition Error.h:26

References OGS_FATAL.

◆ assembleWithJacobian()

void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::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 67 of file HydroMechanicsLocalAssemblerInterface.h.

74 {
75 auto const local_dof_size = local_x_.size();
76
77 _local_u.setZero();
78 for (unsigned i = 0; i < local_dof_size; i++)
79 {
80 _local_u[_dofIndex_to_localIndex[i]] = local_x_[i];
81 }
82 _local_u_prev.setZero();
83 for (unsigned i = 0; i < local_dof_size; i++)
84 {
85 _local_u_prev[_dofIndex_to_localIndex[i]] = local_x_prev_[i];
86 }
87 _local_b.setZero();
88 _local_J.setZero();
89
91 _local_J);
92
93 local_b_data.resize(local_dof_size);
94 for (unsigned i = 0; i < local_dof_size; i++)
95 {
96 local_b_data[i] = _local_b[_dofIndex_to_localIndex[i]];
97 }
98
99 local_Jac_data.resize(local_dof_size * local_dof_size);
100 for (unsigned i = 0; i < local_dof_size; i++)
101 {
102 for (unsigned j = 0; j < local_dof_size; j++)
103 {
104 local_Jac_data[i * local_dof_size + j] = _local_J(
106 }
107 }
108 }
virtual void assembleWithJacobianConcrete(double const t, double const dt, Eigen::VectorXd const &local_u, Eigen::VectorXd const &local_u_prev, Eigen::VectorXd &local_b, Eigen::MatrixXd &local_J)=0

References _dofIndex_to_localIndex, _local_b, _local_J, _local_u, _local_u_prev, and assembleWithJacobianConcrete().

◆ assembleWithJacobianConcrete()

virtual void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::assembleWithJacobianConcrete ( double const t,
double const dt,
Eigen::VectorXd const & local_u,
Eigen::VectorXd const & local_u_prev,
Eigen::VectorXd & local_b,
Eigen::MatrixXd & local_J )
protectedpure virtual

◆ getIntPtDarcyVelocity()

virtual std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::getIntPtDarcyVelocity ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
pure virtual

◆ getIntPtEpsilon()

virtual std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::getIntPtEpsilon ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
pure virtual

◆ getIntPtFractureAperture()

virtual std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::getIntPtFractureAperture ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
pure virtual

◆ getIntPtFracturePermeability()

virtual std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::getIntPtFracturePermeability ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
pure virtual

◆ getIntPtFractureStress()

virtual std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::getIntPtFractureStress ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
pure virtual

◆ getIntPtFractureVelocity()

virtual std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::getIntPtFractureVelocity ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
pure virtual

◆ getIntPtSigma()

virtual std::vector< double > const & ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::getIntPtSigma ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
pure virtual

◆ postTimestepConcrete()

void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::postTimestepConcrete ( Eigen::VectorXd const & local_x_,
Eigen::VectorXd const & ,
const double t,
double const dt,
int const  )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 110 of file HydroMechanicsLocalAssemblerInterface.h.

114 {
115 _local_u.setZero();
116 for (Eigen::Index i = 0; i < local_x_.rows(); i++)
117 {
118 _local_u[_dofIndex_to_localIndex[i]] = local_x_[i];
119 }
120
122 }
virtual void postTimestepConcreteWithVector(double const t, double const dt, Eigen::VectorXd const &local_u)=0

References _dofIndex_to_localIndex, _local_u, and postTimestepConcreteWithVector().

◆ postTimestepConcreteWithVector()

Member Data Documentation

◆ _dofIndex_to_localIndex

std::vector<unsigned> const ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_dofIndex_to_localIndex
private

◆ _element

MeshLib::Element const& ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_element
protected

Definition at line 175 of file HydroMechanicsLocalAssemblerInterface.h.

◆ _is_axially_symmetric

bool const ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_is_axially_symmetric
protected

Definition at line 176 of file HydroMechanicsLocalAssemblerInterface.h.

◆ _local_b

Eigen::VectorXd ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_local_b
private

◆ _local_J

Eigen::MatrixXd ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_local_J
private

◆ _local_u

Eigen::VectorXd ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_local_u
private

◆ _local_u_prev

Eigen::VectorXd ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface::_local_u_prev
private

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