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  IOobject io
52  (
53  IOobject::groupName("phi", this->name()),
54  U.mesh().time().timeName(),
55  U.mesh(),
56  IOobject::MUST_READ,
57  IOobject::AUTO_WRITE,
58  IOobject::REGISTER
59  );
60 
61  if (io.typeHeaderOk<surfaceScalarField>(true))
62  {
63  Info<< "Reading face flux field " << io.name() << endl;
64 
65  return tmp<surfaceScalarField>::New(io, U.mesh());
66  }
67  else
68  {
69  Info<< "Calculating face flux field " << io.name() << endl;
70 
71  io.readOpt(IOobject::NO_READ);
72 
74  (
75  U.boundaryField().size(),
77  );
78 
79  forAll(U.boundaryField(), i)
80  {
81  if
82  (
83  isA<fixedValueFvPatchVectorField>(U.boundaryField()[i])
84  || isA<slipFvPatchVectorField>(U.boundaryField()[i])
85  || isA<partialSlipFvPatchVectorField>(U.boundaryField()[i])
86  )
87  {
88  patchTypes[i] = fixedValueFvPatchScalarField::typeName;
89  }
90  }
91 
93  (
94  io,
95  fvc::flux(U),
97  );
98  }
99 }
100 
101 
102 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
103 
104 template<class BasePhaseModel>
106 (
107  const phaseSystem& fluid,
108  const word& phaseName,
109  const label index
110 )
111 :
112  BasePhaseModel(fluid, phaseName, index),
113  U_
114  (
115  IOobject
116  (
117  IOobject::groupName("U", this->name()),
118  fluid.mesh().time().timeName(),
119  fluid.mesh(),
120  IOobject::MUST_READ,
121  IOobject::AUTO_WRITE,
122  IOobject::REGISTER
123  ),
124  fluid.mesh()
125  ),
126  phi_(phi(U_)),
127  alphaPhi_
128  (
129  IOobject
130  (
131  IOobject::groupName("alphaPhi", this->name()),
132  fluid.mesh().time().timeName(),
133  fluid.mesh()
134  ),
135  fluid.mesh(),
136  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0))
137  ),
138  alphaRhoPhi_
139  (
140  IOobject
141  (
142  IOobject::groupName("alphaRhoPhi", this->name()),
143  fluid.mesh().time().timeName(),
144  fluid.mesh()
145  ),
146  fluid.mesh(),
147  dimensionedScalar(dimensionSet(1, 0, -1, 0, 0))
148  ),
149  DUDt_(nullptr),
150  DUDtf_(nullptr),
151  divU_(nullptr),
152  turbulence_
153  (
155  (
156  *this,
157  this->thermo().rho(),
158  U_,
159  alphaRhoPhi_,
160  phi_,
161  *this
162  )
163  ),
164  continuityErrorFlow_
165  (
166  IOobject
167  (
168  IOobject::groupName("continuityErrorFlow", this->name()),
169  fluid.mesh().time().timeName(),
170  fluid.mesh()
171  ),
172  fluid.mesh(),
174  ),
175  continuityErrorSources_
176  (
177  IOobject
178  (
179  IOobject::groupName("continuityErrorSources", this->name()),
180  fluid.mesh().time().timeName(),
181  fluid.mesh()
182  ),
183  fluid.mesh(),
185  ),
186  K_(nullptr)
187 {
189 
191 }
192 
193 
194 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
195 
196 template<class BasePhaseModel>
198 {
200 
201  this->fluid().MRF().correctBoundaryVelocity(U_);
202 
203  volScalarField& rho = this->thermoRef().rho();
204 
205  continuityErrorFlow_ = fvc::ddt(*this, rho) + fvc::div(alphaRhoPhi_);
206 
207  continuityErrorSources_ = - (this->fluid().fvOptions()(*this, rho)&rho);
208 }
209 
210 
211 template<class BasePhaseModel>
213 {
214  BasePhaseModel::correctKinematics();
215 
216  if (DUDt_.valid())
217  {
218  DUDt_.clear();
219  DUDt();
220  }
221 
222  if (DUDtf_.valid())
223  {
224  DUDtf_.clear();
225  DUDtf();
226  }
227 
228  if (K_.valid())
229  {
230  K_.ref() = 0.5*magSqr(this->U());
231  }
232 }
233 
234 
235 template<class BasePhaseModel>
237 {
239 
240  turbulence_->correct();
241 }
242 
243 
244 template<class BasePhaseModel>
246 {
247  BasePhaseModel::correctEnergyTransport();
248 
249  turbulence_->correctEnergyTransport();
250 }
251 
252 
253 template<class BasePhaseModel>
255 {
256  return false;
257 }
258 
259 
260 template<class BasePhaseModel>
263 {
264  const volScalarField& alpha = *this;
265  const tmp<volScalarField> trho = this->thermo().rho();
266  const volScalarField& rho = trho();
267 
268  return
269  (
270  fvm::ddt(alpha, rho, U_)
271  + fvm::div(alphaRhoPhi_, U_)
272  + fvm::SuSp(- this->continuityError(), U_)
273  + this->fluid().MRF().DDt(alpha*rho, U_)
274  + turbulence_->divDevRhoReff(U_)
275  );
276 }
277 
278 
279 template<class BasePhaseModel>
282 {
283  // As the "normal" U-eqn but without the ddt terms
284 
285  const volScalarField& alpha = *this;
286 
287  return
288  (
289  fvm::div(alphaRhoPhi_, U_)
290  - fvm::Sp(fvc::div(alphaRhoPhi_), U_)
291  + fvm::SuSp(- this->continuityErrorSources(), U_)
292  + this->fluid().MRF().DDt(alpha*this->thermo().rho(), U_)
293  + turbulence_->divDevRhoReff(U_)
294  );
295 }
296 
297 
298 template<class BasePhaseModel>
301 {
302  return U_;
303 }
304 
305 
306 template<class BasePhaseModel>
309 {
310  return U_;
311 }
312 
313 
314 template<class BasePhaseModel>
317 {
318  return phi_;
319 }
320 
321 
322 template<class BasePhaseModel>
325 {
326  return phi_;
327 }
328 
329 
330 template<class BasePhaseModel>
333 {
334  return alphaPhi_;
335 }
336 
337 
338 template<class BasePhaseModel>
341 {
342  return alphaPhi_;
343 }
344 
345 
346 template<class BasePhaseModel>
349 {
350  return alphaRhoPhi_;
351 }
352 
353 
354 template<class BasePhaseModel>
357 {
358  return alphaRhoPhi_;
359 }
360 
361 
362 template<class BasePhaseModel>
365 {
366  if (!DUDt_.valid())
367  {
368  DUDt_ = fvc::ddt(U_) + fvc::div(phi_, U_) - fvc::div(phi_)*U_;
369  }
371  return tmp<volVectorField>(DUDt_());
372 }
373 
374 
375 template<class BasePhaseModel>
378 {
379  if (!DUDtf_.valid())
380  {
381  DUDtf_ = byDt(phi_ - phi_.oldTime());
382  }
384  return tmp<surfaceScalarField>(DUDtf_());
385 }
386 
387 
388 template<class BasePhaseModel>
391 {
392  return continuityErrorFlow_ + continuityErrorSources_;
393 }
394 
395 
396 template<class BasePhaseModel>
399 {
400  return continuityErrorFlow_;
401 }
402 
403 
404 template<class BasePhaseModel>
407 {
408  return continuityErrorSources_;
409 }
410 
411 
412 template<class BasePhaseModel>
415 {
416  if (!K_.valid())
417  {
419  (
420  IOobject::groupName("K", this->name()),
421  0.5*magSqr(this->U())
422  );
423  }
425  return tmp<volScalarField>(K_());
426 }
427 
428 
429 template<class BasePhaseModel>
432 {
433  return divU_.valid() ? tmp<volScalarField>(divU_()) : tmp<volScalarField>();
434 }
435 
436 
437 template<class BasePhaseModel>
439 {
440  divU_ = divU;
441 }
442 
443 
444 template<class BasePhaseModel>
447 {
448  return turbulence_->mut();
449 }
450 
451 
452 template<class BasePhaseModel>
455 {
456  return turbulence_->muEff();
457 }
458 
459 
460 template<class BasePhaseModel>
463 {
464  return turbulence_->nut();
465 }
466 
467 
468 template<class BasePhaseModel>
471 {
472  return turbulence_->nuEff();
473 }
474 
475 
476 template<class BasePhaseModel>
479 {
480  return turbulence_->kappaEff();
481 }
482 
483 
484 template<class BasePhaseModel>
486 Foam::MovingPhaseModel<BasePhaseModel>::kappaEff(const label patchi) const
487 {
488  return turbulence_->kappaEff(patchi);
489 }
490 
491 
492 template<class BasePhaseModel>
495 {
496  return turbulence_->alphaEff();
497 }
498 
499 
500 template<class BasePhaseModel>
502 Foam::MovingPhaseModel<BasePhaseModel>::alphaEff(const label patchi) const
503 {
504  return turbulence_->alphaEff(patchi);
505 }
506 
507 
508 template<class BasePhaseModel>
511 {
512  return turbulence_->k();
513 }
514 
515 
516 template<class BasePhaseModel>
519 {
520  return turbulence_->pPrime();
521 }
522 
523 
524 // ************************************************************************* //
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.
wordList patchTypes(nPatches)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< volScalarField > trho
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:76
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:72
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:402
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:436
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:143
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/v2406/OpenFOAM-v2406/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
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
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:180
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.