]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
including QA histograms and an additional output container to hold them
authormcosenti <mcosenti@cern.ch>
Thu, 13 Feb 2014 15:10:01 +0000 (13:10 -0200)
committermcosenti <mcosenti@cern.ch>
Thu, 13 Feb 2014 15:19:23 +0000 (13:19 -0200)
PWGGA/EMCALTasks/AliAnalysisTaskEMCALIsoPhoton.cxx
PWGGA/EMCALTasks/AliAnalysisTaskEMCALIsoPhoton.h
PWGGA/EMCALTasks/macros/AddTaskEMCALIsoPhoton.C

index 943fe8510a76eadf8cc88ca477f3026d46f70edb..14282db29516fc61baffe6626636e07955aa01e1 100644 (file)
@@ -71,6 +71,8 @@ AliAnalysisTaskEMCALIsoPhoton::AliAnalysisTaskEMCALIsoPhoton() :
   fHigherPtCone(0),
   fImportGeometryFromFile(0),
   fImportGeometryFilePath(""),
+  fMaxPtTrack(0),
+  fMaxEClus(0),
   fESD(0),
   fAOD(0),
   fMCEvent(0),
@@ -96,7 +98,18 @@ AliAnalysisTaskEMCALIsoPhoton::AliAnalysisTaskEMCALIsoPhoton() :
   fAllIsoEtMcGamma(0),
   fAllIsoNoUeEtMcGamma(0),
   fMCDirPhotonPtEtaPhiNoClus(0),
-  fHnOutput(0)
+  fHnOutput(0),
+  fQAList(0),
+  fNTracks(0),     
+  fEmcNCells(0),   
+  fEmcNClus(0),    
+  fEmcNClusCut(0), 
+  fNTracksECut(0), 
+  fEmcNCellsCut(0),
+  fEmcClusE(0),    
+  fEmcClusECut(0), 
+  fTrackPt(0),     
+  fTrackPtCut(0)   
 {
   // Default constructor.
   for(Int_t i = 0; i < 12;    i++)  fGeomMatrix[i] =  0;
@@ -132,6 +145,8 @@ AliAnalysisTaskEMCALIsoPhoton::AliAnalysisTaskEMCALIsoPhoton(const char *name) :
   fHigherPtCone(0),
   fImportGeometryFromFile(0),
   fImportGeometryFilePath(""),
+  fMaxPtTrack(0),
+  fMaxEClus(0),
   fESD(0),
   fAOD(0),
   fMCEvent(0),
@@ -157,7 +172,18 @@ AliAnalysisTaskEMCALIsoPhoton::AliAnalysisTaskEMCALIsoPhoton(const char *name) :
   fAllIsoEtMcGamma(0),
   fAllIsoNoUeEtMcGamma(0),
   fMCDirPhotonPtEtaPhiNoClus(0),
-  fHnOutput(0)
+  fHnOutput(0),
+  fQAList(0),
+  fNTracks(0),     
+  fEmcNCells(0),   
+  fEmcNClus(0),    
+  fEmcNClusCut(0), 
+  fNTracksECut(0), 
+  fEmcNCellsCut(0),
+  fEmcClusE(0),    
+  fEmcClusECut(0), 
+  fTrackPt(0),     
+  fTrackPtCut(0)   
 {
   // Constructor
 
@@ -167,6 +193,7 @@ AliAnalysisTaskEMCALIsoPhoton::AliAnalysisTaskEMCALIsoPhoton(const char *name) :
   // Output slot #0 id reserved by the base class for AOD
   // Output slot #1 writes into a TH1 container
   DefineOutput(1, TList::Class());
+  DefineOutput(2, TList::Class());
 }
 
 //________________________________________________________________________
@@ -260,9 +287,44 @@ void AliAnalysisTaskEMCALIsoPhoton::UserCreateOutputObjects()
   fHnOutput =  new THnSparseF("fHnOutput","Output matrix: E_{T},M02,CeIso,TrIso,AllIso, CeIsoNoUESub, AllIsoNoUESub, d#phi_{trk},d#eta_{trk},#eta_{clus},#phi_{clus},T_{max},mult,mc-p_{T}^{#gamma}", ndims, bins, xmin, xmax);
   fOutputList->Add(fHnOutput);
 
+  //QA outputs
+  fQAList = new TList();
+  fQAList->SetOwner();// Container cleans up all histos (avoids leaks in merging) 
 
+  fNTracks = new TH1F("hNTracks","# of selected tracks;n-tracks;counts",120,-0.5,119.5);
+  fNTracks->Sumw2();
+  fQAList->Add(fNTracks);
+
+  fEmcNCells = new TH1F("fEmcNCells",";n/event;count",120,-0.5,119.5);  
+  fEmcNCells->Sumw2();
+  fQAList->Add(fEmcNCells);
+  fEmcNClus = new TH1F("fEmcNClus",";n/event;count",120,-0.5,119.5);    
+  fEmcNClus->Sumw2();
+  fQAList->Add(fEmcNClus);
+  fEmcNClusCut = new TH1F("fEmcNClusCut",";n/event;count",120,-0.5,119.5); 
+  fEmcNClusCut->Sumw2();
+  fQAList->Add(fEmcNClusCut);
+  fNTracksECut = new TH1F("fNTracksECut",";n/event;count",120,-0.5,119.5); 
+  fNTracksECut->Sumw2();
+  fQAList->Add(fNTracksECut);
+  fEmcNCellsCut = new TH1F("fEmcNCellsCut",";n/event;count",120,-0.5,119.5);
+  fEmcNCellsCut->Sumw2();
+  fQAList->Add(fEmcNCellsCut);
+  fEmcClusE = new TH1F("fEmcClusE",";GeV;",100,-0.25,49.75);    
+  fEmcClusE->Sumw2();
+  fQAList->Add(fEmcClusE);
+  fEmcClusECut = new TH1F("fEmcClusECut",";GeV;",100,-0.25,49.75); 
+  fEmcClusECut->Sumw2();
+  fQAList->Add(fEmcClusECut);
+  fTrackPt = new TH1F("fTrackPt",";GeV;",100,-0.25,49.75);     
+  fTrackPt->Sumw2();
+  fQAList->Add(fTrackPt);
+  fTrackPtCut = new TH1F("fTrackPtCut",";GeV;",100,-0.25,49.75);   
+  fTrackPtCut->Sumw2();
+  fQAList->Add(fTrackPtCut);
 
   PostData(1, fOutputList);
+  PostData(2, fQAList);
 }
 
 //________________________________________________________________________
