ReynoldsStress.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) 2015-2023 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 "fvcSnGrad.H"
28 #include "wallFvPatch.H"
29 
30 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
31 
32 template<class BasicMomentumTransportModel>
34 (
36 ) const
37 {
38  scalar kMin = this->kMin_.value();
39 
40  R.max
41  (
43  (
44  "zero",
45  R.dimensions(),
47  (
48  kMin, -great, -great,
49  kMin, -great,
50  kMin
51  )
52  )
53  );
54 }
55 
56 
57 template<class BasicMomentumTransportModel>
59 (
61 ) const
62 {
63  const fvPatchList& patches = this->mesh_.boundary();
64 
65  volSymmTensorField::Boundary& RBf = R.boundaryFieldRef();
66 
68  {
69  const fvPatch& curPatch = patches[patchi];
70 
71  if (isA<wallFvPatch>(curPatch))
72  {
73  symmTensorField& Rw = RBf[patchi];
74 
75  const scalarField& nutw = this->nut_.boundaryField()[patchi];
76 
77  const vectorField snGradUw
78  (
79  this->U_.boundaryField()[patchi].snGrad()
80  );
81 
82  const vectorField& Sf = this->mesh_.Sf().boundaryField()[patchi];
83 
84  const scalarField& magSf =
85  this->mesh_.magSf().boundaryField()[patchi];
86 
87  forAll(curPatch, facei)
88  {
89  // Calculate near-wall velocity gradient
90  const tensor gradUw
91  = (Sf[facei]/magSf[facei])*snGradUw[facei];
92 
93  // Set the wall Reynolds-stress to the near-wall shear-stress
94  // Note: the spherical part of the normal stress is included in
95  // the pressure
96  Rw[facei] = -nutw[facei]*dev(twoSymm(gradUw));
97  }
98  }
99  }
100 }
101 
102 
103 template<class BasicMomentumTransportModel>
106 {
108  (
110  (
111  this->R_,
112  dimVolume*this->rho_.dimensions()*this->R_.dimensions()/dimTime
113  )
114  );
115 }
116 
117 
118 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
119 
120 template<class BasicMomentumTransportModel>
122 (
123  const word& modelName,
124  const alphaField& alpha,
125  const rhoField& rho,
126  const volVectorField& U,
127  const surfaceScalarField& alphaRhoPhi,
128  const surfaceScalarField& phi,
129  const viscosity& viscosity
130 )
131 :
132  BasicMomentumTransportModel
133  (
134  modelName,
135  alpha,
137  U,
138  alphaRhoPhi,
139  phi,
140  viscosity
141  ),
142 
143  couplingFactor_
144  (
145  dimensioned<scalar>::lookupOrAddToDict
146  (
147  "couplingFactor",
148  this->coeffDict_,
149  1
150  )
151  ),
152 
153  R_
154  (
155  IOobject
156  (
157  this->groupName("R"),
158  this->runTime_.name(),
159  this->mesh_,
160  IOobject::MUST_READ,
161  IOobject::AUTO_WRITE
162  ),
163  this->mesh_
164  ),
165 
166  nut_
167  (
168  IOobject
169  (
170  this->groupName("nut"),
171  this->runTime_.name(),
172  this->mesh_,
173  IOobject::MUST_READ,
174  IOobject::AUTO_WRITE
175  ),
176  this->mesh_
177  )
178 {
179  if (couplingFactor_.value() < 0 || couplingFactor_.value() > 1)
180  {
182  << "couplingFactor = " << couplingFactor_
183  << " is not in range 0 - 1" << nl
184  << exit(FatalError);
185  }
186 }
187 
188 
189 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
190 
191 template<class BasicMomentumTransportModel>
193 {
195 }
196 
197 
198 template<class BasicMomentumTransportModel>
201 {
202  return R_;
203 }
204 
205 
206 template<class BasicMomentumTransportModel>
209 {
210  tmp<Foam::volScalarField> tk(0.5*tr(R_));
211  tk.ref().rename("k");
212  return tk;
213 }
214 
215 
216 template<class BasicMomentumTransportModel>
219 {
221  (
222  this->groupName("devTau"),
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 template<class BasicMomentumTransportModel>
231 template<class RhoFieldType>
234 (
235  const RhoFieldType& rho,
237 ) const
238 {
239  tmp<volTensorField> tgradU = fvc::grad(U);
240  const volTensorField& gradU = tgradU();
241  const surfaceTensorField gradUf(fvc::interpolate(gradU));
242 
243  // Interpolate Reynolds stress to the faces
244  // with either a stress or velocity coupling correction
245  const surfaceVectorField Refff
246  (
247  (this->mesh().Sf() & fvc::interpolate(R_))
248 
249  // Stress coupling
250  + couplingFactor_
251  *(this->mesh().Sf() & fvc::interpolate(this->nut()*gradU))
252 
253  // or velocity gradient coupling
254  // + couplingFactor_
255  // *fvc::interpolate(this->nut())*(this->mesh().Sf() & gradUf)
256 
257  - fvc::interpolate(couplingFactor_*this->nut() + this->nu())
258  *this->mesh().magSf()*fvc::snGrad(U)
259  - fvc::interpolate(this->nu())*(this->mesh().Sf() & dev2(gradUf.T()))
260  );
261 
262  return
263  (
264  fvc::div(fvc::interpolate(this->alpha_*rho)*Refff)
265  - correction(fvm::laplacian(this->alpha_*rho*this->nuEff(), U))
266  );
267 }
268 
269 
270 template<class BasicMomentumTransportModel>
273 (
275 ) const
276 {
277  return DivDevRhoReff(this->rho_, U);
278 }
279 
280 
281 template<class BasicMomentumTransportModel>
284 (
285  const volScalarField& rho,
287 ) const
288 {
289  return DivDevRhoReff(rho, U);
290 }
291 
292 
293 template<class BasicMomentumTransportModel>
295 {
296  correctNut();
297 }
298 
299 
300 template<class BasicMomentumTransportModel>
302 {
304 }
305 
306 
307 // ************************************************************************* //
#define R(A, B, C, D, E, F, K, M)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricBoundaryField class.
Generic GeometricField class.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
BasicMomentumTransportModel::alphaField alphaField
void boundNormalStress(volSymmTensorField &R) const
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual void validate()
Validate the turbulence fields after construction.
ReynoldsStress(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
virtual tmp< volSymmTensorField > devTau() const
Return the effective stress tensor.
virtual bool read()=0
Re-read model coefficients if they have changed.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the source term for the momentum equation.
virtual tmp< volSymmTensorField > sigma() const
Return the Reynolds stress tensor [m^2/s^2].
virtual tmp< fvSymmTensorMatrix > RSource() const
Source term for the R equation.
tmp< fvVectorMatrix > DivDevRhoReff(const RhoFieldType &rho, volVectorField &U) const
Return the source term for the momentum equation.
void correctWallShearStress(volSymmTensorField &R) const
BasicMomentumTransportModel::rhoField rhoField
Generic dimensioned Type class.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const scalar nut
const scalar nuEff
Calculate the snGrad of the given volField.
label patchi
const fvPatchList & patches
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimTime
const dimensionSet dimVolume
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
error FatalError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260