IncompressibleMomentumTransportModel.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 
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 )
42 :
44  <
45  geometricOneField,
46  geometricOneField,
47  incompressibleMomentumTransportModel,
48  TransportModel
49  >
50  (
51  alpha,
52  rho,
53  U,
54  alphaRhoPhi,
55  phi,
56  transport
57  )
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
62 
63 template<class TransportModel>
66 (
67  const volVectorField& U,
68  const surfaceScalarField& phi,
69  const TransportModel& transport
70 )
71 {
73  (
76  <
77  geometricOneField,
78  geometricOneField,
79  incompressibleMomentumTransportModel,
80  TransportModel
81  >::New
82  (
83  geometricOneField(),
84  geometricOneField(),
85  U,
86  phi,
87  phi,
88  transport
89  ).ptr())
90  );
91 }
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
96 template<class TransportModel>
99 {
100  return devTau();
101 }
102 
103 
104 template<class TransportModel>
107 (
108  volVectorField& U
109 ) const
110 {
111  return divDevTau(U);
112 }
113 
114 
115 template<class TransportModel>
118 devTau() const
119 {
121 
122  return devSigma();
123 }
124 
125 
126 template<class TransportModel>
129 divDevTau
130 (
131  volVectorField& U
132 ) const
133 {
135 
136  return divDevSigma(U);
137 }
138 
139 
140 template<class TransportModel>
143 divDevTau
144 (
145  const volScalarField& rho,
146  volVectorField& U
147 ) const
148 {
150 
151  return divDevSigma(U);
152 }
153 
154 
155 // ************************************************************************* //
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the source term for the momentum equation.
static autoPtr< IncompressibleMomentumTransportModel > New(const volVectorField &U, const surfaceScalarField &phi, const TransportModel &transportModel)
Return a reference to the selected turbulence model.
virtual tmp< fvVectorMatrix > divDevSigma(volVectorField &U) const
Return the source term for the momentum equation.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
virtual tmp< volSymmTensorField > devTau() const
Return the effective stress tensor.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
phi
Definition: correctPhi.H:3
IncompressibleMomentumTransportModel(const word &type, const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const TransportModel &transport)
Construct.
Templated abstract base class for single-phase incompressible turbulence models.
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
Templated abstract base class for turbulence models.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:370
virtual tmp< volSymmTensorField > devSigma() const
Return the effective stress tensor.