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-2024 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 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class BasicMomentumTransportModel>
33 (
34  const word& modelName,
35  const alphaField& alpha,
36  const rhoField& rho,
37  const volVectorField& U,
38  const surfaceScalarField& alphaRhoPhi,
39  const surfaceScalarField& phi,
40  const viscosity& viscosity
41 )
42 :
43  eddyViscosity<BasicMomentumTransportModel>
44  (
45  modelName,
46  alpha,
47  rho,
48  U,
49  alphaRhoPhi,
50  phi,
51  viscosity
52  ),
53 
54  nonlinearStress_
55  (
56  IOobject
57  (
58  this->groupName("nonlinearStress"),
59  this->runTime_.name(),
60  this->mesh_
61  ),
62  this->mesh_,
64  (
65  "nonlinearStress",
67  Zero
68  )
69  )
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
75 template<class BasicMomentumTransportModel>
78 {
80  (
82  );
83  tR.ref() += nonlinearStress_;
84  return tR;
85 }
86 
87 
88 template<class BasicMomentumTransportModel>
91 {
93  (
95  );
96 
97  tdevTau.ref() += fvc::dotInterpolate
98  (
99  this->mesh().nf(),
100  this->rho_*nonlinearStress_
101  );
102 
103  return tdevTau;
104 }
105 
106 
107 template<class BasicMomentumTransportModel>
108 template<class RhoFieldType>
111 (
112  const RhoFieldType& rho,
114 ) const
115 {
116  return
117  (
118  fvc::div(rho*nonlinearStress_)
120  );
121 }
122 
123 
124 template<class BasicMomentumTransportModel>
127 (
129 ) const
130 {
131  return DivDevTau(this->rho_, U);
132 }
133 
134 
135 template<class BasicMomentumTransportModel>
138 (
139  const volScalarField& rho,
141 ) const
142 {
143  return DivDevTau(rho, U);
144 }
145 
146 
147 // ************************************************************************* //
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
BasicMomentumTransportModel::alphaField alphaField
Definition: eddyViscosity.H:70
BasicMomentumTransportModel::rhoField rhoField
Definition: eddyViscosity.H:71
Eddy viscosity turbulence model with non-linear correction base class.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor [m^2/s^2].
virtual tmp< surfaceVectorField > devTau() const
Return the effective surface stress.
nonlinearEddyViscosity(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the source term for the momentum equation.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< SurfaceField< typename innerProduct< vector, Type >::type > > dotInterpolate(const surfaceVectorField &Sf, const VolField< Type > &tvf)
Interpolate field onto faces.
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
static const zero Zero
Definition: zero.H:97
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimVelocity