createFields.H
Go to the documentation of this file.
1 Info<< "Reading thermophysical properties\n" << endl;
2 autoPtr<psiReactionThermo> pThermo(psiReactionThermo::New(mesh));
3 psiReactionThermo& thermo = pThermo();
4 thermo.validate(args.executable(), "h", "e");
5 
6 SLGThermo slgThermo(mesh, thermo);
7 
8 basicSpecieMixture& composition = thermo.composition();
9 PtrList<volScalarField>& Y = composition.Y();
10 
11 const word inertSpecie(thermo.get<word>("inertSpecie"));
12 if (!composition.species().found(inertSpecie))
13 {
15  << "Inert specie " << inertSpecie << " not found in available species "
16  << composition.species()
17  << exit(FatalIOError);
18 }
19 
20 Info<< "Creating field rho\n" << endl;
22 (
23  IOobject
24  (
25  "rho",
26  runTime.timeName(),
27  mesh,
28  IOobject::NO_READ,
29  IOobject::AUTO_WRITE
30  ),
31  thermo.rho()
32 );
33 
34 volScalarField& p = thermo.p();
35 
36 Info<< "\nReading field U\n" << endl;
38 (
39  IOobject
40  (
41  "U",
42  runTime.timeName(),
43  mesh,
44  IOobject::MUST_READ,
45  IOobject::AUTO_WRITE
46  ),
47  mesh
48 );
49 
50 #include "compressibleCreatePhi.H"
51 
52 #include "createMRF.H"
53 
54 
55 Info<< "Creating turbulence model\n" << endl;
56 autoPtr<compressible::turbulenceModel> turbulence
57 (
59  (
60  rho,
61  U,
62  phi,
63  thermo
64  )
65 );
66 
67 Info<< "Creating combustion model\n" << endl;
68 autoPtr<CombustionModel<psiReactionThermo>> combustion
69 (
71 );
72 
73 
74 #include "readGravitationalAcceleration.H"
75 #include "readhRef.H"
76 #include "gh.H"
77 #include "readpRef.H"
78 
80 (
81  IOobject
82  (
83  "p_rgh",
84  runTime.timeName(),
85  mesh,
86  IOobject::MUST_READ,
87  IOobject::AUTO_WRITE
88  ),
89  mesh
90 );
91 
92 mesh.setFluxRequired(p_rgh.name());
93 
94 #include "phrghEqn.H"
95 
96 
97 multivariateSurfaceInterpolationScheme<scalar>::fieldTable fields;
98 
99 forAll(Y, i)
100 {
101  fields.add(Y[i]);
102 }
103 fields.add(thermo.he());
104 
106 (
107  IOobject
108  (
109  "Qdot",
110  runTime.timeName(),
111  mesh,
112  IOobject::READ_IF_PRESENT,
113  IOobject::AUTO_WRITE
114  ),
115  mesh,
117 );
118 
119 #include "createDpdt.H"
120 
121 #include "createK.H"
122 
123 #include "createClouds.H"
124 #include "createSurfaceFilmModel.H"
125 #include "createPyrolysisModel.H"
126 #include "createRadiationModel.H"
127 #include "createFvOptions.H"
runTime
engineTime & runTime
Definition: createEngineTime.H:13
gh.H
U
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedVector(dimVelocity, Zero))
p_rgh
p_rgh
Definition: createFields.H:95
forAll
forAll(U.boundaryField(), patchi)
Definition: createFields.H:52
turbulence
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimEnergy
const dimensionSet dimEnergy
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
createPyrolysisModel.H
rho
rho
Definition: createFields.H:81
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
createRadiationModel.H
createFvOptions.H
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Qdot
volScalarField Qdot(IOobject("Qdot", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimEnergy/dimVolume/dimTime, Zero))
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
FatalIOErrorIn
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:467
Y
PtrList< volScalarField > & Y
Definition: createFields.H:9
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
createK.H
p
volScalarField & p
Definition: createFields.H:23
readpRef.H
inertSpecie
const word inertSpecie(thermo.get< word >("inertSpecie"))
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
compressibleCreatePhi.H
Creates and initialises the face-flux field phi.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
createDpdt.H
slgThermo
SLGThermo slgThermo(mesh, thermo)
Foam::argList::executable
const word & executable() const noexcept
Name of executable without the path.
Definition: argListI.H:51
createSurfaceFilmModel.H
readhRef.H
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
createClouds.H
args
Foam::argList args(argc, argv)
phrghEqn.H
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
combustion
Info<< "Creating combustion model\n"<< endl;autoPtr< CombustionModel< psiReactionThermo > > combustion(CombustionModel< psiReactionThermo >::New(thermo, turbulence()))
composition
basicSpecieMixture & composition
Definition: createFields.H:8
pThermo
Info<< "Reading thermophysical properties\n"<< endl;autoPtr< psiReactionThermo > pThermo(psiReactionThermo::New(mesh))