createFluidFields.H
Go to the documentation of this file.
1 // Initialise fluid field pointer lists
2 PtrList<rhoReactionThermo> thermoFluid(fluidRegions.size());
3 PtrList<volScalarField> rhoFluid(fluidRegions.size());
4 PtrList<volVectorField> UFluid(fluidRegions.size());
5 PtrList<surfaceScalarField> phiFluid(fluidRegions.size());
6 PtrList<uniformDimensionedScalarField> hRefFluid(fluidRegions.size());
7 PtrList<volScalarField> ghFluid(fluidRegions.size());
8 PtrList<surfaceScalarField> ghfFluid(fluidRegions.size());
9 PtrList<compressible::turbulenceModel> turbulenceFluid(fluidRegions.size());
10 PtrList<CombustionModel<rhoReactionThermo>> reactionFluid(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 PtrList<multivariateSurfaceInterpolationScheme<scalar>::fieldTable>
16  fieldsFluid(fluidRegions.size());
17 PtrList<volScalarField> QdotFluid(fluidRegions.size());
18 
19 List<scalar> initialMassFluid(fluidRegions.size());
20 List<bool> frozenFlowFluid(fluidRegions.size(), false);
21 
22 PtrList<IOMRFZoneList> MRFfluid(fluidRegions.size());
23 PtrList<fv::options> fluidFvOptions(fluidRegions.size());
24 
25 List<label> pRefCellFluid(fluidRegions.size());
26 List<scalar> pRefValueFluid(fluidRegions.size());
27 
28 PtrList<dimensionedScalar> rhoMinFluid(fluidRegions.size());
29 PtrList<dimensionedScalar> rhoMaxFluid(fluidRegions.size());
30 
31 PtrList<pressureControl> pressureControls(fluidRegions.size());
32 
34 
35 // Populate fluid field pointer lists
37 {
38  Info<< "*** Reading fluid mesh thermophysical properties for region "
39  << fluidRegions[i].name() << nl << endl;
40 
41  Info<< " Adding to thermoFluid\n" << endl;
43 
44  Info<< " Adding to rhoFluid\n" << endl;
45  rhoFluid.set
46  (
47  i,
48  new volScalarField
49  (
50  IOobject
51  (
52  "rho",
53  runTime.timeName(),
54  fluidRegions[i],
55  IOobject::NO_READ,
56  IOobject::AUTO_WRITE
57  ),
58  thermoFluid[i].rho()
59  )
60  );
61 
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  );
79 
80  Info<< " Adding to phiFluid\n" << endl;
81  phiFluid.set
82  (
83  i,
85  (
86  IOobject
87  (
88  "phi",
89  runTime.timeName(),
90  fluidRegions[i],
91  IOobject::READ_IF_PRESENT,
92  IOobject::AUTO_WRITE
93  ),
95  & fluidRegions[i].Sf()
96  )
97  );
98 
99  Info<< " Adding to hRefFluid\n" << endl;
100  hRefFluid.set
101  (
102  i,
104  (
105  IOobject
106  (
107  "hRef",
108  runTime.constant(),
109  fluidRegions[i],
110  IOobject::READ_IF_PRESENT,
111  IOobject::NO_WRITE
112  ),
113  dimensionedScalar("hRef", dimLength, Zero) // uses name
114  )
115  );
116 
117  dimensionedScalar ghRef
118  (
119  mag(g.value()) > SMALL
120  ? g & (cmptMag(g.value())/mag(g.value()))*hRefFluid[i]
121  : dimensionedScalar("ghRef", g.dimensions()*dimLength, 0)
122  );
123 
124  Info<< " Adding to ghFluid\n" << endl;
125  ghFluid.set
126  (
127  i,
128  new volScalarField
129  (
130  "gh",
131  (g & fluidRegions[i].C()) - ghRef
132  )
133  );
134 
135  Info<< " Adding to ghfFluid\n" << endl;
136  ghfFluid.set
137  (
138  i,
140  (
141  "ghf",
142  (g & fluidRegions[i].Cf()) - ghRef
143  )
144  );
145 
146  Info<< " Adding to turbulenceFluid\n" << endl;
147  turbulenceFluid.set
148  (
149  i,
151  (
152  rhoFluid[i],
153  UFluid[i],
154  phiFluid[i],
155  thermoFluid[i]
156  ).ptr()
157  );
158 
159  Info<< " Adding to reactionFluid\n" << endl;
160  reactionFluid.set
161  (
162  i,
164  (
165  thermoFluid[i],
166  turbulenceFluid[i]
167  )
168  );
169 
170  p_rghFluid.set
171  (
172  i,
173  new volScalarField
174  (
175  IOobject
176  (
177  "p_rgh",
178  runTime.timeName(),
179  fluidRegions[i],
180  IOobject::MUST_READ,
181  IOobject::AUTO_WRITE
182  ),
183  fluidRegions[i]
184  )
185  );
186 
187  // Force p_rgh to be consistent with p
188  p_rghFluid[i] = thermoFluid[i].p() - rhoFluid[i]*ghFluid[i];
189 
190  fluidRegions[i].setFluxRequired(p_rghFluid[i].name());
191 
192  Info<< " Adding to radiationFluid\n" << endl;
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],
225  dimensionedScalar(thermoFluid[i].p().dimensions()/dimTime, Zero)
226  )
227  );
228 
229  Info<< " Adding to fieldsFluid\n" << endl;
230  fieldsFluid.set
231  (
232  i,
233  new multivariateSurfaceInterpolationScheme<scalar>::fieldTable
234  );
235  forAll(thermoFluid[i].composition().Y(), j)
236  {
237  fieldsFluid[i].add(thermoFluid[i].composition().Y()[j]);
238  }
239  fieldsFluid[i].add(thermoFluid[i].he());
240 
241  Info<< " Adding to QdotFluid\n" << endl;
242  QdotFluid.set
243  (
244  i,
245  new volScalarField
246  (
247  IOobject
248  (
249  "Qdot",
250  runTime.timeName(),
251  fluidRegions[i],
252  IOobject::READ_IF_PRESENT,
253  IOobject::AUTO_WRITE
254  ),
255  fluidRegions[i],
257  )
258  );
259 
260  const dictionary& pimpleDict =
261  fluidRegions[i].solutionDict().subDict("PIMPLE");
262  pimpleDict.readIfPresent("frozenFlow", frozenFlowFluid[i]);
263 
264  rhoMaxFluid.set
265  (
266  i,
267  new dimensionedScalar("rhoMax", dimDensity, GREAT, pimpleDict)
268  );
269 
270  rhoMinFluid.set
271  (
272  i,
274  );
275 
276  pressureControls.set
277  (
278  i,
279  new pressureControl(thermoFluid[i].p(), rhoFluid[i], pimpleDict, false)
280  );
281 
282  Info<< " Adding MRF\n" << endl;
283  MRFfluid.set
284  (
285  i,
286  new IOMRFZoneList(fluidRegions[i])
287  );
288 
289  Info<< " Adding fvOptions\n" << endl;
290  fluidFvOptions.set
291  (
292  i,
293  new fv::options(fluidRegions[i])
294  );
295 
296  turbulenceFluid[i].validate();
297 
298  pRefCellFluid[i] = -1;
299  pRefValueFluid[i] = 0.0;
300 
301  if (p_rghFluid[i].needReference())
302  {
303  setRefCell
304  (
305  thermoFluid[i].p(),
306  p_rghFluid[i],
307  pimpleDict,
308  pRefCellFluid[i],
309  pRefValueFluid[i]
310  );
311  }
312 
313 }
runTime
engineTime & runTime
Definition: createEngineTime.H:13
p
volScalarField & p
Definition: createFieldRefs.H:8
rhoMaxFluid
PtrList< dimensionedScalar > rhoMaxFluid(fluidRegions.size())
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::dimDensity
const dimensionSet dimDensity
reactionFluid
PtrList< CombustionModel< rhoReactionThermo > > reactionFluid(fluidRegions.size())
QdotFluid
PtrList< volScalarField > QdotFluid(fluidRegions.size())
Foam::fvc::domainIntegrate
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcVolumeIntegrate.C:88
UFluid
PtrList< volVectorField > UFluid(fluidRegions.size())
p_rghFluid
PtrList< volScalarField > p_rghFluid(fluidRegions.size())
frozenFlowFluid
List< bool > frozenFlowFluid(fluidRegions.size(), false)
fieldsFluid
PtrList< multivariateSurfaceInterpolationScheme< scalar >::fieldTable > fieldsFluid(fluidRegions.size())
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
radiation
PtrList< radiation::radiationModel > radiation(fluidRegions.size())
rho
rho
Definition: readInitialConditions.H:88
turbulenceFluid
PtrList< compressible::turbulenceModel > turbulenceFluid(fluidRegions.size())
composition
basicSpecieMixture & composition
Definition: createFieldRefs.H:6
C
volScalarField & C
Definition: readThermalProperties.H:102
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:400
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
pressureControls
PtrList< pressureControl > pressureControls(fluidRegions.size())
rhoFluid
PtrList< volScalarField > rhoFluid(fluidRegions.size())
MRFfluid
PtrList< IOMRFZoneList > MRFfluid(fluidRegions.size())
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
fluidFvOptions
PtrList< fv::options > fluidFvOptions(fluidRegions.size())
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::uniformDimensionedVectorField
UniformDimensionedField< vector > uniformDimensionedVectorField
Definition: uniformDimensionedFields.H:50
Foam::uniformDimensionedScalarField
UniformDimensionedField< scalar > uniformDimensionedScalarField
Definition: uniformDimensionedFields.H:49
Foam::setRefCell
bool 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:34
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
hRefFluid
PtrList< uniformDimensionedScalarField > hRefFluid(fluidRegions.size())
phiFluid
PtrList< surfaceScalarField > phiFluid(fluidRegions.size())
pRefCellFluid
List< label > pRefCellFluid(fluidRegions.size(), -1)
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:24
T
const volScalarField & T
Definition: createFieldRefs.H:2
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
ghFluid
PtrList< volScalarField > ghFluid(fluidRegions.size())
pRefValueFluid
List< scalar > pRefValueFluid(fluidRegions.size(), Zero)
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
forAll
forAll(fluidRegions, i)
Definition: createFluidFields.H:27
Foam::linearInterpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:112
he
volScalarField & he
Definition: YEEqn.H:52
Foam::nl
constexpr char nl
Definition: Ostream.H:385
pimpleDict
const dictionary & pimpleDict
Definition: createRegionControls.H:3
pressureControl
const pressureControl & pressureControl
Definition: setRegionFluidFields.H:67
fluidRegions
PtrList< fvMesh > fluidRegions(fluidNames.size())
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
ghfFluid
PtrList< surfaceScalarField > ghfFluid(fluidRegions.size())
rhoMinFluid
PtrList< dimensionedScalar > rhoMinFluid(fluidRegions.size())
KFluid
PtrList< volScalarField > KFluid(fluidRegions.size())
dpdtFluid
PtrList< volScalarField > dpdtFluid(fluidRegions.size())
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
thermoFluid
PtrList< rhoThermo > thermoFluid(fluidRegions.size())
initialMassFluid
List< scalar > initialMassFluid(fluidRegions.size())