OGS
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim > Class Template Reference

Detailed Description

template<typename ShapeFunctionDisplacement, typename ShapeFunctionPressure, int GlobalDim>
class ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >

Definition at line 29 of file HydroMechanicsLocalAssemblerMatrixNearFracture.h.

#include <HydroMechanicsLocalAssemblerMatrixNearFracture.h>

Inheritance diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >:
[legend]

Public Member Functions

 HydroMechanicsLocalAssemblerMatrixNearFracture (HydroMechanicsLocalAssemblerMatrixNearFracture const &)=delete
 
 HydroMechanicsLocalAssemblerMatrixNearFracture (HydroMechanicsLocalAssemblerMatrixNearFracture &&)=delete
 
 HydroMechanicsLocalAssemblerMatrixNearFracture (MeshLib::Element const &e, std::size_t const n_variables, std::size_t const local_matrix_size, std::vector< unsigned > const &dofIndex_to_localIndex, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, HydroMechanicsProcessData< GlobalDim > &process_data)
 
- Public Member Functions inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >
 HydroMechanicsLocalAssemblerMatrix (HydroMechanicsLocalAssemblerMatrix const &)=delete
 
 HydroMechanicsLocalAssemblerMatrix (HydroMechanicsLocalAssemblerMatrix &&)=delete
 
 HydroMechanicsLocalAssemblerMatrix (MeshLib::Element const &e, std::size_t const n_variables, std::size_t const local_matrix_size, std::vector< unsigned > const &dofIndex_to_localIndex, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, HydroMechanicsProcessData< GlobalDim > &process_data)
 
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 override
 
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 override
 
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 override
 
