]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Including the physics selection and moving to AliAnalysisTaskSE
authorpchrist <pchrist@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sat, 2 Oct 2010 21:48:12 +0000 (21:48 +0000)
committerpchrist <pchrist@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sat, 2 Oct 2010 21:48:12 +0000 (21:48 +0000)
PWG2/EBYE/AliAnalysisTaskBF.cxx
PWG2/EBYE/AliAnalysisTaskBF.h
PWG2/EBYE/AliBalance.cxx
PWG2/EBYE/macros/AddTaskBalanceFunction.C
PWG2/EBYE/macros/runBalanceFunction.C

index dc4838901237720f7aec79970a81722366e17a65..b2b344f75490601383cd7a4934414bbf3afad5c6 100755 (executable)
@@ -1,12 +1,14 @@
 #include "TChain.h"
-#include "TTree.h"
+#include "TList.h"
 #include "TCanvas.h"
 #include "TLorentzVector.h"
 #include "TGraphErrors.h"
+#include "TH1F.h"
 
-#include "AliAnalysisTask.h"
+#include "AliAnalysisTaskSE.h"
 #include "AliAnalysisManager.h"
 
+#include "AliESDVertex.h"
 #include "AliESDEvent.h"
 #include "AliESDInputHandler.h"
 #include "AliAODEvent.h"
@@ -27,66 +29,51 @@ ClassImp(AliAnalysisTaskBF)
 
 //________________________________________________________________________
 AliAnalysisTaskBF::AliAnalysisTaskBF(const char *name) 
