createFluidFields.H
Go to the documentation of this file.
1// Initialise fluid field pointer lists
2PtrList<rhoThermo> thermoFluid(fluidRegions.size());
3PtrList<volScalarField> rhoFluid(fluidRegions.size());
4PtrList<volVectorField> UFluid(fluidRegions.size());
5PtrList<surfaceScalarField> phiFluid(fluidRegions.size());
6PtrList<uniformDimensionedScalarField> hRefFluid(fluidRegions.size());
7PtrList<volScalarField> ghFluid(fluidRegions.size());
8PtrList<surfaceScalarField> ghfFluid(fluidRegions.size());
9PtrList<compressible::turbulenceModel> turbulence(fluidRegions.size());
10PtrList<volScalarField> p_rghFluid(fluidRegions.size());
11PtrList<radiation::radiationModel> radiation(fluidRegions.size());
12
13List<scalar> initialMassFluid(fluidRegions.size());
14List<label> pRefCellFluid(fluidRegions.size(), -1);
15List<scalar> pRefValueFluid(fluidRegions.size(), Zero);
16List<bool> frozenFlowFluid(fluidRegions.size(), false);
17
18PtrList<dimensionedScalar> rhoMax(fluidRegions.size());
19PtrList<dimensionedScalar> rhoMin(fluidRegions.size());
20
21PtrList<IOMRFZoneList> MRFfluid(fluidRegions.size());
22PtrList<fv::options> fluidFvOptions(fluidRegions.size());
23
24PtrList<fvVectorMatrix> UEqFluid(fluidRegions.size());
25
26const uniformDimensionedVectorField& g = meshObjects::gravity::New(runTime);
27
28// Populate fluid field pointer lists
30{
31 Info<< "*** Reading fluid mesh thermophysical properties for region "
32 << fluidRegions[i].name() << nl << endl;
33
34 Info<< " Adding to thermoFluid\n" << endl;
35
36 thermoFluid.set
37 (
38 i,
39 rhoThermo::New(fluidRegions[i]).ptr()
40 );
41
42 Info<< " Adding to rhoFluid\n" << endl;
43 rhoFluid.set
44 (
45 i,
46 new volScalarField
47 (
48 IOobject
49 (
50 "rho",
51 runTime.timeName(),
52 fluidRegions[i],
53 IOobject::NO_READ,
54 IOobject::AUTO_WRITE
55 ),
56 thermoFluid[i].rho()
57 )
58 );
59
60 Info<< " Adding to UFluid\n" << endl;
61 UFluid.set
62 (
63 i,
64 new volVectorField
65 (
66 IOobject
67 (
68 "U",
69 runTime.timeName(),
70 fluidRegions[i],
71 IOobject::MUST_READ,
72 IOobject::AUTO_WRITE
73 ),
75 )
76 );
77
78 Info<< " Adding to phiFluid\n" << endl;
79 phiFluid.set
80 (
81 i,
82 new surfaceScalarField
83 (
84 IOobject
85 (
86 "phi",
87 runTime.timeName(),
88 fluidRegions[i],
89 IOobject::READ_IF_PRESENT,
90 IOobject::AUTO_WRITE
91 ),
92 linearInterpolate(rhoFluid[i]*UFluid[i])
93 & fluidRegions[i].Sf()
94 )
95 );
96
97 Info<< " Adding to hRefFluid\n" << endl;
98 hRefFluid.set
99 (
100 i,
101 new uniformDimensionedScalarField
102 (
103 IOobject
104 (
105 "hRef",
106 runTime.constant(),
107 fluidRegions[i],
108 IOobject::READ_IF_PRESENT,
109 IOobject::NO_WRITE
110 ),
111 dimensionedScalar("hRef", dimLength, Zero) // uses name
112 )
113 );
114
115 dimensionedScalar ghRef
116 (
117 mag(g.value()) > SMALL
118 ? g & (cmptMag(g.value())/mag(g.value()))*hRefFluid[i]
119 : dimensionedScalar("ghRef", g.dimensions()*dimLength, 0)
120 );
121
122 Info<< " Adding to ghFluid\n" << endl;
123 ghFluid.set
124 (
125 i,
126 new volScalarField
127 (
128 "gh",
129 (g & fluidRegions[i].C()) - ghRef
130 )
131 );
132
133 Info<< " Adding to ghfFluid\n" << endl;
134 ghfFluid.set
135 (
136 i,
137 new surfaceScalarField
138 (
139 "ghf",
140 (g & fluidRegions[i].Cf()) - ghRef
141 )
142 );
143
144 Info<< " Adding to turbulence\n" << endl;
145 turbulence.set
146 (
147 i,
148 compressible::turbulenceModel::New
149 (
150 rhoFluid[i],
151 UFluid[i],
152 phiFluid[i],
153 thermoFluid[i]
154 ).ptr()
155 );
156
157 p_rghFluid.set
158 (
159 i,
160 new volScalarField
161 (
162 IOobject
163 (
164 "p_rgh",
165 runTime.timeName(),
166 fluidRegions[i],
167 IOobject::MUST_READ,
168 IOobject::AUTO_WRITE
169 ),
170 fluidRegions[i]
171 )
172 );
173
174 // Force p_rgh to be consistent with p
175 p_rghFluid[i] = thermoFluid[i].p() - rhoFluid[i]*ghFluid[i];
176
177 fluidRegions[i].setFluxRequired(p_rghFluid[i].name());
178
179 radiation.set
180 (
181 i,
182 radiation::radiationModel::New(thermoFluid[i].T())
183 );
184
185 initialMassFluid[i] = fvc::domainIntegrate(rhoFluid[i]).value();
186
187 const dictionary& simpleDict =
188 fluidRegions[i].solutionDict().subDict("SIMPLE");
189
191 (
192 thermoFluid[i].p(),
193 p_rghFluid[i],
194 simpleDict,
195 pRefCellFluid[i],
197 );
198
199 simpleDict.readIfPresent("frozenFlow", frozenFlowFluid[i]);
200
201 rhoMax.set
202 (
203 i,
204 new dimensionedScalar("rhoMax", dimDensity, GREAT, simpleDict)
205 );
206
207 rhoMin.set
208 (
209 i,
210 new dimensionedScalar("rhoMin", dimDensity, Zero, simpleDict)
211 );
212
213 Info<< " Adding MRF\n" << endl;
214 MRFfluid.set
215 (
216 i,
217 new IOMRFZoneList(fluidRegions[i])
218 );
219
220 Info<< " Adding fvOptions\n" << endl;
222 (
223 i,
224 new fv::options(fluidRegions[i])
225 );
226
227 UEqFluid.set
228 (
229 i,
230 new fvVectorMatrix(UFluid[i], dimForce)
231 );
232
233 turbulence[i].validate();
234}
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
PtrList< volScalarField > ghFluid(fluidRegions.size())
const uniformDimensionedVectorField & g
PtrList< rhoThermo > thermoFluid(fluidRegions.size())
PtrList< IOMRFZoneList > MRFfluid(fluidRegions.size())
PtrList< fv::options > fluidFvOptions(fluidRegions.size())
PtrList< uniformDimensionedScalarField > hRefFluid(fluidRegions.size())
List< label > pRefCellFluid(fluidRegions.size(), -1)
PtrList< volVectorField > UFluid(fluidRegions.size())
PtrList< volScalarField > p_rghFluid(fluidRegions.size())
List< scalar > initialMassFluid(fluidRegions.size())
List< bool > frozenFlowFluid(fluidRegions.size(), false)
PtrList< volScalarField > rhoFluid(fluidRegions.size())
PtrList< surfaceScalarField > ghfFluid(fluidRegions.size())
PtrList< radiation::radiationModel > radiation(fluidRegions.size())
PtrList< fvVectorMatrix > UEqFluid(fluidRegions.size())
List< scalar > pRefValueFluid(fluidRegions.size(), Zero)
PtrList< surfaceScalarField > phiFluid(fluidRegions.size())
PtrList< fvMesh > fluidRegions(fluidNames.size())
volScalarField & p
const volScalarField & T
engineTime & runTime
const dimensionedScalar rhoMin
const dimensionedScalar rhoMax
compressible::turbulenceModel & turbulence
setRefCell(p, pimple.dict(), pRefCell, pRefValue)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333