HrenyaSinclairConductivity.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
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace kineticTheoryModels
37{
38namespace conductivityModels
39{
41
43 (
47 );
48}
49}
50}
51
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
56(
57 const dictionary& dict
58)
59:
61 coeffDict_(dict.optionalSubDict(typeName + "Coeffs")),
62 L_("L", dimLength, coeffDict_)
63{}
64
65
66// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67
70{}
71
72
73// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74
77(
79 const volScalarField& Theta,
80 const volScalarField& g0,
81 const volScalarField& rho1,
82 const volScalarField& da,
83 const dimensionedScalar& e
84) const
85{
86 const scalar sqrtPi = sqrt(constant::mathematical::pi);
87
88 volScalarField lamda
89 (
90 scalar(1) + da/(6*sqrt(2.0)*(alpha1 + 1e-5))/L_
91 );
92
93 return rho1*da*sqrt(Theta)*
94 (
95 2*sqr(alpha1)*g0*(1 + e)/sqrtPi
96 + (9.0/8.0)*sqrtPi*g0*0.25*sqr(1 + e)*(2*e - 1)*sqr(alpha1)
97 /(49.0/16.0 - 33*e/16.0)
98 + (15.0/16.0)*sqrtPi*alpha1*(0.5*sqr(e) + 0.25*e - 0.75 + lamda)
99 /((49.0/16.0 - 33*e/16.0)*lamda)
100 + (25.0/64.0)*sqrtPi
101 /((1 + e)*(49.0/16.0 - 33*e/16.0)*lamda*g0)
102 );
103}
104
105
107{
108 coeffDict_ <<= dict_.optionalSubDict(typeName + "Coeffs");
109
110 L_.readIfPresent(coeffDict_);
111
112 return true;
113}
114
115
116// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const volScalarField & alpha1
volScalarField & rho1
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
constexpr scalar pi(M_PI)
Namespace for OpenFOAM.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar sqrt(const dimensionedScalar &ds)
dictionary dict
volScalarField & e
Definition: createFields.H:11