]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
adding cut study of primaries vs. secondaries
authorjgrosseo <jgrosseo@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 11 Jan 2008 07:39:33 +0000 (07:39 +0000)
committerjgrosseo <jgrosseo@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 11 Jan 2008 07:39:33 +0000 (07:39 +0000)
PWG0/esdTrackCuts/AliCutTask.cxx
PWG0/esdTrackCuts/AliCutTask.h
PWG0/esdTrackCuts/run.C

index 6a9112ab9e57a4cffbf1de7148b7cda6c81115c5..bbf4ca2cbefcb9e69c8c7a83cd6fa9145a7966e7 100644 (file)
@@ -4,6 +4,8 @@
 #include "TCanvas.h"
 #include "TList.h"
 #include "TFile.h"
+#include "TParticle.h"
+#include "TParticlePDG.h"
 
 #include "esdTrackCuts/AliESDtrackCuts.h"
 #include "AliPWG0Helper.h"
 #include "AliESDInputHandler.h"
 #include "AliESDVertex.h"
 
+#include "AliMCEventHandler.h"
+#include "AliMCEvent.h"
+#include "AliStack.h"
+#include "AliLog.h"
+
 #include "AliCutTask.h"
 
 // simple task that runs the esd track cuts to evaluate the basic plots created during the cuts
@@ -22,7 +29,8 @@ ClassImp(AliCutTask)
 
 //________________________________________________________________________
 AliCutTask::AliCutTask(const char *name) 
-  : AliAnalysisTask(name, ""), fESD(0), fTrackCuts(0), fVertex(0), fTriggerStats(0), fOutput(0)
+  : AliAnalysisTask(name, ""), fESD(0), fTrackCuts(0), fTrackCutsPrimaries(0),
+    fTrackCutsSecondaries(0), fVertex(0), fTriggerStats(0), fOutput(0)
 {
   // Constructor
 
@@ -77,6 +85,10 @@ void AliCutTask::CreateOutputObjects()
   fOutput->SetOwner();
 
   fOutput->Add(fTrackCuts);
+  if (fTrackCutsPrimaries)
+    fOutput->Add(fTrackCutsPrimaries);
+  if (fTrackCutsSecondaries)
+    fOutput->Add(fTrackCutsSecondaries);
 
   fVertex = new TH1F("fVertex", "fVertex;z vtx (cm);Count", 201, -20, 20);
   fOutput->Add(fVertex);
@@ -88,6 +100,9 @@ void AliCutTask::CreateOutputObjects()
   fTriggerStats->GetXaxis()->SetBinLabel(4, "ITS_SPD_GFO_L0");
   fTriggerStats->GetXaxis()->SetBinLabel(5, "VZERO_OR_LEFT | VZERO_OR_RIGHT");
   fOutput->Add(fTriggerStats);
+
+  // disable info messages of AliMCEvent (per event)
+  AliLog::SetClassDebugLevel("AliMCEvent", AliLog::kWarning - AliLog::kDebug + 1);
 }
 
 //________________________________________________________________________
@@ -115,14 +130,14 @@ void AliCutTask::Exec(Option_t *)
   if (AliPWG0Helper::IsEventTriggered(fESD->GetTriggerMask(), AliPWG0Helper::kMB2))
     fTriggerStats->Fill(2);
 
-  if (fESD->GetTriggerMask() & 32)
+  if (fESD->GetTriggerMask() & (0x1 << 14))
     fTriggerStats->Fill(3);
 
   if (fESD->GetTriggerMask() & 1 || fESD->GetTriggerMask() & 2)
     fTriggerStats->Fill(4);
 
   //if (!AliPWG0Helper::IsEventTriggered(fESD->GetTriggerMask()), AliPWG0Helper::kMB1)
-  if (fESD->GetTriggerMask() & 32 == 0)
+  if (fESD->GetTriggerMask() & (0x1 << 14) == 0)
     return;
 
   if (!AliPWG0Helper::IsVertexReconstructed(fESD->GetVertex()))
@@ -130,6 +145,73 @@ void AliCutTask::Exec(Option_t *)
 
   //Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
   fTrackCuts->CountAcceptedTracks(fESD);
