surfaceArrheniusReactionRateI.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) 2019-2023 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 speciesTable& species,
31  const objectRegistry& ob,
32  const dictionary& dict
33 )
34 :
35  ArrheniusReactionRate(species, dict),
36  AvName_(dict.lookupBackwardsCompatible({"Av", "a"})),
37  ob_(ob),
38  tAv_(nullptr)
39 {}
40 
41 
42 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
43 
45 {
47 
48  tAv_ = ob_.lookupObject<volScalarField::Internal>(AvName_);
49 }
50 
51 
53 {
55 
56  tAv_ = tmp<volScalarField::Internal>(nullptr);
57 }
58 
59 
60 inline Foam::scalar Foam::surfaceArrheniusReactionRate::operator()
61 (
62  const scalar p,
63  const scalar T,
64  const scalarField& c,
65  const label li
66 ) const
67 {
68  return ArrheniusReactionRate::operator()(p, T, c, li)*tAv_()[li];
69 }
70 
71 
73 (
74  const scalar p,
75  const scalar T,
76  const scalarField& c,
77  const label li
78 ) const
79 {
80  return ArrheniusReactionRate::ddT(p, T, c, li)*tAv_()[li];
81 }
82 
83 
85 {
87  writeEntry(os, "a", AvName_);
88 }
89 
90 
91 inline Foam::Ostream& Foam::operator<<
92 (
93  Ostream& os,
95 )
96 {
97  arr.write(os);
98  return os;
99 }
100 
101 
102 // ************************************************************************* //
Arrhenius reaction rate given by:
scalar operator()(const scalar p, const scalar T, const scalarField &c, const label li) const
Return the rate.
void preEvaluate() const
Pre-evaluation hook.
void postEvaluate() const
Post-evaluation hook.
void write(Ostream &os) const
Write to stream.
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
The derivative of the rate w.r.t. temperature.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
A wordList with hashed indices for faster lookup by name.
Registry of regIOobjects.
A modified Arrhenius reaction rate given by:
void preEvaluate() const
Pre-evaluation hook.
void postEvaluate() const
Post-evaluation hook.
void write(Ostream &os) const
Write to stream.
surfaceArrheniusReactionRate(const speciesTable &species, const objectRegistry &ob, const dictionary &dict)
Construct from dictionary.
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
Evaluate the derivative.
A class for managing temporary objects.
Definition: tmp.H:55
const dimensionedScalar c
Speed of light in a vacuum.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p