Loading...
Searching...
No Matches
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-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
34namespace Foam
36namespace LESModels
37{
38
39// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40
41template<class BasicTurbulenceModel>
43(
44 const tmp<volTensorField>& gradU
45)
46{
47 this->nut_ = (flm_/fmm_)*sqr(this->delta())*mag(devSymm(gradU));
48 this->nut_.correctBoundaryConditions();
49 fv::options::New(this->mesh_).correct(this->nut_);
50
51 BasicTurbulenceModel::correctNut();
52}
53
54
55template<class BasicTurbulenceModel>
57{
59}
60
61
62// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63
64template<class BasicTurbulenceModel>
65dynamicLagrangian<BasicTurbulenceModel>::dynamicLagrangian
66(
67 const alphaField& alpha,
68 const rhoField& rho,
69 const volVectorField& U,
70 const surfaceScalarField& alphaRhoPhi,
72 const transportModel& transport,
73 const word& propertiesName,
74 const word& type
75)
76:
77 LESeddyViscosity<BasicTurbulenceModel>
78 (
79 type,
80 alpha,
81 rho,
82 U,
83 alphaRhoPhi,
84 phi,
85 transport,
86 propertiesName
87 ),
88
89 flm_
90 (
92 (
93 IOobject::groupName("flm", this->alphaRhoPhi_.group()),
94 this->runTime_.timeName(),
95 this->mesh_,
96 IOobject::MUST_READ,
97 IOobject::AUTO_WRITE
98 ),
99 this->mesh_
100 ),
101 fmm_
102 (
104 (
105 IOobject::groupName("fmm", this->alphaRhoPhi_.group()),
106 this->runTime_.timeName(),
107 this->mesh_,
108 IOobject::MUST_READ,
109 IOobject::AUTO_WRITE
110 ),
111 this->mesh_
112 ),
113 theta_
114 (
115 dimensioned<scalar>::getOrAddToDict
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);
134}
135
136
137// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
138
139template<class BasicTurbulenceModel>
141{
143 {
144 filter_.read(this->coeffDict());
145 theta_.readIfPresent(this->coeffDict());
146
147 return true;
149
150 return false;
151}
152
153
154template<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
171 tmp<volTensorField> tgradU(fvc::grad(U));
172 const volTensorField& gradU = tgradU();
173
174 volSymmTensorField S(devSymm(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// ************************************************************************* //
scalar delta
#define M(I)
fv::options & fvOptions
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Eddy viscosity LES SGS model base class.
virtual bool read()
Read model coefficients if they have changed.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
virtual void correct()
Correct Eddy-Viscosity and related properties.
void correctNut(const tmp< volTensorField > &gradU)
Update sub-grid eddy-viscosity.
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Read model coefficients if they have changed.
Abstract class for LES filters.
Definition LESfilter.H:54
Generic dimensioned Type class.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition fvMatrix.C:1094
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
Definition fvOptions.H:69
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
A class for managing temporary objects.
Definition tmp.H:75
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > Uf
word timeName
Definition getTimeIndex.H:3
Namespace for LES SGS models.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedSymmTensor sqr(const dimensionedVector &dv)
fvMatrix< scalar > fvScalarMatrix
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition bound.C:29
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
SymmTensor< Cmpt > devSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of the symmetric part of a SymmTensor.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
volScalarField & alpha
const vector L(dict.get< vector >("L"))