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 | Copyright (C) 2011-2015 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 "HrenyaSinclairConductivity.H"
27 #include "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace kineticTheoryModels
35 {
36 namespace conductivityModels
37 {
38  defineTypeNameAndDebug(HrenyaSinclair, 0);
39 
41  (
42  conductivityModel,
43  HrenyaSinclair,
44  dictionary
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const dictionary& dict
56 )
57 :
58  conductivityModel(dict),
59  coeffDict_(dict.subDict(typeName + "Coeffs")),
60  L_("L", dimensionSet(0, 1, 0, 0, 0), coeffDict_)
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
65 
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
75 (
76  const volScalarField& alpha1,
77  const volScalarField& Theta,
78  const volScalarField& g0,
79  const volScalarField& rho1,
80  const volScalarField& da,
81  const dimensionedScalar& e
82 ) const
83 {
84  const scalar sqrtPi = sqrt(constant::mathematical::pi);
85 
86  volScalarField lamda
87  (
88  scalar(1) + da/(6.0*sqrt(2.0)*(alpha1 + scalar(1.0e-5)))/L_
89  );
90 
91  return rho1*da*sqrt(Theta)*
92  (
93  2.0*sqr(alpha1)*g0*(1.0 + e)/sqrtPi
94  + (9.0/8.0)*sqrtPi*g0*0.25*sqr(1.0 + e)*(2.0*e - 1.0)*sqr(alpha1)
95  /(49.0/16.0 - 33.0*e/16.0)
96  + (15.0/16.0)*sqrtPi*alpha1*(0.5*sqr(e) + 0.25*e - 0.75 + lamda)
97  /((49.0/16.0 - 33.0*e/16.0)*lamda)
98  + (25.0/64.0)*sqrtPi
99  /((1.0 + e)*(49.0/16.0 - 33.0*e/16.0)*lamda*g0)
100  );
101 }
102 
103 
105 {
106  coeffDict_ <<= dict_.subDict(typeName + "Coeffs");
107 
108  L_.readIfPresent(coeffDict_);
109 
110  return true;
111 }
112 
113 
114 // ************************************************************************* //
dictionary dict
const double e
Elementary charge.
Definition: doubleFloat.H:78
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
dimensionedScalar sqrt(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
tmp< volScalarField > kappa(const volScalarField &alpha1, const volScalarField &Theta, const volScalarField &g0, const volScalarField &rho1, const volScalarField &da, const dimensionedScalar &e) const
HrenyaSinclair(const dictionary &dict)
Construct from components.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:54
Namespace for OpenFOAM.