EEqn.H
Go to the documentation of this file.
1 {
2  volScalarField& he = thermo.he();
3 
4  const tmp<volScalarField>& tCp = thermo.Cp();
5  const tmp<volScalarField>& tCv = thermo.Cv();
6 
7  const volScalarField& Cp = tCp();
8  const volScalarField& Cv = tCv();
9  const scalar gamma = max(Cp/Cv).value();
10 
11  if (mag(gamma - min(Cp/Cv).value()) > VSMALL)
12  {
13  notImplemented("gamma not constant in space");
14  }
15 
16  const dictionary& thermoDict = thermo.subDict("mixture");
17 
18  const dictionary& eosDict = thermoDict.subDict("equationOfState");
19 
20  bool local = eosDict.getOrDefault("local", false);
21 
22  // Evolve T as:
23  //
24  // T_1 = T_0 \frac{p}{p_0}^{\frac{\gamma - 1}{\gamma}}
25 
26  if (!local)
27  {
28  const scalar T0 = eosDict.get<scalar>("T0");
29  const scalar p0 = eosDict.get<scalar>("p0");
30 
31  he = thermo.he(p, pow(p/p0, (gamma - scalar(1))/gamma)*T0);
32  }
33  else
34  {
35  const volScalarField& T0 = T.oldTime();
36  const volScalarField& p0 = p.oldTime();
37 
38  he = thermo.he(p, pow(p/p0, (gamma - scalar(1))/gamma)*T0);
39  }
40 
41  thermo.correct();
42 
43  psi = 1.0/((Cp - Cv)*T);
44 
45  rho = thermo.rho();
46  rho.relax();
47 
48  rho.writeMinMax(Info);
49 }
p
volScalarField & p
Definition: createFieldRefs.H:8
Cv
const volScalarField & Cv
Definition: EEqn.H:8
tCp
const tmp< volScalarField > & tCp
Definition: EEqn.H:4
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
psi
psi
Definition: EEqn.H:43
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
eosDict
const dictionary & eosDict
Definition: EEqn.H:18
notImplemented
#define notImplemented(functionName)
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:507
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
he
he
Definition: EEqn.H:38
T
const volScalarField & T
Definition: createFieldRefs.H:2
thermoDict
const dictionary & thermoDict
Definition: EEqn.H:16
rho
rho
Definition: EEqn.H:45
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
gamma
const scalar gamma
Definition: EEqn.H:9
Cp
const volScalarField & Cp
Definition: EEqn.H:7
p0
const volScalarField & p0
Definition: EEqn.H:36
T0
scalar T0
Definition: createFields.H:22
tCv
const tmp< volScalarField > & tCv
Definition: EEqn.H:5