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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "ReynoldsStress.H"
29 #include "fvc.H"
30 #include "fvm.H"
31 #include "wallFvPatch.H"
32 
33 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
34 
35 template<class BasicTurbulenceModel>
37 (
39 ) const
40 {
41  scalar kMin = this->kMin_.value();
42 
43  R.max
44  (
46  (
47  "zero",
48  R.dimensions(),
50  (
51  kMin, -GREAT, -GREAT,
52  kMin, -GREAT,
53  kMin
54  )
55  )
56  );
57 }
58 
59 
60 template<class BasicTurbulenceModel>
62 (
64 ) const
65 {
66  const fvPatchList& patches = this->mesh_.boundary();
67 
68  volSymmTensorField::Boundary& RBf = R.boundaryFieldRef();
69 
70  forAll(patches, patchi)
71  {
72  const fvPatch& curPatch = patches[patchi];
73 
74  if (isA<wallFvPatch>(curPatch))
75  {
76  symmTensorField& Rw = RBf[patchi];
77 
78  const scalarField& nutw = this->nut_.boundaryField()[patchi];
79 
80  const vectorField snGradU
81  (
82  this->U_.boundaryField()[patchi].snGrad()
83  );
84 
85  const vectorField& faceAreas
86  = this->mesh_.Sf().boundaryField()[patchi];
87 
88  const scalarField& magFaceAreas
89  = this->mesh_.magSf().boundaryField()[patchi];
90 
91  forAll(curPatch, facei)
92  {
93  // Calculate near-wall velocity gradient
94  const tensor gradUw
95  = (faceAreas[facei]/magFaceAreas[facei])*snGradU[facei];
96 
97  // Set the wall Reynolds-stress to the near-wall shear-stress
98  // Note: the spherical part of the normal stress is included in
99  // the pressure
100  Rw[facei] = -nutw[facei]*2*dev(symm(gradUw));
101  }
102  }
103  }
104 }
105 
106 
107 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
108 
109 template<class BasicTurbulenceModel>
111 (
112  const word& modelName,
113  const alphaField& alpha,
114  const rhoField& rho,
115  const volVectorField& U,
116  const surfaceScalarField& alphaRhoPhi,
117  const surfaceScalarField& phi,
118  const transportModel& transport,
119  const word& propertiesName
120 )
121 :
122  BasicTurbulenceModel
123  (
124  modelName,
125  alpha,
126  rho,
127  U,
128  alphaRhoPhi,
129  phi,
130  transport,
131  propertiesName
132  ),
133 
134  couplingFactor_
135  (
137  (
138  "couplingFactor",
139  this->coeffDict_,
140  0.0
141  )
142  ),
143 
144  R_
145  (
146  IOobject
147  (
148  IOobject::groupName("R", alphaRhoPhi.group()),
149  this->runTime_.timeName(),
150  this->mesh_,
151  IOobject::MUST_READ,
152  IOobject::AUTO_WRITE
153  ),
154  this->mesh_
155  ),
156 
157  nut_
158  (
159  IOobject
160  (
161  IOobject::groupName("nut", alphaRhoPhi.group()),
162  this->runTime_.timeName(),
163  this->mesh_,
164  IOobject::MUST_READ,
165  IOobject::AUTO_WRITE
166  ),
167  this->mesh_
168  )
169 {
170  if (couplingFactor_.value() < 0.0 || couplingFactor_.value() > 1.0)
171  {
173  << "couplingFactor = " << couplingFactor_
174  << " is not in range 0 - 1" << nl
175  << exit(FatalError);
176  }
177 }
178 
179 
180 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
181 
182 template<class BasicTurbulenceModel>
184 {
186 }
187 
188 
189 template<class BasicTurbulenceModel>
192 {
193  return R_;
194 }
195 
196 
197 template<class BasicTurbulenceModel>
200 {
201  tmp<Foam::volScalarField> tk(0.5*tr(R_));
202  tk.ref().rename("k");
203  return tk;
204 }
205 
206 
207 template<class BasicTurbulenceModel>
210 {
212  (
214  (
215  IOobject
216  (
217  IOobject::groupName("devRhoReff", this->alphaRhoPhi_.group()),
218  this->runTime_.timeName(),
219  this->mesh_,
220  IOobject::NO_READ,
221  IOobject::NO_WRITE
222  ),
223  this->alpha_*this->rho_*R_
224  - (this->alpha_*this->rho_*this->nu())
225  *dev(twoSymm(fvc::grad(this->U_)))
226  )
227  );
228 }
229 
230 
231 template<class BasicTurbulenceModel>
232 template<class RhoFieldType>
235 (
236  const RhoFieldType& rho,
238 ) const
239 {
240  if (couplingFactor_.value() > 0.0)
241  {
242  return
243  (
245  (
246  (1.0 - couplingFactor_)*this->alpha_*rho*this->nut(),
247  U,
248  "laplacian(nuEff,U)"
249  )
250  + fvc::div
251  (
252  this->alpha_*rho*R_
253  + couplingFactor_
254  *this->alpha_*rho*this->nut()*fvc::grad(U),
255  "div(devRhoReff)"
256  )
257  - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
258  - fvm::laplacian(this->alpha_*rho*this->nuEff(), U)
259  );
260  }
261  else
262  {
263  return
264  (
266  (
267  this->alpha_*rho*this->nut(),
268  U,
269  "laplacian(nuEff,U)"
270  )
271  + fvc::div(this->alpha_*rho*R_)
272  - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
273  - fvm::laplacian(this->alpha_*rho*this->nuEff(), U)
274  );
275  }
276 }
277 
278 
279 template<class BasicTurbulenceModel>
282 (
284 ) const
285 {
286  return DivDevRhoReff(this->rho_, U);
287 }
288 
289 
290 template<class BasicTurbulenceModel>
293 (
294  const volScalarField& rho,
296 ) const
297 {
298  return DivDevRhoReff(rho, U);
299 }
300 
301 
302 template<class BasicTurbulenceModel>
304 {
305  correctNut();
306 }
307 
308 
309 template<class BasicTurbulenceModel>
311 {
313 }
314 
315 
316 // ************************************************************************* //
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:199
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:59
Foam::ReynoldsStress::read
virtual bool read()=0
Re-read model coefficients if they have changed.
Definition: ReynoldsStress.C:183
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:62
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:96
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
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:258
Foam::ReynoldsStress::divDevRhoReff
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
Definition: ReynoldsStress.C:282
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:111
Foam::Field< symmTensor >
correct
fvOptions correct(rho)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
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.
Definition: symmTensor.H:50
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
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:303
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:209
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:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::ReynoldsStress::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: ReynoldsStress.C:310
Foam::ReynoldsStress::R
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
Definition: ReynoldsStress.C:191
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
fvc.H
Foam::ReynoldsStress::boundNormalStress
void boundNormalStress(volSymmTensorField &R) const
Definition: ReynoldsStress.C:37
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