@@ -376,10 +438,15 @@ void AliAnalysisTaskEMCALIsoPhoton::UserExec(Option_t *)
     AliESDtrack *esdTrack = dynamic_cast<AliESDtrack*>(track);
     if (esdTrack && fPrTrCuts && fPrTrCuts->IsSelected(track)){
       fSelPrimTracks->Add(track);
+      /*if(fTrackMaxPt<track->Pt())
+       fTrackMaxPt = track->Pt();*/
       //printf("pt,eta,phi:%1.1f,%1.1f,%1.1f \n",track->Pt(),track->Eta(), track->Phi());
     }
-    else if(aodTrack)
+    else if(aodTrack){
       fSelPrimTracks->Add(track);
+      /*if(fTrackMaxPt<track->Pt())
+       fTrackMaxPt = track->Pt();*/
+    }
   }
 
   TObjArray *matEMCAL=(TObjArray*)fOADBContainer->GetObject(runnumber,"EmcalMatrices");
@@ -441,12 +508,16 @@ void AliAnalysisTaskEMCALIsoPhoton::UserExec(Option_t *)
   FillMcHists();
   if(fDebug)
     printf("passed calling of FillMcHists\n");
+  FillQA();
+  if(fDebug)
+    printf("passed calling of FillQA\n");
   /*if(fESD)
     fESDClusters->Clear();*/
   fSelPrimTracks->Clear();
   fNClusForDirPho = 0;
 
   PostData(1, fOutputList);
+  PostData(2, fQAList);
 }      
 
 //________________________________________________________________________
@@ -932,6 +1003,40 @@ Float_t AliAnalysisTaskEMCALIsoPhoton::GetMcPtSumInCone(Float_t etaclus, Float_t
   return ptsum;
 }
 //________________________________________________________________________
