]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
AliMCEvent: Access to Transport MC Truth during analysis
authormorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 1 Aug 2007 07:22:54 +0000 (07:22 +0000)
committermorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 1 Aug 2007 07:22:54 +0000 (07:22 +0000)
STEER/AliMCEvent.cxx [new file with mode: 0644]
STEER/AliMCEvent.h [new file with mode: 0644]
STEER/STEERLinkDef.h
STEER/libSTEER.pkg

diff --git a/STEER/AliMCEvent.cxx b/STEER/AliMCEvent.cxx
new file mode 100644 (file)
index 0000000..167d631
--- /dev/null
@@ -0,0 +1,214 @@
+/**************************************************************************
+ * 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$ */
+//---------------------------------------------------------------------------------
+//                          Class AliMCEvent
+// This class gives access to MC truth during the analysis.
+// Monte Carlo truth is containe in the kinematics tree (produced particles) and 
+// the tree of reference hits.
+//      
+// Origin: Andreas Morsch, CERN, andreas.morsch@cern.ch 
+//---------------------------------------------------------------------------------
+
+
+
+#include "AliMCEvent.h"
+#include "AliTrackReference.h"
+#include "AliHeader.h"
+#include "AliStack.h"
+
+#include <TTree.h>
+#include <TFile.h>
+#include <TParticle.h>
+#include <TClonesArray.h>
+#include <TDirectoryFile.h>
+#include <TArrow.h>
+#include <TMarker.h>
+#include <TH2F.h>
+
+
+ClassImp(AliMCEvent)
+
+AliMCEvent::AliMCEvent() :
+    AliVirtualEventHandler(),
+    fFileE(0),
+    fFileK(0),
+    fFileTR(0),
+    fTreeE(0),
+    fTreeK(0),
+    fTreeTR(0),
+    fStack(0),
+    fHeader(0),
+    fTrackReferences(0),
+    fNEvent(-1),
+    fEvent(-1),
+    fNprimaries(-1),
+    fNparticles(-1)
+{
+    // Default constructor
+}
+
+AliMCEvent::AliMCEvent(const char* name, const char* title) :
+    AliVirtualEventHandler(name, title),
+    fFileE(0),
+    fFileK(0),
+    fFileTR(0),
+    fTreeE(0),
+    fTreeK(0),
+    fTreeTR(0),
+    fStack(0),
+    fHeader(new AliHeader()),
+    fTrackReferences(new TClonesArray("AliTrackReference", 200)),
+    fNEvent(-1),
+    fEvent(-1),
+    fNprimaries(-1),
+    fNparticles(-1)
+{
+    // Constructor
+}
+AliMCEvent::~AliMCEvent()
+{ 
+    // Destructor
+    delete fFileE;
+    delete fFileK;
+    delete fFileTR;
+
+    delete fTreeE;
+    delete fTreeK;
+    delete fTreeTR;
+    
+}
+
+Bool_t AliMCEvent::InitIO(Option_t* /*opt*/) 
+{ 
+    // Initialize input
+    fFileE = new TFile("galice.root");
+    fFileE->GetObject("TE", fTreeE);
+    fTreeE->SetBranchAddress("Header", &fHeader);
+    fNEvent = fTreeE->GetEntries();
+    //
+    // Tree K
+    fFileK = new TFile("Kinematics.root");
+    //
+    // Tree TR
+    fFileTR = new TFile("TrackRefs.root");
+    //
+    // Reset the event number
+    fEvent = -1;
+    printf("Number of events in this directory %5d \n", fNEvent);
+    return kTRUE;
+    
+}
+
+Bool_t AliMCEvent::BeginEvent()
+{ 
+    // Read the next event
+    fEvent++;
+    char folder[20];
+    sprintf(folder, "Event%d", fEvent);
+    // TreeE
+    fTreeE->GetEntry(fEvent);
+    fStack = fHeader->Stack();
+    // Tree K
+    TDirectoryFile* dirK  = 0;
+    fFileK->GetObject(folder, dirK);
+    dirK->GetObject("TreeK", fTreeK);
+    fStack->ConnectTree(fTreeK);
+    fStack->GetEvent();
+    //Tree TR 
+    TDirectoryFile* dirTR = 0;
+    fFileTR->GetObject(folder, dirTR);
+    dirTR->GetObject("TreeTR", fTreeTR);
+    fTreeTR->SetBranchAddress("TrackReferences", &fTrackReferences);
+    //
+    fNparticles = fStack->GetNtrack();
+    fNprimaries = fStack->GetNprimary();
+    printf("Event#: %5d Number of particles %5d \n", fEvent, fNparticles);
+    
+    return kTRUE;
+}
+
+Int_t AliMCEvent::GetParticleAndTR(Int_t i, TParticle*& particle, TClonesArray*& trefs)
+{
+    // Retrieve entry i
+    if (i > -1 && i < fNparticles) {
+       fTreeTR->GetEntry(fStack->TreeKEntry(i));
+    } else {
+       printf("AliMCEvent::GetEntry: Index out of range");
+    }
+    particle = fStack->Particle(i);
+    trefs    = fTrackReferences;
+    printf("Returning %5d \n", particle->GetPdgCode());
+    
+    return trefs->GetEntries();
+    
+}
+
+void AliMCEvent::DrawCheck(Int_t i)
+{
+    // Retrieve entry i and draw momentum vector and hits
+    if (i > -1 && i < fNparticles) {
+       fTreeTR->GetEntry(fStack->TreeKEntry(i));
+    } else {
+       printf("AliMCEvent::GetEntry: Index out of range");
+    }
+
+    TParticle* particle = fStack->Particle(i);
+    Int_t nh =  fTrackReferences->GetEntries();
+    
+    TH2F*    h = new TH2F("", "", 100, -500, 500, 100, -500, 500);
+    Float_t x0 = particle->Vx();
+    Float_t y0 = particle->Vy();
+
+    Float_t x1 = particle->Vx() + particle->Px() * 50.;
+    Float_t y1 = particle->Vy() + particle->Py() * 50.;
+    
+    TArrow*  a = new TArrow(x0, y0, x1, y1, 0.01);
+    h->Draw();
+    a->SetLineColor(2);
+    
+    a->Draw();
+    
+    for (Int_t ih = 0; ih < nh; ih++) {
+       AliTrackReference* ref = (AliTrackReference*) fTrackReferences->At(ih);
+       TMarker* m = new TMarker(ref->X(), ref->Y(), 20);
+       m->Draw();
+       m->SetMarkerSize(0.4);
+       
+    }
+}
+
+    
+
+                           
+Bool_t AliMCEvent::FinishEvent()
+{
+    // Dummy 
+    return kTRUE;
+}
+
+Bool_t AliMCEvent::Terminate()
+{ 
+    // Dummy 
+    return kTRUE;
+}
+
+Bool_t AliMCEvent::TerminateIO()
+{ 
+    // Dummy
+    return kTRUE;
+}
+    
diff --git a/STEER/AliMCEvent.h b/STEER/AliMCEvent.h
new file mode 100644 (file)
index 0000000..572d67b
--- /dev/null
@@ -0,0 +1,64 @@
+// -*- mode: C++ -*- 
+#ifndef ALIMCEVENT_H
+#define ALIMCEVENT_H
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+
+/* $Id$ */
+
+//-------------------------------------------------------------------------
+//                          Class AliMCEvent
+// This class gives access to MC truth during the analysis.
+// Monte Carlo truth is containe in the kinematics tree (produced particles) and 
+// the tree of reference hits.
+//      
+// Origin: Andreas Morsch, CERN, andreas.morsch@cern.ch 
+//-------------------------------------------------------------------------
+
+#include "AliVirtualEventHandler.h"
+class TFile;
+class TTree;
+class TParticle;
+class TClonesArray;
+class AliHeader;
+class AliStack;
+
+
+class AliMCEvent : public AliVirtualEventHandler 
+{
+public:
+    AliMCEvent();
+    AliMCEvent(const char* name, const char* title);
+    virtual ~AliMCEvent();
+    virtual void         SetOutputFileName(char* /* fname */) {;}
+    virtual char*        GetOutputFileName() {return 0;}
+    virtual Bool_t       InitIO(Option_t* opt);
+    virtual Bool_t       BeginEvent();
+    virtual Bool_t       FinishEvent();
+    virtual Bool_t       Terminate();
+    virtual Bool_t       TerminateIO();
+    //
+    AliStack* Stack()  {return fStack;}
+    TTree*    TreeTR() {return fTreeTR;}
+    Int_t     GetParticleAndTR(Int_t i, TParticle*& particle, TClonesArray*& trefs);
+    void      DrawCheck(Int_t i);
+           
+private:
+    TFile            *fFileE;            //! File with TreeE
+    TFile            *fFileK;            //! File with TreeK
+    TFile            *fFileTR;           //! File with TreeTR
+    TTree            *fTreeE;            //! TreeE  (Event Headers)
+    TTree            *fTreeK;            //! TreeK  (kinematics tree)
+    TTree            *fTreeTR;           //! TreeTR (track references tree)
+    AliStack         *fStack;            //! Current pointer to stack
+    AliHeader        *fHeader;           //! Current pointer to header
+    TClonesArray     *fTrackReferences;  //! Current list of tarck references
+    Int_t             fNEvent;           //! Number of events
+    Int_t             fEvent;            //! Current event
+    Int_t             fNprimaries;       //! Number of primaries
+    Int_t             fNparticles;       //! Number of particles 
+    ClassDef(AliMCEvent,1)  //MCEvent class 
+};
+#endif 
+
index 2fab564faea7ab91625844f2c4cffb8f4bae2b79..e0252d47d9a666e131c2a42eea371bb5bf745d30 100644 (file)
 
 #pragma link C++ class AliFstream+;
 #pragma link C++ class AliCTPRawData+;
+#pragma link C++ class AliMCEvent+;
 
 #endif
index 710ebaee6f8bf2624772bdff7356d9551b30deb1..a33ee3d9e6e4ae23a34751b1a4ae55e731b11893 100644 (file)
@@ -44,7 +44,8 @@ AliDCSGenDB.cxx \
 AliSurveyObj.cxx AliSurveyPoint.cxx \
 AliCodeTimer.cxx \
 AliFstream.cxx \
-AliCTPRawData.cxx
+AliCTPRawData.cxx \
+AliMCEvent.cxx