CarnahanStarlingRadial.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) 2011-2018 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
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
35namespace kineticTheoryModels
36{
37namespace radialModels
38{
40
42 (
46 );
47}
48}
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
55(
56 const dictionary& dict
57)
58:
60{}
61
62
63// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
64
66{}
67
68
69// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70
73(
74 const volScalarField& alpha,
75 const dimensionedScalar& alphaMinFriction,
77) const
78{
79
80 return
81 1.0/(1 - alpha)
82 + 3*alpha/(2*sqr(1 - alpha))
83 + sqr(alpha)/(2*pow3(1 - alpha));
84}
85
86
89(
90 const volScalarField& alpha,
91 const dimensionedScalar& alphaMinFriction,
93) const
94{
95 return
96 2.5/sqr(1 - alpha)
97 + 4*alpha/pow3(1 - alpha)
98 + 1.5*sqr(alpha)/pow4(1 - alpha);
99}
100
101
102// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
tmp< volScalarField > g0prime(const volScalarField &alpha, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
Derivative of the radial distribution function.
tmp< volScalarField > g0(const volScalarField &alpha, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
Radial distribution function.
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
Namespace for OpenFOAM.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
volScalarField & alpha
dictionary dict
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)