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-2020 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 
66  volSymmTensorField::Boundary& RBf = R.boundaryFieldRef();
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 snGradU
79  (
80  this->U_.boundaryField()[patchi].snGrad()
81  );
82 
83  const vectorField& faceAreas
84  = this->mesh_.Sf().boundaryField()[patchi];
85 
86  const scalarField& magFaceAreas
87  = this->mesh_.magSf().boundaryField()[patchi];
88 
89  forAll(curPatch, facei)
90  {
91  // Calculate near-wall velocity gradient
92  const tensor gradUw
93  = (faceAreas[facei]/magFaceAreas[facei])*snGradU[facei];
94 
95  // Set the wall Reynolds-stress to the near-wall shear-stress
96  // Note: the spherical part of the normal stress is included in
97  // the pressure
98  Rw[facei] = -nutw[facei]*2*dev(symm(gradUw));
99  }
100  }
101  }
102 }
103 
104 
105 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
106 
107 template<class BasicMomentumTransportModel>
109 (
110  const word& modelName,
111  const alphaField& alpha,
112  const rhoField& rho,
113  const volVectorField& U,
114  const surfaceScalarField& alphaRhoPhi,
115  const surfaceScalarField& phi,
116  const transportModel& transport
117 )
118 :
119  BasicMomentumTransportModel
120  (
121  modelName,
122  alpha,
123  rho,
124  U,
125  alphaRhoPhi,
126  phi,
127  transport
128  ),
129 
130  couplingFactor_
131  (
133  (
134  "couplingFactor",
135  this->coeffDict_,
136  0.0
137  )
138  ),
139 
140  R_
141  (
142  IOobject
143  (
144  IOobject::groupName("R", alphaRhoPhi.group()),
145  this->runTime_.timeName(),
146  this->mesh_,
147  IOobject::MUST_READ,
148  IOobject::AUTO_WRITE
149  ),
150  this->mesh_
151  ),
152 
153  nut_
154  (
155  IOobject
156  (
157  IOobject::groupName("nut", alphaRhoPhi.group()),
158  this->runTime_.timeName(),
159  this->mesh_,
160  IOobject::MUST_READ,
161  IOobject::AUTO_WRITE
162  ),
163  this->mesh_
164  )
165 {
166  if (couplingFactor_.value() < 0.0 || couplingFactor_.value() > 1.0)
167  {
169  << "couplingFactor = " << couplingFactor_
170  << " is not in range 0 - 1" << nl
171  << exit(FatalError);
172  }
173 }
174 
175 
176 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
177 
178 template<class BasicMomentumTransportModel>
180 {
182 }
183 
184 
185 template<class BasicMomentumTransportModel>
188 {
189  return R_;
190 }
191 
192 
193 template<class BasicMomentumTransportModel>
196 {
197  tmp<Foam::volScalarField> tk(0.5*tr(R_));
198  tk.ref().rename("k");
199  return tk;
200 }
201 
202 
203 template<class BasicMomentumTransportModel>
206 {
208  (
209  IOobject::groupName("devTau", this->alphaRhoPhi_.group()),
210  this->alpha_*this->rho_*R_
211  - (this->alpha_*this->rho_*this->nu())
212  *dev(twoSymm(fvc::grad(this->U_)))
213  );
214 }
215 
216 
217 template<class BasicMomentumTransportModel>
218 template<class RhoFieldType>
221 (
222  const RhoFieldType& rho,
223  volVectorField& U
224 ) const
225 {
226  if (couplingFactor_.value() > 0.0)
227  {
228  return
229  (
231  (
232  (1.0 - couplingFactor_)*this->alpha_*rho*this->nut(),
233  U,
234  "laplacian(nuEff,U)"
235  )
236  + fvc::div
237  (
238  this->alpha_*rho*R_
239  + couplingFactor_
240  *this->alpha_*rho*this->nut()*fvc::grad(U),
241  "div(devTau)"
242  )
243  - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
244  - fvm::laplacian(this->alpha_*rho*this->nuEff(), U)
245  );
246  }
247  else
248  {
249  return
250  (
252  (
253  this->alpha_*rho*this->nut(),
254  U,
255  "laplacian(nuEff,U)"
256  )
257  + fvc::div(this->alpha_*rho*R_)
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 }
263 
264 
265 template<class BasicMomentumTransportModel>
268 (
269  volVectorField& U
270 ) const
271 {
272  return DivDevRhoReff(this->rho_, U);
273 }
274 
275 
276 template<class BasicMomentumTransportModel>
279 (
280  const volScalarField& rho,
281  volVectorField& U
282 ) const
283 {
284  return DivDevRhoReff(rho, U);
285 }
286 
287 
288 template<class BasicMomentumTransportModel>
290 {
291  correctNut();
292 }
293 
294 
295 template<class BasicMomentumTransportModel>
297 {
299 }
300 
301 
302 // ************************************************************************* //
void correctWallShearStress(volSymmTensorField &R) const
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:144
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
ReynoldsStress(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Construct from components.
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
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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:174
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:61
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Generic dimensioned Type class.
patches[0]
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
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 & dimensions() const
Return dimensions.
virtual bool read()=0
Re-read model coefficients if they have changed.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
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.
void boundNormalStress(volSymmTensorField &R) const
static const char nl
Definition: Ostream.H:260
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/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
BasicMomentumTransportModel::rhoField rhoField
Definition: LESModel.H:99
virtual void rename(const word &newName)
Rename.
Definition: regIOobject.C:400
const volScalarField & T
BasicMomentumTransportModel::transportModel transportModel
Definition: LESModel.H:100
label patchi
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
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.
BasicMomentumTransportModel::alphaField alphaField
Definition: LESModel.H:98
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:92
scalar nut