LESeddyViscosity.C
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-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "LESeddyViscosity.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace LESModels
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 template<class BasicTurbulenceModel>
42 LESeddyViscosity<BasicTurbulenceModel>::LESeddyViscosity
43 (
44  const word& type,
45  const alphaField& alpha,
46  const rhoField& rho,
47  const volVectorField& U,
48  const surfaceScalarField& alphaRhoPhi,
49  const surfaceScalarField& phi,
50  const transportModel& transport,
51  const word& propertiesName
52 )
53 :
55  (
56  type,
57  alpha,
58  rho,
59  U,
60  alphaRhoPhi,
61  phi,
62  transport,
63  propertiesName
64  ),
65 
66  Ce_
67  (
69  (
70  "Ce",
71  this->coeffDict_,
72  1.048
73  )
74  )
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
80 template<class BasicTurbulenceModel>
82 {
84  {
85  Ce_.readIfPresent(this->coeffDict());
86 
87  return true;
88  }
89 
90  return false;
91 }
92 
93 
94 template<class BasicTurbulenceModel>
96 {
97  tmp<volScalarField> tk(this->k());
98 
99  tmp<volScalarField> tepsilon
100  (
101  new volScalarField
102  (
103  IOobject
104  (
105  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
106  this->runTime_.timeName(),
107  this->mesh_,
110  ),
111  Ce_*tk()*sqrt(tk())/this->delta(),
113  )
114  );
115  volScalarField& epsilon = tepsilon.ref();
116  epsilon.correctBoundaryConditions();
117 
118  return tepsilon;
119 }
120 
121 
122 template<class BasicTurbulenceModel>
124 {
125  tmp<volScalarField> tk(this->k());
126  tmp<volScalarField> tepsilon(this->epsilon());
127 
128  auto tomega = tmp<volScalarField>::New
129  (
130  IOobject
131  (
132  IOobject::groupName("omega", this->alphaRhoPhi_.group()),
133  this->runTime_.timeName(),
134  this->mesh_
135  ),
136  tepsilon()/(0.09*tk())
137  );
138  auto& omega = tomega.ref();
139  omega.correctBoundaryConditions();
140 
141  return tomega;
142 }
143 
144 
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 
147 } // End namespace LESModels
148 } // End namespace Foam
149 
150 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
LESeddyViscosity.H
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::LESModels::LESeddyViscosity::read
virtual bool read()
Read model coefficients if they have changed.
Definition: LESeddyViscosity.C:81
Foam::LESModels::LESeddyViscosity::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: LESeddyViscosity.H:82
zeroGradientFvPatchField.H
rho
rho
Definition: readInitialConditions.H:88
Foam::zeroGradientFvPatchField
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Definition: zeroGradientFvPatchField.H:64
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::LESModels::LESeddyViscosity::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: LESeddyViscosity.H:81
Foam::LESModels::LESeddyViscosity::omega
virtual tmp< volScalarField > omega() const
Return sub-grid specific dissipation rate.
Definition: LESeddyViscosity.C:123
U
U
Definition: pEqn.H:72
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::LESModels::LESeddyViscosity::epsilon
virtual tmp< volScalarField > epsilon() const
Return sub-grid dissipation rate.
Definition: LESeddyViscosity.C:95
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:55
Foam::LESModel
Templated abstract base class for LES SGS models.
Definition: LESModel.H:58
epsilon
scalar epsilon
Definition: evaluateNearWall.H:7
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::LESModels::LESeddyViscosity::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: LESeddyViscosity.H:80