linearViscousStress.H
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-2016 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::linearViscousStress
29
30Group
31 grpTurbulence
32
33Description
34 Linear viscous stress turbulence model base class
35
36SourceFiles
37 linearViscousStress.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef linearViscousStress_H
42#define linearViscousStress_H
43
44// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45
46namespace Foam
47{
48
49/*---------------------------------------------------------------------------*\
50 Class linearViscousStress Declaration
51\*---------------------------------------------------------------------------*/
52
53template<class BasicTurbulenceModel>
55:
56 public BasicTurbulenceModel
57{
58
59public:
61 typedef typename BasicTurbulenceModel::alphaField alphaField;
62 typedef typename BasicTurbulenceModel::rhoField rhoField;
63 typedef typename BasicTurbulenceModel::transportModel transportModel;
64
65
66 // Constructors
67
68 //- Construct from components
70 (
71 const word& modelName,
72 const alphaField& alpha,
73 const rhoField& rho,
74 const volVectorField& U,
75 const surfaceScalarField& alphaRhoPhi,
77 const transportModel& transport,
78 const word& propertiesName
79 );
80
81
82 //- Destructor
83 virtual ~linearViscousStress() = default;
84
85
86 // Member Functions
87
88 //- Re-read model coefficients if they have changed
89 virtual bool read() = 0;
90
91 //- Return the effective stress tensor
92 virtual tmp<volSymmTensorField> devRhoReff() const;
93
94 //- Return the effective stress tensor based on a given velocity field
96 (
97 const volVectorField& U
98 ) const;
99
100 //- Return the source term for the momentum equation
102
103 //- Return the source term for the momentum equation
105 (
106 const volScalarField& rho,
108 ) const;
109
110 //- Solve the turbulence equations and correct the turbulence viscosity
111 virtual void correct() = 0;
112};
113
114
115// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
116
117} // End namespace Foam
118
119// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
120
121#ifdef NoRepository
122 #include "linearViscousStress.C"
123#endif
124
125// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
126
127#endif
128
129// ************************************************************************* //
surfaceScalarField & phi
Linear viscous stress turbulence model base class.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
virtual ~linearViscousStress()=default
Destructor.
virtual bool read()=0
Re-read model coefficients if they have changed.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
volScalarField & alpha