HrenyaSinclairViscosity.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 "HrenyaSinclairViscosity.H"
27 #include "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace kineticTheoryModels
35 {
36 namespace viscosityModels
37 {
38  defineTypeNameAndDebug(HrenyaSinclair, 0);
39 
41  (
42  viscosityModel,
43  HrenyaSinclair,
44  dictionary
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const dictionary& dict
56 )
57 :
58  viscosityModel(dict),
59  coeffDict_(dict.optionalSubDict(typeName + "Coeffs")),
60  L_("L", dimensionSet(0, 1, 0, 0, 0), coeffDict_)
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
65 
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
74 (
75  const volScalarField& alpha1,
76  const volScalarField& Theta,
77  const volScalarField& g0,
78  const volScalarField& rho1,
79  const volScalarField& da,
80  const dimensionedScalar& e
81 ) const
82 {
83  const scalar sqrtPi = sqrt(constant::mathematical::pi);
84 
85  const volScalarField lamda(1 + da/(6*sqrt(2.0)*(alpha1 + 1e-5))/L_);
86 
87  return da*sqrt(Theta)*
88  (
89  (4.0/5.0)*sqr(alpha1)*g0*(1 + e)/sqrtPi
90  + (1.0/15.0)*sqrtPi*g0*(1 + e)*(3*e - 1)*sqr(alpha1)/(3 - e)
91  + (1.0/6.0)*sqrtPi*alpha1*(0.5*lamda + 0.25*(3*e - 1))
92  /(0.5*(3 - e)*lamda)
93  + (10.0/96.0)*sqrtPi/((1 + e)*0.5*(3 - e)*g0*lamda)
94  );
95 }
96 
97 
99 {
100  coeffDict_ <<= dict_.optionalSubDict(typeName + "Coeffs");
101 
102  L_.readIfPresent(coeffDict_);
103 
104  return true;
105 }
106 
107 
108 // ************************************************************************* //
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
tmp< volScalarField > nu(const volScalarField &alpha1, const volScalarField &Theta, const volScalarField &g0, const volScalarField &rho1, const volScalarField &da, const dimensionedScalar &e) const
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:759
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.
HrenyaSinclair(const dictionary &dict)
Construct from components.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.