createFields.H
Go to the documentation of this file.
1 #include "createRDeltaT.H"
2 
3 #include "readGravitationalAcceleration.H"
4 
5 Info<< "Creating combustion model\n" << endl;
6 
7 autoPtr<combustionModels::psiCombustionModel> combustion
8 (
10 );
11 
12 psiReactionThermo& thermo = combustion->thermo();
13 thermo.validate(args.executable(), "h", "e");
14 
15 SLGThermo slgThermo(mesh, thermo);
16 
17 basicSpecieMixture& composition = thermo.composition();
18 PtrList<volScalarField>& Y = composition.Y();
19 
20 const word inertSpecie(thermo.lookup("inertSpecie"));
21 if (!composition.species().found(inertSpecie))
22 {
24  << "Inert specie " << inertSpecie << " not found in available species "
25  << composition.species()
26  << exit(FatalIOError);
27 }
28 
29 volScalarField& p = thermo.p();
30 const volScalarField& T = thermo.T();
31 const volScalarField& psi = thermo.psi();
32 
33 multivariateSurfaceInterpolationScheme<scalar>::fieldTable fields;
34 
35 forAll(Y, i)
36 {
37  fields.add(Y[i]);
38 }
39 fields.add(thermo.he());
40 
42 (
43  IOobject
44  (
45  "rho",
46  runTime.timeName(),
47  mesh,
48  IOobject::NO_READ,
49  IOobject::AUTO_WRITE
50  ),
51  thermo.rho()
52 );
53 
54 // lagrangian effective density field - used externally (optional)
56 (
57  IOobject
58  (
59  "rhoEffLagrangian",
60  runTime.timeName(),
61  mesh,
62  IOobject::NO_READ,
63  IOobject::AUTO_WRITE
64  ),
65  mesh,
66  dimensionedScalar("zero", dimDensity, 0.0)
67 );
68 
69 // dynamic pressure field - used externally (optional)
71 (
72  IOobject
73  (
74  "pDyn",
75  runTime.timeName(),
76  mesh,
77  IOobject::NO_READ,
78  IOobject::AUTO_WRITE
79  ),
80  mesh,
81  dimensionedScalar("zero", dimPressure, 0.0)
82 );
83 
84 
85 Info<< "\nReading field U\n" << endl;
87 (
88  IOobject
89  (
90  "U",
91  runTime.timeName(),
92  mesh,
93  IOobject::MUST_READ,
94  IOobject::AUTO_WRITE
95  ),
96  mesh
97 );
98 
99 #include "compressibleCreatePhi.H"
100 
101 mesh.setFluxRequired(p.name());
102 
103 Info<< "Creating turbulence model\n" << endl;
104 autoPtr<compressible::turbulenceModel> turbulence
105 (
107  (
108  rho,
109  U,
110  phi,
111  thermo
112  )
113 );
114 
115 // Set the turbulence into the combustion model
116 combustion->setTurbulence(turbulence());
117 
118 Info<< "Creating field dpdt\n" << endl;
120 (
121  IOobject
122  (
123  "dpdt",
124  runTime.timeName(),
125  mesh
126  ),
127  mesh,
128  dimensionedScalar("dpdt", p.dimensions()/dimTime, 0)
129 );
130 
131 Info<< "Creating field kinetic energy K\n" << endl;
132 volScalarField K("K", 0.5*magSqr(U));
133 
135 (
136  IOobject
137  (
138  "Qdot",
139  runTime.timeName(),
140  mesh,
141  IOobject::READ_IF_PRESENT,
142  IOobject::AUTO_WRITE
143  ),
144  mesh,
146 );
147 
148 #include "createMRF.H"
149 #include "createClouds.H"
150 #include "createRadiationModel.H"
surfaceScalarField & phi
forAll(Y, i)
Definition: createFields.H:106
const word & executable() const
Name of executable without the path.
Definition: argListI.H:30
Info<< "Creating combustion model\"<< endl;autoPtr< combustionModels::psiCombustionModel > combustion(combustionModels::psiCombustionModel::New(mesh))
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
SLGThermo slgThermo(mesh, thermo)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedVector("zero", dimVelocity, Zero))
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
CGAL::Exact_predicates_exact_constructions_kernel K
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:104
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
psiReactionThermo & thermo
Definition: createFields.H:31
rho
Definition: createFields.H:81
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dynamicFvMesh & mesh
volScalarField & dpdt
const dimensionSet dimPressure
dimensioned< scalar > magSqr(const dimensioned< Type > &)
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
Info<< "Reading field U\"<< 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\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
const dimensionSet dimEnergy
const dimensionSet dimDensity
volScalarField rhoEffLagrangian(IOobject("rhoEffLagrangian", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar("zero", dimDensity, 0.0))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
volScalarField pDyn(IOobject("pDyn", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar("zero", dimPressure, 0.0))
Foam::argList args(argc, argv)
const word inertSpecie(thermo.lookup("inertSpecie"))
Creates and initialises the face-flux field phi.
volScalarField Qdot(IOobject("Qdot", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), mesh, dimensionedScalar("Qdot", dimEnergy/dimVolume/dimTime, 0.0))
IOerror FatalIOError