HrenyaSinclairViscosity.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 -------------------------------------------------------------------------------
10 License
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 
29 #include "mathematicalConstants.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace kineticTheoryModels
37 {
38 namespace viscosityModels
39 {
40  defineTypeNameAndDebug(HrenyaSinclair, 0);
41 
43  (
44  viscosityModel,
45  HrenyaSinclair,
46  dictionary
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 
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
76 (
77  const volScalarField& alpha1,
78  const volScalarField& Theta,
79  const volScalarField& g0,
80  const volScalarField& rho1,
81  const volScalarField& da,
82  const dimensionedScalar& e
83 ) const
84 {
85  const scalar sqrtPi = sqrt(constant::mathematical::pi);
86 
87  const volScalarField lamda(1 + da/(6*sqrt(2.0)*(alpha1 + 1e-5))/L_);
88 
89  return da*sqrt(Theta)*
90  (
91  (4.0/5.0)*sqr(alpha1)*g0*(1 + e)/sqrtPi
92  + (1.0/15.0)*sqrtPi*g0*(1 + e)*(3*e - 1)*sqr(alpha1)/(3 - e)
93  + (1.0/6.0)*sqrtPi*alpha1*(0.5*lamda + 0.25*(3*e - 1))
94  /(0.5*(3 - e)*lamda)
95  + (10.0/96.0)*sqrtPi/((1 + e)*0.5*(3 - e)*g0*lamda)
96  );
97 }
98 
99 
101 {
102  coeffDict_ <<= dict_.optionalSubDict(typeName + "Coeffs");
103 
104  L_.readIfPresent(coeffDict_);
105 
106  return true;
107 }
108 
109 
110 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField & rho1
addToRunTimeSelectionTable(viscosityModel, Gidaspow, dictionary)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
constexpr scalar pi(M_PI)
tmp< volScalarField > nu(const volScalarField &alpha1, const volScalarField &Theta, const volScalarField &g0, const volScalarField &rho1, const volScalarField &da, const dimensionedScalar &e) const
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
HrenyaSinclair(const dictionary &dict)
Construct from components.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Namespace for OpenFOAM.
const volScalarField & alpha1