1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2008 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
28 #include "addToRunTimeSelectionTable.H"
29 #include "wallFvPatch.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 namespace incompressible
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 defineTypeNameAndDebug(LRR, 0);
43 addToRunTimeSelectionTable(RASModel, LRR, dictionary);
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 const volVectorField& U,
51 const surfaceScalarField& phi,
52 transportModel& lamTransportModel
55 RASModel(typeName, U, phi, lamTransportModel),
59 dimensioned<scalar>::lookupOrAddToDict
68 dimensioned<scalar>::lookupOrAddToDict
77 dimensioned<scalar>::lookupOrAddToDict
86 dimensioned<scalar>::lookupOrAddToDict
95 dimensioned<scalar>::lookupOrAddToDict
104 dimensioned<scalar>::lookupOrAddToDict
113 dimensioned<scalar>::lookupOrAddToDict
122 dimensioned<scalar>::lookupOrAddToDict
131 dimensioned<scalar>::lookupOrAddToDict
178 nut_(Cmu_*sqr(k_)/(epsilon_ + epsilonSmall_))
180 # include "wallViscosityI.H"
182 if (couplingFactor_.value() < 0.0 || couplingFactor_.value() > 1.0)
187 "(const volVectorField& U, const surfaceScalarField& phi,"
188 "transportModel& lamTransportModel)"
189 ) << "couplingFactor = " << couplingFactor_
190 << " is not in range 0 - 1" << nl
198 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
200 tmp<volSymmTensorField> LRR::devReff() const
202 return tmp<volSymmTensorField>
204 new volSymmTensorField
214 R_ - nu()*dev(twoSymm(fvc::grad(U_)))
220 tmp<fvVectorMatrix> LRR::divDevReff(volVectorField& U) const
222 if (couplingFactor_.value() > 0.0)
226 fvc::div(R_ + couplingFactor_*nut_*fvc::grad(U), "div(R)")
229 (1.0 - couplingFactor_)*nut_,
233 - fvm::laplacian(nuEff(), U)
241 + fvc::laplacian(nut_, U, "laplacian(nuEff,U)")
242 - fvm::laplacian(nuEff(), U)
250 if (RASModel::read())
252 Cmu_.readIfPresent(coeffDict_);
253 Clrr1_.readIfPresent(coeffDict_);
254 Clrr2_.readIfPresent(coeffDict_);
255 C1_.readIfPresent(coeffDict_);
256 C2_.readIfPresent(coeffDict_);
257 Cs_.readIfPresent(coeffDict_);
258 Ceps_.readIfPresent(coeffDict_);
259 alphaEps_.readIfPresent(coeffDict_);
261 couplingFactor_.readIfPresent(coeffDict_);
263 if (couplingFactor_.value() < 0.0 || couplingFactor_.value() > 1.0)
265 FatalErrorIn("LRR::read()")
266 << "couplingFactor = " << couplingFactor_
267 << " is not in range 0 - 1"
282 transportModel_.correct();
291 volSymmTensorField P = -twoSymm(R_ & fvc::grad(U_));
292 volScalarField G = 0.5*tr(P);
294 # include "wallFunctionsI.H"
296 // Dissipation equation
297 tmp<fvScalarMatrix> epsEqn
300 + fvm::div(phi_, epsilon_)
301 //- fvm::laplacian(Ceps*(K/epsilon_)*R, epsilon_)
302 - fvm::laplacian(DepsilonEff(), epsilon_)
305 - fvm::Sp(C2_*epsilon_/k_, epsilon_)
310 # include "wallDissipationI.H"
313 bound(epsilon_, epsilon0_);
316 // Reynolds stress equation
318 const fvPatchList& patches = mesh_.boundary();
320 forAll(patches, patchi)
322 const fvPatch& curPatch = patches[patchi];
324 if (typeid(curPatch) == typeid(wallFvPatch))
326 forAll(curPatch, facei)
328 label faceCelli = curPatch.faceCells()[facei];
330 *= min(G[faceCelli]/(0.5*tr(P[faceCelli]) + SMALL), 1.0);
336 tmp<fvSymmTensorMatrix> REqn
340 //- fvm::laplacian(Cs*(k_/epsilon_)*R_, R_)
341 - fvm::laplacian(DREff(), R_)
342 + fvm::Sp(Clrr1_*epsilon_/k_, R_)
345 - (2.0/3.0*(1 - Clrr1_)*I)*epsilon_
354 dimensionedSymmTensor
360 k0_.value(), -GREAT, -GREAT,
371 // Re-calculate viscosity
372 nut_ = Cmu_*sqr(k_)/epsilon_;
374 # include "wallViscosityI.H"
377 // Correct wall shear stresses
379 forAll(patches, patchi)
381 const fvPatch& curPatch = patches[patchi];
383 if (typeid(curPatch) == typeid(wallFvPatch))
385 symmTensorField& Rw = R_.boundaryField()[patchi];
387 const scalarField& nutw = nut_.boundaryField()[patchi];
389 vectorField snGradU = U_.boundaryField()[patchi].snGrad();
391 const vectorField& faceAreas
392 = mesh_.Sf().boundaryField()[patchi];
394 const scalarField& magFaceAreas
395 = mesh_.magSf().boundaryField()[patchi];
397 forAll(curPatch, facei)
399 // Calculate near-wall velocity gradient
401 = (faceAreas[facei]/magFaceAreas[facei])*snGradU[facei];
403 // Calculate near-wall shear-stress tensor
404 tensor tauw = -nutw[facei]*2*symm(gradUw);
406 // Reset the shear components of the stress tensor
407 Rw[facei].xy() = tauw.xy();
408 Rw[facei].xz() = tauw.xz();
409 Rw[facei].yz() = tauw.yz();
416 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
418 } // End namespace RASModels
419 } // End namespace incompressible
420 } // End namespace Foam
422 // ************************************************************************* //