]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Ensuring that proper units are used when parsing and HepMC file (B Thorsbro)
authormfloris <michele.floris@cern.ch>
Wed, 20 Aug 2014 15:29:57 +0000 (17:29 +0200)
committermfloris <michele.floris@cern.ch>
Wed, 20 Aug 2014 15:29:57 +0000 (17:29 +0200)
The parser asks the library to use GEV and CM and a conversion will
happen if needed to ensure the numbers are correct

EVGEN/AliGenReaderHepMC.cxx
TEvtGen/THepMCParser.cxx
TEvtGen/THepMCParser.h

index 7e25aa200e3a50b300fea84d2be29d7221dc18e9..88636bf657618684e31759b8d51f79f9ec889a63 100644 (file)
@@ -57,7 +57,7 @@ Int_t AliGenReaderHepMC::NextEvent()
    if (fGenEvent) delete fGenEvent;
    // Read the next event
    if ((fGenEvent = fEventsHandle->read_next_event())) {
-      THepMCParser::ParseGenEvent2TCloneArray(fGenEvent,fParticleArray,false);
+      THepMCParser::ParseGenEvent2TCloneArray(fGenEvent,fParticleArray,"GEV","CM",false);
       fParticleIterator->Reset();
       THepMCParser::HeavyIonHeader_t heavyIonHeader;
       THepMCParser::PdfHeader_t pdfHeader;
index f6531e1f0a19ac264c2742f0405aa34466787962..4ff2fb7313f26fa5732b185d355f943acf022d4b 100644 (file)
@@ -9,7 +9,6 @@
 #include <time.h>
 
 #include "THepMCParser.h"
-#include "HepMC/IO_GenEvent.h"
 #include "TObject.h"
 #include "TTree.h"
 #include "TClonesArray.h"
@@ -18,6 +17,7 @@
 #include "TDatabasePDG.h"
 #include "HepMC/IO_GenEvent.h"
 
+
 using namespace std;
 
 
@@ -342,11 +342,12 @@ string THepMCParser::ListReactionChain(TClonesArray * particles, Int_t particleI
 }
 
 
-string THepMCParser::ParseGenEvent2TCloneArray(HepMC::GenEvent * genEvent, TClonesArray * array, bool requireSecondMotherBeforeDaughters)
+string THepMCParser::ParseGenEvent2TCloneArray(HepMC::GenEvent * genEvent, TClonesArray * array, string momUnit, string lenUnit, bool requireSecondMotherBeforeDaughters)
 {
    if (requireSecondMotherBeforeDaughters) {
       return "requireSecondMotherBeforeDaughters not implemented yet!";
    }
+   genEvent->use_units(momUnit, lenUnit);
    array->Clear();
    ostringstream errMsgStream;
    map<int,Int_t> partonMemory; // unordered_map in c++11 - but probably not much performance gain from that: log(n) vs log(1) where constant can be high
index 1d4625800a41f6fcb66bf96ecc79a5413cb762d9..2aabfea1916d35c55d1eaacac3b59a4ad2947fa9 100644 (file)
@@ -89,6 +89,9 @@ public:
 
    // ***** 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
@@ -104,7 +107,7 @@ public:
    //  - 2 = transitory particle
    //  - 4 = beam particle
    // The function is static to enable customized wrappers around it
-   static std::string ParseGenEvent2TCloneArray(HepMC::GenEvent *, TClonesArray *, bool requireSecondMotherBeforeDaughters = false);
+   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