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<radiation::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 PtrList<volScalarField> QdotFluid(fluidRegions.size());
19 
20 List<scalar> initialMassFluid(fluidRegions.size());
21 
22 PtrList<IOMRFZoneList> MRFfluid(fluidRegions.size());
23 PtrList<fv::options> fluidFvOptions(fluidRegions.size());
24 
25 // Populate fluid field pointer lists
27 {
28  Info<< "*** Reading fluid mesh thermophysical properties for region "
29  << fluidRegions[i].name() << nl << endl;
30 
31  Info<< " Adding to thermoFluid\n" << endl;
33 
34  Info<< " Adding to rhoFluid\n" << endl;
35  rhoFluid.set
36  (
37  i,
38  new volScalarField
39  (
40  IOobject
41  (
42  "rho",
43  runTime.timeName(),
44  fluidRegions[i],
45  IOobject::NO_READ,
46  IOobject::AUTO_WRITE
47  ),
48  thermoFluid[i].rho()
49  )
50  );
51 
52  Info<< " Adding to UFluid\n" << endl;
53  UFluid.set
54  (
55  i,
56  new volVectorField
57  (
58  IOobject
59  (
60  "U",
61  runTime.timeName(),
62  fluidRegions[i],
63  IOobject::MUST_READ,
64  IOobject::AUTO_WRITE
65  ),
66  fluidRegions[i]
67  )
68  );
69 
70  Info<< " Adding to phiFluid\n" << endl;
71  phiFluid.set
72  (
73  i,
75  (
76  IOobject
77  (
78  "phi",
79  runTime.timeName(),
80  fluidRegions[i],
81  IOobject::READ_IF_PRESENT,
82  IOobject::AUTO_WRITE
83  ),
85  & fluidRegions[i].Sf()
86  )
87  );
88 
89  Info<< " Adding to gFluid\n" << endl;
90  gFluid.set
91  (
92  i,
94  (
95  IOobject
96  (
97  "g",
98  runTime.constant(),
99  fluidRegions[i],
100  IOobject::MUST_READ,
101  IOobject::NO_WRITE
102  )
103  )
104  );
105 
106  Info<< " Adding to hRefFluid\n" << endl;
107  hRefFluid.set
108  (
109  i,
111  (
112  IOobject
113  (
114  "hRef",
115  runTime.constant(),
116  fluidRegions[i],
117  IOobject::READ_IF_PRESENT,
118  IOobject::NO_WRITE
119  ),
120  dimensionedScalar("hRef", dimLength, 0)
121  )
122  );
123 
124  dimensionedScalar ghRef(- mag(gFluid[i])*hRefFluid[i]);
125 
126  Info<< " Adding to ghFluid\n" << endl;
127  ghFluid.set
128  (
129  i,
130  new volScalarField
131  (
132  "gh",
133  (gFluid[i] & fluidRegions[i].C()) - ghRef
134  )
135  );
136 
137  Info<< " Adding to ghfFluid\n" << endl;
138  ghfFluid.set
139  (
140  i,
142  (
143  "ghf",
144  (gFluid[i] & 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],
228  (
229  "dpdt",
230  thermoFluid[i].p().dimensions()/dimTime,
231  0
232  )
233  )
234  );
235 
236  Info<< " Adding to fieldsFluid\n" << endl;
237  fieldsFluid.set
238  (
239  i,
240  new multivariateSurfaceInterpolationScheme<scalar>::fieldTable
241  );
242  forAll(thermoFluid[i].composition().Y(), j)
243  {
244  fieldsFluid[i].add(thermoFluid[i].composition().Y()[j]);
245  }
246  fieldsFluid[i].add(thermoFluid[i].he());
247 
248  Info<< " Adding to QdotFluid\n" << endl;
249  QdotFluid.set
250  (
251  i,
252  new volScalarField
253  (
254  IOobject
255  (
256  "Qdot",
257  runTime.timeName(),
258  fluidRegions[i],
259  IOobject::READ_IF_PRESENT,
260  IOobject::AUTO_WRITE
261  ),
262  fluidRegions[i],
264  )
265  );
266 
267  Info<< " Adding MRF\n" << endl;
268  MRFfluid.set
269  (
270  i,
271  new IOMRFZoneList(fluidRegions[i])
272  );
273 
274  Info<< " Adding fvOptions\n" << endl;
275  fluidFvOptions.set
276  (
277  i,
278  new fv::options(fluidRegions[i])
279  );
280 
281  turbulenceFluid[i].validate();
282 }
volScalarField & he
Definition: YEEqn.H:51
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:107
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())
PtrList< radiation::radiationModel > radiation(fluidRegions.size())
PtrList< volScalarField > QdotFluid(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())
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
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)
volScalarField & C
const volScalarField & T
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimEnergy
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 > &)
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
PtrList< surfaceScalarField > phiFluid(fluidRegions.size())
PtrList< CombustionModel< rhoReactionThermo > > reactionFluid(fluidRegions.size())
PtrList< volScalarField > rhoFluid(fluidRegions.size())