MovingPhaseModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2018 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 "MovingPhaseModel.H"
27 #include "phaseSystem.H"
28 #include "phaseCompressibleTurbulenceModel.H"
30 #include "slipFvPatchFields.H"
32 
33 #include "fvmDdt.H"
34 #include "fvmDiv.H"
35 #include "fvmSup.H"
36 #include "fvcDdt.H"
37 #include "fvcDiv.H"
38 #include "fvcFlux.H"
39 #include "surfaceInterpolate.H"
40 #include "fvMatrix.H"
41 
42 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
43 
44 template<class BasePhaseModel>
47 {
48  word phiName(IOobject::groupName("phi", this->name()));
49 
50  IOobject phiHeader
51  (
52  phiName,
53  U.mesh().time().timeName(),
54  U.mesh(),
55  IOobject::NO_READ
56  );
57 
58  if (phiHeader.typeHeaderOk<surfaceScalarField>(true))
59  {
60  Info<< "Reading face flux field " << phiName << endl;
61 
62  return tmp<surfaceScalarField>
63  (
65  (
66  IOobject
67  (
68  phiName,
69  U.mesh().time().timeName(),
70  U.mesh(),
71  IOobject::MUST_READ,
72  IOobject::AUTO_WRITE
73  ),
74  U.mesh()
75  )
76  );
77  }
78  else
79  {
80  Info<< "Calculating face flux field " << phiName << endl;
81 
82  wordList phiTypes
83  (
84  U.boundaryField().size(),
85  calculatedFvPatchScalarField::typeName
86  );
87 
88  forAll(U.boundaryField(), i)
89  {
90  if
91  (
92  isA<fixedValueFvPatchVectorField>(U.boundaryField()[i])
93  || isA<slipFvPatchVectorField>(U.boundaryField()[i])
94  || isA<partialSlipFvPatchVectorField>(U.boundaryField()[i])
95  )
96  {
97  phiTypes[i] = fixedValueFvPatchScalarField::typeName;
98  }
99  }
100 
101  return tmp<surfaceScalarField>
102  (
104  (
105  IOobject
106  (
107  phiName,
108  U.mesh().time().timeName(),
109  U.mesh(),
110  IOobject::NO_READ,
111  IOobject::AUTO_WRITE
112  ),
113  fvc::flux(U),
114  phiTypes
115  )
116  );
117  }
118 }
119 
120 
121 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
122 
123 template<class BasePhaseModel>
125 (
126  const phaseSystem& fluid,
127  const word& phaseName,
128  const label index
129 )
130 :
131  BasePhaseModel(fluid, phaseName, index),
132  U_
133  (
134  IOobject
135  (
136  IOobject::groupName("U", this->name()),
137  fluid.mesh().time().timeName(),
138  fluid.mesh(),
139  IOobject::MUST_READ,
140  IOobject::AUTO_WRITE
141  ),
142  fluid.mesh()
143  ),
144  phi_(phi(U_)),
145  alphaPhi_
146  (
147  IOobject
148  (
149  IOobject::groupName("alphaPhi", this->name()),
150  fluid.mesh().time().timeName(),
151  fluid.mesh()
152  ),
153  fluid.mesh(),
154  dimensionedScalar("0", dimensionSet(0, 3, -1, 0, 0), 0)
155  ),
156  alphaRhoPhi_
157  (
158  IOobject
159  (
160  IOobject::groupName("alphaRhoPhi", this->name()),
161  fluid.mesh().time().timeName(),
162  fluid.mesh()
163  ),
164  fluid.mesh(),
165  dimensionedScalar("0", dimensionSet(1, 0, -1, 0, 0), 0)
166  ),
167  DUDt_(nullptr),
168  DUDtf_(nullptr),
169  divU_(nullptr),
170  turbulence_
171  (
173  (
174  *this,
175  this->thermo().rho(),
176  U_,
177  alphaRhoPhi_,
178  phi_,
179  *this
180  )
181  ),
182  continuityErrorFlow_
183  (
184  IOobject
185  (
186  IOobject::groupName("continuityErrorFlow", this->name()),
187  fluid.mesh().time().timeName(),
188  fluid.mesh()
189  ),
190  fluid.mesh(),
192  ),
193  continuityErrorSources_
194  (
195  IOobject
196  (
197  IOobject::groupName("continuityErrorSources", this->name()),
198  fluid.mesh().time().timeName(),
199  fluid.mesh()
200  ),
201  fluid.mesh(),
203  ),
204  K_(nullptr)
205 {
206  phi_.writeOpt() = IOobject::AUTO_WRITE;
207 
208  correctKinematics();
209 }
210 
211 
212 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
213 
214 template<class BasePhaseModel>
216 {}
217 
218 
219 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
220 
221 template<class BasePhaseModel>
223 {
225 
226  this->fluid().MRF().correctBoundaryVelocity(U_);
227 
228  volScalarField& rho = this->thermoRef().rho();
229 
230  continuityErrorFlow_ = fvc::ddt(*this, rho) + fvc::div(alphaRhoPhi_);
231 
232  continuityErrorSources_ = - (this->fluid().fvOptions()(*this, rho)&rho);
233 }
234 
235 
236 template<class BasePhaseModel>
238 {
239  BasePhaseModel::correctKinematics();
240 
241  if (DUDt_.valid())
242  {
243  DUDt_.clear();
244  DUDt();
245  }
246 
247  if (DUDtf_.valid())
248  {
249  DUDtf_.clear();
250  DUDtf();
251  }
252 
253  if (K_.valid())
254  {
255  K_.clear();
256  K();
257  }
258 }
259 
260 
261 template<class BasePhaseModel>
263 {
264  BasePhaseModel::correctTurbulence();
265 
266  turbulence_->correct();
267 }
268 
269 
270 template<class BasePhaseModel>
272 {
273  BasePhaseModel::correctEnergyTransport();
274 
275  turbulence_->correctEnergyTransport();
276 }
277 
278 
279 template<class BasePhaseModel>
281 {
282  return false;
283 }
284 
285 
286 template<class BasePhaseModel>
289 {
290  const volScalarField& alpha = *this;
291  const volScalarField& rho = this->thermo().rho();
292 
293  return
294  (
295  fvm::ddt(alpha, rho, U_)
296  + fvm::div(alphaRhoPhi_, U_)
297  + fvm::SuSp(- this->continuityError(), U_)
298  + this->fluid().MRF().DDt(alpha*rho, U_)
299  + turbulence_->divDevRhoReff(U_)
300  );
301 }
302 
303 
304 template<class BasePhaseModel>
307 {
308  // As the "normal" U-eqn but without the ddt terms
309 
310  const volScalarField& alpha = *this;
311  const volScalarField& rho = this->thermo().rho();
312 
313  return
314  (
315  fvm::div(alphaRhoPhi_, U_)
316  - fvm::Sp(fvc::div(alphaRhoPhi_), U_)
317  + fvm::SuSp(- this->continuityErrorSources(), U_)
318  + this->fluid().MRF().DDt(alpha*rho, U_)
319  + turbulence_->divDevRhoReff(U_)
320  );
321 }
322 
323 
324 template<class BasePhaseModel>
327 {
328  return U_;
329 }
330 
331 
332 template<class BasePhaseModel>
335 {
336  return U_;
337 }
338 
339 
340 template<class BasePhaseModel>
343 {
344  return phi_;
345 }
346 
347 
348 template<class BasePhaseModel>
351 {
352  return phi_;
353 }
354 
355 
356 template<class BasePhaseModel>
359 {
360  return alphaPhi_;
361 }
362 
363 
364 template<class BasePhaseModel>
367 {
368  return alphaPhi_;
369 }
370 
371 
372 template<class BasePhaseModel>
375 {
376  return alphaRhoPhi_;
377 }
378 
379 
380 template<class BasePhaseModel>
383 {
384  return alphaRhoPhi_;
385 }
386 
387 
388 template<class BasePhaseModel>
391 {
392  if (!DUDt_.valid())
393  {
394  DUDt_ = fvc::ddt(U_) + fvc::div(phi_, U_) - fvc::div(phi_)*U_;
395  }
396 
397  return tmp<volVectorField>(DUDt_());
398 }
399 
400 
401 template<class BasePhaseModel>
404 {
405  if (!DUDtf_.valid())
406  {
407  DUDtf_ = byDt(phi_ - phi_.oldTime());
408  }
409 
410  return tmp<surfaceScalarField>(DUDtf_());
411 }
412 
413 
414 template<class BasePhaseModel>
417 {
418  return continuityErrorFlow_ + continuityErrorSources_;
419 }
420 
421 
422 template<class BasePhaseModel>
425 {
426  return continuityErrorFlow_;
427 }
428 
429 
430 template<class BasePhaseModel>
433 {
434  return continuityErrorSources_;
435 }
436 
437 
438 template<class BasePhaseModel>
441 {
442  if (!K_.valid())
443  {
444  K_ =
445  new volScalarField
446  (
447  IOobject::groupName("K", this->name()),
448  0.5*magSqr(this->U())
449  );
450  }
451 
452  return tmp<volScalarField>(K_());
453 }
454 
455 
456 template<class BasePhaseModel>
459 {
460  return divU_.valid() ? tmp<volScalarField>(divU_()) : tmp<volScalarField>();
461 }
462 
463 
464 template<class BasePhaseModel>
465 void Foam::MovingPhaseModel<BasePhaseModel>::divU(tmp<volScalarField> divU)
466 {
467  divU_ = divU;
468 }
469 
470 
471 template<class BasePhaseModel>
474 {
475  return turbulence_->mut();
476 }
477 
478 
479 template<class BasePhaseModel>
482 {
483  return turbulence_->muEff();
484 }
485 
486 
487 template<class BasePhaseModel>
490 {
491  return turbulence_->nut();
492 }
493 
494 
495 template<class BasePhaseModel>
498 {
499  return turbulence_->nuEff();
500 }
501 
502 
503 template<class BasePhaseModel>
506 {
507  return turbulence_->kappaEff();
508 }
509 
510 
511 template<class BasePhaseModel>
514 {
515  return turbulence_->kappaEff(patchi);
516 }
517 
518 
519 template<class BasePhaseModel>
522 {
523  return turbulence_->alphaEff();
524 }
525 
526 
527 template<class BasePhaseModel>
530 {
531  return turbulence_->alphaEff(patchi);
532 }
533 
534 
535 template<class BasePhaseModel>
538 {
539  return turbulence_->k();
540 }
541 
542 
543 template<class BasePhaseModel>
546 {
547  return turbulence_->pPrime();
548 }
549 
550 
551 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual tmp< volScalarField > divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
tmp< GeometricField< Type, fvPatchField, volMesh > > DDt(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &psi)
Definition: fvcDDt.C:45
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
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
surfaceScalarField & phi
multiphaseSystem & fluid
Definition: createFields.H:11
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
IOMRFZoneList & MRF
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual tmp< volScalarField > alphaEff() const
Return the effective thermal diffusivity.
virtual bool stationary() const
Return whether the phase is stationary.
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
virtual tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
virtual tmp< surfaceScalarField > alphaPhi() const
Return the volumetric flux of the phase.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
rhoReactionThermo & thermo
Definition: createFields.H:28
virtual void correctTurbulence()
Correct the turbulence.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
CGAL::Exact_predicates_exact_constructions_kernel K
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Calculate the first temporal derivative.
dynamicFvMesh & mesh
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
virtual ~MovingPhaseModel()
Destructor.
Calculate the face-flux of the given field.
static word groupName(Name name, const word &group)
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calculate the matrix for the first temporal derivative.
virtual tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
word timeName
Definition: getTimeIndex.H:3
virtual tmp< volScalarField > continuityErrorSources() const
Return the continuity error due to any sources.
virtual tmp< volScalarField > continuityErrorFlow() const
Return the continuity error due to the flow field.
virtual volVectorField & URef()
Access the velocity.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
Calculate the divergence of the given field.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
virtual tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
virtual tmp< volScalarField > kappaEff() const
Return the effective thermal conductivity.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
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/ubuntu/OpenFOAM-6/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:72
virtual tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
virtual tmp< volVectorField > U() const
Return the velocity.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
tmp< volScalarField > byDt(const volScalarField &vf)
const dimensionSet dimDensity
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
Calculate the matrix for the divergence of the given field and flux.
List< word > wordList
A List of words.
Definition: fileName.H:54
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
virtual void correctKinematics()
Correct the kinematics.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux.
zeroField divU
Definition: alphaSuSp.H:3
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
Calculate the matrix for implicit and explicit sources.