+void AliAnalysisTaskEMCALIsoPhoton::FillQA() 
+{
+  if(!fSelPrimTracks)
+    return;
+  const int ntracks = fSelPrimTracks->GetEntriesFast();
+  const int ncells = fESDCells->GetNumberOfCells();
+  const Int_t nclus = fESDClusters->GetEntries();
+
+  fNTracks->Fill(ntracks);
+  fEmcNCells->Fill(ncells);
+  fEmcNClus->Fill(nclus);
+  if(fMaxEClus>fECut){
+    fNTracksECut->Fill(ntracks);
+    fEmcNCellsCut->Fill(ncells);
+    fEmcNClusCut->Fill(nclus);
+  }
+  for(int it=0;it<ntracks;it++){
+    AliVTrack *t = (AliVTrack*)fSelPrimTracks->At(it);
+    if(!t)
+      continue;
+    fTrackPt->Fill(t->Pt());
+    if(fMaxEClus>fECut)
+      fTrackPtCut->Fill(t->Pt());
+  }
+  for(int ic=0;ic<nclus;ic++){
+    AliVCluster *c = (AliVCluster*)fESDClusters->At(ic);
+    if(!c)
+      continue;
+    fEmcClusE->Fill(c->E());
+    if(fMaxEClus>fECut)
+      fEmcClusECut->Fill(c->E());
+  }
+}
+//________________________________________________________________________
 void AliAnalysisTaskEMCALIsoPhoton::Terminate(Option_t *) 
 {
   // Called once at the end of the query.
index 9d5e604216351a50c330d377d599ce0c1d010838..4342500796d62edf495f259e715fe9f575e70187 100644 (file)
@@ -41,6 +41,7 @@ class AliAnalysisTaskEMCALIsoPhoton : public AliAnalysisTaskSE {
   void                   GetTrIso(TVector3 vec, Float_t &iso, Float_t &phiband, Float_t &core);
   void                   FillClusHists();
   void                   FillMcHists();
+  void                   FillQA();
   Float_t                GetClusSource(const AliVCluster *cluster);
   void                   FollowGamma();
   void                   GetDaughtersInfo(int firstd, int lastd, int selfid, const char *indputindent);
@@ -88,6 +89,8 @@ class AliAnalysisTaskEMCALIsoPhoton : public AliAnalysisTaskSE {
   Float_t                fHigherPtCone;          // higher pt inside the cone around the candidate
   Bool_t                 fImportGeometryFromFile;  // Import geometry settings in geometry.root file
   TString                fImportGeometryFilePath;  // path fo geometry.root file
+  Double_t               fMaxPtTrack;           //track with highest pt in event
+  Double_t               fMaxEClus;             //cluster with highest energy in event
 
   
  private:
@@ -120,6 +123,20 @@ class AliAnalysisTaskEMCALIsoPhoton : public AliAnalysisTaskSE {
   TH3F        *fMCDirPhotonPtEtaPhiNoClus; //!pt x eta x phi for prompt photons that didn't produce clusters
   THnSparse   *fHnOutput;                  //!Output matrix with 7 dimensions
 
+  //QA histos
+  TList       *fQAList;           //!output list holding QA histos
+  TH1F        *fNTracks;          //!number of tracks from Array->GetEntries()
+  TH1F        *fEmcNCells;        //!number of emcal cells in the event
+  TH1F        *fEmcNClus;         //!# of emcal clusters
+  TH1F        *fEmcNClusCut;      //!# of clusters in an event with at least 1 clus with E > fECut ("triggered event")
+  TH1F        *fNTracksECut;      //!number of tracks from Array->GetEntries() in "triggered event"
+  TH1F        *fEmcNCellsCut;     //!number of emcal cells in a in "triggered event"
+  TH1F        *fEmcClusE;         //!cluster E spectrum
+  TH1F        *fEmcClusECut;      //!cluster E spectrum in "triggered event"
+  TH1F        *fTrackPt;          //!selected tracks pt
+  TH1F        *fTrackPtCut;       //!selected tracks pt in "triggered event"
+
+
   AliAnalysisTaskEMCALIsoPhoton(const AliAnalysisTaskEMCALIsoPhoton&); // not implemented
   AliAnalysisTaskEMCALIsoPhoton& operator=(const AliAnalysisTaskEMCALIsoPhoton&); // not implemented
   
index c5b7f06420786025acd7cbb5f54d1e92a682099c..e5203aa96835e5d29f5bebb85d6c3acff4be7c3c 100644 (file)
@@ -46,9 +46,11 @@ AliAnalysisTaskEMCALIsoPhoton *AddTaskEMCALIsoPhoton(
   ana->SetGeoName(geoname.Data());  
   mgr->AddTask(ana);
   TString containername = "histosEMCALIsoPhoton";
+  TString containernameQA = "histosQA";
   if(pathstrsel != "/"){
     TString dirpth = (TSubString)pathstrsel.operator()(1,1);
     containername += dirpth;
+    containernameQA  += dirpth;
   }
   
   // Create ONLY the output containers for the data produced by the task.
@@ -57,9 +59,13 @@ AliAnalysisTaskEMCALIsoPhoton *AddTaskEMCALIsoPhoton(
   AliAnalysisDataContainer *coutput1 = mgr->CreateContainer(containername.Data(), 
                                                            TList::Class(),AliAnalysisManager::kOutputContainer,
                                                            Form("%s", AliAnalysisManager::GetCommonFileName()));
+  AliAnalysisDataContainer *coutput2 = mgr->CreateContainer(containernameQA.Data(), 
+                                                           TList::Class(),AliAnalysisManager::kOutputContainer,
+                                                           Form("%s", AliAnalysisManager::GetCommonFileName()));
   
   mgr->ConnectInput  (ana, 0, mgr->GetCommonInputContainer());
   mgr->ConnectOutput (ana, 1, coutput1 );
+  mgr->ConnectOutput (ana, 2, coutput2 );
    
   return ana;
 }