PhaseCompressibleTurbulenceModel.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) 2013-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
32template<class TransportModel>
35(
36 const word& type,
37 const volScalarField& alpha,
38 const volScalarField& rho,
39 const volVectorField& U,
40 const surfaceScalarField& alphaRhoPhi,
42 const transportModel& transport,
43 const word& propertiesName
44)
45:
47 <
52 >
53 (
54 alpha,
55 rho,
56 U,
57 alphaRhoPhi,
58 phi,
59 transport,
60 propertiesName
61 )
62{}
63
64
65// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
66
67template<class TransportModel>
70(
71 const volScalarField& alpha,
72 const volScalarField& rho,
73 const volVectorField& U,
74 const surfaceScalarField& alphaRhoPhi,
76 const transportModel& transport,
77 const word& propertiesName
78)
79{
81 (
84 <
89 >::New
90 (
92 rho,
93 U,
94 alphaRhoPhi,
95 phi,
96 transport,
97 propertiesName
98 ).ptr())
99 );
100}
101
102
103// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104
105template<class TransportModel>
108{
110 (
112 (
113 IOobject::groupName("pPrime", this->alphaRhoPhi_.group()),
114 this->runTime_.timeName(),
115 this->mesh_,
118 ),
119 this->mesh_,
121 );
122}
123
124
125template<class TransportModel>
128{
130 (
132 (
133 IOobject::groupName("pPrimef", this->alphaRhoPhi_.group()),
134 this->runTime_.timeName(),
135 this->mesh_,
138 ),
139 this->mesh_,
141 );
142}
143
144
145// ************************************************************************* //
surfaceScalarField & phi
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Templated abstract base class for multiphase compressible turbulence models.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Templated abstract base class for turbulence models.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Abstract base class for turbulence models (RAS, LES and laminar).
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
const dimensionSet dimPressure
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
volScalarField & alpha