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 | 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 a,
33  const scalar b,
34  const scalar c,
35  const scalar d,
36  const scalar e
37 )
38 :
39  a_(a),
40  b_(b),
41  c_(c),
42  d_(d),
43  e_(e)
44 {}
45 
46 
48 :
49  a_(dict.get<scalar>("a")),
50  b_(dict.get<scalar>("b")),
51  c_(dict.get<scalar>("c")),
52  d_(dict.get<scalar>("d")),
53  e_(dict.get<scalar>("e"))
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
58 
59 inline Foam::scalar Foam::SRIFallOffFunction::operator()
60 (
61  const scalar T,
62  const scalar Pr
63 ) const
64 {
65  scalar X = 1.0/(1.0 + sqr(log10(max(Pr, SMALL))));
66  return d_*pow(a_*exp(-b_/T) + exp(-T/c_), X)*pow(T, e_);
67 }
68 
69 
70 inline void Foam::SRIFallOffFunction::write(Ostream& os) const
71 {
72  os.writeEntry("a", a_);
73  os.writeEntry("b", b_);
74  os.writeEntry("c", c_);
75  os.writeEntry("d", d_);
76  os.writeEntry("e", e_);
77 }
78 
79 
80 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
81 
82 inline Foam::Ostream& Foam::operator<<
83 (
85  const Foam::SRIFallOffFunction& srifof
86 )
87 {
88  srifof.write(os);
89  return os;
90 }
91 
92 
93 // ************************************************************************* //
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
The SRI fall-off function.
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)
void write(Ostream &os) const
Write to stream.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
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
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
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.
dimensionedScalar log10(const dimensionedScalar &ds)