phaseSystem.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2015-2017 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "phaseSystem.H"
27 #include "surfaceTensionModel.H"
28 #include "aspectRatioModel.H"
29 #include "surfaceInterpolate.H"
30 #include "fvcDdt.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(phaseSystem, 0);
37 }
38 
39 const Foam::word Foam::phaseSystem::propertiesName("phaseProperties");
40 
41 
42 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43 
45 (
46  const phaseModelList& phaseModels
47 ) const
48 {
49  tmp<surfaceScalarField> tmpPhi
50  (
52  (
53  "phi",
54  fvc::interpolate(phaseModels[0])*phaseModels[0].phi()
55  )
56  );
57 
58  for (label phasei=1; phasei<phaseModels.size(); phasei++)
59  {
60  tmpPhi.ref() +=
61  fvc::interpolate(phaseModels[phasei])*phaseModels[phasei].phi();
62  }
63 
64  return tmpPhi;
65 }
66 
67 
69 (
70  const dictTable& modelDicts
71 )
72 {
73  forAllConstIter(dictTable, modelDicts, iter)
74  {
75  const phasePairKey& key = iter.key();
76 
77  // pair already exists
78  if (phasePairs_.found(key))
79  {}
80 
81  // new ordered pair
82  else if (key.ordered())
83  {
84  phasePairs_.insert
85  (
86  key,
87  autoPtr<phasePair>
88  (
89  new orderedPhasePair
90  (
91  phaseModels_[key.first()],
92  phaseModels_[key.second()]
93  )
94  )
95  );
96  }
97 
98  // new unordered pair
99  else
100  {
101  phasePairs_.insert
102  (
103  key,
104  autoPtr<phasePair>
105  (
106  new phasePair
107  (
108  phaseModels_[key.first()],
109  phaseModels_[key.second()]
110  )
111  )
112  );
113  }
114  }
115 }
116 
117 
118 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
119 
121 (
122  const fvMesh& mesh
123 )
124 :
125  IOdictionary
126  (
127  IOobject
128  (
129  "phaseProperties",
130  mesh.time().constant(),
131  mesh,
132  IOobject::MUST_READ_IF_MODIFIED,
133  IOobject::NO_WRITE
134  )
135  ),
136 
137  mesh_(mesh),
138 
139  phaseModels_(lookup("phases"), phaseModel::iNew(*this)),
140 
141  phi_(calcPhi(phaseModels_)),
142 
143  dpdt_
144  (
145  IOobject
146  (
147  "dpdt",
148  mesh.time().timeName(),
149  mesh
150  ),
151  mesh,
153  ),
154 
155  MRF_(mesh_)
156 {
157  phi_.writeOpt() = IOobject::AUTO_WRITE;
158 
159  // Blending methods
160  forAllConstIter(dictionary, subDict("blending"), iter)
161  {
162  blendingMethods_.insert
163  (
164  iter().dict().dictName(),
166  (
167  iter().dict(),
168  phaseModels_.toc()
169  )
170  );
171  }
172 
173  // Sub-models
174  generatePairsAndSubModels("surfaceTension", surfaceTensionModels_);
175  generatePairsAndSubModels("aspectRatio", aspectRatioModels_);
176 
177  correctKinematics();
178 }
179 
180 
181 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
182 
184 {}
185 
186 
187 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
188 
190 {
191  tmp<volScalarField> tmpRho
192  (
193  phaseModels_[0]*phaseModels_[0].rho()
194  );
195 
196  for (label phasei=1; phasei<phaseModels_.size(); phasei++)
197  {
198  tmpRho.ref() += phaseModels_[phasei]*phaseModels_[phasei].rho();
199  }
200 
201  return tmpRho;
202 }
203 
204 
206 {
207  tmp<volVectorField> tmpU
208  (
209  phaseModels_[0]*phaseModels_[0].U()
210  );
211 
212  for (label phasei=1; phasei<phaseModels_.size(); phasei++)
213  {
214  tmpU.ref() += phaseModels_[phasei]*phaseModels_[phasei].U();
215  }
216 
217  return tmpU;
218 }
219 
220 
222 Foam::phaseSystem::E(const phasePairKey& key) const
223 {
224  if (aspectRatioModels_.found(key))
225  {
226  return aspectRatioModels_[key]->E();
227  }
228  else
229  {
230  return tmp<volScalarField>
231  (
232  new volScalarField
233  (
234  IOobject
235  (
236  aspectRatioModel::typeName + ":E",
237  this->mesh_.time().timeName(),
238  this->mesh_,
241  false
242  ),
243  this->mesh_,
244  dimensionedScalar("zero", dimless, 1)
245  )
246  );
247  }
248 }
249 
250 
252 Foam::phaseSystem::sigma(const phasePairKey& key) const
253 {
254  if (surfaceTensionModels_.found(key))
255  {
256  return surfaceTensionModels_[key]->sigma();
257  }
258  else
259  {
260  return tmp<volScalarField>
261  (
262  new volScalarField
263  (
264  IOobject
265  (
266  surfaceTensionModel::typeName + ":sigma",
267  this->mesh_.time().timeName(),
268  this->mesh_,
271  false
272  ),
273  this->mesh_,
275  )
276  );
277  }
278 }
279 
280 
282 {}
283 
284 
286 {
287  forAll(phaseModels_, phasei)
288  {
289  phaseModels_[phasei].correct();
290  }
291 }
292 
293 
295 {
296  bool updateDpdt = false;
297 
298  forAll(phaseModels_, phasei)
299  {
300  phaseModels_[phasei].correctKinematics();
301 
302  updateDpdt = updateDpdt || phaseModels_[phasei].thermo().dpdt();
303  }
304 
305  // Update the pressure time-derivative if required
306  if (updateDpdt)
307  {
308  dpdt_ = fvc::ddt(phaseModels_.begin()().thermo().p());
309  }
310 }
311 
312 
314 {
315  forAll(phaseModels_, phasei)
316  {
317  phaseModels_[phasei].correctThermo();
318  }
319 }
320 
321 
323 {
324  forAll(phaseModels_, phasei)
325  {
326  phaseModels_[phasei].correctTurbulence();
327  }
328 }
329 
330 
332 {
333  forAll(phaseModels_, phasei)
334  {
335  phaseModels_[phasei].correctEnergyTransport();
336  }
337 }
338 
339 
341 {
342  if (regIOobject::read())
343  {
344  bool readOK = true;
345 
346  forAll(phaseModels_, phasei)
347  {
348  readOK &= phaseModels_[phasei].read();
349  }
350 
351  // models ...
352 
353  return readOK;
354  }
355  else
356  {
357  return false;
358  }
359 }
360 
361 
362 // ************************************************************************* //
surfaceScalarField & phi
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:264
dictionary dict
void generatePairs(const dictTable &modelDicts)
Generate pairs.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
tmp< volScalarField > E(const phasePairKey &key) const
Return the aspect-ratio.
U
Definition: pEqn.H:83
virtual void correctTurbulence()
Correct the turbulence.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
label phasei
Definition: pEqn.H:27
virtual bool read()
Read object.
virtual void solve()
Solve for the phase fractions.
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Calculate the first temporal derivative.
psiReactionThermo & thermo
Definition: createFields.H:31
virtual void correctKinematics()
Correct the kinematics.
stressControl lookup("compactNormalStress") >> compactNormalStress
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
A class for handling words, derived from string.
Definition: word.H:59
tmp< volVectorField > U() const
Return the mixture velocity.
virtual void correctThermo()
Correct the thermodynamics.
word timeName
Definition: getTimeIndex.H:3
const word dictName("particleTrackDict")
const dimensionSet dimPressure
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
static autoPtr< blendingMethod > New(const dictionary &dict, const wordList &phaseNames)
virtual void correct()
Correct the fluid properties other than the thermo and turbulence.
tmp< volScalarField > rho() const
Return the mixture density.
defineTypeNameAndDebug(combustionModel, 0)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
tmp< surfaceScalarField > calcPhi(const phaseModelList &phaseModels) const
Calculate and return the mixture flux.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
static const dimensionSet dimSigma
Surface tension coefficient dimensions.
virtual ~phaseSystem()
Destructor.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool read()
Read base phaseProperties dictionary.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
phasePairKey()
Construct null.
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Namespace for OpenFOAM.