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

Detailed Description

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

Definition at line 32 of file HydroMechanicsLocalAssemblerFracture.h.

#include <HydroMechanicsLocalAssemblerFracture.h>

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

Public Member Functions

 HydroMechanicsLocalAssemblerFracture (HydroMechanicsLocalAssemblerFracture const &)=delete
 
 HydroMechanicsLocalAssemblerFracture (HydroMechanicsLocalAssemblerFracture &&)=delete
 
 HydroMechanicsLocalAssemblerFracture (MeshLib::Element const &e, 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)
 
void preTimestepConcrete (std::vector< double > const &, double const, double const) override
 
void postTimestepConcreteWithVector (const double t, double const dt, Eigen::VectorXd const &local_x) override
 
- 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_, 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, 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...
 

Private Types

using ShapeMatricesTypeDisplacement = ShapeMatrixPolicyType< ShapeFunctionDisplacement, GlobalDim >
 
using HMatricesType = HMatrixPolicyType< ShapeFunctionDisplacement, GlobalDim >
 
using HMatrixType = typename HMatricesType::HMatrixType
 
using ShapeMatricesTypePressure = ShapeMatrixPolicyType< ShapeFunctionPressure, GlobalDim >
 
using IntegrationPointDataType = IntegrationPointDataFracture< HMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, GlobalDim >
 
using GlobalDimVectorType = Eigen::Matrix< double, GlobalDim, 1 >
 

Private Member Functions

void assembleWithJacobianConcrete (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, Eigen::VectorXd &local_b, Eigen::MatrixXd &local_J) override
 
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 &g, Eigen::Ref< const Eigen::VectorXd > const &g_dot, Eigen::Ref< Eigen::VectorXd > rhs_p, Eigen::Ref< Eigen::VectorXd > rhs_g, Eigen::Ref< Eigen::MatrixXd > J_pp, Eigen::Ref< Eigen::MatrixXd > J_pg, Eigen::Ref< Eigen::MatrixXd > J_gg, Eigen::Ref< Eigen::MatrixXd > J_gp)
 

Private Attributes

HydroMechanicsProcessData< GlobalDim > & _process_data
 
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
 

Static Private Attributes

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
 

Additional Inherited Members

- Protected Attributes inherited from ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface
MeshLib::Element const & _element
 
bool const _is_axially_symmetric
 

Member Typedef Documentation

◆ GlobalDimVectorType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::GlobalDimVectorType = Eigen::Matrix<double, GlobalDim, 1>
private

Definition at line 98 of file HydroMechanicsLocalAssemblerFracture.h.

◆ HMatricesType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::HMatricesType = HMatrixPolicyType<ShapeFunctionDisplacement, GlobalDim>
private

Definition at line 83 of file HydroMechanicsLocalAssemblerFracture.h.

◆ HMatrixType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::HMatrixType = typename HMatricesType::HMatrixType
private

Definition at line 85 of file HydroMechanicsLocalAssemblerFracture.h.

◆ IntegrationPointDataType

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::IntegrationPointDataType = IntegrationPointDataFracture<HMatricesType, ShapeMatricesTypeDisplacement, ShapeMatricesTypePressure, GlobalDim>
private

Definition at line 92 of file HydroMechanicsLocalAssemblerFracture.h.

◆ ShapeMatricesTypeDisplacement

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

Definition at line 81 of file HydroMechanicsLocalAssemblerFracture.h.

◆ ShapeMatricesTypePressure

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
using ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::ShapeMatricesTypePressure = ShapeMatrixPolicyType<ShapeFunctionPressure, GlobalDim>
private

Definition at line 88 of file HydroMechanicsLocalAssemblerFracture.h.

Constructor & Destructor Documentation

◆ HydroMechanicsLocalAssemblerFracture() [1/3]

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

◆ HydroMechanicsLocalAssemblerFracture() [2/3]

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

