nonlinearEddyViscosity.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-2018 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 "nonlinearEddyViscosity.H"
27 #include "fvc.H"
28 #include "fvm.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class BasicTurbulenceModel>
34 (
35  const word& modelName,
36  const alphaField& alpha,
37  const rhoField& rho,
38  const volVectorField& U,
39  const surfaceScalarField& alphaRhoPhi,
40  const surfaceScalarField& phi,
41  const transportModel& transport,
42  const word& propertiesName
43 )
44 :
46  (
47  modelName,
48  alpha,
49  rho,
50  U,
51  alphaRhoPhi,
52  phi,
53  transport,
54  propertiesName
55  ),
56 
57  nonlinearStress_
58  (
59  IOobject
60  (
61  IOobject::groupName("nonlinearStress", alphaRhoPhi.group()),
62  this->runTime_.timeName(),
63  this->mesh_
64  ),
65  this->mesh_,
67  (
68  "nonlinearStress",
70  Zero
71  )
72  )
73 {}
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
78 template<class BasicTurbulenceModel>
81 {
83  (
85  );
86  tR.ref() += nonlinearStress_;
87  return tR;
88 }
89 
90 
91 template<class BasicTurbulenceModel>
94 {
95  tmp<volSymmTensorField> tdevRhoReff
96  (
98  );
99  tdevRhoReff.ref() += this->rho_*nonlinearStress_;
100  return tdevRhoReff;
101 }
102 
103 
104 template<class BasicTurbulenceModel>
107 (
108  volVectorField& U
109 ) const
110 {
111  return
112  (
113  fvc::div(this->rho_*nonlinearStress_)
115  );
116 }
117 
118 
119 template<class BasicTurbulenceModel>
122 (
123  const volScalarField& rho,
124  volVectorField& U
125 ) const
126 {
127  return
128  (
129  fvc::div(rho*nonlinearStress_)
131  );
132 }
133 
134 
135 // ************************************************************************* //
surfaceScalarField & phi
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
incompressible::RASModel ::rhoField rhoField
Definition: eddyViscosity.H:71
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
incompressible::RASModel ::transportModel transportModel
Definition: eddyViscosity.H:72
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
A class for handling words, derived from string.
Definition: word.H:59
static const zero Zero
Definition: zero.H:97
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
nonlinearEddyViscosity(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.
incompressible::RASModel ::alphaField alphaField
Definition: eddyViscosity.H:70
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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
const dimensionSet dimVelocity