All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2020 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 BasicMomentumTransportModel>
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 )
52 :
54  (
55  typeName,
56  alpha,
57  rho,
58  U,
59  alphaRhoPhi,
60  phi,
61  transport
62  )
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
68 template<class BasicMomentumTransportModel>
69 const dictionary&
71 {
72  return dictionary::null;
73 }
74 
75 
76 template<class BasicMomentumTransportModel>
78 {
79  return true;
80 }
81 
82 
83 template<class BasicMomentumTransportModel>
86 {
87  return volScalarField::New
88  (
89  IOobject::groupName("nut", this->alphaRhoPhi_.group()),
90  this->mesh_,
92  );
93 }
94 
95 
96 template<class BasicMomentumTransportModel>
99 (
100  const label patchi
101 ) const
102 {
103  return tmp<scalarField>
104  (
105  new scalarField(this->mesh_.boundary()[patchi].size(), 0.0)
106  );
107 }
108 
109 
110 template<class BasicMomentumTransportModel>
113 {
114  return volScalarField::New
115  (
116  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
117  this->nu()
118  );
119 }
120 
121 
122 template<class BasicMomentumTransportModel>
125 (
126  const label patchi
127 ) const
128 {
129  return this->nu(patchi);
130 }
131 
132 
133 template<class BasicMomentumTransportModel>
136 {
137  return volScalarField::New
138  (
139  IOobject::groupName("k", this->alphaRhoPhi_.group()),
140  this->mesh_,
141  dimensionedScalar(sqr(this->U_.dimensions()), 0)
142  );
143 }
144 
145 
146 template<class BasicMomentumTransportModel>
149 {
150  return volScalarField::New
151  (
152  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
153  this->mesh_,
154  dimensionedScalar(sqr(this->U_.dimensions())/dimTime, 0)
155  );
156 }
157 
158 
159 template<class BasicMomentumTransportModel>
162 {
164  (
165  IOobject::groupName("R", this->alphaRhoPhi_.group()),
166  this->mesh_,
167  dimensionedSymmTensor(sqr(this->U_.dimensions()), Zero)
168  );
169 }
170 
171 
172 template<class BasicMomentumTransportModel>
174 {
176 }
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 } // End namespace laminarModels
182 } // End namespace Foam
183 
184 // ************************************************************************* //
virtual tmp< volSymmTensorField > sigma() const
Return the stress tensor [m^2/s^2], i.e. 0 for Stokes flow.
Definition: Stokes.C:161
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:112
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:241
virtual bool read()
Read momentumTransport dictionary.
Definition: Stokes.C:77
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
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:70
BasicMomentumTransportModel::rhoField rhoField
Definition: laminarModel.H:77
BasicMomentumTransportModel::alphaField alphaField
Definition: laminarModel.H:76
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate,.
Definition: Stokes.C:148
BasicMomentumTransportModel::transportModel transportModel
Definition: laminarModel.H:78
phi
Definition: pEqn.H:104
Stokes(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Construct from components.
Definition: Stokes.C:44
Calculate the gradient of the given field.
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.
label patchi
U
Definition: pEqn.H:72
virtual void correct()
Correct the laminar transport.
Definition: laminarModel.C:271
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< volScalarField > k() const
Return the turbulence kinetic energy, i.e. 0 for Stokes flow.
Definition: Stokes.C:135
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:173
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:85
Namespace for OpenFOAM.