MovingPhaseModel.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2018 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "MovingPhaseModel.H"
29 #include "phaseSystem.H"
30 #include "phaseCompressibleTurbulenceModel.H"
32 #include "slipFvPatchFields.H"
34 
35 #include "fvmDdt.H"
36 #include "fvmDiv.H"
37 #include "fvmSup.H"
38 #include "fvcDdt.H"
39 #include "fvcDiv.H"
40 #include "fvcFlux.H"
41 #include "surfaceInterpolate.H"
42 #include "fvMatrix.H"
43 
44 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
45 
46 template<class BasePhaseModel>
49 {
50  word phiName(IOobject::groupName("phi", this->name()));
51 
52  IOobject phiHeader
53  (
54  phiName,
55  U.mesh().time().timeName(),
56  U.mesh(),
57  IOobject::NO_READ
58  );
59 
60  if (phiHeader.typeHeaderOk<surfaceScalarField>(true))
61  {
62  Info<< "Reading face flux field " << phiName << endl;
63 
64  return tmp<surfaceScalarField>
65  (
67  (
68  IOobject
69  (
70  phiName,
71  U.mesh().time().timeName(),
72  U.mesh(),
73  IOobject::MUST_READ,
74  IOobject::AUTO_WRITE
75  ),
76  U.mesh()
77  )
78  );
79  }
80  else
81  {
82  Info<< "Calculating face flux field " << phiName << endl;
83 
84  wordList phiTypes
85  (
86  U.boundaryField().size(),
87  calculatedFvPatchScalarField::typeName
88  );
89 
90  forAll(U.boundaryField(), i)
91  {
92  if
93  (
94  isA<fixedValueFvPatchVectorField>(U.boundaryField()[i])
95  || isA<slipFvPatchVectorField>(U.boundaryField()[i])
96  || isA<partialSlipFvPatchVectorField>(U.boundaryField()[i])
97  )
98  {
99  phiTypes[i] = fixedValueFvPatchScalarField::typeName;
100  }
101  }
102 
103  return tmp<surfaceScalarField>
104  (
106  (
107  IOobject
108  (
109  phiName,
110  U.mesh().time().timeName(),
111  U.mesh(),
112  IOobject::NO_READ,
113  IOobject::AUTO_WRITE
114  ),
115  fvc::flux(U),
116  phiTypes
117  )
118  );
119  }
120 }
121 
122 
123 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
124 
125 template<class BasePhaseModel>
127 (
128  const phaseSystem& fluid,
129  const word& phaseName,
130  const label index
131 )
132 :
133  BasePhaseModel(fluid, phaseName, index),
134  U_
135  (
136  IOobject
137  (
138  IOobject::groupName("U", this->name()),
139  fluid.mesh().time().timeName(),
140  fluid.mesh(),
141  IOobject::MUST_READ,
142  IOobject::AUTO_WRITE
143  ),
144  fluid.mesh()
145  ),
146  phi_(phi(U_)),
147  alphaPhi_
148  (
149  IOobject
150  (
151  IOobject::groupName("alphaPhi", this->name()),
152  fluid.mesh().time().timeName(),
153  fluid.mesh()
154  ),
155  fluid.mesh(),
156  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0))
157  ),
158  alphaRhoPhi_
159  (
160  IOobject
161  (
162  IOobject::groupName("alphaRhoPhi", this->name()),
163  fluid.mesh().time().timeName(),
164  fluid.mesh()
165  ),
166  fluid.mesh(),
167  dimensionedScalar(dimensionSet(1, 0, -1, 0, 0))
168  ),
169  DUDt_(nullptr),
170  DUDtf_(nullptr),
171  divU_(nullptr),
172  turbulence_
173  (
175  (
176  *this,
177  this->thermo().rho(),
178  U_,
179  alphaRhoPhi_,
180  phi_,
181  *this
182  )
183  ),
184  continuityErrorFlow_
185  (
186  IOobject
187  (
188  IOobject::groupName("continuityErrorFlow", this->name()),
189  fluid.mesh().time().timeName(),
190  fluid.mesh()
191  ),
192  fluid.mesh(),
194  ),
195  continuityErrorSources_
196  (
197  IOobject
198  (
199  IOobject::groupName("continuityErrorSources", this->name()),
200  fluid.mesh().time().timeName(),
201  fluid.mesh()
202  ),
203  fluid.mesh(),
205  ),
206  K_(nullptr)
207 {
208  phi_.writeOpt() = IOobject::AUTO_WRITE;
209 
210  correctKinematics();
211 }
212 
213 
214 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
215 
216 template<class BasePhaseModel>
218 {}
219 
220 
221 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
222 
223 template<class BasePhaseModel>
225 {
227 
228  this->fluid().MRF().correctBoundaryVelocity(U_);
229 
230  volScalarField& rho = this->thermoRef().rho();
231 
232  continuityErrorFlow_ = fvc::ddt(*this, rho) + fvc::div(alphaRhoPhi_);
233 
234  continuityErrorSources_ = - (this->fluid().fvOptions()(*this, rho)&rho);
235 }
236 
237 
238 template<class BasePhaseModel>
240 {
241  BasePhaseModel::correctKinematics();
242 
243  if (DUDt_.valid())
244  {
245  DUDt_.clear();
246  DUDt();
247  }
248 
249  if (DUDtf_.valid())
250  {
251  DUDtf_.clear();
252  DUDtf();
253  }
254 
255  if (K_.valid())
256  {
257  K_.clear();
258  K();
259  }
260 }
261 
262 
263 template<class BasePhaseModel>
265 {
267 
268  turbulence_->correct();
269 }
270 
271 
272 template<class BasePhaseModel>
274 {
275  BasePhaseModel::correctEnergyTransport();
276 
277  turbulence_->correctEnergyTransport();
278 }
279 
280 
281 template<class BasePhaseModel>
283 {
284  return false;
285 }
286 
287 
288 template<class BasePhaseModel>
291 {
292  const volScalarField& alpha = *this;
293  const volScalarField& rho = this->thermo().rho();
294 
295  return
296  (
297  fvm::ddt(alpha, rho, U_)
298  + fvm::div(alphaRhoPhi_, U_)
299  + fvm::SuSp(- this->continuityError(), U_)
300  + this->fluid().MRF().DDt(alpha*rho, U_)
301  + turbulence_->divDevRhoReff(U_)
302  );
303 }
304 
305 
306 template<class BasePhaseModel>
309 {
310  // As the "normal" U-eqn but without the ddt terms
311 
312  const volScalarField& alpha = *this;
313  const volScalarField& rho = this->thermo().rho();
314 
315  return
316  (
317  fvm::div(alphaRhoPhi_, U_)
318  - fvm::Sp(fvc::div(alphaRhoPhi_), U_)
319  + fvm::SuSp(- this->continuityErrorSources(), U_)
320  + this->fluid().MRF().DDt(alpha*rho, U_)
321  + turbulence_->divDevRhoReff(U_)
322  );
323 }
324 
325 
326 template<class BasePhaseModel>
329 {
330  return U_;
331 }
332 
333 
334 template<class BasePhaseModel>
337 {
338  return U_;
339 }
340 
341 
342 template<class BasePhaseModel>
345 {
346  return phi_;
347 }
348 
349 
350 template<class BasePhaseModel>
353 {
354  return phi_;
355 }
356 
357 
358 template<class BasePhaseModel>
361 {
362  return alphaPhi_;
363 }
364 
365 
366 template<class BasePhaseModel>
369 {
370  return alphaPhi_;
371 }
372 
373 
374 template<class BasePhaseModel>
377 {
378  return alphaRhoPhi_;
379 }
380 
381 
382 template<class BasePhaseModel>
385 {
386  return alphaRhoPhi_;
387 }
388 
389 
390 template<class BasePhaseModel>
393 {
394  if (!DUDt_.valid())
395  {
396  DUDt_ = fvc::ddt(U_) + fvc::div(phi_, U_) - fvc::div(phi_)*U_;
397  }
398 
399  return tmp<volVectorField>(DUDt_());
400 }
401 
402 
403 template<class BasePhaseModel>
406 {
407  if (!DUDtf_.valid())
408  {
409  DUDtf_ = byDt(phi_ - phi_.oldTime());
410  }
411 
412  return tmp<surfaceScalarField>(DUDtf_());
413 }
414 
415 
416 template<class BasePhaseModel>
419 {
420  return continuityErrorFlow_ + continuityErrorSources_;
421 }
422 
423 
424 template<class BasePhaseModel>
427 {
428  return continuityErrorFlow_;
429 }
430 
431 
432 template<class BasePhaseModel>
435 {
436  return continuityErrorSources_;
437 }
438 
439 
440 template<class BasePhaseModel>
443 {
444  if (!K_.valid())
445  {
447  (
448  IOobject::groupName("K", this->name()),
449  0.5*magSqr(this->U())
450  );
451  }
452 
453  return tmp<volScalarField>(K_());
454 }
455 
456 
457 template<class BasePhaseModel>
460 {
461  return divU_.valid() ? tmp<volScalarField>(divU_()) : tmp<volScalarField>();
462 }
463 
464 
465 template<class BasePhaseModel>
467 {
468  divU_ = divU;
469 }
470 
471 
472 template<class BasePhaseModel>
475 {
476  return turbulence_->mut();
477 }
478 
479 
480 template<class BasePhaseModel>
483 {
484  return turbulence_->muEff();
485 }
486 
487 
488 template<class BasePhaseModel>
491 {
492  return turbulence_->nut();
493 }
494 
495 
496 template<class BasePhaseModel>
499 {
500  return turbulence_->nuEff();
501 }
502 
503 
504 template<class BasePhaseModel>
507 {
508  return turbulence_->kappaEff();
509 }
510 
511 
512 template<class BasePhaseModel>
515 {
516  return turbulence_->kappaEff(patchi);
517 }
518 
519 
520 template<class BasePhaseModel>
523 {
524  return turbulence_->alphaEff();
525 }
526 
527 
528 template<class BasePhaseModel>
531 {
532  return turbulence_->alphaEff(patchi);
533 }
534 
535 
536 template<class BasePhaseModel>
539 {
540  return turbulence_->k();
541 }
542 
543 
544 template<class BasePhaseModel>
547 {
548  return turbulence_->pPrime();
549 }
550 
551 
552 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
partialSlipFvPatchFields.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::MovingPhaseModel::kappaEff
virtual tmp< volScalarField > kappaEff() const
Return the effective thermal conductivity.
Definition: MovingPhaseModel.C:506
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
slipFvPatchFields.H
Foam::MovingPhaseModel::continuityErrorSources
virtual tmp< volScalarField > continuityErrorSources() const
Return the continuity error due to any sources.
Definition: MovingPhaseModel.C:434
Foam::dimDensity
const dimensionSet dimDensity
Foam::MovingPhaseModel::alphaPhi
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access the volumetric flux of the phase.
Definition: MovingPhaseModel.C:97
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::MovingPhaseModel::stationary
virtual bool stationary() const
Return whether the phase is stationary.
Definition: MovingPhaseModel.C:282
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Sp
zeroField Sp
Definition: alphaSuSp.H:2
Foam::MovingPhaseModel::continuityError
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
Definition: MovingPhaseModel.C:418
fvcDiv.H
Calculate the divergence of the given field.
MovingPhaseModel.H
fvMatrix.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
Foam::MovingPhaseModel::mut
virtual tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
Definition: MovingPhaseModel.C:474
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
rho
rho
Definition: readInitialConditions.H:88
Foam::MovingPhaseModel::~MovingPhaseModel
virtual ~MovingPhaseModel()=default
Destructor.
Definition: MovingPhaseModel.C:217
fvmDiv.H
Calculate the matrix for the divergence of the given field and flux.
Foam::MovingPhaseModel::nuEff
virtual tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
Definition: MovingPhaseModel.C:498
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::fam::SuSp
tmp< faMatrix< Type > > SuSp(const areaScalarField &sp, const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famSup.C:151
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::MovingPhaseModel::phiRef
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
Definition: MovingPhaseModel.C:352
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:59
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:22
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::MovingPhaseModel::alphaPhiRef
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
Definition: MovingPhaseModel.C:368
Foam::MovingPhaseModel::phi
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux.
Definition: MovingPhaseModel.C:81
Foam::MovingPhaseModel::correctTurbulence
virtual void correctTurbulence()
Correct the turbulence.
Definition: MovingPhaseModel.C:264
Foam::MovingPhaseModel::k
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
Definition: MovingPhaseModel.C:538
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
correct
fvOptions correct(rho)
Foam::MovingPhaseModel::U
virtual tmp< volVectorField > U() const
Access const reference to U.
Definition: MovingPhaseModel.C:113
Foam::MovingPhaseModel::alphaEff
virtual tmp< volScalarField > alphaEff() const
Return the effective thermal diffusivity.
Definition: MovingPhaseModel.C:522
divU
zeroField divU
Definition: alphaSuSp.H:3
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
correctTurbulence
fluid correctTurbulence()
Foam::MovingPhaseModel::UfEqn
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
Definition: MovingPhaseModel.C:308
Foam::MovingPhaseModel::MovingPhaseModel
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName)
Definition: MovingPhaseModel.C:48
Foam::MovingPhaseModel::UEqn
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
Definition: MovingPhaseModel.C:290
Foam::MovingPhaseModel::DUDtf
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
Definition: MovingPhaseModel.C:405
Foam::MovingPhaseModel::correctKinematics
virtual void correctKinematics()
Correct the kinematics.
Definition: MovingPhaseModel.C:239
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::MovingPhaseModel::URef
virtual volVectorField & URef()
Access the velocity.
Definition: MovingPhaseModel.C:336
Foam::phaseSystem::mesh
const fvMesh & mesh() const
Return mesh.
Definition: phaseSystem.C:979
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
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
U
U
Definition: pEqn.H:72
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::MovingPhaseModel::K
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
Definition: MovingPhaseModel.C:442
Foam::MovingPhaseModel::correct
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
Definition: MovingPhaseModel.C:73
Foam::MovingPhaseModel::alphaRhoPhi
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
Definition: MovingPhaseModel.C:376
Foam::MovingPhaseModel::pPrime
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
Definition: MovingPhaseModel.C:546
Foam::fac::ddt
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:47
fvcFlux.H
Calculate the face-flux of the given field.
fixedValueFvPatchFields.H
Foam::MovingPhaseModel::divU
virtual tmp< volScalarField > divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
Definition: MovingPhaseModel.C:459
Foam::MovingPhaseModel::alphaRhoPhiRef
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
Definition: MovingPhaseModel.C:384
Foam::MovingPhaseModel::nut
virtual tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
Definition: MovingPhaseModel.C:490
Foam::MovingPhaseModel::continuityErrorFlow
virtual tmp< volScalarField > continuityErrorFlow() const
Return the continuity error due to the flow field.
Definition: MovingPhaseModel.C:426
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
fvcDdt.H
Calculate the first temporal derivative.
Foam::MovingPhaseModel::muEff
virtual tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
Definition: MovingPhaseModel.C:482
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:63
Foam::MovingPhaseModel::DUDt
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
Definition: MovingPhaseModel.C:392
Foam::byDt
tmp< volScalarField > byDt(const volScalarField &vf)
Definition: phaseSystem.C:438
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::MovingPhaseModel::correctEnergyTransport
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
Definition: MovingPhaseModel.C:273
fvmDdt.H
Calculate the matrix for the first temporal derivative.
Foam::fvc::DDt
tmp< GeometricField< Type, fvPatchField, volMesh > > DDt(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &psi)
Definition: fvcDDt.C:47