12 #include "MooseMesh.h"
16 #include "libmesh/quadrature.h"
24 InputParameters params = validParams<Material>();
26 params.addRequiredParam<Real>(
27 "mat_porosity",
"The porosity of the material. Should be between 0 and 1. Eg, 0.1");
28 params.addCoupledVar(
"por_change",
30 "An auxillary variable describing porosity changes. "
31 "Porosity = mat_porosity + por_change. If this is not "
32 "provided, zero is used.");
33 params.addRequiredParam<RealTensorValue>(
"mat_permeability",
"The permeability tensor (m^2).");
34 params.addCoupledVar(
"perm_change",
35 "A list of auxillary variable describing permeability "
36 "changes. There must be 9 of these, corresponding to the "
37 "xx, xy, xz, yx, yy, yz, zx, zy, zz components respectively. "
38 " Permeability = mat_permeability*10^(perm_change).");
39 params.addRequiredParam<UserObjectName>(
40 "richardsVarNames_UO",
"The UserObject that holds the list of Richards variable names.");
41 params.addRequiredParam<std::vector<UserObjectName>>(
42 "relperm_UO",
"List of names of user objects that define relative permeability");
43 params.addRequiredParam<std::vector<UserObjectName>>(
45 "List of name of user objects that define effective saturation as a function of "
47 params.addRequiredParam<std::vector<UserObjectName>>(
49 "List of names of user objects that define saturation as a function of effective saturation");
50 params.addRequiredParam<std::vector<UserObjectName>>(
51 "density_UO",
"List of names of user objects that define the fluid density");
52 params.addRequiredParam<std::vector<UserObjectName>>(
53 "SUPG_UO",
"List of names of user objects that define the SUPG");
54 params.addRequiredParam<std::vector<Real>>(
55 "viscosity",
"List of viscosity of fluids (Pa.s). Typical value for water is=1E-3");
56 params.addRequiredParam<RealVectorValue>(
58 "Gravitational acceleration (m/s^2) as a vector pointing downwards. Eg (0,0,-10)");
60 params.addParam<
bool>(
63 "If you are using second-order Lagrange shape functions you need to set this to false.");
68 : Material(parameters),
70 _material_por(getParam<Real>(
"mat_porosity")),
71 _por_change(coupledValue(
"por_change")),
72 _por_change_old(coupledValueOld(
"por_change")),
74 _material_perm(getParam<RealTensorValue>(
"mat_permeability")),
76 _material_gravity(getParam<RealVectorValue>(
"gravity")),
78 _porosity_old(declareProperty<Real>(
"porosity_old")),
79 _porosity(declareProperty<Real>(
"porosity")),
80 _permeability(declareProperty<RealTensorValue>(
"permeability")),
81 _gravity(declareProperty<RealVectorValue>(
"gravity")),
84 _num_p(_richards_name_UO.num_v()),
86 _trace_perm(_material_perm.tr()),
88 _material_viscosity(getParam<std::vector<Real>>(
"viscosity")),
90 _pp_old(declareProperty<std::vector<Real>>(
"porepressure_old")),
91 _pp(declareProperty<std::vector<Real>>(
"porepressure")),
92 _dpp_dv(declareProperty<std::vector<std::vector<Real>>>(
"dporepressure_dv")),
93 _d2pp_dv(declareProperty<std::vector<std::vector<std::vector<Real>>>>(
"d2porepressure_dvdv")),
95 _viscosity(declareProperty<std::vector<Real>>(
"viscosity")),
97 _density_old(declareProperty<std::vector<Real>>(
"density_old")),
98 _density(declareProperty<std::vector<Real>>(
"density")),
99 _ddensity_dv(declareProperty<std::vector<std::vector<Real>>>(
"ddensity_dv")),
101 _seff_old(declareProperty<std::vector<Real>>(
"s_eff_old")),
102 _seff(declareProperty<std::vector<Real>>(
"s_eff")),
103 _dseff_dv(declareProperty<std::vector<std::vector<Real>>>(
"ds_eff_dv")),
104 _d2seff_dv(declareProperty<std::vector<std::vector<std::vector<Real>>>>(
"d2s_eff_dvdv")),
106 _sat_old(declareProperty<std::vector<Real>>(
"sat_old")),
107 _sat(declareProperty<std::vector<Real>>(
"sat")),
108 _dsat_dv(declareProperty<std::vector<std::vector<Real>>>(
"dsat_dv")),
110 _rel_perm(declareProperty<std::vector<Real>>(
"rel_perm")),
111 _drel_perm_dv(declareProperty<std::vector<std::vector<Real>>>(
"drel_perm_dv")),
113 _mass_old(declareProperty<std::vector<Real>>(
"mass_old")),
114 _mass(declareProperty<std::vector<Real>>(
"mass")),
115 _dmass(declareProperty<std::vector<std::vector<Real>>>(
"dmass")),
117 _flux_no_mob(declareProperty<std::vector<RealVectorValue>>(
"flux_no_mob")),
118 _dflux_no_mob_dv(declareProperty<std::vector<std::vector<RealVectorValue>>>(
"dflux_no_mob_dv")),
119 _dflux_no_mob_dgradv(
120 declareProperty<std::vector<std::vector<RealTensorValue>>>(
"dflux_no_mob_dgradv")),
122 _flux(declareProperty<std::vector<RealVectorValue>>(
"flux")),
123 _dflux_dv(declareProperty<std::vector<std::vector<RealVectorValue>>>(
"dflux_dv")),
124 _dflux_dgradv(declareProperty<std::vector<std::vector<RealTensorValue>>>(
"dflux_dgradv")),
126 declareProperty<std::vector<std::vector<std::vector<RealVectorValue>>>>(
"d2flux_dvdv")),
128 declareProperty<std::vector<std::vector<std::vector<RealTensorValue>>>>(
"d2flux_dgradvdv")),
130 declareProperty<std::vector<std::vector<std::vector<RealTensorValue>>>>(
"d2flux_dvdgradv")),
132 _tauvel_SUPG(declareProperty<std::vector<RealVectorValue>>(
"tauvel_SUPG")),
133 _dtauvel_SUPG_dgradp(
134 declareProperty<std::vector<std::vector<RealTensorValue>>>(
"dtauvel_SUPG_dgradv")),
135 _dtauvel_SUPG_dp(declareProperty<std::vector<std::vector<RealVectorValue>>>(
"dtauvel_SUPG_dv"))
142 const std::vector<MooseVariableFEBase *> & coupled_vars =
144 for (
unsigned int i = 0; i < coupled_vars.size(); i++)
145 addMooseVariableDependency(coupled_vars[i]);
148 if (_material_por <= 0 || _material_por >= 1)
149 mooseError(
"Porosity set to ",
_material_por,
" but it must be between 0 and 1");
151 if (isCoupled(
"perm_change") && (coupledComponents(
"perm_change") != LIBMESH_DIM * LIBMESH_DIM))
152 mooseError(LIBMESH_DIM * LIBMESH_DIM,
153 " components of perm_change must be given to a RichardsMaterial. You supplied ",
154 coupledComponents(
"perm_change"),
158 for (
unsigned int i = 0; i < LIBMESH_DIM * LIBMESH_DIM; ++i)
159 _perm_change[i] = (isCoupled(
"perm_change") ? &coupledValue(
"perm_change", i)
165 getParam<std::vector<UserObjectName>>(
"relperm_UO").size() ==
_num_p &&
166 getParam<std::vector<UserObjectName>>(
"seff_UO").size() ==
_num_p &&
167 getParam<std::vector<UserObjectName>>(
"sat_UO").size() ==
_num_p &&
168 getParam<std::vector<UserObjectName>>(
"density_UO").size() ==
_num_p &&
169 getParam<std::vector<UserObjectName>>(
"SUPG_UO").size() ==
_num_p))
170 mooseError(
"There are ",
172 " Richards fluid variables, so you need to specify this "
173 "number of viscosities, relperm_UO, seff_UO, sat_UO, "
174 "density_UO, SUPG_UO");
187 for (
unsigned int i = 0; i <
_num_p; ++i)
198 getParam<std::vector<UserObjectName>>(
"relperm_UO")[i]);
200 &getUserObjectByName<RichardsSeff>(getParam<std::vector<UserObjectName>>(
"seff_UO")[i]);
202 &getUserObjectByName<RichardsSat>(getParam<std::vector<UserObjectName>>(
"sat_UO")[i]);
204 getParam<std::vector<UserObjectName>>(
"density_UO")[i]);
206 &getUserObjectByName<RichardsSUPG>(getParam<std::vector<UserObjectName>>(
"SUPG_UO")[i]);
217 for (
unsigned int i = 0; i <
_num_p; ++i)
225 for (
unsigned int qp = 0; qp < _qrule->n_points(); qp++)
239 for (
unsigned int i = 0; i <
_num_p; ++i)
248 for (
unsigned int j = 0; j <
_num_p; ++j)
258 for (
unsigned int j = 0; j <
_num_p; ++j)
274 for (
unsigned int i = 0; i <
_num_p; ++i)
281 for (
unsigned int i = 0; i <
_num_p; ++i)
286 for (
unsigned int j = 0; j <
_num_p; ++j)
294 for (
unsigned int i = 0; i <
_num_p; ++i)
299 for (
unsigned int j = 0; j <
_num_p; ++j)
306 for (
unsigned int i = 0; i <
_num_p; ++i)
310 for (
unsigned int j = 0; j <
_num_p; ++j)
318 for (
unsigned int i = 0; i <
_num_p; ++i)
323 for (
unsigned int j = 0; j <
_num_p; ++j)
332 for (
unsigned int i = 0; i <
_num_p; ++i)
337 for (
unsigned int j = 0; j <
_num_p; ++j)
341 for (
unsigned int j = 0; j <
_num_p; ++j)
349 for (
unsigned int i = 0; i <
_num_p; ++i)
354 for (
unsigned int j = 0; j <
_num_p; ++j)
364 for (
unsigned int j = 0; j <
_num_p; ++j)
377 for (
unsigned int i = 0; i <
_num_p; ++i)
382 for (
unsigned int j = 0; j <
_num_p; ++j)
396 for (
unsigned int i = 0; i <
_num_p; ++i)
405 for (
unsigned int j = 0; j <
_num_p; ++j)
408 for (
unsigned int k = 0; k <
_num_p; ++k)
417 for (
unsigned int j = 0; j <
_num_p; ++j)
420 for (
unsigned int k = 0; k <
_num_p; ++k)
426 for (
unsigned int j = 0; j <
_num_p; ++j)
428 for (
unsigned int k = 0; k <
_num_p; ++k)
450 for (
unsigned int j = 0; j <
_num_p; ++j)
451 for (
unsigned int k = 0; k <
_num_p; ++k)
454 for (
unsigned int j = 0; j <
_num_p; ++j)
456 for (
unsigned int k = 0; k <
_num_p; ++k)
473 for (
unsigned int i = 0; i <
_num_p; ++i)
486 auto && fe = _assembly.getFE(getParam<bool>(
"linear_shape_fcns") ? FEType(FIRST, LAGRANGE)
487 : FEType(SECOND, LAGRANGE),
488 _current_elem->dim());
491 const std::vector<Real> & dxidx(fe->get_dxidx());
492 const std::vector<Real> & dxidy(fe->get_dxidy());
493 const std::vector<Real> & dxidz(fe->get_dxidz());
494 const std::vector<Real> & detadx(fe->get_detadx());
495 const std::vector<Real> & detady(fe->get_detady());
496 const std::vector<Real> & detadz(fe->get_detadz());
497 const std::vector<Real> & dzetadx(fe->get_dzetadx());
498 const std::vector<Real> & dzetady(fe->get_dzetady());
499 const std::vector<Real> & dzetadz(fe->get_dzetadz());
501 for (
unsigned int qp = 0; qp < _qrule->n_points(); qp++)
505 mooseAssert(qp < dxidx.size(),
"Insufficient data in dxidx array!");
506 mooseAssert(qp < dxidy.size(),
"Insufficient data in dxidy array!");
507 mooseAssert(qp < dxidz.size(),
"Insufficient data in dxidz array!");
508 if (_mesh.dimension() >= 2)
510 mooseAssert(qp < detadx.size(),
"Insufficient data in detadx array!");
511 mooseAssert(qp < detady.size(),
"Insufficient data in detady array!");
512 mooseAssert(qp < detadz.size(),
"Insufficient data in detadz array!");
514 if (_mesh.dimension() >= 3)
516 mooseAssert(qp < dzetadx.size(),
"Insufficient data in dzetadx array!");
517 mooseAssert(qp < dzetady.size(),
"Insufficient data in dzetady array!");
518 mooseAssert(qp < dzetadz.size(),
"Insufficient data in dzetadz array!");
522 RealVectorValue xi_prime(dxidx[qp], dxidy[qp], dxidz[qp]);
523 RealVectorValue eta_prime, zeta_prime;
524 if (_mesh.dimension() >= 2)
526 eta_prime(0) = detadx[qp];
527 eta_prime(1) = detady[qp];
529 if (_mesh.dimension() == 3)
531 eta_prime(2) = detadz[qp];
532 zeta_prime(0) = dzetadx[qp];
533 zeta_prime(1) = dzetady[qp];
534 zeta_prime(2) = dzetadz[qp];
538 for (
unsigned int i = 0; i <
_num_p; ++i)
540 RealVectorValue vel =
544 RealVectorValue dvel_dp =
548 (*
_material_SUPG_UO[i]).bb(vel, _mesh.dimension(), xi_prime, eta_prime, zeta_prime);
549 RealVectorValue dbb2_dgradp =
550 (*
_material_SUPG_UO[i]).dbb2_dgradp(vel, dvel_dgradp, xi_prime, eta_prime, zeta_prime);
551 Real dbb2_dp = (*
_material_SUPG_UO[i]).dbb2_dp(vel, dvel_dp, xi_prime, eta_prime, zeta_prime);
553 RealVectorValue dtau_dgradp =
559 RealTensorValue dtauvel_dgradp = tau * dvel_dgradp;
560 for (
unsigned int j = 0; j < LIBMESH_DIM; ++j)
561 for (
unsigned int k = 0; k < LIBMESH_DIM; ++k)
562 dtauvel_dgradp(j, k) +=
563 dtau_dgradp(j) * vel(k);
564 for (
unsigned int j = 0; j <
_num_p; ++j)
567 RealVectorValue dtauvel_dp = dtau_dp * vel + tau * dvel_dp;
568 for (
unsigned int j = 0; j <
_num_p; ++j)
582 for (
unsigned int qp = 0; qp < _qrule->n_points(); qp++)
588 for (
unsigned int i = 0; i < LIBMESH_DIM; i++)
589 for (
unsigned int j = 0; j < LIBMESH_DIM; j++)
596 for (
unsigned int qp = 0; qp < _qrule->n_points(); qp++)
601 for (
unsigned int qp = 0; qp < _qrule->n_points(); qp++)
605 for (
unsigned int qp = 0; qp < _qrule->n_points(); qp++)
609 bool trivial_supg =
true;
610 for (
unsigned int i = 0; i <
_num_p; ++i)