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

Detailed Description

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

Definition at line 31 of file HydroMechanicsLocalAssemblerMatrixNearFracture.h.

#include <HydroMechanicsLocalAssemblerMatrixNearFracture.h>

Inheritance diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, 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, bool const is_axially_symmetric, unsigned const integration_order, HydroMechanicsProcessData< GlobalDim > &process_data)
 
- Public Member Functions inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, 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, bool const is_axially_symmetric, unsigned const integration_order, HydroMechanicsProcessData< GlobalDim > &process_data)
 
- 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_xdot_, const double, const double, 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_, const double t, double const dt) override
 
- 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 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 assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, const double dxdot_dx, const double dx_dx, 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 std::vector< double > interpolateNodalValuesToIntegrationPoints (std::vector< double > const &)
 
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...
 

Private Types

using Base = HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, GlobalDim >
 
using BMatricesType = BMatrixPolicyType< ShapeFunctionDisplacement, GlobalDim >
 
using ShapeMatricesTypeDisplacement = ShapeMatrixPolicyType< ShapeFunctionDisplacement, GlobalDim >
 

Private Member Functions

void assembleWithJacobianConcrete (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_dot, 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, IntegrationMethod, GlobalDim >
using ShapeMatricesTypeDisplacement = ShapeMatrixPolicyType< ShapeFunctionDisplacement, GlobalDim >
 
using BMatricesType = BMatrixPolicyType< ShapeFunctionDisplacement, GlobalDim >
 
using ShapeMatricesTypePressure = ShapeMatrixPolicyType< ShapeFunctionPressure, GlobalDim >
 
using IntegrationPointDataType = IntegrationPointDataMatrix< BMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, GlobalDim, ShapeFunctionDisplacement::NPOINTS >
 
using GlobalDimVector = Eigen::Matrix< double, GlobalDim, 1 >
 
- Protected Member Functions inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, GlobalDim >
void assembleBlockMatricesWithJacobian (double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &p_dot, Eigen::Ref< const Eigen::VectorXd > const &u, Eigen::Ref< const Eigen::VectorXd > const &u_dot, 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)
 
void setPressureDotOfInactiveNodes (Eigen::Ref< Eigen::VectorXd > p_dot)
 
- Protected Attributes inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrix< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, GlobalDim >
HydroMechanicsProcessData< GlobalDim > & _process_data
 
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_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, IntegrationMethod, 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 , typename IntegrationMethod , int GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, GlobalDim >::Base = HydroMechanicsLocalAssemblerMatrix<ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, GlobalDim>
private

◆ BMatricesType

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

Definition at line 96 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ ShapeMatricesTypeDisplacement

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

Definition at line 94 of file HydroMechanicsLocalAssemblerMatrix.h.

Constructor & Destructor Documentation

◆ HydroMechanicsLocalAssemblerMatrixNearFracture() [1/3]

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

◆ HydroMechanicsLocalAssemblerMatrixNearFracture() [2/3]

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

◆ HydroMechanicsLocalAssemblerMatrixNearFracture() [3/3]

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , typename IntegrationMethod , int GlobalDim>
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, 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,
bool const  is_axially_symmetric,
unsigned const  integration_order,
HydroMechanicsProcessData< GlobalDim > &  process_data 
)

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

34  : HydroMechanicsLocalAssemblerMatrix<ShapeFunctionDisplacement,
35  ShapeFunctionPressure,
36  IntegrationMethod, GlobalDim>(
37  e, n_variables, local_matrix_size, dofIndex_to_localIndex,
38  is_axially_symmetric, integration_order, process_data),
40 {
41  // currently not supporting multiple fractures
42  _fracture_props.push_back(process_data.fracture_property.get());
43  _fracID_to_local.insert({0, 0});
44 }
const T * getCoords() const
Definition: TemplatePoint.h:75
HydroMechanicsLocalAssemblerMatrix(HydroMechanicsLocalAssemblerMatrix const &)=delete
MeshLib::Node getCenterOfGravity(Element const &element)
Calculates the center of gravity for the mesh element.
Definition: Element.cpp:126

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

Member Function Documentation

◆ assembleWithJacobianConcrete()

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

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

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

