OGS
ThermoMechanicalPhaseFieldFEM-impl.h
Go to the documentation of this file.
1 
11 #pragma once
12 
14 
17 
18 namespace ProcessLib
19 {
20 namespace ThermoMechanicalPhaseField
21 {
22 template <typename ShapeFunction, typename IntegrationMethod,
23  int DisplacementDim>
24 void ThermoMechanicalPhaseFieldLocalAssembler<ShapeFunction, IntegrationMethod,
25  DisplacementDim>::
26  assembleWithJacobianForStaggeredScheme(
27  double const t, double const dt, Eigen::VectorXd const& local_x,
28  Eigen::VectorXd const& local_xdot, const double /*dxdot_dx*/,
29  const double /*dx_dx*/, int const process_id,
30  std::vector<double>& /*local_M_data*/,
31  std::vector<double>& /*local_K_data*/,
32  std::vector<double>& local_b_data, std::vector<double>& local_Jac_data)
33 {
34  if (process_id == _phase_field_process_id)
35  {
36  assembleWithJacobianForPhaseFieldEquations(t, local_x, local_b_data,
37  local_Jac_data);
38  return;
39  }
40 
41  if (process_id == _heat_conduction_process_id)
42  {
43  assembleWithJacobianForHeatConductionEquations(
44  t, dt, local_x, local_xdot, local_b_data, local_Jac_data);
45  return;
46  }
47 
48  // For the equations with deformation
49  assembleWithJacobianForDeformationEquations(t, dt, local_x, local_b_data,
50  local_Jac_data);
51 }
52 
53 template <typename ShapeFunction, typename IntegrationMethod,
54  int DisplacementDim>
55 void ThermoMechanicalPhaseFieldLocalAssembler<ShapeFunction, IntegrationMethod,
56  DisplacementDim>::
57  assembleWithJacobianForDeformationEquations(
58  double const t, double const dt, Eigen::VectorXd const& local_x,
59  std::vector<double>& local_b_data, std::vector<double>& local_Jac_data)
60 {
61  assert(local_x.size() ==
62  phasefield_size + displacement_size + temperature_size);
63 
64  auto const d = local_x.template segment<phasefield_size>(phasefield_index);
65  auto const u =
66  local_x.template segment<displacement_size>(displacement_index);
67  auto const T =
68  local_x.template segment<temperature_size>(temperature_index);
69 
70  auto local_Jac = MathLib::createZeroedMatrix<
71  typename ShapeMatricesType::template MatrixType<displacement_size,
72  displacement_size>>(
73  local_Jac_data, displacement_size, displacement_size);
74 
75  auto local_rhs = MathLib::createZeroedVector<
76  typename ShapeMatricesType::template VectorType<displacement_size>>(
77  local_b_data, displacement_size);
78 
80  x_position.setElementID(_element.getID());
81 
82  int const n_integration_points = _integration_method.getNumberOfPoints();
83  for (int ip = 0; ip < n_integration_points; ip++)
84  {
85  x_position.setIntegrationPoint(ip);
86  auto const& w = _ip_data[ip].integration_weight;
87  auto const& N = _ip_data[ip].N;
88  auto const& dNdx = _ip_data[ip].dNdx;
89 
90  auto const x_coord =
91  NumLib::interpolateXCoordinate<ShapeFunction, ShapeMatricesType>(
92  _element, N);
93  auto const& B =
94  LinearBMatrix::computeBMatrix<DisplacementDim,
95  ShapeFunction::NPOINTS,
97  dNdx, N, x_coord, _is_axially_symmetric);
98 
99  auto& eps = _ip_data[ip].eps;
100 
101  auto const& C_tensile = _ip_data[ip].C_tensile;
102  auto const& C_compressive = _ip_data[ip].C_compressive;
103 
104  auto const& sigma = _ip_data[ip].sigma;
105 
106  auto const k = _process_data.residual_stiffness(t, x_position)[0];
107  auto rho_sr = _process_data.solid_density(t, x_position)[0];
108  auto const alpha = _process_data.linear_thermal_expansion_coefficient(
109  t, x_position)[0];
110  double const T0 = _process_data.reference_temperature;
111  auto const& b = _process_data.specific_body_force;
112 
113  double const T_ip = N.dot(T);
114  double const delta_T = T_ip - T0;
115  // calculate real density
116  double const rho_s = rho_sr / (1 + 3 * alpha * delta_T);
117 
118  double const d_ip = N.dot(d);
119  double const degradation = d_ip * d_ip * (1 - k) + k;
120 
121  auto const C_eff = degradation * C_tensile + C_compressive;
122  eps.noalias() = B * u;
123  _ip_data[ip].updateConstitutiveRelation(t, x_position, dt, u, alpha,
124  delta_T, degradation);
125 
126  typename ShapeMatricesType::template MatrixType<DisplacementDim,
127  displacement_size>
128  N_u = ShapeMatricesType::template MatrixType<
129  DisplacementDim, displacement_size>::Zero(DisplacementDim,
130  displacement_size);
131 
132  for (int i = 0; i < DisplacementDim; ++i)
133  {
134  N_u.template block<1, displacement_size / DisplacementDim>(
135  i, i * displacement_size / DisplacementDim)
136  .noalias() = N;
137  }
138 
139  local_Jac.noalias() += B.transpose() * C_eff * B * w;
140 
141  local_rhs.noalias() -=
142  (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
143  }
144 }
145 
146 template <typename ShapeFunction, typename IntegrationMethod,
147  int DisplacementDim>
148 void ThermoMechanicalPhaseFieldLocalAssembler<ShapeFunction, IntegrationMethod,
149  DisplacementDim>::
150  assembleWithJacobianForHeatConductionEquations(
151  double const t, double const dt, Eigen::VectorXd const& local_x,
152  Eigen::VectorXd const& local_xdot, std::vector<double>& local_b_data,
153  std::vector<double>& local_Jac_data)
154 {
155  assert(local_x.size() ==
156  phasefield_size + displacement_size + temperature_size);
157 
158  auto const d = local_x.template segment<phasefield_size>(phasefield_index);
159  auto const T =
160  local_x.template segment<temperature_size>(temperature_index);
161  auto const T_dot =
162  local_xdot.template segment<temperature_size>(temperature_index);
163  auto local_Jac = MathLib::createZeroedMatrix<
164  typename ShapeMatricesType::template MatrixType<temperature_size,
165  temperature_size>>(
166  local_Jac_data, temperature_size, temperature_size);
167 
168  auto local_rhs = MathLib::createZeroedVector<
169  typename ShapeMatricesType::template VectorType<temperature_size>>(
170  local_b_data, temperature_size);
171 
173  x_position.setElementID(_element.getID());
174 
175  int const n_integration_points = _integration_method.getNumberOfPoints();
176  for (int ip = 0; ip < n_integration_points; ip++)
177  {
178  x_position.setIntegrationPoint(ip);
179  auto const& w = _ip_data[ip].integration_weight;
180  auto const& N = _ip_data[ip].N;
181  auto const& dNdx = _ip_data[ip].dNdx;
182 
183  auto& eps_m = _ip_data[ip].eps_m;
184  auto& heatflux = _ip_data[ip].heatflux;
185 
186  auto rho_sr = _process_data.solid_density(t, x_position)[0];
187  auto const alpha = _process_data.linear_thermal_expansion_coefficient(
188  t, x_position)[0];
189  double const c = _process_data.specific_heat_capacity(t, x_position)[0];
190  auto const lambda =
191  _process_data.thermal_conductivity(t, x_position)[0];
192  auto const lambda_res =
193  _process_data.residual_thermal_conductivity(t, x_position)[0];
194  double const T0 = _process_data.reference_temperature;
195 
196  double const d_ip = N.dot(d);
197  double const T_ip = N.dot(T);
198  double const T_dot_ip = N.dot(T_dot);
199  double const delta_T = T_ip - T0;
200  // calculate real density
201  double const rho_s = rho_sr / (1 + 3 * alpha * delta_T);
202  // calculate effective thermal conductivity
203  auto const lambda_eff =
204  d_ip * d_ip * lambda + (1 - d_ip) * (1 - d_ip) * lambda_res;
205 
206  using Invariants = MathLib::KelvinVector::Invariants<
208 
209  double const eps_m_trace = Invariants::trace(eps_m);
210  if (eps_m_trace >= 0)
211  {
212  local_Jac.noalias() += (dNdx.transpose() * lambda_eff * dNdx +
213  N.transpose() * rho_s * c * N / dt) *
214  w;
215 
216  local_rhs.noalias() -= (N.transpose() * rho_s * c * T_dot_ip +
217  dNdx.transpose() * lambda_eff * dNdx * T) *
218  w;
219 
220  heatflux.noalias() = -(lambda_eff * dNdx * T) * w;
221  }
222  else
223  {
224  local_Jac.noalias() += (dNdx.transpose() * lambda * dNdx +
225  N.transpose() * rho_s * c * N / dt) *
226  w;
227 
228  local_rhs.noalias() -= (N.transpose() * rho_s * c * T_dot_ip +
229  dNdx.transpose() * lambda * dNdx * T) *
230  w;
231 
232  heatflux.noalias() = -(lambda * dNdx * T) * w;
233  }
234  }
235 }
236 
237 template <typename ShapeFunction, typename IntegrationMethod,
238  int DisplacementDim>
239 void ThermoMechanicalPhaseFieldLocalAssembler<ShapeFunction, IntegrationMethod,
240  DisplacementDim>::
241  assembleWithJacobianForPhaseFieldEquations(
242  double const t, Eigen::VectorXd const& local_x,
243  std::vector<double>& local_b_data, std::vector<double>& local_Jac_data)
244 {
245  assert(local_x.size() ==
246  phasefield_size + displacement_size + temperature_size);
247 
248  auto const d = local_x.template segment<phasefield_size>(phasefield_index);
249 
250  auto local_Jac = MathLib::createZeroedMatrix<
251  typename ShapeMatricesType::template MatrixType<phasefield_size,
252  phasefield_size>>(
253  local_Jac_data, phasefield_size, phasefield_size);
254 
255  auto local_rhs = MathLib::createZeroedVector<
256  typename ShapeMatricesType::template VectorType<phasefield_size>>(
257  local_b_data, phasefield_size);
258 
260  x_position.setElementID(_element.getID());
261 
262  int const n_integration_points = _integration_method.getNumberOfPoints();
263  for (int ip = 0; ip < n_integration_points; ip++)
264  {
265  x_position.setIntegrationPoint(ip);
266  auto const& w = _ip_data[ip].integration_weight;
267  auto const& N = _ip_data[ip].N;
268  auto const& dNdx = _ip_data[ip].dNdx;
269 
270  auto const& strain_energy_tensile = _ip_data[ip].strain_energy_tensile;
271 
272  auto const gc = _process_data.crack_resistance(t, x_position)[0];
273  auto const ls = _process_data.crack_length_scale(t, x_position)[0];
274 
275  double const d_ip = N.dot(d);
276 
277  local_Jac.noalias() += (dNdx.transpose() * gc * ls * dNdx +
278  N.transpose() * 2 * strain_energy_tensile * N +
279  N.transpose() * gc / ls * N) *
280  w;
281 
282  local_rhs.noalias() -=
283  (dNdx.transpose() * gc * ls * dNdx * d +
284  N.transpose() * d_ip * 2 * strain_energy_tensile -
285  N.transpose() * gc / ls * (1 - d_ip)) *
286  w;
287  }
288 }
289 
290 } // namespace ThermoMechanicalPhaseField
291 } // namespace ProcessLib
void setElementID(std::size_t element_id)
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
Definition: BMatrixPolicy.h:55
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Definition: KelvinVector.h:23
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:32
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
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.
Definition: LinearBMatrix.h:42