SRIFallOffFunctionI.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 | Copyright (C) 2011 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
27 
29 (
30  const scalar a,
31  const scalar b,
32  const scalar c,
33  const scalar d,
34  const scalar e
35 )
36 :
37  a_(a),
38  b_(b),
39  c_(c),
40  d_(d),
41  e_(e)
42 {}
43 
44 
46 :
47  a_(readScalar(is.readBegin("SRIFallOffFunction(Istream&)"))),
48  b_(readScalar(is)),
49  c_(readScalar(is)),
50  d_(readScalar(is)),
51  e_(readScalar(is))
52 {
53  is.readEnd("SRIFallOffFunction(Istream&)");
54 }
55 
56 
58 :
59  a_(readScalar(dict.lookup("a"))),
60  b_(readScalar(dict.lookup("b"))),
61  c_(readScalar(dict.lookup("c"))),
62  d_(readScalar(dict.lookup("d"))),
63  e_(readScalar(dict.lookup("e")))
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
69 inline Foam::scalar Foam::SRIFallOffFunction::operator()
70 (
71  const scalar T,
72  const scalar Pr
73 ) const
74 {
75  scalar X = 1.0/(1.0 + sqr(log10(max(Pr, SMALL))));
76  return d_*pow(a_*exp(-b_/T) + exp(-T/c_), X)*pow(T, e_);
77 }
78 
79 
81 {
82  os.writeKeyword("a") << a_ << token::END_STATEMENT << nl;
83  os.writeKeyword("b") << b_ << token::END_STATEMENT << nl;
84  os.writeKeyword("c") << c_ << token::END_STATEMENT << nl;
85  os.writeKeyword("d") << d_ << token::END_STATEMENT << nl;
86  os.writeKeyword("e") << e_ << token::END_STATEMENT << nl;
87 }
88 
89 
90 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
91 
92 inline Foam::Ostream& Foam::operator<<
93 (
94  Foam::Ostream& os,
95  const Foam::SRIFallOffFunction& srifof
96 )
97 {
98  os << token::BEGIN_LIST
99  << srifof.a_
100  << token::SPACE << srifof.b_
101  << token::SPACE << srifof.c_
102  << token::SPACE << srifof.d_
103  << token::SPACE << srifof.e_
104  << token::END_LIST;
105 
106  return os;
107 }
108 
109 
110 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
dimensionedScalar Pr("Pr", dimless, laminarTransport)
dictionary dict
The SRI fall-off function.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Istream & readEnd(const char *funcName)
Definition: Istream.C:103
stressControl lookup("compactNormalStress") >> compactNormalStress
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:78
SRIFallOffFunction(const scalar a, const scalar b, const scalar c, const scalar d, const scalar e)
Construct from components.
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void write(Ostream &os) const
Write to stream.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar log10(const dimensionedScalar &ds)