solidArrheniusReactionRateI.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-2012 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 Ta,
32  const scalar Tcrit
33  //const scalar nReact
34 )
35 :
36  A_(A),
37  Ta_(Ta),
38  Tcrit_(Tcrit)
39 {}
40 
41 
43 (
44  const speciesTable&,
45  Istream& is
46 )
47 :
48  A_(readScalar(is.readBegin("solidArrheniusReaction(Istream&)"))),
49  Ta_(readScalar(is)),
50  Tcrit_(readScalar(is))
51 {}
52 
53 
55 (
56  const speciesTable&,
57  const dictionary& dict
58 )
59 :
60  A_(readScalar(dict.lookup("A"))),
61  Ta_(readScalar(dict.lookup("Ta"))),
62  Tcrit_(readScalar(dict.lookup("Tcrit")))
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
68 inline Foam::scalar Foam::solidArrheniusReactionRate::operator()
69 (
70  const scalar,
71  const scalar T,
72  const scalarField&
73 ) const
74 {
75  scalar ak = A_;
76 
77  if (T < Tcrit_)
78  {
79  ak *= 0.0;
80  }
81  else
82  {
83  ak *= exp(-Ta_/T);
84  }
85 
86  return ak;
87 }
88 
89 
91 {
92  os.writeKeyword("A") << A_ << token::END_STATEMENT << nl;
93  os.writeKeyword("Ta") << Ta_ << token::END_STATEMENT << nl;
94  os.writeKeyword("Tcrit") << Tcrit_ << token::END_STATEMENT << nl;
95 }
96 
97 
98 inline Foam::Ostream& Foam::operator<<
99 (
100  Ostream& os,
101  const solidArrheniusReactionRate& arr
102 )
103 {
104  os << token::BEGIN_LIST
105  << arr.A_ << token::SPACE << arr.Ta_ << token::SPACE << arr.Tcrit_
106  << token::END_LIST;
107  return os;
108 }
109 
110 
111 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
Istream & readBegin(const char *funcName)
Definition: Istream.C:86
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
dimensionedScalar exp(const dimensionedScalar &ds)
void write(Ostream &os) const
Write to stream.
Arrhenius reaction rate for solids.
solidArrheniusReactionRate(const scalar A, const scalar Ta, const scalar Tcrit)
Construct from components.
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
const volScalarField & T
A wordList with hashed indices for faster lookup by name.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451