solidArrheniusReactionRateI.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-2018 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  const dictionary& dict
46 )
47 :
48  A_(readScalar(dict.lookup("A"))),
49  Ta_(readScalar(dict.lookup("Ta"))),
50  Tcrit_(readScalar(dict.lookup("Tcrit")))
51 {}
52 
53 
54 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
55 
56 inline Foam::scalar Foam::solidArrheniusReactionRate::operator()
57 (
58  const scalar,
59  const scalar T,
60  const scalarField&
61 ) const
62 {
63  if (T < Tcrit_)
64  {
65  return 0;
66  }
67  else
68  {
69  return A_*exp(-Ta_/T);
70  }
71 }
72 
73 
74 inline Foam::scalar Foam::solidArrheniusReactionRate::ddT
75 (
76  const scalar p,
77  const scalar T,
78  const scalarField&
79 ) const
80 {
81  if (T < Tcrit_)
82  {
83  return 0;
84  }
85  else
86  {
87  return A_*exp(-Ta_/T)*Ta_/sqr(T);
88  }
89 }
90 
91 
94 {
95  return NullObjectRef<List<Tuple2<label, scalar>>>();
96 }
97 
98 
100 (
101  const scalar p,
102  const scalar T,
103  const scalarField& c,
104  scalarField& dcidc
105 ) const
106 {}
107 
108 
109 inline Foam::scalar Foam::solidArrheniusReactionRate::dcidT
110 (
111  const scalar p,
112  const scalar T,
113  const scalarField& c
114 ) const
115 {
116  return 0;
117 }
118 
119 
121 {
122  os.writeKeyword("A") << A_ << token::END_STATEMENT << nl;
123  os.writeKeyword("Ta") << Ta_ << token::END_STATEMENT << nl;
124  os.writeKeyword("Tcrit") << Tcrit_ << token::END_STATEMENT << nl;
125 }
126 
127 
128 inline Foam::Ostream& Foam::operator<<
129 (
130  Ostream& os,
131  const solidArrheniusReactionRate& arr
132 )
133 {
134  arr.write(os);
135  return os;
136 }
137 
138 
139 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
const List< Tuple2< label, scalar > > & beta() const
Third-body efficiencies (beta = 1-alpha)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
scalar dcidT(const scalar p, const scalar T, const scalarField &c) const
Temperature derivative of the pressure dependent term.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void write(Ostream &os) const
Write to stream.
dimensionedScalar exp(const dimensionedScalar &ds)
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
scalar ddT(const scalar p, const scalar T, const scalarField &c) const
static const char nl
Definition: Ostream.H:265
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.
void dcidc(const scalar p, const scalar T, const scalarField &c, scalarField &dcidc) const
Species concentration derivative of the pressure dependent term.
virtual Ostream & write(const token &)=0
Write next token to stream.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576