Stokes.H
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) 2013-2018 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 Class
25  Foam::laminarModels::Stokes
26 
27 Description
28  Turbulence model for Stokes flow.
29 
30 SourceFiles
31  Stokes.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef Stokes_H
36 #define Stokes_H
37 
38 #include "laminarModel.H"
39 #include "linearViscousStress.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace laminarModels
46 {
47 
48 /*---------------------------------------------------------------------------* \
49  Class Stokes Declaration
50 \*---------------------------------------------------------------------------*/
51 
52 template<class BasicTurbulenceModel>
53 class Stokes
54 :
55  public linearViscousStress<laminarModel<BasicTurbulenceModel>>
56 {
57 
58 public:
59 
60  typedef typename BasicTurbulenceModel::alphaField alphaField;
61  typedef typename BasicTurbulenceModel::rhoField rhoField;
62  typedef typename BasicTurbulenceModel::transportModel transportModel;
63 
64 
65  //- Runtime type information
66  TypeName("Stokes");
67 
68 
69  // Constructors
70 
71  //- Construct from components
72  Stokes
73  (
74  const alphaField& alpha,
75  const rhoField& rho,
76  const volVectorField& U,
77  const surfaceScalarField& alphaRhoPhi,
78  const surfaceScalarField& phi,
79  const transportModel& transport,
80  const word& propertiesName = turbulenceModel::propertiesName
81  );
82 
83 
84  // Selectors
85 
86  //- Return a reference to the selected turbulence model
87  static autoPtr<Stokes> New
88  (
89  const alphaField& alpha,
90  const rhoField& rho,
91  const volVectorField& U,
92  const surfaceScalarField& alphaRhoPhi,
93  const surfaceScalarField& phi,
94  const transportModel& transport,
95  const word& propertiesName = turbulenceModel::propertiesName
96  );
97 
98 
99  //- Destructor
100  virtual ~Stokes()
101  {}
102 
103 
104  // Member Functions
105 
106  //- Const access to the coefficients dictionary
107  virtual const dictionary& coeffDict() const;
108 
109  //- Read turbulenceProperties dictionary
110  virtual bool read();
111 
112  //- Return the turbulence viscosity, i.e. 0 for Stokes flow
113  virtual tmp<volScalarField> nut() const;
114 
115  //- Return the turbulence viscosity on patch
116  virtual tmp<scalarField> nut(const label patchi) const;
117 
118  //- Return the effective viscosity, i.e. the Stokes viscosity
119  virtual tmp<volScalarField> nuEff() const;
120 
121  //- Return the effective viscosity on patch
122  virtual tmp<scalarField> nuEff(const label patchi) const;
123 
124  //- Return the turbulence kinetic energy, i.e. 0 for Stokes flow
125  virtual tmp<volScalarField> k() const;
126 
127  //- Return the turbulence kinetic energy dissipation rate,
128  // i.e. 0 for Stokes flow
129  virtual tmp<volScalarField> epsilon() const;
130 
131  //- Return the Reynolds stress tensor, i.e. 0 for Stokes flow
132  virtual tmp<volSymmTensorField> R() const;
133 
134  //- Correct the Stokes viscosity
135  virtual void correct();
136 };
137 
138 
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140 
141 } // End namespace laminarModels
142 } // End namespace Foam
143 
144 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
145 
146 #ifdef NoRepository
147  #include "Stokes.C"
148 #endif
149 
150 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
151 
152 #endif
153 
154 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity, i.e. the Stokes viscosity.
Definition: Stokes.C:125
surfaceScalarField & phi
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual bool read()
Read turbulenceProperties dictionary.
Definition: Stokes.C:79
Linear viscous stress turbulence model base class.
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: Stokes.C:72
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate,.
Definition: Stokes.C:174
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
TypeName("Stokes")
Runtime type information.
virtual ~Stokes()
Destructor.
Definition: Stokes.H:99
Stokes(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName)
Construct from components.
Definition: Stokes.C:44
label patchi
U
Definition: pEqn.H:72
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.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor, i.e. 0 for Stokes flow.
Definition: Stokes.C:201
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy, i.e. 0 for Stokes flow.
Definition: Stokes.C:150
BasicTurbulenceModel::rhoField rhoField
Definition: Stokes.H:60
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Turbulence model for Stokes flow.
Definition: Stokes.H:52
virtual void correct()
Correct the Stokes viscosity.
Definition: Stokes.C:227
BasicTurbulenceModel::alphaField alphaField
Definition: Stokes.H:59
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity, i.e. 0 for Stokes flow.
Definition: Stokes.C:87
BasicTurbulenceModel::transportModel transportModel
Definition: Stokes.H:61
Namespace for OpenFOAM.