-  : AliAnalysisTask(name, ""), fESD(0), fAOD(0), fMC(0), fBalance(0) {
+  : AliAnalysisTaskSE(name), 
+    fBalance(0),
+    fList(0),
+    fHistEventStats(0) {
   // Constructor
 
   // Define input and output slots here
   // Input slot #0 works with a TChain
   DefineInput(0, TChain::Class());
   // Output slot #0 writes into a TH1 container
-  DefineOutput(0, AliBalance::Class());
+  DefineOutput(1, AliBalance::Class());
+  DefineOutput(2, TList::Class());
 }
 
 //________________________________________________________________________
-void AliAnalysisTaskBF::ConnectInputData(Option_t *) {
-  // Connect ESD or AOD here
+void AliAnalysisTaskBF::UserCreateOutputObjects() {
+  // Create histograms
   // Called once
-  TString gAnalysisLevel = fBalance->GetAnalysisLevel();
-
-  TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
-  if (!tree) {
-    Printf("ERROR: Could not read chain from input slot 0");
-  } 
-  else {
-    if(gAnalysisLevel == "ESD") {
-      AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
-      
-      if (!esdH) {
-       Printf("ERROR: Could not get ESDInputHandler");
-      } 
-      else
-       fESD = esdH->GetEvent();
-    }//ESD
-    else if(gAnalysisLevel == "AOD") {
-      AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
-      
-      if (!aodH) {
-       Printf("ERROR: Could not get AODInputHandler");
-      } else
-       fAOD = aodH->GetEvent();
-    }//AOD
-    else if(gAnalysisLevel == "MC") {
-      AliMCEventHandler* mcH = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
-      if (!mcH) {
-       Printf("ERROR: Could not retrieve MC event handler");
-      }
-      else
-       fMC = mcH->MCEvent();
-    }//MC
-    else
-      Printf("Wrong analysis type: Only ESD, AOD and MC types are allowed!");
+  if(!fBalance) {
+    fBalance = new AliBalance();
+    fBalance->SetAnalysisLevel("ESD");
+    fBalance->SetAnalysisType(1);
+    fBalance->SetNumberOfBins(18);
+    fBalance->SetInterval(-0.9,0.9);
   }
-}
 
-//________________________________________________________________________
-void AliAnalysisTaskBF::CreateOutputObjects() {
-  // Create histograms
-  // Called once
-  
+  fList = new TList();
+  fList->SetName("listQA");
+
+  TString gCutName[4] = {"Total","Offline trigger",
+                         "Vertex","Analyzed"};
+  fHistEventStats = new TH1F("fHistEventStats",
+                             "Event statistics;;N_{events}",
+                             4,0.5,4.5);
+  for(Int_t i = 1; i <= 4; i++)
+    fHistEventStats->GetXaxis()->SetBinLabel(i,gCutName[i-1].Data());
+  fList->Add(fHistEventStats);
+
+  // Post output data.
+  PostData(1, fBalance);
+  PostData(2, fList);
 }
 
 //________________________________________________________________________
-void AliAnalysisTaskBF::Exec(Option_t *) {
+void AliAnalysisTaskBF::UserExec(Option_t *) {
   // Main loop
   // Called for each event
   TString gAnalysisLevel = fBalance->GetAnalysisLevel();
@@ -95,31 +82,45 @@ void AliAnalysisTaskBF::Exec(Option_t *) {
   
   //ESD analysis
   if(gAnalysisLevel == "ESD") {
-    if (!fESD) {
-      Printf("ERROR: fESD not available");
+    AliESDEvent* gESD = dynamic_cast<AliESDEvent*>(InputEvent()); // from TaskSE
+    if (!gESD) {
+      Printf("ERROR: gESD not available");
       return;
     }
-    
-    Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
-    for (Int_t iTracks = 0; iTracks < fESD->GetNumberOfTracks(); iTracks++) {
-      AliESDtrack* track = fESD->GetTrack(iTracks);
-      if (!track) {
-       Printf("ERROR: Could not receive track %d", iTracks);
-       continue;
-      }
-      array->Add(track);
-    } //track loop
+
+    fHistEventStats->Fill(1); //all events
+    Bool_t isSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
+    if(isSelected) {
+      fHistEventStats->Fill(2); //triggered events
+
+      const AliESDVertex *vertex = gESD->GetPrimaryVertex();
+      if(vertex) {
+       fHistEventStats->Fill(3); //events with a proper vertex
+       fHistEventStats->Fill(4); //analayzed events
+       
+       Printf("There are %d tracks in this event", gESD->GetNumberOfTracks());
+       for (Int_t iTracks = 0; iTracks < gESD->GetNumberOfTracks(); iTracks++) {
+         AliESDtrack* track = gESD->GetTrack(iTracks);
+         if (!track) {
+           Printf("ERROR: Could not receive track %d", iTracks);
+           continue;
+         }
+         array->Add(track);
+       } //track loop
+      }//vertex object valid
+    }//triggered event 
   }//ESD analysis
   //AOD analysis
   else if(gAnalysisLevel == "AOD") {
-    if (!fAOD) {
-      Printf("ERROR: fAOD not available");
+    AliAODEvent* gAOD = dynamic_cast<AliAODEvent*>(InputEvent()); // from TaskSE
+    if(!gAOD) {
+      Printf("ERROR: gAOD not available");
       return;
     }
-    
-    Printf("There are %d tracks in this event", fAOD->GetNumberOfTracks());
-    for (Int_t iTracks = 0; iTracks < fAOD->GetNumberOfTracks(); iTracks++) {
-      AliAODTrack* track = fAOD->GetTrack(iTracks);
+
+    Printf("There are %d tracks in this event", gAOD->GetNumberOfTracks());
+    for (Int_t iTracks = 0; iTracks < gAOD->GetNumberOfTracks(); iTracks++) {
+      AliAODTrack* track = gAOD->GetTrack(iTracks);
       if (!track) {
        Printf("ERROR: Could not receive track %d", iTracks);
        continue;
@@ -129,14 +130,16 @@ void AliAnalysisTaskBF::Exec(Option_t *) {
   }//AOD analysis
   //MC analysis
   else if(gAnalysisLevel == "MC") {
-    if (!fMC) {
-      Printf("ERROR: fMC not available");
+
+    AliMCEvent*  mcEvent = MCEvent(); 
+    if (!mcEvent) {
+      Printf("ERROR: mcEvent not available");
       return;
     }
     
-    Printf("There are %d tracks in this event", fMC->GetNumberOfPrimaries());
-    for (Int_t iTracks = 0; iTracks < fMC->GetNumberOfPrimaries(); iTracks++) {
-      AliMCParticle* track = dynamic_cast<AliMCParticle *>(fMC->GetTrack(iTracks));
+    Printf("There are %d tracks in this event", mcEvent->GetNumberOfPrimaries());
+    for (Int_t iTracks = 0; iTracks < mcEvent->GetNumberOfPrimaries(); iTracks++) {
+      AliMCParticle* track = dynamic_cast<AliMCParticle *>(mcEvent->GetTrack(iTracks));
       if (!track) {
        Printf("ERROR: Could not receive particle %d", iTracks);
        continue;
@@ -148,17 +151,14 @@ void AliAnalysisTaskBF::Exec(Option_t *) {
   fBalance->CalculateBalance(array);
   
   delete array;
-
-  // Post output data.
-  PostData(0, fBalance);
+  
 }      
 
 //________________________________________________________________________
 void AliAnalysisTaskBF::Terminate(Option_t *) {
   // Draw result to the screen
   // Called once at the end of the query
-
-  fBalance = dynamic_cast<AliBalance*> (GetOutputData(0));
+  fBalance = dynamic_cast<AliBalance*> (GetOutputData(1));
   if (!fBalance) {
     Printf("ERROR: fBalance not available");
     return;
index d921edc84e996b36dd78894e1a475937ce6cf107..130b9e5169029ec5cd23fbc5f096ac43f4bb7403 100755 (executable)
@@ -4,37 +4,37 @@
 // Analysis task for the BF code\r
 // Authors: Panos Cristakoglou@cern.ch\r
 \r
+class TList;\r
+class TH1F;\r
+\r
 class AliBalance;\r
 class AliESDEvent;\r
 class AliAODEvent;\r
 class AliMCEvent;\r
 \r
-#include "AliAnalysisTask.h"\r
+#include "AliAnalysisTaskSE.h"\r
 \r
-class AliAnalysisTaskBF : public AliAnalysisTask {\r
+class AliAnalysisTaskBF : public AliAnalysisTaskSE {\r
  public:\r
   AliAnalysisTaskBF(const char *name = "AliAnalysisTaskBF");\r
   virtual ~AliAnalysisTaskBF() {}\r
   \r
-  virtual void   ConnectInputData(Option_t *);\r
-  virtual void   CreateOutputObjects();\r
-  virtual void   Exec(Option_t *option);\r
+  virtual void   UserCreateOutputObjects();\r
+  virtual void   UserExec(Option_t *option);\r
   virtual void   Terminate(Option_t *);\r
 \r
   void SetAnalysisObject(AliBalance *const analysis) {\r
     fBalance = analysis;}\r
   \r
  private:\r
-  AliESDEvent *fESD;    //ESD object\r
-  AliAODEvent *fAOD;    //AOD object\r
-  AliMCEvent  *fMC;     //MC object\r
-\r
   AliBalance *fBalance; //BF object\r
-   \r
+  TList *fList; //fList object\r
+  TH1F *fHistEventStats; //event stats\r
+\r
   AliAnalysisTaskBF(const AliAnalysisTaskBF&); // not implemented\r
   AliAnalysisTaskBF& operator=(const AliAnalysisTaskBF&); // not implemented\r
   \r
-  ClassDef(AliAnalysisTaskBF, 1); // example of analysis\r
+  ClassDef(AliAnalysisTaskBF, 2); // example of analysis\r
 };\r
 \r
 #endif\r
index 3c0882b90e07e3c9999782e72719d04515599800..e21fd2edcff6ce62ae8082a10f94193737c189c5 100644 (file)
@@ -224,6 +224,7 @@ void AliBalance::CalculateBalance(TObjArray *gTrackArray) {
   AliVParticle* track1 = 0;
   AliVParticle* track2 = 0;
     
+  //Printf("(AliBalance) Number of tracks: %d",gTrackArray->GetEntries());
   Int_t gNtrack = gTrackArray->GetEntries();
   for(i = 0; i < gNtrack; i++) {
     if(fAnalysisLevel == "ESD")
index 67f1133dbfd90b4f6e171a4b90612d76f9e65953..15123938eb0561305a226a060e604c1d89911be0 100644 (file)
@@ -32,9 +32,11 @@ AliAnalysisTaskBF *AddTaskBalanceFunction() {
   //==============================================================================
   TString outputFileName = AliAnalysisManager::GetCommonFileName();
   outputFileName += ":PWG2EbyE.outputBalanceFunctionAnalysis.root";
-  AliAnalysisDataContainer *cout_bf = mgr->CreateContainer("bfOutput", AliBalance::Class(),AliAnalysisManager::kOutputContainer,outputFileName.Data());
+  AliAnalysisDataContainer *coutBF = mgr->CreateContainer("bfOutput", AliBalance::Class(),AliAnalysisManager::kOutputContainer,outputFileName.Data());
+  AliAnalysisDataContainer *coutQA = mgr->CreateContainer("listQA", TList::Class(),AliAnalysisManager::kOutputContainer,outputFileName.Data());
   mgr->ConnectInput(taskBF, 0, mgr->GetCommonInputContainer());
-  mgr->ConnectOutput(taskBF, 0, cout_bf);
+  mgr->ConnectOutput(taskBF, 1, coutBF);
+  mgr->ConnectOutput(taskBF, 2, coutQA);
 
   // Return task pointer at the end
   return taskBF;
index 7b4495cfdab3d636e91e421c895bbb2a0a1e9ab6..0f4ac14473bcd9e82d9299da91cd407b727df752 100755 (executable)
@@ -22,7 +22,7 @@ void runBalanceFunction(Int_t analysisMode = kLocal,
   if (analysisMode==kLocal || analysisMode == kLocalPAR) {
     TChain *chain = new TChain("esdTree");
     chain->Add("Set1/AliESDs.root");
-    chain->Add("Set2/AliESDs.root");
+    //chain->Add("Set2/AliESDs.root");
   }
 
   //___________________________________________________//
@@ -70,7 +70,7 @@ void runBalanceFunction(Int_t analysisMode = kLocal,
   timer.Print();
 }
 
-void runLocal(const char* mode) {
+/*void runLocal(const char* mode) {
   //____________________________________________________//
   //_____________Setting up the par files_______________//
   //____________________________________________________//
@@ -108,7 +108,7 @@ void runLocal(const char* mode) {
   // Add task
   mgr->AddTask(taskBF);
 
-  // Create containers for input/output
+  // Create conainers for input/output
   AliAnalysisDataContainer *cinput = mgr->GetCommonInputContainer();
   AliAnalysisDataContainer *coutput = mgr->CreateContainer("contBF", AliBalance::Class(), AliAnalysisManager::kOutputContainer, "BF.ESD.root");
 
@@ -125,7 +125,7 @@ void runLocal(const char* mode) {
   mgr->PrintStatus();
 
   mgr->StartAnalysis("local", chain);
-}
+  }*/
 
 //__________________________________________________________//
 Int_t setupPar(const char* pararchivename) {