Transferred copyright to the OpenFOAM Foundation
[OpenFOAM-2.0.x.git] / applications / solvers / combustion / fireFoam / createFields.H
blobf1e152dc320b47f428e326350f66aae48c947338
1     Info<< "Reading thermophysical properties\n" << endl;
3     autoPtr<hsCombustionThermo> pThermo
4     (
5         hsCombustionThermo::New(mesh)
6     );
7     hsCombustionThermo& thermo = pThermo();
9     SLGThermo slgThermo(mesh, thermo);
11     basicMultiComponentMixture& composition = thermo.composition();
12     PtrList<volScalarField>& Y = composition.Y();
14     const word inertSpecie(thermo.lookup("inertSpecie"));
16     Info<< "Creating field rho\n" << endl;
17     volScalarField rho
18     (
19         IOobject
20         (
21             "rho",
22             runTime.timeName(),
23             mesh,
24             IOobject::NO_READ,
25             IOobject::AUTO_WRITE
26         ),
27         thermo.rho()
28     );
30     volScalarField& p = thermo.p();
31     volScalarField& hs = thermo.hs();
32     const volScalarField& T = thermo.T();
33     const volScalarField& psi = thermo.psi();
35     Info<< "\nReading field U\n" << endl;
36     volVectorField U
37     (
38         IOobject
39         (
40             "U",
41             runTime.timeName(),
42             mesh,
43             IOobject::MUST_READ,
44             IOobject::AUTO_WRITE
45         ),
46         mesh
47     );
49     #include "compressibleCreatePhi.H"
51     Info<< "Creating turbulence model\n" << endl;
52     autoPtr<compressible::turbulenceModel> turbulence
53     (
54         compressible::turbulenceModel::New
55         (
56             rho,
57             U,
58             phi,
59             thermo
60         )
61     );
63     IOdictionary combustionProperties
64     (
65         IOobject
66         (
67             "combustionProperties",
68             runTime.constant(),
69             mesh,
70             IOobject::MUST_READ_IF_MODIFIED,
71             IOobject::NO_WRITE
72         )
73     );
75     Info<< "Creating combustion model\n" << endl;
76     autoPtr<combustionModel> combustion
77     (
78         combustionModel::combustionModel::New
79         (
80             combustionProperties,
81             thermo,
82             turbulence(),
83             phi,
84             rho
85         )
86     );
88     volScalarField dQ
89     (
90         IOobject
91         (
92             "dQ",
93             runTime.timeName(),
94             mesh,
95             IOobject::NO_READ,
96             IOobject::AUTO_WRITE
97         ),
98         mesh,
99         dimensionedScalar("dQ", dimMass/pow3(dimTime)/dimLength, 0.0)
100     );
102     Info<< "Creating field DpDt\n" << endl;
103     volScalarField DpDt
104     (
105         "DpDt",
106         fvc::DDt(surfaceScalarField("phiU", phi/fvc::interpolate(rho)), p)
107     );
110     Info<< "Calculating field g.h\n" << endl;
111     volScalarField gh("gh", g & mesh.C());
113     surfaceScalarField ghf("ghf", g & mesh.Cf());
115     volScalarField p_rgh
116     (
117         IOobject
118         (
119             "p_rgh",
120             runTime.timeName(),
121             mesh,
122             IOobject::MUST_READ,
123             IOobject::AUTO_WRITE
124         ),
125         mesh
126     );
128     // Force p_rgh to be consistent with p
129     p_rgh = p - rho*gh;
131     multivariateSurfaceInterpolationScheme<scalar>::fieldTable fields;
133     forAll(Y, i)
134     {
135         fields.add(Y[i]);
136     }
137     fields.add(hs);
139     IOdictionary additionalControlsDict
140     (
141         IOobject
142         (
143             "additionalControls",
144             runTime.constant(),
145             mesh,
146             IOobject::MUST_READ_IF_MODIFIED,
147             IOobject::NO_WRITE
148         )
149     );
151     Switch solvePrimaryRegion
152     (
153         additionalControlsDict.lookup("solvePrimaryRegion")
154     );