1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "SyamlalConductivity.H"
#include "mathematicalConstants.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace kineticTheoryModels
{
namespace conductivityModels
{
defineTypeNameAndDebug(Syamlal, 0);
addToRunTimeSelectionTable
(
conductivityModel,
Syamlal,
dictionary
);
}
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::kineticTheoryModels::conductivityModels::Syamlal::Syamlal
(
const dictionary& dict
)
:
conductivityModel(dict)
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::kineticTheoryModels::conductivityModels::Syamlal::~Syamlal()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::tmp<Foam::volScalarField>
Foam::kineticTheoryModels::conductivityModels::Syamlal::kappa
(
const volScalarField& alpha1,
const volScalarField& Theta,
const volScalarField& g0,
const volScalarField& rho1,
const volScalarField& da,
const dimensionedScalar& e
) const
{
const scalar sqrtPi = sqrt(constant::mathematical::pi);
return rho1*da*sqrt(Theta)*
(
2.0*sqr(alpha1)*g0*(1.0 + e)/sqrtPi
+ (9.0/8.0)*sqrtPi*g0*0.25*sqr(1.0 + e)*(2.0*e - 1.0)*sqr(alpha1)
/(49.0/16.0 - 33.0*e/16.0)
+ (15.0/32.0)*sqrtPi*alpha1/(49.0/16.0 - 33.0*e/16.0)
);
}
// ************************************************************************* //
|