223 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
224 updateConstitutiveRelations(
225 Eigen::Ref<Eigen::VectorXd const>
const local_x,
226 Eigen::Ref<Eigen::VectorXd const>
const local_x_prev,
228 double const dt,
IpData& ip_data,
231 assert(local_x.size() ==
232 pressure_size + displacement_size + temperature_size);
234 auto const [T, p, u] = localDOF(local_x);
235 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
237 auto const& solid_material =
239 _process_data.solid_materials, _process_data.material_ids,
242 auto const& medium = _process_data.media_map.getMedium(_element.getID());
243 auto const& liquid_phase = medium->phase(
"AqueousLiquid");
244 auto const& solid_phase = medium->phase(
"Solid");
245 auto*
const frozen_liquid_phase = medium->hasPhase(
"FrozenLiquid")
246 ? &medium->phase(
"FrozenLiquid")
250 auto const& N_u = ip_data.
N_u;
251 auto const& dNdx_u = ip_data.
dNdx_u;
253 auto const& N = ip_data.
N;
254 auto const& dNdx = ip_data.
dNdx;
256 auto const T_int_pt = N.dot(T);
257 auto const T_prev_int_pt = N.dot(T_prev);
258 double const dT_int_pt = T_int_pt - T_prev_int_pt;
266 ShapeFunctionDisplacement::NPOINTS,
268 dNdx_u, N_u, x_coord, _is_axially_symmetric);
272 auto& eps = ip_data.
eps;
273 eps.noalias() = B * u;
278 double const p_int_pt = N.dot(p);
283 auto const solid_density =
285 .template value<double>(vars, x_position, t, dt);
287 auto const porosity =
289 .template value<double>(vars, x_position, t, dt);
295 .template value<double>(vars, x_position, t, dt);
296 auto const& alpha = crv.alpha_biot;
299 t, x_position, dt,
static_cast<double>(T_int_pt));
300 auto const solid_skeleton_compressibility =
301 1 / solid_material.getBulkModulus(t, x_position, &C_el);
303 crv.beta_SR = (1 - alpha) * solid_skeleton_compressibility;
308 auto const& identity2 = Invariants::identity2;
309 auto const sigma_total =
310 (ip_data.
sigma_eff - alpha * p_int_pt * identity2).eval();
319 auto const intrinsic_permeability =
322 .value(vars, x_position, t, dt));
324 auto const fluid_density =
326 .template value<double>(vars, x_position, t, dt);
332 .template dValue<double>(
336 crv.fluid_compressibility = 1 / fluid_density * drho_dp;
338 double const fluid_volumetric_thermal_expansion_coefficient =
340 liquid_phase, vars, fluid_density, x_position, t, dt);
345 .template value<double>(vars, x_position, t, dt);
346 crv.K_over_mu = intrinsic_permeability / ip_data_output.
viscosity;
348 auto const& b = _process_data.specific_body_force;
351 crv.solid_linear_thermal_expansion_coefficient =
356 .value(vars, x_position, t, dt));
360 crv.solid_linear_thermal_expansion_coefficient * dT_int_pt;
362 crv.K_pT_thermal_osmosis =
363 (solid_phase.hasProperty(
368 thermal_osmosis_coefficient)
369 .value(vars, x_position, t, dt))
370 : Eigen::MatrixXd::Zero(DisplacementDim, DisplacementDim));
373 crv.K_pT_thermal_osmosis * dNdx * T +
374 crv.K_over_mu * fluid_density * b;
380 auto& eps_m = ip_data.
eps_m;
382 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
390 crv.rho = solid_density * (1 - porosity) + porosity * fluid_density;
393 porosity * fluid_volumetric_thermal_expansion_coefficient +
395 Invariants::trace(crv.solid_linear_thermal_expansion_coefficient);
408 .template value<double>(vars, x_position, t, dt);
409 crv.effective_thermal_conductivity =
414 .value(vars, x_position, t, dt));
418 crv.effective_thermal_conductivity.noalias() +=
419 fluid_density * crv.c_f *
421 _process_data.stabilizer, _element.getID(),
422 GlobalDimMatrixType::Zero(DisplacementDim, DisplacementDim),
429 .template value<double>(vars, x_position, t, dt);
432 crv.effective_volumetric_heat_capacity =
433 porosity * fluid_density * crv.c_f +
434 (1.0 - porosity) * solid_density * c_s;
436 if (frozen_liquid_phase)
439 double const phi_fr =
441 .template value<double>(vars, x_position, t, dt);
444 auto const frozen_liquid_value =
447 return (*frozen_liquid_phase)[p].template value<double>(
448 vars, x_position, t, dt);
451 double const rho_fr =
454 double const c_fr = frozen_liquid_value(
457 double const l_fr = frozen_liquid_value(
460 double const dphi_fr_dT =
462 .template dValue<double>(
466 double const phi_fr_prev = [&]()
471 .template value<double>(vars_prev, x_position, t, dt);
477 DisplacementDim>
const ice_linear_thermal_expansion_coefficient =
482 .value(vars, x_position, t, dt));
485 dthermal_strain_ice =
486 ice_linear_thermal_expansion_coefficient * dT_int_pt;
491 phase_change_expansion_coefficient =
495 phase_change_expansivity)
496 .value(vars, x_position, t, dt));
499 dphase_change_strain = phase_change_expansion_coefficient *
500 (phi_fr - phi_fr_prev) / porosity;
504 auto& eps0 = ip_data.
eps0;
505 auto const& eps0_prev = ip_data.
eps0_prev;
511 eps_m_ice.noalias() = eps_m_ice_prev + eps - eps_prev -
512 (eps0 - eps0_prev) - dthermal_strain_ice -
513 dphase_change_strain;
519 *_process_data.ice_constitutive_relation, vars_ice, t, x_position,
521 crv.effective_volumetric_heat_capacity +=
522 -phi_fr * fluid_density * crv.c_f + phi_fr * rho_fr * c_fr -
523 l_fr * rho_fr * dphi_fr_dT;
526 double const d2phi_fr_dT2 =
528 .template d2Value<double>(
533 crv.J_uu_fr = phi_fr * C_IR;
536 crv.r_u_fr = phi_fr * sigma_eff_ice;
538 crv.J_uT_fr = phi_fr * C_IR * ice_linear_thermal_expansion_coefficient;
540 crv.J_TT_fr = ((rho_fr * c_fr - fluid_density * crv.c_f) * dphi_fr_dT +
541 l_fr * rho_fr * d2phi_fr_dT2) *
553 DisplacementDim>::assembleWithJacobian(
double const t,
double const dt,
554 std::vector<double>
const& local_x,
555 std::vector<double>
const&
557 std::vector<double>& local_rhs_data,
558 std::vector<double>& local_Jac_data)
560 assert(local_x.size() ==
561 pressure_size + displacement_size + temperature_size);
564 Eigen::Map<Eigen::VectorXd const>(local_x.data(), local_x.size());
565 auto const x_prev = Eigen::Map<Eigen::VectorXd const>(local_x_prev.data(),
566 local_x_prev.size());
568 auto const [T, p, u] = localDOF(local_x);
569 auto const [T_prev, p_prev, u_prev] = localDOF(local_x_prev);
572 typename ShapeMatricesTypeDisplacement::template MatrixType<
573 temperature_size + displacement_size + pressure_size,
574 temperature_size + displacement_size + pressure_size>>(
575 local_Jac_data, displacement_size + pressure_size + temperature_size,
576 displacement_size + pressure_size + temperature_size);
579 typename ShapeMatricesTypeDisplacement::template VectorType<
580 displacement_size + pressure_size + temperature_size>>(
581 local_rhs_data, displacement_size + pressure_size + temperature_size);
584 MTT.setZero(temperature_size, temperature_size);
587 KTT.setZero(temperature_size, temperature_size);
590 KTp.setZero(temperature_size, pressure_size);
593 dKTT_dp.setZero(temperature_size, pressure_size);
596 laplace_p.setZero(pressure_size, pressure_size);
599 laplace_T.setZero(pressure_size, temperature_size);
602 storage_p.setZero(pressure_size, pressure_size);
605 storage_T.setZero(pressure_size, temperature_size);
607 typename ShapeMatricesTypeDisplacement::template MatrixType<
608 displacement_size, pressure_size>
610 Kup.setZero(displacement_size, pressure_size);
612 auto const& medium = _process_data.media_map.getMedium(_element.getID());
613 bool const has_frozen_liquid_phase = medium->hasPhase(
"FrozenLiquid");
615 unsigned const n_integration_points =
616 _integration_method.getNumberOfPoints();
618 std::vector<GlobalDimVectorType> ip_flux_vector;
619 double average_velocity_norm = 0.0;
620 ip_flux_vector.reserve(n_integration_points);
622 for (
unsigned ip = 0; ip < n_integration_points; ip++)
624 auto const& N_u = _ip_data[ip].N_u;
626 std::nullopt, _element.getID(), ip,
632 auto const crv = updateConstitutiveRelations(
633 x, x_prev, x_position, t, dt, _ip_data[ip], _ip_data_output[ip]);
635 auto const& w = _ip_data[ip].integration_weight;
637 auto const& dNdx_u = _ip_data[ip].dNdx_u;
639 auto const& N = _ip_data[ip].N;
640 auto const& dNdx = _ip_data[ip].dNdx;
642 auto const T_int_pt = N.dot(T);
650 ShapeFunctionDisplacement::NPOINTS,
652 dNdx_u, N_u, x_coord, _is_axially_symmetric);
654 auto const& b = _process_data.specific_body_force;
655 auto const velocity = _ip_data_output[ip].velocity;
661 if (has_frozen_liquid_phase)
664 .template block<displacement_size, displacement_size>(
665 displacement_index, displacement_index)
666 .noalias() += B.transpose() * crv.J_uu_fr * B * w;
668 local_rhs.template segment<displacement_size>(displacement_index)
669 .noalias() -= B.transpose() * crv.r_u_fr * w;
672 .template block<displacement_size, temperature_size>(
673 displacement_index, temperature_index)
674 .noalias() -= B.transpose() * crv.J_uT_fr * N * w;
678 .template block<displacement_size, displacement_size>(
679 displacement_index, displacement_index)
680 .noalias() += B.transpose() * crv.C * B * w;
683 .template block<displacement_size, temperature_size>(
684 displacement_index, temperature_index)
685 .noalias() -= B.transpose() * crv.C *
686 crv.solid_linear_thermal_expansion_coefficient * N *
689 local_rhs.template segment<displacement_size>(displacement_index)
690 .noalias() -= (B.transpose() * _ip_data[ip].sigma_eff -
691 N_u_op(N_u).transpose() * crv.rho * b) *
698 B.transpose() * crv.alpha_biot * Invariants::identity2 * N * w;
703 laplace_p.noalias() += dNdx.transpose() * crv.K_over_mu * dNdx * w;
705 storage_p.noalias() +=
707 (_ip_data[ip].porosity * crv.fluid_compressibility +
708 (crv.alpha_biot - _ip_data[ip].porosity) * crv.beta_SR) *
711 laplace_T.noalias() +=
712 dNdx.transpose() * crv.K_pT_thermal_osmosis * dNdx * w;
716 double const fluid_density = _ip_data_output[ip].fluid_density;
717 local_rhs.template segment<pressure_size>(pressure_index).noalias() +=
718 dNdx.transpose() * fluid_density * crv.K_over_mu * b * w;
722 storage_T.noalias() += N.transpose() * crv.beta * N * w;
733 dNdx.transpose() * crv.effective_thermal_conductivity * dNdx * w;
735 ip_flux_vector.emplace_back(velocity * fluid_density * crv.c_f);
736 average_velocity_norm += velocity.norm();
738 if (has_frozen_liquid_phase)
741 .template block<temperature_size, temperature_size>(
742 temperature_index, temperature_index)
743 .noalias() -= N.transpose() * crv.J_TT_fr * N * w;
747 N.transpose() * crv.effective_volumetric_heat_capacity * N * w;
753 dNdx.transpose() * T_int_pt * crv.K_pT_thermal_osmosis * dNdx * w;
756 dKTT_dp.noalias() -= fluid_density * crv.c_f * N.transpose() *
757 (dNdx * T).transpose() * crv.K_over_mu * dNdx * w;
803 _process_data.stabilizer, _ip_data, ip_flux_vector,
804 average_velocity_norm /
static_cast<double>(n_integration_points), KTT);
808 .template block<temperature_size, temperature_size>(temperature_index,
810 .noalias() += KTT + MTT / dt;
814 .template block<temperature_size, pressure_size>(temperature_index,
816 .noalias() += KTp + dKTT_dp;
820 .template block<displacement_size, pressure_size>(displacement_index,
826 .template block<pressure_size, temperature_size>(pressure_index,
828 .noalias() -= storage_T / dt - laplace_T;
832 .template block<pressure_size, pressure_size>(pressure_index,
834 .noalias() += laplace_p + storage_p / dt;
838 .template block<pressure_size, displacement_size>(pressure_index,
840 .noalias() += Kup.transpose() / dt;
843 local_rhs.template segment<pressure_size>(pressure_index).noalias() -=
844 laplace_p * p + laplace_T * T + storage_p * (p - p_prev) / dt -
845 storage_T * (T - T_prev) / dt + Kup.transpose() * (u - u_prev) / dt;
848 local_rhs.template segment<displacement_size>(displacement_index)
849 .noalias() += Kup * p;
852 local_rhs.template segment<temperature_size>(temperature_index).noalias() -=
853 KTT * T + MTT * (T - T_prev) / dt;
855 local_rhs.template segment<temperature_size>(temperature_index).noalias() -=
918 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
919 computeSecondaryVariableConcrete(
double const ,
double const ,
920 Eigen::VectorXd
const& local_x,
921 Eigen::VectorXd
const& )
923 auto const p = local_x.template segment<pressure_size>(pressure_index);
925 local_x.template segment<temperature_size>(temperature_index);
927 unsigned const n_integration_points =
928 _integration_method.getNumberOfPoints();
930 double fluid_density_avg = 0;
931 double viscosity_avg = 0;
934 KV sigma_avg = KV::Zero();
936 for (
unsigned ip = 0; ip < n_integration_points; ip++)
938 auto& ip_data = _ip_data[ip];
940 fluid_density_avg += _ip_data_output[ip].fluid_density;
941 viscosity_avg += _ip_data_output[ip].viscosity;
942 sigma_avg += ip_data.sigma_eff;
945 fluid_density_avg /= n_integration_points;
946 viscosity_avg /= n_integration_points;
947 sigma_avg /= n_integration_points;
949 (*_process_data.element_fluid_density)[_element.getID()] =
951 (*_process_data.element_viscosity)[_element.getID()] = viscosity_avg;
953 Eigen::Map<KV>(&(*_process_data.element_stresses)[_element.getID() *
954 KV::RowsAtCompileTime]) =
958 ShapeFunctionPressure,
typename ShapeFunctionDisplacement::MeshElement,
959 DisplacementDim>(_element, _is_axially_symmetric, p,
960 *_process_data.pressure_interpolated);
963 ShapeFunctionPressure,
typename ShapeFunctionDisplacement::MeshElement,
964 DisplacementDim>(_element, _is_axially_symmetric, T,
965 *_process_data.temperature_interpolated);