◆ HydroMechanicsLocalAssemblerFracture() [3/3]

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::HydroMechanicsLocalAssemblerFracture ( MeshLib::Element const &  e,
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 37 of file HydroMechanicsLocalAssemblerFracture-impl.h.

46  e, is_axially_symmetric,
47  ShapeFunctionDisplacement::NPOINTS * GlobalDim +
48  ShapeFunctionPressure::NPOINTS,
49  dofIndex_to_localIndex),
50  _process_data(process_data)
51 {
52  assert(e.getDimension() == GlobalDim - 1);
53 
54  unsigned const n_integration_points =
55  integration_method.getNumberOfPoints();
56 
57  _ip_data.reserve(n_integration_points);
58 
59  auto const shape_matrices_u =
60  NumLib::initShapeMatrices<ShapeFunctionDisplacement,
62  e, is_axially_symmetric, integration_method);
63 
64  auto const shape_matrices_p =
65  NumLib::initShapeMatrices<ShapeFunctionPressure,
66  ShapeMatricesTypePressure, GlobalDim>(
67  e, is_axially_symmetric, integration_method);
68 
69  auto const& frac_prop = *_process_data.fracture_property;
70 
71  // Get element nodes for aperture0 interpolation from nodes to integration
72  // point. The aperture0 parameter is time-independent.
74  aperture0_node_values = frac_prop.aperture0.getNodalValuesOnElement(
75  e, /*time independent*/ 0);
76 
78  x_position.setElementID(e.getID());
79  for (unsigned ip = 0; ip < n_integration_points; ip++)
80  {
81  x_position.setIntegrationPoint(ip);
82 
83  _ip_data.emplace_back(*_process_data.fracture_model);
84  auto const& sm_u = shape_matrices_u[ip];
85  auto const& sm_p = shape_matrices_p[ip];
86  auto& ip_data = _ip_data[ip];
87  ip_data.integration_weight =
88  sm_u.detJ * sm_u.integralMeasure *
89  integration_method.getWeightedPoint(ip).getWeight();
90 
91  ip_data.H_u.setZero(GlobalDim,
92  ShapeFunctionDisplacement::NPOINTS * GlobalDim);
94  GlobalDim, ShapeFunctionDisplacement::NPOINTS,
96  HMatrixType>(sm_u.N, ip_data.H_u);
97  ip_data.N_p = sm_p.N;
98  ip_data.dNdx_p = sm_p.dNdx;
99 
100  // Initialize current time step values
101  ip_data.w.setZero(GlobalDim);
102  ip_data.sigma_eff.setZero(GlobalDim);
103 
104  // Previous time step values are not initialized and are set later.
105  ip_data.w_prev.resize(GlobalDim);
106  ip_data.sigma_eff_prev.resize(GlobalDim);
107 
108  ip_data.C.resize(GlobalDim, GlobalDim);
109 
110  ip_data.aperture0 = aperture0_node_values.dot(sm_u.N);
111  ip_data.aperture = ip_data.aperture0;
112 
113  ip_data.permeability_state =
114  frac_prop.permeability_model->getNewState();
115 
116  auto const initial_effective_stress =
117  _process_data.initial_fracture_effective_stress(0, x_position);
118  for (int i = 0; i < GlobalDim; i++)
119  {
120  ip_data.sigma_eff[i] = initial_effective_stress[i];
121  ip_data.sigma_eff_prev[i] = initial_effective_stress[i];
122  }
123  }
124 }
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
ShapeMatrixPolicyType< ShapeFunctionPressure, GlobalDim > ShapeMatricesTypePressure
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
ShapeMatrixPolicyType< ShapeFunctionDisplacement, GlobalDim > ShapeMatricesTypeDisplacement
HydroMechanicsLocalAssemblerInterface(MeshLib::Element const &element, bool const is_axially_symmetric, std::size_t n_local_size, std::vector< unsigned > dofIndex_to_localIndex)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
void computeHMatrix(N_Type const &N, HMatrixType &H)
Fills a H-matrix based on given shape function.
Definition: HMatrixUtils.h:53
VectorType< ShapeFunction::NPOINTS > NodalVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType

References ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_ip_data, ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::_process_data, ProcessLib::computeHMatrix(), MeshLib::Element::getDimension(), MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), ParameterLib::SpatialPosition::setElementID(), and ParameterLib::SpatialPosition::setIntegrationPoint().

Member Function Documentation

◆ assembleBlockMatricesWithJacobian()

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
void ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::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 &  g,
Eigen::Ref< const Eigen::VectorXd > const &  g_dot,
Eigen::Ref< Eigen::VectorXd >  rhs_p,
Eigen::Ref< Eigen::VectorXd >  rhs_g,
Eigen::Ref< Eigen::MatrixXd >  J_pp,
Eigen::Ref< Eigen::MatrixXd >  J_pg,
Eigen::Ref< Eigen::MatrixXd >  J_gg,
Eigen::Ref< Eigen::MatrixXd >  J_gp 
)
private

Definition at line 160 of file HydroMechanicsLocalAssemblerFracture-impl.h.

