dynamicLagrangian.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-2018 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 "dynamicLagrangian.H"
27 #include "fvOptions.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace LESModels
34 {
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
38 template<class BasicTurbulenceModel>
40 (
41  const tmp<volTensorField>& gradU
42 )
43 {
44  this->nut_ = (flm_/fmm_)*sqr(this->delta())*mag(dev(symm(gradU)));
45  this->nut_.correctBoundaryConditions();
46  fv::options::New(this->mesh_).correct(this->nut_);
47 
48  BasicTurbulenceModel::correctNut();
49 }
50 
51 
52 template<class BasicTurbulenceModel>
54 {
55  correctNut(fvc::grad(this->U_));
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
61 template<class BasicTurbulenceModel>
63 (
64  const alphaField& alpha,
65  const rhoField& rho,
66  const volVectorField& U,
67  const surfaceScalarField& alphaRhoPhi,
68  const surfaceScalarField& phi,
69  const transportModel& transport,
70  const word& propertiesName,
71  const word& type
72 )
73 :
75  (
76  type,
77  alpha,
78  rho,
79  U,
80  alphaRhoPhi,
81  phi,
82  transport,
83  propertiesName
84  ),
85 
86  flm_
87  (
88  IOobject
89  (
90  IOobject::groupName("flm", this->alphaRhoPhi_.group()),
91  this->runTime_.timeName(),
92  this->mesh_,
95  ),
96  this->mesh_
97  ),
98  fmm_
99  (
100  IOobject
101  (
102  IOobject::groupName("fmm", this->alphaRhoPhi_.group()),
103  this->runTime_.timeName(),
104  this->mesh_,
107  ),
108  this->mesh_
109  ),
110  theta_
111  (
113  (
114  "theta",
115  this->coeffDict_,
116  1.5
117  )
118  ),
119 
120  simpleFilter_(U.mesh()),
121  filterPtr_(LESfilter::New(U.mesh(), this->coeffDict())),
122  filter_(filterPtr_()),
123 
124  flm0_("flm0", flm_.dimensions(), 0.0),
125  fmm0_("fmm0", fmm_.dimensions(), vSmall)
126 {
127  if (type == typeName)
128  {
129  this->printCoeffs(type);
130  }
131 }
132 
133 
134 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
135 
136 template<class BasicTurbulenceModel>
138 {
140  {
141  filter_.read(this->coeffDict());
142  theta_.readIfPresent(this->coeffDict());
143 
144  return true;
145  }
146  else
147  {
148  return false;
149  }
150 }
151 
152 
153 template<class BasicTurbulenceModel>
155 {
156  if (!this->turbulence_)
157  {
158  return;
159  }
160 
161  // Local references
162  const alphaField& alpha = this->alpha_;
163  const rhoField& rho = this->rho_;
164  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
165  const volVectorField& U = this->U_;
166  fv::options& fvOptions(fv::options::New(this->mesh_));
167 
169 
170  tmp<volTensorField> tgradU(fvc::grad(U));
171  const volTensorField& gradU = tgradU();
172 
173  volSymmTensorField S(dev(symm(gradU)));
174  volScalarField magS(mag(S));
175 
176  volVectorField Uf(filter_(U));
178  volScalarField magSf(mag(Sf));
179 
180  volSymmTensorField L(dev(filter_(sqr(U)) - (sqr(filter_(U)))));
182  (
183  2.0*sqr(this->delta())*(filter_(magS*S) - 4.0*magSf*Sf)
184  );
185 
186  volScalarField invT
187  (
188  alpha*rho*(1.0/(theta_.value()*this->delta()))*pow(flm_*fmm_, 1.0/8.0)
189  );
190 
191  volScalarField LM(L && M);
192 
193  fvScalarMatrix flmEqn
194  (
195  fvm::ddt(alpha, rho, flm_)
196  + fvm::div(alphaRhoPhi, flm_)
197  ==
198  invT*LM
199  - fvm::Sp(invT, flm_)
200  + fvOptions(alpha, rho, flm_)
201  );
202 
203  flmEqn.relax();
204  fvOptions.constrain(flmEqn);
205  flmEqn.solve();
206  fvOptions.correct(flm_);
207  bound(flm_, flm0_);
208 
209  volScalarField MM(M && M);
210 
211  fvScalarMatrix fmmEqn
212  (
213  fvm::ddt(alpha, rho, fmm_)
214  + fvm::div(alphaRhoPhi, fmm_)
215  ==
216  invT*MM
217  - fvm::Sp(invT, fmm_)
218  + fvOptions(alpha, rho, fmm_)
219  );
220 
221  fmmEqn.relax();
222  fvOptions.constrain(fmmEqn);
223  fmmEqn.solve();
224  fvOptions.correct(fmm_);
225  bound(fmm_, fmm0_);
226 
227  correctNut(gradU);
228 }
229 
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 } // End namespace LESModels
234 } // End namespace Foam
235 
236 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
scalar delta
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fv::options & fvOptions
surfaceScalarField & phi
Dynamic SGS model with Lagrangian averaging.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Finite-volume options.
Definition: fvOptions.H:52
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
autoPtr< surfaceVectorField > Uf
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
virtual bool read()
Read model coefficients if they have changed.
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:57
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:523
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
BasicTurbulenceModel::rhoField rhoField
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
virtual void correct()
Correct Eddy-Viscosity and related properties.
BasicTurbulenceModel::alphaField alphaField
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
static autoPtr< LESfilter > New(const fvMesh &, const dictionary &, const word &filterDictName="filter")
Return a reference to the selected LES filter.
Definition: LESfilter.C:41
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
BasicTurbulenceModel::transportModel transportModel
#define M(I)
Namespace for OpenFOAM.