34 std::size_t
const n_variables,
36 std::vector<unsigned>
const& dofIndex_to_localIndex,
38 bool const is_axially_symmetric,
41 n_variables * ShapeFunction::NPOINTS * DisplacementDim,
42 dofIndex_to_localIndex),
43 _process_data(process_data),
44 _integration_method(integration_method),
47 DisplacementDim>(e, is_axially_symmetric,
48 _integration_method)),
53 unsigned const n_integration_points =
56 _ip_data.reserve(n_integration_points);
60 auto frac_id =
_process_data.map_materialID_to_fractureID[mat_id];
69 ranges::views::transform(
72 ranges::to<std::vector>;
76 for (
unsigned ip = 0; ip < n_integration_points; ip++)
81 ip_data.integration_weight =
83 sm.integralMeasure * sm.detJ;
84 ip_data.h_matrices.setZero(DisplacementDim,
85 ShapeFunction::NPOINTS * DisplacementDim);
92 ip_data.w.setZero(DisplacementDim);
93 ip_data.sigma.setZero(DisplacementDim);
96 ip_data.sigma_prev.resize(DisplacementDim);
97 ip_data.w_prev.resize(DisplacementDim);
99 ip_data.C.resize(DisplacementDim, DisplacementDim);
102 ip_data.aperture_prev = ip_data.aperture0;
111 Eigen::VectorXd
const& local_u,
112 Eigen::VectorXd& local_b, Eigen::MatrixXd& local_J)
114 auto const N_DOF_PER_VAR = ShapeFunction::NPOINTS * DisplacementDim;
115 auto const n_fractures = _fracture_props.size();
116 auto const n_junctions = _junction_props.size();
117 auto const n_enrich_var = n_fractures + n_junctions;
131 using BlockVectorType =
132 typename Eigen::VectorXd::FixedSegmentReturnType<N_DOF_PER_VAR>::Type;
133 using BlockMatrixType =
134 Eigen::Block<Eigen::MatrixXd, N_DOF_PER_VAR, N_DOF_PER_VAR>;
136 std::vector<BlockVectorType> vec_local_b_g;
137 for (
unsigned i = 0; i < n_enrich_var; i++)
139 vec_local_b_g.push_back(
140 local_b.segment<N_DOF_PER_VAR>(N_DOF_PER_VAR * i));
142 std::vector<std::vector<BlockMatrixType>> vec_local_J_gg(n_enrich_var);
143 for (
unsigned i = 0; i < n_enrich_var; i++)
145 for (
unsigned j = 0; j < n_enrich_var; j++)
147 auto sub_gg = local_J.block<N_DOF_PER_VAR, N_DOF_PER_VAR>(
148 N_DOF_PER_VAR * i, N_DOF_PER_VAR * j);
149 vec_local_J_gg[i].push_back(sub_gg);
153 std::vector<Eigen::VectorXd> vec_nodal_g;
154 for (
unsigned i = 0; i < n_enrich_var; i++)
156 auto sub =
const_cast<Eigen::VectorXd&
>(local_u).segment<N_DOF_PER_VAR>(
158 vec_nodal_g.push_back(sub);
166 int const index_normal = DisplacementDim - 1;
167 auto const& R = _fracture_property->R;
169 unsigned const n_integration_points =
170 _integration_method.getNumberOfPoints();
175 for (
unsigned ip = 0; ip < n_integration_points; ip++)
177 auto& ip_data = _ip_data[ip];
178 auto const& integration_weight = ip_data.integration_weight;
179 auto const& H = ip_data.h_matrices;
180 auto& mat = ip_data.fracture_material;
181 auto& sigma = ip_data.sigma;
182 auto const& sigma_prev = ip_data.sigma_prev;
184 auto const& w_prev = ip_data.w_prev;
186 auto& state = *ip_data.material_state_variables;
187 auto const& N = _secondary_data.N[ip];
191 _fracture_property->fracture_id, _fracture_props, _junction_props,
192 _fracID_to_local, ip_physical_coords));
196 Eigen::VectorXd nodal_gap(N_DOF_PER_VAR);
198 for (
unsigned i = 0; i < n_enrich_var; i++)
200 nodal_gap += levelsets[i] * vec_nodal_g[i];
204 w.noalias() = R * H * nodal_gap;
207 ip_data.aperture = ip_data.aperture0 + w[index_normal];
210 mat.computeConstitutiveRelation(
211 t, x_position, ip_data.aperture0,
212 Eigen::Matrix<double, DisplacementDim, 1>::Zero(),
216 w_prev, w, sigma_prev, sigma, C, state);
219 for (
unsigned i = 0; i < n_enrich_var; i++)
221 vec_local_b_g[i].noalias() -= levelsets[i] * H.transpose() *
222 R.transpose() * sigma *
227 for (
unsigned i = 0; i < n_enrich_var; i++)
229 for (
unsigned j = 0; j < n_enrich_var; j++)
232 vec_local_J_gg[i][j].noalias() +=
233 (levelsets[i] * H.transpose() * R.transpose()) * C *
234 (levelsets[j] * R * H) * integration_weight;
243 Eigen::VectorXd
const& local_u)
245 auto const N_DOF_PER_VAR = ShapeFunction::NPOINTS * DisplacementDim;
246 auto const n_fractures = _fracture_props.size();
247 auto const n_junctions = _junction_props.size();
248 auto const n_enrich_var = n_fractures + n_junctions;
250 auto const& R = _fracture_property->R;
254 int const index_normal = DisplacementDim - 1;
256 unsigned const n_integration_points =
257 _integration_method.getNumberOfPoints();
260 auto const e_id = _element.getID();
263 std::vector<Eigen::VectorXd> vec_nodal_g;
264 for (
unsigned i = 0; i < n_enrich_var; i++)
266 auto sub =
const_cast<Eigen::VectorXd&
>(local_u).segment<N_DOF_PER_VAR>(
268 vec_nodal_g.push_back(sub);
271 for (
unsigned ip = 0; ip < n_integration_points; ip++)
273 auto& ip_data = _ip_data[ip];
274 auto const& H = ip_data.h_matrices;
275 auto& mat = ip_data.fracture_material;
276 auto& sigma = ip_data.sigma;
277 auto const& sigma_prev = ip_data.sigma_prev;
279 auto const& w_prev = ip_data.w_prev;
281 auto& state = *ip_data.material_state_variables;
282 auto& b_m = ip_data.aperture;
283 auto const& N = _secondary_data.N[ip];
287 _fracture_property->fracture_id, _fracture_props, _junction_props,
288 _fracID_to_local, ip_physical_coords));
292 Eigen::VectorXd nodal_gap(N_DOF_PER_VAR);
294 for (
unsigned i = 0; i < n_enrich_var; i++)
296 nodal_gap += levelsets[i] * vec_nodal_g[i];
300 w.noalias() = R * H * nodal_gap;
303 b_m = ip_data.aperture0 + w[index_normal];
307 "Element {:d}, gp {:d}: Fracture aperture is {:g}, but it must "
309 _element.getID(), ip, b_m);
313 mat.computeConstitutiveRelation(
314 t, x_position, ip_data.aperture0,
315 Eigen::Matrix<double, DisplacementDim, 1>::Zero(),
319 w_prev, w, sigma_prev, sigma, C, state);
326 for (
unsigned ip = 0; ip < n_integration_points; ip++)
328 auto const& ip_data = _ip_data[ip];
329 ele_b += ip_data.aperture;
331 ele_sigma += ip_data.sigma;
333 ele_b /= n_integration_points;
334 ele_w /= n_integration_points;
335 ele_sigma /= n_integration_points;
336 (*_process_data.mesh_prop_b)[_element.getID()] = ele_b;
338 Eigen::Map<GlobalDimVectorType>(
339 &(*_process_data.element_fracture_stresses)[e_id * DisplacementDim]) =
342 Eigen::Map<GlobalDimVectorType>(
343 &(*_process_data.element_local_jumps)[e_id * DisplacementDim]) = ele_w;