170 {
171  auto const& frac_prop = *_process_data.fracture_property;
172  auto const& R = frac_prop.R;
173 
174  // the index of a normal (normal to a fracture plane) component
175  // in a displacement vector
176  auto constexpr index_normal = GlobalDim - 1;
177 
178  typename ShapeMatricesTypePressure::NodalMatrixType laplace_p =
179  ShapeMatricesTypePressure::NodalMatrixType::Zero(pressure_size,
180  pressure_size);
181 
182  typename ShapeMatricesTypePressure::NodalMatrixType storage_p =
183  ShapeMatricesTypePressure::NodalMatrixType::Zero(pressure_size,
184  pressure_size);
185 
186  typename ShapeMatricesTypeDisplacement::template MatrixType<
188  Kgp = ShapeMatricesTypeDisplacement::template MatrixType<
190  pressure_size);
191 
192  // Projection of the body force vector at the element.
193  Eigen::MatrixXd const global2local_rotation =
194  R.template topLeftCorner<ShapeFunctionPressure::DIM, GlobalDim>();
195 
196  GlobalDimVectorType const gravity_vec = global2local_rotation.transpose() *
197  global2local_rotation *
198  _process_data.specific_body_force;
199 
201  x_position.setElementID(_element.getID());
202 
203  unsigned const n_integration_points = _ip_data.size();
204  for (unsigned ip = 0; ip < n_integration_points; ip++)
205  {
206  x_position.setIntegrationPoint(ip);
207 
208  auto& ip_data = _ip_data[ip];
209  auto const& ip_w = ip_data.integration_weight;
210  auto const& N_p = ip_data.N_p;
211  auto const& dNdx_p = ip_data.dNdx_p;
212  auto const& H_g = ip_data.H_u;
213  auto const& identity2 =
215 
216  auto& mat = ip_data.fracture_material;
217  auto& effective_stress = ip_data.sigma_eff;
218  auto const& effective_stress_prev = ip_data.sigma_eff_prev;
219  auto& w = ip_data.w;
220  auto const& w_prev = ip_data.w_prev;
221  auto& C = ip_data.C;
222  auto& state = *ip_data.material_state_variables;
223  auto& b_m = ip_data.aperture;
224 
225  double const S = frac_prop.specific_storage(t, x_position)[0];
226  double const mu = _process_data.fluid_viscosity(t, x_position)[0];
227  auto const alpha = frac_prop.biot_coefficient(t, x_position)[0];
228  auto const rho_fr = _process_data.fluid_density(t, x_position)[0];
229 
230  // displacement jumps in local coordinates
231  w.noalias() = R * H_g * g;
232 
233  // aperture
234  b_m = ip_data.aperture0 + w[index_normal];
235  if (b_m < 0.0)
236  {
237  DBUG(
238  "Element {:d}, gp {:d}: Fracture aperture is {:g}, but it must "
239  "be "
240  "non-negative. Setting it to zero.",
241  _element.getID(), ip, b_m);
242  b_m = 0;
243  }
244 
245  auto const initial_effective_stress =
246  _process_data.initial_fracture_effective_stress(0, x_position);
247 
248  Eigen::Map<typename HMatricesType::ForceVectorType const> const stress0(
249  initial_effective_stress.data(), initial_effective_stress.size());
250 
251  // local C, local stress
252  mat.computeConstitutiveRelation(
253  t, x_position, ip_data.aperture0, stress0, w_prev, w,
254  effective_stress_prev, effective_stress, C, state);
255 
256  auto& k = ip_data.permeability;
257  k = frac_prop.permeability_model->permeability(
258  ip_data.permeability_state.get(), ip_data.aperture0, b_m);
259 
260  // derivative of permeability respect to aperture
261  double const dk_db =
262  frac_prop.permeability_model->dpermeability_daperture(
263  ip_data.permeability_state.get(), ip_data.aperture0, b_m);
264 
265  //
266  // displacement equation, displacement jump part
267  //
268  rhs_g.noalias() -=
269  H_g.transpose() * R.transpose() * effective_stress * ip_w;
270  J_gg.noalias() += H_g.transpose() * R.transpose() * C * R * H_g * ip_w;
271 
272  //
273  // displacement equation, pressure part
274  //
275  Kgp.noalias() +=
276  H_g.transpose() * R.transpose() * alpha * identity2 * N_p * ip_w;
277 
278  //
279  // pressure equation, pressure part.
280  //
281  storage_p.noalias() += N_p.transpose() * b_m * S * N_p * ip_w;
282  laplace_p.noalias() +=
283  dNdx_p.transpose() * b_m * k / mu * dNdx_p * ip_w;
284  rhs_p.noalias() +=
285  dNdx_p.transpose() * b_m * k / mu * rho_fr * gravity_vec * ip_w;
286 
287  //
288  // pressure equation, displacement jump part.
289  //
290  GlobalDimVectorType const grad_head_over_mu =
291  (dNdx_p * p + rho_fr * gravity_vec) / mu;
292  Eigen::Matrix<double, 1, displacement_size> const mT_R_Hg =
293  identity2.transpose() * R * H_g;
294  // velocity in global coordinates
295  ip_data.darcy_velocity = -k * grad_head_over_mu;
296  J_pg.noalias() += N_p.transpose() * S * N_p * p_dot * mT_R_Hg * ip_w;
297  J_pg.noalias() +=
298  dNdx_p.transpose() * k * grad_head_over_mu * mT_R_Hg * ip_w;
299  J_pg.noalias() += dNdx_p.transpose() * b_m * dk_db * grad_head_over_mu *
300  mT_R_Hg * ip_w;
301  }
302 
303  // displacement equation, pressure part
304  J_gp.noalias() -= Kgp;
305 
306  // pressure equation, pressure part.
307  J_pp.noalias() += laplace_p + storage_p / dt;
308 
309  // pressure equation, displacement jump part.
310  J_pg.noalias() += Kgp.transpose() / dt;
311 
312  // pressure equation
313  rhs_p.noalias() -=
314  laplace_p * p + storage_p * p_dot + Kgp.transpose() * g_dot;
315 
316  // displacement equation
317  rhs_g.noalias() -= -Kgp * p;
318 }
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition: Logging.h:29
virtual std::size_t getID() const final
Returns the ID of the element.
Definition: Element.h:89
static const double t
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType

