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-2024 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 dimensionSet& dims,
38  const dictionary& dict
39 )
40 :
41  ArrheniusReactionRate(species, dims, dict),
42  phaseName_(dict.lookup("phase")),
43  ob_(ob),
44  tAv_(nullptr)
45 {}
46 
47 
48 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
49 
51 {
53 
54  const phaseModel& phase =
55  ob_.lookupObject<phaseModel>(IOobject::groupName("alpha", phaseName_));
56 
57  tAv_ = phase.diameter().Av();
58 }
59 
60 
62 {
64 
65  tAv_.clear();
66 }
67 
68 
69 inline Foam::scalar Foam::phaseSurfaceArrheniusReactionRate::operator()
70 (
71  const scalar p,
72  const scalar T,
73  const scalarField& c,
74  const label li
75 ) const
76 {
77  return ArrheniusReactionRate::operator()(p, T, c, li)*tAv_()[li];
78 }
79 
80 
82 (
83  const scalar p,
84  const scalar T,
85  const scalarField& c,
86  const label li
87 ) const
88 {
89  return ArrheniusReactionRate::ddT(p, T, c, li)*tAv_()[li];
90 }
91 
92 
94 {
96  writeEntry(os, "phase", phaseName_);
97 }
98 
99 
100 inline Foam::Ostream& Foam::operator<<
101 (
102  Ostream& os,
104 )
105 {
106  arr.write(os);
107  return os;
108 }
109 
110 
111 // ************************************************************************* //
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.
virtual tmp< volScalarField > Av() const =0
Return the surface area per unit volume.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
A wordList with hashed indices for faster lookup by name.
Registry of regIOobjects.
const diameterModel & diameter() const
Return a reference to the diameterModel of the phase.
Definition: phaseModel.C:187
A modified Arrhenius reaction rate given by:
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.
phaseSurfaceArrheniusReactionRate(const speciesTable &species, const objectRegistry &ob, const dimensionSet &dims, const dictionary &dict)
Construct from dictionary.
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