PhaseCompressibleMomentumTransportModel.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 volScalarField& alpha,
36  const volScalarField& rho,
37  const volVectorField& U,
38  const surfaceScalarField& alphaRhoPhi,
39  const surfaceScalarField& phi,
40  const transportModel& transport
41 )
42 :
44  <
47  compressibleMomentumTransportModel,
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 volScalarField& alpha,
68  const volScalarField& rho,
69  const volVectorField& U,
70  const surfaceScalarField& alphaRhoPhi,
71  const surfaceScalarField& phi,
72  const transportModel& transport
73 )
74 {
76  (
79  <
82  compressibleMomentumTransportModel,
84  >::New
85  (
86  alpha,
87  rho,
88  U,
89  alphaRhoPhi,
90  phi,
91  transport
92  ).ptr())
93  );
94 }
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
99 template<class TransportModel>
102 {
103  return volScalarField::New
104  (
105  IOobject::groupName("pPrime", this->alphaRhoPhi_.group()),
106  this->mesh_,
108  );
109 }
110 
111 
112 template<class TransportModel>
115 {
117  (
118  IOobject::groupName("pPrimef", this->alphaRhoPhi_.group()),
119  this->mesh_,
121  );
122 }
123 
124 
125 // ************************************************************************* //
Templated abstract base class for multiphase compressible turbulence models.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
static autoPtr< PhaseCompressibleMomentumTransportModel > New(const alphaField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transportModel)
Return a reference to the selected turbulence model.
phi
Definition: pEqn.H:104
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
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
const dimensionSet dimPressure
PhaseCompressibleMomentumTransportModel(const word &type, const alphaField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Construct.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure&#39;.
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
Templated abstract base class for turbulence models.