Stokes.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 | Copyright (C) 2013-2016 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 "Stokes.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "fvcGrad.H"
30 #include "fvcDiv.H"
31 #include "fvmLaplacian.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace laminarModels
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class BasicTurbulenceModel>
44 (
45  const alphaField& alpha,
46  const rhoField& rho,
47  const volVectorField& U,
48  const surfaceScalarField& alphaRhoPhi,
49  const surfaceScalarField& phi,
50  const transportModel& transport,
51  const word& propertiesName
52 )
53 :
55  (
56  typeName,
57  alpha,
58  rho,
59  U,
60  alphaRhoPhi,
61  phi,
62  transport,
63  propertiesName
64  )
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
70 template<class BasicTurbulenceModel>
71 const dictionary&
73 {
74  return dictionary::null;
75 }
76 
77 
78 template<class BasicTurbulenceModel>
80 {
81  return true;
82 }
83 
84 
85 template<class BasicTurbulenceModel>
88 {
89  return tmp<volScalarField>
90  (
91  new volScalarField
92  (
93  IOobject
94  (
95  IOobject::groupName("nut", this->U_.group()),
96  this->runTime_.timeName(),
97  this->mesh_,
100  false
101  ),
102  this->mesh_,
103  dimensionedScalar("nut", dimViscosity, 0.0)
104  )
105  );
106 }
107 
108 
109 template<class BasicTurbulenceModel>
112 (
113  const label patchi
114 ) const
115 {
116  return tmp<scalarField>
117  (
118  new scalarField(this->mesh_.boundary()[patchi].size(), 0.0)
119  );
120 }
121 
122 
123 template<class BasicTurbulenceModel>
126 {
127  return tmp<volScalarField>
128  (
129  new volScalarField
130  (
131  IOobject::groupName("nuEff", this->U_.group()), this->nu()
132  )
133  );
134 }
135 
136 
137 template<class BasicTurbulenceModel>
140 (
141  const label patchi
142 ) const
143 {
144  return this->nu(patchi);
145 }
146 
147 
148 template<class BasicTurbulenceModel>
151 {
152  return tmp<volScalarField>
153  (
154  new volScalarField
155  (
156  IOobject
157  (
158  IOobject::groupName("k", this->U_.group()),
159  this->runTime_.timeName(),
160  this->mesh_,
163  false
164  ),
165  this->mesh_,
166  dimensionedScalar("k", sqr(this->U_.dimensions()), 0.0)
167  )
168  );
169 }
170 
171 
172 template<class BasicTurbulenceModel>
175 {
176  return tmp<volScalarField>
177  (
178  new volScalarField
179  (
180  IOobject
181  (
182  IOobject::groupName("epsilon", this->U_.group()),
183  this->runTime_.timeName(),
184  this->mesh_,
187  false
188  ),
189  this->mesh_,
191  (
192  "epsilon", sqr(this->U_.dimensions())/dimTime, 0.0
193  )
194  )
195  );
196 }
197 
198 
199 template<class BasicTurbulenceModel>
202 {
204  (
206  (
207  IOobject
208  (
209  IOobject::groupName("R", this->U_.group()),
210  this->runTime_.timeName(),
211  this->mesh_,
214  false
215  ),
216  this->mesh_,
218  (
219  "R", sqr(this->U_.dimensions()), Zero
220  )
221  )
222  );
223 }
224 
225 
226 template<class BasicTurbulenceModel>
228 {
230 }
231 
232 
233 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234 
235 } // End namespace laminarModels
236 } // End namespace Foam
237 
238 // ************************************************************************* //
Foam::surfaceFields.
surfaceScalarField & phi
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
U
Definition: pEqn.H:83
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static const dictionary null
Null dictionary.
Definition: dictionary.H:202
virtual bool read()
Read turbulenceProperties dictionary.
Definition: Stokes.C:79
const dimensionSet dimViscosity
Calculate the matrix for the laplacian of the field.
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
BasicTurbulenceModel::transportModel transportModel
Definition: laminarModel.H:89
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
BasicTurbulenceModel::rhoField rhoField
Definition: laminarModel.H:88
BasicTurbulenceModel::alphaField alphaField
Definition: laminarModel.H:87
Calculate the gradient of the given field.
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
static const zero Zero
Definition: zero.H:91
Calculate the divergence of the given field.
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
virtual void correct()
Correct the laminar transport.
Definition: laminarModel.C:326
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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
virtual void correct()
Correct the Stokes viscosity.
Definition: Stokes.C:227
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
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
volScalarField & nu
Namespace for OpenFOAM.