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