+  
+  if (fTrackCutsPrimaries && fTrackCutsSecondaries)
+  {
+    AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
+    if (!eventHandler) {
+      Printf("ERROR: Could not retrieve MC event handler");
+      return;
+    }
+
+    AliMCEvent* mcEvent = eventHandler->MCEvent();
+    if (!mcEvent) {
+      Printf("ERROR: Could not retrieve MC event");
+      return;
+    }
+
+    AliStack* stack = mcEvent->Stack();
+    if (!stack)
+    {
+      Printf("ERROR: Stack not available");
+      return;
+    }
+    
+    for (Int_t trackId = 0; trackId < fESD->GetNumberOfTracks(); trackId++)
+    {
+      AliESDtrack* track = fESD->GetTrack(trackId);
+      if (!track)
+      {
+        Printf("UNEXPECTED: Could not get track with id %d", trackId);
+        continue;
+      }
+      
+      if (track->GetLabel() == 0)
+      {
+        Printf("Track %d has no label. Skipping", trackId);
+        continue;
+      }
+      
+      Int_t label = TMath::Abs(track->GetLabel());
+      if (stack->IsPhysicalPrimary(label) == kTRUE)
+      {
+        fTrackCutsPrimaries->AcceptTrack(track);
+      }
+      else
+      {
+        if (!stack->Particle(label)) {
+          Printf("ERROR: Could not get particle with label %d", label);
+          continue;
+        }
+        
+        // Deuteron is ok, but missing in Pdg particles in root
+        if (stack->Particle(label)->GetPdgCode() != 10010020) {
+          if (!stack->Particle(label)->GetPDG()) {
+            Printf("ERROR: Could not get PDG particle of particle with label %d (pdg code is %d)", label, stack->Particle(label)->GetPdgCode());
+            stack->Particle(label)->Print();
+            continue;
+          }
+          
+          if (stack->Particle(label)->GetPDG()->Charge() == 0) {
+            Printf("Particle %d has 0 charge. Skipping.", label);
+            continue;
+          }
+        }
+          
+        fTrackCutsSecondaries->AcceptTrack(track);
+      }
+    }
+  }
 
   // get the ESD vertex
   fVertex->Fill(fESD->GetVertex()->GetZv());
@@ -165,13 +247,26 @@ void AliCutTask::Terminate(Option_t *)
     return;
   }
 
+  fTrackCutsPrimaries = dynamic_cast<AliESDtrackCuts*> (fOutput->FindObject("fTrackCutsPrimaries"));
+  fTrackCutsSecondaries = dynamic_cast<AliESDtrackCuts*> (fOutput->FindObject("fTrackCutsSecondaries"));
+  
   TFile* file = TFile::Open("trackCuts.root", "RECREATE");
 
   fTrackCuts->SaveHistograms();
   fVertex->Write();
   fTriggerStats->Write();
+  
+  if (fTrackCutsPrimaries)
+    fTrackCutsPrimaries->SaveHistograms();
 
+  if (fTrackCutsSecondaries)
+    fTrackCutsSecondaries->SaveHistograms();
+  
   file->Close();
+  
+  Printf("Writting results to trackCuts.root.");
+  
+  return;
 
        fTrackCuts->DrawHistograms();
 
@@ -181,3 +276,18 @@ void AliCutTask::Terminate(Option_t *)
   new TCanvas;
   fTriggerStats->Draw();
 }
+
+void AliCutTask::EnableSecondaryStudy()
+{ 
+  // 
+  // clones the fTrackCuts
+  //
+  
+  if (!fTrackCuts) {
+    Printf("ERROR: fTrackCuts 0. Do not call this function before SetTrackCuts");
+    return;
+  }
+  
+  fTrackCutsPrimaries = dynamic_cast<AliESDtrackCuts*> (fTrackCuts->Clone("fTrackCutsPrimaries"));
+  fTrackCutsSecondaries = dynamic_cast<AliESDtrackCuts*> (fTrackCuts->Clone("fTrackCutsSecondaries"));
+}
index a7c889adffa81d6ba2e99ac30319b576125fe76a..238225e147029043c226c61922530c622934688c 100644 (file)
@@ -21,10 +21,14 @@ class AliCutTask : public AliAnalysisTask {
   virtual void   Terminate(Option_t *);
 
   void SetTrackCuts(AliESDtrackCuts* cuts) { fTrackCuts = cuts; }
+  void EnableSecondaryStudy();
   
  private:
   AliESDEvent *fESD;           //! ESD object
   AliESDtrackCuts* fTrackCuts; // track cuts
+  
+  AliESDtrackCuts* fTrackCutsPrimaries; // cuts for tracks from primary particles
+  AliESDtrackCuts* fTrackCutsSecondaries; // cuts for tracks from secondary particles
 
   TH1F* fVertex;   //! event z vertex distribution
   TH1F* fTriggerStats;  //! triggers
index 80f6920b5a265b37d1b32c3c43b93661660caf14..8d94bc12dccbc60105d68f0582c500c839f28d57 100644 (file)
@@ -1,4 +1,4 @@
-void run(Char_t* data, Int_t nRuns=20, Int_t offset=0, Bool_t aDebug = kFALSE, Bool_t aProof = kFALSE, const char* option = "")
+void run(Char_t* data, Int_t nRuns=20, Int_t offset=0, Bool_t aDebug = kFALSE, Bool_t aProof = kFALSE, Bool_t mc = kFALSE, const char* option = "")
 {
   if (aProof)
   {
@@ -57,6 +57,14 @@ void run(Char_t* data, Int_t nRuns=20, Int_t offset=0, Bool_t aDebug = kFALSE, B
   // Add ESD handler
   AliESDInputHandler* esdH = new AliESDInputHandler;
   mgr->SetInputEventHandler(esdH);
+  
+  if (mc) {
+    task->EnableSecondaryStudy();
+    // Enable MC event handler
+    AliMCEventHandler* handler = new AliMCEventHandler;
+    handler->SetReadTR(kFALSE);
+    mgr->SetMCtruthEventHandler(handler);
+  }
 
   // Attach input
   cInput  = mgr->CreateContainer("cInput", TChain::Class(), AliAnalysisManager::kInputContainer);