powerSeriesReactionRateI.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-------------------------------------------------------------------------------
10License
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 beta,
34 const scalar Ta,
36)
37:
38 A_(A),
39 beta_(beta),
40 Ta_(Ta),
41 coeffs_(coeffs)
42{}
43
44
46(
47 const speciesTable&,
48 const dictionary& dict
49)
50:
51 A_(dict.get<scalar>("A")),
52 beta_(dict.get<scalar>("beta")),
53 Ta_(dict.get<scalar>("Ta")),
54 coeffs_(dict.lookup("coeffs"))
55{}
56
57
58// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
59
61(
62 const scalar p,
63 const scalar T,
64 const scalarField&
65) const
66{
67 scalar lta = A_;
68
69 if (mag(beta_) > VSMALL)
70 {
71 lta *= pow(T, beta_);
72 }
73
74 scalar expArg = 0.0;
75
76 forAll(coeffs_, n)
77 {
78 expArg += coeffs_[n]/pow(T, n + 1);
79 }
80
81 lta *= exp(expArg);
82
83 return lta;
84}
85
86
88{
89 os.writeEntry("A", A_);
90 os.writeEntry("beta", beta_);
91 os.writeEntry("Ta", Ta_);
92 os.writeEntry("coeffs", coeffs_);
93}
94
95
96inline Foam::Ostream& Foam::operator<<
97(
98 Ostream& os,
99 const powerSeriesReactionRate& psrr
100)
101{
102 psrr.write(os);
103 return os;
104}
105
106
107// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
label n
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: FixedList.H:81
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
virtual bool write()
Write the output fields.
A wordList with hashed named lookup, which can be faster in some situations than using the normal lis...
Power series reaction rate.
Lookup type of boundary radiation properties.
Definition: lookup.H:66
volScalarField & p
const volScalarField & T
OBJstream os(runTime.globalPath()/outputName)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dictionary dict
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333