createFields.H
Go to the documentation of this file.
1#include "createRDeltaT.H"
2
3Info<< "Reading thermophysical properties\n" << endl;
4autoPtr<psiReactionThermo> pThermo(psiReactionThermo::New(mesh));
5psiReactionThermo& thermo = pThermo();
6thermo.validate(args.executable(), "h", "e");
7
8basicSpecieMixture& composition = thermo.composition();
9PtrList<volScalarField>& Y = composition.Y();
10
11const word inertSpecie(thermo.get<word>("inertSpecie"));
12if (!composition.species().found(inertSpecie))
13{
15 << "Inert specie " << inertSpecie << " not found in available species "
16 << composition.species() << exit(FatalIOError);
17}
18
19volScalarField rho
20(
21 IOobject
22 (
23 "rho",
24 runTime.timeName(),
25 mesh
26 ),
27 thermo.rho()
28);
29
30Info<< "Reading field U\n" << endl;
31volVectorField U
32(
33 IOobject
34 (
35 "U",
36 runTime.timeName(),
37 mesh,
38 IOobject::MUST_READ,
39 IOobject::AUTO_WRITE
40 ),
41 mesh
42);
43
44volScalarField& p = thermo.p();
45
47
49
50mesh.setFluxRequired(p.name());
51
52Info << "Creating turbulence model.\n" << nl;
53autoPtr<compressible::turbulenceModel> turbulence
54(
55 compressible::turbulenceModel::New
56 (
57 rho,
58 U,
59 phi,
60 thermo
61 )
62);
63
64Info<< "Creating reaction model\n" << endl;
65autoPtr<CombustionModel<psiReactionThermo>> reaction
66(
67 CombustionModel<psiReactionThermo>::New(thermo, turbulence())
68);
69
70multivariateSurfaceInterpolationScheme<scalar>::fieldTable fields;
71
73{
74 fields.add(Y[i]);
75}
76fields.add(thermo.he());
77
78volScalarField Qdot
79(
80 IOobject
81 (
82 "Qdot",
83 runTime.timeName(),
84 mesh,
85 IOobject::READ_IF_PRESENT,
86 IOobject::AUTO_WRITE
87 ),
88 mesh,
89 dimensionedScalar(dimEnergy/dimVolume/dimTime, Zero)
90);
91
92#include "createDpdt.H"
93
94#include "createK.H"
95
96#include "createMRF.H"
97#include "createFvOptions.H"
surfaceScalarField & phi
pimpleControl & pimple
const word & executable() const noexcept
Name of executable without the path.
Definition: argListI.H:51
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
U
Definition: pEqn.H:72
volScalarField & p
Creates and initialises the face-flux field phi.
dynamicFvMesh & mesh
engineTime & runTime
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:467
const pressureControl & pressureControl
compressible::turbulenceModel & turbulence
CombustionModel< rhoReactionThermo > & reaction
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
Foam::argList args(argc, argv)
scalar Qdot
Definition: solveChemistry.H:2
Info<< "Reading thermophysical properties\n"<< endl;autoPtr< psiReactionThermo > pThermo(psiReactionThermo::New(mesh))
const word inertSpecie(thermo.get< word >("inertSpecie"))
basicSpecieMixture & composition
Definition: createFields.H:8
PtrList< volScalarField > & Y
Definition: createFields.H:9
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333