]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Add control histograms for time stamp, move histogram filling to separate method...
authorgconesab <gconesab@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 10 Dec 2012 09:57:03 +0000 (09:57 +0000)
committergconesab <gconesab@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 10 Dec 2012 09:57:03 +0000 (09:57 +0000)
PWG/CaloTrackCorrBase/AliAnaCaloTrackCorrMaker.cxx
PWG/CaloTrackCorrBase/AliAnaCaloTrackCorrMaker.h
PWG/CaloTrackCorrBase/AliCaloTrackReader.cxx
PWG/CaloTrackCorrBase/AliCaloTrackReader.h

index 7003cf7c13238fc303f2e83c04fec7d2751a36ae..9fa3d662440e97f0404c5c3e610707e94528935a 100755 (executable)
@@ -52,13 +52,12 @@ fScaleFactor(-1),
 fhNEvents(0),                 fhNPileUpEvents(0),
 fhZVertex(0),                 
 fhPileUpClusterMult(0),       fhPileUpClusterMultAndSPDPileUp(0),
-fh2PileUpClusterMult(0),      fh2PileUpClusterMultAndSPDPileUp(0),
 fhTrackMult(0),
 fhCentrality(0),              fhEventPlaneAngle(0),
 fhNMergedFiles(0),            fhScaleFactor(0),
 fhEMCalBCEvent(0),            fhEMCalBCEventCut(0),
 fhTrackBCEvent(0),            fhTrackBCEventCut(0),
-fhPrimaryVertexBC(0)
+fhPrimaryVertexBC(0),         fhTimeStampFraction(0)
 {
   //Default Ctor
   if(fAnaDebug > 1 ) printf("*** Analysis Maker Constructor *** \n");
@@ -81,8 +80,6 @@ fhNPileUpEvents(maker.fhNPileUpEvents),
 fhZVertex(maker.fhZVertex),    
 fhPileUpClusterMult(maker.fhPileUpClusterMult),
 fhPileUpClusterMultAndSPDPileUp(maker.fhPileUpClusterMultAndSPDPileUp),
-fh2PileUpClusterMult(maker.fh2PileUpClusterMult),
-fh2PileUpClusterMultAndSPDPileUp(maker.fh2PileUpClusterMultAndSPDPileUp),
 fhTrackMult(maker.fhTrackMult),
 fhCentrality(maker.fhCentrality),
 fhEventPlaneAngle(maker.fhEventPlaneAngle),
@@ -161,6 +158,77 @@ TList * AliAnaCaloTrackCorrMaker::FillAndGetAODBranchList()
        
 }
 
