-
Notifications
You must be signed in to change notification settings - Fork 496
Expand file tree
/
Copy pathAODToHepMC.h
More file actions
610 lines (599 loc) · 21.9 KB
/
AODToHepMC.h
File metadata and controls
610 lines (599 loc) · 21.9 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
// Copyright 2023-2099 CERN and copyright holders of ALICE O2.
// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
// All rights not expressly granted are reserved.
//
// This software is distributed under the terms of the GNU General Public
// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
//
// In applying this license CERN does not waive the privileges and immunities
// granted to it by virtue of its status as an Intergovernmental Organization
// or submit itself to any jurisdiction.
/// @author Christian Holm Christensen <cholm@nbi.dk>
#ifndef ALICEO2_EVENTGEN_AODTOHEPMC_H_
#define ALICEO2_EVENTGEN_AODTOHEPMC_H_
#include <Framework/AnalysisDataModel.h>
#include <Framework/AnalysisManagers.h>
#include <Framework/Configurable.h>
#include <HepMC3/GenEvent.h>
#include <HepMC3/GenParticle.h>
#include <HepMC3/GenVertex.h>
#include <HepMC3/GenPdfInfo.h>
#include <HepMC3/GenHeavyIon.h>
#include <HepMC3/GenCrossSection.h>
#include <HepMC3/WriterAscii.h>
#include <fstream>
#include <list>
namespace o2
{
namespace eventgen
{
/**
* Convert AOD tables of MC information into a HepMC event structure.
*
* The event structure is kept in memory.
*
* The conventions used here are
*
* - A @e track is a kinematics particle stored in the @c
* o2::aod::McParticles table (aliased as AODToHepMC::Tracks), and
* correspond to a track used during simulation transport. These
* are of type @c o2::aod::McParticle aliased to AODToHepMC::Track.
*
* - A @e generated track is a particle originally made by the event
* generator. These can be recognised by @c
* Track::producedByGenerator()
*
* - A @e particle is a particle stored in the HepMC event structure
* of the class @c HepMC3::GenParticle.
*
* - A @e vertex is where a particle is either produced or disappears.
* A vertex with incoming particles will have a number of outgoing
* particles. Thus a vertex is the @a production vertex for
* out-going particles of that vertex, and the @e end vertex of
* incoming particles of that vertex. Vertexes are of the type @c
* HepMC3::GenVertex.
*
* - The relationship between mother and daughter tracks in
* AODToHepMC::Tracks is encoded by indexes.
*
* - At most two mothers can be stored per track. If a track has
* no mothers, then both indexes are -1. If a track has a
* single mother track, then both indexes are the same.
*
* - A track can have any number of daughters. The first daughter
* track is identified by the first index and the last daughter
* by the second stored index. All tracks indexes with in the
* range from the first to the second index (both inclusive) are
* daughter tracks. If a track has no daughters, then both
* indexes are -1. If a track has one daughter, then both
* indexes are equal. The number of daughters can be obtained
* by
*
* last - first + 1
*
* if both last and first are zero or greater.
*
* - An event header is the information stored in a row of @c
* o2::aod::McCollisions table (aliased to AODToHepMC::Headers).
* Each event has one such header (aliased as AODToHepMC::Header).
*
* The header stores information on the event such as event number,
* weight, primary interaction point (IP), and so on.
*
* In addition, auxiliary information (impact parameter,
* @f$N_{\mathrm{part}}@f$, and so on) may be stored separate tables.
*
* - The table @c o2::aod::HepMCXSections (aliased to
* AODToHepMC::XSections) stores the event cross section and weight.
*
* - The table @c o2::aod::HepMCPdfInfos (aliased to
* AODToHepMC::PdfInfos) stores the parton distribution function
* parameters used in the event.
*
* - The table @c o2::aod::HepMCHeavyIons (aliased to
* AODToHepMC::PdfHeavyIons) stores the heavy-ion auxiliary
* information, such as the event impact parameter @f$b@f$,
* @f$N_{\mathrm{part}}@f$, @f$N_{\mathrm{coll}}@f$, and so on.
*
* - A HepMC event has a simple header which contains the event
* number, number of particles, and number of vertexes in the event.
*
* Other event information is stored in specific structures.
*
* - The event cross-section(s) and weight(s) is stored in a @c
* HepMC3::GenCrossSection (aliased to
* AODToHepMC::CrossSectionPtr) object. The information to fill
* into that is read from AODToHepMC::XSections table
*
* - The event parton distribution function parameters are stored in
* a @c HepMC3::GenPdfInfo (aliased to AODToHepMC::PdfInfoPtr)
* object. The information to fill into that is read from
* AODToHepMC::PdfInfos table.
*
* - The event heavy-ion parameters are stored in
* a @c HepMC3::GenHeavyIon (aliased to AODToHepMC::HeavyIonPtr)
* object. The information to fill into that is read from
* AODToHepMC::HeavyIons table.
*
* The conversion is done as follows:
*
* - For all MC tracks, create the correspond particle
* (AODToHepMC::makeParticleRecursive)
*
* - Check if we have already created the corresponding particle
* (AODToHepMC::getParticle). If so, then we go on to the next
* track.
*
* - If we are asked to only keep generated tracks
* (AODToHepMC::mOnlyGen set by option @c --hepmc-only-generated),
* i.e., tracks that were posted by the event generator, and this
* track is not such a track (AODToHepMC::isIgnored), nor forced
* (argument @c force) then return nothing. Note that mothers of
* particles that we keep are @e always made (@c force is set to
* true).
*
* - If we do not have a particle yet for the track (look-up in
* AODToHepMC::mParticles) and it is not excluded, then create a
* particle that corresponds to the track, and store the mapping
* from track to particle (in AODToHepMC::mParticles).
*
* - For all mother tracks of the current track, create the
* corresponding particle, following this algorithm (recursion,
* AODToHepMC::makeParticleRecursive).
*
* - If a mother particle has an end vertex, set that vertex as the
* production vertex of the current particle.
*
* - If no mother particle has an end vertex, and this particle is
* not a beam particle and it does have mothers, then create a
* vertex (AODToHepMC::makeVertex) corresponding to the track
* production position, and add this particle as an outgoing
* particle of that vertex.
*
* - In this case, if some mother does not have an outgoing vertex,
* add that mother as an incoming particle to the created
* vertex.
*
* - If the particle is a beam particle (status=4) then store this
* particle as a beam particle.
*
* - If not a beam particle, and the particle has no mothers, mark
* this particle as an orphan.
*
* - Once we have made particle for all tracks, we flesh-out all
* particles. For all tracks (AODToHepMC::fleshOutParticle)
*
* - If this track is ignored (AODToHepMC::isIgnored), or we have no
* particle (AODToHepMC::getParticle) corresponding to this track,
* go on to the next track.
*
* - Get the end vertex of the particle. If any. Set the candidate
* end vertex to this vertex, whether it exists or not.
*
* - Then for each daughter track of the current track, check if it is
* consistent with the end vertex of the particle.
*
* - Check if the dauther is ignored (AODToHepMC::isIgnored). If so,
* move on to the next daughter track.
*
* - Get the particle corresponding to the daughter track
* (AODToHepMC::getParticle.). If no particle is found, move on
* to the next daughter track.
*
* - Check that the daughter particle has an end vertex. If it
* doesn't, mark it as a @a head-less particle.
*
* - If the production vertex of the daughter doesn't match the end
* vertex of the current particle, or the current candidate end
* vertex, then issue a warning, and move on to the next daughter.
*
* - Update the candidate end vertex to the daughter end vertex.
*
* - After processing all daughter particles, and if we have no end
* vertex for the current particle, then
*
* - if we found no candiate end vertex, and the particle is either a
* beam (status=4) or decayed (status=2) particle, issue a
* warning.
*
* - if we do have a candidate end vertex, set that as the end vertex
* of the current particle.
*
* - If, after this, the current particle does have an end vertex,
* loop over all daughters previsouly as head-less and set their
* production vertex to this end vertex.
*
* - At this point, we should have the particles in a proper HepMC
* event structure.
*
* - During simulation transport, the interaction point vertex (IP)
* may not be at (0,0,0,0). Since some consumers of the the HepMC
* event structure may expect the IP to be at (0,0,0,0), we can
* recenter (AODMCToHepMC::recenter) all vertexes of the event.
* This is governed by the member AODMCToHepMC::mRecenter set by the
* option @c --hepmc-recenter
*
* - We can then fill in remaing information into the HepMC event header.
*
* - The event number and weight is set from event header
* (AODToHepMC::makeHeader).
*
* - The event cross section(s and weight(s) is set from
* AODToHepMC::CrossSections table.
* (AODToHepMC::makeXSection). If no AODToHepMC::CrossSections row
* is passed, then dummy values are filled in.
*
* - The event parton distribution function parameters are set from
* AODToHepMC::PdfInfos table. (AODToHepMC::makePdfInfo). If no
* AODToHepMC::PdfInfos row is passed, then dummy values are
* filled in.
*
* - The event heavy-ion parameters are set from
* AODToHepMC::HeavyIons table. (AODToHepMC::makeHeavyIon). If no
* AODToHepMC::HeavyIons row is passed, then dummy values are
* filled in.
*
* - Once all the above is done, we have a complete HepMC event
* (AODToHepMC::mEvent). This event structure is kept in memory.
*
* - Optionally (option @c --hepmc-dumb @e filename) we may write
* the events to a file (AODToHepMC::mOutput). The event is still
* kept in memory.
*
* - Clients of this class (e.g., o2::pwgmm::RivetWrapper) may access
* the event structure for further processing.
*
* The utility @c o2-aod-mc-to-hepmc will read in AODs and write out a
* HepMC event file (plain ASCII).
*
*/
struct AODToHepMC {
/**
* Group of configurables which will be added to an program that
* uses this class. Note that it is really the specialisation of
* framework::OptionManager<AODToHepMC> that propagates the options
* to the program.
*/
struct {
/** Option for dumping HepMC event structures to disk. Takes one
* argument - the name of the file to write to. */
std::string dump{""};
/** Option for only storing particles from the event generator.
* Note, if a particle is stored down, then its mothers will also
* be stored. */
bool onlyGen{false};
/** Use HepMC's tree parsing for building event structure */
bool useTree{false};
/** Floating point precision used when writing to disk */
int precision{8};
/** Recenter event at IP=(0,0,0,0). */
bool recenter{false};
} configs;
/**
* @{
* @name The containers we subscribe to
*/
/** Alias of MC collisions table */
using Headers = o2::aod::McCollisions;
/** Alias of MC collisions table row */
using Header = o2::aod::McCollision;
/** Alias MC particles (tracks) table */
using Tracks = o2::aod::McParticles;
/** Alias MC particles (tracks) table row */
using Track = typename Tracks::iterator;
/** Alias auxiliary MC table of cross-sections */
using XSections = o2::aod::HepMCXSections;
/** Alias auxiliary MC table of parton distribution function
* parameters */
using PdfInfos = o2::aod::HepMCPdfInfos;
/** Alias auxiliary MC table of heavy-ion parameters */
using HeavyIons = o2::aod::HepMCHeavyIons;
/** Alias row of auxiliary MC table of cross-sections */
using XSection = typename XSections::iterator;
/** Alias row of auxiliary MC table of parton distribution function
* parameters */
using PdfInfo = typename PdfInfos::iterator;
/** Alias row of auxiliary MC table of heavy-ion parameters */
using HeavyIon = typename HeavyIons::iterator;
/** @} */
/**
* @{
* @name Types from HepMC3
*/
/** Alias HepMC four-vector */
using FourVector = HepMC3::FourVector;
/** Alias (smart-)pointer to HepMC particle */
using ParticlePtr = HepMC3::GenParticlePtr;
/** Alias (smart-)pointer to HepMC vertex */
using VertexPtr = HepMC3::GenVertexPtr;
/** Alias HepMC eventt structure */
using Event = HepMC3::GenEvent;
/** Alias (smart-)pointer to HepMC heavy-ion object */
using HeavyIonPtr = HepMC3::GenHeavyIonPtr;
/** Alias (smart-)pointer to HepMC cross section object */
using CrossSectionPtr = HepMC3::GenCrossSectionPtr;
/** Alias (smart-)pointer to HepMC parton distribution function
* object */
using PdfInfoPtr = HepMC3::GenPdfInfoPtr;
/** Type used to map tracks to HepMC particles */
using ParticleMap = std::unordered_map<long, ParticlePtr>;
/** A container of pointers to particles */
using ParticleVector = std::vector<ParticlePtr>;
/** A container of pointers to particles */
using ParticleList = std::list<ParticlePtr>;
/** A container of pointers to vertexes */
using VertexVector = std::vector<VertexPtr>;
/** A container of pointers to vertexes */
using VertexList = std::list<VertexPtr>;
/** Alias of HepMC writer class */
using WriterAscii = HepMC3::WriterAscii;
/** The of pointer to HepMC writer class */
using WriterAsciiPtr = std::shared_ptr<WriterAscii>;
/** @} */
/**
* @{
* @name HepMC3 objects
*/
/** The result of processing */
Event mEvent;
/** Pointer to cross section-ion information */
CrossSectionPtr mCrossSec = nullptr;
/** Pointer to heavy-ion information */
HeavyIonPtr mIon = nullptr;
/** Pointer to parton distribution function information */
PdfInfoPtr mPdf = nullptr;
/** @} */
/**
* @{
* @name Containers etc.
*/
/** Maps tracks to particles */
ParticleMap mParticles; //! Cache of particles
/** List of vertexes made */
VertexList mVertices; //! Cache of vertices
/** List of beam particles */
ParticleVector mBeams; //! Cache of beam particles
/** Particles without a mother */
ParticleList mOrphans; //! Cache of particles w/o mothers
/** @} */
/**
* @{
* @name Options and such
*/
/** Output writer, if enabled */
WriterAsciiPtr mWriter = nullptr;
/** Current sequential event number */
int mEventNo = 0;
/** The last bunch crossing identifier */
int mLastBC = -1;
/** If true, only store particles from the generator */
bool mOnlyGen = false;
/** If true, use HepMC tree parser */
bool mUseTree = true;
/** Output stream if enabled */
std::ofstream* mOutput = nullptr;
/** Precision used on the output stream */
int mPrecision = 16;
/** If true, recenter IP to (0,0,0,0) */
bool mRecenter = false;
/** @} */
/**
* @{
* @name Interface member functions
*/
/**
* Initialize the converter. Sets internal parameters based on the
* configurables.
*/
virtual void init();
/**
* Call before starting to process an event. This clears the
* current event and internal data structures.
*/
virtual void startEvent();
/**
* Process the collision header and tracks
*
* @param collision Header information
* @param tracks Particle tracks
*/
virtual void process(Header const& collision, Tracks const& tracks);
/**
* Process collision header and HepMC auxiliary information
*
* @param collision Header information
* @param xsections Cross-section table (possible no rows)
* @param pdfs Parton-distribution function table (possible no rows)
* @param heavyions Heavy ion collision table (possible no rows)
*/
virtual void process(Header const& collision,
XSections const& xsections,
PdfInfos const& pdfs,
HeavyIons const& heavyions);
/**
* Call after process an. Thisf finalises the event and optionally
* outputs to dump.
*/
virtual void endEvent();
/**
* End of run - closes output file if enabled. This is called via
* specialisation of o2::framework::OutputManager<AODToHepMC>.
*/
virtual bool postRun()
{
enableDump("");
return true;
}
/** @} */
protected:
/**
* @{
* @name Actual worker member functions
*/
/**
* Generate the final event, including fleshing out the vertexes,
* and so on
*
* @param collision Header information
* @param tracks Particle tracks
*/
virtual void makeEvent(Header const& collision,
Tracks const& tracks);
/**
* Set the various fields in the header of the HepMC3::GenEvent
* object
*
* @param header Header object */
virtual void makeHeader(Header const& header);
/**
* Make cross-section information. If no entry in the table,
* then make dummy information
*/
virtual void makeXSection(XSections const& xsection);
/**
* Make parton-distribition function information. If no entry
* in the table, then make dummy information
*/
virtual void makePdfInfo(PdfInfos const& pdf);
/**
* Make heavy-ion collision information. If no header given,
* then fill in other reasonable values
*/
virtual void makeHeavyIon(HeavyIons const& heavyion,
Header const& header);
/**
* This is supposed to make the beam particles from the information
* available. However, it seems like we really don't have enough
* information, so for now we will do nothing here. Perhaps the
* user will be forced to provide that information - either via the
* analyser configurables or somehow from somewhere else.
*/
virtual void makeBeams(Header const& header, const VertexPtr ip) {}
/**
* Make all particles. We loop through the MC particles, and for
* each of them create that particle and any possible mother
* particles (recursively). This allows us to traverse the data
* rather straight-forwardly.
*
* @param tracks The MC tracks
*/
virtual void makeParticles(Tracks const& tracks);
/**
* Get particle corresponding to track @a no from particle cache
*
* @param ref Track reference
* @return Pointer to HepMC3::GenParticle or null
*/
virtual ParticlePtr getParticle(Track const& ref) const;
/**
* Check if we are ignoring this track
*
* @param track Track to check */
virtual bool isIgnored(Track const& track) const;
/**
* Truely make a particle, and its mother particles if any. We add
* vertexes as needed here too.
*
* Note that this traverses the tree from the bottom up. That is,
* we check if a particle has any mothers, and if it does, we create
* those.
*
* However, this can be a bit problematic, since the Kinematic tree
* (and thus the McParticles table) only allows for at most 2
* mothers. In the case a vertex has 3 or more incoming particles,
* then some of the intermediate particles will be lost.
*
* We remedy that problem by traversing the tree again, but this
* time from the bottom up - that is, we look for daughters of all
* particles, and if they are not registered with the out-going
* vertex, then we reattch the parent to the incoming vertex of the
* daughters. For this to work, we need to map from track index to
* HepMC3::GenParticle.
*
* @param track Current track
* @param tracks MC tracks
* @param force If true, do make the particle even if it would
* otherwise be skipped.
*
* @return Pointer to created particle.
*/
virtual ParticlePtr makeParticleRecursive(Track const& track,
Tracks const& tracks,
bool force = false);
/**
* Generate a HepMC particle from a track. Note that the job here
* is simply to make the object. The more complicated job of adding
* the track to the tree is done above in makeParticleRecursive
*
* @param track MC track
* @param mst Mother status code - updated on return
* @param force Force generation of particle even if onlyGen
*
* @returns Shared pointer to new HepMC particle
*/
virtual ParticlePtr makeParticle(const Track& track,
Int_t& motherStatus,
bool force) const;
/**
* Generate vertex from production vertex of track
*
* @param track MC track
*
* @returns Shared pointer to new HepMC vertex
*/
virtual VertexPtr makeVertex(const Track& track);
/**
* Visit all tracks, but this time look for daughters and
* attach mothers as incoming particles if not already done.
*/
virtual void fleshOut(Tracks const& tracks);
/**
* Flesh out a single particle
*
* @param track The track for which we should flesh out the
* corresponding particle.
*/
virtual void fleshOutParticle(Track const& track, Tracks const& tracks);
/**
* Recenter event to (0,0,0,0). This will use the vertex
* information from the event header to translate all vertexes in
* the event.
*
* @param header Event header
*/
virtual void recenter(Header const& collision);
/** @} */
/**
* Open dump output, or close if an empty string was given.
*
* @param dump
*/
void enableDump(const std::string& dump);
}; /** class Generator **/
} // namespace eventgen
namespace framework
{
struct AODToHepMCPostRun {
static AODToHepMCPostRun& instance()
{
static AODToHepMCPostRun inst{};
return inst;
}
AODToHepMCPostRun(eventgen::AODToHepMC* ptr_ = nullptr)
: ptr{ptr_}
{
}
void endOfStream()
{
if (ptr != nullptr) {
ptr->postRun();
}
}
eventgen::AODToHepMC* ptr = nullptr;
};
} // namespace framework
} // namespace o2
#endif /* ALICEO2_EVENTGEN_GENERATOR_H_ */
// Local Variables:
// mode: C++
// End: