ENH: Time: access to libs
[OpenFOAM-2.0.x.git] / applications / solvers / heatTransfer / chtMultiRegionFoam / chtMultiRegionSimpleFoam / fluid / createFluidFields.H
blob3d8531e1bb6bb98ac532724d71f3a62337950744
1     // Initialise fluid field pointer lists
2     PtrList<basicRhoThermo> thermoFluid(fluidRegions.size());
3     PtrList<volScalarField> rhoFluid(fluidRegions.size());
4     PtrList<volScalarField> KFluid(fluidRegions.size());
5     PtrList<volVectorField> UFluid(fluidRegions.size());
6     PtrList<surfaceScalarField> phiFluid(fluidRegions.size());
7     PtrList<uniformDimensionedVectorField> gFluid(fluidRegions.size());
8     PtrList<compressible::turbulenceModel> turbulence(fluidRegions.size());
9     PtrList<volScalarField> p_rghFluid(fluidRegions.size());
10     PtrList<volScalarField> ghFluid(fluidRegions.size());
11     PtrList<surfaceScalarField> ghfFluid(fluidRegions.size());
12     PtrList<radiation::radiationModel> radiation(fluidRegions.size());
14     List<scalar> initialMassFluid(fluidRegions.size());
15     List<label> pRefCellFluid(fluidRegions.size(),0);
16     List<scalar> pRefValueFluid(fluidRegions.size(),0.0);
18     PtrList<dimensionedScalar> rhoMax(fluidRegions.size());
19     PtrList<dimensionedScalar> rhoMin(fluidRegions.size());
21     // Populate fluid field pointer lists
22     forAll(fluidRegions, i)
23     {
24         Info<< "*** Reading fluid mesh thermophysical properties for region "
25             << fluidRegions[i].name() << nl << endl;
27         Info<< "    Adding to thermoFluid\n" << endl;
29         thermoFluid.set
30         (
31             i,
32             basicRhoThermo::New(fluidRegions[i]).ptr()
33         );
35         Info<< "    Adding to rhoFluid\n" << endl;
36         rhoFluid.set
37         (
38             i,
39             new volScalarField
40             (
41                 IOobject
42                 (
43                     "rho",
44                     runTime.timeName(),
45                     fluidRegions[i],
46                     IOobject::NO_READ,
47                     IOobject::AUTO_WRITE
48                 ),
49                 thermoFluid[i].rho()
50             )
51         );
53         Info<< "    Adding to KFluid\n" << endl;
54         KFluid.set
55         (
56             i,
57             new volScalarField
58             (
59                 IOobject
60                 (
61                     "K",
62                     runTime.timeName(),
63                     fluidRegions[i],
64                     IOobject::NO_READ,
65                     IOobject::NO_WRITE
66                 ),
67                 thermoFluid[i].Cp()*thermoFluid[i].alpha()
68             )
69         );
71         Info<< "    Adding to UFluid\n" << endl;
72         UFluid.set
73         (
74             i,
75             new volVectorField
76             (
77                 IOobject
78                 (
79                     "U",
80                     runTime.timeName(),
81                     fluidRegions[i],
82                     IOobject::MUST_READ,
83                     IOobject::AUTO_WRITE
84                 ),
85                 fluidRegions[i]
86             )
87         );
89         Info<< "    Adding to phiFluid\n" << endl;
90         phiFluid.set
91         (
92             i,
93             new surfaceScalarField
94             (
95                 IOobject
96                 (
97                     "phi",
98                     runTime.timeName(),
99                     fluidRegions[i],
100                     IOobject::READ_IF_PRESENT,
101                     IOobject::AUTO_WRITE
102                 ),
103                 linearInterpolate(rhoFluid[i]*UFluid[i])
104                     & fluidRegions[i].Sf()
105             )
106         );
108         Info<< "    Adding to gFluid\n" << endl;
109         gFluid.set
110         (
111             i,
112             new uniformDimensionedVectorField
113             (
114                 IOobject
115                 (
116                     "g",
117                     runTime.constant(),
118                     fluidRegions[i],
119                     IOobject::MUST_READ,
120                     IOobject::NO_WRITE
121                 )
122             )
123         );
125         Info<< "    Adding to turbulence\n" << endl;
126         turbulence.set
127         (
128             i,
129             compressible::turbulenceModel::New
130             (
131                 rhoFluid[i],
132                 UFluid[i],
133                 phiFluid[i],
134                 thermoFluid[i]
135             ).ptr()
136         );
138         Info<< "    Adding to ghFluid\n" << endl;
139         ghFluid.set
140         (
141             i,
142             new volScalarField("gh", gFluid[i] & fluidRegions[i].C())
143         );
145         Info<< "    Adding to ghfFluid\n" << endl;
146         ghfFluid.set
147         (
148             i,
149             new surfaceScalarField("ghf", gFluid[i] & fluidRegions[i].Cf())
150         );
152         p_rghFluid.set
153         (
154             i,
155             new volScalarField
156             (
157                 IOobject
158                 (
159                     "p_rgh",
160                     runTime.timeName(),
161                     fluidRegions[i],
162                     IOobject::MUST_READ,
163                     IOobject::AUTO_WRITE
164                 ),
165                 fluidRegions[i]
166             )
167         );
169         // Force p_rgh to be consistent with p
170         p_rghFluid[i] = thermoFluid[i].p() - rhoFluid[i]*ghFluid[i];
172         radiation.set
173         (
174             i,
175             radiation::radiationModel::New(thermoFluid[i].T())
176         );
178         initialMassFluid[i] = fvc::domainIntegrate(rhoFluid[i]).value();
180         setRefCell
181         (
182             thermoFluid[i].p(),
183             p_rghFluid[i],
184             fluidRegions[i].solutionDict().subDict("SIMPLE"),
185             pRefCellFluid[i],
186             pRefValueFluid[i]
187         );
189         rhoMax.set
190         (
191             i,
192             new dimensionedScalar
193             (
194                 fluidRegions[i].solutionDict().subDict("SIMPLE").lookup
195                 (
196                     "rhoMax"
197                 )
198             )
199         );
201         rhoMin.set
202         (
203             i,
204             new dimensionedScalar
205             (
206                 fluidRegions[i].solutionDict().subDict("SIMPLE").lookup
207                 (
208                     "rhoMin"
209                 )
210             )
211         );
212     }