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-2020 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"
29 #include "slipFvPatchFields.H"
31 
32 #include "fvmDdt.H"
33 #include "fvmDiv.H"
34 #include "fvmSup.H"
35 #include "fvcDdt.H"
36 #include "fvcDiv.H"
37 #include "fvcFlux.H"
38 
39 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
40 
41 template<class BasePhaseModel>
44 {
45  word phiName(IOobject::groupName("phi", this->name()));
46 
47  IOobject phiHeader
48  (
49  phiName,
50  U.mesh().time().timeName(),
51  U.mesh(),
52  IOobject::NO_READ
53  );
54 
55  if (phiHeader.typeHeaderOk<surfaceScalarField>(true))
56  {
57  Info<< "Reading face flux field " << phiName << endl;
58 
59  return tmp<surfaceScalarField>
60  (
62  (
63  IOobject
64  (
65  phiName,
66  U.mesh().time().timeName(),
67  U.mesh(),
68  IOobject::MUST_READ,
69  IOobject::AUTO_WRITE
70  ),
71  U.mesh()
72  )
73  );
74  }
75  else
76  {
77  Info<< "Calculating face flux field " << phiName << endl;
78 
79  wordList phiTypes
80  (
81  U.boundaryField().size(),
82  calculatedFvPatchScalarField::typeName
83  );
84 
85  forAll(U.boundaryField(), patchi)
86  {
87  if (!U.boundaryField()[patchi].assignable())
88  {
89  phiTypes[patchi] = fixedValueFvPatchScalarField::typeName;
90  }
91  }
92 
93  return tmp<surfaceScalarField>
94  (
96  (
97  IOobject
98  (
99  phiName,
100  U.mesh().time().timeName(),
101  U.mesh(),
102  IOobject::NO_READ,
103  IOobject::AUTO_WRITE
104  ),
105  fvc::flux(U),
106  phiTypes
107  )
108  );
109  }
110 }
111 
112 
113 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
114 
115 template<class BasePhaseModel>
117 (
118  const phaseSystem& fluid,
119  const word& phaseName,
120  const bool referencePhase,
121  const label index
122 )
123 :
124  BasePhaseModel(fluid, phaseName, referencePhase, index),
125  U_
126  (
127  IOobject
128  (
129  IOobject::groupName("U", this->name()),
130  fluid.mesh().time().timeName(),
131  fluid.mesh(),
132  IOobject::MUST_READ,
133  IOobject::AUTO_WRITE
134  ),
135  fluid.mesh()
136  ),
137  phi_(phi(U_)),
138  alphaPhi_
139  (
140  IOobject
141  (
142  IOobject::groupName("alphaPhi", this->name()),
143  fluid.mesh().time().timeName(),
144  fluid.mesh()
145  ),
146  fluid.mesh(),
147  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), 0)
148  ),
149  alphaRhoPhi_
150  (
151  IOobject
152  (
153  IOobject::groupName("alphaRhoPhi", this->name()),
154  fluid.mesh().time().timeName(),
155  fluid.mesh()
156  ),
157  fluid.mesh(),
158  dimensionedScalar(dimensionSet(1, 0, -1, 0, 0), 0)
159  ),
160  DUDt_(nullptr),
161  DUDtf_(nullptr),
162  divU_(nullptr),
163  turbulence_
164  (
166  (
167  *this,
168  this->thermo().rho(),
169  U_,
170  alphaRhoPhi_,
171  phi_,
172  *this
173  )
174  ),
175  thermophysicalTransport_
176  (
177  PhaseThermophysicalTransportModel
178  <
180  typename BasePhaseModel::thermoModel
181  >::New(turbulence_, this->thermo_)
182  ),
183  continuityError_
184  (
185  IOobject
186  (
187  IOobject::groupName("continuityError", this->name()),
188  fluid.mesh().time().timeName(),
189  fluid.mesh()
190  ),
191  fluid.mesh(),
193  ),
194  K_(nullptr)
195 {
196  phi_.writeOpt() = IOobject::AUTO_WRITE;
197 
198  correctKinematics();
199 }
200 
201 
202 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
203 
204 template<class BasePhaseModel>
206 {}
207 
208 
209 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
210 
211 template<class BasePhaseModel>
213 (
214  const volScalarField& source
215 )
216 {
217  volScalarField& rho = this->thermoRef().rho();
218 
219  continuityError_ = fvc::ddt(*this, rho) + fvc::div(alphaRhoPhi_) - source;
220 }
221 
222 
223 template<class BasePhaseModel>
225 {
227  this->fluid().MRF().correctBoundaryVelocity(U_);
228 }
229 
230 
231 template<class BasePhaseModel>
233 {
234  BasePhaseModel::correctKinematics();
235 
236  if (DUDt_.valid())
237  {
238  DUDt_.clear();
239  DUDt();
240  }
241 
242  if (DUDtf_.valid())
243  {
244  DUDtf_.clear();
245  DUDtf();
246  }
247 
248  if (K_.valid())
249  {
250  K_.clear();
251  K();
252  }
253 }
254 
255 
256 template<class BasePhaseModel>
258 {
259  BasePhaseModel::correctTurbulence();
260 
261  turbulence_->correct();
262 }
263 
264 
265 template<class BasePhaseModel>
267 {
268  BasePhaseModel::correctEnergyTransport();
269  thermophysicalTransport_->correct();
270 }
271 
272 
273 template<class BasePhaseModel>
275 {
276  return false;
277 }
278 
279 
280 template<class BasePhaseModel>
283 {
284  const volScalarField& alpha = *this;
285  const volScalarField& rho = this->thermo().rho();
286 
287  return
288  (
289  fvm::ddt(alpha, rho, U_)
290  + fvm::div(alphaRhoPhi_, U_)
291  + fvm::SuSp(-this->continuityError(), U_)
292  + this->fluid().MRF().DDt(alpha*rho, U_)
293  + turbulence_->divDevTau(U_)
294  );
295 }
296 
297 
298 template<class BasePhaseModel>
301 {
302  // As the "normal" U-eqn but without the ddt terms
303 
304  const volScalarField& alpha = *this;
305  const volScalarField& rho = this->thermo().rho();
306 
307  return
308  (
309  fvm::div(alphaRhoPhi_, U_)
310  + fvm::SuSp(fvc::ddt(*this, rho) - this->continuityError(), U_)
311  + this->fluid().MRF().DDt(alpha*rho, U_)
312  + turbulence_->divDevTau(U_)
313  );
314 }
315 
316 
317 template<class BasePhaseModel>
320 {
321  return U_;
322 }
323 
324 
325 template<class BasePhaseModel>
328 {
329  return U_;
330 }
331 
332 
333 template<class BasePhaseModel>
336 {
337  return phi_;
338 }
339 
340 
341 template<class BasePhaseModel>
344 {
345  return phi_;
346 }
347 
348 
349 template<class BasePhaseModel>
352 {
353  return alphaPhi_;
354 }
355 
356 
357 template<class BasePhaseModel>
360 {
361  return alphaPhi_;
362 }
363 
364 
365 template<class BasePhaseModel>
368 {
369  return alphaRhoPhi_;
370 }
371 
372 
373 template<class BasePhaseModel>
376 {
377  return alphaRhoPhi_;
378 }
379 
380 
381 template<class BasePhaseModel>
384 {
385  if (!DUDt_.valid())
386  {
387  DUDt_ = fvc::ddt(U_) + fvc::div(phi_, U_) - fvc::div(phi_)*U_;
388  }
389 
390  return tmp<volVectorField>(DUDt_());
391 }
392 
393 
394 template<class BasePhaseModel>
397 {
398  if (!DUDtf_.valid())
399  {
400  DUDtf_ = byDt(phi_ - phi_.oldTime());
401  }
402 
403  return tmp<surfaceScalarField>(DUDtf_());
404 }
405 
406 
407 template<class BasePhaseModel>
410 {
411  return continuityError_;
412 }
413 
414 
415 template<class BasePhaseModel>
418 {
419  if (!K_.valid())
420  {
422  (
423  IOobject::groupName("K", this->name()),
424  0.5*magSqr(this->U())
425  );
426  }
427 
428  return tmp<volScalarField>(K_());
429 }
430 
431 
432 template<class BasePhaseModel>
435 {
436  return divU_.valid() ? tmp<volScalarField>(divU_()) : tmp<volScalarField>();
437 }
438 
439 
440 template<class BasePhaseModel>
441 void Foam::MovingPhaseModel<BasePhaseModel>::divU(tmp<volScalarField> divU)
442 {
443  divU_ = divU;
444 }
445 
446 
447 template<class BasePhaseModel>
450 {
451  return thermophysicalTransport_->kappaEff(patchi);
452 }
453 
454 
455 template<class BasePhaseModel>
458 {
459  return turbulence_->k();
460 }
461 
462 
463 template<class BasePhaseModel>
466 {
467  return turbulence_->pPrime();
468 }
469 
470 
471 template<class BasePhaseModel>
474 {
475  return thermophysicalTransport_->divq(he);
476 }
477 
478 
479 template<class BasePhaseModel>
482 {
483  return thermophysicalTransport_->divj(Yi);
484 }
485 
486 
487 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
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
PhaseCompressibleMomentumTransportModel< phaseModel > phaseCompressibleMomentumTransportModel
Typedef for phaseCompressibleMomentumTransportModel.
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
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual bool stationary() const
Return whether the phase is stationary.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
virtual tmp< surfaceScalarField > alphaPhi() const
Return the volumetric flux of the phase.
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction.
rhoReactionThermo & thermo
Definition: createFields.H:28
virtual void correctTurbulence()
Correct the turbulence.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
CGAL::Exact_predicates_exact_constructions_kernel K
phi
Definition: pEqn.H:104
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:57
Calculate the first temporal derivative.
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
dynamicFvMesh & mesh
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
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.
phaseSystem & fluid
Definition: createFields.H:11
word timeName
Definition: getTimeIndex.H:3
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 > &)
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(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/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:68
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.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
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.
label patchi
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.
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
virtual tmp< scalarField > kappaEff(const label patchi) const
Return the effective thermal conductivity on a patch.
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.