createFluidFields.H
Go to the documentation of this file.
1 // Initialise fluid field pointer lists
2 PtrList<rhoThermo> thermoFluid(fluidRegions.size());
3 PtrList<volScalarField> rhoFluid(fluidRegions.size());
4 PtrList<volVectorField> UFluid(fluidRegions.size());
5 PtrList<surfaceScalarField> phiFluid(fluidRegions.size());
6 PtrList<uniformDimensionedVectorField> gFluid(fluidRegions.size());
7 PtrList<uniformDimensionedScalarField> hRefFluid(fluidRegions.size());
8 PtrList<volScalarField> ghFluid(fluidRegions.size());
9 PtrList<surfaceScalarField> ghfFluid(fluidRegions.size());
10 PtrList<compressible::turbulenceModel> turbulence(fluidRegions.size());
11 PtrList<volScalarField> p_rghFluid(fluidRegions.size());
12 PtrList<radiation::radiationModel> radiation(fluidRegions.size());
13 
14 List<scalar> initialMassFluid(fluidRegions.size());
15 List<label> pRefCellFluid(fluidRegions.size(), 0);
16 List<scalar> pRefValueFluid(fluidRegions.size(), 0.0);
17 List<bool> frozenFlowFluid(fluidRegions.size(), false);
18 
19 PtrList<dimensionedScalar> rhoMax(fluidRegions.size());
20 PtrList<dimensionedScalar> rhoMin(fluidRegions.size());
21 
22 PtrList<IOMRFZoneList> MRFfluid(fluidRegions.size());
23 PtrList<fv::options> fluidFvOptions(fluidRegions.size());
24 
25 // Populate fluid field pointer lists
27 {
28  Info<< "*** Reading fluid mesh thermophysical properties for region "
29  << fluidRegions[i].name() << nl << endl;
30 
31  Info<< " Adding to thermoFluid\n" << endl;
32 
33  thermoFluid.set
34  (
35  i,
37  );
38 
39  Info<< " Adding to rhoFluid\n" << endl;
40  rhoFluid.set
41  (
42  i,
43  new volScalarField
44  (
45  IOobject
46  (
47  "rho",
48  runTime.timeName(),
49  fluidRegions[i],
50  IOobject::NO_READ,
51  IOobject::AUTO_WRITE
52  ),
53  thermoFluid[i].rho()
54  )
55  );
56 
57  Info<< " Adding to UFluid\n" << endl;
58  UFluid.set
59  (
60  i,
61  new volVectorField
62  (
63  IOobject
64  (
65  "U",
66  runTime.timeName(),
67  fluidRegions[i],
68  IOobject::MUST_READ,
69  IOobject::AUTO_WRITE
70  ),
71  fluidRegions[i]
72  )
73  );
74 
75  Info<< " Adding to phiFluid\n" << endl;
76  phiFluid.set
77  (
78  i,
80  (
81  IOobject
82  (
83  "phi",
84  runTime.timeName(),
85  fluidRegions[i],
86  IOobject::READ_IF_PRESENT,
87  IOobject::AUTO_WRITE
88  ),
90  & fluidRegions[i].Sf()
91  )
92  );
93 
94  Info<< " Adding to gFluid\n" << endl;
95  gFluid.set
96  (
97  i,
99  (
100  IOobject
101  (
102  "g",
103  runTime.constant(),
104  fluidRegions[i],
105  IOobject::MUST_READ,
106  IOobject::NO_WRITE
107  )
108  )
109  );
110 
111  Info<< " Adding to hRefFluid\n" << endl;
112  hRefFluid.set
113  (
114  i,
116  (
117  IOobject
118  (
119  "hRef",
120  runTime.constant(),
121  fluidRegions[i],
122  IOobject::READ_IF_PRESENT,
123  IOobject::NO_WRITE
124  ),
125  dimensionedScalar("hRef", dimLength, 0)
126  )
127  );
128 
129  dimensionedScalar ghRef
130  (
131  mag(gFluid[i].value()) > SMALL
132  ? gFluid[i]
133  & (cmptMag(gFluid[i].value())/mag(gFluid[i].value()))*hRefFluid[i]
134  : dimensionedScalar("ghRef", gFluid[i].dimensions()*dimLength, 0)
135  );
136 
137  Info<< " Adding to ghFluid\n" << endl;
138  ghFluid.set
139  (
140  i,
141  new volScalarField
142  (
143  "gh",
144  (gFluid[i] & fluidRegions[i].C()) - ghRef
145  )
146  );
147 
148  Info<< " Adding to ghfFluid\n" << endl;
149  ghfFluid.set
150  (
151  i,
153  (
154  "ghf",
155  (gFluid[i] & fluidRegions[i].Cf()) - ghRef
156  )
157  );
158 
159  Info<< " Adding to turbulence\n" << endl;
160  turbulence.set
161  (
162  i,
164  (
165  rhoFluid[i],
166  UFluid[i],
167  phiFluid[i],
168  thermoFluid[i]
169  ).ptr()
170  );
171 
172  p_rghFluid.set
173  (
174  i,
175  new volScalarField
176  (
177  IOobject
178  (
179  "p_rgh",
180  runTime.timeName(),
181  fluidRegions[i],
182  IOobject::MUST_READ,
183  IOobject::AUTO_WRITE
184  ),
185  fluidRegions[i]
186  )
187  );
188 
189  // Force p_rgh to be consistent with p
190  p_rghFluid[i] = thermoFluid[i].p() - rhoFluid[i]*ghFluid[i];
191 
192  fluidRegions[i].setFluxRequired(p_rghFluid[i].name());
193 
194  radiation.set
195  (
196  i,
198  );
199 
201 
202  const dictionary& simpleDict =
203  fluidRegions[i].solutionDict().subDict("SIMPLE");
204 
205  setRefCell
206  (
207  thermoFluid[i].p(),
208  p_rghFluid[i],
209  simpleDict,
210  pRefCellFluid[i],
211  pRefValueFluid[i]
212  );
213 
214  simpleDict.readIfPresent("frozenFlow", frozenFlowFluid[i]);
215 
216  rhoMax.set
217  (
218  i,
220  (
221  dimensionedScalar::lookupOrDefault
222  (
223  "rhoMax",
224  simpleDict,
225  dimDensity,
226  GREAT
227  )
228  )
229  );
230 
231  rhoMin.set
232  (
233  i,
235  (
236  dimensionedScalar::lookupOrDefault
237  (
238  "rhoMin",
239  simpleDict,
240  dimDensity,
241  0
242  )
243  )
244  );
245 
246  Info<< " Adding MRF\n" << endl;
247  MRFfluid.set
248  (
249  i,
250  new IOMRFZoneList(fluidRegions[i])
251  );
252 
253  Info<< " Adding fvOptions\n" << endl;
254  fluidFvOptions.set
255  (
256  i,
257  new fv::options(fluidRegions[i])
258  );
259 
260  turbulence[i].validate();
261 }
PtrList< dimensionedScalar > rhoMax(fluidRegions.size())
UniformDimensionedField< vector > uniformDimensionedVectorField
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:107
PtrList< dimensionedScalar > rhoMin(fluidRegions.size())
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
PtrList< IOMRFZoneList > MRFfluid(fluidRegions.size())
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
UniformDimensionedField< scalar > uniformDimensionedScalarField
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
PtrList< compressible::turbulenceModel > turbulence(fluidRegions.size())
PtrList< volScalarField > rhoFluid(fluidRegions.size())
forAll(fluidRegions, i)
List< bool > frozenFlowFluid(fluidRegions.size(), false)
PtrList< radiation::radiationModel > radiation(fluidRegions.size())
PtrList< volScalarField > p_rghFluid(fluidRegions.size())
PtrList< volScalarField > ghFluid(fluidRegions.size())
PtrList< volVectorField > UFluid(fluidRegions.size())
PtrList< rhoThermo > thermoFluid(fluidRegions.size())
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
static const char nl
Definition: Ostream.H:262
volScalarField & C
List< scalar > pRefValueFluid(fluidRegions.size(), 0.0)
const volScalarField & T
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimDensity
PtrList< fv::options > fluidFvOptions(fluidRegions.size())
PtrList< uniformDimensionedVectorField > gFluid(fluidRegions.size())
PtrList< surfaceScalarField > ghfFluid(fluidRegions.size())
List< scalar > initialMassFluid(fluidRegions.size())
PtrList< surfaceScalarField > phiFluid(fluidRegions.size())
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< uniformDimensionedScalarField > hRefFluid(fluidRegions.size())
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
PtrList< fvMesh > fluidRegions(fluidNames.size())
void setRefCell(const volScalarField &field, const volScalarField &fieldRef, const dictionary &dict, label &refCelli, scalar &refValue, const bool forceReference=false)
If the field fieldRef needs referencing find the reference cell nearest.
Definition: findRefCell.C:31
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
List< label > pRefCellFluid(fluidRegions.size(), 0)