55 {
56  auto p = const_cast<Eigen::VectorXd&>(local_x).segment(pressure_index,
58  auto p_dot = const_cast<Eigen::VectorXd&>(local_x_dot)
59  .segment(pressure_index, pressure_size);
60  if (_process_data.deactivate_matrix_in_flow)
61  {
64  }
65  auto const u = local_x.segment(displacement_index, displacement_size);
66  auto const u_dot =
67  local_x_dot.segment(displacement_index, displacement_size);
68 
69  auto rhs_p = local_b.segment(pressure_index, pressure_size);
70  auto rhs_u = local_b.segment(displacement_index, displacement_size);
71 
72  auto J_pp = local_J.block(pressure_index, pressure_index, pressure_size,
74  auto J_pu = local_J.block(pressure_index, displacement_index, pressure_size,
76  auto J_up = local_J.block(displacement_index, pressure_index,
78  auto J_uu = local_J.block(displacement_index, displacement_index,
80 
81  // levelset value of the element
82  // remark: this assumes the levelset function is uniform within an element
83  std::vector<double> levelsets = uGlobalEnrichments(
85  double const ele_levelset = levelsets[0]; // single fracture
86 
87  if (ele_levelset == 0)
88  {
89  // no DoF exists for displacement jumps. do the normal assembly
91  t, dt, p, p_dot, u, u_dot, rhs_p, rhs_u, J_pp, J_pu, J_uu, J_up);
92  return;
93  }
94 
95  // Displacement jumps should be taken into account
96 
97  // compute true displacements
98  auto const g = local_x.segment(displacement_jump_index, displacement_size);
99  auto const g_dot =
100  local_x_dot.segment(displacement_jump_index, displacement_size);
101  Eigen::VectorXd const total_u = u + ele_levelset * g;
102  Eigen::VectorXd const total_u_dot = u_dot + ele_levelset * g_dot;
103 
104  // evaluate residuals and Jacobians for pressure and displacements
105  Base::assembleBlockMatricesWithJacobian(t, dt, p, p_dot, total_u,
106  total_u_dot, rhs_p, rhs_u, J_pp,
107  J_pu, J_uu, J_up);
108 
109  // compute residuals and Jacobians for displacement jumps
110  auto rhs_g = local_b.segment(displacement_jump_index, displacement_size);
111  auto J_pg = local_J.block(pressure_index, displacement_jump_index,
113  auto J_ug = local_J.block(displacement_index, displacement_jump_index,
115  auto J_gp = local_J.block(displacement_jump_index, pressure_index,
117  auto J_gu = local_J.block(displacement_jump_index, displacement_index,
119  auto J_gg = local_J.block(displacement_jump_index, displacement_jump_index,
121 
122  rhs_g = ele_levelset * rhs_u;
123  J_pg = ele_levelset * J_pu;
124  J_ug = ele_levelset * J_uu;
125  J_gp = ele_levelset * J_up;
126  J_gu = ele_levelset * J_uu;
127  J_gg = ele_levelset * ele_levelset * J_uu;
128 }
void assembleBlockMatricesWithJacobian(double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &p_dot, Eigen::Ref< const Eigen::VectorXd > const &u, Eigen::Ref< const Eigen::VectorXd > const &u_dot, 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 , typename IntegrationMethod , int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, GlobalDim >::postTimestepConcreteWithVector ( double const  t,
double const  dt,
Eigen::VectorXd const &  local_x 
)
overrideprivatevirtual

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

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

136 {
137  auto p = const_cast<Eigen::VectorXd&>(local_x).segment(pressure_index,
138  pressure_size);
139  if (_process_data.deactivate_matrix_in_flow)
140  {
142  }
143  auto u = local_x.segment(displacement_index, displacement_size);
144 
145  // levelset value of the element
146  // remark: this assumes the levelset function is uniform within an element
147  std::vector<double> levelsets = uGlobalEnrichments(
149  double const ele_levelset = levelsets[0]; // single fracture
150 
151  if (ele_levelset == 0)
152  {
153  // no DoF exists for displacement jumps. do the normal assembly
155  return;
156  }
157 
158  // Displacement jumps should be taken into account
159 
160  // compute true displacements
161  auto const g = local_x.segment(displacement_jump_index, displacement_size);
162  Eigen::VectorXd const total_u = u + ele_levelset * g;
163 
164  // evaluate residuals and Jacobians for pressure and displacements
166 }
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 , typename IntegrationMethod , int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, IntegrationMethod, GlobalDim >::preTimestepConcrete ( std::vector< double > const &  ,
double const  ,
double const   
)
inlineoverrideprivatevirtual

Member Data Documentation

◆ _e_center_coords

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

◆ _fracID_to_local

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

◆ _fracture_props

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

◆ _ip_data

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

◆ _junction_props

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

◆ _process_data

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

Definition at line 112 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ displacement_index

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

Definition at line 120 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ displacement_jump_index

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

◆ displacement_size

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

Definition at line 121 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ kelvin_vector_size

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

Definition at line 123 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ pressure_index

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

Definition at line 118 of file HydroMechanicsLocalAssemblerMatrix.h.

◆ pressure_size

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

Definition at line 119 of file HydroMechanicsLocalAssemblerMatrix.h.


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