WALE.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) 2015-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2021 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 "WALE.H"
30 #include "fvOptions.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace LESModels
37 {
38 
39 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40 
41 template<class BasicTurbulenceModel>
42 tmp<volSymmTensorField> WALE<BasicTurbulenceModel>::Sd
43 (
44  const volTensorField& gradU
45 ) const
46 {
47  return dev(symm(gradU & gradU));
48 }
49 
50 
51 template<class BasicTurbulenceModel>
53 (
54  const volTensorField& gradU
55 ) const
56 {
57  volScalarField magSqrSd(magSqr(Sd(gradU)));
58 
59  return tmp<volScalarField>
60  (
61  new volScalarField
62  (
63  IOobject
64  (
65  IOobject::groupName("k", this->alphaRhoPhi_.group()),
66  this->runTime_.timeName(),
67  this->mesh_
68  ),
69  sqr(sqr(Cw_)*this->delta()/Ck_)*
70  (
71  pow3(magSqrSd)
72  /(
73  sqr
74  (
75  pow(magSqr(symm(gradU)), 5.0/2.0)
76  + pow(magSqrSd, 5.0/4.0)
77  )
79  (
80  "SMALL",
81  dimensionSet(0, 0, -10, 0, 0),
82  SMALL
83  )
84  )
85  )
86  )
87  );
88 }
89 
90 
91 template<class BasicTurbulenceModel>
93 {
94  this->nut_ = Ck_*this->delta()*sqrt(this->k(fvc::grad(this->U_)));
95  this->nut_.correctBoundaryConditions();
96  fv::options::New(this->mesh_).correct(this->nut_);
97 
98  BasicTurbulenceModel::correctNut();
99 }
100 
101 
102 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
103 
104 template<class BasicTurbulenceModel>
106 (
107  const alphaField& alpha,
108  const rhoField& rho,
109  const volVectorField& U,
110  const surfaceScalarField& alphaRhoPhi,
111  const surfaceScalarField& phi,
112  const transportModel& transport,
113  const word& propertiesName,
114  const word& type
115 )
116 :
118  (
119  type,
120  alpha,
121  rho,
122  U,
123  alphaRhoPhi,
124  phi,
125  transport,
126  propertiesName
127  ),
128 
129  Ck_
130  (
132  (
133  "Ck",
134  this->coeffDict_,
135  0.094
136  )
137  ),
138 
139  Cw_
140  (
142  (
143  "Cw",
144  this->coeffDict_,
145  0.325
146  )
147  )
148 {
149  if (type == typeName)
150  {
151  this->printCoeffs(type);
152  }
153 }
154 
155 
156 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
157 
158 template<class BasicTurbulenceModel>
160 {
162  {
163  Ck_.readIfPresent(this->coeffDict());
164  Cw_.readIfPresent(this->coeffDict());
165 
166  return true;
167  }
168 
169  return false;
170 }
171 
172 
173 template<class BasicTurbulenceModel>
175 {
177  correctNut();
178 }
179 
180 
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 
183 } // End namespace LESModels
184 } // End namespace Foam
185 
186 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:355
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::LESModels::WALE::correctNut
virtual void correctNut()
Update the SGS eddy-viscosity.
Definition: WALE.C:92
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:103
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:108
Foam::LESModels::LESeddyViscosity::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: LESeddyViscosity.H:75
rho
rho
Definition: readInitialConditions.H:88
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::LESModels::WALE::Sd
tmp< volSymmTensorField > Sd(const volTensorField &gradU) const
Return the deviatoric symmetric part of the square of the given.
Definition: WALE.C:43
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::LESModels::WALE::k
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: WALE.H:151
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::LESModels::LESeddyViscosity
Eddy viscosity LES SGS model base class.
Definition: LESeddyViscosity.H:58
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::LESModels::LESeddyViscosity::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: LESeddyViscosity.H:74
U
U
Definition: pEqn.H:72
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::LESModels::WALE::read
virtual bool read()
Read model coefficients if they have changed.
Definition: WALE.C:159
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::LESModels::WALE
The Wall-adapting local eddy-viscosity (WALE) SGS model.
Definition: WALE.H:79
Foam::LESModels::WALE::correct
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: WALE.C:174
WALE.H
Foam::eddyViscosity< LESModel< BasicTurbulenceModel > >::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: eddyViscosity.C:134
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::GeometricField< tensor, fvPatchField, volMesh >
Foam::LESModels::LESeddyViscosity::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: LESeddyViscosity.H:73
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106