Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
47template<class BasePhaseModel>
49Foam::MovingPhaseModel<BasePhaseModel>::phi(const volVectorField& U) const
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(),
76 fvsPatchFieldBase::calculatedType()
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
92 return tmp<surfaceScalarField>::New
93 (
94 io,
95 fvc::flux(U),
97 );
98 }
99}
100
101
102// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
103
104template<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 (
154 phaseCompressibleTurbulenceModel::New
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(),
173 dimensionedScalar(dimDensity/dimTime)
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(),
184 dimensionedScalar(dimDensity/dimTime)
185 ),
186 K_(nullptr)
187{
189
191}
192
193
194// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
195
196template<class BasePhaseModel>
198{
199 BasePhaseModel::correct();
200
201 this->fluid().MRF().correctBoundaryVelocity(U_);
202
203 volScalarField& rho = this->thermoRef().rho();
204
206
207 continuityErrorSources_ = - (this->fluid().fvOptions()(*this, rho)&rho);
208}
209
210
211template<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())
230 K_.ref() = 0.5*magSqr(this->U());
231 }
232}
233
234
235template<class BasePhaseModel>
237{
238 BasePhaseModel::correctTurbulence();
239
240 turbulence_->correct();
241}
242
243
244template<class BasePhaseModel>
246{
247 BasePhaseModel::correctEnergyTransport();
248
249 turbulence_->correctEnergyTransport();
250}
251
252
253template<class BasePhaseModel>
256 return false;
257}
258
259
260template<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
279template<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
298template<class BasePhaseModel>
302 return U_;
303}
304
305
306template<class BasePhaseModel>
309{
310 return U_;
311}
312
313
314template<class BasePhaseModel>
318 return phi_;
319}
320
321
322template<class BasePhaseModel>
325{
326 return phi_;
327}
328
329
330template<class BasePhaseModel>
334 return alphaPhi_;
335}
336
337
338template<class BasePhaseModel>
342 return alphaPhi_;
343}
344
345
346template<class BasePhaseModel>
350 return alphaRhoPhi_;
351}
352
353
354template<class BasePhaseModel>
358 return alphaRhoPhi_;
359}
360
361
362template<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
375template<class BasePhaseModel>
378{
379 if (!DUDtf_.valid())
380 {
381 DUDtf_ = byDt(phi_ - phi_.oldTime());
382 }
385}
386
387
388template<class BasePhaseModel>
395
396template<class BasePhaseModel>
401}
402
403
404template<class BasePhaseModel>
409}
410
411
412template<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
429template<class BasePhaseModel>
435
436
437template<class BasePhaseModel>
440 divU_ = divU;
441}
442
443
444template<class BasePhaseModel>
448 return turbulence_->mut();
449}
450
451
452template<class BasePhaseModel>
456 return turbulence_->muEff();
457}
458
459
460template<class BasePhaseModel>
464 return turbulence_->nut();
465}
466
467
468template<class BasePhaseModel>
472 return turbulence_->nuEff();
473}
474
475
476template<class BasePhaseModel>
480 return turbulence_->kappaEff();
481}
482
483
484template<class BasePhaseModel>
488 return turbulence_->kappaEff(patchi);
489}
490
491
492template<class BasePhaseModel>
496 return turbulence_->alphaEff();
497}
498
499
500template<class BasePhaseModel>
504 return turbulence_->alphaEff(patchi);
505}
506
507
508template<class BasePhaseModel>
512 return turbulence_->k();
513}
514
515
516template<class BasePhaseModel>
519{
520 return turbulence_->pPrime();
521}
522
523
524// ************************************************************************* //
IOMRFZoneList & MRF
twoPhaseSystem & fluid
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())
writeOption writeOpt() const noexcept
Get the write option.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
virtual void correctKinematics()
Correct the kinematics.
virtual tmp< volScalarField > continuityErrorSources() const
Return the continuity error due to any sources.
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
virtual tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
tmp< volScalarField > divU_
Dilatation rate.
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
virtual tmp< volVectorField > U() const
Access const reference to U.
tmp< surfaceScalarField > DUDtf_
Lagrangian acceleration field on the faces (needed for virtual-mass).
virtual tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
volScalarField continuityErrorFlow_
Continuity error due to the flow.
autoPtr< phaseCompressibleTurbulenceModel > turbulence_
Turbulence model.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
surfaceScalarField alphaRhoPhi_
Mass flux.
virtual void correctTurbulence()
Correct the turbulence.
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
virtual tmp< volScalarField > alphaEff() const
Return the effective thermal diffusivity.
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access 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.
virtual tmp< volScalarField > continuityErrorFlow() const
Return the continuity error due to the flow field.
virtual tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
MovingPhaseModel(const multiphaseInterSystem &fluid, const word &phaseName)
virtual tmp< volScalarField > divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi)).
tmp< volVectorField > DUDt_
Lagrangian acceleration field (needed for virtual-mass).
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
Constant access the volumetric flux.
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
virtual tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
virtual volVectorField & URef()
Access the velocity.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
virtual tmp< volScalarField > kappaEff() const
Return the effective thermal conductivity.
virtual bool stationary() const
Return whether the phase is stationary.
tmp< volScalarField > K_
Kinetic Energy.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
volScalarField continuityErrorSources_
Continuity error due to any sources.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Class to represent a system of phases and model interfacial transfers between them.
Definition phaseSystem.H:72
const IOMRFZoneList & MRF() const
Return MRF zones.
fv::options & fvOptions() const
Access the fvOptions.
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...
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
dynamicFvMesh & mesh
const auto & io
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 finiteVolume matrix for implicit and explicit sources.
auto & name
zeroField divU
Definition alphaSuSp.H:3
word timeName
Definition getTimeIndex.H:3
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition fvcDdt.C:40
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
List< word > wordList
List of word.
Definition fileName.H:60
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.
GeometricField< vector, fvPatchField, volMesh > volVectorField
tmp< volScalarField > byDt(const volScalarField &vf)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
wordList patchTypes(nPatches)
volScalarField & alpha
tmp< volScalarField > trho
psiReactionThermo & thermo
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299