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 -------------------------------------------------------------------------------
10 License
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 
34 template<class BasicTurbulenceModel>
36 (
37  const word& modelName,
38  const alphaField& alpha,
39  const rhoField& rho,
40  const volVectorField& U,
41  const surfaceScalarField& alphaRhoPhi,
42  const surfaceScalarField& phi,
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 
63 template<class BasicTurbulenceModel>
65 {
67 }
68 
69 
70 template<class BasicTurbulenceModel>
73 {
74  return devRhoReff(this->U_);
75 }
76 
77 
78 template<class BasicTurbulenceModel>
81 (
82  const volVectorField& U
83 ) const
84 {
86  (
88  (
89  IOobject
90  (
91  IOobject::groupName("devRhoReff", this->alphaRhoPhi_.group()),
92  this->runTime_.timeName(),
93  this->mesh_,
94  IOobject::NO_READ,
95  IOobject::NO_WRITE
96  ),
97  (-(this->alpha_*this->rho_*this->nuEff()))
99  )
100  );
101 }
102 
103 
104 template<class BasicTurbulenceModel>
107 (
109 ) const
110 {
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 
119 template<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 
135 template<class BasicTurbulenceModel>
137 {
139 }
140 
141 
142 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::fac::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:56
Foam::linearViscousStress::read
virtual bool read()=0
Re-read model coefficients if they have changed.
Definition: linearViscousStress.C:64
Foam::dev2
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:117
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
rho
rho
Definition: readInitialConditions.H:88
Foam::blockMeshTools::read
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:57
Foam::linearViscousStress::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: linearViscousStress.C:136
correct
fvOptions correct(rho)
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::linearViscousStress::linearViscousStress
linearViscousStress(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
Construct from components.
Definition: linearViscousStress.C:36
Foam::linearViscousStress::devRhoReff
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
Definition: linearViscousStress.C:72
fvm.H
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:98
Foam::linearViscousStress::divDevRhoReff
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
Definition: linearViscousStress.C:107
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:97
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:99
U
U
Definition: pEqn.H:72
fvc.H
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
linearViscousStress.H
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106