]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Add histogram counting events with tracks or EMCal clusters in different bunch crossings
authorgconesab <gconesab@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 19 Nov 2012 15:26:10 +0000 (15:26 +0000)
committergconesab <gconesab@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 19 Nov 2012 15:26:10 +0000 (15:26 +0000)
PWG/CaloTrackCorrBase/AliAnaCaloTrackCorrMaker.cxx
PWG/CaloTrackCorrBase/AliAnaCaloTrackCorrMaker.h
PWG/CaloTrackCorrBase/AliCaloTrackReader.cxx
PWG/CaloTrackCorrBase/AliCaloTrackReader.h

index ee92da33f01295951b34e83bc06e0941106d5ed3..967a0655a7934e695ab7717073358fcfd22df51d 100755 (executable)
@@ -53,7 +53,10 @@ fhPileUpClusterMult(0),       fhPileUpClusterMultAndSPDPileUp(0),
 fh2PileUpClusterMult(0),      fh2PileUpClusterMultAndSPDPileUp(0),
 fhTrackMult(0),
 fhCentrality(0),              fhEventPlaneAngle(0),
-fhNMergedFiles(0),            fhScaleFactor(0)
+fhNMergedFiles(0),            fhScaleFactor(0),
+fhEMCalBCEvent(0),            fhEMCalBCEventCut(0),
+fhTrackBCEvent(0),            fhTrackBCEventCut(0)
+
 {
   //Default Ctor
   if(fAnaDebug > 1 ) printf("*** Analysis Maker Constructor *** \n");
@@ -82,7 +85,12 @@ fhTrackMult(maker.fhTrackMult),
 fhCentrality(maker.fhCentrality),
 fhEventPlaneAngle(maker.fhEventPlaneAngle),
 fhNMergedFiles(maker.fhNMergedFiles),          
-fhScaleFactor(maker.fhScaleFactor)
+fhScaleFactor(maker.fhScaleFactor),
+fhEMCalBCEvent(maker.fhEMCalBCEvent),
+fhEMCalBCEventCut(maker.fhEMCalBCEventCut),
+fhTrackBCEvent(maker.fhTrackBCEvent),
+fhTrackBCEventCut(maker.fhTrackBCEventCut)
+
 {
   // cpy ctor
 }
@@ -195,9 +203,37 @@ TList *AliAnaCaloTrackCorrMaker::GetOutputContainer()
   fhNPileUpEvents->GetXaxis()->SetBinLabel(6 ,"!EMCal && SPD");
   fhNPileUpEvents->GetXaxis()->SetBinLabel(7 ,"EMCal && !SPD");
   fhNPileUpEvents->GetXaxis()->SetBinLabel(8 ,"!EMCal && !SPD");
-
   fOutputContainer->Add(fhNPileUpEvents);
   
+  fhTrackBCEvent      = new TH1I("hTrackBCEvent",   "Number of events with at least 1 track in a bunch crossing ", 19 , 0 , 19 ) ;
+  fhTrackBCEvent->SetYTitle("# events");
+  fhTrackBCEvent->SetXTitle("Bunch crossing");
+  for(Int_t i = 1; i < 20; i++)
+    fhTrackBCEvent->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
+  fOutputContainer->Add(fhTrackBCEvent);
+  
+  fhTrackBCEventCut      = new TH1I("hTrackBCEventCut",   "Number of events with at least 1 track in a bunch crossing ", 19 , 0 , 19 ) ;
+  fhTrackBCEventCut->SetYTitle("# events");
+  fhTrackBCEventCut->SetXTitle("Bunch crossing");
+  for(Int_t i = 1; i < 20; i++)
+    fhTrackBCEventCut->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
+  fOutputContainer->Add(fhTrackBCEventCut);
+
+  
+  fhEMCalBCEvent      = new TH1I("hEMCalBCEvent",   "Number of events with at least 1 cluster in a bunch crossing ", 19 , 0 , 19 ) ;
+  fhEMCalBCEvent->SetYTitle("# events");
+  fhEMCalBCEvent->SetXTitle("Bunch crossing");
+  for(Int_t i = 1; i < 20; i++)
+    fhEMCalBCEvent->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
+  fOutputContainer->Add(fhEMCalBCEvent);
+  
+  fhEMCalBCEventCut      = new TH1I("hEMCalBCEventCut",   "Number of events with at least 1 cluster in a bunch crossing", 19 , 0 , 19 ) ;
+  fhEMCalBCEventCut->SetYTitle("# events");
+  fhEMCalBCEventCut->SetXTitle("Bunch crossing");
+  for(Int_t i = 1; i < 20; i++)
+    fhEMCalBCEventCut->GetXaxis()->SetBinLabel(i ,Form("%d",i-10));
+  fOutputContainer->Add(fhEMCalBCEventCut);
+  
   fhZVertex      = new TH1F("hZVertex", " Z vertex distribution"   , 200 , -50 , 50  ) ;
   fhZVertex->SetXTitle("v_{z} (cm)");
   fOutputContainer->Add(fhZVertex);
@@ -478,12 +514,22 @@ void AliAnaCaloTrackCorrMaker::ProcessEvent(const Int_t iEntry,
     fhPileUpClusterMultAndSPDPileUp ->Fill(fReader->GetNPileUpClusters());
     fh2PileUpClusterMultAndSPDPileUp->Fill(fReader->GetNPileUpClusters(),fReader->GetNNonPileUpClusters());
   }
+  
   fhPileUpClusterMult ->Fill(fReader->GetNPileUpClusters  ());
   fh2PileUpClusterMult->Fill(fReader->GetNPileUpClusters  (),fReader->GetNNonPileUpClusters());
   fhTrackMult         ->Fill(fReader->GetTrackMultiplicity());
   fhCentrality        ->Fill(fReader->GetEventCentrality  ());
   fhEventPlaneAngle   ->Fill(fReader->GetEventPlaneAngle  ());
 
+  
+  for(Int_t i = 0; i < 19; i++)
+  {
+    if(fReader->GetTrackEventBC(i))   fhTrackBCEvent   ->Fill(i);
+    if(fReader->GetTrackEventBCcut(i))fhTrackBCEventCut->Fill(i);
+    if(fReader->GetEMCalEventBC(i))   fhEMCalBCEvent   ->Fill(i);
+    if(fReader->GetEMCalEventBCcut(i))fhEMCalBCEventCut->Fill(i);
+  }
+  
   Double_t v[3];
   fReader->GetInputEvent()->GetPrimaryVertex()->GetXYZ(v) ;
   fhZVertex->Fill(v[2]);
index dfaa008f1e45b28b61f2458a0dee0bb0494e7ca7..cfb4120e939597fbe6810ba98daf2864f18772f7 100755 (executable)
@@ -104,10 +104,14 @@ class AliAnaCaloTrackCorrMaker : public TObject {
   TH1F *   fhEventPlaneAngle;   //! Histogram with Event plane angle
   TH1I *   fhNMergedFiles;      //! Number of files merged
   TH1F *   fhScaleFactor;       //! Factor to scale histograms
+  TH1I *   fhEMCalBCEvent;      //! N events depending on the existance of a cluster in a given bunch crossing
+  TH1I *   fhEMCalBCEventCut;   //! N events depending on the existance of a cluster above acceptance and E cut in a given bunch crossing
+  TH1I *   fhTrackBCEvent;      //! N events depending on the existance of a track in a given bunch crossing
+  TH1I *   fhTrackBCEventCut;   //! N events depending on the existance of a track above acceptance and pt cut in a given bunch crossing
 
   AliAnaCaloTrackCorrMaker & operator = (const AliAnaCaloTrackCorrMaker & ) ; // cpy assignment
   
-  ClassDef(AliAnaCaloTrackCorrMaker,13)
+  ClassDef(AliAnaCaloTrackCorrMaker,14)
 } ;
  
 
index 7c1f2973106eaa17523d4de085db50dc9c7a9ccd..9f39297e275bf1c64fddee255f8778d7a2b0f795 100755 (executable)
@@ -69,6 +69,7 @@ fCTSPtMax(0),                fEMCALPtMax(0),                  fPHOSPtMax(0),
 fEMCALTimeCutMin(-10000),    fEMCALTimeCutMax(10000),
 fEMCALParamTimeCutMin(),     fEMCALParamTimeCutMax(),
 fUseParamTimeCut(kFALSE),
+fTrackTimeCutMin(-10000),    fTrackTimeCutMax(10000),
 fAODBranchList(0x0),
 fCTSTracks(0x0),             fEMCALClusters(0x0),             fPHOSClusters(0x0),
 fEMCALCells(0x0),            fPHOSCells(0x0),
@@ -1024,13 +1025,19 @@ void AliCaloTrackReader::FillInputCTS()
   Int_t nTracks = fInputEvent->GetNumberOfTracks() ;
   fTrackMult    = 0;
   Int_t nstatus = 0;
-  
-  for (Int_t itrack =  0; itrack <  nTracks; itrack++) 
+  for(Int_t i = 0; i < 19; i++)
+  {
+    fTrackBCEvent   [i] = 0;
+    fTrackBCEventCut[i] = 0;
+  }
+  for (Int_t itrack =  0; itrack <  nTracks; itrack++)
   {////////////// track loop
     AliVTrack * track = (AliVTrack*)fInputEvent->GetTrack(itrack) ; // retrieve track from esd
     
     //Select tracks under certain conditions, TPCrefit, ITSrefit ... check the set bits
-    if (fTrackStatus && !((track->GetStatus() & fTrackStatus) == fTrackStatus)) 
+    ULong_t status = track->GetStatus();
+
+    if (fTrackStatus && !((status & fTrackStatus) == fTrackStatus))
       continue ;
     
     nstatus++;
@@ -1105,26 +1112,48 @@ void AliCaloTrackReader::FillInputCTS()
     
     TLorentzVector momentum(pTrack[0],pTrack[1],pTrack[2],0);
     
-    if(fCTSPtMin < momentum.Pt() && fCTSPtMax > momentum.Pt())
+    Bool_t okTOF = ( (status & AliVTrack::kTOFout) == AliVTrack::kTOFout ) && ( (status & AliVTrack::kTIME) == AliVTrack::kTIME );
+    Double_t tof = -1000;
+    Int_t    bc  = -1000;
+    if(okTOF)
     {
-      if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"CTS")) 
-        continue;
+      tof = track->GetTOFsignal()*1e-3;
       
-      if(fDebug > 2 && momentum.Pt() > 0.1) 
-        printf("AliCaloTrackReader::FillInputCTS() - Selected tracks E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
-               momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
+      //printf("track TOF %e\n",tof);
+      bc = TMath::Nint((tof-25)/50) + 9;
+      //printf("track pt %f, tof %2.2f, bc=%d\n",track->Pt(),tof,bc);
       
-      if (fMixedEvent) 
-      {
-        track->SetID(itrack);
-      }
+      SetTrackEventBC(bc);
       
-      fCTSTracks->Add(track);        
+    }
+    
+    if(fCTSPtMin > momentum.Pt() || fCTSPtMax < momentum.Pt()) continue ;
+        
+    if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"CTS")) continue;
+    
+    if(okTOF)
+    {
+       SetTrackEventBCcut(bc);
       
-    }//Pt and Fiducial cut passed. 
+      //In any case, the time should to be larger than the fixed window ...
+      if( tof < fTrackTimeCutMin  || tof > fTrackTimeCutMax )
+      {
+        //printf("Remove track time %f\n",tof);
+        continue ;
+      }
+    }
+    
+    if(fDebug > 2 && momentum.Pt() > 0.1)
+      printf("AliCaloTrackReader::FillInputCTS() - Selected tracks E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
+             momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
+    
+    if (fMixedEvent)  track->SetID(itrack);
+
+    fCTSTracks->Add(track);
+    
   }// track loop
        
