]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Reading friends in analysis framework inside HLT
authorzampolli <chiara.zampolli@cern.ch>
Wed, 30 Jul 2014 11:29:31 +0000 (13:29 +0200)
committerzampolli <chiara.zampolli@cern.ch>
Wed, 30 Jul 2014 11:29:31 +0000 (13:29 +0200)
ANALYSIS/AliAnalysisManager.cxx
ANALYSIS/AliAnalysisManager.h
ANALYSIS/AliAnalysisTaskPt.cxx
ANALYSIS/AliAnalysisTaskPt.h
HLT/global/physics/AliHLTAnaManagerComponent.cxx
STEER/STEER/AliHLTTestInputHandler.cxx
STEER/STEER/AliHLTTestInputHandler.h
STEER/STEERBase/AliVEventHandler.h

index 5f025abdfab2042023b905339f06a3d59bbfd8af..281adc03afd2c6d653c976d772e6829a11eaa916 100644 (file)
@@ -2950,7 +2950,7 @@ void AliAnalysisManager::Changed()
 }
 
 //______________________________________________________________________________
-void AliAnalysisManager::InitInpuData(AliVEvent* esdEvent)
+void AliAnalysisManager::InitInputData(AliVEvent* esdEvent, AliESDfriend* esdFriend)
 {
 
 // Method to propagte to all the connected tasks the HLT event.
@@ -2961,7 +2961,7 @@ void AliAnalysisManager::InitInpuData(AliVEvent* esdEvent)
     TString classInputHandler = fInputEventHandler->ClassName();
     if (classInputHandler.Contains("HLT")){
       TObjArray* arrTasks = GetTasks();
-      fInputEventHandler->InitTaskInputData(esdEvent, arrTasks);
+      fInputEventHandler->InitTaskInputData(esdEvent, esdFriend, arrTasks);
     }
     else {
       Fatal("PropagateHLTEvent", "Input Handler not of type HLT, we cannot use this method!");
index 323337a5fa727d2d996c04cd533b6ef17eb51084..38179c0e794a81dc86d0be0f983e49e0c7cc4561 100644 (file)
@@ -37,6 +37,7 @@ class AliVEventPool;
 class AliAnalysisGrid;
 class AliAnalysisStatistics;
 class AliVEvent;
+class AliESDfriend;
 
 class AliAnalysisManager : public TNamed {
 
@@ -229,7 +230,7 @@ enum EAliAnalysisFlags {
    void                 Lock();
    void                 UnLock();
    void                 Changed();
-   void                 InitInpuData(AliVEvent* esdEvent);
+   void                 InitInputData(AliVEvent* esdEvent, AliESDfriend* esdFriend);
 protected:
    void                 CreateReadCache();
    void                 ImportWrappers(TList *source);
index 72fc2a332be196976f56ee190a70b6e59d4b293c..3f67553d30b56e9ba63bf793cb9cc29fdd0335e0 100644 (file)
@@ -12,6 +12,8 @@
 #include "AliVVtrack.h"
 #include "AliESDtrackCuts.h"
 #include "AliVEventHandler.h"
+#include "../TPC/Rec/AliTPCseed.h"
+#include "../TPC/Rec/AliTPCclusterMI.h"
 
 #include "AliAnalysisTaskPt.h"
 
@@ -22,7 +24,7 @@ ClassImp(AliAnalysisTaskPt)
 
 //________________________________________________________________________
 AliAnalysisTaskPt::AliAnalysisTaskPt(const char *name) 
-: AliAnalysisTask(name, ""), fESD(0), fHistPt(0), fCuts(0), fEv(0)
+: AliAnalysisTask(name, ""), fESD(0), fESDfriend(0), fHistPt(0), fCuts(0), fEv(0), fHistQ(0), fListOut(0)
 {
   // Constructor
 
@@ -30,7 +32,7 @@ AliAnalysisTaskPt::AliAnalysisTaskPt(const char *name)
   // Input slot #0 works with a TChain
   DefineInput(0, TChain::Class());
   // Output slot #0 writes into a TH1 container
-  DefineOutput(0, TH1F::Class());
+  DefineOutput(0, TList::Class());
 }
 
 //________________________________________________________________________
@@ -59,13 +61,17 @@ void AliAnalysisTaskPt::ConnectInputData(Option_t *)
       Printf("----> AliAnalysisTaskPt::ConnectInputData Getting the event from handler %p", esdH);
       //fESD = dynamic_cast<AliESDEvent*>(esdH->GetEvent());
       fESD = esdH->GetEvent();
+      fESDfriend = esdH->GetFriendEvent();
     }
     if (!fESD) {
-      Printf("ERROR, the dynamic cast did not work");
+      Printf("ERROR, no ESD event");
+    }
+    if (!fESDfriend) {
+      Printf("ERROR, no ESD friend");
     }
   }
 
-  Printf("fESD = %p", fESD);
+  Printf("fESD = %p, fESDfriend = %p", fESD, fESDfriend);
   printf("<---- AliAnalysisTaskPt::ConnectInputData\n");
 }
 
@@ -75,14 +81,26 @@ void AliAnalysisTaskPt::CreateOutputObjects()
   // Create histograms
   // Called once
 
+  fListOut = new TList();
+  fListOut->SetOwner();
+  fListOut->SetName("listHistos");
 
   fHistPt = new TH1F("fHistPt", "P_{T} distribution", 15, 0.1, 3.1);
   fHistPt->GetXaxis()->SetTitle("P_{T} (GeV/c)");
   fHistPt->GetYaxis()->SetTitle("dN/dP_{T} (c/GeV)");
   fHistPt->SetMarkerStyle(kFullCircle);
 
+  fHistQ = new TH1F("fHistQ", "TPC clusters: Q distribution", 1000, 0, 10000);
+  fHistQ->GetXaxis()->SetTitle("Q");
+  fHistQ->GetYaxis()->SetTitle("dN/dQ");
+  fHistQ->SetMarkerStyle(kFullCircle);
+
+  fListOut->Add(fHistPt);
+  fListOut->Add(fHistQ);
+
+  PostData(0, fListOut);
+
   fCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(1);
-  PostData(0, fHistPt);
 }
 
 //________________________________________________________________________
