TroeFallOffFunctionI.H
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-2017 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 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 (
32  const scalar alpha,
33  const scalar Tsss,
34  const scalar Ts,
35  const scalar Tss
36 )
37 :
38  alpha_(alpha),
39  Tsss_(Tsss),
40  Ts_(Ts),
41  Tss_(Tss)
42 {}
43 
44 
46 :
47  alpha_(dict.get<scalar>("alpha")),
48  Tsss_(dict.get<scalar>("Tsss")),
49  Ts_(dict.get<scalar>("Ts")),
50  Tss_(dict.get<scalar>("Tss"))
51 {}
52 
53 
54 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
55 
56 inline Foam::scalar Foam::TroeFallOffFunction::operator()
57 (
58  const scalar T,
59  const scalar Pr
60 ) const
61 {
62  scalar logFcent = log10
63  (
64  max
65  (
66  (1 - alpha_)*exp(-T/Tsss_) + alpha_*exp(-T/Ts_) + exp(-Tss_/T),
67  SMALL
68  )
69  );
70 
71  scalar c = -0.4 - 0.67*logFcent;
72  static const scalar d = 0.14;
73  scalar n = 0.75 - 1.27*logFcent;
74 
75  scalar logPr = log10(max(Pr, SMALL));
76  return pow(10.0, logFcent/(1.0 + sqr((logPr + c)/(n - d*(logPr + c)))));
77 }
78 
79 
81 {
82  os.writeEntry("alpha", alpha_);
83  os.writeEntry("Tsss", Tsss_);
84  os.writeEntry("Ts", Ts_);
85  os.writeEntry("Tss", Tss_);
86 }
87 
88 
89 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
90 
91 inline Foam::Ostream& Foam::operator<<
92 (
94  const Foam::TroeFallOffFunction& tfof
95 )
96 {
97  tfof.write(os);
98  return os;
99 }
100 
101 
102 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
The Troe fall-off function.
void write(Ostream &os) const
Write to stream.
dimensionedScalar exp(const dimensionedScalar &ds)
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
TroeFallOffFunction(const scalar alpha, const scalar Tsss, const scalar Ts, const scalar Tss)
Construct from components.
OBJstream os(runTime.globalPath()/outputName)
const volScalarField & T
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionedScalar c
Speed of light in a vacuum.
label n
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
dimensionedScalar log10(const dimensionedScalar &ds)