linearViscousStress.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2013-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "linearViscousStress.H"
29#include "fvc.H"
30#include "fvm.H"
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
34template<class BasicTurbulenceModel>
36(
37 const word& modelName,
38 const alphaField& alpha,
39 const rhoField& rho,
40 const volVectorField& U,
41 const surfaceScalarField& alphaRhoPhi,
43 const transportModel& transport,
44 const word& propertiesName
45)
46:
47 BasicTurbulenceModel
48 (
49 modelName,
50 alpha,
51 rho,
52 U,
53 alphaRhoPhi,
54 phi,
55 transport,
56 propertiesName
57 )
58{}
59
60
61// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62
63template<class BasicTurbulenceModel>
65{
66 return BasicTurbulenceModel::read();
67}
69
70template<class BasicTurbulenceModel>
73{
74 return devRhoReff(this->U_);
75}
76
77
78template<class BasicTurbulenceModel>
81(
82 const volVectorField& U
83) const
84{
86 (
88 (
90 (
91 IOobject::groupName("devRhoReff", this->alphaRhoPhi_.group()),
92 this->runTime_.timeName(),
93 this->mesh_,
96 ),
97 (-(this->alpha_*this->rho_*this->nuEff()))
99 )
100 );
101}
102
104template<class BasicTurbulenceModel>
107(
109) const
111 return
112 (
113 - fvc::div((this->alpha_*this->rho_*this->nuEff())*dev2(T(fvc::grad(U))))
114 - fvm::laplacian(this->alpha_*this->rho_*this->nuEff(), U)
115 );
116}
117
118
119template<class BasicTurbulenceModel>
122(
123 const volScalarField& rho,
125) const
126{
127 return
128 (
129 - fvc::div((this->alpha_*rho*this->nuEff())*dev2(T(fvc::grad(U))))
130 - fvm::laplacian(this->alpha_*rho*this->nuEff(), U)
131 );
132}
133
134
135template<class BasicTurbulenceModel>
137{
138 BasicTurbulenceModel::correct();
139}
140
141
142// ************************************************************************* //
surfaceScalarField & phi
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
tmp< volVectorField > divDevRhoReff()
Return the effective viscous stress (laminar + turbulent).
Definition: momentumError.C:52
Linear viscous stress turbulence model base class.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
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.
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
const volScalarField & T
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
volScalarField & alpha