DPMIncompressibleTurbulenceModel.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) 2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class TransportModel>
35 (
36  const word& type,
37  const volScalarField& alpha,
38  const geometricOneField& rho,
39  const volVectorField& U,
40  const surfaceScalarField& alphaRhoPhi,
41  const surfaceScalarField& phi,
42  const TransportModel& transportModel,
43  const word& propertiesName
44 )
45 :
46  TurbulenceModel
47  <
49  geometricOneField,
50  incompressibleTurbulenceModel,
51  TransportModel
52  >
53  (
54  alpha,
55  rho,
56  U,
57  alphaRhoPhi,
58  phi,
59  transportModel,
60  propertiesName
61  )
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
66 
67 template<class TransportModel>
70 (
71  const volScalarField& alpha,
72  const volVectorField& U,
73  const surfaceScalarField& alphaRhoPhi,
74  const surfaceScalarField& phi,
75  const TransportModel& transportModel,
76  const word& propertiesName
77 )
78 {
79  return autoPtr<DPMIncompressibleTurbulenceModel>
80  (
81  static_cast<DPMIncompressibleTurbulenceModel*>(
82  TurbulenceModel
83  <
85  geometricOneField,
86  incompressibleTurbulenceModel,
87  TransportModel
88  >::New
89  (
90  alpha,
91  geometricOneField(),
92  U,
93  alphaRhoPhi,
94  phi,
95  transportModel,
96  propertiesName
97  ).ptr())
98  );
99 }
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
104 template<class TransportModel>
107 {
109  (
110  IOobject
111  (
112  IOobject::groupName("pPrime", this->alphaRhoPhi_.group()),
113  this->runTime_.timeName(),
114  this->mesh_,
115  IOobject::NO_READ,
116  IOobject::NO_WRITE
117  ),
118  this->mesh_,
120  );
121 }
122 
123 
124 template<class TransportModel>
127 {
129  (
130  IOobject
131  (
132  IOobject::groupName("pPrimef", this->alphaRhoPhi_.group()),
133  this->runTime_.timeName(),
134  this->mesh_,
135  IOobject::NO_READ,
136  IOobject::NO_WRITE
137  ),
138  this->mesh_,
140  );
141 }
142 
143 
144 template<class TransportModel>
147 {
148  return devRhoReff();
149 }
150 
151 
152 template<class TransportModel>
155 (
156  const volVectorField& U
157 ) const
158 {
159  return devRhoReff(U);
160 }
161 
162 
163 template<class TransportModel>
166 (
168 ) const
169 {
170  return divDevRhoReff(U);
171 }
172 
173 
174 template<class TransportModel>
177 {
179 
180  return devReff();
181 }
182 
183 
184 template<class TransportModel>
187 (
188  const volVectorField& U
189 ) const
190 {
192 
193  return nullptr;
194 }
195 
196 
197 template<class TransportModel>
200 (
202 ) const
203 {
205 
206  return divDevReff(U);
207 }
208 
209 
210 // ************************************************************************* //
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const
Return the source term for the momentum equation.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure&#39;.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
type
Types of root.
Definition: Roots.H:52
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
const dimensionSet dimPressure
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor.
DPMIncompressibleTurbulenceModel(const word &type, const alphaField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const TransportModel &transportModel, const word &propertiesName)
Construct.
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static autoPtr< DPMIncompressibleTurbulenceModel > New(const alphaField &alpha, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const TransportModel &transportModel, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected turbulence model.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:666
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133