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<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> turbulenceFluid(fluidRegions.size());
11 PtrList<CombustionModel<rhoReactionThermo>> reactionFluid(fluidRegions.size());
12 PtrList<volScalarField> p_rghFluid(fluidRegions.size());
13 PtrList<radiationModel> radiation(fluidRegions.size());
14 PtrList<volScalarField> KFluid(fluidRegions.size());
15 PtrList<volScalarField> dpdtFluid(fluidRegions.size());
16 PtrList<multivariateSurfaceInterpolationScheme<scalar>::fieldTable>
17  fieldsFluid(fluidRegions.size());
18 
19 List<scalar> initialMassFluid(fluidRegions.size());
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;
32 
33  Info<< " Adding to rhoFluid\n" << endl;
34  rhoFluid.set
35  (
36  i,
37  new volScalarField
38  (
39  IOobject
40  (
41  "rho",
42  runTime.timeName(),
43  fluidRegions[i],
44  IOobject::NO_READ,
45  IOobject::AUTO_WRITE
46  ),
47  thermoFluid[i].rho()
48  )
49  );
50 
51  Info<< " Adding to UFluid\n" << endl;
52  UFluid.set
53  (
54  i,
55  new volVectorField
56  (
57  IOobject
58  (
59  "U",
60  runTime.timeName(),
61  fluidRegions[i],
62  IOobject::MUST_READ,
63  IOobject::AUTO_WRITE
64  ),
65  fluidRegions[i]
66  )
67  );
68 
69  Info<< " Adding to phiFluid\n" << endl;
70  phiFluid.set
71  (
72  i,
74  (
75  IOobject
76  (
77  "phi",
78  runTime.timeName(),
79  fluidRegions[i],
80  IOobject::READ_IF_PRESENT,
81  IOobject::AUTO_WRITE
82  ),
84  & fluidRegions[i].Sf()
85  )
86  );
87 
88  Info<< " Adding to gFluid\n" << endl;
89  gFluid.set
90  (
91  i,
93  (
94  IOobject
95  (
96  "g",
97  runTime.constant(),
98  fluidRegions[i],
99  IOobject::MUST_READ,
100  IOobject::NO_WRITE
101  )
102  )
103  );
104 
105  Info<< " Adding to hRefFluid\n" << endl;
106  hRefFluid.set
107  (
108  i,
110  (
111  IOobject
112  (
113  "hRef",
114  runTime.constant(),
115  fluidRegions[i],
116  IOobject::READ_IF_PRESENT,
117  IOobject::NO_WRITE
118  ),
120  )
121  );
122 
123  dimensionedScalar ghRef(- mag(gFluid[i])*hRefFluid[i]);
124 
125  Info<< " Adding to ghFluid\n" << endl;
126  ghFluid.set
127  (
128  i,
129  new volScalarField
130  (
131  "gh",
132  (gFluid[i] & fluidRegions[i].C()) - ghRef
133  )
134  );
135 
136  Info<< " Adding to ghfFluid\n" << endl;
137  ghfFluid.set
138  (
139  i,
141  (
142  "ghf",
143  (gFluid[i] & fluidRegions[i].Cf()) - ghRef
144  )
145  );
146 
147  Info<< " Adding to turbulenceFluid\n" << endl;
148  turbulenceFluid.set
149  (
150  i,
152  (
153  rhoFluid[i],
154  UFluid[i],
155  phiFluid[i],
156  thermoFluid[i]
157  ).ptr()
158  );
159 
160  Info<< " Adding to reactionFluid\n" << endl;
161  reactionFluid.set
162  (
163  i,
165  (
166  thermoFluid[i],
167  turbulenceFluid[i]
168  )
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  Info<< " Adding to radiationFluid\n" << endl;
194  radiation.set
195  (
196  i,
198  );
199 
201 
202  Info<< " Adding to KFluid\n" << endl;
203  KFluid.set
204  (
205  i,
206  new volScalarField
207  (
208  "K",
209  0.5*magSqr(UFluid[i])
210  )
211  );
212 
213  Info<< " Adding to dpdtFluid\n" << endl;
214  dpdtFluid.set
215  (
216  i,
217  new volScalarField
218  (
219  IOobject
220  (
221  "dpdt",
222  runTime.timeName(),
223  fluidRegions[i]
224  ),
225  fluidRegions[i],
227  (
228  thermoFluid[i].p().dimensions()/dimTime,
229  0
230  )
231  )
232  );
233 
234  Info<< " Adding to fieldsFluid\n" << endl;
235  fieldsFluid.set
236  (
237  i,
238  new multivariateSurfaceInterpolationScheme<scalar>::fieldTable
239  );
240  forAll(thermoFluid[i].composition().Y(), j)
241  {
242  fieldsFluid[i].add(thermoFluid[i].composition().Y()[j]);
243  }
244  fieldsFluid[i].add(thermoFluid[i].he());
245 
246  Info<< " Adding MRF\n" << endl;
247  MRFfluid.set
248  (
249  i,
250  new IOMRFZoneList(fluidRegions[i])
251  );
252 
253  Info<< " Adding fvOptions\n" << endl;
254  fluidFvOptions.set
255  (
256  i,
257  new fv::options(fluidRegions[i])
258  );
259 
260  turbulenceFluid[i].validate();
261 }
volScalarField & he
Definition: YEEqn.H:50
PtrList< uniformDimensionedVectorField > gFluid(fluidRegions.size())
PtrList< volVectorField > UFluid(fluidRegions.size())
PtrList< fv::options > fluidFvOptions(fluidRegions.size())
PtrList< IOMRFZoneList > MRFfluid(fluidRegions.size())
UniformDimensionedField< vector > uniformDimensionedVectorField
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:108
basicSpecieMixture & composition
PtrList< volScalarField > p_rghFluid(fluidRegions.size())
PtrList< uniformDimensionedScalarField > hRefFluid(fluidRegions.size())
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
UniformDimensionedField< scalar > uniformDimensionedScalarField
PtrList< compressible::turbulenceModel > turbulenceFluid(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
List< scalar > initialMassFluid(fluidRegions.size())
PtrList< volScalarField > ghFluid(fluidRegions.size())
PtrList< rhoReactionThermo > thermoFluid(fluidRegions.size())
PtrList< volScalarField > KFluid(fluidRegions.size())
PtrList< volScalarField > dpdtFluid(fluidRegions.size())
PtrList< multivariateSurfaceInterpolationScheme< scalar >::fieldTable > fieldsFluid(fluidRegions.size())
PtrList< fvMesh > fluidRegions(fluidNames.size())
dimensioned< scalar > magSqr(const dimensioned< Type > &)
PtrList< surfaceScalarField > ghfFluid(fluidRegions.size())
static const char nl
Definition: Ostream.H:265
forAll(fluidRegions, i)
const volScalarField & T
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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< volScalarField > & Y
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
PtrList< radiationModel > radiation(fluidRegions.size())
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
PtrList< surfaceScalarField > phiFluid(fluidRegions.size())
PtrList< CombustionModel< rhoReactionThermo > > reactionFluid(fluidRegions.size())
PtrList< volScalarField > rhoFluid(fluidRegions.size())