ReynoldsStress.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) 2015-2017 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "ReynoldsStress.H"
30 #include "fvc.H"
31 #include "fvm.H"
32 #include "wallFvPatch.H"
33 
34 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
35 
36 template<class BasicTurbulenceModel>
38 (
40 ) const
41 {
42  scalar kMin = this->kMin_.value();
43 
44  R.max
45  (
47  (
48  "zero",
49  R.dimensions(),
51  (
52  kMin, -GREAT, -GREAT,
53  kMin, -GREAT,
54  kMin
55  )
56  )
57  );
58 }
59 
60 
61 template<class BasicTurbulenceModel>
63 (
65 ) const
66 {
67  const fvPatchList& patches = this->mesh_.boundary();
68 
69  volSymmTensorField::Boundary& RBf = R.boundaryFieldRef();
70 
71  forAll(patches, patchi)
72  {
73  const fvPatch& curPatch = patches[patchi];
74 
75  if (isA<wallFvPatch>(curPatch))
76  {
77  symmTensorField& Rw = RBf[patchi];
78 
79  const scalarField& nutw = this->nut_.boundaryField()[patchi];
80 
81  const vectorField snGradU
82  (
83  this->U_.boundaryField()[patchi].snGrad()
84  );
85 
86  const vectorField& faceAreas
87  = this->mesh_.Sf().boundaryField()[patchi];
88 
89  const scalarField& magFaceAreas
90  = this->mesh_.magSf().boundaryField()[patchi];
91 
92  forAll(curPatch, facei)
93  {
94  // Calculate near-wall velocity gradient
95  const tensor gradUw
96  = (faceAreas[facei]/magFaceAreas[facei])*snGradU[facei];
97 
98  // Set the wall Reynolds-stress to the near-wall shear-stress
99  // Note: the spherical part of the normal stress is included in
100  // the pressure
101  Rw[facei] = -nutw[facei]*2*dev(symm(gradUw));
102  }
103  }
104  }
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
109 
110 template<class BasicTurbulenceModel>
112 (
113  const word& modelName,
114  const alphaField& alpha,
115  const rhoField& rho,
116  const volVectorField& U,
117  const surfaceScalarField& alphaRhoPhi,
118  const surfaceScalarField& phi,
119  const transportModel& transport,
120  const word& propertiesName
121 )
122 :
123  BasicTurbulenceModel
124  (
125  modelName,
126  alpha,
127  rho,
128  U,
129  alphaRhoPhi,
130  phi,
131  transport,
132  propertiesName
133  ),
134 
135  couplingFactor_
136  (
138  (
139  "couplingFactor",
140  this->coeffDict_,
141  0.0
142  )
143  ),
144 
145  R_
146  (
147  IOobject
148  (
149  IOobject::groupName("R", alphaRhoPhi.group()),
150  this->runTime_.timeName(),
151  this->mesh_,
152  IOobject::MUST_READ,
153  IOobject::AUTO_WRITE
154  ),
155  this->mesh_
156  ),
157 
158  nut_
159  (
160  IOobject
161  (
162  IOobject::groupName("nut", alphaRhoPhi.group()),
163  this->runTime_.timeName(),
164  this->mesh_,
165  IOobject::MUST_READ,
166  IOobject::AUTO_WRITE
167  ),
168  this->mesh_
169  )
170 {
171  if (couplingFactor_.value() < 0.0 || couplingFactor_.value() > 1.0)
172  {
174  << "couplingFactor = " << couplingFactor_
175  << " is not in range 0 - 1" << nl
176  << exit(FatalError);
177  }
178 }
179 
180 
181 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
182 
183 template<class BasicTurbulenceModel>
185 {
187 }
188 
189 
190 template<class BasicTurbulenceModel>
193 {
194  return R_;
195 }
196 
197 
198 template<class BasicTurbulenceModel>
201 {
202  tmp<Foam::volScalarField> tk(0.5*tr(R_));
203  tk.ref().rename("k");
204  return tk;
205 }
206 
207 
208 template<class BasicTurbulenceModel>
211 {
213  (
215  (
216  IOobject
217  (
218  IOobject::groupName("devRhoReff", this->alphaRhoPhi_.group()),
219  this->runTime_.timeName(),
220  this->mesh_,
221  IOobject::NO_READ,
222  IOobject::NO_WRITE
223  ),
224  this->alpha_*this->rho_*R_
225  - (this->alpha_*this->rho_*this->nu())
226  *dev(twoSymm(fvc::grad(this->U_)))
227  )
228  );
229 }
230 
231 
232 template<class BasicTurbulenceModel>
233 template<class RhoFieldType>
236 (
237  const RhoFieldType& rho,
239 ) const
240 {
241  if (couplingFactor_.value() > 0.0)
242  {
243  return
244  (
246  (
247  (1.0 - couplingFactor_)*this->alpha_*rho*this->nut(),
248  U,
249  "laplacian(nuEff,U)"
250  )
251  + fvc::div
252  (
253  this->alpha_*rho*R_
254  + couplingFactor_
255  *this->alpha_*rho*this->nut()*fvc::grad(U),
256  "div(devRhoReff)"
257  )
258  - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
259  - fvm::laplacian(this->alpha_*rho*this->nuEff(), U)
260  );
261  }
262  else
263  {
264  return
265  (
267  (
268  this->alpha_*rho*this->nut(),
269  U,
270  "laplacian(nuEff,U)"
271  )
272  + fvc::div(this->alpha_*rho*R_)
273  - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
274  - fvm::laplacian(this->alpha_*rho*this->nuEff(), U)
275  );
276  }
277 }
278 
279 
280 template<class BasicTurbulenceModel>
283 (
285 ) const
286 {
287  return DivDevRhoReff(this->rho_, U);
288 }
289 
290 
291 template<class BasicTurbulenceModel>
294 (
295  const volScalarField& rho,
297 ) const
298 {
299  return DivDevRhoReff(rho, U);
300 }
301 
302 
303 template<class BasicTurbulenceModel>
305 {
306  correctNut();
307 }
308 
309 
310 template<class BasicTurbulenceModel>
312 {
314 }
315 
316 
317 // ************************************************************************* //
Foam::Tensor< scalar >
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
Foam::ReynoldsStress::k
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: ReynoldsStress.C:200
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::ReynoldsStress::read
virtual bool read()=0
Re-read model coefficients if they have changed.
Definition: ReynoldsStress.C:184
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::ReynoldsStress::correctWallShearStress
void correctWallShearStress(volSymmTensorField &R) const
Definition: ReynoldsStress.C:63
wallFvPatch.H
Foam::fac::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:56
Foam::dev2
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:117
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
rho
rho
Definition: readInitialConditions.H:88
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
ReynoldsStress.H
R
#define R(A, B, C, D, E, F, K, M)
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::ReynoldsStress::divDevRhoReff
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
Definition: ReynoldsStress.C:283
Foam::ReynoldsStress::ReynoldsStress
ReynoldsStress(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
Construct from components.
Definition: ReynoldsStress.C:112
Foam::Field< symmTensor >
correct
fvOptions correct(rho)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::ReynoldsStress::DivDevRhoReff
tmp< fvVectorMatrix > DivDevRhoReff(const RhoFieldType &rho, volVectorField &U) const
Return the source term for the momentum equation.
Foam::symmTensor
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition: symmTensor.H:59
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::blockMeshTools::read
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:33
fvm.H
Foam::FatalError
error FatalError
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:97
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam::ReynoldsStress::validate
virtual void validate()
Validate the turbulence fields after construction.
Definition: ReynoldsStress.C:304
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:96
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::ReynoldsStress::devRhoReff
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
Definition: ReynoldsStress.C:210
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:98
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::ReynoldsStress::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: ReynoldsStress.C:311
Foam::ReynoldsStress::R
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
Definition: ReynoldsStress.C:192
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
fvc.H
Foam::ReynoldsStress::boundNormalStress
void boundNormalStress(volSymmTensorField &R) const
Definition: ReynoldsStress.C:38
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:51
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
Foam::GeometricField< symmTensor, fvPatchField, volMesh >
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106