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-2017 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(dict.lookup("a"))),
48  b_(readScalar(dict.lookup("b"))),
49  c_(readScalar(dict.lookup("c"))),
50  d_(readScalar(dict.lookup("d"))),
51  e_(readScalar(dict.lookup("e")))
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
56 
57 inline Foam::scalar Foam::SRIFallOffFunction::operator()
58 (
59  const scalar T,
60  const scalar Pr
61 ) const
62 {
63  scalar X = 1.0/(1.0 + sqr(log10(max(Pr, SMALL))));
64  return d_*pow(a_*exp(-b_/T) + exp(-T/c_), X)*pow(T, e_);
65 }
66 
67 
69 {
70  os.writeKeyword("a") << a_ << token::END_STATEMENT << nl;
71  os.writeKeyword("b") << b_ << token::END_STATEMENT << nl;
72  os.writeKeyword("c") << c_ << token::END_STATEMENT << nl;
73  os.writeKeyword("d") << d_ << token::END_STATEMENT << nl;
74  os.writeKeyword("e") << e_ << token::END_STATEMENT << nl;
75 }
76 
77 
78 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
79 
80 inline Foam::Ostream& Foam::operator<<
81 (
82  Foam::Ostream& os,
83  const Foam::SRIFallOffFunction& srifof
84 )
85 {
86  srifof.write(os);
87  return os;
88 }
89 
90 
91 // ************************************************************************* //
#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)
void write(Ostream &os) const
Write to stream.
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
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)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionedScalar c
Speed of light in a vacuum.
virtual Ostream & write(const token &)=0
Write next token to stream.
dimensionedScalar log10(const dimensionedScalar &ds)