OGS
PhaseFieldFEM-impl.h
Go to the documentation of this file.
1
11#pragma once
12
15
16namespace ProcessLib
17{
18namespace PhaseField
19{
20template <typename ShapeFunction, int DisplacementDim>
23 double const t, double const dt, Eigen::VectorXd const& local_x,
24 Eigen::VectorXd const& /*local_x_prev*/, int const process_id,
25 std::vector<double>& /*local_M_data*/,
26 std::vector<double>& /*local_K_data*/,
27 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data)
28{
29 // For the equations with phase field.
30 if (process_id == phase_process_id)
31 {
32 assembleWithJacobianPhaseFieldEquations(t, dt, local_x, local_b_data,
33 local_Jac_data);
34 return;
35 }
36
37 // For the equations with deformation
38 assembleWithJacobianForDeformationEquations(t, dt, local_x, local_b_data,
39 local_Jac_data);
40}
41
42template <typename ShapeFunction, int DisplacementDim>
45 double const t, double const dt, Eigen::VectorXd const& local_x,
46 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data)
47{
48 using DeformationMatrix =
49 typename ShapeMatricesType::template MatrixType<displacement_size,
50 displacement_size>;
51 auto const d = local_x.template segment<phasefield_size>(phasefield_index);
52 auto const u =
53 local_x.template segment<displacement_size>(displacement_index);
54
55 auto local_Jac = MathLib::createZeroedMatrix<DeformationMatrix>(
56 local_Jac_data, displacement_size, displacement_size);
57
58 auto local_rhs = MathLib::createZeroedVector<DeformationVector>(
59 local_b_data, displacement_size);
60
62 x_position.setElementID(_element.getID());
63
64 auto local_pressure = 0.0;
65 if (_process_data.pressurized_crack)
66 {
67 local_pressure = _process_data.unity_pressure;
68 }
69
70 int const n_integration_points = _integration_method.getNumberOfPoints();
71 for (int ip = 0; ip < n_integration_points; ip++)
72 {
73 x_position.setIntegrationPoint(ip);
74 auto const& w = _ip_data[ip].integration_weight;
75 auto const& N = _ip_data[ip].N;
76 auto const& dNdx = _ip_data[ip].dNdx;
77
78 auto const x_coord =
79 NumLib::interpolateXCoordinate<ShapeFunction, ShapeMatricesType>(
80 _element, N);
81
82 auto const& B =
83 LinearBMatrix::computeBMatrix<DisplacementDim,
84 ShapeFunction::NPOINTS,
86 dNdx, N, x_coord, _is_axially_symmetric);
87
88 auto& eps = _ip_data[ip].eps;
89 eps.noalias() = B * u;
90 double const k = _process_data.residual_stiffness(t, x_position)[0];
91 double const ls = _process_data.crack_length_scale(t, x_position)[0];
92 double const d_ip = N.dot(d);
93 double const degradation =
94 _process_data.degradation_derivative->degradation(d_ip, k, ls);
95 _ip_data[ip].updateConstitutiveRelation(
96 t, x_position, dt, u, degradation,
97 _process_data.energy_split_model);
98
99 auto& sigma = _ip_data[ip].sigma;
100 auto const& D = _ip_data[ip].D;
101
102 typename ShapeMatricesType::template MatrixType<DisplacementDim,
103 displacement_size>
104 N_u = ShapeMatricesType::template MatrixType<
105 DisplacementDim, displacement_size>::Zero(DisplacementDim,
106 displacement_size);
107
108 for (int i = 0; i < DisplacementDim; ++i)
109 {
110 N_u.template block<1, displacement_size / DisplacementDim>(
111 i, i * displacement_size / DisplacementDim)
112 .noalias() = N;
113 }
114
115 auto const rho_sr = _process_data.solid_density(t, x_position)[0];
116 auto const& b = _process_data.specific_body_force;
117
118 local_rhs.noalias() -=
119 (B.transpose() * sigma - N_u.transpose() * rho_sr * b -
120 local_pressure * N_u.transpose() * dNdx * d) *
121 w;
122 local_Jac.noalias() += B.transpose() * D * B * w;
123 }
124}
125
126template <typename ShapeFunction, int DisplacementDim>
128 assembleWithJacobianPhaseFieldEquations(double const t, double const dt,
129 Eigen::VectorXd const& local_x,
130 std::vector<double>& local_b_data,
131 std::vector<double>& local_Jac_data)
132{
133 auto const d = local_x.template segment<phasefield_size>(phasefield_index);
134 auto const u =
135 local_x.template segment<displacement_size>(displacement_index);
136
137 auto local_Jac = MathLib::createZeroedMatrix<PhaseFieldMatrix>(
138 local_Jac_data, phasefield_size, phasefield_size);
139 auto local_rhs = MathLib::createZeroedVector<PhaseFieldVector>(
140 local_b_data, phasefield_size);
141
143 x_position.setElementID(_element.getID());
144
145 auto local_pressure = 0.0;
146 if (_process_data.static_pressurized_crack)
147 {
148 local_pressure = _process_data.unity_pressure;
149 }
150 else if (_process_data.propagating_pressurized_crack)
151 {
152 local_pressure = _process_data.pressure;
153 }
154
155 double const k = _process_data.residual_stiffness(t, x_position)[0];
156 double const ls = _process_data.crack_length_scale(t, x_position)[0];
157 double const gc = _process_data.crack_resistance(t, x_position)[0];
158
159 int const n_integration_points = _integration_method.getNumberOfPoints();
160 for (int ip = 0; ip < n_integration_points; ip++)
161 {
162 x_position.setIntegrationPoint(ip);
163 auto const& w = _ip_data[ip].integration_weight;
164 auto const& N = _ip_data[ip].N;
165 auto const& dNdx = _ip_data[ip].dNdx;
166
167 double const d_ip = N.dot(d);
168 double const degradation =
169 _process_data.degradation_derivative->degradation(d_ip, k, ls);
170 double const degradation_df1 =
171 _process_data.degradation_derivative->degradation_df1(d_ip, ls);
172 double const degradation_df2 =
173 _process_data.degradation_derivative->degradation_df2(d_ip, ls);
174 auto const x_coord =
175 NumLib::interpolateXCoordinate<ShapeFunction, ShapeMatricesType>(
176 _element, N);
177 auto const& B =
178 LinearBMatrix::computeBMatrix<DisplacementDim,
179 ShapeFunction::NPOINTS,
181 dNdx, N, x_coord, _is_axially_symmetric);
182
183 auto& eps = _ip_data[ip].eps;
184 eps.noalias() = B * u;
185 _ip_data[ip].updateConstitutiveRelation(
186 t, x_position, dt, u, degradation,
187 _process_data.energy_split_model);
188
189 auto const& strain_energy_tensile = _ip_data[ip].strain_energy_tensile;
190
191 auto& ip_data = _ip_data[ip];
192 ip_data.strain_energy_tensile = strain_energy_tensile;
193
194 typename ShapeMatricesType::template MatrixType<DisplacementDim,
195 displacement_size>
196 N_u = ShapeMatricesType::template MatrixType<
197 DisplacementDim, displacement_size>::Zero(DisplacementDim,
198 displacement_size);
199
200 for (int i = 0; i < DisplacementDim; ++i)
201 {
202 N_u.template block<1, displacement_size / DisplacementDim>(
203 i, i * displacement_size / DisplacementDim)
204 .noalias() = N;
205 }
206
207 local_Jac.noalias() +=
208 (N.transpose() * N * degradation_df2 * strain_energy_tensile) * w;
209
210 local_rhs.noalias() -=
211 (N.transpose() * degradation_df1 * strain_energy_tensile -
212 local_pressure * dNdx.transpose() * N_u * u) *
213 w;
214
215 switch (_process_data.phasefield_model)
216 {
218 {
219 auto const local_Jac_AT1 =
220 (gc * 0.75 * dNdx.transpose() * dNdx * ls * w).eval();
221 local_Jac.noalias() += local_Jac_AT1;
222
223 local_rhs.noalias() -=
224 gc * (-0.375 * N.transpose() / ls) * w + local_Jac_AT1 * d;
225 break;
226 }
228 {
229 auto const local_Jac_AT2 =
230 (gc *
231 (N.transpose() * N / ls + dNdx.transpose() * dNdx * ls) *
232 w)
233 .eval();
234 local_Jac.noalias() += local_Jac_AT2;
235
236 local_rhs.noalias() -=
237 local_Jac_AT2 * d - gc * (N.transpose() / ls) * w;
238 break;
239 }
241 {
242 auto const local_Jac_COHESIVE =
243 (2.0 / boost::math::double_constants::pi * gc *
244 (-N.transpose() * N / ls + dNdx.transpose() * dNdx * ls) *
245 w)
246 .eval();
247
248 local_Jac.noalias() += local_Jac_COHESIVE;
249
250 local_rhs.noalias() -= local_Jac_COHESIVE * d;
251 break;
252 }
253 }
254 }
255}
256
257template <typename ShapeFunction, int DisplacementDim>
260 std::size_t mesh_item_id,
261 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
262 std::vector<GlobalVector*> const& x, double const /*t*/,
263 double& crack_volume)
264{
265 std::vector<std::vector<GlobalIndexType>> indices_of_processes;
266 indices_of_processes.reserve(dof_tables.size());
267 std::transform(dof_tables.begin(), dof_tables.end(),
268 std::back_inserter(indices_of_processes),
269 [&](auto const dof_table)
270 { return NumLib::getIndices(mesh_item_id, *dof_table); });
271
272 auto local_coupled_xs = getCoupledLocalSolutions(x, indices_of_processes);
273 assert(local_coupled_xs.size() == displacement_size + phasefield_size);
274
275 auto const d = Eigen::Map<PhaseFieldVector const>(
276 &local_coupled_xs[phasefield_index], phasefield_size);
277 auto const u = Eigen::Map<DeformationVector const>(
278 &local_coupled_xs[displacement_index], displacement_size);
279
281 x_position.setElementID(_element.getID());
282
283 int const n_integration_points = _integration_method.getNumberOfPoints();
284 for (int ip = 0; ip < n_integration_points; ip++)
285 {
286 x_position.setIntegrationPoint(ip);
287 auto const& w = _ip_data[ip].integration_weight;
288 auto const& N = _ip_data[ip].N;
289 auto const& dNdx = _ip_data[ip].dNdx;
290
291 typename ShapeMatricesType::template MatrixType<DisplacementDim,
292 displacement_size>
293 N_u = ShapeMatricesType::template MatrixType<
294 DisplacementDim, displacement_size>::Zero(DisplacementDim,
295 displacement_size);
296
297 for (int i = 0; i < DisplacementDim; ++i)
298 {
299 N_u.template block<1, displacement_size / DisplacementDim>(
300 i, i * displacement_size / DisplacementDim)
301 .noalias() = N;
302 }
303
304 crack_volume += (N_u * u).dot(dNdx * d) * w;
305 }
306}
307
308template <typename ShapeFunction, int DisplacementDim>
310 std::size_t mesh_item_id,
311 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
312 std::vector<GlobalVector*> const& x, double const t, double& elastic_energy,
313 double& surface_energy, double& pressure_work)
314{
315 std::vector<std::vector<GlobalIndexType>> indices_of_processes;
316 indices_of_processes.reserve(dof_tables.size());
317 std::transform(dof_tables.begin(), dof_tables.end(),
318 std::back_inserter(indices_of_processes),
319 [&](auto const dof_table)
320 { return NumLib::getIndices(mesh_item_id, *dof_table); });
321
322 auto const local_coupled_xs =
323 getCoupledLocalSolutions(x, indices_of_processes);
324 assert(local_coupled_xs.size() == displacement_size + phasefield_size);
325
326 auto const d = Eigen::Map<PhaseFieldVector const>(
327 &local_coupled_xs[phasefield_index], phasefield_size);
328 auto const u = Eigen::Map<DeformationVector const>(
329 &local_coupled_xs[displacement_index], displacement_size);
330
332 x_position.setElementID(_element.getID());
333
334 double element_elastic_energy = 0.0;
335 double element_surface_energy = 0.0;
336 double element_pressure_work = 0.0;
337
338 double const gc = _process_data.crack_resistance(t, x_position)[0];
339 double const ls = _process_data.crack_length_scale(t, x_position)[0];
340 int const n_integration_points = _integration_method.getNumberOfPoints();
341 for (int ip = 0; ip < n_integration_points; ip++)
342 {
343 x_position.setIntegrationPoint(ip);
344 auto const& w = _ip_data[ip].integration_weight;
345 auto const& N = _ip_data[ip].N;
346 auto const& dNdx = _ip_data[ip].dNdx;
347 auto pressure_ip = _process_data.pressure;
348
349 typename ShapeMatricesType::template MatrixType<DisplacementDim,
350 displacement_size>
351 N_u = ShapeMatricesType::template MatrixType<
352 DisplacementDim, displacement_size>::Zero(DisplacementDim,
353 displacement_size);
354
355 for (int i = 0; i < DisplacementDim; ++i)
356 {
357 N_u.template block<1, displacement_size / DisplacementDim>(
358 i, i * displacement_size / DisplacementDim)
359 .noalias() = N;
360 }
361
362 element_elastic_energy += _ip_data[ip].elastic_energy * w;
363
364 double const d_ip = N.dot(d);
365 switch (_process_data.phasefield_model)
366 {
368 {
369 element_surface_energy +=
370 gc * 0.375 *
371 ((1 - d_ip) / ls + (dNdx * d).dot((dNdx * d)) * ls) * w;
372
373 break;
374 }
376 {
377 element_surface_energy += 0.5 * gc *
378 ((1 - d_ip) * (1 - d_ip) / ls +
379 (dNdx * d).dot((dNdx * d)) * ls) *
380 w;
381 break;
382 }
384 {
385 element_surface_energy +=
386 gc / boost::math::double_constants::pi *
387 ((1 - d_ip * d_ip) / ls + (dNdx * d).dot((dNdx * d)) * ls) *
388 w;
389 break;
390 }
391 }
392
393 if (_process_data.pressurized_crack)
394 {
395 element_pressure_work += pressure_ip * (N_u * u).dot(dNdx * d) * w;
396 }
397 }
398
399#ifdef USE_PETSC
400 int const n_all_nodes = indices_of_processes[1].size();
401 int const n_regular_nodes = std::count_if(
402 begin(indices_of_processes[1]), end(indices_of_processes[1]),
403 [](GlobalIndexType const& index) { return index >= 0; });
404 if (n_all_nodes != n_regular_nodes)
405 {
406 element_elastic_energy *=
407 static_cast<double>(n_regular_nodes) / n_all_nodes;
408 element_surface_energy *=
409 static_cast<double>(n_regular_nodes) / n_all_nodes;
410 element_pressure_work *=
411 static_cast<double>(n_regular_nodes) / n_all_nodes;
412 }
413#endif // USE_PETSC
414 elastic_energy += element_elastic_energy;
415 surface_energy += element_surface_energy;
416 pressure_work += element_pressure_work;
417}
418} // namespace PhaseField
419} // namespace ProcessLib
GlobalMatrix::IndexType GlobalIndexType
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
void assembleWithJacobianForDeformationEquations(double const t, double const dt, Eigen::VectorXd const &local_x, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
void assembleWithJacobianPhaseFieldEquations(double const t, double const dt, Eigen::VectorXd const &local_x, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
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) override
void computeCrackIntegral(std::size_t mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double &crack_volume) override
void computeEnergy(std::size_t mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double &elastic_energy, double &surface_energy, double &pressure_work) override
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, const double radius, const bool is_axially_symmetric)
Fills a B-matrix based on given shape function dN/dx values.
std::vector< double > getCoupledLocalSolutions(std::vector< GlobalVector * > const &global_solutions, std::vector< std::vector< GlobalIndexType > > const &indices)