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-2022 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 :
50  (
51  type,
52  alpha,
53  rho,
54  U,
55  alphaRhoPhi,
56  phi,
57  viscosity
58  ),
59 
60  Ck_
61  (
63  (
64  "Ck",
65  this->coeffDict_,
66  0.094
67  )
68  ),
69 
70  Ce_
71  (
73  (
74  "Ce",
75  this->coeffDict_,
76  1.048
77  )
78  )
79 {}
80 
81 
82 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
83 
84 template<class BasicMomentumTransportModel>
86 {
88  {
89  Ck_.readIfPresent(this->coeffDict());
90  Ce_.readIfPresent(this->coeffDict());
91 
92  return true;
93  }
94  else
95  {
96  return false;
97  }
98 }
99 
100 
101 template<class BasicMomentumTransportModel>
104 {
105  tmp<volScalarField> tk(this->k());
106 
107  return volScalarField::New
108  (
109  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
110  Ce_*tk()*sqrt(tk())/this->delta()
111  );
112 }
113 
114 
115 template<class BasicMomentumTransportModel>
118 {
119  tmp<volScalarField> tk(this->k());
120 
121  return volScalarField::New
122  (
123  IOobject::groupName("omega", this->alphaRhoPhi_.group()),
124  (Ce_/Ck_)*sqrt(tk())/this->delta()
125  );
126 }
127 
128 
129 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
130 
131 } // End namespace LESModels
132 } // End namespace Foam
133 
134 // ************************************************************************* //
scalar delta
BasicMomentumTransportModel::rhoField rhoField
Templated abstract base class for LES SGS models.
Definition: LESModel.H:60
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
bool read(const char *, int32_t &)
Definition: int32IO.C:85
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
Info<< "Reading field p\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, simple.dict());mesh.schemes().setFluxRequired(p.name());Info<< "Reading field pa\"<< endl;volScalarField pa(IOobject("pa", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Ua\"<< endl;volVectorField Ua(IOobject("Ua", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);# 65 "/home/ubuntu/OpenFOAM-10/applications/solvers/incompressible/adjointShapeOptimisationFoam/createFields.H" 2label paRefCell=0;scalar paRefValue=0.0;setRefCell(pa, simple.dict(), paRefCell, paRefValue);mesh.schemes().setFluxRequired(pa.name());autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
phi
Definition: correctPhi.H:3
virtual bool read()
Read model coefficients if they have changed.
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.
BasicMomentumTransportModel::alphaField alphaField
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.