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