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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2017 OpenFOAM Foundation
9  Copyright (C) 2020-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "Stokes.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "fvcGrad.H"
33 #include "fvcDiv.H"
34 #include "fvmLaplacian.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace laminarModels
41 {
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 template<class BasicTurbulenceModel>
47 (
48  const alphaField& alpha,
49  const rhoField& rho,
50  const volVectorField& U,
51  const surfaceScalarField& alphaRhoPhi,
52  const surfaceScalarField& phi,
53  const transportModel& transport,
54  const word& propertiesName
55 )
56 :
58  (
59  typeName,
60  alpha,
61  rho,
62  U,
63  alphaRhoPhi,
64  phi,
65  transport,
66  propertiesName
67  )
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
73 template<class BasicTurbulenceModel>
74 const dictionary&
76 {
77  return dictionary::null;
78 }
79 
80 
81 template<class BasicTurbulenceModel>
83 {
84  return true;
85 }
86 
87 
88 template<class BasicTurbulenceModel>
91 {
92  return volScalarField::New
93  (
94  IOobject::groupName("nut", this->alphaRhoPhi_.group()),
96  this->mesh_,
98  );
99 }
100 
101 
102 template<class BasicTurbulenceModel>
105 (
106  const label patchi
107 ) const
108 {
109  return tmp<scalarField>::New(this->mesh_.boundary()[patchi].size(), Zero);
110 }
111 
112 
113 template<class BasicTurbulenceModel>
116 {
117  return volScalarField::New
118  (
119  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
121  this->nu()
122  );
123 }
124 
125 
126 template<class BasicTurbulenceModel>
129 (
130  const label patchi
131 ) const
132 {
133  return this->nu(patchi);
134 }
135 
136 
137 template<class BasicTurbulenceModel>
139 {
141 }
142 
143 
144 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
145 
146 } // End namespace laminarModels
147 } // End namespace Foam
148 
149 // ************************************************************************* //
Foam::surfaceFields.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity, i.e. the Stokes viscosity.
Definition: Stokes.C:108
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual bool read()
Read turbulenceProperties dictionary.
Definition: Stokes.C:75
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:68
BasicTurbulenceModel::transportModel transportModel
Definition: laminarModel.H:95
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
BasicTurbulenceModel::rhoField rhoField
Definition: laminarModel.H:94
BasicTurbulenceModel::alphaField alphaField
Definition: laminarModel.H:93
Calculate the gradient of the given field.
A class for handling words, derived from Foam::string.
Definition: word.H:63
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:486
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
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:40
U
Definition: pEqn.H:72
virtual void correct()
Correct the laminar transport.
Definition: laminarModel.C:287
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual void correct()
Correct the Stokes viscosity.
Definition: Stokes.C:131
A class for managing temporary objects.
Definition: HashPtrTable.H:50
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:83
volScalarField & nu
Do not request registration (bool: false)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127