16template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
19 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
20 HydroMechanicsLocalAssemblerMatrixNearFracture(
22 std::size_t
const n_variables,
23 std::size_t
const local_matrix_size,
24 std::vector<unsigned>
const& dofIndex_to_localIndex,
26 bool const is_axially_symmetric,
29 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>(
30 e, n_variables, local_matrix_size, dofIndex_to_localIndex,
31 integration_method, is_axially_symmetric, process_data),
36 _fracID_to_local.insert({fid, _fracture_props.size()});
43template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
45void HydroMechanicsLocalAssemblerMatrixNearFracture<
46 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
47 assembleWithJacobianConcrete(
double const t,
double const dt,
48 Eigen::VectorXd
const& local_x,
49 Eigen::VectorXd
const& local_x_prev,
50 Eigen::VectorXd& local_b,
51 Eigen::MatrixXd& local_J)
53 auto p =
const_cast<Eigen::VectorXd&
>(local_x).segment(
pressure_index,
55 auto p_prev =
const_cast<Eigen::VectorXd&
>(local_x_prev)
81 double const ele_levelset = levelsets[0];
83 if (ele_levelset == 0)
87 t, dt, p, p_prev, u, u_prev, rhs_p, rhs_u, J_pp, J_pu, J_uu, J_up);
97 Eigen::VectorXd
const total_u = u + ele_levelset * g;
98 Eigen::VectorXd
const total_u_prev = u_prev + ele_levelset * g_prev;
102 total_u_prev, rhs_p, rhs_u, J_pp,
118 rhs_g = ele_levelset * rhs_u;
119 J_pg = ele_levelset * J_pu;
120 J_ug = ele_levelset * J_uu;
121 J_gp = ele_levelset * J_up;
122 J_gu = ele_levelset * J_uu;
123 J_gg = ele_levelset * ele_levelset * J_uu;
126template <
typename ShapeFunctionDisplacement,
typename ShapeFunctionPressure,
129 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
130 postTimestepConcreteWithVector(
double const t,
double const dt,
131 Eigen::VectorXd
const& local_x)
133 auto p =
const_cast<Eigen::VectorXd&
>(local_x).segment(
pressure_index,
145 double const ele_levelset = levelsets[0];
147 if (ele_levelset == 0)
158 Eigen::VectorXd
const total_u = u + ele_levelset * g;
std::size_t getID() const
Returns the ID of the element.
HydroMechanicsProcessData< DisplacementDim > & _process_data
static const int displacement_jump_index
std::vector< FractureProperty * > _fracture_props
std::vector< JunctionProperty * > _junction_props
static const int displacement_size
static const int pressure_size
static const int displacement_index
Eigen::Vector3d _e_center_coords
std::unordered_map< int, int > _fracID_to_local
static const int pressure_index
void postTimestepConcreteWithBlockVectors(double const t, double const dt, Eigen::Ref< const Eigen::VectorXd > const &p, Eigen::Ref< const Eigen::VectorXd > const &u)
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 setPressureOfInactiveNodes(double const t, Eigen::Ref< Eigen::VectorXd > p)
HydroMechanicsLocalAssemblerMatrix(HydroMechanicsLocalAssemblerMatrix const &)=delete
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)
std::vector< std::vector< int > > vec_ele_connected_fractureIDs