]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Added AliGenHepMCEventHeader class and put it in STEERBase library. Updated AliGenRea...
authorbthorsbr <brian.peter.thorsbro@cern.ch>
Mon, 14 Jul 2014 15:26:57 +0000 (17:26 +0200)
committerhristov <Peter.Hristov@cern.ch>
Fri, 8 Aug 2014 14:37:53 +0000 (16:37 +0200)
EVGEN/AliGenReaderHepMC.cxx
STEER/CMakelibSTEERBase.pkg
STEER/STEERBase/AliGenHepMCEventHeader.cxx [new file with mode: 0644]
STEER/STEERBase/AliGenHepMCEventHeader.h [new file with mode: 0644]
STEER/STEERBaseLinkDef.h

index 0bd03677a180bde8cf03e8ba279f49eea7a05832..628808f5fa817afb13cbe6c83240d6f1f59659d8 100644 (file)
@@ -5,8 +5,7 @@
 #include "AliGenReaderHepMC.h"
 #include "AliRun.h"
 #include "AliStack.h"
-
-
+#include "AliGenHepMCEventHeader.h"
 
 
 ClassImp(AliGenReaderHepMC)
@@ -33,9 +32,6 @@ void AliGenReaderHepMC::Init()
    fEventsHandle = new HepMC::IO_GenEvent(fFileName, std::ios::in);
    fParticleArray = new TClonesArray("TParticle");
    fParticleIterator = new TIter(fParticleArray);
-
-
-
 }
 
 Int_t AliGenReaderHepMC::NextEvent()
@@ -44,10 +40,33 @@ Int_t AliGenReaderHepMC::NextEvent()
    if ((fGenEvent = fEventsHandle->read_next_event())) {
       THepMCParser::ParseGenEvent2TCloneArray(fGenEvent,fParticleArray,false);
       fParticleIterator->Reset();
-
-
-      // implement header... somewhere
-
+      THepMCParser::HeavyIonHeader_t heavyIonHeader;
+      THepMCParser::PdfHeader_t pdfHeader;
+      THepMCParser::ParseGenEvent2HeaderStructs(fGenEvent,heavyIonHeader,pdfHeader,true,true);
+      fGenEventHeader = new AliGenHepMCEventHeader(
+            heavyIonHeader.Ncoll_hard,
+            heavyIonHeader.Npart_proj,
+            heavyIonHeader.Npart_targ,
+            heavyIonHeader.Ncoll,
+            heavyIonHeader.spectator_neutrons,
+            heavyIonHeader.spectator_protons,
+            heavyIonHeader.N_Nwounded_collisions,
+            heavyIonHeader.Nwounded_N_collisions,
+            heavyIonHeader.Nwounded_Nwounded_collisions,
+            heavyIonHeader.impact_parameter,
+            heavyIonHeader.event_plane_angle,
+            heavyIonHeader.eccentricity,
+            heavyIonHeader.sigma_inel_NN,
+            pdfHeader.id1,
+            pdfHeader.id2,
+            pdfHeader.pdf_id1,
+            pdfHeader.pdf_id2,
+            pdfHeader.x1,
+            pdfHeader.x2,
+            pdfHeader.scalePDF,
+            pdfHeader.pdf1,
+            pdfHeader.pdf2
+      );
       return fGenEvent->particles_size();
    }
    return 0;
index 384911c936c49430c1a853d74924592859eae25b..3731a820ab071a978c1061ee72c5b6a5407a89e3 100644 (file)
@@ -97,6 +97,7 @@ set ( SRCS
     STEERBase/AliVTOFMatch.cxx
     STEERBase/AliVTOFcluster.cxx
     STEERBase/AliVMultiplicity.cxx
+    STEERBase/AliGenHepMCEventHeader.cxx
   )
 
 string(REPLACE ".cxx" ".h" HDRS  "${SRCS}")
