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-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::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
104 
107 
108 
109  // Protected Member Functions
110 
112  (
115  );
116 
117  void setCdAv
118  (
121  );
122 
124 
125  //- Bound epsilon and return Cmu*sqr(k) for nut
127 
128  //- Correct the eddy-viscosity nut
129  virtual void correctNut();
130 
132  (
133  const volScalarField::Internal& magU,
134  const volScalarField::Internal& magU3
135  ) const;
136 
138  (
139  const volScalarField::Internal& magU,
140  const volScalarField::Internal& magU3
141  ) const;
142 
143 
144 public:
145 
146  typedef typename BasicMomentumTransportModel::alphaField alphaField;
147  typedef typename BasicMomentumTransportModel::rhoField rhoField;
148 
149 
150  //- Runtime type information
151  TypeName("kEpsilonLopesdaCosta");
152 
153 
154  // Constructors
155 
156  //- Construct from components
158  (
159  const alphaField& alpha,
160  const rhoField& rho,
161  const volVectorField& U,
162  const surfaceScalarField& alphaRhoPhi,
163  const surfaceScalarField& phi,
164  const viscosity& viscosity,
165  const word& type = typeName
166  );
167 
168  //- Disallow default bitwise copy construction
170 
171 
172  //- Destructor
173  virtual ~kEpsilonLopesdaCosta()
174  {}
175 
176 
177  // Member Functions
178 
179  //- Re-read model coefficients if they have changed
180  virtual bool read();
181 
182  //- Return the effective diffusivity for k
183  tmp<volScalarField> DkEff() const
184  {
185  return volScalarField::New
186  (
187  "DkEff",
188  (this->nut_/sigmak_ + this->nu())
189  );
190  }
191 
192  //- Return the effective diffusivity for epsilon
194  {
195  return volScalarField::New
196  (
197  "DepsilonEff",
198  (this->nut_/sigmaEps_ + this->nu())
199  );
200  }
201 
202  //- Return the turbulence kinetic energy
203  virtual tmp<volScalarField> k() const
204  {
205  return k_;
206  }
207 
208  //- Return the turbulence kinetic energy dissipation rate
209  virtual tmp<volScalarField> epsilon() const
210  {
211  return epsilon_;
212  }
213 
214  //- Return the turbulence specific dissipation rate
215  virtual tmp<volScalarField> omega() const
216  {
217  return volScalarField::New
218  (
219  "omega",
220  epsilon_/(Cmu_*k_)
221  );
222  }
223 
224  //- Solve the turbulence equations and correct the turbulence viscosity
225  virtual void correct();
226 
227 
228  // Member Operators
229 
230  //- Disallow default bitwise assignment
231  void operator=(const kEpsilonLopesdaCosta&) = delete;
232 };
233 
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 } // End namespace RASModels
238 } // End namespace Foam
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242 #ifdef NoRepository
243  #include "kEpsilonLopesdaCosta.C"
244 #endif
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 #endif
249 
250 // ************************************************************************* //
Graphite solid properties.
Definition: C.H:51
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
Variant of the standard k-epsilon turbulence model with additional source terms to handle the changes...
BasicMomentumTransportModel::alphaField alphaField
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
virtual tmp< fvScalarMatrix > kSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
TypeName("kEpsilonLopesdaCosta")
Runtime type information.
kEpsilonLopesdaCosta(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
void operator=(const kEpsilonLopesdaCosta &)=delete
Disallow default bitwise assignment.
void setCdAv(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
virtual void correctNut()
Correct the eddy-viscosity nut.
virtual tmp< fvScalarMatrix > epsilonSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
void setPorosityCoefficient(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
virtual bool read()
Re-read model coefficients if they have changed.
BasicMomentumTransportModel::rhoField rhoField
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Variant of the power law porosity model with spatially varying drag coefficient.
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
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.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488