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 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 BasicTurbulenceModel>
79 :
80  public eddyViscosity<RASModel<BasicTurbulenceModel>>
81 {
82  // Private Member Functions
83 
84  // Disallow default bitwise copy construct and assignment
86  void operator=(const kEpsilonLopesdaCosta&);
87 
88 
89 protected:
90 
91  // Protected data
92 
93  // Standard model coefficients
94 
100 
101  // Lopes da Costa porosity coefficients
108 
109 
110  // Fields
114 
115 
116  // Protected Member Functions
117 
119  (
122  );
123 
124  void setCdSigma
125  (
128  );
129 
131 
132  virtual void correctNut();
133 
135  (
136  const volScalarField::Internal& magU,
137  const volScalarField::Internal& magU3
138  ) const;
139 
141  (
142  const volScalarField::Internal& magU,
143  const volScalarField::Internal& magU3
144  ) const;
145 
146 
147 public:
149  typedef typename BasicTurbulenceModel::alphaField alphaField;
150  typedef typename BasicTurbulenceModel::rhoField rhoField;
151  typedef typename BasicTurbulenceModel::transportModel transportModel;
152 
153 
154  //- Runtime type information
155  TypeName("kEpsilonLopesdaCosta");
156 
157 
158  // Constructors
159 
160  //- Construct from components
162  (
163  const alphaField& alpha,
164  const rhoField& rho,
165  const volVectorField& U,
166  const surfaceScalarField& alphaRhoPhi,
167  const surfaceScalarField& phi,
168  const transportModel& transport,
169  const word& propertiesName = turbulenceModel::propertiesName,
170  const word& type = typeName
171  );
172 
173 
174  //- Destructor
175  virtual ~kEpsilonLopesdaCosta()
176  {}
177 
178 
179  // Member Functions
180 
181  //- Re-read model coefficients if they have changed
182  virtual bool read();
183 
184  //- Return the effective diffusivity for k
185  tmp<volScalarField> DkEff() const
186  {
187  return tmp<volScalarField>
188  (
189  new volScalarField
190  (
191  "DkEff",
192  (this->nut_/sigmak_ + this->nu())
193  )
194  );
195  }
196 
197  //- Return the effective diffusivity for epsilon
199  {
200  return tmp<volScalarField>
201  (
202  new volScalarField
203  (
204  "DepsilonEff",
205  (this->nut_/sigmaEps_ + this->nu())
206  )
207  );
208  }
209 
210  //- Return the turbulence kinetic energy
211  virtual tmp<volScalarField> k() const
212  {
213  return k_;
214  }
215 
216  //- Return the turbulence kinetic energy dissipation rate
217  virtual tmp<volScalarField> epsilon() const
218  {
219  return epsilon_;
220  }
221 
222  //- Solve the turbulence equations and correct the turbulence viscosity
223  virtual void correct();
224 };
225 
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 } // End namespace RASModels
230 } // End namespace Foam
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 #ifdef NoRepository
235  #include "kEpsilonLopesdaCosta.C"
236 #endif
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 #endif
241 
242 // ************************************************************************* //
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
surfaceScalarField & phi
Variant of the power law porosity model with spatially varying drag coefficient.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
BasicTurbulenceModel::rhoField rhoField
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.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual tmp< fvScalarMatrix > epsilonSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
void setCdSigma(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
U
Definition: pEqn.H:72
Variant of the standard k-epsilon turbulence model with additional source terms to handle the changes...
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
BasicTurbulenceModel::alphaField alphaField
volScalarField & nu
Namespace for OpenFOAM.
TypeName("kEpsilonLopesdaCosta")
Runtime type information.