PhaseIncompressibleTurbulenceModel.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 volScalarField& alpha,
36  const geometricOneField& rho,
37  const volVectorField& U,
38  const surfaceScalarField& alphaRhoPhi,
39  const surfaceScalarField& phi,
40  const TransportModel& transportModel,
41  const word& propertiesName
42 )
43 :
45  <
47  geometricOneField,
48  incompressibleTurbulenceModel,
49  TransportModel
50  >
51  (
52  alpha,
53  rho,
54  U,
55  alphaRhoPhi,
56  phi,
57  transportModel,
58  propertiesName
59  )
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
64 
65 template<class TransportModel>
68 (
69  const volScalarField& alpha,
70  const volVectorField& U,
71  const surfaceScalarField& alphaRhoPhi,
72  const surfaceScalarField& phi,
73  const TransportModel& transportModel,
74  const word& propertiesName
75 )
76 {
78  (
81  <
83  geometricOneField,
84  incompressibleTurbulenceModel,
85  TransportModel
86  >::New
87  (
88  alpha,
89  geometricOneField(),
90  U,
91  alphaRhoPhi,
92  phi,
93  transportModel,
94  propertiesName
95  ).ptr())
96  );
97 }
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
102 template<class TransportModel>
105 {
106  return tmp<volScalarField>
107  (
108  new volScalarField
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_,
119  dimensionedScalar("pPrimef", dimPressure, 0.0)
120  )
121  );
122 }
123 
124 
125 template<class TransportModel>
128 {
130  (
132  (
133  IOobject
134  (
135  IOobject::groupName("pPrimef", this->alphaRhoPhi_.group()),
136  this->runTime_.timeName(),
137  this->mesh_,
138  IOobject::NO_READ,
139  IOobject::NO_WRITE
140  ),
141  this->mesh_,
142  dimensionedScalar("pPrimef", dimPressure, 0.0)
143  )
144  );
145 }
146 
147 
148 template<class TransportModel>
151 {
152  return devRhoReff();
153 }
154 
155 
156 template<class TransportModel>
159 (
160  volVectorField& U
161 ) const
162 {
163  return divDevRhoReff(U);
164 }
165 
166 
167 template<class TransportModel>
170 devRhoReff() const
171 {
173 
174  return devReff();
175 }
176 
177 
178 template<class TransportModel>
182 (
183  volVectorField& U
184 ) const
185 {
187 
188  return divDevReff(U);
189 }
190 
191 
192 // ************************************************************************* //
surfaceScalarField & phi
Templated abstract base class for turbulence models.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
static autoPtr< PhaseIncompressibleTurbulenceModel > New(const alphaField &alpha, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const TransportModel &trasportModel, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected turbulence model.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Templated abstract base class for multiphase incompressible turbulence models.
A class for handling words, derived from string.
Definition: word.H:59
PhaseIncompressibleTurbulenceModel(const word &type, const alphaField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const TransportModel &trasportModel, const word &propertiesName)
Construct.
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const
Return the source term for the momentum equation.
const dimensionSet dimPressure
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor.
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
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:366
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure&#39;.