23 Eigen::VectorXd
const& local_x,
24 Eigen::VectorXd
const& local_x_prev,
26 std::vector<double>& local_b_data,
27 std::vector<double>& local_Jac_data)
29 if (process_id == _phase_field_process_id)
31 assembleWithJacobianForPhaseFieldEquations(t, local_x, local_b_data,
36 if (process_id == _heat_conduction_process_id)
38 assembleWithJacobianForHeatConductionEquations(
39 t, dt, local_x, local_x_prev, local_b_data, local_Jac_data);
44 assembleWithJacobianForDeformationEquations(t, dt, local_x, local_b_data,
51 double const t,
double const dt, Eigen::VectorXd
const& local_x,
52 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data)
54 assert(local_x.size() ==
55 phasefield_size + displacement_size + temperature_size);
57 auto const d = local_x.template segment<phasefield_size>(phasefield_index);
59 local_x.template segment<displacement_size>(displacement_index);
61 local_x.template segment<temperature_size>(temperature_index);
64 typename ShapeMatricesType::template MatrixType<displacement_size,
66 local_Jac_data, displacement_size, displacement_size);
69 typename ShapeMatricesType::template VectorType<displacement_size>>(
70 local_b_data, displacement_size);
75 int const n_integration_points = _integration_method.getNumberOfPoints();
76 for (
int ip = 0; ip < n_integration_points; ip++)
78 auto const& w = _ip_data[ip].integration_weight;
79 auto const& N = _ip_data[ip].N;
80 auto const& dNdx = _ip_data[ip].dNdx;
87 ShapeFunction::NPOINTS,
89 dNdx, N, x_coord, _is_axially_symmetric);
91 auto& eps = _ip_data[ip].eps;
93 auto const& D = _ip_data[ip].D;
94 auto const& sigma = _ip_data[ip].sigma;
96 auto const k = _process_data.residual_stiffness(t, x_position)[0];
97 auto rho_sr = _process_data.solid_density(t, x_position)[0];
98 auto const alpha = _process_data.linear_thermal_expansion_coefficient(
100 double const T0 = _process_data.reference_temperature;
101 auto const& b = _process_data.specific_body_force;
103 double const T_ip = N.dot(T);
104 double const delta_T = T_ip - T0;
106 double const rho_s = rho_sr / (1 + 3 * alpha * delta_T);
108 double const d_ip = N.dot(d);
109 double const degradation = d_ip * d_ip * (1 - k) + k;
111 eps.noalias() = B * u;
112 _ip_data[ip].updateConstitutiveRelation(t, x_position, dt, u, alpha,
113 delta_T, degradation);
115 typename ShapeMatricesType::template MatrixType<DisplacementDim,
117 N_u = ShapeMatricesType::template MatrixType<
118 DisplacementDim, displacement_size>::Zero(DisplacementDim,
121 for (
int i = 0; i < DisplacementDim; ++i)
123 N_u.template block<1, displacement_size / DisplacementDim>(
124 i, i * displacement_size / DisplacementDim)
128 local_Jac.noalias() += B.transpose() * D * B * w;
130 local_rhs.noalias() -=
131 (B.transpose() * sigma - N_u.transpose() * rho_s * b) * w;
138 double const t,
double const dt, Eigen::VectorXd
const& local_x,
139 Eigen::VectorXd
const& local_x_prev, std::vector<double>& local_b_data,
140 std::vector<double>& local_Jac_data)
142 assert(local_x.size() ==
143 phasefield_size + displacement_size + temperature_size);
145 auto const d = local_x.template segment<phasefield_size>(phasefield_index);
147 local_x.template segment<temperature_size>(temperature_index);
149 local_x_prev.template segment<temperature_size>(temperature_index);
151 typename ShapeMatricesType::template MatrixType<temperature_size,
153 local_Jac_data, temperature_size, temperature_size);
156 typename ShapeMatricesType::template VectorType<temperature_size>>(
157 local_b_data, temperature_size);
162 int const n_integration_points = _integration_method.getNumberOfPoints();
163 for (
int ip = 0; ip < n_integration_points; ip++)
165 auto const& w = _ip_data[ip].integration_weight;
166 auto const& N = _ip_data[ip].N;
167 auto const& dNdx = _ip_data[ip].dNdx;
169 auto& eps_m = _ip_data[ip].eps_m;
170 auto& heatflux = _ip_data[ip].heatflux;
172 auto rho_sr = _process_data.solid_density(t, x_position)[0];
173 auto const alpha = _process_data.linear_thermal_expansion_coefficient(
175 double const c = _process_data.specific_heat_capacity(t, x_position)[0];
177 _process_data.thermal_conductivity(t, x_position)[0];
178 auto const lambda_res =
179 _process_data.residual_thermal_conductivity(t, x_position)[0];
180 double const T0 = _process_data.reference_temperature;
182 double const d_ip = N.dot(d);
183 double const T_ip = N.dot(T);
184 double const T_dot_ip = (T_ip - N.dot(T_prev)) / dt;
185 double const delta_T = T_ip - T0;
187 double const rho_s = rho_sr / (1 + 3 * alpha * delta_T);
189 auto const lambda_eff =
190 d_ip * d_ip * lambda + (1 - d_ip) * (1 - d_ip) * lambda_res;
195 double const eps_m_trace = Invariants::trace(eps_m);
196 if (eps_m_trace >= 0)
198 local_Jac.noalias() += (dNdx.transpose() * lambda_eff * dNdx +
199 N.transpose() * rho_s * c * N / dt) *
202 local_rhs.noalias() -= (N.transpose() * rho_s * c * T_dot_ip +
203 dNdx.transpose() * lambda_eff * dNdx * T) *
206 heatflux.noalias() = -(lambda_eff * dNdx * T) * w;
210 local_Jac.noalias() += (dNdx.transpose() * lambda * dNdx +
211 N.transpose() * rho_s * c * N / dt) *
214 local_rhs.noalias() -= (N.transpose() * rho_s * c * T_dot_ip +
215 dNdx.transpose() * lambda * dNdx * T) *
218 heatflux.noalias() = -(lambda * dNdx * T) * w;
226 double const t, Eigen::VectorXd
const& local_x,
227 std::vector<double>& local_b_data, std::vector<double>& local_Jac_data)
229 assert(local_x.size() ==
230 phasefield_size + displacement_size + temperature_size);
232 auto const d = local_x.template segment<phasefield_size>(phasefield_index);
235 typename ShapeMatricesType::template MatrixType<phasefield_size,
237 local_Jac_data, phasefield_size, phasefield_size);
240 typename ShapeMatricesType::template VectorType<phasefield_size>>(
241 local_b_data, phasefield_size);
246 int const n_integration_points = _integration_method.getNumberOfPoints();
247 for (
int ip = 0; ip < n_integration_points; ip++)
249 auto const& w = _ip_data[ip].integration_weight;
250 auto const& N = _ip_data[ip].N;
251 auto const& dNdx = _ip_data[ip].dNdx;
253 auto const& strain_energy_tensile = _ip_data[ip].strain_energy_tensile;
255 auto const gc = _process_data.crack_resistance(t, x_position)[0];
256 auto const ls = _process_data.crack_length_scale(t, x_position)[0];
258 double const d_ip = N.dot(d);
260 local_Jac.noalias() += (dNdx.transpose() * gc * ls * dNdx +
261 N.transpose() * 2 * strain_energy_tensile * N +
262 N.transpose() * gc / ls * N) *
265 local_rhs.noalias() -=
266 (dNdx.transpose() * gc * ls * dNdx * d +
267 N.transpose() * d_ip * 2 * strain_energy_tensile -
268 N.transpose() * gc / ls * (1 - d_ip)) *
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.