dynamicLagrangian.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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2020 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 "dynamicLagrangian.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>
43 (
44  const tmp<volTensorField>& gradU
45 )
46 {
47  this->nut_ = (flm_/fmm_)*sqr(this->delta())*mag(dev(symm(gradU)));
48  this->nut_.correctBoundaryConditions();
49  fv::options::New(this->mesh_).correct(this->nut_);
50 
51  BasicTurbulenceModel::correctNut();
52 }
53 
54 
55 template<class BasicTurbulenceModel>
57 {
58  correctNut(fvc::grad(this->U_));
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63 
64 template<class BasicTurbulenceModel>
66 (
67  const alphaField& alpha,
68  const rhoField& rho,
69  const volVectorField& U,
70  const surfaceScalarField& alphaRhoPhi,
71  const surfaceScalarField& phi,
72  const transportModel& transport,
73  const word& propertiesName,
74  const word& type
75 )
76 :
78  (
79  type,
80  alpha,
81  rho,
82  U,
83  alphaRhoPhi,
84  phi,
85  transport,
86  propertiesName
87  ),
88 
89  flm_
90  (
91  IOobject
92  (
93  IOobject::groupName("flm", this->alphaRhoPhi_.group()),
94  this->runTime_.timeName(),
95  this->mesh_,
98  ),
99  this->mesh_
100  ),
101  fmm_
102  (
103  IOobject
104  (
105  IOobject::groupName("fmm", this->alphaRhoPhi_.group()),
106  this->runTime_.timeName(),
107  this->mesh_,
110  ),
111  this->mesh_
112  ),
113  theta_
114  (
116  (
117  "theta",
118  this->coeffDict_,
119  1.5
120  )
121  ),
122 
123  simpleFilter_(U.mesh()),
124  filterPtr_(LESfilter::New(U.mesh(), this->coeffDict())),
125  filter_(filterPtr_()),
126 
127  flm0_("flm0", flm_.dimensions(), Zero),
128  fmm0_("fmm0", fmm_.dimensions(), VSMALL)
129 {
130  if (type == typeName)
131  {
132  this->printCoeffs(type);
133  }
134 }
135 
136 
137 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
138 
139 template<class BasicTurbulenceModel>
141 {
143  {
144  filter_.read(this->coeffDict());
145  theta_.readIfPresent(this->coeffDict());
146 
147  return true;
148  }
149 
150  return false;
151 }
152 
153 
154 template<class BasicTurbulenceModel>
156 {
157  if (!this->turbulence_)
158  {
159  return;
160  }
161 
162  // Local references
163  const alphaField& alpha = this->alpha_;
164  const rhoField& rho = this->rho_;
165  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
166  const volVectorField& U = this->U_;
167  fv::options& fvOptions(fv::options::New(this->mesh_));
168 
170 
172  const volTensorField& gradU = tgradU();
173 
174  volSymmTensorField S(dev(symm(gradU)));
175  volScalarField magS(mag(S));
176 
177  volVectorField Uf(filter_(U));
179  volScalarField magSf(mag(Sf));
180 
181  volSymmTensorField L(dev(filter_(sqr(U)) - (sqr(filter_(U)))));
183  (
184  2.0*sqr(this->delta())*(filter_(magS*S) - 4.0*magSf*Sf)
185  );
186 
187  volScalarField invT
188  (
189  alpha*rho*(1.0/(theta_.value()*this->delta()))*pow(flm_*fmm_, 1.0/8.0)
190  );
191 
192  volScalarField LM(L && M);
193 
194  fvScalarMatrix flmEqn
195  (
196  fvm::ddt(alpha, rho, flm_)
197  + fvm::div(alphaRhoPhi, flm_)
198  ==
199  invT*LM
200  - fvm::Sp(invT, flm_)
201  + fvOptions(alpha, rho, flm_)
202  );
203 
204  flmEqn.relax();
205  fvOptions.constrain(flmEqn);
206  flmEqn.solve();
207  fvOptions.correct(flm_);
208  bound(flm_, flm0_);
209 
210  volScalarField MM(M && M);
211 
212  fvScalarMatrix fmmEqn
213  (
214  fvm::ddt(alpha, rho, fmm_)
215  + fvm::div(alphaRhoPhi, fmm_)
216  ==
217  invT*MM
218  - fvm::Sp(invT, fmm_)
219  + fvOptions(alpha, rho, fmm_)
220  );
221 
222  fmmEqn.relax();
223  fvOptions.constrain(fmmEqn);
224  fmmEqn.solve();
225  fvOptions.correct(fmm_);
226  bound(fmm_, fmm0_);
227 
228  correctNut(gradU);
229 }
230 
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 } // End namespace LESModels
235 } // End namespace Foam
236 
237 // ************************************************************************* //
Foam::LESModels::dynamicLagrangian::correctNut
virtual void correctNut()
Definition: dynamicLagrangian.C:56
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
L
const vector L(dict.get< vector >("L"))
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
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::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
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::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
Foam::LESModels::dynamicLagrangian::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: dynamicLagrangian.H:106
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:33
rho
rho
Definition: readInitialConditions.H:88
Foam::LESModels::dynamicLagrangian::read
virtual bool read()
Read model coefficients if they have changed.
Definition: dynamicLagrangian.C:140
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1183
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
Foam::LESModels::LESeddyViscosity
Eddy viscosity LES SGS model base class.
Definition: LESeddyViscosity.H:58
Foam::LESModels::dynamicLagrangian
Dynamic SGS model with Lagrangian averaging.
Definition: dynamicLagrangian.H:65
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::LESModels::dynamicLagrangian::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: dynamicLagrangian.H:105
Foam::LESModels::dynamicLagrangian::correct
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: dynamicLagrangian.C:155
U
U
Definition: pEqn.H:72
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::LESfilter::New
static autoPtr< LESfilter > New(const fvMesh &, const dictionary &, const word &filterDictName="filter")
Return a reference to the selected LES filter.
Definition: LESfilter.C:44
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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
M
#define M(I)
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:319
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::LESModels::dynamicLagrangian::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: dynamicLagrangian.H:107
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
dynamicLagrangian.H
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106