initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / solvers / heatTransfer / chtMultiRegionFoam / fluid / createFluidFields.H
blob1826a77217e86ac89941941c2d8e452d783e94fb
1     // Initialise fluid field pointer lists
2     PtrList<basicPsiThermo> 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> DpDtFluid(fluidRegions.size());
11     List<scalar> initialMassFluid(fluidRegions.size());
13     // Populate fluid field pointer lists
14     forAll(fluidRegions, i)
15     {
16         Info<< "*** Reading fluid mesh thermophysical properties for region "
17             << fluidRegions[i].name() << nl << endl;
19         Info<< "    Adding to thermoFluid\n" << endl;
20         thermoFluid.set
21         (
22             i,
23             basicPsiThermo::New(fluidRegions[i]).ptr()
24         );
26         Info<< "    Adding to rhoFluid\n" << endl;
27         rhoFluid.set
28         (
29             i,
30             new volScalarField
31             (
32                 IOobject
33                 (
34                     "rho",
35                     runTime.timeName(),
36                     fluidRegions[i],
37                     IOobject::NO_READ,
38                     IOobject::AUTO_WRITE
39                 ),
40                 thermoFluid[i].rho()
41             )
42         );
44         Info<< "    Adding to KFluid\n" << endl;
45         KFluid.set
46         (
47             i,
48             new volScalarField
49             (
50                 IOobject
51                 (
52                     "K",
53                     runTime.timeName(),
54                     fluidRegions[i],
55                     IOobject::NO_READ,
56                     IOobject::NO_WRITE
57                 ),
58                 thermoFluid[i].Cp()*thermoFluid[i].alpha()
59             )
60         );
62         Info<< "    Adding to UFluid\n" << endl;
63         UFluid.set
64         (
65             i,
66             new volVectorField
67             (
68                 IOobject
69                 (
70                     "U",
71                     runTime.timeName(),
72                     fluidRegions[i],
73                     IOobject::MUST_READ,
74                     IOobject::AUTO_WRITE
75                 ),
76                 fluidRegions[i]
77             )
78         );
80         Info<< "    Adding to phiFluid\n" << endl;
81         phiFluid.set
82         (
83             i,
84             new surfaceScalarField
85             (
86                 IOobject
87                 (
88                     "phi",
89                     runTime.timeName(),
90                     fluidRegions[i],
91                     IOobject::READ_IF_PRESENT,
92                     IOobject::AUTO_WRITE
93                 ),
94                 linearInterpolate(rhoFluid[i]*UFluid[i])
95                     & fluidRegions[i].Sf()
96             )
97         );
99         Info<< "    Adding to gFluid\n" << endl;
100         gFluid.set
101         (
102             i,
103             new uniformDimensionedVectorField
104             (
105                 IOobject
106                 (
107                     "g",
108                     runTime.constant(),
109                     fluidRegions[i],
110                     IOobject::MUST_READ,
111                     IOobject::NO_WRITE
112                 )
113             )
114         );
116         Info<< "    Adding to turbulence\n" << endl;
117         turbulence.set
118         (
119             i,
120             autoPtr<compressible::turbulenceModel>
121             (
122                 compressible::turbulenceModel::New
123                 (
124                     rhoFluid[i],
125                     UFluid[i],
126                     phiFluid[i],
127                     thermoFluid[i]
128                 )
129             ).ptr()
130         );
132         Info<< "    Adding to DpDtFluid\n" << endl;
133         DpDtFluid.set
134         (
135             i,
136             new volScalarField
137             (
138                 "DpDt",
139                 fvc::DDt
140                 (
141                     surfaceScalarField
142                     (
143                         "phiU",
144                         phiFluid[i]/fvc::interpolate(rhoFluid[i])
145                     ),
146                     thermoFluid[i].p()
147                 )
148             )
149         );
151         initialMassFluid[i] = fvc::domainIntegrate(rhoFluid[i]).value();
152     }