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  DUDt_(nullptr),
166  DUDtf_(nullptr),
167  divU_(nullptr),
168  momentumTransport_
169  (
170  phaseCompressible::momentumTransportModel::New
171  (
172  *this,
173  this->rho(),
174  U_,
175  alphaRhoPhi_,
176  phi_,
177  *this
178  )
179  ),
180  thermophysicalTransport_
181  (
183  <
184  phaseCompressible::momentumTransportModel,
186  >::New(momentumTransport_, this->thermo_)
187  ),
188  continuityError_
189  (
190  IOobject
191  (
192  IOobject::groupName("continuityError", this->name()),
193  fluid.mesh().time().name(),
194  fluid.mesh()
195  ),
196  fluid.mesh(),
198  ),
199  K_(nullptr)
200 {
202 
203  if (fluid.mesh().dynamic() || this->fluid().MRF().size())
204  {
205  Uf_ = new surfaceVectorField
206  (
207  IOobject
208  (
209  IOobject::groupName("Uf", this->name()),
210  fluid.mesh().time().name(),
211  fluid.mesh(),
214  ),
216  );
217  }
218 
220 }
221 
222 
223 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
224 
225 template<class BasePhaseModel>
227 {}
228 
229 
230 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
231 
232 template<class BasePhaseModel>
234 (
235  const volScalarField& source
236 )
237 {
238  volScalarField& rho = this->rho();
239  continuityError_ = fvc::ddt(*this, rho) + fvc::div(alphaRhoPhi_) - source;
240 }
241 
242 
243 template<class BasePhaseModel>
245 {
247 }
248 
249 
250 template<class BasePhaseModel>
252 {
253  BasePhaseModel::correctKinematics();
254 
255  if (DUDt_.valid())
256  {
257  DUDt_.clear();
258  DUDt();
259  }
260 
261  if (DUDtf_.valid())
262  {
263  DUDtf_.clear();
264  DUDtf();
265  }
266 
267  if (K_.valid())
268  {
269  K_.ref() = 0.5*magSqr(this->U());
270  }
271 }
272 
273 
274 template<class BasePhaseModel>
276 {
277  BasePhaseModel::predictMomentumTransport();
278  momentumTransport_->predict();
279 }
280 
281 
282 template<class BasePhaseModel>
284 {
285  BasePhaseModel::predictThermophysicalTransport();
286  thermophysicalTransport_->predict();
287 }
288 
289 
290 template<class BasePhaseModel>
292 {
293  BasePhaseModel::correctMomentumTransport();
294  momentumTransport_->correct();
295 }
296 
297 
298 template<class BasePhaseModel>
300 {
301  BasePhaseModel::correctThermophysicalTransport();
302  thermophysicalTransport_->correct();
303 }
304 
305 
306 template<class BasePhaseModel>
308 {
309  const fvMesh& mesh = this->fluid().mesh();
310 
311  if (Uf_.valid())
312  {
313  Uf_() = fvc::interpolate(U_);
314  surfaceVectorField n(mesh.Sf()/mesh.magSf());
315  Uf_() +=
316  n*(
317  this->fluid().MRF().absolute(fvc::absolute(phi_, U_))
318  /mesh.magSf()
319  - (n & Uf_())
320  );
321  }
322 }
323 
324 
325 template<class BasePhaseModel>
327 {
328  return false;
329 }
330 
331 
332 template<class BasePhaseModel>
335 {
336  const volScalarField& alpha = *this;
337  const volScalarField& rho = this->rho();
338 
339  return
340  (
341  fvm::ddt(alpha, rho, U_)
342  + fvm::div(alphaRhoPhi_, U_)
343  + fvm::SuSp(-this->continuityError(), U_)
344  + this->fluid().MRF().DDt(alpha*rho, U_)
345  + momentumTransport_->divDevTau(U_)
346  );
347 }
348 
349 
350 template<class BasePhaseModel>
353 {
354  // As the "normal" U-eqn but without the ddt terms
355 
356  const volScalarField& alpha = *this;
357  const volScalarField& rho = this->rho();
358 
359  return
360  (
361  fvm::div(alphaRhoPhi_, U_)
362  + fvm::SuSp(fvc::ddt(*this, rho) - this->continuityError(), U_)
363  + this->fluid().MRF().DDt(alpha*rho, U_)
364  + momentumTransport_->divDevTau(U_)
365  );
366 }
367 
368 
369 template<class BasePhaseModel>
372 {
373  return U_;
374 }
375 
376 
377 template<class BasePhaseModel>
380 {
381  return U_;
382 }
383 
384 
385 template<class BasePhaseModel>
388 {
389  return U_;
390 }
391 
392 
393 template<class BasePhaseModel>
396 {
397  return phi_;
398 }
399 
400 
401 template<class BasePhaseModel>
404 {
405  return phi_;
406 }
407 
408 
409 template<class BasePhaseModel>
412 {
413  return phi_;
414 }
415 
416 
417 template<class BasePhaseModel>
420 {
421  return Uf_;
422 }
423 
424 
425 template<class BasePhaseModel>
428 {
429  if (Uf_.valid())
430  {
431  return Uf_();
432  }
433  else
434  {
436  << "Uf has not been allocated."
437  << exit(FatalError);
438 
439  return const_cast<surfaceVectorField&>(surfaceVectorField::null());
440  }
441 }
442 
443 
444 template<class BasePhaseModel>
447 {
448  if (Uf_.valid())
449  {
450  return Uf_();
451  }
452  else
453  {
455  << "Uf has not been allocated."
456  << exit(FatalError);
457 
458  return const_cast<surfaceVectorField&>(surfaceVectorField::null());
459  }
460 }
461 
462 
463 template<class BasePhaseModel>
466 {
467  return alphaPhi_;
468 }
469 
470 
471 template<class BasePhaseModel>
474 {
475  return alphaPhi_;
476 }
477 
478 
479 template<class BasePhaseModel>
482 {
483  return alphaPhi_;
484 }
485 
486 
487 template<class BasePhaseModel>
490 {
491  return alphaRhoPhi_;
492 }
493 
494 
495 template<class BasePhaseModel>
498 {
499  return alphaRhoPhi_;
500 }
501 
502 
503 template<class BasePhaseModel>
506 {
507  return alphaRhoPhi_;
508 }
509 
510 
511 template<class BasePhaseModel>
514 {
515  if (!DUDt_.valid())
516  {
517  const tmp<surfaceScalarField> taphi(fvc::absolute(phi_, U_));
518  const surfaceScalarField& aphi(taphi());
519  DUDt_ =
520  new volVectorField
521  (
522  IOobject::groupName("DUDt", this->name()),
523  fvc::ddt(U_) + fvc::div(aphi, U_) - fvc::div(aphi)*U_
524  );
525  }
526 
527  return tmp<volVectorField>(DUDt_());
528 }
529 
530 
531 template<class BasePhaseModel>
534 {
535  if (!DUDtf_.valid())
536  {
537  DUDtf_ =
539  (
540  IOobject::groupName("DUDtf", this->name()),
541  byDt(phi_ - phi_.oldTime())
542  );
543  }
544 
545  return tmp<surfaceScalarField>(DUDtf_());
546 }
547 
548 
549 template<class BasePhaseModel>
552 {
553  return continuityError_;
554 }
555 
556 
557 template<class BasePhaseModel>
560 {
561  if (!K_.valid())
562  {
563  K_ =
564  new volScalarField
565  (
566  IOobject::groupName("K", this->name()),
567  0.5*magSqr(this->U())
568  );
569  }
570 
571  return tmp<volScalarField>(K_());
572 }
573 
574 
575 template<class BasePhaseModel>
578 {
579  return divU_;
580 }
581 
582 
583 template<class BasePhaseModel>
585 {
586  if (!divU_.valid())
587  {
588  divU_ = divU.ptr();
589  divU_().rename(IOobject::groupName("divU", this->name()));
590  divU_().checkIn();
591  }
592  else
593  {
594  divU_() = divU;
595  }
596 }
597 
598 
599 template<class BasePhaseModel>
602 {
603  return momentumTransport_->k();
604 }
605 
606 
607 template<class BasePhaseModel>
610 {
611  return momentumTransport_->pPrime();
612 }
613 
614 
615 template<class BasePhaseModel>
618 {
619  return thermophysicalTransport_->kappaEff(patchi);
620 }
621 
622 
623 template<class BasePhaseModel>
626 {
627  return thermophysicalTransport_->divq(he);
628 }
629 
630 
631 template<class BasePhaseModel>
634 {
635  return thermophysicalTransport_->divj(Yi);
636 }
637 
638 
639 // ************************************************************************* //
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.
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
MovingPhaseModelTransportThermoModel< typename BasePhaseModel::thermoModel >::type transportThermoModel
Thermo type for the thermophysical transport model.
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 void predictThermophysicalTransport()
Predict the energy transport e.g. alphat.
virtual ~MovingPhaseModel()
Destructor.
virtual tmp< volVectorField > U() const
Return the velocity.
virtual void correctMomentumTransport()
Correct the momentumTransport.
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction.
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 void correctThermophysicalTransport()
Correct the energy transport e.g. alphat.
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.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
surfaceScalarField phi_
Flux.
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
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 tmp< volVectorField > DUDt() const
Return the substantive acceleration.
virtual bool stationary() const
Return whether the phase is stationary.
virtual void predictMomentumTransport()
Predict the momentumTransport.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
virtual tmp< scalarField > kappaEff(const label patchi) const
Return the effective thermal conductivity on a patch.
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.
autoPtr< surfaceVectorField > Uf_
Face velocity field.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Templated base class for multiphase thermophysical transport models.
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:122
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:406
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:641
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:28
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:306
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< 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 > > DDt(const surfaceScalarField &phi, const VolField< Type > &psi)
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 > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
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:62
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:251
SurfaceField< scalar > surfaceScalarField
messageStream Info
tmp< volScalarField > byDt(const volScalarField &vf)
Definition: phaseSystem.C:853
const dimensionSet dimTime
const dimensionSet dimDensity
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
error FatalError
SurfaceField< vector > surfaceVectorField
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensioned< scalar > magSqr(const dimensioned< Type > &)
thermo he()