+//_________________________________________________________
+void AliAnaCaloTrackCorrMaker::FillControlHistograms()
+{
+  // Event control histograms
+  
+  fhNEvents        ->Fill(0); // Number of events analyzed
+  
+  if( fReader->IsPileUpFromSPD())
+    fhNPileUpEvents->Fill(0.5);
+  if( fReader->GetInputEvent()->IsPileupFromSPDInMultBins())
+    fhNPileUpEvents->Fill(1.5);
+  if( fReader->IsPileUpFromEMCal())
+    fhNPileUpEvents->Fill(2.5);
+  if( fReader->IsPileUpFromSPDOrEMCal() )
+    fhNPileUpEvents->Fill(3.5);
+  if( fReader->IsPileUpFromSPDAndEMCal() )
+    fhNPileUpEvents->Fill(4.5);
+  if( fReader->IsPileUpFromSPDAndNotEMCal() )
+    fhNPileUpEvents->Fill(5.5);
+  if( fReader->IsPileUpFromEMCalAndNotSPD() )
+    fhNPileUpEvents->Fill(6.5);
+  if( fReader->IsPileUpFromNotSPDAndNotEMCal() )
+    fhNPileUpEvents->Fill(7.5);
+  
+  if(fReader->IsPileUpFromSPD())
+    fhPileUpClusterMultAndSPDPileUp ->Fill(fReader->GetNPileUpClusters());
+  
+  fhPileUpClusterMult ->Fill(fReader->GetNPileUpClusters  ());
+  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]);
+  
+  Int_t primaryBC = -1000;
+  AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (fReader->GetInputEvent());
+  AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (fReader->GetInputEvent());
+  
+  if     (esdevent)
+    primaryBC = esdevent->GetPrimaryVertex()->GetBC();
+  else if(aodevent)
+    primaryBC = aodevent->GetPrimaryVertex()->GetBC();
+  
+  fhPrimaryVertexBC->Fill(primaryBC);
+  
+  // Time stamp
+  if(fReader->IsSelectEventTimeStampOn() && esdevent)
+  {
+    Int_t timeStamp = esdevent->GetTimeStamp();
+    Float_t timeStampFrac = 1.*(timeStamp-fReader->GetRunTimeStampMin()) /
+                               (fReader->GetRunTimeStampMax()-fReader->GetRunTimeStampMin());
+    
+    //printf("stamp %d, min %d, max %d, frac %f\n", timeStamp, fReader->GetRunTimeStampMin(), fReader->GetRunTimeStampMax(), timeStampFrac);
+
+    fhTimeStampFraction->Fill(timeStampFrac);
+  }
+  
+  
+}
+
 //_______________________________________________________
 TList * AliAnaCaloTrackCorrMaker::GetListOfAnalysisCuts()
 { 
@@ -256,16 +324,6 @@ TList *AliAnaCaloTrackCorrMaker::GetOutputContainer()
   fhPileUpClusterMultAndSPDPileUp->SetXTitle("# clusters");
   fOutputContainer->Add(fhPileUpClusterMultAndSPDPileUp);
   
-  fh2PileUpClusterMult = new TH2F("h2PileUpClusterMult", "Number of clusters per event with large time (|t| > 20 ns)" , 100 , 0 , 100 , 100 , 0 , 100 ) ;
-  fh2PileUpClusterMult->SetXTitle("# clusters (large t)");
-  fh2PileUpClusterMult->SetYTitle("# clusters (small t)");
-  fOutputContainer->Add(fh2PileUpClusterMult);
-  
-  fh2PileUpClusterMultAndSPDPileUp = new TH2F("h2PileUpClusterMultAndSPDPileUp", "Number of clusters per event with large time (|t| > 20 ns, events tagged as pile-up by SPD)" , 100 , 0 , 100 , 100 , 0 , 100) ;
-  fh2PileUpClusterMultAndSPDPileUp->SetXTitle("# clusters (large t)");
-  fh2PileUpClusterMultAndSPDPileUp->SetYTitle("# clusters (small t)");
-  fOutputContainer->Add(fh2PileUpClusterMultAndSPDPileUp);
-  
   fhCentrality   = new TH1F("hCentrality","Number of events in centrality bin",100,0.,100) ;
   fhCentrality->SetXTitle("Centrality bin");
   fOutputContainer->Add(fhCentrality) ;  
@@ -274,6 +332,13 @@ TList *AliAnaCaloTrackCorrMaker::GetOutputContainer()
   fhEventPlaneAngle->SetXTitle("EP angle (rad)");
   fOutputContainer->Add(fhEventPlaneAngle) ;
   
+  if(fReader->IsSelectEventTimeStampOn())
+  {
+    fhTimeStampFraction = new TH1F("hTimeStampFraction","Fraction of events within a given time stamp range",150, -1, 2) ;
+    fhTimeStampFraction->SetXTitle("fraction");
+    fOutputContainer->Add(fhTimeStampFraction) ;
+  }
+  
   if(fScaleFactor > 0)
   {
     fhNMergedFiles = new TH1F("hNMergedFiles",   "Number of merged output files"     , 1 , 0 , 1  ) ;
@@ -495,63 +560,8 @@ void AliAnaCaloTrackCorrMaker::ProcessEvent(const Int_t iEntry,
     }
   }
   
-  // Event control histograms
-  fhNEvents        ->Fill(0); // Number of events analyzed
-  
-  if( fReader->IsPileUpFromSPD()) 
-    fhNPileUpEvents->Fill(0.5);
-  if( fReader->GetInputEvent()->IsPileupFromSPDInMultBins()) 
-    fhNPileUpEvents->Fill(1.5);
-  if( fReader->IsPileUpFromEMCal())
-    fhNPileUpEvents->Fill(2.5);
-  if( fReader->IsPileUpFromSPDOrEMCal() )
-    fhNPileUpEvents->Fill(3.5);
-  if( fReader->IsPileUpFromSPDAndEMCal() )
-    fhNPileUpEvents->Fill(4.5);
-  if( fReader->IsPileUpFromSPDAndNotEMCal() )
-    fhNPileUpEvents->Fill(5.5);
-  if( fReader->IsPileUpFromEMCalAndNotSPD() )
-    fhNPileUpEvents->Fill(6.5);
-  if( fReader->IsPileUpFromNotSPDAndNotEMCal() )
-    fhNPileUpEvents->Fill(7.5);
-  
-  if(fReader->IsPileUpFromSPD())
-  {
-    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  ());
-
+  FillControlHistograms();
   
-  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]);
-  
-  Int_t primaryBC = -1000;
-  AliESDEvent* esdevent = dynamic_cast<AliESDEvent*> (fReader->GetInputEvent());
-  AliAODEvent* aodevent = dynamic_cast<AliAODEvent*> (fReader->GetInputEvent());
-
-  if     (esdevent)
-    primaryBC = esdevent->GetPrimaryVertex()->GetBC();
-  else if(aodevent)
-    primaryBC = aodevent->GetPrimaryVertex()->GetBC();
-
-  fhPrimaryVertexBC->Fill(primaryBC);
-
   //printf(">>>>>>>>>> AFTER >>>>>>>>>>>\n");
   //gObjectTable->Print();
        