-  if(fDebug > 1) 
+  if(fDebug > 1)
     printf("AliCaloTrackReader::FillInputCTS()   - aod entries %d, input tracks %d, pass status %d, multipliticy %d\n", fCTSTracks->GetEntriesFast(), nTracks, nstatus, fTrackMult);//fCTSTracksNormalInputEntries);
   
 }
@@ -1221,24 +1250,31 @@ void AliCaloTrackReader::FillInputEMCALAlgorithm(AliVCluster * clus,
     clus->SetE(rdmEnergy);
   }
   
+  Double_t tof = clus->GetTOF()*1e9;
+
+  Int_t bc = TMath::Nint(tof/50) + 9;
+  //printf("tof %2.2f, bc+5=%d\n",tof,bc);
+  
+  SetEMCalEventBC(bc);
+  
+  if(fEMCALPtMin > clus->E() || fEMCALPtMax < clus->E()) return ;
+
   TLorentzVector momentum ;
   
-  clus->GetMomentum(momentum, fVertex[vindex]);      
+  clus->GetMomentum(momentum, fVertex[vindex]);
   
   if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"EMCAL")) return ;
   
-  if(fEMCALPtMin > momentum.E() || fEMCALPtMax < momentum.E())         return ;
-  
-  Double_t tof = clus->GetTOF()*1e9;
+  SetEMCalEventBCcut(bc);
 