References MaterialPropertyLib::alpha, DBUG(), MaterialPropertyLib::effective_stress, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), and MathLib::t.

◆ assembleWithJacobianConcrete()

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

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 130 of file HydroMechanicsLocalAssemblerFracture-impl.h.

135 {
136  auto const p = local_x.segment(pressure_index, pressure_size);
137  auto const p_dot = local_xdot.segment(pressure_index, pressure_size);
138  auto const g = local_x.segment(displacement_index, displacement_size);
139  auto const g_dot =
140  local_xdot.segment(displacement_index, displacement_size);
141 
142  auto rhs_p = local_b.segment(pressure_index, pressure_size);
143  auto rhs_g = local_b.segment(displacement_index, displacement_size);
144  auto J_pp = local_J.block(pressure_index, pressure_index, pressure_size,
145  pressure_size);
146  auto J_pg = local_J.block(pressure_index, displacement_index, pressure_size,
148  auto J_gp = local_J.block(displacement_index, pressure_index,
150  auto J_gg = local_J.block(displacement_index, displacement_index,
152 
153  assembleBlockMatricesWithJacobian(t, dt, p, p_dot, g, g_dot, rhs_p, rhs_g,
154  J_pp, J_pg, J_gg, J_gp);
155 }
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 &g, Eigen::Ref< const Eigen::VectorXd > const &g_dot, Eigen::Ref< Eigen::VectorXd > rhs_p, Eigen::Ref< Eigen::VectorXd > rhs_g, Eigen::Ref< Eigen::MatrixXd > J_pp, Eigen::Ref< Eigen::MatrixXd > J_pg, Eigen::Ref< Eigen::MatrixXd > J_gg, Eigen::Ref< Eigen::MatrixXd > J_gp)

References MathLib::t.

◆ postTimestepConcreteWithVector()

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

Implements ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerInterface.

Definition at line 324 of file HydroMechanicsLocalAssemblerFracture-impl.h.

327 {
328  auto const nodal_g = local_x.segment(displacement_index, displacement_size);
329 
330  auto const& frac_prop = *_process_data.fracture_property;
331  auto const& R = frac_prop.R;
332  // the index of a normal (normal to a fracture plane) component
333  // in a displacement vector
334  auto constexpr index_normal = GlobalDim - 1;
335 
337  x_position.setElementID(_element.getID());
338 
339  unsigned const n_integration_points = _ip_data.size();
340  for (unsigned ip = 0; ip < n_integration_points; ip++)
341  {
342  x_position.setIntegrationPoint(ip);
343 
344  auto& ip_data = _ip_data[ip];
345  auto const& H_g = ip_data.H_u;
346  auto& mat = ip_data.fracture_material;
347  auto& effective_stress = ip_data.sigma_eff;
348  auto const& effective_stress_prev = ip_data.sigma_eff_prev;
349  auto& w = ip_data.w;
350  auto const& w_prev = ip_data.w_prev;
351  auto& C = ip_data.C;
352  auto& state = *ip_data.material_state_variables;
353  auto& b_m = ip_data.aperture;
354 
355  // displacement jumps in local coordinates
356  w.noalias() = R * H_g * nodal_g;
357 
358  // aperture
359  b_m = ip_data.aperture0 + w[index_normal];
360  if (b_m < 0.0)
361  {
362  DBUG(
363  "Element {:d}, gp {:d}: Fracture aperture is {:g}, but it is "
364  "expected to be non-negative.",
365  _element.getID(), ip, b_m);
366  }
367 
368  auto const initial_effective_stress =
369  _process_data.initial_fracture_effective_stress(0, x_position);
370 
371  Eigen::Map<typename HMatricesType::ForceVectorType const> const stress0(
372  initial_effective_stress.data(), initial_effective_stress.size());
373 
374  // local C, local stress
375  mat.computeConstitutiveRelation(
376  t, x_position, ip_data.aperture0, stress0, w_prev, w,
377  effective_stress_prev, effective_stress, C, state);
378  }
379 
380  double ele_b = 0;
381  double ele_k = 0;
382  typename HMatricesType::ForceVectorType ele_sigma_eff =
383  HMatricesType::ForceVectorType::Zero(GlobalDim);
384  typename HMatricesType::ForceVectorType ele_w =
385  HMatricesType::ForceVectorType::Zero(GlobalDim);
386  double ele_Fs = -std::numeric_limits<double>::max();
387  GlobalDimVectorType ele_velocity = GlobalDimVectorType::Zero();
388  for (auto const& ip : _ip_data)
389  {
390  ele_b += ip.aperture;
391  ele_k += ip.permeability;
392  ele_w += ip.w;
393  ele_sigma_eff += ip.sigma_eff;
394  ele_Fs = std::max(
395  ele_Fs, ip.material_state_variables->getShearYieldFunctionValue());
396  ele_velocity += ip.darcy_velocity;
397  }
398  ele_b /= static_cast<double>(n_integration_points);
399  ele_k /= static_cast<double>(n_integration_points);
400  ele_w /= static_cast<double>(n_integration_points);
401  ele_sigma_eff /= static_cast<double>(n_integration_points);
402  ele_velocity /= static_cast<double>(n_integration_points);
403  auto const element_id = _element.getID();
404  (*_process_data.mesh_prop_b)[element_id] = ele_b;
405  (*_process_data.mesh_prop_k_f)[element_id] = ele_k;
406  (*_process_data.mesh_prop_w_n)[element_id] = ele_w[index_normal];
407  (*_process_data.mesh_prop_w_s)[element_id] = ele_w[0];
408  (*_process_data.mesh_prop_fracture_stress_normal)[element_id] =
409  ele_sigma_eff[index_normal];
410  (*_process_data.mesh_prop_fracture_stress_shear)[element_id] =
411  ele_sigma_eff[0];
412  (*_process_data.mesh_prop_fracture_shear_failure)[element_id] = ele_Fs;
413 
414  if (GlobalDim == 3)
415  {
416  (*_process_data.mesh_prop_w_s2)[element_id] = ele_w[1];
417  (*_process_data.mesh_prop_fracture_stress_shear2)[element_id] =
418  ele_sigma_eff[1];
419  }
420 
421  for (unsigned i = 0; i < GlobalDim; i++)
422  {
423  (*_process_data.mesh_prop_velocity)[element_id * GlobalDim + i] =
424  ele_velocity[i];
425  }
426 }
VectorType< DisplacementDim > ForceVectorType
Definition: HMatrixUtils.h:45

References DBUG(), MaterialPropertyLib::effective_stress, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), and MathLib::t.

◆ preTimestepConcrete()

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

Member Data Documentation

◆ _ip_data

◆ _process_data

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

◆ displacement_index

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

Definition at line 108 of file HydroMechanicsLocalAssemblerFracture.h.

◆ displacement_size

template<typename ShapeFunctionDisplacement , typename ShapeFunctionPressure , int GlobalDim>
const int ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::displacement_size
staticprivate
Initial value:
=
ShapeFunctionDisplacement::NPOINTS * GlobalDim

Definition at line 109 of file HydroMechanicsLocalAssemblerFracture.h.

◆ pressure_index

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

Definition at line 106 of file HydroMechanicsLocalAssemblerFracture.h.

◆ pressure_size

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

Definition at line 107 of file HydroMechanicsLocalAssemblerFracture.h.


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