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 PtrList<volScalarField> KFluid(fluidRegions.size());
14 PtrList<volScalarField> dpdtFluid(fluidRegions.size());
15 
16 List<scalar> initialMassFluid(fluidRegions.size());
17 List<bool> frozenFlowFluid(fluidRegions.size(), false);
18 List<bool> residualReachedFluid(fluidRegions.size(), true);
19 List<bool> residualControlUsedFluid(fluidRegions.size(), false);
20 
21 PtrList<IOMRFZoneList> MRFfluid(fluidRegions.size());
22 PtrList<fv::options> fluidFvOptions(fluidRegions.size());
23 
24 // Populate fluid field pointer lists
26 {
27  Info<< "*** Reading fluid mesh thermophysical properties for region "
28  << fluidRegions[i].name() << nl << endl;
29 
30  Info<< " Adding to thermoFluid\n" << endl;
31 
32  thermoFluid.set
33  (
34  i,
36  );
37 
38  Info<< " Adding to rhoFluid\n" << endl;
39  rhoFluid.set
40  (
41  i,
42  new volScalarField
43  (
44  IOobject
45  (
46  "rho",
47  runTime.timeName(),
48  fluidRegions[i],
49  IOobject::NO_READ,
50  IOobject::AUTO_WRITE
51  ),
52  thermoFluid[i].rho()
53  )
54  );
55 
56  Info<< " Adding to UFluid\n" << endl;
57  UFluid.set
58  (
59  i,
60  new volVectorField
61  (
62  IOobject
63  (
64  "U",
65  runTime.timeName(),
66  fluidRegions[i],
67  IOobject::MUST_READ,
68  IOobject::AUTO_WRITE
69  ),
70  fluidRegions[i]
71  )
72  );
73 
74  Info<< " Adding to phiFluid\n" << endl;
75  phiFluid.set
76  (
77  i,
79  (
80  IOobject
81  (
82  "phi",
83  runTime.timeName(),
84  fluidRegions[i],
85  IOobject::READ_IF_PRESENT,
86  IOobject::AUTO_WRITE
87  ),
89  & fluidRegions[i].Sf()
90  )
91  );
92 
93  Info<< " Adding to gFluid\n" << endl;
94  gFluid.set
95  (
96  i,
98  (
99  IOobject
100  (
101  "g",
102  runTime.constant(),
103  fluidRegions[i],
104  IOobject::MUST_READ,
105  IOobject::NO_WRITE
106  )
107  )
108  );
109 
110  Info<< " Adding to hRefFluid\n" << endl;
111  hRefFluid.set
112  (
113  i,
115  (
116  IOobject
117  (
118  "hRef",
119  runTime.constant(),
120  fluidRegions[i],
121  IOobject::READ_IF_PRESENT,
122  IOobject::NO_WRITE
123  ),
124  dimensionedScalar("hRef", dimLength, 0)
125  )
126  );
127 
128  dimensionedScalar ghRef
129  (
130  mag(gFluid[i].value()) > SMALL
131  ? gFluid[i]
132  & (cmptMag(gFluid[i].value())/mag(gFluid[i].value()))*hRefFluid[i]
133  : dimensionedScalar("ghRef", gFluid[i].dimensions()*dimLength, 0)
134  );
135 
136  Info<< " Adding to ghFluid\n" << endl;
137  ghFluid.set
138  (
139  i,
140  new volScalarField
141  (
142  "gh",
143  (gFluid[i] & fluidRegions[i].C()) - ghRef
144  )
145  );
146 
147  Info<< " Adding to ghfFluid\n" << endl;
148  ghfFluid.set
149  (
150  i,
152  (
153  "ghf",
154  (gFluid[i] & fluidRegions[i].Cf()) - ghRef
155  )
156  );
157 
158  Info<< " Adding to turbulence\n" << endl;
159  turbulence.set
160  (
161  i,
163  (
164  rhoFluid[i],
165  UFluid[i],
166  phiFluid[i],
167  thermoFluid[i]
168  ).ptr()
169  );
170 
171  p_rghFluid.set
172  (
173  i,
174  new volScalarField
175  (
176  IOobject
177  (
178  "p_rgh",
179  runTime.timeName(),
180  fluidRegions[i],
181  IOobject::MUST_READ,
182  IOobject::AUTO_WRITE
183  ),
184  fluidRegions[i]
185  )
186  );
187 
188  // Force p_rgh to be consistent with p
189  p_rghFluid[i] = thermoFluid[i].p() - rhoFluid[i]*ghFluid[i];
190 
191  fluidRegions[i].setFluxRequired(p_rghFluid[i].name());
192 
193  radiation.set
194  (
195  i,
197  );
198 
200 
201  Info<< " Adding to KFluid\n" << endl;
202  KFluid.set
203  (
204  i,
205  new volScalarField
206  (
207  "K",
208  0.5*magSqr(UFluid[i])
209  )
210  );
211 
212  Info<< " Adding to dpdtFluid\n" << endl;
213  dpdtFluid.set
214  (
215  i,
216  new volScalarField
217  (
218  IOobject
219  (
220  "dpdt",
221  runTime.timeName(),
222  fluidRegions[i]
223  ),
224  fluidRegions[i],
226  (
227  "dpdt",
228  thermoFluid[i].p().dimensions()/dimTime,
229  0
230  )
231  )
232  );
233 
234  const dictionary& pimpleDict =
235  fluidRegions[i].solutionDict().subDict("PIMPLE");
236  pimpleDict.readIfPresent("frozenFlow", frozenFlowFluid[i]);
237 
238  Info<< " Adding MRF\n" << endl;
239  MRFfluid.set
240  (
241  i,
242  new IOMRFZoneList(fluidRegions[i])
243  );
244 
245  Info<< " Adding fvOptions\n" << endl;
246  fluidFvOptions.set
247  (
248  i,
249  new fv::options(fluidRegions[i])
250  );
251 
252  turbulence[i].validate();
253 }
PtrList< volScalarField > dpdtFluid(fluidRegions.size())
UniformDimensionedField< vector > uniformDimensionedVectorField
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:107
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
PtrList< volScalarField > KFluid(fluidRegions.size())
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())
List< bool > residualControlUsedFluid(fluidRegions.size(), false)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const char nl
Definition: Ostream.H:262
volScalarField & C
const volScalarField & T
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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())
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const dictionary & pimpleDict
Definition: setRDeltaT.H:29
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
PtrList< fvMesh > fluidRegions(fluidNames.size())
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
List< bool > residualReachedFluid(fluidRegions.size(), true)