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-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 "nonlinearEddyViscosity.H"
27 #include "fvc.H"
28 #include "fvm.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class BasicMomentumTransportModel>
35 (
36  const word& modelName,
37  const alphaField& alpha,
38  const rhoField& rho,
39  const volVectorField& U,
40  const surfaceScalarField& alphaRhoPhi,
41  const surfaceScalarField& phi,
42  const transportModel& transport
43 )
44 :
46  (
47  modelName,
48  alpha,
49  rho,
50  U,
51  alphaRhoPhi,
52  phi,
53  transport
54  ),
55 
56  nonlinearStress_
57  (
58  IOobject
59  (
60  IOobject::groupName("nonlinearStress", alphaRhoPhi.group()),
61  this->runTime_.timeName(),
62  this->mesh_
63  ),
64  this->mesh_,
66  (
67  "nonlinearStress",
69  Zero
70  )
71  )
72 {}
73 
74 
75 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76 
77 template<class BasicMomentumTransportModel>
80 {
82  (
84  );
85  tR.ref() += nonlinearStress_;
86  return tR;
87 }
88 
89 
90 template<class BasicMomentumTransportModel>
93 {
95  (
97  );
98  tdevTau.ref() += this->rho_*nonlinearStress_;
99  return tdevTau;
100 }
101 
102 
103 template<class BasicMomentumTransportModel>
106 (
107  volVectorField& U
108 ) const
109 {
110  return
111  (
112  fvc::div(this->rho_*nonlinearStress_)
114  );
115 }
116 
117 
118 template<class BasicMomentumTransportModel>
121 (
122  const volScalarField& rho,
123  volVectorField& U
124 ) const
125 {
126  return
127  (
128  fvc::div(rho*nonlinearStress_)
130  );
131 }
132 
133 
134 // ************************************************************************* //
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 ::rhoField rhoField
Definition: eddyViscosity.H:71
phi
Definition: pEqn.H:104
virtual tmp< volSymmTensorField > devTau() const
Return the effective stress tensor.
incompressible::RASModel ::transportModel transportModel
Definition: eddyViscosity.H:72
A class for handling words, derived from string.
Definition: word.H:59
incompressible::RASModel ::alphaField alphaField
Definition: eddyViscosity.H:70
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the source term for the momentum equation.
static const zero Zero
Definition: zero.H:97
virtual tmp< volSymmTensorField > sigma() const
Return the Reynolds stress tensor [m^2/s^2].
U
Definition: pEqn.H:72
nonlinearEddyViscosity(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Construct from components.
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