https://mooseframework.inl.gov
FlowModel1PhaseUtils.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 
13 #include "THMIndicesVACE.h"
14 #include "MooseVariable.h"
15 
16 #include "libmesh/elem.h"
17 
19 {
20 
27 template <bool is_ad>
28 std::vector<GenericReal<is_ad>>
31 {
32  const auto & rhoA = U[THMVACE1D::RHOA];
33  const auto & rhouA = U[THMVACE1D::RHOUA];
34  const auto & rhoEA = U[THMVACE1D::RHOEA];
35  const auto & A = U[THMVACE1D::AREA];
36 
37  const auto rho = rhoA / A;
38  const auto vel = rhouA / rhoA;
39  const auto v = 1.0 / rho;
40  const auto e = rhoEA / rhoA - 0.5 * vel * vel;
41  const auto p = fp.p_from_v_e(v, e);
42  const auto T = fp.T_from_v_e(v, e);
43 
44  std::vector<GenericReal<is_ad>> W(THMVACE1D::N_PRIM_VARS);
45  W[THMVACE1D::PRESSURE] = fp.p_from_v_e(v, e);
46  W[THMVACE1D::VELOCITY] = vel;
47  W[THMVACE1D::TEMPERATURE] = fp.T_from_v_e(v, e);
48 
49  return W;
50 }
51 
59 template <bool is_ad>
60 std::vector<GenericReal<is_ad>>
62  const GenericReal<is_ad> & A,
64 {
65  const auto & p = W[THMVACE1D::PRESSURE];
66  const auto & T = W[THMVACE1D::TEMPERATURE];
67  const auto & vel = W[THMVACE1D::VELOCITY];
68 
69  const ADReal rho = fp.rho_from_p_T(p, T);
70  const ADReal e = fp.e_from_p_rho(p, rho);
71  const ADReal E = e + 0.5 * vel * vel;
72 
73  std::vector<GenericReal<is_ad>> U(THMVACE1D::N_FLUX_INPUTS);
74  U[THMVACE1D::RHOA] = rho * A;
75  U[THMVACE1D::RHOUA] = U[THMVACE1D::RHOA] * vel;
77  U[THMVACE1D::AREA] = A;
78 
79  return U;
80 }
81 
89 template <bool is_ad>
90 std::vector<GenericReal<is_ad>>
92  const GenericReal<is_ad> & A,
94 {
95  const auto & p = W[THMVACE1D::PRESSURE];
96  const auto & T = W[THMVACE1D::TEMPERATURE];
97  const auto & vel = W[THMVACE1D::VELOCITY];
98 
99  const auto rho = fp.rho_from_p_T(p, T);
100  const auto e = fp.e_from_p_rho(p, rho);
101  const auto E = e + 0.5 * vel * vel;
102 
103  std::vector<ADReal> F(THMVACE1D::N_FLUX_OUTPUTS, 0.0);
104  F[THMVACE1D::MASS] = rho * vel * A;
105  F[THMVACE1D::MOMENTUM] = (rho * vel * vel + p) * A;
106  F[THMVACE1D::ENERGY] = vel * (rho * E + p) * A;
107 
108  return F;
109 }
110 
118 template <bool is_ad>
119 std::vector<GenericReal<is_ad>>
120 getElementalSolutionVector(const Elem * elem,
121  const std::vector<MooseVariable *> & U_vars,
122  bool is_implicit)
123 {
124  mooseAssert(elem, "The supplied element is a nullptr.");
125 
126  std::vector<GenericReal<is_ad>> U(THMVACE1D::N_FLUX_INPUTS, 0.0);
127 
128  if (is_implicit)
129  {
130  for (unsigned int i = 0; i < THMVACE1D::N_FLUX_INPUTS; i++)
131  {
132  mooseAssert(U_vars[i], "The supplied variable is a nullptr.");
133  U[i] = U_vars[i]->getElementalValue(elem);
134  }
135 
136  std::vector<dof_id_type> dof_indices;
137 
138  const std::vector<unsigned int> ind = {THMVACE1D::RHOA, THMVACE1D::RHOUA, THMVACE1D::RHOEA};
139  for (unsigned int j = 0; j < ind.size(); j++)
140  {
141  const auto i = ind[j];
142  U_vars[i]->dofMap().dof_indices(elem, dof_indices, U_vars[i]->number());
143  Moose::derivInsert(U[i].derivatives(), dof_indices[0], 1.0);
144  }
145  }
146  else
147  {
148  for (unsigned int i = 0; i < THMVACE1D::N_FLUX_INPUTS; i++)
149  U[i] = U_vars[i]->getElementalValueOld(elem);
150  }
151 
152  return U;
153 }
154 }
Moose::GenericType< Real, is_ad > GenericReal
std::vector< GenericReal< is_ad > > computeConservativeSolutionVector(const std::vector< GenericReal< is_ad >> &W, const GenericReal< is_ad > &A, const SinglePhaseFluidProperties &fp)
Computes the conservative solution vector from the primitive solution vector.
std::vector< GenericReal< is_ad > > computeFluxFromPrimitive(const std::vector< GenericReal< is_ad >> &W, const GenericReal< is_ad > &A, const SinglePhaseFluidProperties &fp)
Computes the numerical flux vector from the primitive solution vector.
static const std::string F
Definition: NS.h:165
std::vector< GenericReal< is_ad > > computePrimitiveSolutionVector(const std::vector< GenericReal< is_ad >> &U, const SinglePhaseFluidProperties &fp)
Computes the primitive solution vector from the conservative solution vector.
Common class for single phase fluid properties.
static const unsigned int N_FLUX_INPUTS
Number of numerical flux function inputs for 1D.
std::vector< GenericReal< is_ad > > getElementalSolutionVector(const Elem *elem, const std::vector< MooseVariable *> &U_vars, bool is_implicit)
Gets the elemental conservative solution vector.
static const unsigned int N_FLUX_OUTPUTS
Number of numerical flux function outputs for 1D.
static const unsigned int N_PRIM_VARS
static const std::string v
Definition: NS.h:84
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")