phaseSurfaceArrheniusReactionRateI.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) 2021-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 
27 #include "phaseSystem.H"
28 #include "diameterModel.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
34 (
35  const speciesTable& species,
36  const objectRegistry& ob,
37  const dictionary& dict
38 )
39 :
40  ArrheniusReactionRate(species, dict),
41  phaseName_(dict.lookup("phase")),
42  ob_(ob),
43  tAv_(nullptr)
44 {}
45 
46 
47 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
48 
50 {
52 
53  const phaseModel& phase =
54  ob_.lookupObject<phaseModel>(IOobject::groupName("alpha", phaseName_));
55 
56  tAv_ = phase.dPtr()->Av();
57 }
58 
59 
61 {
63 
64  tAv_.clear();
65 }
66 
67 
68 inline Foam::scalar Foam::phaseSurfaceArrheniusReactionRate::operator()
69 (
70  const scalar p,
71  const scalar T,
72  const scalarField& c,
73  const label li
74 ) const
75 {
76  return ArrheniusReactionRate::operator()(p, T, c, li)*tAv_()[li];
77 }
78 
79 
81 (
82  const scalar p,
83  const scalar T,
84  const scalarField& c,
85  const label li
86 ) const
87 {
88  return ArrheniusReactionRate::ddT(p, T, c, li)*tAv_()[li];
89 }
90 
91 
93 {
95  writeEntry(os, "phase", phaseName_);
96 }
97 
98 
99 inline Foam::Ostream& Foam::operator<<
100 (
101  Ostream& os,
103 )
104 {
105  arr.write(os);
106  return os;
107 }
108 
109 
110 // ************************************************************************* //
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.
static word groupName(Name name, const word &group)
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.
const autoPtr< diameterModel > & dPtr() const
Return const-reference to diameterModel of the phase.
Definition: phaseModel.C:151
A modified Arrhenius reaction rate given by:
phaseSurfaceArrheniusReactionRate(const speciesTable &species, const objectRegistry &ob, const dictionary &dict)
Construct from dictionary.
void write(Ostream &os) const
Write to stream.
scalar ddT(const scalar p, const scalar T, const scalarField &c, const label li) const
Evaluate the derivative.
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