-  if(!IsInTimeWindow(tof,momentum.E()))
+  if(!IsInTimeWindow(tof,clus->E()))
   {
     fNPileUpClusters++ ;
     return ;
   }
   else
     fNNonPileUpClusters++;
-
+  
   if(fDebug > 2 && momentum.E() > 0.1) 
     printf("AliCaloTrackReader::FillInputEMCAL() - Selected clusters E %3.2f, pt %3.2f, phi %3.2f, eta %3.2f\n",
            momentum.E(),momentum.Pt(),momentum.Phi()*TMath::RadToDeg(),momentum.Eta());
@@ -1265,6 +1301,11 @@ void AliCaloTrackReader::FillInputEMCAL()
   
   fNPileUpClusters    = 0; // Init counter
   fNNonPileUpClusters = 0; // Init counter
+  for(Int_t i = 0; i < 19; i++)
+  {
+    fEMCalBCEvent   [i] = 0;
+    fEMCalBCEventCut[i] = 0;
+  }
   
   //Loop to select clusters in fiducial cut and fill container with aodClusters
   if(fEMCALClustersListName=="")
@@ -1319,6 +1360,11 @@ void AliCaloTrackReader::FillInputEMCAL()
 
     fNPileUpClusters    = 0; // Init counter
     fNNonPileUpClusters = 0; // Init counter