index b0b303e00ab9a27c873e232ea431f1ad7719bd3e..2febc65fc7bd7116c838cd8e86c201c4481eda6a 100755 (executable)
@@ -34,6 +34,8 @@ class AliAnaCaloTrackCorrMaker : public TObject {
   
   void    AddAnalysis(TObject* ana, Int_t n) ;
 
+  void    FillControlHistograms();
+  
   TList * GetListOfAnalysisContainers() { return fAnalysisContainer ; }
   TList * GetListOfAnalysisCuts();
   TList * GetOutputContainer() ;
@@ -97,8 +99,6 @@ class AliAnaCaloTrackCorrMaker : public TObject {
   TH1F *   fhZVertex;           //! Vertex of accepted event
   TH1F *   fhPileUpClusterMult; //! N clusters with high time
   TH1F *   fhPileUpClusterMultAndSPDPileUp; //! N clusters with high time in events tagged as pile-up by SPD
-  TH2F *   fh2PileUpClusterMult; //! N clusters with high time vs N clusterd with small time
-  TH2F *   fh2PileUpClusterMultAndSPDPileUp; //! N clusters with high time vs N clusterd with small time in events tagged as pile-up by SPD
   TH1F *   fhTrackMult;         //! Number of tracks per event histogram
   TH1F *   fhCentrality;        //! Histogram with centrality bins
   TH1F *   fhEventPlaneAngle;   //! Histogram with Event plane angle
@@ -109,10 +109,11 @@ class AliAnaCaloTrackCorrMaker : public TObject {
   TH1F *   fhTrackBCEvent;      //! N events depending on the existance of a track in a given bunch crossing
   TH1F *   fhTrackBCEventCut;   //! N events depending on the existance of a track above acceptance and pt cut in a given bunch crossing
   TH1F *   fhPrimaryVertexBC;   //! Primary vertex BC
-
+  TH1F *   fhTimeStampFraction; //! event fraction depending on Time Stamp, only if activated on reader
+  
   AliAnaCaloTrackCorrMaker & operator = (const AliAnaCaloTrackCorrMaker & ) ; // cpy assignment
   
-  ClassDef(AliAnaCaloTrackCorrMaker,14)
+  ClassDef(AliAnaCaloTrackCorrMaker,15)
 } ;
  
 
index 68aaa8c3b69c6ca2d4359e8a33b4bf30de9e3f08..829403d6f5e426229295df71b274efd419373b50 100755 (executable)
@@ -455,7 +455,13 @@ void AliCaloTrackReader::InitParameters()
 
   // Parametrized time cut (LHC11c)
   //fEMCALParamTimeCutMin[0] =-5;   fEMCALParamTimeCutMin[1] =-1 ; fEMCALParamTimeCutMin[2] = 1.87; fEMCALParamTimeCutMin[3] = 0.4;   
-  //fEMCALParamTimeCutMax[0] = 3.5; fEMCALParamTimeCutMax[1] = 50; fEMCALParamTimeCutMax[2] = 0.15; fEMCALParamTimeCutMax[3] = 1.6;   
+  //fEMCALParamTimeCutMax[0] = 3.5; fEMCALParamTimeCutMax[1] = 50; fEMCALParamTimeCutMax[2] = 0.15; fEMCALParamTimeCutMax[3] = 1.6;
+  
+  fTimeStampRunMin = -1;
+  fTimeStampRunMax = 1e12;
+  fTimeStampEventFracMin = -1;
+  fTimeStampEventFracMax = 2;
+
 }
 
 //___________________________________________________________
@@ -696,7 +702,7 @@ Bool_t AliCaloTrackReader::FillInputEvent(const Int_t iEntry,
     Int_t timeStamp = esd->GetTimeStamp();
     Float_t timeStampFrac = 1.*(timeStamp-fTimeStampRunMin) / (fTimeStampRunMax-fTimeStampRunMin);
   
-    //printf("stamp0 %d, max0 %d, frac %f\n", timeStamp-fTimeStampRunMin,fTimeStampRunMax-fTimeStampRunMin, //timeStampFrac);
+    //printf("stamp0 %d, max0 %d, frac %f\n", timeStamp-fTimeStampRunMin,fTimeStampRunMax-fTimeStampRunMin, timeStampFrac);
     
     if(timeStampFrac < fTimeStampEventFracMin || timeStampFrac > fTimeStampEventFracMax) return kFALSE;
   
index 0b44272c70c330d7c70c14d57a57dc6b0ce4766d..db148da055edfdf3c8f4b81a753db2d3bc499953 100755 (executable)
@@ -284,6 +284,7 @@ public:
   void             SwitchOnSelectEventTimeStamp()          { fTimeStampEventSelect = kTRUE   ; }
   void             SwitchOffSelectEventTimeStamp()         { fTimeStampEventSelect = kFALSE  ; }
   
+  Bool_t           IsSelectEventTimeStampOn()              {return  fTimeStampEventSelect    ; }
   
   Bool_t           IsPileUpFromSPD()               const ;
   Bool_t           IsPileUpFromEMCal()             const ;