@@ -95,8 +113,13 @@ void AliAnalysisTaskPt::Exec(Option_t *)
     Printf("ERROR: fESD not available");
     return;
   }
+  if (!fESDfriend) {
+    Printf("ERROR: fESDfriend not available");
+    return;
+  }
 
   Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
+  Printf("... and there are %d friends in this event", fESDfriend->GetNumberOfTracks());
 
   // Track loop to fill a pT spectrum
   for (Int_t iTracks = 0; iTracks < fESD->GetNumberOfTracks(); iTracks++) {
@@ -109,8 +132,38 @@ void AliAnalysisTaskPt::Exec(Option_t *)
     fHistPt->Fill(track->Pt());
   } //track loop 
 
+
+  // Friend Track loop
+  for (Int_t iFriend = 0; iFriend < fESDfriend->GetNumberOfTracks(); iFriend++) {
+    AliESDfriendTrack* friendTrack = fESDfriend->GetTrack(iFriend);
+    if (!friendTrack) {
+      Printf("ERROR: Could not receive track %d", iFriend);
+      continue;
+    } 
+    TObject* calibObject;
+    AliTPCseed* seed = NULL;
+    for (Int_t idx = 0; (calibObject = friendTrack->GetCalibObject(idx)); ++idx) {
+      Printf(" |Cal %d = %p", idx, calibObject); 
+      if ((seed = dynamic_cast<AliTPCseed*>(calibObject))) {
+       Printf("Found TPC seed");
+       for (Int_t irow = 0; irow < 160; irow++){
+         AliTPCclusterMI* cluMI = seed->GetClusterPointer(irow);
+         if (cluMI){
+           Printf("Found cluster at row %d", irow);
+           Float_t q = cluMI->GetQ();
+           Printf("Q = %f", q);
+           fHistQ->Fill(q);
+         }
+         else {
+           Printf("Row %d does not contain clusters", irow);
+         }
+       }        
+      }
+    }    
+  }
+
   // Post output data.
-  PostData(0, fHistPt);
+  PostData(0, fListOut);
   fEv++;
 }      
 
@@ -122,13 +175,21 @@ void AliAnalysisTaskPt::Terminate(Option_t *)
 
   Printf("Terminate called: fESD = %p", fESD);
 
-  fHistPt = dynamic_cast<TH1F*> (GetOutputData(0));
-  if (!fHistPt) {
-    Printf("ERROR: fHistPt not available");
-    return;
-  }
+  fListOut = dynamic_cast<TList*> (GetOutputData(0)); 
+
+  if (fListOut) {
+    fHistPt = dynamic_cast<TH1F*>(fListOut->FindObject("fHistPt")); 
+    if (!fHistPt) {
+      Printf("ERROR: fHistPt not available");
+      return;
+    }
    
-  TCanvas *c1 = new TCanvas("AliAnalysisTaskPt","Pt",10,10,510,510);
-  c1->cd(1)->SetLogy();
-  fHistPt->DrawCopy("E");
+    TCanvas *c1 = new TCanvas("AliAnalysisTaskPt","Pt",10,10,510,510);
+    c1->cd(1)->SetLogy();
+    fHistPt->DrawCopy("E");
+  }
+  else {
+    Printf("In Terminate: no TList found");
+  }
+
 }
