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