]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TEvtGen/THepMCParser.h
addition to previous fix: remove unnecessary if in deleting the pointer
[u/mrichter/AliRoot.git] / TEvtGen / THepMCParser.h
1 // module identifier line...
2 // Author: Brian Thorsbro, 24/6-2014
3
4
5 #ifndef THEPMCPARSER_H
6 #define THEPMCPARSER_H
7
8 #include <string>
9 #include <list>
10 #include <set>
11 #include "TTree.h"
12 #include "TClonesArray.h"
13 #include "TParticle.h"
14
15 namespace HepMC {
16   class IO_BaseClass;
17   class GenVertex;
18   class GenEvent;
19 }
20
21 class THepMCParser {
22
23 private:
24
25    TTree * fTree;
26    void init(HepMC::IO_BaseClass *);
27
28    // 'name(pdgcode,invmass)' if db entry or 'pdgcode()' if no db entry
29    static std::string GetParticleName(TParticle *);
30
31    // Fold out the vertex tree in a list, such that daughters are always last
32    static void ExploreVertex(HepMC::GenVertex *, std::list<HepMC::GenVertex*> &, std::set<int> &, bool);
33
34
35 public:
36
37    // Header struct
38    static const char * fgHeavyIonHeaderBranchString; // "Ncoll_hard/I,Npart_proj,Npart_targ,Ncoll,spectator_neutrons,spectator_protons,N_Nwounded_collisions,Nwounded_N_collisions,Nwounded_Nwounded_collisions,impact_parameter/F,event_plane_angle,eccentricity,sigma_inel_NN";
39    struct HeavyIonHeader_t {
40       Int_t   Ncoll_hard;                    // Number of hard scatterings
41       Int_t   Npart_proj;                    // Number of projectile participants
42       Int_t   Npart_targ;                    // Number of target participants
43       Int_t   Ncoll;                         // Number of NN (nucleon-nucleon) collisions
44       Int_t   spectator_neutrons;            // Number of spectator neutrons
45       Int_t   spectator_protons;             // Number of spectator protons
46       Int_t   N_Nwounded_collisions;         // Number of N-Nwounded collisions
47       Int_t   Nwounded_N_collisions;         // Number of Nwounded-N collisons
48       Int_t   Nwounded_Nwounded_collisions;  // Number of Nwounded-Nwounded collisions
49       Float_t impact_parameter;              // Impact Parameter(in fm) of collision
50       Float_t event_plane_angle;             // Azimuthal angle of event plane
51       Float_t eccentricity;                  // eccentricity of participating nucleons in the transverse plane (as in phobos nucl-ex/0510031)
52       Float_t sigma_inel_NN;                 // nucleon-nucleon inelastic (including diffractive) cross-section
53    };
54    static const char * fgPdfHeaderBranchString; // "id1/I,id2,pdf_id1,pdf_id2,x1/D,x2,scalePDF,pdf1,pdf2";
55    struct PdfHeader_t {
56       Int_t    id1;        // flavour code of first parton
57       Int_t    id2;        // flavour code of second parton
58       Int_t    pdf_id1;    // LHAPDF set id of first parton
59       Int_t    pdf_id2;    // LHAPDF set id of second parton
60       Double_t x1;         // fraction of beam momentum carried by first parton ("beam side")
61       Double_t x2;         // fraction of beam momentum carried by second parton ("target side")
62       Double_t scalePDF;   // Q-scale used in evaluation of PDF's   (in GeV)
63       Double_t pdf1;       // PDF (id1, x1, Q) - x*f(x)
64       Double_t pdf2;       // PDF (id2, x2, Q) - x*f(x)
65    };
66
67    // Default constructor/destructor stuff, don't inherit from this class unless you handle the tree pointer
68    inline THepMCParser() : fTree(0) {;} // nullptr in c++11
69    inline virtual ~THepMCParser() {;} // should be a memory management of the TTree...
70
71    // The actual useful constructors, either take:
72    //  - a file name for a file with HepMC data or
73    //  - a HepMC event data structure
74    THepMCParser(const char *);
75    THepMCParser(HepMC::IO_BaseClass *);
76
77    // Optional validators, set the argument to true for verbose output to STDERR
78    // WARNING: including status code 2 may produce invalid flag when it is in fact valid, not recommended to use
79    bool IsValidMotherDaughtersConsitency(bool useStdErr = false, bool requireSecondMotherBeforeDaughters = false);
80    bool IsValidParticleInvariantMass(bool useStdErr = false, bool includeStatusCode2Particles = false);
81    bool IsValidVertexInvariantMass(bool useStdErr = false, bool includeStatusCode2Particles = false);
82
83    // Show the decay chain by point of view of some particle
84    static std::string ListReactionChain(TClonesArray *, Int_t);
85
86    // Access the TTree generated or write it to a file
87    TTree * GetTTree();
88    void WriteTTreeToFile(const char *);
89
90    // ***** The constructors will parse all the events with this function for you *****
91    // The work horse of the parser, takes one event and generates the corresponding TClonesArray of TParticles,
92    // requireSecondMotherBeforeDaughters defaults to false
93    //  - Success if length of return string is 0, the provided TClonesArray has been filled
94    //  - Failure otherwise and the return string then contains the error message
95    //
96    // Note:
97    // The TClonesArray must be initialized before hand
98    // the array will be cleared and filled with TParticles
99    // The capacity of the array will be set to the number of particles parsed
100    // First mother will be before the daughters, but second mother may be after
101    // The two first particles in the array will be the beam particles
102    // The status code set on the particle is copied from HepMC, i.e.
103    //  - 1 = final particle
104    //  - 2 = transitory particle
105    //  - 4 = beam particle
106    // The function is static to enable customized wrappers around it
107    static std::string ParseGenEvent2TCloneArray(HepMC::GenEvent *, TClonesArray *, bool requireSecondMotherBeforeDaughters = false);
108
109
110    // Depending on the implementation of HepMC::IO_BaseClass there may be information on
111    // heavy ions or parton distribution functions available. This function will pull them
112    // out if available.
113    // The caller must supply allocated structures which will then be filled.
114    static std::string ParseGenEvent2HeaderStructs(HepMC::GenEvent *, HeavyIonHeader_t &, PdfHeader_t &, bool fillZeroOnMissingHeavyIon = true, bool fillZeroOnMissingPdf = true);
115
116
117 };
118
119 #endif