https://mooseframework.inl.gov
Functions
FlowModelGasMixUtils Namespace Reference

Functions

template<bool is_ad>
std::vector< GenericReal< is_ad > > computePrimitiveSolution (const std::vector< GenericReal< is_ad >> &U, const VaporMixtureFluidProperties &fp)
 Computes the primitive solution vector from the conservative solution vector. More...
 
template<bool is_ad>
std::vector< GenericReal< is_ad > > computeConservativeSolution (const std::vector< GenericReal< is_ad >> &W, const GenericReal< is_ad > &A, const VaporMixtureFluidProperties &fp)
 Computes the conservative solution vector from the primitive solution vector. More...
 
template<bool is_ad>
std::vector< GenericReal< is_ad > > computeFluxFromPrimitive (const std::vector< GenericReal< is_ad >> &W, const GenericReal< is_ad > &A, const VaporMixtureFluidProperties &fp)
 Computes the numerical flux vector from the primitive solution vector. More...
 
template<bool is_ad>
std::vector< GenericReal< is_ad > > getElementalSolutionVector (const Elem *elem, const std::vector< MooseVariable *> &U_vars, bool is_implicit)
 Gets the elemental conservative solution vector. More...
 
template<bool is_ad>
GenericReal< is_ad > computeSecondaryMoleFraction (const GenericReal< is_ad > &xi_secondary, const VaporMixtureFluidProperties &fp)
 
ADReal computeSecondaryMoleFraction (const ADReal &xi_secondary, const VaporMixtureFluidProperties &fp)
 

Function Documentation

◆ computeConservativeSolution()

template<bool is_ad>
std::vector<GenericReal<is_ad> > FlowModelGasMixUtils::computeConservativeSolution ( const std::vector< GenericReal< is_ad >> &  W,
const GenericReal< is_ad > &  A,
const VaporMixtureFluidProperties fp 
)

Computes the conservative solution vector from the primitive solution vector.

Parameters
[in]WPrimitive solution vector
[in]ACross-sectional area
[in]fpFluid properties object

Definition at line 64 of file FlowModelGasMixUtils.h.

67 {
68  const auto & xi = W[THMGasMix1D::MASS_FRACTION];
69  const auto & p = W[THMGasMix1D::PRESSURE];
70  const auto & T = W[THMGasMix1D::TEMPERATURE];
71  const auto & vel = W[THMGasMix1D::VELOCITY];
72 
73  const auto rho = fp.rho_from_p_T(p, T, {xi});
74  const auto rhoA = rho * A;
75  const auto e = fp.e_from_p_rho(p, rho, {xi});
76  const auto E = e + 0.5 * vel * vel;
77 
78  std::vector<GenericReal<is_ad>> U(THMGasMix1D::N_FLUX_INPUTS);
79  U[THMGasMix1D::XIRHOA] = xi * rhoA;
80  U[THMGasMix1D::RHOA] = rhoA;
81  U[THMGasMix1D::RHOUA] = rhoA * vel;
82  U[THMGasMix1D::RHOEA] = rhoA * E;
83  U[THMGasMix1D::AREA] = A;
84 
85  return U;
86 }
static const unsigned int N_FLUX_INPUTS
Number of numerical flux function inputs.

◆ computeFluxFromPrimitive()

template<bool is_ad>
std::vector<GenericReal<is_ad> > FlowModelGasMixUtils::computeFluxFromPrimitive ( const std::vector< GenericReal< is_ad >> &  W,
const GenericReal< is_ad > &  A,
const VaporMixtureFluidProperties fp 
)

Computes the numerical flux vector from the primitive solution vector.

Parameters
[in]WPrimitive solution vector
[in]ACross-sectional area
[in]fpFluid properties object

Definition at line 97 of file FlowModelGasMixUtils.h.

100 {
101  const auto & xi = W[THMGasMix1D::MASS_FRACTION];
102  const auto & p = W[THMGasMix1D::PRESSURE];
103  const auto & T = W[THMGasMix1D::TEMPERATURE];
104  const auto & vel = W[THMGasMix1D::VELOCITY];
105 
106  const auto rho = fp.rho_from_p_T(p, T, {xi});
107  const auto e = fp.e_from_p_rho(p, rho, {xi});
108  const auto E = e + 0.5 * vel * vel;
109 
110  std::vector<ADReal> F(THMGasMix1D::N_FLUX_OUTPUTS, 0.0);
111  F[THMGasMix1D::SPECIES] = xi * rho * vel * A;
112  F[THMGasMix1D::MASS] = rho * vel * A;
113  F[THMGasMix1D::MOMENTUM] = (rho * vel * vel + p) * A;
114  F[THMGasMix1D::ENERGY] = vel * (rho * E + p) * A;
115 
116  return F;
117 }
static const std::string F
Definition: NS.h:165
static const unsigned int N_FLUX_OUTPUTS
Number of numerical flux function outputs.

◆ computePrimitiveSolution()

template<bool is_ad>
std::vector<GenericReal<is_ad> > FlowModelGasMixUtils::computePrimitiveSolution ( const std::vector< GenericReal< is_ad >> &  U,
const VaporMixtureFluidProperties fp 
)

Computes the primitive solution vector from the conservative solution vector.

Parameters
[in]UConservative solution vector
[in]fpFluid properties object

Definition at line 29 of file FlowModelGasMixUtils.h.

