kEpsilonLopesdaCosta.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) 2018-2020 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::kEpsilonLopesdaCosta
26 
27 Description
28  Variant of the standard k-epsilon turbulence model with additional source
29  terms to handle the changes in turbulence in porous regions represented by
30  the powerLawLopesdaCosta porosity model.
31 
32  Reference:
33  \verbatim
34  Costa, J. C. P. L. D. (2007).
35  Atmospheric flow over forested and non-forested complex terrain.
36  \endverbatim
37 
38  The default model coefficients are
39  \verbatim
40  kEpsilonLopesdaCostaCoeffs
41  {
42  Cmu 0.09;
43  C1 1.44;
44  C2 1.92;
45  sigmak 1.0;
46  sigmaEps 1.3;
47  }
48  \endverbatim
49 
50 See also
51  Foam::RASModels::kEpsilon
52  Foam::porosityModels::powerLawLopesdaCosta
53 
54 SourceFiles
55  kEpsilonLopesdaCosta.C
56 
57 \*---------------------------------------------------------------------------*/
58 
59 #ifndef kEpsilonLopesdaCosta_H
60 #define kEpsilonLopesdaCosta_H
61 
62 #include "RASModel.H"
63 #include "eddyViscosity.H"
64 #include "powerLawLopesdaCosta.H"
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 namespace Foam
69 {
70 namespace RASModels
71 {
72 
73 /*---------------------------------------------------------------------------*\
74  Class kEpsilonLopesdaCosta Declaration
75 \*---------------------------------------------------------------------------*/
76 
77 template<class BasicMomentumTransportModel>
79 :
80  public eddyViscosity<RASModel<BasicMomentumTransportModel>>
81 {
82 protected:
83 
84  // Protected data
85 
86  // Standard model coefficients
87 
93 
94  // Lopes da Costa porosity coefficients
95 
101 
102 
103  // Fields
107 
108 
109  // Protected Member Functions
110 
112  (
115  );
116 
117  void setCdSigma
118  (
121  );
122 
124 
125  virtual void correctNut();
126 
128  (
129  const volScalarField::Internal& magU,
130  const volScalarField::Internal& magU3
131  ) const;
132 
134  (
135  const volScalarField::Internal& magU,
136  const volScalarField::Internal& magU3
137  ) const;
138 
139 
140 public:
142  typedef typename BasicMomentumTransportModel::alphaField alphaField;
143  typedef typename BasicMomentumTransportModel::rhoField rhoField;
144  typedef typename BasicMomentumTransportModel::transportModel transportModel;
145 
146 
147  //- Runtime type information
148  TypeName("kEpsilonLopesdaCosta");
149 
150 
151  // Constructors
152 
153  //- Construct from components
155  (
156  const alphaField& alpha,
157  const rhoField& rho,
158  const volVectorField& U,
159  const surfaceScalarField& alphaRhoPhi,
160  const surfaceScalarField& phi,
161  const transportModel& transport,
162  const word& type = typeName
163  );
164 
165  //- Disallow default bitwise copy construction
167 
168 
169  //- Destructor
170  virtual ~kEpsilonLopesdaCosta()
171  {}
172 
173 
174  // Member Functions
175 
176  //- Re-read model coefficients if they have changed
177  virtual bool read();
178 
179  //- Return the effective diffusivity for k
180  tmp<volScalarField> DkEff() const
181  {
182  return volScalarField::New
183  (
184  "DkEff",
185  (this->nut_/sigmak_ + this->nu())
186  );
187  }
188 
189  //- Return the effective diffusivity for epsilon
191  {
192  return volScalarField::New
193  (
194  "DepsilonEff",
195  (this->nut_/sigmaEps_ + this->nu())
196  );
197  }
198 
199  //- Return the turbulence kinetic energy
200  virtual tmp<volScalarField> k() const
201  {
202  return k_;
203  }
204 
205  //- Return the turbulence kinetic energy dissipation rate
206  virtual tmp<volScalarField> epsilon() const
207  {
208  return epsilon_;
209  }
210 
211  //- Solve the turbulence equations and correct the turbulence viscosity
212  virtual void correct();
213 
214 
215  // Member Operators
216 
217  //- Disallow default bitwise assignment
218  void operator=(const kEpsilonLopesdaCosta&) = delete;
219 };
220 
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 } // End namespace RASModels
225 } // End namespace Foam
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 #ifdef NoRepository
230  #include "kEpsilonLopesdaCosta.C"
231 #endif
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 #endif
236 
237 // ************************************************************************* //
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
void setPorosityCoefficient(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
Graphite solid properties.
Definition: C.H:48
Variant of the power law porosity model with spatially varying drag coefficient.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual tmp< fvScalarMatrix > kSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
virtual bool read()
Re-read model coefficients if they have changed.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
phi
Definition: pEqn.H:104
virtual tmp< fvScalarMatrix > epsilonSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
BasicMomentumTransportModel::transportModel transportModel
A class for handling words, derived from string.
Definition: word.H:59
void setCdSigma(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
void operator=(const kEpsilonLopesdaCosta &)=delete
Disallow default bitwise assignment.
BasicMomentumTransportModel::rhoField rhoField
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
U
Definition: pEqn.H:72
kEpsilonLopesdaCosta(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Variant of the standard k-epsilon turbulence model with additional source terms to handle the changes...
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
BasicMomentumTransportModel::alphaField alphaField
Namespace for OpenFOAM.
TypeName("kEpsilonLopesdaCosta")
Runtime type information.