Stokes.C
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 \*---------------------------------------------------------------------------*/
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 volScalarField::New
90  (
91  IOobject::groupName("nut", this->alphaRhoPhi_.group()),
92  this->mesh_,
94  );
95 }
96 
97 
98 template<class BasicTurbulenceModel>
101 (
102  const label patchi
103 ) const
104 {
105  return tmp<scalarField>
106  (
107  new scalarField(this->mesh_.boundary()[patchi].size(), 0.0)
108  );
109 }
110 
111 
112 template<class BasicTurbulenceModel>
115 {
116  return volScalarField::New
117  (
118  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
119  this->nu()
120  );
121 }
122 
123 
124 template<class BasicTurbulenceModel>
127 (
128  const label patchi
129 ) const
130 {
131  return this->nu(patchi);
132 }
133 
134 
135 template<class BasicTurbulenceModel>
138 {
139  return volScalarField::New
140  (
141  IOobject::groupName("k", this->alphaRhoPhi_.group()),
142  this->mesh_,
143  dimensionedScalar(sqr(this->U_.dimensions()), 0)
144  );
145 }
146 
147 
148 template<class BasicTurbulenceModel>
151 {
152  return volScalarField::New
153  (
154  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
155  this->mesh_,
156  dimensionedScalar(sqr(this->U_.dimensions())/dimTime, 0)
157  );
158 }
159 
160 
161 template<class BasicTurbulenceModel>
164 {
166  (
167  IOobject::groupName("R", this->alphaRhoPhi_.group()),
168  this->mesh_,
169  dimensionedSymmTensor(sqr(this->U_.dimensions()), Zero)
170  );
171 }
172 
173 
174 template<class BasicTurbulenceModel>
176 {
178 }
179 
180 
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 
183 } // End namespace laminarModels
184 } // End namespace Foam
185 
186 // ************************************************************************* //
Foam::surfaceFields.
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:114
surfaceScalarField & phi
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static const dictionary null
Null dictionary.
Definition: dictionary.H:223
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:150
BasicTurbulenceModel::transportModel transportModel
Definition: laminarModel.H:78
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
BasicTurbulenceModel::rhoField rhoField
Definition: laminarModel.H:77
BasicTurbulenceModel::alphaField alphaField
Definition: laminarModel.H:76
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:97
Calculate the divergence of the given field.
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
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
virtual void correct()
Correct the laminar transport.
Definition: laminarModel.C:274
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:163
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy, i.e. 0 for Stokes flow.
Definition: Stokes.C:137
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
virtual void correct()
Correct the Stokes viscosity.
Definition: Stokes.C:175
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity, i.e. 0 for Stokes flow.
Definition: Stokes.C:87
volScalarField & nu
Namespace for OpenFOAM.