ENH: directMappedFixedValue: use different tag
[OpenFOAM-2.0.x.git] / applications / solvers / combustion / XiFoam / createFields.H
blob7a7e796b62c1df021b6d19b2439babc15c985bdb
1     Info<< "Reading thermophysical properties\n" << endl;
3     autoPtr<hhuCombustionThermo> pThermo
4     (
5         hhuCombustionThermo::New(mesh)
6     );
7     hhuCombustionThermo& thermo = pThermo();
8     basicMultiComponentMixture& composition = thermo.composition();
10     volScalarField rho
11     (
12         IOobject
13         (
14             "rho",
15             runTime.timeName(),
16             mesh,
17             IOobject::NO_READ,
18             IOobject::AUTO_WRITE
19         ),
20         thermo.rho()
21     );
23     volScalarField& p = thermo.p();
24     const volScalarField& psi = thermo.psi();
25     volScalarField& h = thermo.h();
26     volScalarField& hu = thermo.hu();
28     volScalarField& b = composition.Y("b");
29     Info<< "min(b) = " << min(b).value() << endl;
31     const volScalarField& T = thermo.T();
34     Info<< "\nReading field U\n" << endl;
35     volVectorField U
36     (
37         IOobject
38         (
39             "U",
40             runTime.timeName(),
41             mesh,
42             IOobject::MUST_READ,
43             IOobject::AUTO_WRITE
44         ),
45         mesh
46     );
48 #   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     Info<< "Creating field DpDt\n" << endl;
64     volScalarField DpDt
65     (
66         fvc::DDt(surfaceScalarField("phiU", phi/fvc::interpolate(rho)), p)
67     );
70     Info<< "Creating field Xi\n" << endl;
71     volScalarField Xi
72     (
73         IOobject
74         (
75             "Xi",
76             runTime.timeName(),
77             mesh,
78             IOobject::MUST_READ,
79             IOobject::AUTO_WRITE
80         ),
81         mesh
82     );
85     Info<< "Creating the unstrained laminar flame speed\n" << endl;
86     autoPtr<laminarFlameSpeed> unstrainedLaminarFlameSpeed
87     (
88         laminarFlameSpeed::New(thermo)
89     );
92     Info<< "Reading strained laminar flame speed field Su\n" << endl;
93     volScalarField Su
94     (
95         IOobject
96         (
97             "Su",
98             runTime.timeName(),
99             mesh,
100             IOobject::MUST_READ,
101             IOobject::AUTO_WRITE
102         ),
103         mesh
104     );
106     dimensionedScalar SuMin = 0.01*Su.average();
107     dimensionedScalar SuMax = 4*Su.average();
109     Info<< "Calculating turbulent flame speed field St\n" << endl;
110     volScalarField St
111     (
112         IOobject
113         (
114             "St",
115             runTime.timeName(),
116             mesh,
117             IOobject::NO_READ,
118             IOobject::AUTO_WRITE
119         ),
120         Xi*Su
121     );
124     multivariateSurfaceInterpolationScheme<scalar>::fieldTable fields;
126     if (composition.contains("ft"))
127     {
128         fields.add(composition.Y("ft"));
129     }
131     fields.add(b);
132     fields.add(h);
133     fields.add(hu);