227 ShapeFunctionDisplacement, ShapeFunctionPressure, DisplacementDim>::
228 updateConstitutiveRelations(
229 Eigen::Ref<Eigen::VectorXd const>
const local_x,
230 Eigen::Ref<Eigen::VectorXd const>
const local_x_prev,
232 double const dt,
IpData& ip_data,
235 assert(local_x.size() ==
238 auto const [T, p, u] =
localDOF(local_x);
239 auto const [T_prev, p_prev, u_prev] =
localDOF(local_x_prev);
241 auto const& solid_material =
247 auto const& liquid_phase = medium->phase(
"AqueousLiquid");
248 auto const& solid_phase = medium->phase(
"Solid");
249 auto*
const frozen_liquid_phase = medium->hasPhase(
"FrozenLiquid")
250 ? &medium->phase(
"FrozenLiquid")
254 auto const& N_u = ip_data.
N_u;
255 auto const& dNdx_u = ip_data.
dNdx_u;
257 auto const& N = ip_data.
N;
258 auto const& dNdx = ip_data.
dNdx;
260 auto const T_int_pt = N.dot(T);
261 auto const T_prev_int_pt = N.dot(T_prev);
262 double const dT_int_pt = T_int_pt - T_prev_int_pt;
268 ShapeFunctionDisplacement::NPOINTS,
274 auto& eps = ip_data.
eps;
275 eps.noalias() = B * u;
280 double const p_int_pt = N.dot(p);
282 double const p_prev_int_pt = N.dot(p_prev);
283 double const dp_int_pt = p_int_pt - p_prev_int_pt;
287 auto const solid_density =
289 .template value<double>(vars, x_position, t, dt);
291 auto const drho_SR_dT =
293 .template dValue<double>(vars,
297 auto const porosity =
299 .template value<double>(vars, x_position, t, dt);
305 .template value<double>(vars, x_position, t, dt);
306 auto const& alpha = crv.alpha_biot;
309 t, x_position, dt,
static_cast<double>(T_int_pt));
310 auto const solid_skeleton_compressibility =
311 1 / solid_material.getBulkModulus(t, x_position, &C_el);
313 crv.beta_SR = (1 - alpha) * solid_skeleton_compressibility;
319 auto const sigma_total =
320 (ip_data.
sigma_eff - alpha * p_int_pt * identity2).eval();
329 auto const intrinsic_permeability =
332 .value(vars, x_position, t, dt));
334 auto const fluid_density =
336 .template value<double>(vars, x_position, t, dt);
342 .template dValue<double>(
346 crv.fluid_compressibility = 1 / fluid_density * drho_dp;
350 .template dValue<double>(vars,
354 double const fluid_volumetric_thermal_expansion_coefficient =
356 liquid_phase, vars, fluid_density, x_position, t, dt);
361 .template value<double>(vars, x_position, t, dt);
362 crv.K_over_mu = intrinsic_permeability / ip_data_output.
viscosity;
366 if (frozen_liquid_phase)
370 .template value<double>(vars, x_position, t, dt);
372 double const dphi_fr_dT =
374 .template dValue<double>(
387 .template value<double>(vars, x_position, t, dt);
391 auto const dk_rel_dphi_ice =
395 .template dValue<double>(
398 crv.dk_rel_dT = dk_rel_dphi_ice * dphi_fr_dT / porosity;
404 crv.solid_linear_thermal_expansion_coefficient =
409 .value(vars, x_position, t, dt));
413 crv.solid_linear_thermal_expansion_coefficient * dT_int_pt;
415 crv.K_pT_thermal_osmosis =
416 (solid_phase.hasProperty(
421 thermal_osmosis_coefficient)
422 .value(vars, x_position, t, dt))
423 : Eigen::MatrixXd::Zero(DisplacementDim, DisplacementDim));
426 -crv.k_rel * crv.K_over_mu * dNdx * p -
427 crv.K_pT_thermal_osmosis * dNdx * T +
428 (fluid_density * crv.k_rel) * crv.K_over_mu * b;
431 -crv.dk_rel_dT * crv.K_over_mu * dNdx * p +
433 (crv.dk_rel_dT * fluid_density + crv.k_rel * crv.drho_LR_dT) *
439 auto& eps_m = ip_data.
eps_m;
441 eps_m.noalias() = eps_m_prev + eps - eps_prev - dthermal_strain;
449 crv.rho = solid_density * (1 - porosity) + porosity * fluid_density;
452 porosity * fluid_volumetric_thermal_expansion_coefficient +
467 .template value<double>(vars, x_position, t, dt);
468 crv.effective_thermal_conductivity =
473 .value(vars, x_position, t, dt));
483 crv.effective_thermal_conductivity.noalias() +=
484 fluid_density * crv.c_f *
487 GlobalDimMatrixType::Zero(DisplacementDim, DisplacementDim),
494 .template value<double>(vars, x_position, t, dt);
497 crv.effective_volumetric_heat_capacity =
498 porosity * fluid_density * crv.c_f +
499 (1.0 - porosity) * solid_density * c_s;
500 double dC_eff_dT = porosity * crv.drho_LR_dT * crv.c_f +
501 (1.0 - porosity) * drho_SR_dT * c_s;
503 if (frozen_liquid_phase)
506 double const phi_fr = ip_data.
phi_fr;
508 auto const frozen_liquid_value =
511 return (*frozen_liquid_phase)[p].template value<double>(
512 vars, x_position, t, dt);
515 double const c_fr = frozen_liquid_value(
518 double const l_fr = frozen_liquid_value(
521 double const dphi_fr_dT =
523 .template dValue<double>(
527 double const d2phi_fr_dT2 =
529 .template d2Value<double>(
534 double const phi_fr_prev = [&]()
539 .template value<double>(vars_prev, x_position, t, dt);
543 double const rho_fr =
545 ip_data_output.
rho_fr = rho_fr;
547 crv.rho += ip_data.
phi_fr * rho_fr - ip_data.
phi_fr * fluid_density;
549 -dphi_fr_dT * porosity * (1. - rho_fr / fluid_density);
550 double const dmass_exchange_dT =
551 -d2phi_fr_dT2 * porosity * (1. - rho_fr / fluid_density) +
552 dphi_fr_dT * porosity * rho_fr * crv.drho_LR_dT /
553 (fluid_density * fluid_density);
557 DisplacementDim>
const ice_linear_thermal_expansion_coefficient =
562 .value(vars, x_position, t, dt));
565 dthermal_strain_ice =
566 ice_linear_thermal_expansion_coefficient * dT_int_pt;
573 crv.solid_linear_thermal_expansion_coefficient);
578 phase_change_expansion_coefficient =
582 phase_change_expansivity)
583 .value(vars, x_position, t, dt));
586 dphase_change_strain = phase_change_expansion_coefficient *
587 (phi_fr - phi_fr_prev) / porosity;
591 auto& eps0 = ip_data.
eps0;
592 auto const& eps0_prev = ip_data.
eps0_prev;
598 eps_m_ice.noalias() = eps_m_ice_prev + eps - eps_prev -
599 (eps0 - eps0_prev) - dthermal_strain_ice -
600 dphase_change_strain;
606 *
_process_data.ice_constitutive_relation, vars_ice, t, x_position,
611 static_cast<double>(T_int_pt));
613 1. /
_process_data.ice_constitutive_relation->getBulkModulus(
614 t, x_position, &C_el_ice);
616 crv.effective_volumetric_heat_capacity +=
617 -phi_fr * fluid_density * crv.c_f + phi_fr * rho_fr * c_fr -
618 l_fr * rho_fr * dphi_fr_dT;
620 crv.J_uu_fr = phi_fr * C_IR;
623 crv.r_u_fr = phi_fr * sigma_eff_ice;
625 crv.J_uT_fr = phi_fr * C_IR * ice_linear_thermal_expansion_coefficient;
628 dC_eff_dT += -dphi_fr_dT * fluid_density * crv.c_f -
629 phi_fr * crv.drho_LR_dT * crv.c_f +
630 dphi_fr_dT * rho_fr * c_fr - l_fr * rho_fr * d2phi_fr_dT2;
631 double const storage_p_fr_coeff =
632 (porosity * crv.beta_IR + (alpha - porosity) * crv.beta_SR) *
633 rho_fr / fluid_density -
634 (porosity * crv.fluid_compressibility +
635 (alpha - porosity) * crv.beta_SR);
636 crv.storage_p_fr = phi_fr / porosity * storage_p_fr_coeff;
638 double const dstorage_p_fr_coeff_dT =
639 (porosity * crv.beta_IR + (alpha - porosity) * crv.beta_SR) *
640 crv.drho_LR_dT / (fluid_density * fluid_density);
642 crv.J_pT_fr = (dphi_fr_dT * storage_p_fr_coeff +
643 phi_fr * dstorage_p_fr_coeff_dT) /
644 porosity * dp_int_pt / dt;
648 (crv.beta_T_SI * rho_fr / fluid_density - crv.beta) -
650 double const dstorage_T_fr_dT =
651 dphi_fr_dT / porosity *
652 (crv.beta_T_SI * rho_fr / fluid_density - crv.beta) +
653 phi_fr / porosity * crv.beta_T_SI * rho_fr * crv.drho_LR_dT /
654 (fluid_density * fluid_density) -
656 crv.J_pT_fr += dstorage_T_fr_dT * dT_int_pt / dt;
658 crv.J_TT = dC_eff_dT * dT_int_pt / dt;
669 std::vector<double>
const& local_x,
670 std::vector<double>
const&
672 std::vector<double>& local_rhs_data,
673 std::vector<double>& local_Jac_data)
675 assert(local_x.size() ==
679 Eigen::Map<Eigen::VectorXd const>(local_x.data(), local_x.size());
680 auto const x_prev = Eigen::Map<Eigen::VectorXd const>(local_x_prev.data(),
681 local_x_prev.size());
683 auto const [T, p, u] =
localDOF(local_x);
684 auto const [T_prev, p_prev, u_prev] =
localDOF(local_x_prev);
687 typename ShapeMatricesTypeDisplacement::template MatrixType<
694 typename ShapeMatricesTypeDisplacement::template VectorType<
725 typename ShapeMatricesTypeDisplacement::template MatrixType<
731 bool const has_frozen_liquid_phase = medium->hasPhase(
"FrozenLiquid");
733 unsigned const n_integration_points =
736 std::vector<GlobalDimVectorType> ip_flux_vector;
737 double average_velocity_norm = 0.0;
738 ip_flux_vector.reserve(n_integration_points);
740 for (
unsigned ip = 0; ip < n_integration_points; ip++)
743 auto const& N_u = ip_data.N_u;
754 auto const& w = ip_data.integration_weight;
756 auto const& dNdx_u = ip_data.dNdx_u;
758 auto const& N = ip_data.N;
759 auto const& dNdx = ip_data.dNdx;
761 auto const T_int_pt = N.dot(T);
767 ShapeFunctionDisplacement::NPOINTS,
778 if (has_frozen_liquid_phase)
781 .template block<displacement_size, displacement_size>(
783 .noalias() += B.transpose() * crv.J_uu_fr * B * w;
786 .noalias() -= B.transpose() * crv.r_u_fr * w;
789 .template block<displacement_size, temperature_size>(
791 .noalias() -= B.transpose() * crv.J_uT_fr * N * w;
795 .template block<displacement_size, displacement_size>(
797 .noalias() += B.transpose() * crv.C * B * w;
800 .template block<displacement_size, temperature_size>(
802 .noalias() -= B.transpose() * crv.C *
803 crv.solid_linear_thermal_expansion_coefficient * N *
807 .noalias() -= (B.transpose() * ip_data.sigma_eff -
808 N_u_op(N_u).transpose() * crv.rho * b) *
818 if (has_frozen_liquid_phase)
822 (crv.alpha_biot * ip_data.phi_fr / ip_data.porosity *
830 laplace_p.noalias() +=
831 dNdx.transpose() * crv.K_over_mu * dNdx * (crv.k_rel * w);
835 .noalias() += dNdx.transpose() * crv.K_over_mu * (dNdx * p) * N *
838 storage_p.noalias() +=
840 ((ip_data.porosity * crv.fluid_compressibility +
841 (crv.alpha_biot - ip_data.porosity) * crv.beta_SR) *
844 if (has_frozen_liquid_phase)
846 storage_p.noalias() += N.transpose() * crv.storage_p_fr * N * w;
849 .template block<pressure_size, temperature_size>(
851 .noalias() += N.transpose() * crv.J_pT_fr * N * w;
854 laplace_T.noalias() +=
855 dNdx.transpose() * crv.K_pT_thermal_osmosis * dNdx * w;
859 local_rhs.template segment<pressure_size>(
pressure_index).noalias() +=
860 dNdx.transpose() * crv.K_over_mu * b *
861 (fluid_density * crv.k_rel * w);
866 dNdx.transpose() * crv.K_over_mu * b * N *
867 ((fluid_density * crv.dk_rel_dT + crv.drho_LR_dT * crv.k_rel) * w);
871 storage_T.noalias() += N.transpose() * crv.beta * N * w;
873 if (has_frozen_liquid_phase)
875 storage_T.noalias() += N.transpose() * crv.storage_T_fr * N * w;
887 dNdx.transpose() * crv.effective_thermal_conductivity * dNdx * w;
888 dKTT_dT_T.noalias() +=
889 dNdx.transpose() * crv.dlambda_eff_dT * dNdx * T * N * w;
891 ip_flux_vector.emplace_back(velocity * fluid_density * crv.c_f);
896 crv.dvelocity_dT * fluid_density * crv.c_f +
897 velocity * crv.drho_LR_dT * crv.c_f;
898 dKTT_dT_T.noalias() +=
899 N.transpose() * dip_flux_vector_dT.transpose() * dNdx * T * N * w;
900 average_velocity_norm += velocity.norm();
903 N.transpose() * crv.effective_volumetric_heat_capacity * N * w;
905 .template block<temperature_size, temperature_size>(
907 .noalias() += N.transpose() * crv.J_TT * N * w;
913 dNdx.transpose() * crv.K_pT_thermal_osmosis * dNdx * (T_int_pt * w);
916 dKTT_dp_T.noalias() -= N.transpose() * (dNdx * T).transpose() *
917 crv.K_over_mu * dNdx *
918 (fluid_density * crv.c_f * crv.k_rel * w);
965 average_velocity_norm /
static_cast<double>(n_integration_points), KTT);
971 .noalias() += KTT + dKTT_dT_T + MTT / dt;
977 .noalias() += KTp + dKTT_dp_T;
989 .noalias() += -storage_T / dt + laplace_T;
995 .noalias() += laplace_p + storage_p / dt;
1001 .noalias() += Kup.transpose() / dt;
1004 local_rhs.template segment<pressure_size>(
pressure_index).noalias() -=
1005 laplace_p * p + laplace_T * T + storage_p * (p - p_prev) / dt -
1006 storage_T * (T - T_prev) / dt + Kup.transpose() * (u - u_prev) / dt;
1010 .noalias() += Kup * p;
1014 KTT * T + MTT * (T - T_prev) / dt;