Loading...
Searching...
No Matches
LandauTellerReactionRateI.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,
35 const scalar B,
36 const scalar C
37)
38:
39 A_(A),
40 beta_(beta),
41 Ta_(Ta),
42 B_(B),
43 C_(C)
44{}
45
46
48(
49 const speciesTable&,
50 const dictionary& dict
51)
52:
53 A_(dict.get<scalar>("A")),
54 beta_(dict.get<scalar>("beta")),
55 Ta_(dict.get<scalar>("Ta")),
56 B_(dict.get<scalar>("B")),
57 C_(dict.get<scalar>("C"))
58{}
59
60
61// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62
63inline Foam::scalar Foam::LandauTellerReactionRate::operator()
64(
65 const scalar p,
66 const scalar T,
67 const scalarField&
68) const
69{
70 scalar lta = A_;
71
72 if (mag(beta_) > VSMALL)
73 {
74 lta *= pow(T, beta_);
75 }
76
77 scalar expArg = 0.0;
78
79 if (mag(Ta_) > VSMALL)
80 {
81 expArg -= Ta_/T;
82 }
83
84 if (mag(B_) > VSMALL)
85 {
86 expArg += B_/cbrt(T);
87 }
88
89 if (mag(C_) > VSMALL)
90 {
91 expArg += C_/pow(T, 2.0/3.0);
92 }
93
94 if (mag(expArg) > VSMALL)
95 {
96 lta *= exp(expArg);
97 }
98
99 return lta;
100}
101
102
103inline void Foam::LandauTellerReactionRate::write(Ostream& os) const
104{
105 os.writeEntry("A", A_);
106 os.writeEntry("beta", beta_);
107 os.writeEntry("Ta", Ta_);
108 os.writeEntry("B", B_);
109 os.writeEntry("C", C_);
110}
111
112
113inline Foam::Ostream& Foam::operator<<
114(
115 Ostream& os,
116 const LandauTellerReactionRate& arr
117)
118{
119 arr.write(os);
120 return os;
121}
122
123
124// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Graphite solid properties.
Definition C.H:49
void write(Ostream &os) const
Write to stream.
LandauTellerReactionRate(const scalar A, const scalar beta, const scalar Ta, const scalar B, const scalar C)
Construct from components.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
volScalarField & p
const volScalarField & T
OBJstream os(runTime.globalPath()/outputName)
dimensionedScalar exp(const dimensionedScalar &ds)
hashedWordList speciesTable
A table of species as a hashedWordList.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar cbrt(const dimensionedScalar &ds)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)