LESeddyViscosity.C
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-2024 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 \*---------------------------------------------------------------------------*/
25 
26 #include "LESeddyViscosity.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace LESModels
33 {
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
37 template<class BasicMomentumTransportModel>
39 (
40  const word& type,
41  const alphaField& alpha,
42  const rhoField& rho,
43  const volVectorField& U,
44  const surfaceScalarField& alphaRhoPhi,
45  const surfaceScalarField& phi,
46  const viscosity& viscosity
47 )
48 :
49  eddyViscosity<LESModel<BasicMomentumTransportModel>>
50  (
51  type,
52  alpha,
53  rho,
54  U,
55  alphaRhoPhi,
56  phi,
57  viscosity
58  ),
59 
60  Ck_("Ck", this->coeffDict(), 0.094),
61  Ce_("Ce", this->coeffDict(), 1.048)
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
67 template<class BasicMomentumTransportModel>
69 {
71  {
72  Ck_.readIfPresent(this->coeffDict());
73  Ce_.readIfPresent(this->coeffDict());
74 
75  return true;
76  }
77  else
78  {
79  return false;
80  }
81 }
82 
83 
84 template<class BasicMomentumTransportModel>
87 {
88  tmp<volScalarField> tk(this->k());
89 
90  return volScalarField::New
91  (
92  this->groupName("epsilon"),
93  Ce_*tk()*sqrt(tk())/this->delta()
94  );
95 }
96 
97 
98 template<class BasicMomentumTransportModel>
101 {
102  tmp<volScalarField> tk(this->k());
103 
104  return volScalarField::New
105  (
106  this->groupName("omega"),
107  (Ce_/Ck_)*sqrt(tk())/this->delta()
108  );
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
113 
114 } // End namespace LESModels
115 } // End namespace Foam
116 
117 // ************************************************************************* //
label k
scalar delta
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
Templated abstract base class for LES SGS models.
Definition: LESModel.H:63
BasicMomentumTransportModel::alphaField alphaField
LESeddyViscosity(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
virtual bool read()
Read model coefficients if they have changed.
BasicMomentumTransportModel::rhoField rhoField
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
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.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488