+    for(Int_t i = 0; i < 19; i++)
+    {
+      fEMCalBCEvent   [i] = 0;
+      fEMCalBCEventCut[i] = 0;
+    }
     
     for (Int_t iclus =  0; iclus < fInputEvent->GetNumberOfCaloClusters(); iclus++)
     {
@@ -1329,8 +1375,6 @@ void AliCaloTrackReader::FillInputEMCAL()
         if (IsEMCALCluster(clus))
         {
           
-          if(fEMCALPtMin > clus->E() || fEMCALPtMax < clus->E()) continue ;
-
           Float_t  frac     =-1;
           Int_t    absIdMax = GetCaloUtils()->GetMaxEnergyCell(fEMCALCells, clus,frac);
           Double_t tof = clus->GetTOF();
@@ -1342,12 +1386,24 @@ void AliCaloTrackReader::FillInputEMCAL()
           //Reject clusters with bad channels, close to borders and exotic;
           if(!GetCaloUtils()->GetEMCALRecoUtils()->IsGoodCluster(clus,GetCaloUtils()->GetEMCALGeometry(),GetEMCALCells(),fInputEvent->GetBunchCrossNumber()))  continue;
 
+          Int_t bc = TMath::Nint(tof/50) + 9;
+          SetEMCalEventBC(bc);
+          
+          if(fEMCALPtMin > clus->E() || fEMCALPtMax < clus->E()) continue ;
+
+          TLorentzVector momentum ;
+          
+          clus->GetMomentum(momentum, fVertex[0]);
+          
+          if(fCheckFidCut && !fFiducialCut->IsInFiducialCut(momentum,"EMCAL")) return ;
+
+          SetEMCalEventBCcut(bc);
+
           if(!IsInTimeWindow(tof,clus->E()))
-          {
             fNPileUpClusters++ ;
-          }
           else
             fNNonPileUpClusters++;
+          
         }
       }
     }
index 9e11d56bc4627b65dd1dd70461a9fbfe60e93a1c..7f4bbf2fffef84aa23940e64f5e58e360a65ed49 100755 (executable)
@@ -120,6 +120,13 @@ public:
   
   //Time cut
   
