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-2023 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  typeIOobject<surfaceScalarField> phiHeader
48  (
49  phiName,
50  U.mesh().time().name(),
51  U.mesh(),
52  IOobject::NO_READ
53  );
54 
55  if (phiHeader.headerOk())
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().name(),
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().name(),
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().name(),
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().name(),
144  fluid.mesh(),
145  IOobject::READ_IF_PRESENT,
146  IOobject::NO_WRITE
147  ),
148  fluid.mesh(),
149  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), 0)
150  ),
151  alphaRhoPhi_
152  (
153  IOobject
154  (
155  IOobject::groupName("alphaRhoPhi", this->name()),
156  fluid.mesh().time().name(),
157  fluid.mesh(),
158  IOobject::READ_IF_PRESENT,
159  IOobject::AUTO_WRITE
160  ),
161  fluid.mesh(),
162  dimensionedScalar(dimensionSet(1, 0, -1, 0, 0), 0)
163  ),
164  Uf_(nullptr),
165  divU_(nullptr),
166  momentumTransport_
167  (
168  phaseCompressible::momentumTransportModel::New
169  (
170  *this,
171  this->rho(),
172  U_,
173  alphaRhoPhi_,
174  phi_,
175  *this
176  )
177  ),
178  continuityError_
179  (
180  IOobject
181  (
182  IOobject::groupName("continuityError", this->name()),
183  fluid.mesh().time().name(),
184  fluid.mesh()
185  ),
186  fluid.mesh(),
188  ),
189  K_(nullptr)
190 {
192 
193  if (fluid.mesh().dynamic() || this->fluid().MRF().size())
194  {
195  Uf_ = new surfaceVectorField
196  (
197  IOobject
198  (
199  IOobject::groupName("Uf", this->name()),
200  fluid.mesh().time().name(),
201  fluid.mesh(),
204  ),
206  );
207  }
208 
210 }
211 
212 
213 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
214 
215 template<class BasePhaseModel>
217 {}
218 
219 
220 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221 
222 template<class BasePhaseModel>
224 (
225  const volScalarField& source
226 )
227 {
228  volScalarField& rho = this->rho();
229  continuityError_ = fvc::ddt(*this, rho) + fvc::div(alphaRhoPhi_) - source;
230 }
231 
232 
233 template<class BasePhaseModel>
235 {
237 }
238 
239 
240 template<class BasePhaseModel>
242 {
243  BasePhaseModel::correctKinematics();
244 
245  if (K_.valid())
246  {
247  K_.ref() = 0.5*magSqr(this->U());
248  }
249 }
250 
251 
252 template<class BasePhaseModel>
254 {
255  BasePhaseModel::predictMomentumTransport();
256  momentumTransport_->predict();
257 }
258 
259 
260 template<class BasePhaseModel>
262 {
263  BasePhaseModel::correctMomentumTransport();
264  momentumTransport_->correct();
265 }
266 
267 
268 template<class BasePhaseModel>
270 {
271  const fvMesh& mesh = this->fluid().mesh();
272 
273  if (Uf_.valid())
274  {
275  Uf_() = fvc::interpolate(U_);
276  surfaceVectorField n(mesh.Sf()/mesh.magSf());
277  Uf_() +=
278  n*(
279  this->fluid().MRF().absolute(fvc::absolute(phi_, U_))
280  /mesh.magSf()
281  - (n & Uf_())
282  );
283  }
284 }
285 
286 
287 template<class BasePhaseModel>
289 {
290  return false;
291 }
292 
293 
294 template<class BasePhaseModel>
297 {
298  const volScalarField& alpha = *this;
299  const volScalarField& rho = this->rho();
300 
301  return
302  (
303  fvm::ddt(alpha, rho, U_)
304  + fvm::div(alphaRhoPhi_, U_)
305  + fvm::SuSp(-this->continuityError(), U_)
306  + this->fluid().MRF().DDt(alpha*rho, U_)
307  + momentumTransport_->divDevTau(U_)
308  );
309 }
310 
311 
312 template<class BasePhaseModel>
315 {
316  // As the "normal" U-eqn but without the ddt terms
317 
318  const volScalarField& alpha = *this;
319  const volScalarField& rho = this->rho();
320 
321  return
322  (
323  fvm::div(alphaRhoPhi_, U_)
324  + fvm::SuSp(fvc::ddt(*this, rho) - this->continuityError(), U_)
325  + this->fluid().MRF().DDt(alpha*rho, U_)
326  + momentumTransport_->divDevTau(U_)
327  );
328 }
329 
330 
331 template<class BasePhaseModel>
334 {
335  return U_;
336 }
337 
338 
339 template<class BasePhaseModel>
342 {
343  return U_;
344 }
345 
346 
347 template<class BasePhaseModel>
350 {
351  return U_;
352 }
353 
354 
355 template<class BasePhaseModel>
358 {
359  return phi_;
360 }
361 
362 
363 template<class BasePhaseModel>
366 {
367  return phi_;
368 }
369 
370 
371 template<class BasePhaseModel>
374 {
375  return phi_;
376 }
377 
378 
379 template<class BasePhaseModel>
382 {
383  return Uf_;
384 }
385 
386 
387 template<class BasePhaseModel>
390 {
391  if (Uf_.valid())
392  {
393  return Uf_();
394  }
395  else
396  {
398  << "Uf has not been allocated."
399  << exit(FatalError);
400 
401  return const_cast<surfaceVectorField&>(surfaceVectorField::null());
402  }
403 }
404 
405 
406 template<class BasePhaseModel>
409 {
410  if (Uf_.valid())
411  {
412  return Uf_();
413  }
414  else
415  {
417  << "Uf has not been allocated."
418  << exit(FatalError);
419 
420  return const_cast<surfaceVectorField&>(surfaceVectorField::null());
421  }
422 }
423 
424 
425 template<class BasePhaseModel>
428 {
429  return alphaPhi_;
430 }
431 
432 
433 template<class BasePhaseModel>
436 {
437  return alphaPhi_;
438 }
439 
440 
441 template<class BasePhaseModel>
444 {
445  return alphaPhi_;
446 }
447 
448 
449 template<class BasePhaseModel>
452 {
453  return alphaRhoPhi_;
454 }
455 
456 
457 template<class BasePhaseModel>
460 {
461  return alphaRhoPhi_;
462 }
463 
464 
465 template<class BasePhaseModel>
468 {
469  return alphaRhoPhi_;
470 }
471 
472 
473 template<class BasePhaseModel>
476 {
477  const tmp<surfaceScalarField> taphi(fvc::absolute(phi_, U_));
478  const surfaceScalarField& aphi(taphi());
479 
480  return
481  fvm::div(aphi, U_) - fvm::Sp(fvc::div(aphi), U_)
482  + this->fluid().MRF().DDt(U_);
483 }
484 
485 
486 template<class BasePhaseModel>
489 {
490  return fvm::ddt(U_) + UgradU();
491 }
492 
493 
494 template<class BasePhaseModel>
497 {
498  return continuityError_;
499 }
500 
501 
502 template<class BasePhaseModel>
505 {
506  if (!K_.valid())
507  {
508  K_ =
509  new volScalarField
510  (
511  IOobject::groupName("K", this->name()),
512  0.5*magSqr(this->U())
513  );
514  }
515 
516  return tmp<volScalarField>(K_());
517 }
518 
519 
520 template<class BasePhaseModel>
523 {
524  return divU_;
525 }
526 
527 
528 template<class BasePhaseModel>
530 {
531  if (!divU_.valid())
532  {
533  divU_ = divU.ptr();
534  divU_().rename(IOobject::groupName("divU", this->name()));
535  divU_().checkIn();
536  }
537  else
538  {
539  divU_() = divU;
540  }
541 }
542 
543 
544 template<class BasePhaseModel>
547 {
548  return momentumTransport_->k();
549 }
550 
551 
552 template<class BasePhaseModel>
555 {
556  return momentumTransport_->pPrimef();
557 }
558 
559 
560 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
writeOption writeOpt() const
Definition: IOobject.H:370
static word groupName(Name name, const word &group)
virtual void correctUf()
Correct the face velocity for moving meshes.
virtual void correctKinematics()
Correct the kinematics.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
volVectorField U_
Velocity field.
virtual void correct()
Correct the phase properties other than the thermo.
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
virtual tmp< fvVectorMatrix > DUDt() const
Return the substantive acceleration matrix.
virtual ~MovingPhaseModel()
Destructor.
virtual tmp< volVectorField > U() const
Return the velocity.
virtual void correctMomentumTransport()
Correct the momentumTransport.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
virtual const autoPtr< volScalarField > & divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
virtual tmp< surfaceScalarField > alphaPhi() const
Return the volumetric flux of the phase.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
surfaceScalarField phi_
Flux.
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
virtual const autoPtr< surfaceVectorField > & Uf() const
Return the face velocity.
virtual surfaceVectorField & UfRef()
Access the face velocity.
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux.
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
virtual volVectorField & URef()
Access the velocity.
virtual bool stationary() const
Return whether the phase is stationary.
virtual void predictMomentumTransport()
Predict the momentumTransport.
virtual tmp< fvVectorMatrix > UgradU() const
Return the velocity transport matrix.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.
autoPtr< surfaceVectorField > Uf_
Face velocity field.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Dimension set for the base types.
Definition: dimensionSet.H:125
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
const surfaceVectorField & Sf() const
Return cell face area vectors.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
bool dynamic() const
Is this mesh dynamic?
Definition: fvMesh.C:610
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:441
Abstract base class for turbulence models (RAS, LES and laminar).
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:89
A class for managing temporary objects.
Definition: tmp.H:55
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:205
A class for handling words, derived from string.
Definition: word.H:62
IOMRFZoneList MRF(mesh)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for implicit and explicit sources.
label patchi
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
tmp< VolField< Type > > DDt(const surfaceScalarField &phi, const VolField< Type > &psi)
Definition: fvcDDt.C:45
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< word > wordList
A List of words.
Definition: fileName.H:54
VolField< vector > volVectorField
Definition: volFieldsFwd.H:65
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
SurfaceField< scalar > surfaceScalarField
messageStream Info
const dimensionSet dimTime
const dimensionSet dimDensity
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
error FatalError
SurfaceField< vector > surfaceVectorField
dimensioned< scalar > magSqr(const dimensioned< Type > &)