IncompressibleTurbulenceModel.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 
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class TransportModel>
33 (
34  const word& type,
35  const geometricOneField& alpha,
36  const geometricOneField& rho,
37  const volVectorField& U,
38  const surfaceScalarField& alphaRhoPhi,
39  const surfaceScalarField& phi,
40  const TransportModel& transport,
41  const word& propertiesName
42 )
43 :
45  <
46  geometricOneField,
47  geometricOneField,
48  incompressibleTurbulenceModel,
49  TransportModel
50  >
51  (
52  alpha,
53  rho,
54  U,
55  alphaRhoPhi,
56  phi,
57  transport,
58  propertiesName
59  )
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
64 
65 template<class TransportModel>
68 (
69  const volVectorField& U,
70  const surfaceScalarField& phi,
71  const TransportModel& transport,
72  const word& propertiesName
73 )
74 {
76  (
77  static_cast<IncompressibleTurbulenceModel*>(
79  <
80  geometricOneField,
81  geometricOneField,
82  incompressibleTurbulenceModel,
83  TransportModel
84  >::New
85  (
86  geometricOneField(),
87  geometricOneField(),
88  U,
89  phi,
90  phi,
91  transport,
92  propertiesName
93  ).ptr())
94  );
95 }
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
100 template<class TransportModel>
103 {
104  return devRhoReff();
105 }
106 
107 
108 template<class TransportModel>
111 (
112  volVectorField& U
113 ) const
114 {
115  return divDevRhoReff(U);
116 }
117 
118 
119 template<class TransportModel>
122 devRhoReff() const
123 {
125 
126  return devReff();
127 }
128 
129 
130 template<class TransportModel>
134 (
135  volVectorField& U
136 ) const
137 {
139 
140  return divDevReff(U);
141 }
142 
143 
144 template<class TransportModel>
148 (
149  const volScalarField& rho,
150  volVectorField& U
151 ) const
152 {
154 
155  return divDevReff(U);
156 }
157 
158 
159 // ************************************************************************* //
surfaceScalarField & phi
IncompressibleTurbulenceModel(const word &type, const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const TransportModel &transport, const word &propertiesName)
Construct.
Templated abstract base class for turbulence models.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Templated abstract base class for single-phase incompressible turbulence models.
A class for handling words, derived from string.
Definition: word.H:59
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const
Return the source term for the momentum equation.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor.
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
static autoPtr< IncompressibleTurbulenceModel > New(const volVectorField &U, const surfaceScalarField &phi, const TransportModel &trasportModel, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected turbulence model.
U
Definition: pEqn.H:72
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366