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 | Copyright (C) 2015 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 "ReynoldsStress.H"
27 #include "fvc.H"
28 #include "fvm.H"
29 #include "wallFvPatch.H"
30 
31 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
32 
33 template<class BasicTurbulenceModel>
35 (
37 ) const
38 {
39  scalar kMin = this->kMin_.value();
40 
41  R.max
42  (
44  (
45  "zero",
46  R.dimensions(),
48  (
49  kMin, -GREAT, -GREAT,
50  kMin, -GREAT,
51  kMin
52  )
53  )
54  );
55 }
56 
57 
58 template<class BasicTurbulenceModel>
60 (
62 ) const
63 {
64  const fvPatchList& patches = this->mesh_.boundary();
65 
67  {
68  const fvPatch& curPatch = patches[patchi];
69 
70  if (isA<wallFvPatch>(curPatch))
71  {
72  symmTensorField& Rw = R.boundaryField()[patchi];
73 
74  const scalarField& nutw = this->nut_.boundaryField()[patchi];
75 
76  const vectorField snGradU
77  (
78  this->U_.boundaryField()[patchi].snGrad()
79  );
80 
81  const vectorField& faceAreas
82  = this->mesh_.Sf().boundaryField()[patchi];
83 
84  const scalarField& magFaceAreas
85  = this->mesh_.magSf().boundaryField()[patchi];
86 
87  forAll(curPatch, facei)
88  {
89  // Calculate near-wall velocity gradient
90  tensor gradUw
91  = (faceAreas[facei]/magFaceAreas[facei])*snGradU[facei];
92 
93  // Calculate near-wall shear-stress tensor
94  tensor tauw = -nutw[facei]*2*dev(symm(gradUw));
95 
96  // Reset the shear components of the stress tensor
97  Rw[facei].xy() = tauw.xy();
98  Rw[facei].xz() = tauw.xz();
99  Rw[facei].yz() = tauw.yz();
100  }
101  }
102  }
103 }
104 
105 
106 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
107 
108 template<class BasicTurbulenceModel>
110 (
111  const word& modelName,
112  const alphaField& alpha,
113  const rhoField& rho,
114  const volVectorField& U,
115  const surfaceScalarField& alphaRhoPhi,
116  const surfaceScalarField& phi,
117  const transportModel& transport,
118  const word& propertiesName
119 )
120 :
121  BasicTurbulenceModel
122  (
123  modelName,
124  alpha,
125  rho,
126  U,
127  alphaRhoPhi,
128  phi,
129  transport,
130  propertiesName
131  ),
132 
133  couplingFactor_
134  (
136  (
137  "couplingFactor",
138  this->coeffDict_,
139  0.0
140  )
141  ),
142 
143  R_
144  (
145  IOobject
146  (
147  IOobject::groupName("R", U.group()),
148  this->runTime_.timeName(),
149  this->mesh_,
150  IOobject::MUST_READ,
151  IOobject::AUTO_WRITE
152  ),
153  this->mesh_
154  ),
155 
156  nut_
157  (
158  IOobject
159  (
160  IOobject::groupName("nut", U.group()),
161  this->runTime_.timeName(),
162  this->mesh_,
163  IOobject::MUST_READ,
164  IOobject::AUTO_WRITE
165  ),
166  this->mesh_
167  )
168 {
169  if (couplingFactor_.value() < 0.0 || couplingFactor_.value() > 1.0)
170  {
172  << "couplingFactor = " << couplingFactor_
173  << " is not in range 0 - 1" << nl
174  << exit(FatalError);
175  }
176 }
177 
178 
179 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
180 
181 template<class BasicTurbulenceModel>
183 {
185 }
186 
187 
188 template<class BasicTurbulenceModel>
191 {
192  return R_;
193 }
194 
195 
196 template<class BasicTurbulenceModel>
199 {
200  tmp<Foam::volScalarField> tk(0.5*tr(R_));
201  tk().rename("k");
202  return tk;
203 }
204 
205 
206 template<class BasicTurbulenceModel>
209 {
211  (
213  (
214  IOobject
215  (
216  IOobject::groupName("devRhoReff", this->U_.group()),
217  this->runTime_.timeName(),
218  this->mesh_,
219  IOobject::NO_READ,
220  IOobject::NO_WRITE
221  ),
222  this->alpha_*this->rho_*R_
223  - (this->alpha_*this->rho_*this->nu())
224  *dev(twoSymm(fvc::grad(this->U_)))
225  )
226  );
227 }
228 
229 
230 template<class BasicTurbulenceModel>
233 (
235 ) const
236 {
237  if (couplingFactor_.value() > 0.0)
238  {
239  return
240  (
241  fvc::div
242  (
243  this->alpha_*this->rho_*R_
244  + couplingFactor_
245  *this->alpha_*this->rho_*this->nut()*fvc::grad(U),
246  "div(devRhoReff)"
247  )
249  (
250  (1.0 - couplingFactor_)*this->alpha_*this->rho_*this->nut(),
251  U,
252  "laplacian(nuEff,U)"
253  )
254  - fvm::laplacian(this->alpha_*this->rho_*this->nuEff(), U)
255  - fvc::div(this->alpha_*this->rho_*this->nu()*dev2(T(fvc::grad(U))))
256  );
257  }
258  else
259  {
260  return
261  (
262  fvc::div(this->alpha_*this->rho_*R_)
264  (
265  this->alpha_*this->rho_*this->nut(),
266  U,
267  "laplacian(nuEff,U)"
268  )
269  - fvm::laplacian(this->alpha_*this->rho_*this->nuEff(), U)
270  - fvc::div(this->alpha_*this->rho_*this->nu()*dev2(T(fvc::grad(U))))
271  );
272  }
273 
274  return
275  (
276  - fvm::laplacian(this->alpha_*this->rho_*this->nuEff(), U)
277  - fvc::div((this->alpha_*this->rho_*this->nuEff())*dev2(T(fvc::grad(U))))
278  );
279 }
280 
281 
282 template<class BasicTurbulenceModel>
285 (
286  const volScalarField& rho,
288 ) const
289 {
290  return
291  (
292  - fvm::laplacian(this->alpha_*rho*this->nuEff(), U)
293  - fvc::div((this->alpha_*rho*this->nuEff())*dev2(T(fvc::grad(U))))
294  );
295 }
296 
297 
298 template<class BasicTurbulenceModel>
300 {
301  correctNut();
302 }
303 
304 
305 template<class BasicTurbulenceModel>
307 {
309 }
310 
311 
312 // ************************************************************************* //
Foam::Tensor
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
Definition: complexI.H:224
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:82
Foam::ReynoldsStress::k
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: ReynoldsStress.C:198
Foam::SymmTensor< scalar >
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::ReynoldsStress::read
virtual bool read()=0
Re-read model coefficients if they have changed.
Definition: ReynoldsStress.C:182
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::ReynoldsStress::correctWallShearStress
void correctWallShearStress(volSymmTensorField &R) const
Definition: ReynoldsStress.C:60
wallFvPatch.H
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::dev2
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:115
U
U
Definition: pEqn.H:46
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
ReynoldsStress.H
R
#define R(A, B, C, D, E, F, K, M)
Foam::ReynoldsStress::divDevRhoReff
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
Definition: ReynoldsStress.C:233
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:110
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
nut
nut
Definition: createTDFields.H:71
Foam::nl
static const char nl
Definition: Ostream.H:260
correct
fvOptions correct(rho)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
tauw
scalar tauw
Definition: evaluateNearWall.H:12
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
fvm.H
Foam::FatalError
error FatalError
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:100
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::ReynoldsStress::validate
virtual void validate()
Validate the turbulence fields after construction.
Definition: ReynoldsStress.C:299
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::ReynoldsStress::devRhoReff
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
Definition: ReynoldsStress.C:208
rho
rho
Definition: pEqn.H:3
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
T
const volScalarField & T
Definition: createFields.H:25
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::ReynoldsStress::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: ReynoldsStress.C:306
Foam::ReynoldsStress::R
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
Definition: ReynoldsStress.C:190
fvc.H
Foam::ReynoldsStress::boundNormalStress
void boundNormalStress(volSymmTensorField &R) const
Definition: ReynoldsStress.C:35
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:49
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:93
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104