+  Double_t         GetTrackTimeCutMin()              const { return fTrackTimeCutMin       ; }
+  Double_t         GetTrackTimeCutMax()              const { return fTrackTimeCutMax       ; }
+  
+  
+  void             SetTrackTimeCut(Double_t a, Double_t b) { fTrackTimeCutMin = a ;
+                                                             fTrackTimeCutMax = b          ; } // ns
+  
   Double_t         GetEMCALTimeCutMin()              const { return fEMCALTimeCutMin       ; }
   Double_t         GetEMCALTimeCutMax()              const { return fEMCALTimeCutMax       ; } 
 
@@ -258,7 +265,18 @@ public:
   
   Int_t            GetNPileUpClusters()                    { return  fNPileUpClusters     ; }
   Int_t            GetNNonPileUpClusters()                 { return  fNNonPileUpClusters  ; }
+  
+  Int_t            GetEMCalEventBC(Int_t bc)     const     { if(bc >=0 && bc < 19) return  fEMCalBCEvent   [bc] ; else return 0 ; }
+  Int_t            GetTrackEventBC(Int_t bc)     const     { if(bc >=0 && bc < 19) return  fTrackBCEvent   [bc] ; else return 0 ; }
+  Int_t            GetEMCalEventBCcut(Int_t bc)  const     { if(bc >=0 && bc < 19) return  fEMCalBCEventCut[bc] ; else return 0 ; }
+  Int_t            GetTrackEventBCcut(Int_t bc)  const     { if(bc >=0 && bc < 19) return  fTrackBCEventCut[bc] ; else return 0 ; }
+
+  void             SetEMCalEventBC(Int_t bc)               { if(bc >=0 && bc < 19) fEMCalBCEvent   [bc] = 1 ; }
+  void             SetTrackEventBC(Int_t bc)               { if(bc >=0 && bc < 19) fTrackBCEvent   [bc] = 1 ; }
+  void             SetEMCalEventBCcut(Int_t bc)            { if(bc >=0 && bc < 19) fEMCalBCEventCut[bc] = 1 ; }
+  void             SetTrackEventBCcut(Int_t bc)            { if(bc >=0 && bc < 19) fTrackBCEventCut[bc] = 1 ; }
 
+  
   // Track selection
   ULong_t          GetTrackStatus()                  const { return fTrackStatus          ; }
   void             SetTrackStatus(ULong_t bit)             { fTrackStatus = bit           ; }          
@@ -450,8 +468,10 @@ public:
   Float_t          fEMCALParamTimeCutMin[4];// Remove clusters/cells with time smaller than parametrized value, in ns
   Double_t         fEMCALParamTimeCutMax[4];// Remove clusters/cells with time larger than parametrized value, in ns
   Bool_t           fUseParamTimeCut;        // Use simple or parametrized time cut
+  Double_t         fTrackTimeCutMin;        // Remove tracks with time smaller than this value, in ns
+  Double_t         fTrackTimeCutMax;        // Remove tracks with time larger than this value, in ns
   
-  TList          * fAODBranchList ;         //-> List with AOD branches created and needed in analysis  
+  TList          * fAODBranchList ;         //-> List with AOD branches created and needed in analysis
   TObjArray      * fCTSTracks ;             //-> temporal array with tracks
   TObjArray      * fEMCALClusters ;         //-> temporal array with EMCAL CaloClusters
   TObjArray      * fPHOSClusters ;          //-> temporal array with PHOS  CaloClusters
@@ -529,6 +549,10 @@ public:
   Int_t            fNPileUpClusters;             // Number of clusters with time avobe 20 ns
   Int_t            fNNonPileUpClusters;          // Number of clusters with time below 20 ns
   Int_t            fNPileUpClustersCut;          // Cut to select event as pile-up
+  Int_t            fEMCalBCEvent[19];            // Fill one entry per event if there is a cluster in a given BC
+  Int_t            fEMCalBCEventCut[19];         // Fill one entry per event if there is a cluster in a given BC, depend on cluster E, acceptance cut
+  Int_t            fTrackBCEvent[19];            // Fill one entry per event if there is a track in a given BC
+  Int_t            fTrackBCEventCut[19];         // Fill one entry per event if there is a track in a given BC, depend on track pT, acceptance cut
 
   //Centrality/Event plane
   TString          fCentralityClass;        // Name of selected centrality class