CarnahanStarlingRadial.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) 2011-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 
26 #include "CarnahanStarlingRadial.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace kineticTheoryModels
34 {
35 namespace radialModels
36 {
37  defineTypeNameAndDebug(CarnahanStarling, 0);
38 
40  (
41  radialModel,
42  CarnahanStarling,
43  dictionary
44  );
45 }
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const dictionary& dict
55 )
56 :
57  radialModel(dict)
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
71 (
72  const volScalarField& alpha,
73  const dimensionedScalar& alphaMinFriction,
74  const dimensionedScalar& alphaMax
75 ) const
76 {
77 
78  return
79  1.0/(1 - alpha)
80  + 3*alpha/(2*sqr(1 - alpha))
81  + sqr(alpha)/(2*pow3(1 - alpha));
82 }
83 
84 
87 (
88  const volScalarField& alpha,
89  const dimensionedScalar& alphaMinFriction,
90  const dimensionedScalar& alphaMax
91 ) const
92 {
93  return
94  2.5/sqr(1 - alpha)
95  + 4*alpha/pow3(1 - alpha)
96  + 1.5*sqr(alpha)/pow4(1 - alpha);
97 }
98 
99 
100 // ************************************************************************* //
dictionary dict
tmp< volScalarField > g0(const volScalarField &alpha, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
Radial distribution function.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
tmp< volScalarField > g0prime(const volScalarField &alpha, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
Derivative of the radial distribution function.
CarnahanStarling(const dictionary &dict)
Construct from components.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pow4(const dimensionedScalar &ds)
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Namespace for OpenFOAM.