phasePressureModel.H
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) 2013-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 Class
25  Foam::RASModels::phasePressureModel
26 
27 Description
28  Particle-particle phase-pressure RAS model
29 
30  The derivative of the phase-pressure with respect to the phase-fraction
31  is evaluated as
32 
33  g0*min(exp(preAlphaExp*(alpha - alphaMax)), expMax)
34 
35  The default model coefficients correspond to the following:
36  \verbatim
37  phasePressureCoeffs
38  {
39  preAlphaExp 500;
40  expMax 1000;
41  g0 1000;
42  }
43  \endverbatim
44 
45 SourceFiles
46  phasePressureModel.C
47 
48 \*---------------------------------------------------------------------------*/
49 
50 #ifndef phasePressureModel_H
51 #define phasePressureModel_H
52 
53 #include "RASModel.H"
54 #include "eddyViscosity.H"
56 #include "phaseModel.H"
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 namespace Foam
61 {
62 namespace RASModels
63 {
64 
65 /*---------------------------------------------------------------------------*\
66  Class phasePressureModel Declaration
67 \*---------------------------------------------------------------------------*/
68 
70 :
71  public eddyViscosity<RASModel<phaseCompressible::momentumTransportModel>>
72 {
73  // Private Data
74 
75  const phaseModel& phase_;
76 
77  // Phase pressure coefficients
78 
79  //- Pre-exponential factor
80  scalar preAlphaExp_;
81 
82  //- Maximum limit of the exponential
83  scalar expMax_;
84 
85  //- g0
87 
88 
89  // Private Member Functions
90 
91  void correctNut()
92  {}
93 
94  //- Return the phase-pressure'
95  // (derivative of phase-pressure w.r.t. phase-fraction)
96  tmp<volScalarField> pPrime() const;
97 
98 
99 public:
100 
101  //- Runtime type information
102  TypeName("phasePressure");
103 
104 
105  // Constructors
106 
107  //- Construct from components
109  (
110  const volScalarField& alpha,
111  const volScalarField& rho,
112  const volVectorField& U,
113  const surfaceScalarField& alphaRhoPhi,
114  const surfaceScalarField& phi,
115  const viscosity& viscosity,
116  const word& type = typeName
117  );
118 
119  //- Disallow default bitwise copy construction
120  phasePressureModel(const phasePressureModel&) = delete;
121 
122 
123  //- Destructor
124  virtual ~phasePressureModel();
125 
126 
127  // Member Functions
128 
129  //- Re-read model coefficients if they have changed
130  virtual bool read();
131 
132  //- Return the effective viscosity
133  virtual tmp<volScalarField> nuEff() const
134  {
135  return this->nut();
136  }
137 
138  //- Return the effective viscosity on patch
139  virtual tmp<scalarField> nuEff(const label patchi) const
140  {
141  return this->nut(patchi);
142  }
143 
144  //- Return the turbulence kinetic energy
145  virtual tmp<volScalarField> k() const;
146 
147  //- Return the turbulence kinetic energy dissipation rate
148  virtual tmp<volScalarField> epsilon() const;
149 
150  //- Return the turbulence specific dissipation rate
151  virtual tmp<volScalarField> omega() const;
152 
153  //- Return the stress tensor [m^2/s^2]
154  virtual tmp<volSymmTensorField> sigma() const;
155 
156  //- Return the face-phase-pressure'
157  // (derivative of phase-pressure w.r.t. phase-fraction)
158  virtual tmp<surfaceScalarField> pPrimef() const;
159 
160  //- Return the effective stress tensor
161  virtual tmp<volSymmTensorField> devTau() const;
162 
163  //- Return the source term for the momentum equation
165 
166  //- Solve the kinetic theory equations and correct the viscosity
167  virtual void correct();
168 
169 
170  // Member Operators
171 
172  //- Disallow default bitwise assignment
173  void operator=(const phasePressureModel&) = delete;
174 };
175 
176 
177 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 
179 } // End namespace RASModels
180 } // End namespace Foam
181 
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 
184 #endif
185 
186 // ************************************************************************* //
Generic GeometricField class.
Particle-particle phase-pressure RAS model.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
TypeName("phasePressure")
Runtime type information.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual tmp< volSymmTensorField > devTau() const
Return the effective stress tensor.
phasePressureModel(const volScalarField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the source term for the momentum equation.
virtual tmp< volSymmTensorField > sigma() const
Return the stress tensor [m^2/s^2].
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
void operator=(const phasePressureModel &)=delete
Disallow default bitwise assignment.
virtual bool read()
Re-read model coefficients if they have changed.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity.
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
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))
Namespace for OpenFOAM.
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488