DESModel.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) 2012-2015 OpenFOAM Foundation
9  Copyright (C) 2022 Upstream CFD GmbH
10  Copyright (C) 2022 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
30 #include "DESModel.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace LESModels
37 {
38 
39 // * * * * * * * * * * * * * * * Constructor * * * * * * * * * * * * * * * * //
40 
41 template<class BasicTurbulenceModel>
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  Ctrans_(dimless, Zero)
66 {
67  // Note: Ctrans is model-specific and initialised in derived classes
68 }
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
73 template<class BasicTurbulenceModel>
75 {
77  {
78  Ctrans_.readIfPresent(this->coeffDict());
79 
80  return true;
81  }
82 
83  return false;
84 }
85 
86 
87 template<class BasicTurbulenceModel>
89 (
90  const volTensorField& gradU,
91  const dimensionedScalar& coeff
92 )
93 {
94  // Limiter
95  const dimensionedScalar eps0(dimless, SMALL);
96  const dimensionedScalar eps2(dimless/sqr(dimTime), SMALL);
97  const dimensionedScalar eps4(dimless/pow4(dimTime), SMALL);
98  const dimensionedScalar max2(dimless/sqr(dimTime), GREAT);
99  const dimensionedTensor maxTen2
100  (
103  );
104  const dimensionedTensor minTen2
105  (
108  );
109 
110  const volTensorField G(max(min(gradU.T() & gradU, maxTen2), minTen2));
111 
112  // Tensor invariants
113  const volScalarField I1(tr(G));
114  const volScalarField I2(0.5*(sqr(I1) - tr(G & G)));
115  tmp<volScalarField> tI3 = det(G);
116 
117  const volScalarField alpha1(max(sqr(I1)/9.0 - I2/3.0, eps4));
118 
119  tmp<volScalarField> talpha2 =
120  pow3(min(I1, max2))/27.0 - I1*I2/6.0 + 0.5*tI3;
121 
122  const volScalarField alpha3
123  (
124  1.0/3.0
125  *acos
126  (
127  max
128  (
129  scalar(-1) + eps0,
130  min(scalar(1) - eps0, talpha2/pow(alpha1, 3.0/2.0))
131  )
132  )
133  );
134 
135  const scalar piBy3 = constant::mathematical::pi/3.0;
136  const volScalarField sigma1
137  (
138  sqrt(max(I1/3.0 + 2.0*sqrt(alpha1)*cos(alpha3), eps2))
139  );
140  const volScalarField sigma2
141  (
142  sqrt(max(I1/3.0 - 2.0*sqrt(alpha1)*cos(piBy3 + alpha3), eps2))
143  );
144  const volScalarField sigma3
145  (
146  sqrt(max(I1/3.0 - 2.0*sqrt(alpha1)*cos(piBy3 - alpha3), eps2))
147  );
148 
149  return
150  coeff
151  *sigma3
152  *(sigma1 - sigma2)
153  *(sigma2 - sigma3)
154  /max(sqr(sigma1), eps2);
155 }
156 
157 
158 template<class BasicTurbulenceModel>
160 (
161  const volTensorField& gradU
162 ) const
163 {
164  return Ssigma(gradU, Ctrans_);
165 }
166 
167 
168 template<class BasicTurbulenceModel>
170 {
172  (
173  IOobject
174  (
175  "fd",
176  this->mesh_.time().timeName(),
177  this->mesh_,
180  ),
181  this->mesh_,
183  );
184 }
185 
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 } // End namespace LESModels
190 } // End namespace Foam
191 
192 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
dimensionedScalar acos(const dimensionedScalar &ds)
const dimensionedScalar G
Newtonian constant of gravitation.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
BasicTurbulenceModel::rhoField rhoField
Definition: DESModel.H:82
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar det(const dimensionedSphericalTensor &dt)
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
dimensionedScalar cos(const dimensionedScalar &ds)
A class for handling words, derived from Foam::string.
Definition: word.H:63
static tmp< volScalarField > Ssigma(const volTensorField &gradU, const dimensionedScalar &coeff)
Return modified strain rate.
Definition: DESModel.C:82
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
constexpr scalar pi(M_PI)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
BasicTurbulenceModel::transportModel transportModel
Definition: DESModel.H:83
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pow4(const dimensionedScalar &ds)
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual tmp< volScalarField > fd() const
Return the shielding function.
Definition: DESModel.C:162
BasicTurbulenceModel::alphaField alphaField
Definition: DESModel.H:81
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
virtual bool read()
Re-read model coefficients if they have changed.
Definition: DESModel.C:67
Templated abstract base class for DES models.
Definition: DESModel.H:51
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const volScalarField & alpha1