PhaseCompressibleTurbulenceModel.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) 2019 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::PhaseCompressibleTurbulenceModel
29
30Description
31 Templated abstract base class for multiphase compressible
32 turbulence models.
33
34SourceFiles
35 PhaseCompressibleTurbulenceModel.C
36
37\*---------------------------------------------------------------------------*/
38
39#ifndef PhaseCompressibleTurbulenceModel_H
40#define PhaseCompressibleTurbulenceModel_H
41
42#include "TurbulenceModel.H"
44
45// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46
47namespace Foam
48{
49
50/*---------------------------------------------------------------------------*\
51 Class PhaseCompressibleTurbulenceModel Declaration
52\*---------------------------------------------------------------------------*/
53
54template<class TransportModel>
56:
57 public TurbulenceModel
58 <
59 volScalarField,
60 volScalarField,
61 compressibleTurbulenceModel,
62 TransportModel
63 >
64{
65
66public:
70 typedef TransportModel transportModel;
71
72
73 // Constructors
74
75 //- Construct
77 (
78 const word& type,
79 const alphaField& alpha,
80 const volScalarField& rho,
81 const volVectorField& U,
85 const word& propertiesName
86 );
87
88
89 // Selectors
90
91 //- Return a reference to the selected turbulence model
93 (
94 const alphaField& alpha,
95 const volScalarField& rho,
96 const volVectorField& U,
101 );
102
103
104 //- Destructor
105 virtual ~PhaseCompressibleTurbulenceModel() = default;
106
107
108 // Member Functions
109
110 //- Return the laminar dynamic viscosity
111 virtual tmp<volScalarField> mu() const
112 {
113 return this->transport_.mu();
114 }
115
116 //- Return the laminar dynamic viscosity on patch
117 virtual tmp<scalarField> mu(const label patchi) const
118 {
119 return this->transport_.mu(patchi);
120 }
121
122 //- Return the turbulence dynamic viscosity
123 virtual tmp<volScalarField> mut() const
124 {
125 return this->rho_*this->nut();
126 }
127
128 //- Return the turbulence dynamic viscosity on patch
129 virtual tmp<scalarField> mut(const label patchi) const
130 {
131 return this->rho_.boundaryField()[patchi]*this->nut(patchi);
132 }
133
134 //- Return the effective dynamic viscosity
135 virtual tmp<volScalarField> muEff() const
136 {
137 return mut() + mu();
138 }
139
140 //- Return the effective dynamic viscosity on patch
141 virtual tmp<scalarField> muEff(const label patchi) const
142 {
143 return mut(patchi) + mu(patchi);
144 }
145
146 //- Return the phase-pressure'
147 // (derivative of phase-pressure w.r.t. phase-fraction)
148 virtual tmp<volScalarField> pPrime() const;
149
150 //- Return the face-phase-pressure'
151 // (derivative of phase-pressure w.r.t. phase-fraction)
152 virtual tmp<surfaceScalarField> pPrimef() const;
153};
154
155
156// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157
158} // End namespace Foam
159
160// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161
162#ifdef NoRepository
164#endif
165
166// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167
168#endif
169
170// ************************************************************************* //
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Templated abstract base class for multiphase compressible turbulence models.
virtual tmp< scalarField > muEff(const label patchi) const
Return the effective dynamic viscosity on patch.
virtual tmp< scalarField > mu(const label patchi) const
Return the laminar dynamic viscosity on patch.
virtual ~PhaseCompressibleTurbulenceModel()=default
Destructor.
static autoPtr< PhaseCompressibleTurbulenceModel > New(const alphaField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transportModel, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected turbulence model.
virtual tmp< volScalarField > mu() const
Return the laminar dynamic viscosity.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
virtual tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
virtual tmp< scalarField > mut(const label patchi) const
Return the turbulence dynamic viscosity on patch.
virtual tmp< volScalarField > mut() const
Return the turbulence dynamic viscosity.
Templated abstract base class for turbulence models.
const transportModel & transport() const
Access function to incompressible transport model.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
const volScalarField & rho() const
Return the density field.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux field.
A class for managing temporary objects.
Definition: tmp.H:65
const volVectorField & U() const
Access function to velocity field.
static const word propertiesName
Default name of the turbulence properties dictionary.
const surfaceScalarField & alphaRhoPhi() const
Access function to phase flux field.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
A class for handling words, derived from Foam::string.
Definition: word.H:68
Namespace for OpenFOAM.
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