Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / applications / solvers / lagrangian / reactingParcelFilmFoam / createFields.H
blobe2db33f1f7f64960cd07f9898b039970a63c9b76
1     Info<< "Reading thermophysical properties\n" << endl;
3     autoPtr<psiChemistryModel> pChemistry
4     (
5         psiChemistryModel::New(mesh)
6     );
7     psiChemistryModel& chemistry = pChemistry();
9     hsCombustionThermo& thermo = chemistry.thermo();
11     SLGThermo slgThermo(mesh, thermo);
13     basicMultiComponentMixture& composition = thermo.composition();
14     PtrList<volScalarField>& Y = composition.Y();
16     const word inertSpecie(thermo.lookup("inertSpecie"));
18     Info<< "Creating field rho\n" << endl;
19     volScalarField rho
20     (
21         IOobject
22         (
23             "rho",
24             runTime.timeName(),
25             mesh,
26             IOobject::NO_READ,
27             IOobject::AUTO_WRITE
28         ),
29         thermo.rho()
30     );
32     volScalarField& p = thermo.p();
33     volScalarField& hs = thermo.hs();
34     const volScalarField& T = thermo.T();
35     const volScalarField& psi = thermo.psi();
37     Info<< "\nReading field U\n" << endl;
38     volVectorField U
39     (
40         IOobject
41         (
42             "U",
43             runTime.timeName(),
44             mesh,
45             IOobject::MUST_READ,
46             IOobject::AUTO_WRITE
47         ),
48         mesh
49     );
51     #include "compressibleCreatePhi.H"
53     Info<< "Creating field kappa\n" << endl;
54     DimensionedField<scalar, volMesh> kappa
55     (
56         IOobject
57         (
58             "kappa",
59             runTime.timeName(),
60             mesh,
61             IOobject::NO_READ,
62             IOobject::AUTO_WRITE
63         ),
64         mesh,
65         dimensionedScalar("zero", dimless, 0.0)
66     );
68     Info<< "Creating turbulence model\n" << endl;
69     autoPtr<compressible::turbulenceModel> turbulence
70     (
71         compressible::turbulenceModel::New
72         (
73             rho,
74             U,
75             phi,
76             thermo
77         )
78     );
80     Info<< "Creating field DpDt\n" << endl;
81     volScalarField DpDt
82     (
83         "DpDt",
84         fvc::DDt(surfaceScalarField("phiU", phi/fvc::interpolate(rho)), p)
85     );
88     Info<< "Calculating field g.h\n" << endl;
89     volScalarField gh("gh", g & mesh.C());
91     surfaceScalarField ghf("ghf", g & mesh.Cf());
93     volScalarField p_rgh
94     (
95         IOobject
96         (
97             "p_rgh",
98             runTime.timeName(),
99             mesh,
100             IOobject::MUST_READ,
101             IOobject::AUTO_WRITE
102         ),
103         mesh
104     );
106     // Force p_rgh to be consistent with p
107     p_rgh = p - rho*gh;
109     multivariateSurfaceInterpolationScheme<scalar>::fieldTable fields;
111     forAll(Y, i)
112     {
113         fields.add(Y[i]);
114     }
115     fields.add(hs);
117     IOdictionary additionalControlsDict
118     (
119         IOobject
120         (
121             "additionalControls",
122             runTime.constant(),
123             mesh,
124             IOobject::MUST_READ_IF_MODIFIED,
125             IOobject::NO_WRITE
126         )
127     );
129     Switch solvePrimaryRegion
130     (
131         additionalControlsDict.lookup("solvePrimaryRegion")
132     );
134     DimensionedField<scalar, volMesh> chemistrySh
135     (
136         IOobject
137         (
138             "chemistry::Sh",
139             runTime.timeName(),
140             mesh,
141             IOobject::NO_READ,
142             IOobject::NO_WRITE
143         ),
144         mesh,
145         dimensionedScalar("chemistrySh", dimEnergy/dimTime/dimVolume, 0.0)
146     );