ROOT logo
// module identifier line...
// Author: Brian Thorsbro, 24/6-2014


#ifndef THEPMCPARSER_H
#define THEPMCPARSER_H

#include <string>
#include <list>
#include <set>
#include "TTree.h"
#include "TClonesArray.h"
#include "TParticle.h"

namespace HepMC {
  class IO_BaseClass;
  class GenVertex;
  class GenEvent;
}

class THepMCParser {

private:

   TTree * fTree;
   void init(HepMC::IO_BaseClass *);

   // 'name(pdgcode,invmass)' if db entry or 'pdgcode()' if no db entry
   static std::string GetParticleName(TParticle *);

   // Fold out the vertex tree in a list, such that daughters are always last
   static void ExploreVertex(HepMC::GenVertex *, std::list<HepMC::GenVertex*> &, std::set<int> &, bool);


public:

   // Header struct
   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";
   struct HeavyIonHeader_t {
      Int_t   Ncoll_hard;                    // Number of hard scatterings
      Int_t   Npart_proj;                    // Number of projectile participants
      Int_t   Npart_targ;                    // Number of target participants
      Int_t   Ncoll;                         // Number of NN (nucleon-nucleon) collisions
      Int_t   spectator_neutrons;            // Number of spectator neutrons
      Int_t   spectator_protons;             // Number of spectator protons
      Int_t   N_Nwounded_collisions;         // Number of N-Nwounded collisions
      Int_t   Nwounded_N_collisions;         // Number of Nwounded-N collisons
      Int_t   Nwounded_Nwounded_collisions;  // Number of Nwounded-Nwounded collisions
      Float_t impact_parameter;              // Impact Parameter(in fm) of collision
      Float_t event_plane_angle;             // Azimuthal angle of event plane
      Float_t eccentricity;                  // eccentricity of participating nucleons in the transverse plane (as in phobos nucl-ex/0510031)
      Float_t sigma_inel_NN;                 // nucleon-nucleon inelastic (including diffractive) cross-section
   };
   static const char * fgPdfHeaderBranchString; // "id1/I,id2,pdf_id1,pdf_id2,x1/D,x2,scalePDF,pdf1,pdf2";
   struct PdfHeader_t {
      Int_t    id1;        // flavour code of first parton
      Int_t    id2;        // flavour code of second parton
      Int_t    pdf_id1;    // LHAPDF set id of first parton
      Int_t    pdf_id2;    // LHAPDF set id of second parton
      Double_t x1;         // fraction of beam momentum carried by first parton ("beam side")
      Double_t x2;         // fraction of beam momentum carried by second parton ("target side")
      Double_t scalePDF;   // Q-scale used in evaluation of PDF's   (in GeV)
      Double_t pdf1;       // PDF (id1, x1, Q) - x*f(x)
      Double_t pdf2;       // PDF (id2, x2, Q) - x*f(x)
   };

   // Default constructor/destructor stuff, don't inherit from this class unless you handle the tree pointer
   inline THepMCParser() : fTree(0) {;} // nullptr in c++11
   inline virtual ~THepMCParser() {;} // should be a memory management of the TTree...

   // The actual useful constructors, either take:
   //  - a file name for a file with HepMC data or
   //  - a HepMC event data structure
   THepMCParser(const char *);
   THepMCParser(HepMC::IO_BaseClass *);

   // Optional validators, set the argument to true for verbose output to STDERR
   // WARNING: including status code 2 may produce invalid flag when it is in fact valid, not recommended to use
   bool IsValidMotherDaughtersConsitency(bool useStdErr = false, bool requireSecondMotherBeforeDaughters = false);
   bool IsValidParticleInvariantMass(bool useStdErr = false, bool includeStatusCode2Particles = false);
   bool IsValidVertexInvariantMass(bool useStdErr = false, bool includeStatusCode2Particles = false);

   // Show the decay chain by point of view of some particle
   static std::string ListReactionChain(TClonesArray *, Int_t);

   // Access the TTree generated or write it to a file
   TTree * GetTTree();
   void WriteTTreeToFile(const char *);

   // ***** The constructors will parse all the events with this function for you *****
   // The work horse of the parser, takes one event and generates the corresponding TClonesArray of TParticles,
   // Units links directly to HepMC library which governs which units are available
   // The default momentum unit is GEV, other option is MEV
   // The default length unit is CM, other option is MM
   // requireSecondMotherBeforeDaughters defaults to false
   //  - Success if length of return string is 0, the provided TClonesArray has been filled
   //  - Failure otherwise and the return string then contains the error message
   //
   // Note:
   // The TClonesArray must be initialized before hand
   // the array will be cleared and filled with TParticles
   // The capacity of the array will be set to the number of particles parsed
   // First mother will be before the daughters, but second mother may be after
   // The two first particles in the array will be the beam particles
   // The status code set on the particle is copied from HepMC, i.e.
   //  - 1 = final particle
   //  - 2 = transitory particle
   //  - 4 = beam particle
   // The function is static to enable customized wrappers around it
   static std::string ParseGenEvent2TCloneArray(HepMC::GenEvent *, TClonesArray *, std::string momUnit = "GEV", std::string lenUnit = "CM", bool requireSecondMotherBeforeDaughters = false);


