WALE.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) 2015-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 \*---------------------------------------------------------------------------*/
25 
26 #include "WALE.H"
27 #include "fvModels.H"
28 #include "fvConstraints.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace LESModels
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
39 template<class BasicMomentumTransportModel>
41 (
42  const volTensorField& gradU
43 ) const
44 {
45  return dev(symm(gradU & gradU));
46 }
47 
48 
49 template<class BasicMomentumTransportModel>
51 (
52  const volTensorField& gradU
53 ) const
54 {
55  volScalarField magSqrSd(magSqr(Sd(gradU)));
56 
57  return volScalarField::New
58  (
59  this->groupName("k"),
60  sqr(sqr(Cw_)*this->delta()/Ck_)*
61  (
62  pow3(magSqrSd)
63  /(
64  sqr
65  (
66  pow(magSqr(symm(gradU)), 5.0/2.0)
67  + pow(magSqrSd, 5.0/4.0)
68  )
70  (
71  "small",
72  dimensionSet(0, 0, -10, 0, 0),
73  small
74  )
75  )
76  )
77  );
78 }
79 
80 
81 template<class BasicMomentumTransportModel>
83 {
84  this->nut_ = Ck_*this->delta()*sqrt(this->k(fvc::grad(this->U_)));
85  this->nut_.correctBoundaryConditions();
86  fvConstraints::New(this->mesh_).constrain(this->nut_);
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
91 
92 template<class BasicMomentumTransportModel>
94 (
95  const alphaField& alpha,
96  const rhoField& rho,
97  const volVectorField& U,
98  const surfaceScalarField& alphaRhoPhi,
99  const surfaceScalarField& phi,
100  const viscosity& viscosity,
101  const word& type
102 )
103 :
104  LESeddyViscosity<BasicMomentumTransportModel>
105  (
106  type,
107  alpha,
108  rho,
109  U,
110  alphaRhoPhi,
111  phi,
112  viscosity
113  ),
114 
115  Ck_
116  (
117  dimensioned<scalar>::lookupOrAddToDict
118  (
119  "Ck",
120  this->coeffDict_,
121  0.094
122  )
123  ),
124 
125  Cw_
126  (
127  dimensioned<scalar>::lookupOrAddToDict
128  (
129  "Cw",
130  this->coeffDict_,
131  0.325
132  )
133  )
134 {
135  if (type == typeName)
136  {
137  this->printCoeffs(type);
138  }
139 }
140 
141 
142 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
143 
144 template<class BasicMomentumTransportModel>
146 {
148  {
149  Ck_.readIfPresent(this->coeffDict());
150  Cw_.readIfPresent(this->coeffDict());
151 
152  return true;
153  }
154  else
155  {
156  return false;
157  }
158 }
159 
160 
161 template<class BasicMomentumTransportModel>
163 {
164  volScalarField k(this->k(fvc::grad(this->U_)));
165 
166  return volScalarField::New
167  (
168  this->groupName("epsilon"),
169  this->Ce_*k*sqrt(k)/this->delta()
170  );
171 }
172 
173 
174 template<class BasicMomentumTransportModel>
176 {
178  correctNut();
179 }
180 
181 
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 
184 } // End namespace LESModels
185 } // End namespace Foam
186 
187 // ************************************************************************* //
label k
scalar delta
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,.
Eddy viscosity LES SGS model base class.
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: WALE.H:138
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: WALE.C:175
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: WALE.C:162
tmp< volSymmTensorField > Sd(const volTensorField &gradU) const
Return the deviatoric symmetric part of the square of the given.
Definition: WALE.C:41
WALE(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.
Definition: WALE.C:94
virtual void correctNut()
Update the SGS eddy-viscosity.
Definition: WALE.C:82
virtual bool read()
Read model coefficients if they have changed.
Definition: WALE.C:145
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Dimension set for the base types.
Definition: dimensionSet.H:125
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
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))
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488