GidaspowConductivity.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 "GidaspowConductivity.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(Gidaspow, 0);
39 
41  (
42  conductivityModel,
43  Gidaspow,
44  dictionary
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const dictionary& dict
56 )
57 :
58  conductivityModel(dict)
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
72 (
73  const volScalarField& alpha1,
74  const volScalarField& Theta,
75  const volScalarField& g0,
76  const volScalarField& rho1,
77  const volScalarField& da,
78  const dimensionedScalar& e
79 ) const
80 {
81  const scalar sqrtPi = sqrt(constant::mathematical::pi);
82 
83  return rho1*da*sqrt(Theta)*
84  (
85  2.0*sqr(alpha1)*g0*(1.0 + e)/sqrtPi
86  + (9.0/8.0)*sqrtPi*g0*0.5*(1.0 + e)*sqr(alpha1)
87  + (15.0/16.0)*sqrtPi*alpha1
88  + (25.0/64.0)*sqrtPi/((1.0 + e)*g0)
89  );
90 }
91 
92 
93 // ************************************************************************* //
dictionary dict
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Gidaspow(const dictionary &dict)
Construct from components.
const volScalarField & alpha1
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar & rho1
Definition: createFields.H:39
Namespace for OpenFOAM.
tmp< volScalarField > kappa(const volScalarField &alpha1, const volScalarField &Theta, const volScalarField &g0, const volScalarField &rho1, const volScalarField &da, const dimensionedScalar &e) const