   // Depending on the implementation of HepMC::IO_BaseClass there may be information on
   // heavy ions or parton distribution functions available. This function will pull them
   // out if available.
   // The caller must supply allocated structures which will then be filled.
   static std::string ParseGenEvent2HeaderStructs(HepMC::GenEvent *, HeavyIonHeader_t &, PdfHeader_t &, bool fillZeroOnMissingHeavyIon = true, bool fillZeroOnMissingPdf = true);


};

#endif
 THepMCParser.h:1
 THepMCParser.h:2
 THepMCParser.h:3
 THepMCParser.h:4
 THepMCParser.h:5
 THepMCParser.h:6
 THepMCParser.h:7
 THepMCParser.h:8
 THepMCParser.h:9
 THepMCParser.h:10
 THepMCParser.h:11
 THepMCParser.h:12
 THepMCParser.h:13
 THepMCParser.h:14
 THepMCParser.h:15
 THepMCParser.h:16
 THepMCParser.h:17
 THepMCParser.h:18
 THepMCParser.h:19
 THepMCParser.h:20
 THepMCParser.h:21
 THepMCParser.h:22
 THepMCParser.h:23
 THepMCParser.h:24
 THepMCParser.h:25
 THepMCParser.h:26
 THepMCParser.h:27
 THepMCParser.h:28
 THepMCParser.h:29
 THepMCParser.h:30
 THepMCParser.h:31
 THepMCParser.h:32
 THepMCParser.h:33
 THepMCParser.h:34
 THepMCParser.h:35
 THepMCParser.h:36
 THepMCParser.h:37
 THepMCParser.h:38
 THepMCParser.h:39
 THepMCParser.h:40
 THepMCParser.h:41
 THepMCParser.h:42
 THepMCParser.h:43
 THepMCParser.h:44
 THepMCParser.h:45
 THepMCParser.h:46
 THepMCParser.h:47
 THepMCParser.h:48
 THepMCParser.h:49
 THepMCParser.h:50
 THepMCParser.h:51
 THepMCParser.h:52
 THepMCParser.h:53
 THepMCParser.h:54
 THepMCParser.h:55
 THepMCParser.h:56
 THepMCParser.h:57
 THepMCParser.h:58
 THepMCParser.h:59
 THepMCParser.h:60
 THepMCParser.h:61
 THepMCParser.h:62
 THepMCParser.h:63
 THepMCParser.h:64
 THepMCParser.h:65
 THepMCParser.h:66
 THepMCParser.h:67
 THepMCParser.h:68
 THepMCParser.h:69
 THepMCParser.h:70
 THepMCParser.h:71
 THepMCParser.h:72
 THepMCParser.h:73
 THepMCParser.h:74
 THepMCParser.h:75
 THepMCParser.h:76
 THepMCParser.h:77
 THepMCParser.h:78
 THepMCParser.h:79
 THepMCParser.h:80
 THepMCParser.h:81
 THepMCParser.h:82
 THepMCParser.h:83
 THepMCParser.h:84
 THepMCParser.h:85
 THepMCParser.h:86
 THepMCParser.h:87
 THepMCParser.h:88
 THepMCParser.h:89
 THepMCParser.h:90
 THepMCParser.h:91
 THepMCParser.h:92
 THepMCParser.h:93
 THepMCParser.h:94
 THepMCParser.h:95
 THepMCParser.h:96
 THepMCParser.h:97
 THepMCParser.h:98
 THepMCParser.h:99
 THepMCParser.h:100
 THepMCParser.h:101
 THepMCParser.h:102
 THepMCParser.h:103
 THepMCParser.h:104
 THepMCParser.h:105
 THepMCParser.h:106
 THepMCParser.h:107
 THepMCParser.h:108
 THepMCParser.h:109
 THepMCParser.h:110
 THepMCParser.h:111
 THepMCParser.h:112
 THepMCParser.h:113
 THepMCParser.h:114
 THepMCParser.h:115
 THepMCParser.h:116
 THepMCParser.h:117
 THepMCParser.h:118
 THepMCParser.h:119
 THepMCParser.h:120
 THepMCParser.h:121
 THepMCParser.h:122