diff --git a/STEER/STEERBase/AliGenHepMCEventHeader.cxx b/STEER/STEERBase/AliGenHepMCEventHeader.cxx
new file mode 100644 (file)
index 0000000..7dd7c94
--- /dev/null
@@ -0,0 +1,125 @@
+/**************************************************************************
+ * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Author: The ALICE Off-line Project.                                    *
+ * Contributors are mentioned in the code where appropriate.              *
+ *                                                                        *
+ * Permission to use, copy, modify and distribute this software and its   *
+ * documentation strictly for non-commercial purposes is hereby granted   *
+ * without fee, provided that the above copyright notice appears in all   *
+ * copies and that both the copyright notice and this permission notice   *
+ * appear in the supporting documentation. The authors make no claims     *
+ * about the suitability of this software for any purpose. It is          *
+ * provided "as is" without express or implied warranty.                  *
+ **************************************************************************/
+
+/* $Id$ */
+
+#include "AliGenHepMCEventHeader.h"
+ClassImp(AliGenHepMCEventHeader)
+
+
+AliGenHepMCEventHeader::AliGenHepMCEventHeader():
+   fNcoll_hard(0),
+   fNpart_proj(0),
+   fNpart_targ(0),
+   fNcoll(0),
+   fspectator_neutrons(0),
+   fspectator_protons(0),
+   fN_Nwounded_collisions(0),
+   fNwounded_N_collisions(0),
+   fNwounded_Nwounded_collisions(0),
+   fimpact_parameter(0.0),
+   fevent_plane_angle(0.0),
+   feccentricity(0.0),
+   fsigma_inel_NN(0.0),
+   fid1(0),
+   fid2(0),
+   fpdf_id1(0),
+   fpdf_id2(0),
+   fx1(0.0),
+   fx2(0.0),
+   fscalePDF(0.0),
+   fpdf1(0.0),
+   fpdf2(0.0)
+{
+   // Default Constructor
+}
+
+AliGenHepMCEventHeader::AliGenHepMCEventHeader(const char* name):
+   AliGenEventHeader(name),
+   fNcoll_hard(0),
+   fNpart_proj(0),
+   fNpart_targ(0),
+   fNcoll(0),
+   fspectator_neutrons(0),
+   fspectator_protons(0),
+   fN_Nwounded_collisions(0),
+   fNwounded_N_collisions(0),
+   fNwounded_Nwounded_collisions(0),
+   fimpact_parameter(0.0),
+   fevent_plane_angle(0.0),
+   feccentricity(0.0),
+   fsigma_inel_NN(0.0),
+   fid1(0),
+   fid2(0),
+   fpdf_id1(0),
+   fpdf_id2(0),
+   fx1(0.0),
+   fx2(0.0),
+   fscalePDF(0.0),
+   fpdf1(0.0),
+   fpdf2(0.0)
+{
+   // Constructor
+}
+
+AliGenHepMCEventHeader::AliGenHepMCEventHeader(
+      Int_t    var_Ncoll_hard,                    // Number of hard scatterings
+      Int_t    var_Npart_proj,                    // Number of projectile participants
+      Int_t    var_Npart_targ,                    // Number of target participants
+      Int_t    var_Ncoll,                         // Number of NN (nucleon-nucleon) collisions
+      Int_t    var_spectator_neutrons,            // Number of spectator neutrons
+      Int_t    var_spectator_protons,             // Number of spectator protons
+      Int_t    var_N_Nwounded_collisions,         // Number of N-Nwounded collisions
+      Int_t    var_Nwounded_N_collisions,         // Number of Nwounded-N collisons
+      Int_t    var_Nwounded_Nwounded_collisions,  // Number of Nwounded-Nwounded collisions
+      Float_t  var_impact_parameter,              // Impact Parameter(in fm) of collision
+      Float_t  var_event_plane_angle,             // Azimuthal angle of event plane
+      Float_t  var_eccentricity,                  // eccentricity of participating nucleons in the transverse plane (as in phobos nucl-ex/0510031)
+      Float_t  var_sigma_inel_NN,                 // nucleon-nucleon inelastic (including diffractive) cross-section
+      Int_t    var_id1,        // flavour code of first parton
+      Int_t    var_id2,        // flavour code of second parton
+      Int_t    var_pdf_id1,    // LHAPDF set id of first parton
+      Int_t    var_pdf_id2,    // LHAPDF set id of second parton
+      Double_t var_x1,         // fraction of beam momentum carried by first parton ("beam side")
+      Double_t var_x2,         // fraction of beam momentum carried by second parton ("target side")
+      Double_t var_scalePDF,   // Q-scale used in evaluation of PDF's   (in GeV)
+      Double_t var_pdf1,       // PDF (id1, x1, Q) - x*f(x)
+      Double_t var_pdf2        // PDF (id2, x2, Q) - x*f(x)
+):
+   fNcoll_hard(var_Ncoll_hard),
+   fNpart_proj(var_Npart_proj),
+   fNpart_targ(var_Npart_targ),
+   fNcoll(var_Ncoll),
+   fspectator_neutrons(var_spectator_neutrons),
+   fspectator_protons(var_spectator_protons),
+   fN_Nwounded_collisions(var_N_Nwounded_collisions),
+   fNwounded_N_collisions(var_Nwounded_N_collisions),
+   fNwounded_Nwounded_collisions(var_Nwounded_Nwounded_collisions),
+   fimpact_parameter(var_impact_parameter),
+   fevent_plane_angle(var_event_plane_angle),
+   feccentricity(var_eccentricity),
+   fsigma_inel_NN(var_sigma_inel_NN),
+   fid1(var_id1),
+   fid2(var_id2),
+   fpdf_id1(var_pdf_id1),
+   fpdf_id2(var_pdf_id2),
+   fx1(var_x1),
+   fx2(var_x2),
+   fscalePDF(var_scalePDF),
+   fpdf1(var_pdf1),
+   fpdf2(var_pdf2)
+{
+   // The Constructor
+}
diff --git a/STEER/STEERBase/AliGenHepMCEventHeader.h b/STEER/STEERBase/AliGenHepMCEventHeader.h
new file mode 100644 (file)
index 0000000..6ccdf8d
--- /dev/null
@@ -0,0 +1,99 @@
+#ifndef ALIGENHEPMCEVENTHEADER_H
+#define ALIGENHEPMCEVENTHEADER_H
+
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+/* $Id$ */
+
+#include "AliGenEventHeader.h"
+
+
+class AliGenHepMCEventHeader : public AliGenEventHeader
+{
+public:
+   AliGenHepMCEventHeader();
+   AliGenHepMCEventHeader(const char* name);
+   AliGenHepMCEventHeader(
+         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
+         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)
+   );
+   virtual ~AliGenHepMCEventHeader() {}
+
+
+   Int_t    Ncoll_hard() const {return fNcoll_hard;} // Number of hard scatterings
+   Int_t    Npart_proj() const {return fNpart_proj;} // Number of projectile participants
+   Int_t    Npart_targ() const {return fNpart_targ;} // Number of target participants
+   Int_t    Ncoll() const {return fNcoll;} // Number of NN (nucleon-nucleon) collisions
+   Int_t    spectator_neutrons() const {return fspectator_neutrons;} // Number of spectator neutrons
+   Int_t    spectator_protons() const {return fspectator_protons;} // Number of spectator protons
+   Int_t    N_Nwounded_collisions() const {return fN_Nwounded_collisions;} // Number of N-Nwounded collisions
+   Int_t    Nwounded_N_collisions() const {return fNwounded_N_collisions;} // Number of Nwounded-N collisons
+   Int_t    Nwounded_Nwounded_collisions() const {return fNwounded_Nwounded_collisions;} // Number of Nwounded-Nwounded collisions
+   Float_t  impact_parameter() const {return fimpact_parameter;} // Impact Parameter(in fm) of collision
+   Float_t  event_plane_angle() const {return fevent_plane_angle;} // Azimuthal angle of event plane
+   Float_t  eccentricity() const {return feccentricity;} // eccentricity of participating nucleons in the transverse plane (as in phobos nucl-ex/0510031)
+   Float_t  sigma_inel_NN() const {return fsigma_inel_NN;} // nucleon-nucleon inelastic (including diffractive) cross-section
+
+   Int_t    id1() const {return fid1;} // flavour code of first parton
+   Int_t    id2() const {return fid2;} // flavour code of second parton
+   Int_t    pdf_id1() const {return fpdf_id1;} // LHAPDF set id of first parton
+   Int_t    pdf_id2() const {return fpdf_id2;} // LHAPDF set id of second parton
+   Double_t x1() const {return fx1;} // fraction of beam momentum carried by first parton ("beam side")
+   Double_t x2() const {return fx2;} // fraction of beam momentum carried by second parton ("target side")
+   Double_t scalePDF() const {return fscalePDF;} // Q-scale used in evaluation of PDF's   (in GeV)
+   Double_t pdf1() const {return fpdf1;} // PDF (id1, x1, Q) - x*f(x)
+   Double_t pdf2() const {return fpdf2;} // PDF (id2, x2, Q) - x*f(x)
+
+protected:
+
+   Int_t    fNcoll_hard;                    // Number of hard scatterings
+   Int_t    fNpart_proj;                    // Number of projectile participants
+   Int_t    fNpart_targ;                    // Number of target participants
+   Int_t    fNcoll;                         // Number of NN (nucleon-nucleon) collisions
+   Int_t    fspectator_neutrons;            // Number of spectator neutrons
+   Int_t    fspectator_protons;             // Number of spectator protons
+   Int_t    fN_Nwounded_collisions;         // Number of N-Nwounded collisions
+   Int_t    fNwounded_N_collisions;         // Number of Nwounded-N collisons
+   Int_t    fNwounded_Nwounded_collisions;  // Number of Nwounded-Nwounded collisions
+   Float_t  fimpact_parameter;              // Impact Parameter(in fm) of collision
+   Float_t  fevent_plane_angle;             // Azimuthal angle of event plane
+   Float_t  feccentricity;                  // eccentricity of participating nucleons in the transverse plane (as in phobos nucl-ex/0510031)
+   Float_t  fsigma_inel_NN;                 // nucleon-nucleon inelastic (including diffractive) cross-section
+
+   Int_t    fid1;        // flavour code of first parton
+   Int_t    fid2;        // flavour code of second parton
+   Int_t    fpdf_id1;    // LHAPDF set id of first parton
+   Int_t    fpdf_id2;    // LHAPDF set id of second parton
+   Double_t fx1;         // fraction of beam momentum carried by first parton ("beam side")
+   Double_t fx2;         // fraction of beam momentum carried by second parton ("target side")
+   Double_t fscalePDF;   // Q-scale used in evaluation of PDF's   (in GeV)
+   Double_t fpdf1;       // PDF (id1, x1, Q) - x*f(x)
+   Double_t fpdf2;       // PDF (id2, x2, Q) - x*f(x)
+
+   ClassDef(AliGenHepMCEventHeader, 2)  // Event header for HepMC event
+};
+
+
+
+#endif
index 942ea0b7ffc73072a765a3ccbb2a7ff2658e2e3f..34c917e56667d9f7f3a02d1bd39f927627476cc4 100644 (file)
 #pragma link C++ class  AliVTOFMatch+;
 #pragma link C++ class AliVTOFcluster+;
 #pragma link C++ class AliVMultiplicity+;
+#pragma link C++ class AliGenHepMCEventHeader+;
 
 #endif