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-2022 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 BasicMomentumTransportModel>
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 BasicMomentumTransportModel>
60 (
62 ) const
63 {
64  const fvPatchList& patches = this->mesh_.boundary();
65 
67 
68  forAll(patches, patchi)
69  {
70  const fvPatch& curPatch = patches[patchi];
71 
72  if (isA<wallFvPatch>(curPatch))
73  {
74  symmTensorField& Rw = RBf[patchi];
75 
76  const scalarField& nutw = this->nut_.boundaryField()[patchi];
77 
78  const vectorField snGradUw
79  (
80  this->U_.boundaryField()[patchi].snGrad()
81  );
82 
83  const vectorField& Sf = this->mesh_.Sf().boundaryField()[patchi];
84 
85  const scalarField& magSf =
86  this->mesh_.magSf().boundaryField()[patchi];
87 
88  forAll(curPatch, facei)
89  {
90  // Calculate near-wall velocity gradient
91  const tensor gradUw
92  = (Sf[facei]/magSf[facei])*snGradUw[facei];
93 
94  // Set the wall Reynolds-stress to the near-wall shear-stress
95  // Note: the spherical part of the normal stress is included in
96  // the pressure
97  Rw[facei] = -nutw[facei]*dev(twoSymm(gradUw));
98  }
99  }
100  }
101 }
102 
103 
104 template<class BasicMomentumTransportModel>
107 {
109  (
111  (
112  this->R_,
113  dimVolume*this->rho_.dimensions()*this->R_.dimensions()/dimTime
114  )
115  );
116 }
117 
118 
119 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
120 
121 template<class BasicMomentumTransportModel>
123 (
124  const word& modelName,
125  const alphaField& alpha,
126  const rhoField& rho,
127  const volVectorField& U,
128  const surfaceScalarField& alphaRhoPhi,
129  const surfaceScalarField& phi,
130  const viscosity& viscosity
131 )
132 :
133  BasicMomentumTransportModel
134  (
135  modelName,
136  alpha,
137  rho,
138  U,
139  alphaRhoPhi,
140  phi,
141  viscosity
142  ),
143 
144  couplingFactor_
145  (
147  (
148  "couplingFactor",
149  this->coeffDict_,
150  1
151  )
152  ),
153 
154  R_
155  (
156  IOobject
157  (
158  IOobject::groupName("R", alphaRhoPhi.group()),
159  this->runTime_.timeName(),
160  this->mesh_,
161  IOobject::MUST_READ,
162  IOobject::AUTO_WRITE
163  ),
164  this->mesh_
165  ),
166 
167  nut_
168  (
169  IOobject
170  (
171  IOobject::groupName("nut", alphaRhoPhi.group()),
172  this->runTime_.timeName(),
173  this->mesh_,
174  IOobject::MUST_READ,
175  IOobject::AUTO_WRITE
176  ),
177  this->mesh_
178  )
179 {
180  if (couplingFactor_.value() < 0 || couplingFactor_.value() > 1)
181  {
183  << "couplingFactor = " << couplingFactor_
184  << " is not in range 0 - 1" << nl
185  << exit(FatalError);
186  }
187 }
188 
189 
190 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
191 
192 template<class BasicMomentumTransportModel>
194 {
196 }
197 
198 
199 template<class BasicMomentumTransportModel>
202 {
203  return R_;
204 }
205 
206 
207 template<class BasicMomentumTransportModel>
210 {
211  tmp<Foam::volScalarField> tk(0.5*tr(R_));
212  tk.ref().rename("k");
213  return tk;
214 }
215 
216 
217 template<class BasicMomentumTransportModel>
220 {
222  (
223  IOobject::groupName("devTau", this->alphaRhoPhi_.group()),
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 template<class BasicMomentumTransportModel>
232 template<class RhoFieldType>
235 (
236  const RhoFieldType& rho,
237  volVectorField& U
238 ) const
239 {
240  tmp<volTensorField> tgradU = fvc::grad(U);
241  const volTensorField& gradU = tgradU();
242  const surfaceTensorField gradUf(fvc::interpolate(gradU));
243 
244  // Interpolate Reynolds stress to the faces
245  // with either a stress or velocity coupling correction
246  const surfaceVectorField Refff
247  (
248  (this->mesh().Sf() & fvc::interpolate(R_))
249 
250  // Stress coupling
251  + couplingFactor_
252  *(this->mesh().Sf() & fvc::interpolate(this->nut()*gradU))
253 
254  // or velocity gradient coupling
255  // + couplingFactor_
256  // *fvc::interpolate(this->nut())*(this->mesh().Sf() & gradUf)
257 
258  - fvc::interpolate(couplingFactor_*this->nut() + this->nu())
259  *this->mesh().magSf()*fvc::snGrad(U)
260  - fvc::interpolate(this->nu())*(this->mesh().Sf() & dev2(gradUf.T()))
261  );
262 
263  return
264  (
265  fvc::div(fvc::interpolate(this->alpha_*rho)*Refff)
266  - correction(fvm::laplacian(this->alpha_*rho*this->nuEff(), U))
267  );
268 }
269 
270 
271 template<class BasicMomentumTransportModel>
274 (
275  volVectorField& U
276 ) const
277 {
278  return DivDevRhoReff(this->rho_, U);
279 }
280 
281 
282 template<class BasicMomentumTransportModel>
285 (
286  const volScalarField& rho,
287  volVectorField& U
288 ) const
289 {
290  return DivDevRhoReff(rho, U);
291 }
292 
293 
294 template<class BasicMomentumTransportModel>
296 {
297  correctNut();
298 }
299 
300 
301 template<class BasicMomentumTransportModel>
303 {
305 }
306 
307 
308 // ************************************************************************* //
void correctWallShearStress(volSymmTensorField &R) const
const fvPatchList & patches
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:134
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the source term for the momentum equation.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
virtual tmp< volSymmTensorField > sigma() const
Return the Reynolds stress tensor [m^2/s^2].
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
U
Definition: pEqn.H:72
error FatalError
virtual tmp< fvSymmTensorMatrix > RSource() const
Source term for the R equation.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
tmp< fvVectorMatrix > DivDevRhoReff(const RhoFieldType &rho, volVectorField &U) const
Return the source term for the momentum equation.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
fvMatrix< symmTensor > fvSymmTensorMatrix
Definition: fvMatricesFwd.H:47
Generic dimensioned Type class.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
fvMesh & mesh
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimTime
const dimensionSet & dimensions() const
Return dimensions.
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 bool read()=0
Re-read model coefficients if they have changed.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
SymmTensor< scalar > symmTensor
SymmTensor of scalars.
Definition: symmTensor.H:48
A class for handling words, derived from string.
Definition: word.H:59
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
virtual void validate()
Validate the turbulence fields after construction.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
void boundNormalStress(volSymmTensorField &R) const
static const char nl
Definition: Ostream.H:260
BasicMomentumTransportModel::rhoField rhoField
Definition: LESModel.H:107
virtual void rename(const word &newName)
Rename.
Definition: regIOobject.C:417
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
label patchi
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
virtual tmp< volSymmTensorField > devTau() const
Return the effective stress tensor.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
void max(const dimensioned< Type > &)
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
const dimensionSet dimVolume
BasicMomentumTransportModel::alphaField alphaField
Definition: LESModel.H:106
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
scalar nut