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-2021 OpenFOAM Foundation
9  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "MovingPhaseModel.H"
30 #include "phaseSystem.H"
31 #include "phaseCompressibleTurbulenceModel.H"
33 #include "slipFvPatchFields.H"
35 
36 #include "fvmDdt.H"
37 #include "fvmDiv.H"
38 #include "fvmSup.H"
39 #include "fvcDdt.H"
40 #include "fvcDiv.H"
41 #include "fvcFlux.H"
42 #include "surfaceInterpolate.H"
43 #include "fvMatrix.H"
44 
45 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
46 
47 template<class BasePhaseModel>
50 {
51  word phiName(IOobject::groupName("phi", this->name()));
52 
53  IOobject phiHeader
54  (
55  phiName,
56  U.mesh().time().timeName(),
57  U.mesh(),
58  IOobject::NO_READ
59  );
60 
61  if (phiHeader.typeHeaderOk<surfaceScalarField>(true))
62  {
63  Info<< "Reading face flux field " << phiName << endl;
64 
65  return tmp<surfaceScalarField>
66  (
68  (
69  IOobject
70  (
71  phiName,
72  U.mesh().time().timeName(),
73  U.mesh(),
74  IOobject::MUST_READ,
75  IOobject::AUTO_WRITE
76  ),
77  U.mesh()
78  )
79  );
80  }
81  else
82  {
83  Info<< "Calculating face flux field " << phiName << endl;
84 
85  wordList phiTypes
86  (
87  U.boundaryField().size(),
89  );
90 
91  forAll(U.boundaryField(), i)
92  {
93  if
94  (
95  isA<fixedValueFvPatchVectorField>(U.boundaryField()[i])
96  || isA<slipFvPatchVectorField>(U.boundaryField()[i])
97  || isA<partialSlipFvPatchVectorField>(U.boundaryField()[i])
98  )
99  {
100  phiTypes[i] = fixedValueFvPatchScalarField::typeName;
101  }
102  }
103 
104  return tmp<surfaceScalarField>
105  (
107  (
108  IOobject
109  (
110  phiName,
111  U.mesh().time().timeName(),
112  U.mesh(),
113  IOobject::NO_READ,
114  IOobject::AUTO_WRITE
115  ),
116  fvc::flux(U),
117  phiTypes
118  )
119  );
120  }
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
125 
126 template<class BasePhaseModel>
128 (
129  const phaseSystem& fluid,
130  const word& phaseName,
131  const label index
132 )
133 :
134  BasePhaseModel(fluid, phaseName, index),
135  U_
136  (
137  IOobject
138  (
139  IOobject::groupName("U", this->name()),
140  fluid.mesh().time().timeName(),
141  fluid.mesh(),
142  IOobject::MUST_READ,
143  IOobject::AUTO_WRITE
144  ),
145  fluid.mesh()
146  ),
147  phi_(phi(U_)),
148  alphaPhi_
149  (
150  IOobject
151  (
152  IOobject::groupName("alphaPhi", this->name()),
153  fluid.mesh().time().timeName(),
154  fluid.mesh()
155  ),
156  fluid.mesh(),
157  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0))
158  ),
159  alphaRhoPhi_
160  (
161  IOobject
162  (
163  IOobject::groupName("alphaRhoPhi", this->name()),
164  fluid.mesh().time().timeName(),
165  fluid.mesh()
166  ),
167  fluid.mesh(),
168  dimensionedScalar(dimensionSet(1, 0, -1, 0, 0))
169  ),
170  DUDt_(nullptr),
171  DUDtf_(nullptr),
172  divU_(nullptr),
173  turbulence_
174  (
176  (
177  *this,
178  this->thermo().rho(),
179  U_,
180  alphaRhoPhi_,
181  phi_,
182  *this
183  )
184  ),
185  continuityErrorFlow_
186  (
187  IOobject
188  (
189  IOobject::groupName("continuityErrorFlow", this->name()),
190  fluid.mesh().time().timeName(),
191  fluid.mesh()
192  ),
193  fluid.mesh(),
195  ),
196  continuityErrorSources_
197  (
198  IOobject
199  (
200  IOobject::groupName("continuityErrorSources", this->name()),
201  fluid.mesh().time().timeName(),
202  fluid.mesh()
203  ),
204  fluid.mesh(),
206  ),
207  K_(nullptr)
208 {
210 
212 }
213 
214 
215 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
216 
217 template<class BasePhaseModel>
219 {
221 
222  this->fluid().MRF().correctBoundaryVelocity(U_);
223 
224  volScalarField& rho = this->thermoRef().rho();
225 
226  continuityErrorFlow_ = fvc::ddt(*this, rho) + fvc::div(alphaRhoPhi_);
227 
228  continuityErrorSources_ = - (this->fluid().fvOptions()(*this, rho)&rho);
229 }
230 
231 
232 template<class BasePhaseModel>
234 {
235  BasePhaseModel::correctKinematics();
236 
237  if (DUDt_.valid())
238  {
239  DUDt_.clear();
240  DUDt();
241  }
242 
243  if (DUDtf_.valid())
244  {
245  DUDtf_.clear();
246  DUDtf();
247  }
248 
249  if (K_.valid())
250  {
251  K_.ref() = 0.5*magSqr(this->U());
252  }
253 }
254 
255 
256 template<class BasePhaseModel>
258 {
260 
261  turbulence_->correct();
262 }
263 
264 
265 template<class BasePhaseModel>
267 {
268  BasePhaseModel::correctEnergyTransport();
269 
270  turbulence_->correctEnergyTransport();
271 }
272 
273 
274 template<class BasePhaseModel>
276 {
277  return false;
278 }
279 
280 
281 template<class BasePhaseModel>
284 {
285  const volScalarField& alpha = *this;
286  const volScalarField& rho = this->thermo().rho();
287 
288  return
289  (
290  fvm::ddt(alpha, rho, U_)
291  + fvm::div(alphaRhoPhi_, U_)
292  + fvm::SuSp(- this->continuityError(), U_)
293  + this->fluid().MRF().DDt(alpha*rho, U_)
294  + turbulence_->divDevRhoReff(U_)
295  );
296 }
297 
298 
299 template<class BasePhaseModel>
302 {
303  // As the "normal" U-eqn but without the ddt terms
304 
305  const volScalarField& alpha = *this;
306  const volScalarField& rho = this->thermo().rho();
307 
308  return
309  (
310  fvm::div(alphaRhoPhi_, U_)
311  - fvm::Sp(fvc::div(alphaRhoPhi_), U_)
312  + fvm::SuSp(- this->continuityErrorSources(), U_)
313  + this->fluid().MRF().DDt(alpha*rho, U_)
314  + turbulence_->divDevRhoReff(U_)
315  );
316 }
317 
318 
319 template<class BasePhaseModel>
322 {
323  return U_;
324 }
325 
326 
327 template<class BasePhaseModel>
330 {
331  return U_;
332 }
333 
334 
335 template<class BasePhaseModel>
338 {
339  return phi_;
340 }
341 
342 
343 template<class BasePhaseModel>
346 {
347  return phi_;
348 }
349 
350 
351 template<class BasePhaseModel>
354 {
355  return alphaPhi_;
356 }
357 
358 
359 template<class BasePhaseModel>
362 {
363  return alphaPhi_;
364 }
365 
366 
367 template<class BasePhaseModel>
370 {
371  return alphaRhoPhi_;
372 }
373 
374 
375 template<class BasePhaseModel>
378 {
379  return alphaRhoPhi_;
380 }
381 
382 
383 template<class BasePhaseModel>
386 {
387  if (!DUDt_.valid())
388  {
389  DUDt_ = fvc::ddt(U_) + fvc::div(phi_, U_) - fvc::div(phi_)*U_;
390  }
392  return tmp<volVectorField>(DUDt_());
393 }
394 
395 
396 template<class BasePhaseModel>
399 {
400  if (!DUDtf_.valid())
401  {
402  DUDtf_ = byDt(phi_ - phi_.oldTime());
403  }
405  return tmp<surfaceScalarField>(DUDtf_());
406 }
407 
408 
409 template<class BasePhaseModel>
412 {
413  return continuityErrorFlow_ + continuityErrorSources_;
414 }
415 
416 
417 template<class BasePhaseModel>
420 {
421  return continuityErrorFlow_;
422 }
423 
424 
425 template<class BasePhaseModel>
428 {
429  return continuityErrorSources_;
430 }
431 
432 
433 template<class BasePhaseModel>
436 {
437  if (!K_.valid())
438  {
440  (
441  IOobject::groupName("K", this->name()),
442  0.5*magSqr(this->U())
443  );
444  }
446  return tmp<volScalarField>(K_());
447 }
448 
449 
450 template<class BasePhaseModel>
453 {
454  return divU_.valid() ? tmp<volScalarField>(divU_()) : tmp<volScalarField>();
455 }
456 
457 
458 template<class BasePhaseModel>
460 {
461  divU_ = divU;
462 }
463 
464 
465 template<class BasePhaseModel>
468 {
469  return turbulence_->mut();
470 }
471 
472 
473 template<class BasePhaseModel>
476 {
477  return turbulence_->muEff();
478 }
479 
480 
481 template<class BasePhaseModel>
484 {
485  return turbulence_->nut();
486 }
487 
488 
489 template<class BasePhaseModel>
492 {
493  return turbulence_->nuEff();
494 }
495 
496 
497 template<class BasePhaseModel>
500 {
501  return turbulence_->kappaEff();
502 }
503 
504 
505 template<class BasePhaseModel>
507 Foam::MovingPhaseModel<BasePhaseModel>::kappaEff(const label patchi) const
508 {
509  return turbulence_->kappaEff(patchi);
510 }
511 
512 
513 template<class BasePhaseModel>
516 {
517  return turbulence_->alphaEff();
518 }
519 
520 
521 template<class BasePhaseModel>
523 Foam::MovingPhaseModel<BasePhaseModel>::alphaEff(const label patchi) const
524 {
525  return turbulence_->alphaEff(patchi);
526 }
527 
528 
529 template<class BasePhaseModel>
532 {
533  return turbulence_->k();
534 }
535 
536 
537 template<class BasePhaseModel>
540 {
541  return turbulence_->pPrime();
542 }
543 
544 
545 // ************************************************************************* //
virtual tmp< volVectorField > U() const
Access const reference to U.
twoPhaseSystem & fluid
fv::options & fvOptions() const
Access the fvOptions.
Definition: phaseSystemI.H:136
tmp< GeometricField< Type, fvPatchField, volMesh > > DDt(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &psi)
Definition: fvcDDt.C:40
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
virtual tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
virtual tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
MovingPhaseModel(const multiphaseInterSystem &fluid, const word &phaseName)
virtual void correctKinematics()
Correct the kinematics.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
fluid correctTurbulence()
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
virtual tmp< volScalarField > continuityErrorFlow() const
Return the continuity error due to the flow field.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:40
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
Calculate the first temporal derivative.
psiReactionThermo & thermo
Definition: createFields.H:28
const word calculatedType
A calculated patch field type.
word timeName
Definition: getTimeIndex.H:3
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
writeOption writeOpt() const noexcept
Get the write option.
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZoneList.C:394
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual volVectorField & URef()
Access the velocity.
IOMRFZoneList & MRF
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:69
A class for handling words, derived from Foam::string.
Definition: word.H:63
Calculate the face-flux of the given field.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
Calculate the matrix for the first temporal derivative.
const IOMRFZoneList & MRF() const
Return MRF zones.
Definition: phaseSystemI.H:130
tmp< volScalarField > byDt(const volScalarField &vf)
Definition: phaseSystem.C:435
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
virtual tmp< volScalarField > kappaEff() const
Return the effective thermal conductivity.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
virtual tmp< volScalarField > rho() const
Density [kg/m^3] - uses current value of pressure.
Definition: psiThermo.C:139
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
Calculate the divergence of the given field.
virtual tmp< volScalarField > divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2312/OpenFOAM-v2312/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
virtual tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
virtual tmp< volScalarField > alphaEff() const
Return the effective thermal diffusivity.
const dimensionSet dimDensity
Calculate the matrix for the divergence of the given field and flux.
List< word > wordList
List of word.
Definition: fileName.H:59
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
U
Definition: pEqn.H:72
Automatically write from objectRegistry::writeObject()
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
zeroField divU
Definition: alphaSuSp.H:3
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual void correctTurbulence()
Correct the turbulence.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
virtual tmp< volScalarField > continuityErrorSources() const
Return the continuity error due to any sources.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Calculate the finiteVolume matrix for implicit and explicit sources.
virtual bool stationary() const
Return whether the phase is stationary.
virtual tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux.
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access the volumetric flux of the phase.