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<uniformDimensionedScalarField> hRefFluid(fluidRegions.size());
7 PtrList<volScalarField> ghFluid(fluidRegions.size());
8 PtrList<surfaceScalarField> ghfFluid(fluidRegions.size());
9 PtrList<compressible::turbulenceModel> turbulence(fluidRegions.size());
10 PtrList<volScalarField> p_rghFluid(fluidRegions.size());
11 PtrList<radiation::radiationModel> radiation(fluidRegions.size());
12 
13 List<scalar> initialMassFluid(fluidRegions.size());
14 List<label> pRefCellFluid(fluidRegions.size(), -1);
15 List<scalar> pRefValueFluid(fluidRegions.size(), Zero);
16 List<bool> frozenFlowFluid(fluidRegions.size(), false);
17 
18 PtrList<dimensionedScalar> rhoMax(fluidRegions.size());
19 PtrList<dimensionedScalar> rhoMin(fluidRegions.size());
20 
21 PtrList<IOMRFZoneList> MRFfluid(fluidRegions.size());
22 PtrList<fv::options> fluidFvOptions(fluidRegions.size());
23 
25 
26 // Populate fluid field pointer lists
28 {
29  Info<< "*** Reading fluid mesh thermophysical properties for region "
30  << fluidRegions[i].name() << nl << endl;
31 
32  Info<< " Adding to thermoFluid\n" << endl;
33 
34  thermoFluid.set
35  (
36  i,
38  );
39 
40  Info<< " Adding to rhoFluid\n" << endl;
41  rhoFluid.set
42  (
43  i,
44  new volScalarField
45  (
46  IOobject
47  (
48  "rho",
49  runTime.timeName(),
50  fluidRegions[i],
51  IOobject::NO_READ,
52  IOobject::AUTO_WRITE
53  ),
54  thermoFluid[i].rho()
55  )
56  );
57 
58  Info<< " Adding to UFluid\n" << endl;
59  UFluid.set
60  (
61  i,
62  new volVectorField
63  (
64  IOobject
65  (
66  "U",
67  runTime.timeName(),
68  fluidRegions[i],
69  IOobject::MUST_READ,
70  IOobject::AUTO_WRITE
71  ),
72  fluidRegions[i]
73  )
74  );
75 
76  Info<< " Adding to phiFluid\n" << endl;
77  phiFluid.set
78  (
79  i,
81  (
82  IOobject
83  (
84  "phi",
85  runTime.timeName(),
86  fluidRegions[i],
87  IOobject::READ_IF_PRESENT,
88  IOobject::AUTO_WRITE
89  ),
91  & fluidRegions[i].Sf()
92  )
93  );
94 
95  Info<< " Adding to hRefFluid\n" << endl;
96  hRefFluid.set
97  (
98  i,
100  (
101  IOobject
102  (
103  "hRef",
104  runTime.constant(),
105  fluidRegions[i],
106  IOobject::READ_IF_PRESENT,
107  IOobject::NO_WRITE
108  ),
109  dimensionedScalar("hRef", dimLength, Zero) // uses name
110  )
111  );
112 
113  dimensionedScalar ghRef
114  (
115  mag(g.value()) > SMALL
116  ? g & (cmptMag(g.value())/mag(g.value()))*hRefFluid[i]
117  : dimensionedScalar("ghRef", g.dimensions()*dimLength, 0)
118  );
119 
120  Info<< " Adding to ghFluid\n" << endl;
121  ghFluid.set
122  (
123  i,
124  new volScalarField
125  (
126  "gh",
127  (g & fluidRegions[i].C()) - ghRef
128  )
129  );
130 
131  Info<< " Adding to ghfFluid\n" << endl;
132  ghfFluid.set
133  (
134  i,
136  (
137  "ghf",
138  (g & fluidRegions[i].Cf()) - ghRef
139  )
140  );
141 
142  Info<< " Adding to turbulence\n" << endl;
143  turbulence.set
144  (
145  i,
147  (
148  rhoFluid[i],
149  UFluid[i],
150  phiFluid[i],
151  thermoFluid[i]
152  ).ptr()
153  );
154 
155  p_rghFluid.set
156  (
157  i,
158  new volScalarField
159  (
160  IOobject
161  (
162  "p_rgh",
163  runTime.timeName(),
164  fluidRegions[i],
165  IOobject::MUST_READ,
166  IOobject::AUTO_WRITE
167  ),
168  fluidRegions[i]
169  )
170  );
171 
172  // Force p_rgh to be consistent with p
173  p_rghFluid[i] = thermoFluid[i].p() - rhoFluid[i]*ghFluid[i];
174 
175  fluidRegions[i].setFluxRequired(p_rghFluid[i].name());
176 
177  radiation.set
178  (
179  i,
181  );
182 
184 
185  const dictionary& simpleDict =
186  fluidRegions[i].solutionDict().subDict("SIMPLE");
187 
188  setRefCell
189  (
190  thermoFluid[i].p(),
191  p_rghFluid[i],
192  simpleDict,
193  pRefCellFluid[i],
194  pRefValueFluid[i]
195  );
196 
197  simpleDict.readIfPresent("frozenFlow", frozenFlowFluid[i]);
198 
199  rhoMax.set
200  (
201  i,
202  new dimensionedScalar("rhoMax", dimDensity, GREAT, simpleDict)
203  );
204 
205  rhoMin.set
206  (
207  i,
208  new dimensionedScalar("rhoMin", dimDensity, Zero, simpleDict)
209  );
210 
211  Info<< " Adding MRF\n" << endl;
212  MRFfluid.set
213  (
214  i,
215  new IOMRFZoneList(fluidRegions[i])
216  );
217 
218  Info<< " Adding fvOptions\n" << endl;
219  fluidFvOptions.set
220  (
221  i,
222  new fv::options(fluidRegions[i])
223  );
224 
225  turbulence[i].validate();
226 }
runTime
engineTime & runTime
Definition: createEngineTime.H:13
p
volScalarField & p
Definition: createFieldRefs.H:8
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::dimDensity
const dimensionSet dimDensity
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)
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
C
volScalarField & C
Definition: readThermalProperties.H:102
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:400
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
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
Foam::nl
constexpr char nl
Definition: Ostream.H:385
fluidRegions
PtrList< fvMesh > fluidRegions(fluidNames.size())
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
turbulence
PtrList< compressible::turbulenceModel > turbulence(fluidRegions.size())
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
ghfFluid
PtrList< surfaceScalarField > ghfFluid(fluidRegions.size())
rhoMax
PtrList< dimensionedScalar > rhoMax(fluidRegions.size())
thermoFluid
PtrList< rhoThermo > thermoFluid(fluidRegions.size())
initialMassFluid
List< scalar > initialMassFluid(fluidRegions.size())
rhoMin
PtrList< dimensionedScalar > rhoMin(fluidRegions.size())