phasePressureModel.H
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) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 2020 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
27Class
28 Foam::phasePressureModel
29
30Description
31 Particle-particle phase-pressure RAS model
32
33 The derivative of the phase-pressure with respect to the phase-fraction
34 is evaluated as
35
36 g0*min(exp(preAlphaExp*(alpha - alphaMax)), expMax)
37
38 The default model coefficients correspond to the following:
39 \verbatim
40 phasePressureCoeffs
41 {
42 preAlphaExp 500;
43 expMax 1000;
44 alphaMax 0.62;
45 g0 1000;
46 }
47 \endverbatim
48
49SourceFiles
50 phasePressureModel.C
51
52\*---------------------------------------------------------------------------*/
53
54#ifndef phasePressureModel_H
55#define phasePressureModel_H
56
57#include "RASModel.H"
58#include "eddyViscosity.H"
60#include "ThermalDiffusivity.H"
61#include "EddyDiffusivity.H"
62#include "phaseModel.H"
63
64// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65
66namespace Foam
67{
68namespace RASModels
69{
70
71/*---------------------------------------------------------------------------*\
72 Class phasePressureModel Declaration
73\*---------------------------------------------------------------------------*/
74
76:
77 public eddyViscosity
78 <
79 RASModel<EddyDiffusivity<ThermalDiffusivity
80 <
81 PhaseCompressibleTurbulenceModel<phaseModel>
82 >>>
83 >
84{
85 // Private Data
86
87 // Kinetic Theory Model coefficients
88
89 //- Maximum packing phase-fraction
90 scalar alphaMax_;
91
92 //- Pre-exponential factor
93 scalar preAlphaExp_;
94
95 //- Maximum limit of the exponential
96 scalar expMax_;
97
98 //- g0
100
101
102 // Private Member Functions
103
104 //- Disabled
105 void correctNut()
106 {}
107
108 //- No copy construct
109 phasePressureModel(const phasePressureModel&) = delete;
110
111 //- No copy assignment
112 void operator=(const phasePressureModel&) = delete;
113
114
115public:
116
117 //- Runtime type information
118 TypeName("phasePressure");
119
120
121 // Constructors
122
123 //- Construct from components
125 (
126 const volScalarField& alpha,
127 const volScalarField& rho,
128 const volVectorField& U,
129 const surfaceScalarField& alphaRhoPhi,
130 const surfaceScalarField& phi,
131 const phaseModel& transport,
132 const word& propertiesName = turbulenceModel::propertiesName,
133 const word& type = typeName
134 );
135
136
137 //- Destructor
138 virtual ~phasePressureModel() = default;
139
140
141 // Member Functions
142
143 //- Re-read model coefficients if they have changed
144 virtual bool read();
145
146 //- Return the effective viscosity
147 virtual tmp<volScalarField> nuEff() const
148 {
149 return this->nut();
150 }
151
152 //- Return the effective viscosity on patch
153 virtual tmp<scalarField> nuEff(const label patchi) const
154 {
155 return this->nut(patchi);
156 }
157
158 //- Return the turbulence kinetic energy
159 virtual tmp<volScalarField> k() const;
160
161 //- Return the turbulence kinetic energy dissipation rate
162 virtual tmp<volScalarField> epsilon() const;
163
164 //- Return the specific dissipation rate
165 virtual tmp<volScalarField> omega() const;
166
167 //- Return the Reynolds stress tensor
168 virtual tmp<volSymmTensorField> R() const;
169
170 //- Return the phase-pressure'
171 // (derivative of phase-pressure w.r.t. phase-fraction)
172 virtual tmp<volScalarField> pPrime() const;
173
174 //- Return the face-phase-pressure'
175 // (derivative of phase-pressure w.r.t. phase-fraction)
176 virtual tmp<surfaceScalarField> pPrimef() const;
177
178 //- Return the effective stress tensor
179 virtual tmp<volSymmTensorField> devRhoReff() const;
180
181 //- Return the effective stress tensor based on a given velocity field
183 (
184 const volVectorField& U
185 ) const;
186
187 //- Return the source term for the momentum equation
189
190 //- Solve the kinetic theory equations and correct the viscosity
191 virtual void correct();
192};
193
194
195// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196
197} // End namespace RASModels
198} // End namespace Foam
199
200// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201
202#endif
203
204// ************************************************************************* //
surfaceScalarField & phi
Particle-particle phase-pressure RAS model.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
TypeName("phasePressure")
Runtime type information.
virtual tmp< scalarField > nuEff(const label patchi) const
Return the effective viscosity on patch.
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual ~phasePressureModel()=default
Destructor.
virtual tmp< volSymmTensorField > devRhoReff(const volVectorField &U) const
Return the effective stress tensor based on a given velocity field.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
phasePressureModel(const volScalarField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const phaseModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
virtual bool read()
Re-read model coefficients if they have changed.
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
virtual bool read()
Re-read model coefficients if they have changed.
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
Particle-particle phase-pressure RAS model.
A class for managing temporary objects.
Definition: tmp.H:65
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
volScalarField & alpha
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73