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-2024 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 "fvModels.H"
28 #include "fvConstraints.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace LESModels
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
39 template<class BasicMomentumTransportModel>
41 (
42  const tmp<volTensorField>& gradU
43 )
44 {
45  this->nut_ = (flm_/fmm_)*sqr(this->delta())*mag(dev(symm(gradU)));
46  this->nut_.correctBoundaryConditions();
47  fvConstraints::New(this->mesh_).constrain(this->nut_);
48 }
49 
50 
51 template<class BasicMomentumTransportModel>
53 {
54  correctNut(fvc::grad(this->U_));
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
60 template<class BasicMomentumTransportModel>
62 (
63  const alphaField& alpha,
64  const rhoField& rho,
65  const volVectorField& U,
66  const surfaceScalarField& alphaRhoPhi,
67  const surfaceScalarField& phi,
68  const viscosity& viscosity,
69  const word& type
70 )
71 :
72  LESeddyViscosity<BasicMomentumTransportModel>
73  (
74  type,
75  alpha,
76  rho,
77  U,
78  alphaRhoPhi,
79  phi,
80  viscosity
81  ),
82 
83  flm_
84  (
85  IOobject
86  (
87  this->groupName("flm"),
88  this->runTime_.name(),
89  this->mesh_,
90  IOobject::MUST_READ,
91  IOobject::AUTO_WRITE
92  ),
93  this->mesh_
94  ),
95  fmm_
96  (
97  IOobject
98  (
99  this->groupName("fmm"),
100  this->runTime_.name(),
101  this->mesh_,
102  IOobject::MUST_READ,
103  IOobject::AUTO_WRITE
104  ),
105  this->mesh_
106  ),
107  theta_("theta", this->coeffDict(), 1.5),
108 
109  simpleFilter_(U.mesh()),
110  filterPtr_(LESfilter::New(U.mesh(), this->coeffDict())),
111  filter_(filterPtr_()),
112 
113  flm0_("flm0", flm_.dimensions(), 0.0),
114  fmm0_("fmm0", fmm_.dimensions(), vSmall)
115 {}
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
120 template<class BasicMomentumTransportModel>
122 {
124  {
125  filter_.read(this->coeffDict());
126  theta_.readIfPresent(this->coeffDict());
127 
128  return true;
129  }
130  else
131  {
132  return false;
133  }
134 }
135 
136 
137 template<class BasicMomentumTransportModel>
139 {
140  if (!this->turbulence_)
141  {
142  return;
143  }
144 
145  // Local references
146  const alphaField& alpha = this->alpha_;
147  const rhoField& rho = this->rho_;
148  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
149  const volVectorField& U = this->U_;
150  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
152  (
153  Foam::fvConstraints::New(this->mesh_)
154  );
155 
157 
159  const volTensorField& gradU = tgradU();
160 
161  volSymmTensorField S(dev(symm(gradU)));
162  volScalarField magS(mag(S));
163 
164  volVectorField Uf(filter_(U));
166  volScalarField magSf(mag(Sf));
167 
168  volSymmTensorField L(dev(filter_(sqr(U)) - (sqr(filter_(U)))));
170  (
171  2.0*sqr(this->delta())*(filter_(magS*S) - 4.0*magSf*Sf)
172  );
173 
174  volScalarField invT
175  (
176  alpha*rho*(1.0/(theta_.value()*this->delta()))*pow(flm_*fmm_, 1.0/8.0)
177  );
178 
179  volScalarField LM(L && M);
180 
181  fvScalarMatrix flmEqn
182  (
183  fvm::ddt(alpha, rho, flm_)
184  + fvm::div(alphaRhoPhi, flm_)
185  ==
186  invT*LM
187  - fvm::Sp(invT, flm_)
188  + fvModels.source(alpha, rho, flm_)
189  );
190 
191  flmEqn.relax();
192  fvConstraints.constrain(flmEqn);
193  flmEqn.solve();
194  fvConstraints.constrain(flm_);
195  bound(flm_, flm0_);
196 
197  volScalarField MM(M && M);
198 
199  fvScalarMatrix fmmEqn
200  (
201  fvm::ddt(alpha, rho, fmm_)
202  + fvm::div(alphaRhoPhi, fmm_)
203  ==
204  invT*MM
205  - fvm::Sp(invT, fmm_)
206  + fvModels.source(alpha, rho, fmm_)
207  );
208 
209  fmmEqn.relax();
210  fvConstraints.constrain(fmmEqn);
211  fmmEqn.solve();
212  fvConstraints.constrain(fmm_);
213  bound(fmm_, fmm0_);
214 
215  correctNut(gradU);
216 }
217 
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 } // End namespace LESModels
222 } // End namespace Foam
223 
224 // ************************************************************************* //
scalar delta
#define M(I)
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Eddy viscosity LES SGS model base class.
BasicMomentumTransportModel::alphaField alphaField
virtual void correct()
Correct Eddy-Viscosity and related properties.
dynamicLagrangian(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
virtual bool read()
Read model coefficients if they have changed.
BasicMomentumTransportModel::rhoField rhoField
Abstract class for LES filters.
Definition: LESfilter.H:55
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Finite volume constraints.
Definition: fvConstraints.H:67
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:603
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
Finite volume models.
Definition: fvModels.H:65
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488