SRIFallOffFunctionI.H
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-2024 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 "SRIFallOffFunction.H"
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.lookup<scalar>("a")),
50  b_(dict.lookup<scalar>("b")),
51  c_(dict.lookup<scalar>("c")),
52  d_(dict.lookup<scalar>("d")),
53  e_(dict.lookup<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  const scalar logPr = log10(max(Pr, small));
66 
67  const scalar X = 1/(1 + sqr(logPr));
68 
69  const scalar psi = a_*exp(-b_/T) + exp(-T/c_);
70 
71  return d_*pow(psi, X)*pow(T, e_);
72 }
73 
74 
75 inline Foam::scalar Foam::SRIFallOffFunction::ddT
76 (
77  const scalar T,
78  const scalar Pr,
79  const scalar F
80 ) const
81 {
82  const scalar logPr = log10(max(Pr, small));
83 
84  const scalar X = 1/(1 + sqr(logPr));
85 
86  const scalar psi = a_*exp(-b_/T) + exp(-T/c_);
87  const scalar dpsidT = a_*b_/sqr(T)*exp(-b_/T) - 1/c_*exp(-T/c_);
88 
89  return F*(X/psi*dpsidT + e_/T);
90 }
91 
92 
93 inline Foam::scalar Foam::SRIFallOffFunction::ddPr
94 (
95  const scalar T,
96  const scalar Pr,
97  const scalar F
98 ) const
99 {
100  static const scalar logTen = log(scalar(10));
101 
102  const scalar logPr = log10(max(Pr, small));
103  const scalar dlogPrdPr = Pr >= small ? 1/(logTen*Pr) : 0;
104 
105  const scalar X = 1/(1 + sqr(logPr));
106  const scalar dXdPr = -sqr(X)*2*logPr*dlogPrdPr;
107 
108  const scalar psi = a_*exp(-b_/T) + exp(-T/c_);
109 
110  return F*log(psi)*dXdPr;
111 }
112 
113 
115 {
116  writeEntry(os, "a", a_);
117  writeEntry(os, "b", b_);
118  writeEntry(os, "c", c_);
119  writeEntry(os, "d", d_);
120  writeEntry(os, "e", e_);
121 }
122 
123 
124 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
125 
126 inline Foam::Ostream& Foam::operator<<
127 (
128  Foam::Ostream& os,
129  const Foam::SRIFallOffFunction& srifof
130 )
131 {
132  srifof.write(os);
133  return os;
134 }
135 
136 
137 // ************************************************************************* //
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual Ostream & write(const char)=0
Write character.
The SRI fall-off function.
scalar ddT(const scalar T, const scalar Pr, const scalar F) const
void write(Ostream &os) const
Write to stream.
SRIFallOffFunction(const scalar a, const scalar b, const scalar c, const scalar d, const scalar e)
Construct from components.
scalar ddPr(const scalar T, const scalar Pr, const scalar F) const
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const volScalarField & psi
volScalarField & b
Definition: createFields.H:25
const dimensionedScalar e
Elementary charge.
const dimensionedScalar F
Faraday constant: default SI units: [C/kmol].
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar log10(const dimensionedScalar &ds)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict