https://mooseframework.inl.gov
Numerics.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include "libmesh/libmesh_common.h"
13 #include "libmesh/vector_value.h"
14 #include "libmesh/dense_vector.h"
15 
16 #include "ADReal.h"
17 #include "NavierStokesMethods.h"
18 #include "HeatTransferUtils.h"
19 
20 using namespace libMesh;
21 
22 namespace THM
23 {
24 
25 // Default value for magnitude of acceleration due to gravity
26 static const Real gravity_const = 9.81;
27 
28 // Default value for gravitational acceleration vector
30 
31 // Stefan-Boltzman constant, in [W/m^2-K]
32 static const Real Stefan_Boltzman_const = 5.670e-8;
33 
39 template <typename T>
40 int
41 sgn(T val)
42 {
43  return (T(0) < val) - (val < T(0));
44 }
45 
54  const RealVectorValue & b,
56 
65  const RealVectorValue & b,
67 
76  const RealVectorValue & b,
78 
88 Real
89 applyQuotientRule(const Real & num, const Real & den, const Real & dnum_dy, const Real & dden_dy);
90 
101  const Real & den,
102  const DenseVector<Real> & dnum_dy,
103  const DenseVector<Real> & dden_dy);
104 
116 template <typename T1, typename T2, typename T3, typename T4, typename T5>
117 auto
118 Reynolds(const T1 & volume_fraction, const T2 & rho, const T3 & vel, const T4 & D_h, const T5 & mu)
119 {
120  return volume_fraction * HeatTransferUtils::reynolds(rho, vel, D_h, mu);
121 }
122 
131 template <typename T1, typename T2, typename T3>
132 auto
133 Prandtl(const T1 & cp, const T2 & mu, const T3 & k)
134 {
135  return HeatTransferUtils::prandtl(cp, mu, k);
136 }
137 
151 template <typename T1, typename T2, typename T3, typename T4, typename T5, typename T6>
152 auto
153 Peclet(const T1 & volume_fraction,
154  const T2 & cp,
155  const T3 & rho,
156  const T4 & vel,
157  const T5 & D_h,
158  const T6 & k)
159 {
160  const auto diffusivity = HeatTransferUtils::thermalDiffusivity(k, rho, cp);
161  return volume_fraction * HeatTransferUtils::peclet(vel, D_h, diffusivity);
162 }
163 
176 template <typename T1, typename T2, typename T3, typename T4, typename T5>
177 auto
178 Grashof(const T1 & beta,
179  const T2 & dT,
180  const T3 & D_h,
181  const T4 & rho_liquid,
182  const T5 & mu_liquid,
183  const Real & gravity_magnitude)
184 {
185  return gravity_magnitude * beta * dT * std::pow(D_h, 3) * (rho_liquid * rho_liquid) /
186  (mu_liquid * mu_liquid);
187 }
188 
198 template <typename T1, typename T2>
199 auto
200 Laplace(const T1 & surf_tension, const T2 & delta_rho, const Real & gravity_magnitude)
201 {
202  return std::sqrt(surf_tension / (gravity_magnitude * delta_rho));
203 }
204 
216 template <typename T1, typename T2, typename T3, typename T4>
217 auto
218 viscosityNumber(const T1 & viscosity,
219  const T2 & surf_tension,
220  const T3 & rho_k,
221  const T4 & delta_rho,
222  const Real & gravity_magnitude)
223 {
224  return viscosity /
225  std::sqrt(rho_k * surf_tension * std::sqrt(surf_tension / gravity_magnitude / delta_rho));
226 }
227 
229 
237 template <typename T1, typename T2>
238 auto
239 Dean(const T1 & Re, const T2 & doD)
240 {
241  return Re * std::sqrt(doD);
242 }
243 
253 void
255 
264 
272 Real dvel_darhoA(Real arhoA, Real arhouA);
273 
280 Real dvel_darhouA(Real arhoA);
281 
293  Real arhoA, Real alpha, Real A, Real & rho, Real & drho_darhoA, Real & drho_dalpha);
294 
304 
312 void v_from_rhoA_A(Real rhoA, Real A, Real & v, Real & dv_drhoA);
313 
322 
332 void
333 v_from_arhoA_alpha_A(Real arhoA, Real alpha, Real A, Real & v, Real & dv_darhoA, Real & dv_dalpha);
334 
344 
352 void v_from_rho(Real rho, Real & v, Real & dv_drho);
353 
362 Real dv_dalpha_liquid(Real area, Real arhoA, bool is_liquid);
363 
370 Real dv_darhoA(Real area, Real arhoA);
371 
385  Real arhouA,
386  Real arhoEA,
387  Real & e,
388  Real & de_darhoA,
389  Real & de_darhouA,
390  Real & de_darhoEA);
391 
392 ADReal e_from_arhoA_arhouA_arhoEA(ADReal arhoA, ADReal arhouA, ADReal arhoEA);
393 
403 void e_from_E_vel(Real E, Real vel, Real & e, Real & de_dE, Real & de_dvel);
404 
413 
421 Real de_darhoA(Real arhoA, Real arhouA, Real arhoEA);
422 
429 Real de_darhouA(Real arhoA, Real arhouA);
430 
436 Real de_darhoEA(Real arhoA);
437 
447 void E_from_arhoA_arhoEA(Real arhoA, Real arhoEA, Real & E, Real & dE_darhoA, Real & dE_darhoEA);
448 
456 ADReal E_from_arhoA_arhoEA(ADReal arhoA, ADReal arhoEA);
457 
467 void E_from_e_vel(Real e, Real vel, Real & E, Real & dE_de, Real & dE_dvel);
468 
481 void h_from_e_p_rho(Real e, Real p, Real rho, Real & h, Real & dh_de, Real & dh_dp, Real & dh_drho);
482 
484 
492 bool isInlet(Real vel, Real normal);
493 bool isInlet(ADReal vel, Real normal);
494 
502 bool isOutlet(Real vel, Real normal);
503 bool isOutlet(ADReal vel, Real normal);
504 }
ADReal v_from_rhoA_A(ADReal rhoA, ADReal A)
Computes specific volume and its derivatives from rho*A, and area.
Definition: Numerics.C:107
ADReal rho_from_arhoA_alpha_A(ADReal arhoA, ADReal alpha, ADReal A)
Computes density from alpha*rho*A, alpha, and area.
Definition: Numerics.C:94
Real de_darhoA(Real arhoA, Real arhouA, Real arhoEA)
Derivative of specific internal energy wrt density of the phase (rhoA or arhoA)
Definition: Numerics.C:181
Real de_darhoEA(Real arhoA)
Derivative of specific internal energy wrt total energy of the phase (rhoEA or arhoEA) ...
Definition: Numerics.C:193
ADReal e_from_E_vel(ADReal E, ADReal vel)
Computes specific internal energy from specific total energy and velocity.
Definition: Numerics.C:175
static constexpr Real TOLERANCE
bool isOutlet(ADReal vel, Real normal)
Definition: Numerics.C:253
const double tol
auto Prandtl(const T1 &cp, const T2 &mu, const T3 &k)
Compute Prandtl number.
Definition: Numerics.h:133
ADReal vel_from_arhoA_arhouA(ADReal arhoA, ADReal arhouA)
Computes velocity from alpha*rho*A and alpha*rho*u*A.
Definition: Numerics.C:66
ADReal v_from_arhoA_alpha_A(ADReal arhoA, ADReal alpha, ADReal A)
Computes specific volume and its derivatives from alpha*rho*A, volume fraction, and area...
Definition: Numerics.C:121
auto Reynolds(const T1 &volume_fraction, const T2 &rho, const T3 &vel, const T4 &D_h, const T5 &mu)
Compute Reynolds number.
Definition: Numerics.h:118
int sgn(T val)
The sign function.
Definition: Numerics.h:41
Real dv_darhoA(Real area, Real arhoA)
Derivative of specific volume wrt density equation solution variable.
Definition: Numerics.C:141
auto Grashof(const T1 &beta, const T2 &dT, const T3 &D_h, const T4 &rho_liquid, const T5 &mu_liquid, const Real &gravity_magnitude)
Compute Grashof number.
Definition: Numerics.h:178
auto peclet(const T1 &vel, const T2 &L, const T3 &diffusivity)
Compute Peclet number.
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
bool areParallelVectors(const RealVectorValue &a, const RealVectorValue &b, const Real &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
Tests if two real-valued vectors are parallel within some absolute tolerance.
Definition: Numerics.C:26
auto thermalDiffusivity(const T1 &k, const T2 &rho, const T3 &cp)
Compute thermal diffusivity.
auto viscosityNumber(const T1 &viscosity, const T2 &surf_tension, const T3 &rho_k, const T4 &delta_rho, const Real &gravity_magnitude)
Compute viscosity number (or coefficient)
Definition: Numerics.h:218
static const std::string cp
Definition: NS.h:121
ADReal e_from_arhoA_arhouA_arhoEA(ADReal arhoA, ADReal arhouA, ADReal arhoEA)
Definition: Numerics.C:162
auto Laplace(const T1 &surf_tension, const T2 &delta_rho, const Real &gravity_magnitude)
Compute Laplace number (or coefficient)
Definition: Numerics.h:200
static const std::string mu
Definition: NS.h:123
Real de_darhouA(Real arhoA, Real arhouA)
Derivative of specific internal energy wrt momentum of the phase (rhouA or arhouA) ...
Definition: Numerics.C:187
ADReal h_from_e_p_rho(ADReal e, ADReal p, ADReal rho)
Definition: Numerics.C:229
auto Peclet(const T1 &volume_fraction, const T2 &cp, const T3 &rho, const T4 &vel, const T5 &D_h, const T6 &k)
Compute Peclet number.
Definition: Numerics.h:153
Real dvel_darhouA(Real arhoA)
Derivative of velocity w.r.t.
Definition: Numerics.C:78
auto wallHeatTransferCoefficient(const T1 &Nu, const T2 &k, const T3 &D_h)
Compute wall heat transfer coefficient.
auto Dean(const T1 &Re, const T2 &doD)
Compute Dean number.
Definition: Numerics.h:239
auto reynolds(const T1 &rho, const T2 &vel, const T3 &L, const T4 &mu)
Compute Reynolds number.
auto prandtl(const T1 &cp, const T2 &mu, const T3 &k)
Compute Prandtl number.
bool haveSameDirection(const RealVectorValue &a, const RealVectorValue &b, const Real &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
Tests if two real-valued vectors are in the same direction.
Definition: Numerics.C:33
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const std::string v
Definition: NS.h:84
static const std::string alpha
Definition: NS.h:134
void v_from_rho(Real rho, Real &v, Real &dv_drho)
Computes specific volume and its derivative with respect to density.
Definition: Numerics.C:127
void E_from_e_vel(Real e, Real vel, Real &E, Real &dE_de, Real &dE_dvel)
Computes specific total energy and its derivatives from specific internal energy and velocity...
Definition: Numerics.C:212
Real dv_dalpha_liquid(Real area, Real arhoA, bool is_liquid)
Derivative of specific volume wrt alpha_liquid.
Definition: Numerics.C:134
ADReal E_from_arhoA_arhoEA(ADReal arhoA, ADReal arhoEA)
Computes specific total energy from alpha*rho*A and alpha*rho*E*A.
Definition: Numerics.C:206
Real dvel_darhoA(Real arhoA, Real arhouA)
Derivative of velocity w.r.t.
Definition: Numerics.C:72
static const Real gravity_const
Definition: Numerics.h:26
MooseUnits pow(const MooseUnits &, int)
static const std::string k
Definition: NS.h:130
static VectorValue< Real > default_gravity_vector
Definition: Numerics.h:29
DenseVector< Real > applyQuotientRule(const Real &num, const Real &den, const DenseVector< Real > &dnum_dy, const DenseVector< Real > &dden_dy)
Computes a derivative of a fraction using quotient rule for a derivative w.r.t.
bool isInlet(ADReal vel, Real normal)
Definition: Numerics.C:241
static const Real Stefan_Boltzman_const
Definition: Numerics.h:32
bool absoluteFuzzyEqualVectors(const RealVectorValue &a, const RealVectorValue &b, const Real &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
Tests if two real-valued vectors are equal within some absolute tolerance.
Definition: Numerics.C:18