Loading...
Searching...
No Matches
nonlinearEddyViscosity.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) 2015-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
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 eddyViscosity<BasicTurbulenceModel>
48 (
49 modelName,
50 alpha,
51 rho,
52 U,
53 alphaRhoPhi,
54 phi,
55 transport,
56 propertiesName
57 ),
58
60 (
62 (
63 IOobject::groupName("nonlinearStress", alphaRhoPhi.group()),
64 this->runTime_.timeName(),
65 this->mesh_
66 ),
67 this->mesh_,
69 )
70{}
71
72
73// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74
75template<class BasicTurbulenceModel>
78{
80 (
82 );
83 tR.ref() += nonlinearStress_;
84 return tR;
85}
86
87
88template<class BasicTurbulenceModel>
92 return devRhoReff(this->U_);
93}
94
95
96template<class BasicTurbulenceModel>
99(
100 const volVectorField& U
101) const
102{
103 tmp<volSymmTensorField> tdevRhoReff
104 (
106 );
107 tdevRhoReff.ref() += this->rho_*nonlinearStress_;
109 return tdevRhoReff;
110}
112
113template<class BasicTurbulenceModel>
116(
118) const
119{
120 return
121 (
122 fvc::div(this->rho_*nonlinearStress_)
124 );
125}
126
127
128template<class BasicTurbulenceModel>
131(
132 const volScalarField& rho,
134) const
135{
136 return
137 (
138 fvc::div(rho*nonlinearStress_)
140 );
141}
142
143
144// ************************************************************************* //
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
eddyViscosity(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.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
BasicTurbulenceModel::alphaField alphaField
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.
BasicTurbulenceModel::rhoField rhoField
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
word timeName
Definition getTimeIndex.H:3
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
volScalarField & alpha