SpalartAllmarasDDES.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) 2011-2021 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 "SpalartAllmarasDDES.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace LESModels
33 {
34 
35 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
36 
37 template<class BasicMomentumTransportModel>
38 tmp<volScalarField::Internal>
39 SpalartAllmarasDDES<BasicMomentumTransportModel>::rd
40 (
41  const volScalarField::Internal& magGradU
42 ) const
43 {
45  (
46  modelName("rd"),
47  min
48  (
49  this->nuEff()()
50  /(
51  max
52  (
53  magGradU,
54  dimensionedScalar(magGradU.dimensions(), small)
55  )
56  *sqr(this->kappa_*this->y_())
57  ),
58  scalar(10)
59  )
60  );
61 }
62 
63 
64 template<class BasicMomentumTransportModel>
65 tmp<volScalarField::Internal>
66 SpalartAllmarasDDES<BasicMomentumTransportModel>::fd
67 (
68  const volScalarField::Internal& magGradU
69 ) const
70 {
72  (
73  modelName("fd"),
74  1 - tanh(pow3(8*rd(magGradU)))
75  );
76 }
77 
78 
79 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
80 
81 template<class BasicMomentumTransportModel>
82 tmp<volScalarField::Internal>
84 (
85  const volScalarField::Internal& chi,
86  const volScalarField::Internal& fv1,
87  const volTensorField::Internal& gradU
88 ) const
89 {
91  (
92  modelName("dTilda"),
93  max
94  (
95  this->y_
96  - fd(mag(gradU))
97  *max
98  (
99  this->y_() - this->CDES_*this->delta()(),
101  ),
103  )
104  );
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
109 
110 template<class BasicMomentumTransportModel>
112 (
113  const alphaField& alpha,
114  const rhoField& rho,
115  const volVectorField& U,
116  const surfaceScalarField& alphaRhoPhi,
117  const surfaceScalarField& phi,
118  const viscosity& viscosity,
119  const word& type
120 )
121 :
123  (
124  alpha,
125  rho,
126  U,
127  alphaRhoPhi,
128  phi,
129  viscosity
130  )
131 {}
132 
133 
134 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
135 
136 } // End namespace LESModels
137 } // End namespace Foam
138 
139 // ************************************************************************* //
scalar delta
dimensionedScalar tanh(const dimensionedScalar &ds)
BasicMomentumTransportModel::rhoField rhoField
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual tmp< volScalarField::Internal > dTilda(const volScalarField::Internal &chi, const volScalarField::Internal &fv1, const volTensorField::Internal &gradU) const
Length scale.
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
U
Definition: pEqn.H:72
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
const dimensionSet dimLength
SpalartAllmarasDDES(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.
A class for handling words, derived from string.
Definition: word.H:59
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))
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
phi
Definition: correctPhi.H:3
Abstract base class for all fluid physical properties.
Definition: viscosity.H:49
dimensionedScalar pow3(const dimensionedScalar &ds)
BasicMomentumTransportModel::alphaField alphaField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic 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
dimensioned< scalar > mag(const dimensioned< Type > &)
Namespace for OpenFOAM.