31 {
32  const auto & xirhoA = U[THMGasMix1D::XIRHOA];
33  const auto & rhoA = U[THMGasMix1D::RHOA];
34  const auto & rhouA = U[THMGasMix1D::RHOUA];
35  const auto & rhoEA = U[THMGasMix1D::RHOEA];
36  const auto & A = U[THMGasMix1D::AREA];
37 
38  const auto xi = xirhoA / rhoA;
39  const auto rho = rhoA / A;
40  const auto vel = rhouA / rhoA;
41  const auto v = 1.0 / rho;
42  const auto e = rhoEA / rhoA - 0.5 * vel * vel;
43  const auto p = fp.p_from_v_e(v, e, {xi});
44  const auto T = fp.T_from_v_e(v, e, {xi});
45 
46  std::vector<GenericReal<is_ad>> W(THMGasMix1D::N_PRIM_VARS);
48  W[THMGasMix1D::PRESSURE] = p;
49  W[THMGasMix1D::VELOCITY] = vel;
51 
52  return W;
53 }
static const unsigned int N_PRIM_VARS
static const std::string v
Definition: NS.h:84

◆ computeSecondaryMoleFraction() [1/2]

ADReal FlowModelGasMixUtils::computeSecondaryMoleFraction ( const ADReal xi_secondary,
const VaporMixtureFluidProperties fp 
)

Definition at line 18 of file FlowModelGasMixUtils.C.

19 {
20  mooseAssert(fp.getNumberOfSecondaryVapors() == 1,
21  "This function assumes there is a single secondary fluid.");
22  const SinglePhaseFluidProperties & fp_primary = fp.getPrimaryFluidProperties();
23  const SinglePhaseFluidProperties & fp_secondary = fp.getSecondaryFluidProperties();
24 
25  const ADReal xi_primary = 1 - xi_secondary;
26 
27  const ADReal moles_primary = xi_primary / fp_primary.molarMass();
28  const ADReal moles_secondary = xi_secondary / fp_secondary.molarMass();
29 
30  return moles_secondary / (moles_primary + moles_secondary);
31 }
virtual Real molarMass() const
Molar mass [kg/mol].
DualNumber< Real, DNDerivativeType, true > ADReal
Common class for single phase fluid properties.

◆ computeSecondaryMoleFraction() [2/2]

template<bool is_ad>
GenericReal<is_ad> FlowModelGasMixUtils::computeSecondaryMoleFraction ( const GenericReal< is_ad > &  xi_secondary,
const VaporMixtureFluidProperties fp 
)

Definition at line 166 of file FlowModelGasMixUtils.h.

168 {
169  mooseAssert(fp.getNumberOfSecondaryVapors() == 1,
170  "This function assumes there is a single secondary fluid.");
171  const SinglePhaseFluidProperties & fp_primary = fp.getPrimaryFluidProperties();
172  const SinglePhaseFluidProperties & fp_secondary = fp.getSecondaryFluidProperties();
173 
174  const GenericReal<is_ad> xi_primary = 1 - xi_secondary;
175 
176  const GenericReal<is_ad> moles_primary = xi_primary / fp_primary.molarMass();
177  const GenericReal<is_ad> moles_secondary = xi_secondary / fp_secondary.molarMass();
178 
179  return moles_secondary / (moles_primary + moles_secondary);
180 }
Moose::GenericType< Real, is_ad > GenericReal
virtual Real molarMass() const
Molar mass [kg/mol].
Common class for single phase fluid properties.

◆ getElementalSolutionVector()

template<bool is_ad>
std::vector<GenericReal<is_ad> > FlowModelGasMixUtils::getElementalSolutionVector ( const Elem *  elem,
const std::vector< MooseVariable *> &  U_vars,
bool  is_implicit 
)

Gets the elemental conservative solution vector.

Parameters
[in]elemElement
[in]U_varsVector of conservative variable pointers
[in]is_implicitIs implicit?

Definition at line 128 of file FlowModelGasMixUtils.h.

131 {
132  mooseAssert(elem, "The supplied element is a nullptr.");
133 
134  std::vector<GenericReal<is_ad>> U(THMGasMix1D::N_FLUX_INPUTS, 0.0);
135 
136  if (is_implicit)
137  {
138  for (const auto i : make_range(THMGasMix1D::N_FLUX_INPUTS))
139  {
140  mooseAssert(U_vars[i], "The supplied variable is a nullptr.");
141  U[i] = U_vars[i]->getElementalValue(elem);
142  }
143 
144  std::vector<dof_id_type> dof_indices;
145 
146  const std::vector<unsigned int> ind = {
148  for (const auto j : index_range(ind))
149  {
150  const auto i = ind[j];
151  U_vars[i]->dofMap().dof_indices(elem, dof_indices, U_vars[i]->number());
152  Moose::derivInsert(U[i].derivatives(), dof_indices[0], 1.0);
153  }
154  }
155  else
156  {
157  for (const auto i : make_range(THMGasMix1D::N_FLUX_INPUTS))
158  U[i] = U_vars[i]->getElementalValueOld(elem);
159  }
160 
161  return U;
162 }
IntRange< T > make_range(T beg, T end)
void derivInsert(SemiDynamicSparseNumberArray< Real, libMesh::dof_id_type, NWrapper< N >> &derivs, libMesh::dof_id_type index, Real value)
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
static const unsigned int N_FLUX_INPUTS
Number of numerical flux function inputs.
auto index_range(const T &sizable)