std::vector< double > const & getIntPtFractureVelocity (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtFractureStress (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtFractureAperture (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtFracturePermeability (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) 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::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)
 
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
 
- 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 ~ExtrapolatableElement ()=default
 

Private Types

using Base
 
using BMatricesType
 
using IntegrationPointDataType
 
using ShapeMatricesTypeDisplacement
 

Private Member Functions

void assembleWithJacobianConcrete (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, Eigen::VectorXd &local_b, Eigen::MatrixXd &local_J) override
 
void preTimestepConcrete (std::vector< double > const &, double const, double const) override
 
void postTimestepConcreteWithVector (double const t, double const dt, Eigen::VectorXd const &local_x) override
 

Private Attributes

std::vector< FractureProperty * > _fracture_props
 
std::vector< JunctionProperty * > _junction_props
 
std::unordered_map< int, int > _fracID_to_local
 
Eigen::Vector3d _e_center_coords
 
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
 
HydroMechanicsProcessData< GlobalDim > & _process_data
 

Static Private Attributes

static const int displacement_jump_index
 
static const int displacement_index
 
static const int displacement_size
 
static const int kelvin_vector_size
 
static const int pressure_index
 
static const int pressure_size
 

Additional Inherited Members

- Protected Types inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >
using ShapeMatricesTypeDisplacement
 
using BMatricesType
 
using ShapeMatricesTypePressure
 
using IntegrationPointDataType
 
using GlobalDimVector = Eigen::Matrix<double, GlobalDim, 1>
 
- Protected Member Functions inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >
void assembleBlockMatricesWithJacobian (double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &p_prev, Eigen::Ref< const Eigen::VectorXd > const &u, Eigen::Ref< const Eigen::VectorXd > const &u_prev, Eigen::Ref< Eigen::VectorXd > rhs_p, Eigen::Ref< Eigen::VectorXd > rhs_u, Eigen::Ref< Eigen::MatrixXd > J_pp, Eigen::Ref< Eigen::MatrixXd > J_pu, Eigen::Ref< Eigen::MatrixXd > J_uu, Eigen::Ref< Eigen::MatrixXd > J_up)
 
void postTimestepConcreteWithBlockVectors (double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &u)
 
void setPressureOfInactiveNodes (double const t, Eigen::Ref< Eigen::VectorXd > p)
 
- Protected Attributes inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >
HydroMechanicsProcessData< GlobalDim > & _process_data
 
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
 
SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType > _secondary_data
 
- Protected Attributes inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface
MeshLib::Element const & _element
 
bool const _is_axially_symmetric
 
- Static Protected Attributes inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >
static const int pressure_index = 0
 
static const int pressure_size = ShapeFunctionPressure::NPOINTS
 
static const int displacement_index = ShapeFunctionPressure::NPOINTS
 
static const int displacement_size
 
static const int kelvin_vector_size
 

Member Typedef Documentation

◆ Base

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::Base
private
Initial value:
HydroMechanicsLocalAssemblerMatrix<ShapeFunctionDisplacement,
ShapeFunctionPressure,
GlobalDim>
HydroMechanicsLocalAssemblerMatrix(HydroMechanicsLocalAssemblerMatrix const &)=delete

Definition at line 35 of file HydroMechanicsLocalAssemblerMatrixNearFracture.h.

◆ BMatricesType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::BMatricesType
private

Definition at line 159 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ IntegrationPointDataType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::IntegrationPointDataType
private

Definition at line 166 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ ShapeMatricesTypeDisplacement

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::ShapeMatricesTypeDisplacement
private

Definition at line 157 of file HydroMechanicsLocalAssemblerMatrix.h.

Constructor & Destructor Documentation

◆ HydroMechanicsLocalAssemblerMatrixNearFracture() [1/3]

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::HydroMechanicsLocalAssemblerMatrixNearFracture ( HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim > const & )
delete

◆ HydroMechanicsLocalAssemblerMatrixNearFracture() [2/3]

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::HydroMechanicsLocalAssemblerMatrixNearFracture ( HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim > && )
delete

◆ HydroMechanicsLocalAssemblerMatrixNearFracture() [3/3]

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::HydroMechanicsLocalAssemblerMatrixNearFracture ( MeshLib::Element const & e,
std::size_t const n_variables,
std::size_t const local_matrix_size,
std::vector< unsigned > const & dofIndex_to_localIndex,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
HydroMechanicsProcessData< GlobalDim > & process_data )

Definition at line 25 of file HydroMechanicsLocalAssemblerMatrixNearFracture-impl.h.

34 : HydroMechanicsLocalAssemblerMatrix<ShapeFunctionDisplacement,
35 ShapeFunctionPressure, GlobalDim>(
36 e, n_variables, local_matrix_size, dofIndex_to_localIndex,
37 integration_method, is_axially_symmetric, process_data),
39{
40 // currently not supporting multiple fractures
41 _fracture_props.push_back(process_data.fracture_property.get());
42 _fracID_to_local.insert({0, 0});
43}
const double * data() const
Definition Point3d.h:59
MathLib::Point3d getCenterOfGravity(Element const &element)
Calculates the center of gravity for the mesh element.
Definition Element.cpp:124

References ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_fracID_to_local, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_fracture_props, and ProcessLib::LIE::HydroMechanics::HydroMechanicsProcessData< GlobalDim >::fracture_property.

Member Function Documentation

◆ assembleWithJacobianConcrete()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::assembleWithJacobianConcrete ( double const t,
double const dt,
Eigen::VectorXd const & local_x,
Eigen::VectorXd const & local_x_prev,
Eigen::VectorXd & local_b,
Eigen::MatrixXd & local_J )
overrideprivatevirtual

Reimplemented from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >.

Definition at line 49 of file HydroMechanicsLocalAssemblerMatrixNearFracture-impl.h.

55{
56 auto p = const_cast<Eigen::VectorXd&>(local_x).segment(pressure_index,
58 auto p_prev = const_cast<Eigen::VectorXd&>(local_x_prev)
60 if (_process_data.deactivate_matrix_in_flow)
61 {
63 }
64 auto const u = local_x.segment(displacement_index, displacement_size);
65 auto const u_prev =
66 local_x_prev.segment(displacement_index, displacement_size);
67
68 auto rhs_p = local_b.segment(pressure_index, pressure_size);
69 auto rhs_u = local_b.segment(displacement_index, displacement_size);
70
71 auto J_pp = local_J.block(pressure_index, pressure_index, pressure_size,
73 auto J_pu = local_J.block(pressure_index, displacement_index, pressure_size,
75 auto J_up = local_J.block(displacement_index, pressure_index,
77 auto J_uu = local_J.block(displacement_index, displacement_index,
79
80 // levelset value of the element
81 // remark: this assumes the levelset function is uniform within an element
82 std::vector<double> levelsets = uGlobalEnrichments(
84 double const ele_levelset = levelsets[0]; // single fracture
85
86 if (ele_levelset == 0)
87 {
88 // no DoF exists for displacement jumps. do the normal assembly
90 t, dt, p, p_prev, u, u_prev, rhs_p, rhs_u, J_pp, J_pu, J_uu, J_up);
91 return;
92 }
93
94 // Displacement jumps should be taken into account
95
96 // compute true displacements
97 auto const g = local_x.segment(displacement_jump_index, displacement_size);
98 auto const g_prev =
99 local_x_prev.segment(displacement_jump_index, displacement_size);
100 Eigen::VectorXd const total_u = u + ele_levelset * g;
101 Eigen::VectorXd const total_u_prev = u_prev + ele_levelset * g_prev;
102
103 // evaluate residuals and Jacobians for pressure and displacements
104 Base::assembleBlockMatricesWithJacobian(t, dt, p, p_prev, total_u,
105 total_u_prev, rhs_p, rhs_u, J_pp,
106 J_pu, J_uu, J_up);
107
108 // compute residuals and Jacobians for displacement jumps
109 auto rhs_g = local_b.segment(displacement_jump_index, displacement_size);
110 auto J_pg = local_J.block(pressure_index, displacement_jump_index,
112 auto J_ug = local_J.block(displacement_index, displacement_jump_index,
114 auto J_gp = local_J.block(displacement_jump_index, pressure_index,
116 auto J_gu = local_J.block(displacement_jump_index, displacement_index,
118 auto J_gg = local_J.block(displacement_jump_index, displacement_jump_index,
120
121 rhs_g = ele_levelset * rhs_u;
122 J_pg = ele_levelset * J_pu;
123 J_ug = ele_levelset * J_uu;
124 J_gp = ele_levelset * J_up;
125 J_gu = ele_levelset * J_uu;
126 J_gg = ele_levelset * ele_levelset * J_uu;
127}
void assembleBlockMatricesWithJacobian(double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &p_prev, Eigen::Ref< const Eigen::VectorXd > const &u, Eigen::Ref< const Eigen::VectorXd > const &u_prev, Eigen::Ref< Eigen::VectorXd > rhs_p, Eigen::Ref< Eigen::VectorXd > rhs_u, Eigen::Ref< Eigen::MatrixXd > J_pp, Eigen::Ref< Eigen::MatrixXd > J_pu, Eigen::Ref< Eigen::MatrixXd > J_uu, Eigen::Ref< Eigen::MatrixXd > J_up)
std::vector< double > uGlobalEnrichments(std::vector< FractureProperty * > const &frac_props, std::vector< JunctionProperty * > const &junction_props, std::unordered_map< int, int > const &fracID_to_local, Eigen::Vector3d const &x)

References ProcessLib::LIE::uGlobalEnrichments().

◆ postTimestepConcreteWithVector()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::postTimestepConcreteWithVector ( double const t,
double const dt,
Eigen::VectorXd const & local_x )
overrideprivatevirtual

Reimplemented from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >.

Definition at line 133 of file HydroMechanicsLocalAssemblerMatrixNearFracture-impl.h.

135{
136 auto p = const_cast<Eigen::VectorXd&>(local_x).segment(pressure_index,
138 if (_process_data.deactivate_matrix_in_flow)
139 {
141 }
142 auto u = local_x.segment(displacement_index, displacement_size);
143
144 // levelset value of the element
145 // remark: this assumes the levelset function is uniform within an element
146 std::vector<double> levelsets = uGlobalEnrichments(
148 double const ele_levelset = levelsets[0]; // single fracture
149
150 if (ele_levelset == 0)
151 {
152 // no DoF exists for displacement jumps. do the normal assembly
154 return;
155 }
156
157 // Displacement jumps should be taken into account
158
159 // compute true displacements
160 auto const g = local_x.segment(displacement_jump_index, displacement_size);
161 Eigen::VectorXd const total_u = u + ele_levelset * g;
162
163 // evaluate residuals and Jacobians for pressure and displacements
165}
void postTimestepConcreteWithBlockVectors(double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &u)

References ProcessLib::LIE::uGlobalEnrichments().

◆ preTimestepConcrete()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::preTimestepConcrete ( std::vector< double > const & ,
double const ,
double const  )
inlineoverrideprivatevirtual

Member Data Documentation

◆ _e_center_coords

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
Eigen::Vector3d ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_e_center_coords
private

◆ _fracID_to_local

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
std::unordered_map<int, int> ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_fracID_to_local
private

◆ _fracture_props

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
std::vector<FractureProperty*> ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_fracture_props
private

◆ _ip_data

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
std::vector<IntegrationPointDataType, Eigen::aligned_allocator<IntegrationPointDataType> > ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_ip_data
private

◆ _junction_props

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
std::vector<JunctionProperty*> ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_junction_props
private

◆ _process_data

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
HydroMechanicsProcessData<GlobalDim>& ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_process_data
private

Definition at line 173 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ displacement_index

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::displacement_index
staticprivate

Definition at line 181 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ displacement_jump_index

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::displacement_jump_index
staticprivate

◆ displacement_size

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::displacement_size
staticprivate

Definition at line 182 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ kelvin_vector_size

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::kelvin_vector_size
staticprivate

Definition at line 184 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ pressure_index

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::pressure_index
staticprivate

Definition at line 179 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ pressure_size

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::pressure_size
staticprivate

Definition at line 180 of file HydroMechanicsLocalAssemblerMatrix.h.


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