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-2019 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  alphaMax 0.62;
42  g0 1000;
43  }
44  \endverbatim
45 
46 SourceFiles
47  phasePressureModel.C
48 
49 \*---------------------------------------------------------------------------*/
50 
51 #ifndef phasePressureModel_H
52 #define phasePressureModel_H
53 
54 #include "RASModel.H"
55 #include "eddyViscosity.H"
57 #include "ThermalDiffusivity.H"
58 #include "EddyDiffusivity.H"
59 #include "phaseModel.H"
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 namespace Foam
64 {
65 namespace RASModels
66 {
67 
68 /*---------------------------------------------------------------------------*\
69  Class phasePressureModel Declaration
70 \*---------------------------------------------------------------------------*/
71 
72 class phasePressureModel
73 :
74  public eddyViscosity
75  <
76  RASModel<EddyDiffusivity<ThermalDiffusivity
77  <
78  PhaseCompressibleTurbulenceModel<phaseModel>
79  >>>
80  >
81 {
82  // Private Data
83 
84  // Input Fields
85 
86  const phaseModel& phase_;
87 
88 
89  // Kinetic Theory Model coefficients
90 
91  //- Maximum packing phase-fraction
92  scalar alphaMax_;
93 
94  //- Pre-exponential factor
95  scalar preAlphaExp_;
96 
97  //- Maximum limit of the exponential
98  scalar expMax_;
99 
100  //- g0
101  dimensionedScalar g0_;
102 
103 
104  // Private Member Functions
105 
106  void correctNut()
107  {}
108 
109 
110 public:
111 
112  //- Runtime type information
113  TypeName("phasePressure");
114 
115 
116  // Constructors
117 
118  //- Construct from components
120  (
121  const volScalarField& alpha,
122  const volScalarField& rho,
123  const volVectorField& U,
125  const surfaceScalarField& phi,
126  const phaseModel& transport,
128  const word& type = typeName
129  );
130 
131  //- Disallow default bitwise copy construction
132  phasePressureModel(const phasePressureModel&) = delete;
133 
134 
135  //- Destructor
136  virtual ~phasePressureModel();
137 
138 
139  // Member Functions
140 
141  //- Re-read model coefficients if they have changed
142  virtual bool read();
143 
144  //- Return the effective viscosity
145  virtual tmp<volScalarField> nuEff() const
146  {
147  return this->nut();
148  }
149 
150  //- Return the effective viscosity on patch
151  virtual tmp<scalarField> nuEff(const label patchi) const
152  {
153  return this->nut(patchi);
154  }
155 
156  //- Return the turbulence kinetic energy
157  virtual tmp<volScalarField> k() const;
158 
159  //- Return the turbulence kinetic energy dissipation rate
160  virtual tmp<volScalarField> epsilon() const;
161 
162  //- Return the Reynolds stress tensor
163  virtual tmp<volSymmTensorField> R() const;
164 
165  //- Return the phase-pressure'
166  // (derivative of phase-pressure w.r.t. phase-fraction)
167  virtual tmp<volScalarField> pPrime() const;
168 
169  //- Return the face-phase-pressure'
170  // (derivative of phase-pressure w.r.t. phase-fraction)
171  virtual tmp<surfaceScalarField> pPrimef() const;
172 
173  //- Return the effective stress tensor
174  virtual tmp<volSymmTensorField> devRhoReff() const;
175 
176  //- Return the source term for the momentum equation
178 
179  //- Solve the kinetic theory equations and correct the viscosity
180  virtual void correct();
181 
182 
183  // Member Operators
184 
185  //- Disallow default bitwise assignment
186  void operator=(const phasePressureModel&) = delete;
187 };
188 
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 } // End namespace RASModels
193 } // End namespace Foam
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 #endif
198 
199 // ************************************************************************* //
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
const volVectorField & U() const
Access function to velocity field.
Particle-particle phase-pressure RAS model.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
const transportModel & transport() const
Access function to incompressible transport model.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
const surfaceScalarField & alphaRhoPhi() const
Access function to phase flux field.
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
const volScalarField & rho() const
Return the density field.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual bool read()
Re-read model coefficients if they have changed.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
virtual ~phasePressureModel()
Destructor.
label patchi
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.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void operator=(const phasePressureModel &)=delete
Disallow default bitwise assignment.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux field.
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
virtual tmp< volScalarField > alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
TypeName("phasePressure")
Runtime type information.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
Namespace for OpenFOAM.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure&#39;.