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-2021 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 :
73  (
74  type,
75  alpha,
76  rho,
77  U,
78  alphaRhoPhi,
79  phi,
80  viscosity
81  ),
82 
83  flm_
84  (
85  IOobject
86  (
87  IOobject::groupName("flm", this->alphaRhoPhi_.group()),
88  this->runTime_.timeName(),
89  this->mesh_,
92  ),
93  this->mesh_
94  ),
95  fmm_
96  (
97  IOobject
98  (
99  IOobject::groupName("fmm", this->alphaRhoPhi_.group()),
100  this->runTime_.timeName(),
101  this->mesh_,
104  ),
105  this->mesh_
106  ),
107  theta_
108  (
110  (
111  "theta",
112  this->coeffDict_,
113  1.5
114  )
115  ),
116 
117  simpleFilter_(U.mesh()),
118  filterPtr_(LESfilter::New(U.mesh(), this->coeffDict())),
119  filter_(filterPtr_()),
120 
121  flm0_("flm0", flm_.dimensions(), 0.0),
122  fmm0_("fmm0", fmm_.dimensions(), vSmall)
123 {
124  if (type == typeName)
125  {
126  this->printCoeffs(type);
127  }
128 }
129 
130 
131 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
132 
133 template<class BasicMomentumTransportModel>
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 BasicMomentumTransportModel>
152 {
153  if (!this->turbulence_)
154  {
155  return;
156  }
157 
158  // Local references
159  const alphaField& alpha = this->alpha_;
160  const rhoField& rho = this->rho_;
161  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
162  const volVectorField& U = this->U_;
163  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
165  (
166  Foam::fvConstraints::New(this->mesh_)
167  );
168 
170 
171  tmp<volTensorField> tgradU(fvc::grad(U));
172  const volTensorField& gradU = tgradU();
173 
174  volSymmTensorField S(dev(symm(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  + fvModels.source(alpha, rho, flm_)
202  );
203 
204  flmEqn.relax();
205  fvConstraints.constrain(flmEqn);
206  flmEqn.solve();
207  fvConstraints.constrain(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  + fvModels.source(alpha, rho, fmm_)
220  );
221 
222  fmmEqn.relax();
223  fvConstraints.constrain(fmmEqn);
224  fmmEqn.solve();
225  fvConstraints.constrain(fmm_);
226  bound(fmm_, fmm0_);
227 
228  correctNut(gradU);
229 }
230 
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 } // End namespace LESModels
235 } // End namespace Foam
236 
237 // ************************************************************************* //
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
U
Definition: pEqn.H:72
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
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 void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Foam::fvConstraints & fvConstraints
A class for handling words, derived from string.
Definition: word.H:59
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const GeometricField< Type, fvPatchField, volMesh > &)
static word groupName(Name name, const word &group)
Finite volume models.
Definition: fvModels.H:60
Info<< "Reading field p\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, simple.dict());mesh.schemes().setFluxRequired(p.name());Info<< "Reading field pa\"<< endl;volScalarField pa(IOobject("pa", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Ua\"<< endl;volVectorField Ua(IOobject("Ua", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);# 65 "/home/ubuntu/OpenFOAM-10/applications/solvers/incompressible/adjointShapeOptimisationFoam/createFields.H" 2label paRefCell=0;scalar paRefValue=0.0;setRefCell(pa, simple.dict(), paRefCell, paRefValue);mesh.schemes().setFluxRequired(pa.name());autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
phi
Definition: correctPhi.H:3
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
virtual bool read()
Read model coefficients if they have changed.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:604
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
BasicMomentumTransportModel::rhoField rhoField
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Finite volume constraints.
Definition: fvConstraints.H:57
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void correct()
Correct Eddy-Viscosity and related properties.
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:98
BasicMomentumTransportModel::alphaField alphaField
#define M(I)
Namespace for OpenFOAM.