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  this->printCoeffs(type);
126  }
127 }
128 
129 
130 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 
132 template<class BasicTurbulenceModel>
134 {
136  {
137  filter_.read(this->coeffDict());
138  theta_.readIfPresent(this->coeffDict());
139 
140  return true;
141  }
142  else
143  {
144  return false;
145  }
146 }
147 
148 
149 template<class BasicTurbulenceModel>
151 {
152  if (!this->turbulence_)
153  {
154  return;
155  }
156 
157  // Local references
158  const alphaField& alpha = this->alpha_;
159  const rhoField& rho = this->rho_;
160  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
161  const volVectorField& U = this->U_;
162 
164 
166  const volTensorField& gradU = tgradU();
167 
168  volSymmTensorField S(dev(symm(gradU)));
169  volScalarField magS(mag(S));
170 
171  volVectorField Uf(filter_(U));
173  volScalarField magSf(mag(Sf));
174 
175  volSymmTensorField L(dev(filter_(sqr(U)) - (sqr(filter_(U)))));
177  (
178  2.0*sqr(this->delta())*(filter_(magS*S) - 4.0*magSf*Sf)
179  );
180 
181  volScalarField invT
182  (
183  (1.0/(theta_.value()*this->delta()))*pow(flm_*fmm_, 1.0/8.0)
184  );
185 
186  volScalarField LM(L && M);
187 
188  fvScalarMatrix flmEqn
189  (
190  fvm::ddt(alpha, rho, flm_)
191  + fvm::div(alphaRhoPhi, flm_)
192  ==
193  rho*invT*LM
194  - fvm::Sp(rho*invT, flm_)
195  );
196 
197  flmEqn.relax();
198  flmEqn.solve();
199 
200  bound(flm_, flm0_);
201 
202  volScalarField MM(M && M);
203 
204  fvScalarMatrix fmmEqn
205  (
206  fvm::ddt(alpha, rho, fmm_)
207  + fvm::div(alphaRhoPhi, fmm_)
208  ==
209  rho*invT*MM
210  - fvm::Sp(rho*invT, fmm_)
211  );
212 
213  fmmEqn.relax();
214  fmmEqn.solve();
215 
216  bound(fmm_, fmm0_);
217 
218  correctNut(gradU);
219 }
220 
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 } // End namespace LESModels
225 } // End namespace Foam
226 
227 // ************************************************************************* //
Foam::LESModels::dynamicLagrangian::correctNut
virtual void correctNut()
Definition: dynamicLagrangian.C:49
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:82
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Foam::LESModels::dynamicLagrangian::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: dynamicLagrangian.H:101
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::LESModels::dynamicLagrangian::read
virtual bool read()
Read model coefficients if they have changed.
Definition: dynamicLagrangian.C:133
U
U
Definition: pEqn.H:46
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:519
Foam::LESModels::dynamicLagrangian::dynamicLagrangian
dynamicLagrangian(const dynamicLagrangian &)
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::LESModels::LESeddyViscosity
Eddy viscosity LES SGS model base class.
Definition: LESeddyViscosity.H:55
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::LESModels::dynamicLagrangian::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: dynamicLagrangian.H:100
Foam::LESModels::dynamicLagrangian::correct
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: dynamicLagrangian.C:150
rho
rho
Definition: pEqn.H:3
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
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:41
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::eddyViscosity< LESModel< BasicTurbulenceModel > >::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: eddyViscosity.C:132
Foam::LESModels::dynamicLagrangian::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: dynamicLagrangian.H:102
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Uf
Uf
Definition: pEqn.H:78
dynamicLagrangian.H
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104