reactionsSensitivityAnalysis.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2016-2020 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29#include "dictionary.H"
30
31// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32
33template<class chemistryType>
36{
37 if (writeToFile() && !prodFilePtr_)
38 {
39 prodFilePtr_ = createFile("production");
40 writeHeader(prodFilePtr_(), "production");
41 writeFileHeader(prodFilePtr_());
42
43 consFilePtr_ = createFile("consumption");
44 writeHeader(consFilePtr_(), "consumption");
45 writeFileHeader(consFilePtr_());
46
47 prodIntFilePtr_ = createFile("productionInt");
48 writeHeader(prodIntFilePtr_(), "productionInt");
49 writeFileHeader(prodIntFilePtr_());
50
51 consIntFilePtr_ = createFile("consumptionInt");
52 writeHeader(consIntFilePtr_(), "consumptionInt");
53 writeFileHeader(consIntFilePtr_());
54 }
55}
56
57
58template<class chemistryType>
61(
62 OFstream& os
63)
64{
65 writeCommented(os, "Reaction");
66
67 forAll(speciesNames_, k)
68 {
69 os << tab << speciesNames_[k] << tab;
70 }
71
72 os << nl << endl;
73}
74
75
76template<class chemistryType>
79(
80 const basicChemistryModel& basicChemistry
81)
82{
83 auto RRt = tmp<DimensionedField<scalar, volMesh>>::New
84 (
85 IOobject
86 (
87 "RR",
88 time_.timeName(),
89 mesh_,
90 IOobject::NO_READ,
91 IOobject::NO_WRITE
92 ),
93 mesh_,
94 dimensionedScalar(dimMass/dimVolume/dimTime, Zero)
95 );
96 auto& RR = RRt.ref();
97
98 scalar dt = time_.deltaT().value();
99
100 endTime_ += dt;
101
102 forAll(production_, speciei)
103 {
104 forAll(production_[speciei], reactioni)
105 {
106 RR = basicChemistry.calculateRR(reactioni, speciei);
107
108 if (RR[0] > 0.0)
109 {
110 production_[speciei][reactioni] = RR[0];
111 productionInt_[speciei][reactioni] += dt*RR[0];
112 }
113 else if (RR[0] < 0.0)
114 {
115 consumption_[speciei][reactioni] = RR[0];
116 consumptionInt_[speciei][reactioni] += dt*RR[0];
117 }
118 else
119 {
120 production_[speciei][reactioni] = 0.0;
121 consumption_[speciei][reactioni] = 0.0;
122 }
123 }
124 }
125}
126
127
128template<class chemistryType>
131{
132
133 consFilePtr_() << "time : " << mesh_.time().value() << tab << nl;
134 consFilePtr_() << "delta T : "<< mesh_.time().deltaT().value() << nl << nl;
135 prodFilePtr_() << "time : " << mesh_.time().value() << tab << nl;
136 prodFilePtr_() << "delta T : "<< mesh_.time().deltaT().value() << nl << nl;
137
138 consIntFilePtr_() << "start time : " << startTime_ << tab
139 << "end time :" << endTime_ << nl;
140
141 prodIntFilePtr_() << "start time : " << startTime_ << tab
142 << "end time :" << endTime_ << nl;
143
144 for (label reactioni = 0; reactioni < nReactions_; ++reactioni)
145 {
146 consFilePtr_() << reactioni << tab;
147 consIntFilePtr_() << reactioni << tab;
148 prodFilePtr_() << reactioni << tab;
149 prodIntFilePtr_() << reactioni << tab;
150
151 forAll(speciesNames_, i)
152 {
153 prodFilePtr_() << production_[i][reactioni] << tab;
154 consFilePtr_() << consumption_[i][reactioni] << tab;
155 prodIntFilePtr_() << productionInt_[i][reactioni] << tab;
156 consIntFilePtr_() << consumptionInt_[i][reactioni] << tab;
157 consumptionInt_[i][reactioni] = 0.0;
158 productionInt_[i][reactioni] = 0.0;
159 }
160 consFilePtr_() << nl;
161 consIntFilePtr_() << nl;
162 prodFilePtr_() << nl;
163 prodIntFilePtr_() << nl;
164 }
165 consFilePtr_() << nl << nl;
166 consIntFilePtr_() << nl << nl;
167 prodFilePtr_() << nl << nl;
168 prodIntFilePtr_() << nl << nl;
169}
170
171
172// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
173
174template<class chemistryType>
177(
178 const word& name,
179 const Time& runTime,
180 const dictionary& dict
181)
182:
184 writeFile(mesh_, name),
185 nReactions_(0),
186 startTime_(0),
187 endTime_(0),
188 production_(0),
189 consumption_(0),
190 productionInt_(0),
191 consumptionInt_(0),
192 speciesNames_(),
193 prodFilePtr_(),
194 consFilePtr_(),
195 prodIntFilePtr_(),
196 consIntFilePtr_()
197{
198 read(dict);
199
200 if (mesh_.nCells() != 1)
201 {
203 << "Function object only applicable to single cell cases"
204 << abort(FatalError);
205 }
206
207 if (foundObject<basicChemistryModel>("chemistryProperties"))
208 {
209 const chemistryType& chemistry = refCast<const chemistryType>
210 (
211 lookupObject<basicChemistryModel>("chemistryProperties")
212 );
213
214 speciesNames_.setSize
215 (
216 chemistry.thermo().composition().species().size()
217 );
218
219 forAll(speciesNames_, i)
220 {
221 speciesNames_[i] = chemistry.thermo().composition().species()[i];
222 }
223
224 nReactions_ = chemistry.nReaction();
225
226 if (production_.size() == 0)
227 {
228 production_.setSize(speciesNames_.size());
229 consumption_.setSize(production_.size());
230 productionInt_.setSize(production_.size());
231 consumptionInt_.setSize(production_.size());
232
233 forAll(production_, i)
234 {
235 production_[i].setSize(nReactions_, 0.0);
236 consumption_[i].setSize(nReactions_, 0.0);
237 productionInt_[i].setSize(nReactions_, 0.0);
238 consumptionInt_[i].setSize(nReactions_, 0.0);
239 }
240 }
241 }
242 else
243 {
245 << " No chemistry model found. "
246 << " Objects available are : " << mesh_.names()
247 << exit(FatalError);
248 }
249}
250
251
252// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
253
254template<class chemistryType>
256(
257 const dictionary& dict
258)
259{
262
263 return true;
264}
265
266
267template<class chemistryType>
269execute()
270{
271 createFileNames();
272
274 lookupObject<basicChemistryModel>("chemistryProperties");
275 calculateSpeciesRR(chemistry);
276
277 return true;
278}
279
280
281template<class chemistryType>
283write()
284{
285 if (Pstream::master())
286 {
287 writeSpeciesRR();
288
289 startTime_ = endTime_;
290 }
291
292 return true;
293}
294
295
296// ************************************************************************* //
label k
void setSize(const label n)
Alias for resize()
Definition: List.H:218
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Base class for chemistry models.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const Type & value() const
Return const reference to value.
Watches for presence of the named trigger file in the case directory and signals a simulation stop (o...
Definition: abort.H:128
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
Computes indicators for reaction rates of creation or destruction of species in each reaction.
virtual bool write()
Calculate the reactionsSensitivityAnalysis and write.
virtual bool read(const dictionary &)
Read the reactionsSensitivityAnalysis data.
const Time & time() const
Return time database.
Base class for writing single files from the function objects.
Definition: writeFile.H:120
wordList names() const
The unsorted names of all objects.
label nCells() const noexcept
Number of mesh cells.
splitCell * master() const
Definition: splitCell.H:113
A class for handling words, derived from Foam::string.
Definition: word.H:68
BasicChemistryModel< psiReactionThermo > & chemistry
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const scalar RR
Universal gas constant: default in [J/(kmol K)].
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
constexpr char tab
The tab '\t' character(0x09)
Definition: Ostream.H:52
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333