Stokes.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-2021 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::laminarModels::Stokes
29
30Group
31 grpLaminar
32
33Description
34 Turbulence model for Stokes flow.
35
36SourceFiles
37 Stokes.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef Stokes_H
42#define Stokes_H
43
44#include "laminarModel.H"
45#include "linearViscousStress.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49namespace Foam
50{
51namespace laminarModels
52{
53
54/*---------------------------------------------------------------------------* \
55 Class Stokes Declaration
56\*---------------------------------------------------------------------------*/
57
58template<class BasicTurbulenceModel>
59class Stokes
60:
61 public linearViscousStress<laminarModel<BasicTurbulenceModel>>
62{
63
64public:
66 typedef typename BasicTurbulenceModel::alphaField alphaField;
67 typedef typename BasicTurbulenceModel::rhoField rhoField;
68 typedef typename BasicTurbulenceModel::transportModel transportModel;
69
70
71 //- Runtime type information
72 TypeName("Stokes");
73
74
75 // Constructors
76
77 //- Construct from components
78 Stokes
79 (
80 const alphaField& alpha,
81 const rhoField& rho,
82 const volVectorField& U,
83 const surfaceScalarField& alphaRhoPhi,
85 const transportModel& transport,
86 const word& propertiesName = turbulenceModel::propertiesName
87 );
88
89
90 // Selectors
91
92 //- Return a reference to the selected turbulence model
93 static autoPtr<Stokes> New
94 (
95 const alphaField& alpha,
96 const rhoField& rho,
97 const volVectorField& U,
98 const surfaceScalarField& alphaRhoPhi,
100 const transportModel& transport,
101 const word& propertiesName = turbulenceModel::propertiesName
102 );
103
104
105 //- Destructor
106 virtual ~Stokes() = default;
107
108
109 // Member Functions
110
111 //- Const access to the coefficients dictionary
112 virtual const dictionary& coeffDict() const;
113
114 //- Read turbulenceProperties dictionary
115 virtual bool read();
116
117 //- Return the turbulence viscosity, i.e. 0 for Stokes flow
118 virtual tmp<volScalarField> nut() const;
119
120 //- Return the turbulence viscosity on patch
121 virtual tmp<scalarField> nut(const label patchi) const;
122
123 //- Return the effective viscosity, i.e. the Stokes viscosity
124 virtual tmp<volScalarField> nuEff() const;
125
126 //- Return the effective viscosity on patch
127 virtual tmp<scalarField> nuEff(const label patchi) const;
128
129 //- Correct the Stokes viscosity
130 virtual void correct();
131};
132
133
134// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
135
136} // End namespace laminarModels
137} // End namespace Foam
138
139// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140
141#ifdef NoRepository
142 #include "Stokes.C"
143#endif
144
145// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146
147#endif
148
149// ************************************************************************* //
surfaceScalarField & phi
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Turbulence model for Stokes flow.
Definition: Stokes.H:61
BasicTurbulenceModel::alphaField alphaField
Definition: Stokes.H:65
BasicTurbulenceModel::rhoField rhoField
Definition: Stokes.H:66
virtual void correct()
Correct the Stokes viscosity.
Definition: Stokes.C:149
static autoPtr< Stokes > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected turbulence model.
TypeName("Stokes")
Runtime type information.
virtual ~Stokes()=default
Destructor.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity, i.e. the Stokes viscosity.
Definition: Stokes.C:125
BasicTurbulenceModel::transportModel transportModel
Definition: Stokes.H:67
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity, i.e. 0 for Stokes flow.
Definition: Stokes.C:90
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: Stokes.C:75
virtual bool read()
Read turbulenceProperties dictionary.
Definition: Stokes.C:82
Linear viscous stress turbulence model base class.
A class for managing temporary objects.
Definition: tmp.H:65
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
volScalarField & alpha
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73