index 00c5a111cb620ea7300ab0d5360b547fa2d0b0e9..153b1a22ffd8bc847fedef28c414d2775d977bd7 100644 (file)
@@ -6,14 +6,16 @@
 
 class TH1F;
 class AliESDEvent;
+class AliESDfriend;
 class AliVVevent;
 class AliESDtrackCuts;
+class TList;
 
 #include "AliAnalysisTask.h"
 
 class AliAnalysisTaskPt : public AliAnalysisTask {
  public:
- AliAnalysisTaskPt() : AliAnalysisTask(), fESD(0), fHistPt(0), fCuts(0), fEv(0) {}
+ AliAnalysisTaskPt() : AliAnalysisTask(), fESD(0), fESDfriend(0), fHistPt(0), fCuts(0), fEv(0), fHistQ(0), fListOut(0) {}
   AliAnalysisTaskPt(const char *name);
   virtual ~AliAnalysisTaskPt() {}
   
@@ -23,11 +25,14 @@ class AliAnalysisTaskPt : public AliAnalysisTask {
   virtual void   Terminate(Option_t *);
 
  private:
-  AliVVevent *fESD;    //ESD object
-  TH1F        *fHistPt; //Pt spectrum
-  AliESDtrackCuts* fCuts;
+  AliVVevent*       fESD;          // ESD object
+  AliESDfriend*     fESDfriend;    // ESD friend object
+  TH1F*             fHistPt;       // Pt spectrum
+  AliESDtrackCuts*  fCuts;         // cuts
   Int_t fEv;
-   
+  TH1F*             fHistQ;        // TPC clusters Q spectrum
+  TList*            fListOut;      // output list
+
   AliAnalysisTaskPt(const AliAnalysisTaskPt&); // not implemented
   AliAnalysisTaskPt& operator=(const AliAnalysisTaskPt&); // not implemented
   
index fce8f64f1ebb8c5d098c60700c92b1ccd80fb023..98d5c0d8b03900cc46987730f7bb2ab2b2fb9213 100644 (file)
@@ -91,6 +91,7 @@ void AliHLTAnaManagerComponent::GetInputDataTypes( vector<AliHLTComponentDataTyp
   list.push_back(kAliHLTDataTypeESDObject|kAliHLTDataOriginAny);
   list.push_back(kAliHLTDataTypeClusters|kAliHLTDataOriginITSSPD);
   list.push_back(kAliHLTDataTypeESDContent|kAliHLTDataOriginVZERO);
+  list.push_back(kAliHLTDataTypeESDfriendObject|kAliHLTDataOriginAny);
 }
 
 // #################################################################################
@@ -111,11 +112,11 @@ void AliHLTAnaManagerComponent::GetOCDBObjectDescription( TMap* const targetMap)
   // see header file for class documentation
 
   if (!targetMap) return;
-  targetMap->Add(new TObjString("HLT/ConfigGlobal/MultiplicityCorrelations"),
+  /*  targetMap->Add(new TObjString("HLT/ConfigGlobal/MultiplicityCorrelations"),
                 new TObjString("configuration object"));
   targetMap->Add(new TObjString("HLT/ConfigGlobal/MultiplicityCorrelationsCentrality"),
                 new TObjString("centrality configuration object"));
-
+  */
   return;
 }
 
@@ -150,7 +151,7 @@ Int_t AliHLTAnaManagerComponent::DoInit( Int_t /*argc*/, const Char_t** /*argv*/
   fAnalysisManager->AddTask(task);
   AliAnalysisDataContainer *cinput  = fAnalysisManager->GetCommonInputContainer();
   Printf("Defining output file");
-  AliAnalysisDataContainer *coutput1 = fAnalysisManager->CreateContainer("pt", TH1::Class(),
+  AliAnalysisDataContainer *coutput1 = fAnalysisManager->CreateContainer("pt", TList::Class(),
       AliAnalysisManager::kOutputContainer, "Pt.ESD.root");
 
   //           connect containers
@@ -207,6 +208,7 @@ Int_t AliHLTAnaManagerComponent::DoEvent(const AliHLTComponentEventData& evtData
   // -- Get ESD object
   // -------------------
   AliESDEvent *esdEvent = NULL;
+  AliESDfriend *esdFriend = NULL;
   for ( const TObject *iter = GetFirstInputObject(kAliHLTDataTypeESDObject); iter != NULL; iter = GetNextInputObject() ) {
     esdEvent = dynamic_cast<AliESDEvent*>(const_cast<TObject*>( iter ) );
     if( !esdEvent ){ 
@@ -217,8 +219,17 @@ Int_t AliHLTAnaManagerComponent::DoEvent(const AliHLTComponentEventData& evtData
     esdEvent->GetStdContent();
   }
   printf("----> ESDEvent %p has %d tracks: \n", esdEvent, esdEvent->GetNumberOfTracks());
+  for ( const TObject *iter = GetFirstInputObject(kAliHLTDataTypeESDfriendObject); iter != NULL; iter = GetNextInputObject() ) {
+    esdFriend = dynamic_cast<AliESDfriend*>(const_cast<TObject*>( iter ) );
+    if( !esdFriend ){ 
+      HLTWarning("Wrong ESDFriend object received");
+      iResult = -1;
+      continue;
+    }
+  }
+  printf("----> ESDFriend %p has %d tracks: \n", esdFriend, esdFriend->GetNumberOfTracks());
 
-  fAnalysisManager->InitInpuData(esdEvent);
+  fAnalysisManager->InitInputData(esdEvent, esdFriend);
   //  fInputHandler->BeginEvent(0);
   fAnalysisManager->ExecAnalysis();
   fInputHandler->FinishEvent();
index c06a509f8fd74215fcb6f1044c65dcfdf421511a..eecc8e8bd6fc32b75b979323fdeb2587295781e4 100644 (file)
@@ -54,13 +54,14 @@ Bool_t AliHLTTestInputHandler::BeginEvent(Long64_t)
 }     
 
 //______________________________________________________________________________
-Bool_t AliHLTTestInputHandler::InitTaskInputData(AliVEvent* esdEvent, TObjArray* arrTasks) {
+Bool_t AliHLTTestInputHandler::InitTaskInputData(AliVEvent* esdEvent, AliESDfriend* friendEvent, TObjArray* arrTasks) {
 
 // Method to propagte to all the connected tasks the HLT event.
 // The method gets the list of tasks from the manager
 
   Printf("----> AliHLTTestInputHandler::InitTaskInpuData: Setting the event...");
   SetEvent(esdEvent);
+  SetFriendEvent(friendEvent);
   // set transient pointer to event inside tracks
   fEvent->ConnectTracks();
   Printf("----> AliHLTTestInputHandler::InitTaskInpuData: ...Event set");
index 83853b7850e107035262218794cc1bfef609a85d..2eea8f7a92018fc92dfc10cc06737391043f1e9e 100644 (file)
@@ -14,6 +14,7 @@
 
 class TObjArray;
 class AliVEvent;
+class AliESDfriend;
 
 class AliHLTTestInputHandler : public AliVEventHandler {
 
@@ -37,16 +38,20 @@ class AliHLTTestInputHandler : public AliVEventHandler {
     virtual Bool_t TerminateIO() {return kTRUE;}
 
     // Especially needed for HLT
-    Bool_t InitTaskInputData(AliVEvent* /*esdEvent*/, TObjArray* /*arrTasks*/);
+    Bool_t InitTaskInputData(AliVEvent* /*esdEvent*/, AliESDfriend* /*friendEvent*/, TObjArray* /*arrTasks*/);
 
     AliVEvent* GetEvent() const {return fEvent;}
     void  SetEvent(AliVEvent *event) {fEvent = event;}
+
+    AliESDfriend* GetFriendEvent() const {return fFriendEvent;}
+    void  SetFriendEvent(AliESDfriend *friendEvent) {fFriendEvent = friendEvent;}
       
  private:
     AliHLTTestInputHandler(const AliVEventHandler& handler);             
     AliHLTTestInputHandler& operator=(const AliVEventHandler& handler);  
     
-    AliVEvent *fEvent;    //! Pointer to the event
+    AliVEvent    *fEvent;          //! Pointer to the event
+    AliESDfriend *fFriendEvent;    //! Pointer to the friend event
 
     ClassDef(AliHLTTestInputHandler, 1);
 };
index 8b23f56bb5ab51bcc3c3a51ddddd9b550d656400..1cce6cd9173db4c9f9f453184d951e7766b3b5aa 100644 (file)
@@ -15,6 +15,7 @@
 class TTree;
 class TObjArray;
 class AliVEvent;
+class AliESDfriend;
 
 class AliVEventHandler : public TNamed {
 
@@ -60,10 +61,11 @@ enum EEventHandlerFlags {
     void                 Changed();
     virtual void         SetCacheSize(Long64_t) {}
     virtual TList        *GetUserInfo() const {return 0x0;};
-    // HLT
-    virtual Bool_t       InitTaskInputData(AliVEvent* /*event*/, TObjArray* /*arrTasks*/) {printf("OOOOPS!!!\n"); return kTRUE;};
-    virtual AliVEvent*   GetEvent() const {return 0x0;};
 
+    // HLT
+    virtual Bool_t          InitTaskInputData(AliVEvent* /*event*/, AliESDfriend* /*esdFriend*/, TObjArray* /*arrTasks*/) {printf("OOOOPS!!!\n"); return kTRUE;};
+    virtual AliVEvent*      GetEvent() const {return 0x0;};
+    virtual AliESDfriend*   GetFriendEvent() const {return 0x0;};
 
  private :
   ClassDef(AliVEventHandler, 1);