createFields.H
Go to the documentation of this file.
1 #include "createRDeltaT.H"
2 
3 #include "readGravitationalAcceleration.H"
4 
5 Info<< "Reading thermophysical properties\n" << endl;
6 autoPtr<psiReactionThermo> pThermo(psiReactionThermo::New(mesh));
7 psiReactionThermo& thermo = pThermo();
8 thermo.validate(args.executable(), "h", "e");
9 
10 SLGThermo slgThermo(mesh, thermo);
11 
12 basicSpecieMixture& composition = thermo.composition();
13 PtrList<volScalarField>& Y = composition.Y();
14 
15 const word inertSpecie(thermo.get<word>("inertSpecie"));
16 if (!composition.species().found(inertSpecie))
17 {
19  << "Inert specie " << inertSpecie << " not found in available species "
20  << composition.species()
21  << exit(FatalIOError);
22 }
23 
25 
26 multivariateSurfaceInterpolationScheme<scalar>::fieldTable fields;
27 
28 forAll(Y, i)
29 {
30  fields.add(Y[i]);
31 }
32 fields.add(thermo.he());
33 
35 (
36  IOobject
37  (
38  "rho",
39  runTime.timeName(),
40  mesh,
41  IOobject::NO_READ,
42  IOobject::AUTO_WRITE
43  ),
44  thermo.rho()
45 );
46 
47 // lagrangian effective density field - used externally (optional)
49 (
50  IOobject
51  (
52  "rhoEffLagrangian",
53  runTime.timeName(),
54  mesh,
55  IOobject::NO_READ,
56  IOobject::AUTO_WRITE
57  ),
58  mesh,
60 );
61 
62 // dynamic pressure field - used externally (optional)
64 (
65  IOobject
66  (
67  "pDyn",
68  runTime.timeName(),
69  mesh,
70  IOobject::NO_READ,
71  IOobject::AUTO_WRITE
72  ),
73  mesh,
75 );
76 
77 
78 Info<< "\nReading field U\n" << endl;
80 (
81  IOobject
82  (
83  "U",
84  runTime.timeName(),
85  mesh,
86  IOobject::MUST_READ,
87  IOobject::AUTO_WRITE
88  ),
89  mesh
90 );
91 
92 #include "compressibleCreatePhi.H"
93 
94 mesh.setFluxRequired(p.name());
95 
96 Info<< "Creating turbulence model\n" << endl;
97 autoPtr<compressible::turbulenceModel> turbulence
98 (
100  (
101  rho,
102  U,
103  phi,
104  thermo
105  )
106 );
107 
108 Info<< "Creating combustion model\n" << endl;
109 autoPtr<CombustionModel<psiReactionThermo>> combustion
110 (
112 );
113 
115 (
116  IOobject
117  (
118  "Qdot",
119  runTime.timeName(),
120  mesh,
121  IOobject::READ_IF_PRESENT,
122  IOobject::AUTO_WRITE
123  ),
124  mesh,
126 );
127 
128 #include "createDpdt.H"
129 
130 #include "createK.H"
131 
132 #include "createMRF.H"
133 #include "createClouds.H"
134 #include "createRadiationModel.H"
135 #include "createFvOptions.H"
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::dimPressure
const dimensionSet dimPressure
U
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedVector(dimVelocity, Zero))
createClouds.H
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
pDyn
volScalarField pDyn(IOobject("pDyn", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimPressure, Zero))
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::dimDensity
const dimensionSet dimDensity
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...
createRDeltaT.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
rhoEffLagrangian
volScalarField rhoEffLagrangian(IOobject("rhoEffLagrangian", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimDensity, Zero))
createK.H
p
volScalarField & p
Definition: createFields.H:23
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
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
args
Foam::argList args(argc, argv)
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))