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 PtrList<fvVectorMatrix> UEqFluid(fluidRegions.size());
20 
21 List<scalar> initialMassFluid(fluidRegions.size());
22 List<bool> frozenFlowFluid(fluidRegions.size(), false);
23 
24 PtrList<IOMRFZoneList> MRFfluid(fluidRegions.size());
25 PtrList<fv::options> fluidFvOptions(fluidRegions.size());
26 
27 List<label> pRefCellFluid(fluidRegions.size());
28 List<scalar> pRefValueFluid(fluidRegions.size());
29 
30 PtrList<dimensionedScalar> rhoMinFluid(fluidRegions.size());
31 PtrList<dimensionedScalar> rhoMaxFluid(fluidRegions.size());
32 
33 PtrList<pressureControl> pressureControls(fluidRegions.size());
34 
36 
37 // Populate fluid field pointer lists
39 {
40  Info<< "*** Reading fluid mesh thermophysical properties for region "
41  << fluidRegions[i].name() << nl << endl;
42 
43  Info<< " Adding to thermoFluid\n" << endl;
45 
46  Info<< " Adding to rhoFluid\n" << endl;
47  rhoFluid.set
48  (
49  i,
50  new volScalarField
51  (
52  IOobject
53  (
54  "rho",
55  runTime.timeName(),
56  fluidRegions[i],
57  IOobject::NO_READ,
58  IOobject::AUTO_WRITE
59  ),
60  thermoFluid[i].rho()
61  )
62  );
63 
64  Info<< " Adding to UFluid\n" << endl;
65  UFluid.set
66  (
67  i,
68  new volVectorField
69  (
70  IOobject
71  (
72  "U",
73  runTime.timeName(),
74  fluidRegions[i],
75  IOobject::MUST_READ,
76  IOobject::AUTO_WRITE
77  ),
78  fluidRegions[i]
79  )
80  );
81 
82  Info<< " Adding to phiFluid\n" << endl;
83  phiFluid.set
84  (
85  i,
87  (
88  IOobject
89  (
90  "phi",
91  runTime.timeName(),
92  fluidRegions[i],
93  IOobject::READ_IF_PRESENT,
94  IOobject::AUTO_WRITE
95  ),
97  & fluidRegions[i].Sf()
98  )
99  );
100 
101  Info<< " Adding to hRefFluid\n" << endl;
102  hRefFluid.set
103  (
104  i,
106  (
107  IOobject
108  (
109  "hRef",
110  runTime.constant(),
111  fluidRegions[i],
112  IOobject::READ_IF_PRESENT,
113  IOobject::NO_WRITE
114  ),
115  dimensionedScalar("hRef", dimLength, Zero) // uses name
116  )
117  );
118 
119  dimensionedScalar ghRef
120  (
121  mag(g.value()) > SMALL
122  ? g & (cmptMag(g.value())/mag(g.value()))*hRefFluid[i]
123  : dimensionedScalar("ghRef", g.dimensions()*dimLength, 0)
124  );
125 
126  Info<< " Adding to ghFluid\n" << endl;
127  ghFluid.set
128  (
129  i,
130  new volScalarField
131  (
132  "gh",
133  (g & fluidRegions[i].C()) - ghRef
134  )
135  );
136 
137  Info<< " Adding to ghfFluid\n" << endl;
138  ghfFluid.set
139  (
140  i,
142  (
143  "ghf",
144  (g & fluidRegions[i].Cf()) - ghRef
145  )
146  );
147 
148  Info<< " Adding to turbulenceFluid\n" << endl;
149  turbulenceFluid.set
150  (
151  i,
153  (
154  rhoFluid[i],
155  UFluid[i],
156  phiFluid[i],
157  thermoFluid[i]
158  ).ptr()
159  );
160 
161  Info<< " Adding to reactionFluid\n" << endl;
162  reactionFluid.set
163  (
164  i,
166  (
167  thermoFluid[i],
168  turbulenceFluid[i]
169  )
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  Info<< " Adding to radiationFluid\n" << endl;
195  radiation.set
196  (
197  i,
199  );
200 
202 
203  Info<< " Adding to KFluid\n" << endl;
204  KFluid.set
205  (
206  i,
207  new volScalarField
208  (
209  "K",
210  0.5*magSqr(UFluid[i])
211  )
212  );
213 
214  Info<< " Adding to dpdtFluid\n" << endl;
215  dpdtFluid.set
216  (
217  i,
218  new volScalarField
219  (
220  IOobject
221  (
222  "dpdt",
223  runTime.timeName(),
224  fluidRegions[i]
225  ),
226  fluidRegions[i],
227  dimensionedScalar(thermoFluid[i].p().dimensions()/dimTime, Zero)
228  )
229  );
230 
231  Info<< " Adding to fieldsFluid\n" << endl;
232  fieldsFluid.set
233  (
234  i,
235  new multivariateSurfaceInterpolationScheme<scalar>::fieldTable
236  );
237  forAll(thermoFluid[i].composition().Y(), j)
238  {
239  fieldsFluid[i].add(thermoFluid[i].composition().Y()[j]);
240  }
241  fieldsFluid[i].add(thermoFluid[i].he());
242 
243  Info<< " Adding to QdotFluid\n" << endl;
244  QdotFluid.set
245  (
246  i,
247  new volScalarField
248  (
249  IOobject
250  (
251  "Qdot",
252  runTime.timeName(),
253  fluidRegions[i],
254  IOobject::READ_IF_PRESENT,
255  IOobject::AUTO_WRITE
256  ),
257  fluidRegions[i],
259  )
260  );
261 
262  const dictionary& pimpleDict =
263  fluidRegions[i].solutionDict().subDict("PIMPLE");
264  pimpleDict.readIfPresent("frozenFlow", frozenFlowFluid[i]);
265 
266  rhoMaxFluid.set
267  (
268  i,
269  new dimensionedScalar("rhoMax", dimDensity, GREAT, pimpleDict)
270  );
271 
272  rhoMinFluid.set
273  (
274  i,
276  );
277 
278  pressureControls.set
279  (
280  i,
281  new pressureControl(thermoFluid[i].p(), rhoFluid[i], pimpleDict, false)
282  );
283 
284  Info<< " Adding MRF\n" << endl;
285  MRFfluid.set
286  (
287  i,
288  new IOMRFZoneList(fluidRegions[i])
289  );
290 
291  Info<< " Adding fvOptions\n" << endl;
292  fluidFvOptions.set
293  (
294  i,
295  new fv::options(fluidRegions[i])
296  );
297 
298  UEqFluid.set
299  (
300  i,
302  );
303 
304  turbulenceFluid[i].validate();
305 
306  pRefCellFluid[i] = -1;
307  pRefValueFluid[i] = 0.0;
308 
309  if (p_rghFluid[i].needReference())
310  {
311  setRefCell
312  (
313  thermoFluid[i].p(),
314  p_rghFluid[i],
315  pimpleDict,
316  pRefCellFluid[i],
317  pRefValueFluid[i]
318  );
319  }
320 
321 }
runTime
engineTime & runTime
Definition: createEngineTime.H:13
p
volScalarField & p
Definition: createFieldRefs.H:8
UEqFluid
PtrList< fvVectorMatrix > UEqFluid(fluidRegions.size())
rhoMaxFluid
PtrList< dimensionedScalar > rhoMaxFluid(fluidRegions.size())
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
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:369
Foam::dimForce
const dimensionSet dimForce
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:53
pressureControls
PtrList< pressureControl > pressureControls(fluidRegions.size())
rhoFluid
PtrList< volScalarField > rhoFluid(fluidRegions.size())
Foam::fvVectorMatrix
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:47
MRFfluid
PtrList< IOMRFZoneList > MRFfluid(fluidRegions.size())
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
fluidFvOptions
PtrList< fv::options > fluidFvOptions(fluidRegions.size())
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:42
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:26
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:29
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:404
pimpleDict
const dictionary & pimpleDict
Definition: createRegionControls.H:3
pressureControl
const pressureControl & pressureControl
Definition: setRegionFluidFields.H:69
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())
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
KFluid
PtrList< volScalarField > KFluid(fluidRegions.size())
dpdtFluid
PtrList< volScalarField > dpdtFluid(fluidRegions.size())
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
thermoFluid
PtrList< rhoThermo > thermoFluid(fluidRegions.size())
initialMassFluid
List< scalar > initialMassFluid(fluidRegions.size())