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 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() +=
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  // do nothing ...
81  }
82 
83  // new ordered pair
84  else if (key.ordered())
85  {
86  phasePairs_.insert
87  (
88  key,
89  autoPtr<phasePair>
90  (
91  new orderedPhasePair
92  (
93  phaseModels_[key.first()],
94  phaseModels_[key.second()]
95  )
96  )
97  );
98  }
99 
100  // new unordered pair
101  else
102  {
103  phasePairs_.insert
104  (
105  key,
106  autoPtr<phasePair>
107  (
108  new phasePair
109  (
110  phaseModels_[key.first()],
111  phaseModels_[key.second()]
112  )
113  )
114  );
115  }
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
123 (
124  const fvMesh& mesh
125 )
126 :
127  IOdictionary
128  (
129  IOobject
130  (
131  "phaseProperties",
132  mesh.time().constant(),
133  mesh,
134  IOobject::MUST_READ_IF_MODIFIED,
135  IOobject::NO_WRITE
136  )
137  ),
138 
139  mesh_(mesh),
140 
141  phaseModels_(lookup("phases"), phaseModel::iNew(*this)),
142 
143  phi_(calcPhi(phaseModels_)),
144 
145  dpdt_
146  (
147  IOobject
148  (
149  "dpdt",
150  mesh.time().timeName(),
151  mesh
152  ),
153  mesh,
155  ),
156 
157  MRF_(mesh_),
158  fvOptions_(mesh_)
159 {
160  phi_.writeOpt() = IOobject::AUTO_WRITE;
161 
162  // Blending methods
163  forAllConstIter(dictionary, subDict("blending"), iter)
164  {
165  blendingMethods_.insert
166  (
167  iter().dict().dictName(),
169  (
170  iter().dict(),
171  phaseModels_.toc()
172  )
173  );
174  }
175 
176  // Sub-models
177  generatePairsAndSubModels("surfaceTension", surfaceTensionModels_);
178  generatePairsAndSubModels("aspectRatio", aspectRatioModels_);
179 
180  correctKinematics();
181 }
182 
183 
184 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
185 
187 {}
188 
189 
190 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
191 
193 {
194  tmp<volScalarField> tmpRho
195  (
196  phaseModels_[0]*phaseModels_[0].rho()
197  );
198 
199  for (label phasei=1; phasei<phaseModels_.size(); phasei++)
200  {
201  tmpRho() += phaseModels_[phasei]*phaseModels_[phasei].rho();
202  }
203 
204  return tmpRho;
205 }
206 
207 
209 {
210  tmp<volVectorField> tmpU
211  (
212  phaseModels_[0]*phaseModels_[0].U()
213  );
214 
215  for (label phasei=1; phasei<phaseModels_.size(); phasei++)
216  {
217  tmpU() += phaseModels_[phasei]*phaseModels_[phasei].U();
218  }
219 
220  return tmpU;
221 }
222 
223 
225 Foam::phaseSystem::E(const phasePairKey& key) const
226 {
227  if (aspectRatioModels_.found(key))
228  {
229  return aspectRatioModels_[key]->E();
230  }
231  else
232  {
233  return tmp<volScalarField>
234  (
235  new volScalarField
236  (
237  IOobject
238  (
239  aspectRatioModel::typeName + ":E",
240  this->mesh_.time().timeName(),
241  this->mesh_,
244  false
245  ),
246  this->mesh_,
247  dimensionedScalar("zero", dimless, 1)
248  )
249  );
250  }
251 }
252 
253 
255 Foam::phaseSystem::sigma(const phasePairKey& key) const
256 {
257  if (surfaceTensionModels_.found(key))
258  {
259  return surfaceTensionModels_[key]->sigma();
260  }
261  else
262  {
263  return tmp<volScalarField>
264  (
265  new volScalarField
266  (
267  IOobject
268  (
269  surfaceTensionModel::typeName + ":sigma",
270  this->mesh_.time().timeName(),
271  this->mesh_,
274  false
275  ),
276  this->mesh_,
278  )
279  );
280  }
281 }
282 
283 
285 {}
286 
287 
289 {
290  forAll(phaseModels_, phasei)
291  {
292  phaseModels_[phasei].correct();
293  }
294 }
295 
296 
298 {
299  bool updateDpdt = false;
300 
301  forAll(phaseModels_, phasei)
302  {
303  phaseModels_[phasei].correctKinematics();
304 
305  updateDpdt = updateDpdt || phaseModels_[phasei].thermo().dpdt();
306  }
307 
308  // Update the pressure time-derivative if required
309  if (updateDpdt)
310  {
311  dpdt_ = fvc::ddt(phaseModels_.begin()().thermo().p());
312  }
313 }
314 
315 
317 {
318  forAll(phaseModels_, phasei)
319  {
320  phaseModels_[phasei].correctThermo();
321  }
322 }
323 
324 
326 {
327  forAll(phaseModels_, phasei)
328  {
329  phaseModels_[phasei].correctTurbulence();
330  }
331 }
332 
333 
335 {
336  forAll(phaseModels_, phasei)
337  {
338  phaseModels_[phasei].correctEnergyTransport();
339  }
340 }
341 
342 
344 {
345  if (regIOobject::read())
346  {
347  bool readOK = true;
348 
349  forAll(phaseModels_, phasei)
350  {
351  readOK &= phaseModels_[phasei].read();
352  }
353 
354  // models ...
355 
356  return readOK;
357  }
358  else
359  {
360  return false;
361  }
362 }
363 
364 
365 // ************************************************************************* //
Foam::dimPressure
const dimensionSet dimPressure
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::surfaceTensionModel::dimSigma
static const dimensionSet dimSigma
Coefficient dimensions.
Definition: surfaceTensionModel.H:88
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::phaseSystem::correctKinematics
virtual void correctKinematics()
Correct the kinematics.
Foam::blendingMethod::New
static autoPtr< blendingMethod > New(const dictionary &dict, const wordList &phaseNames)
Foam::phaseSystem::calcPhi
tmp< surfaceScalarField > calcPhi(const phaseModelList &phaseModels) const
Calculate and return the mixture flux.
Foam::fvc::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Definition: surfaceInterpolate.C:110
Foam::phaseSystem::phaseSystem
phaseSystem(const fvMesh &mesh)
Construct from fvMesh.
Foam::regIOobject::read
virtual bool read()
Read object.
Definition: regIOobjectRead.C:171
Foam::phaseSystem::read
virtual bool read()
Read base phaseProperties dictionary.
phaseSystem.H
Foam::phaseSystem::propertiesName
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:267
Foam::phaseSystem::rho
tmp< volScalarField > rho() const
Return the mixture density.
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
U
U
Definition: pEqn.H:46
constant
Constant dispersed-phase particle diameter model.
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::phaseSystem::correct
virtual void correct()
Correct the fluid properties other than the thermo and turbulence.
dictName
const word dictName("particleTrackDict")
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::phaseSystem::solve
virtual void solve()
Solve for the phase fractions.
Foam::label
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
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::phaseSystem::correctEnergyTransport
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
Foam::phaseSystem::~phaseSystem
virtual ~phaseSystem()
Destructor.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::phaseSystem::correctTurbulence
virtual void correctTurbulence()
Correct the turbulence.
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::phaseSystem::E
tmp< volScalarField > E(const phasePairKey &key) const
Return the aspect-ratio.
surfaceInterpolate.H
Surface Interpolation.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
rho
rho
Definition: pEqn.H:3
Foam::phaseSystem::sigma
tmp< volScalarField > sigma(const phasePairKey &key) const
Return the surface tension coefficient.
phasei
label phasei
Definition: pEqn.H:37
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
Foam::phaseSystem::U
tmp< volVectorField > U() const
Return the mixture velocity.
Foam::phaseSystem::correctThermo
virtual void correctThermo()
Correct the thermodynamics.
surfaceTensionModel.H
timeName
word timeName
Definition: getTimeIndex.H:3
fvcDdt.H
Calculate the first temporal derivative.
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::phaseSystem::generatePairs
void generatePairs(const dictTable &modelDicts)
Generate pairs.
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress
Foam::phaseSystem::phi
const surfaceScalarField & phi() const
Constant access the mixture flux.
Definition: phaseSystemI.H:55