]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
first draft of physics event selection class
authorjgrosseo <jgrosseo@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 11 Dec 2009 17:31:11 +0000 (17:31 +0000)
committerjgrosseo <jgrosseo@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 11 Dec 2009 17:31:11 +0000 (17:31 +0000)
PWG0/AliPhysicsSelection.cxx [new file with mode: 0644]
PWG0/AliPhysicsSelection.h [new file with mode: 0644]
PWG0/AliTriggerAnalysis.cxx
PWG0/AliTriggerAnalysis.h
PWG0/PWG0baseLinkDef.h
PWG0/libPWG0base.pkg
PWG0/trigger/AliTriggerTask.cxx
PWG0/trigger/run.C

diff --git a/PWG0/AliPhysicsSelection.cxx b/PWG0/AliPhysicsSelection.cxx
new file mode 100644 (file)
index 0000000..a6267b5
--- /dev/null
@@ -0,0 +1,306 @@
+/* $Id: AliPhysicsSelection.cxx 35782 2009-10-22 11:54:31Z jgrosseo $ */
+
+/**************************************************************************
+ * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Author: The ALICE Off-line Project.                                    *
+ * Contributors are mentioned in the code where appropriate.              *
+ *                                                                        *
+ * Permission to use, copy, modify and distribute this software and its   *
+ * documentation strictly for non-commercial purposes is hereby granted   *
+ * without fee, provided that the above copyright notice appears in all   *
+ * copies and that both the copyright notice and this permission notice   *
+ * appear in the supporting documentation. The authors make no claims     *
+ * about the suitability of this software for any purpose. It is          *
+ * provided "as is" without express or implied warranty.                  *
+ **************************************************************************/
+
+//-------------------------------------------------------------------------
+//                      Implementation of   Class AliPhysicsSelection
+//   This class provides function to check if events have been triggered based on the data in the ESD
+//   The trigger bits, trigger class inputs and only the data (offline trigger) can be used
+//   Origin: Jan Fiete Grosse-Oetringhaus, CERN
+//-------------------------------------------------------------------------
+
+#include <Riostream.h>
+#include <TH1F.h>
+#include <TH2F.h>
+#include <TList.h>
+#include <TIterator.h>
+#include <TDirectory.h>
+
+#include <AliPhysicsSelection.h>
+
+#include <AliTriggerAnalysis.h>
+#include <AliLog.h>
+
+#include <AliESDEvent.h>
+
+ClassImp(AliPhysicsSelection)
+
+AliPhysicsSelection::AliPhysicsSelection() :
+  fTriggerAnalysis(0),
+  fHistStatistics(0),
+  fHistBunchCrossing(0)
+{
+  // constructor
+  
+  AliLog::SetClassDebugLevel("AliPhysicsSelection", AliLog::kWarning);
+}
+    
+AliPhysicsSelection::~AliPhysicsSelection()
+{
+  // destructor
+  
+  if (fTriggerAnalysis)
+  {
+    delete fTriggerAnalysis;
+    fTriggerAnalysis = 0;
+  }
+
+  if (fHistStatistics)
+  {
+    delete fHistStatistics;
+    fHistStatistics = 0;
+  }
+  
+  if (fHistBunchCrossing)
+  {
+    delete fHistBunchCrossing;
+    fHistBunchCrossing = 0;
+  }
+}
+    
+Bool_t AliPhysicsSelection::IsCollisionCandidate(const AliESDEvent* aEsd)
+{
+  // checks if the given event is a collision candidate
+  
+  if (!fTriggerAnalysis)
+    AliFatal("Not initialized!");
+    
+  const AliESDHeader* esdHeader = aEsd->GetHeader();
+  if (!esdHeader)
+  {
+    AliError("ESD Header could not be retrieved");
+    return kFALSE;
+  }
+  
+  // check event type (should be PHYSICS = 7)
+  if (esdHeader->GetEventType() != 7)
+    return kFALSE;  
+  
+  fHistStatistics->Fill(1);
+  
+  fTriggerAnalysis->FillTriggerClasses(aEsd);
+    
+  for (Int_t i=0; i < fRequTrigClasses.GetEntries(); i++)
+  {
+    const char* triggerClass = ((TObjString*) fRequTrigClasses.At(i))->String();
+    if (!aEsd->IsTriggerClassFired(triggerClass))
+    {
+      AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is not present", triggerClass));
+      return kFALSE;
+    }
+  }      
+  
+  for (Int_t i=0; i < fRejTrigClasses.GetEntries(); i++)
+  {
+    const char* triggerClass = ((TObjString*) fRejTrigClasses.At(i))->String();
+    if (aEsd->IsTriggerClassFired(triggerClass))
+    {
+      AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is present", triggerClass));
+      return kFALSE;
+    }
+  }
+  
+  fTriggerAnalysis->FillHistograms(aEsd);
+  
+  fHistStatistics->Fill(2);
+    
+  Bool_t fastOR = fTriggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kSPDGFO);
+  Bool_t v0BB = fTriggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0A) || fTriggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0C);
+  Bool_t v0BG = fTriggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0ABG) || fTriggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0CBG);
+  if (fastOR)
+    fHistStatistics->Fill(3);
+  if (v0BB)
+    fHistStatistics->Fill(4);
+  if (v0BG)
+    fHistStatistics->Fill(5);
+    
+  if (fastOR || v0BB)
+    fHistStatistics->Fill(6);
+    
+  if (!fastOR && !v0BB)
+  {
+    AliDebug(AliLog::kDebug, "Rejecting event because neither FO nor V0 has triggered");
+    return kFALSE;
+  }
+  
+  if (v0BG)
+  {
+    AliDebug(AliLog::kDebug, "Rejecting event because of V0 BG flag");
+    return kFALSE;
+  }
+      
+  fHistStatistics->Fill(7);
+  
+  // TODO additional background identification
+  
+  fHistStatistics->Fill(9);
+  
+  fHistBunchCrossing->Fill(aEsd->GetBunchCrossNumber());
+  
+  AliDebug(AliLog::kDebug, "Accepted event");
+  
+  return kTRUE;
+}
+    
+Bool_t AliPhysicsSelection::Initialize(UInt_t runNumber)
+{
+  // initializes the object for the given run
+  // TODO having the run number here and parameters hardcoded is clearly temporary, a way needs to be found to have a CDB-like configuration also for analysis
+  
+  AliInfo(Form("Initializing for run %d", runNumber));
+  
+  fRequTrigClasses.Clear();
+  fRejTrigClasses.Clear();
+  
+  fRequTrigClasses.Add(new TObjString("CINT1B-ABCE-NOPF-ALL"));
+  
+  if (!fTriggerAnalysis)
+  {
+    fTriggerAnalysis = new AliTriggerAnalysis;
+    fTriggerAnalysis->EnableHistograms();
+  }
+    
+  fTriggerAnalysis->SetSPDGFOThreshhold(1);
+  fTriggerAnalysis->SetV0TimeOffset(0);
+  
+  if (runNumber == 104321)
+    fTriggerAnalysis->SetV0TimeOffset(7.5);
+  
+  if (fHistStatistics)
+  {
+    fHistStatistics->Reset();
+  }
+  else
+  {
+    fHistStatistics = new TH1F("fHistStatistics", "fHistStatistics;;event count", 10, 0.5, 10.5);
+    
+    Int_t n = 1;
+    fHistStatistics->GetXaxis()->SetBinLabel(n++, "Total");
+    fHistStatistics->GetXaxis()->SetBinLabel(n++, "Correct trigger class(es)");
+    fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO");
+    fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0 BB");
+    fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0 BG");
+    fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO | V0 BB");
+    fHistStatistics->GetXaxis()->SetBinLabel(n++, "(FO | V0 BB) & !V0 BG");
+    fHistStatistics->GetXaxis()->SetBinLabel(n++, "Background identification");
+    fHistStatistics->GetXaxis()->SetBinLabel(n++, "Accepted");
+  }
+  
+  if (fHistBunchCrossing)
+  {
+    fHistBunchCrossing->Reset();
+  }
+  else
+    fHistBunchCrossing = new TH1F("fHistBunchCrossing", "fHistBunchCrossing;bunch crossing number;accepted events", 4000, -0.5, 3999.5);
+    
+  return kTRUE;
+}
+
+void AliPhysicsSelection::Print(Option_t* /* option */) const
+{
+  // print the configuration
+  
+  AliInfo("Initialized with:");
+  
+  TString str("Required trigger classes: ");
+  for (Int_t i=0; i < fRequTrigClasses.GetEntries(); i++)
+    str += ((TObjString*) fRequTrigClasses.At(i))->String() + " ";
+  AliInfo(str);
+  
+  str = "Rejected trigger classes: ";
+  for (Int_t i=0; i < fRejTrigClasses.GetEntries(); i++)
+    str += ((TObjString*) fRejTrigClasses.At(i))->String() + " ";
+  AliInfo(str);
+  
+  AliInfo(Form("Requiring %d FO chips (offline) or V0A or V0C and no V0 BG flag", fTriggerAnalysis->GetSPDGFOThreshhold()));
+  
+  if (fTriggerAnalysis->GetV0TimeOffset() > 0)
+    AliInfo(Form("V0 time offset active: %.2f ns", fTriggerAnalysis->GetV0TimeOffset()));
+    
+  fTriggerAnalysis->PrintTriggerClasses();
+}
+
+Long64_t AliPhysicsSelection::Merge(TCollection* list)
+{
+  // Merge a list of AliMultiplicityCorrection objects with this (needed for
+  // PROOF).
+  // Returns the number of merged objects (including this).
+
+  if (!list)
+    return 0;
+
+  if (list->IsEmpty())
+    return 1;
+
+  TIterator* iter = list->MakeIterator();
+  TObject* obj;
+
+  // collections of all histograms
+  const Int_t nHists = 9;
+  TList collections[nHists];
+
+  Int_t count = 0;
+  while ((obj = iter->Next())) {
+
+    AliPhysicsSelection* entry = dynamic_cast<AliPhysicsSelection*> (obj);
+    if (entry == 0) 
+      continue;
+
+    Int_t n = 0;
+    collections[n++].Add(entry->fTriggerAnalysis);
+    collections[n++].Add(entry->fHistStatistics);
+    collections[n++].Add(entry->fHistBunchCrossing);
+
+    count++;
+  }
+
+  Int_t n = 0;
+  fTriggerAnalysis->Merge(&collections[n++]);
+  fHistStatistics->Merge(&collections[n++]);
+  fHistBunchCrossing->Merge(&collections[n++]);
+  
+  delete iter;
+
+  return count+1;
+}
+
+void AliPhysicsSelection::SaveHistograms(const char* folder) const
+{
+  // write histograms to current directory
+  
+  if (!fHistStatistics)
+    return;
+    
+  if (folder)
+  {
+    gDirectory->mkdir(folder);
+    gDirectory->cd(folder);
+  }
+  
+  fHistStatistics->Write();
+  fHistBunchCrossing->Write();
+  
+  gDirectory->mkdir("trigger_histograms");
+  gDirectory->cd("trigger_histograms");
+  
+  fTriggerAnalysis->SaveHistograms();
+  
+  gDirectory->cd("..");
+  
+  if (folder)
+    gDirectory->cd("..");
+}
diff --git a/PWG0/AliPhysicsSelection.h b/PWG0/AliPhysicsSelection.h
new file mode 100644 (file)
index 0000000..1e3c643
--- /dev/null
@@ -0,0 +1,59 @@
+/* $Id: AliPhysicsSelection.h 35782 2009-10-22 11:54:31Z jgrosseo $ */
+
+#ifndef ALIPHYSICSSELECTION_H
+#define ALIPHYSICSSELECTION_H
+
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+//-------------------------------------------------------------------------
+//                      Implementation of   Class AliPhysicsSelection
+// 
+// This class selects collision candidates from data runs, applying selection cuts on triggers 
+// and background rejection based on the content of the ESD
+//
+//   Origin: Jan Fiete Grosse-Oetringhaus, CERN
+//-------------------------------------------------------------------------
+
+#include <TObject.h>
+#include <TList.h>
+
+class AliESDEvent;
+class TH1;
+class TH2;
+class TCollection;
+class AliTriggerAnalysis;
+
+class AliPhysicsSelection : public TObject
+{
+  public:
+    AliPhysicsSelection();
+    virtual ~AliPhysicsSelection();
+    
+    Bool_t IsCollisionCandidate(const AliESDEvent* aEsd);
+    Bool_t Initialize(UInt_t runNumber);
+    
+    virtual void Print(Option_t* option = "") const;
+    virtual Long64_t Merge(TCollection* list);
+    void SaveHistograms(const char* folder = 0) const;
+    
+    TList* GetRequiredTriggerClasses() { return &fRequTrigClasses; }
+    TList* GetRejectedTriggerClasses() { return &fRejTrigClasses; }
+    
+  protected:
+    AliTriggerAnalysis* fTriggerAnalysis; // offline trigger object
+  
+    TList fRequTrigClasses; // list of trigger class required for collision candidates
+    TList fRejTrigClasses;  // list of trigger classes not allowed for collision candidates
+    
+    TH1* fHistStatistics;      // how many events are cut away why
+    TH1* fHistBunchCrossing;   // histograms of accepted bunch crossing numbers
+    
+    ClassDef(AliPhysicsSelection, 1)
+    
+  private:
+    AliPhysicsSelection(const AliPhysicsSelection&);
+    AliPhysicsSelection& operator=(const AliPhysicsSelection&);
+};
+
+#endif
index e34dc075f9570bd3dc57245b84c9bbf29f27c46b..650e36429902fa1e608f225d104c51b6cfdf2be2 100644 (file)
@@ -45,14 +45,13 @@ ClassImp(AliTriggerAnalysis)
 
 AliTriggerAnalysis::AliTriggerAnalysis() :
   fSPDGFOThreshold(2),
-  fV0AThreshold(1),
-  fV0CThreshold(1),
+  fV0TimeOffset(0),
   fFMDLowCut(0.2),
   fFMDHitCut(0.5),
   fHistBitsSPD(0),
   fHistFiredBitsSPD(0),
   fHistV0A(0),       
-  fHistV0C(0),    
+  fHistV0C(0),
   fHistZDC(0),    
   fHistFMDA(0),    
   fHistFMDC(0),   
@@ -60,6 +59,72 @@ AliTriggerAnalysis::AliTriggerAnalysis() :
   fHistFMDSum(0),
   fTriggerClasses(0)
 {
+  // constructor
+}
+
+AliTriggerAnalysis::~AliTriggerAnalysis()
+{
+  // destructor
+  
+  if (fHistBitsSPD)
+  {
+    delete fHistBitsSPD;
+    fHistBitsSPD = 0;
+  }
+
+  if (fHistFiredBitsSPD)
+  {
+    delete fHistFiredBitsSPD;
+    fHistFiredBitsSPD = 0;
+  }
+
+  if (fHistV0A)
+  {
+    delete fHistV0A;
+    fHistV0A = 0;
+  }
+
+  if (fHistV0C)
+  {
+    delete fHistV0C;
+    fHistV0C = 0;
+  }
+
+  if (fHistZDC)
+  {
+    delete fHistZDC;
+    fHistZDC = 0;
+  }
+
+  if (fHistFMDA)
+  {
+    delete fHistFMDA;
+    fHistFMDA = 0;
+  }
+
+  if (fHistFMDC)
+  {
+    delete fHistFMDC;
+    fHistFMDC = 0;
+  }
+
+  if (fHistFMDSingle)
+  {
+    delete fHistFMDSingle;
+    fHistFMDSingle = 0;
+  }
+
+  if (fHistFMDSum)
+  {
+    delete fHistFMDSum;
+    fHistFMDSum = 0;
+  }
+
+  if (fTriggerClasses)
+  {
+    delete fTriggerClasses;
+    fTriggerClasses = 0;
+  }
 }
 
 void AliTriggerAnalysis::EnableHistograms()
@@ -79,6 +144,7 @@ void AliTriggerAnalysis::EnableHistograms()
   fHistFMDSum = new TH1F("fHistFMDSum", "FMD sum;multiplicity value;counts", 1000, 0, 10);
   
   fTriggerClasses = new TMap;
+  fTriggerClasses->SetOwner();
 }
 
 //____________________________________________________________________
@@ -99,8 +165,10 @@ const char* AliTriggerAnalysis::GetTriggerName(Trigger trigger)
     case kMB3 : str = "MB3"; break;
     case kSPDGFO : str = "SPD GFO"; break;
     case kSPDGFOBits : str = "SPD GFO Bits"; break;
-    case kV0A : str = "V0 A"; break;
-    case kV0C : str = "V0 C"; break;
+    case kV0A : str = "V0 A BB"; break;
+    case kV0C : str = "V0 C BB"; break;
+    case kV0ABG : str = "V0 A BG"; break;
+    case kV0CBG : str = "V0 C BG"; break;
     case kZDC : str = "ZDC"; break;
     case kZDCA : str = "ZDC A"; break;
     case kZDCC : str = "ZDC C"; break;
@@ -116,7 +184,7 @@ const char* AliTriggerAnalysis::GetTriggerName(Trigger trigger)
   return str;
 }
 
-Bool_t AliTriggerAnalysis::IsTriggerFired(const AliESDEvent* aEsd, Trigger trigger) const
+Bool_t AliTriggerAnalysis::IsTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
 {
   // checks if an event has been triggered
 
@@ -191,7 +259,7 @@ Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, ULong64_t
   return (trigmask & (1ull << (tclass-1)));
 }
 
-Bool_t AliTriggerAnalysis::IsOfflineTriggerFired(const AliESDEvent* aEsd, Trigger trigger) const
+Bool_t AliTriggerAnalysis::IsOfflineTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
 {
   // checks if an event has been triggered "offline"
 
@@ -246,6 +314,18 @@ Bool_t AliTriggerAnalysis::IsOfflineTriggerFired(const AliESDEvent* aEsd, Trigge
         return kTRUE;
       break;
     }
+    case kV0ABG:
+    {
+      if (V0Trigger(aEsd, kASide) == kV0BG)
+        return kTRUE;
+      break;
+    }
+    case kV0CBG:
+    {
+      if (V0Trigger(aEsd, kCSide) == kV0BG)
+        return kTRUE;
+      break;
+    }
     case kZDC:
     {
       if (ZDCTrigger(aEsd, kASide) || ZDCTrigger(aEsd, kCentralBarrel) || ZDCTrigger(aEsd, kCSide))
@@ -419,7 +499,7 @@ void AliTriggerAnalysis::FillTriggerClasses(const AliESDEvent* aEsd)
   // TODO add first and last orbit number here
 }
 
-Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHist) const
+Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHists)
 {
   // returns the number of fired chips in the SPD
   //
@@ -443,7 +523,7 @@ Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, B
       if (mult->TestFastOrFiredChips(i) == kTRUE)
       {
         nChips++;
-        if (fillHist)
+        if (fillHists)
           fHistFiredBitsSPD->Fill(i);
       }
     return nChips;
@@ -452,7 +532,7 @@ Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, B
   return -1;
 }
 
-Bool_t AliTriggerAnalysis::SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin) const
+Bool_t AliTriggerAnalysis::SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin)
 {
   // Returns if the SPD gave a global Fast OR trigger
   
@@ -463,7 +543,7 @@ Bool_t AliTriggerAnalysis::SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin)
   return kFALSE;
 }
 
-AliTriggerAnalysis::V0Decision AliTriggerAnalysis::V0Trigger(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists) const
+AliTriggerAnalysis::V0Decision AliTriggerAnalysis::V0Trigger(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists)
 {
   // Returns the V0 trigger decision in V0A | V0C
   //
@@ -495,10 +575,12 @@ AliTriggerAnalysis::V0Decision AliTriggerAnalysis::V0Trigger(const AliESDEvent*
   Float_t time = 0;
   Int_t ntime = 0;
   for (Int_t i = begin; i < end; ++i) {
-    if (esdV0->GetTime(i) > 1e-6) {
+    if (esdV0->GetTime(i) > 1e-6 && esdV0->GetAdc(i) > 6.0) {
       Float_t correctedTime = V0CorrectLeadingTime(i, esdV0->GetTime(i), esdV0->GetAdc(i));
     
       time += correctedTime;
+      time += fV0TimeOffset;
+      
       ntime++;
     }
   }
@@ -598,7 +680,7 @@ Bool_t AliTriggerAnalysis::ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) c
   return kFALSE;
 }
 
-Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHistograms) const
+Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists)
 {
   // returns number of hit combinations agove threshold
   //
@@ -628,7 +710,7 @@ Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide
           Float_t mult = fmdData->Multiplicity(det,ring,sec,strip);
           if (mult == AliESDFMD::kInvalidMult) continue;
           
-          if (fillHistograms)
+          if (fillHists)
             fHistFMDSingle->Fill(mult);
           
           if (mult > fFMDLowCut)
@@ -638,7 +720,7 @@ Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide
             if (totalMult > fFMDHitCut)
               triggers++;
               
-            if (fillHistograms)
+            if (fillHists)
               fHistFMDSum->Fill(totalMult);
               
             totalMult = 0;
@@ -651,7 +733,7 @@ Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide
   return triggers;
 }
 
-Bool_t AliTriggerAnalysis::FMDTrigger(const AliESDEvent* aEsd, AliceSide side) const
+Bool_t AliTriggerAnalysis::FMDTrigger(const AliESDEvent* aEsd, AliceSide side)
 {
   // Returns if the FMD triggered
   //
@@ -742,7 +824,7 @@ Long64_t AliTriggerAnalysis::Merge(TCollection* list)
   return count+1;
 }
 
-void AliTriggerAnalysis::WriteHistograms() const
+void AliTriggerAnalysis::SaveHistograms() const
 {
   // write histograms to current directory
   
@@ -750,8 +832,8 @@ void AliTriggerAnalysis::WriteHistograms() const
     return;
     
   fHistBitsSPD->Write();
-  fHistBitsSPD->ProjectionX()->Write();
-  fHistBitsSPD->ProjectionY()->Write();
+  fHistBitsSPD->ProjectionX();
+  fHistBitsSPD->ProjectionY();
   fHistFiredBitsSPD->Write();
   fHistV0A->Write();
   fHistV0C->Write();
index 903c216de80d56340f31a5cf53a8cb63c1e6bc0a..5d869b6e1d17630a2414521c31843d4967a15af4 100644 (file)
@@ -29,11 +29,11 @@ class AliTriggerAnalysis : public TObject
     enum V0Decision { kV0Invalid = -1, kV0Empty = 0, kV0BB, kV0BG };
     
     AliTriggerAnalysis();
-    virtual ~AliTriggerAnalysis() {}
+    virtual ~AliTriggerAnalysis();
     
     void EnableHistograms();
-
-    Bool_t IsTriggerFired(const AliESDEvent* aEsd, Trigger trigger) const;
+    
+    Bool_t IsTriggerFired(const AliESDEvent* aEsd, Trigger trigger);
     
     // using trigger bits in ESD
     Bool_t IsTriggerBitFired(const AliESDEvent* aEsd, Trigger trigger) const;
@@ -41,16 +41,16 @@ class AliTriggerAnalysis : public TObject
     Bool_t IsTriggerBitFired(const AliESDEvent* aEsd, ULong64_t tclass) const;
     
     // using ESD data from detectors
-    Bool_t IsOfflineTriggerFired(const AliESDEvent* aEsd, Trigger trigger) const;
+    Bool_t IsOfflineTriggerFired(const AliESDEvent* aEsd, Trigger trigger);
 
     // using trigger classes in ESD
     Bool_t IsTriggerClassFired(const AliESDEvent* aEsd, const Char_t* tclass) const;
     
     // some "raw" trigger functions
-    Bool_t SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin) const;
-    V0Decision V0Trigger(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists = kFALSE) const;
+    Bool_t SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin);
+    V0Decision V0Trigger(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists = kFALSE);
     Bool_t ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) const;
-    Bool_t FMDTrigger(const AliESDEvent* aEsd, AliceSide side) const;
+    Bool_t FMDTrigger(const AliESDEvent* aEsd, AliceSide side);
     
     static const char* GetTriggerName(Trigger trigger);
     
@@ -58,17 +58,16 @@ class AliTriggerAnalysis : public TObject
     void FillTriggerClasses(const AliESDEvent* aEsd);
     
     void SetSPDGFOThreshhold(Int_t t) { fSPDGFOThreshold = t; }
-    void SetV0Threshhold(Int_t aSide, Int_t cSide) { fV0AThreshold = aSide; fV0CThreshold = cSide; }
+    void SetV0TimeOffset(Float_t offset) { fV0TimeOffset = offset; }
     void SetFMDThreshold(Float_t low, Float_t hit) { fFMDLowCut = low; fFMDHitCut = hit; }
     
     Int_t GetSPDGFOThreshhold() const { return fSPDGFOThreshold; }
-    Int_t GetV0AThreshold() const { return fV0AThreshold; }
-    Int_t GetV0CThreshold() const { return fV0CThreshold; }
+    Float_t GetV0TimeOffset() const { return fV0TimeOffset; }
     Float_t GetFMDLowThreshold() const { return fFMDLowCut; }
     Float_t GetFMDHitThreshold() const { return fFMDHitCut; }
     
     virtual Long64_t Merge(TCollection* list);
-    void WriteHistograms() const;
+    void SaveHistograms() const;
     
     void PrintTriggerClasses() const;
 
@@ -78,15 +77,14 @@ class AliTriggerAnalysis : public TObject
     Bool_t IsL2InputFired(const AliESDEvent* aEsd, UInt_t input) const;
     Bool_t IsInputFired(const AliESDEvent* aEsd, Char_t level, UInt_t input) const;
     
-    Int_t SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHists = kFALSE) const;
+    Int_t SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHists = kFALSE);
     
     Float_t V0CorrectLeadingTime(Int_t i, Float_t time, Float_t adc) const;
     
-    Int_t FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHistograms) const;
+    Int_t FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists = kFALSE);
 
     Int_t fSPDGFOThreshold;         // number of chips to accept a SPD GF0 trigger
-    Int_t fV0AThreshold;            // threshold for number of BB triggers in V0A
-    Int_t fV0CThreshold;            // threshold for number of BB triggers in V0C
+    Float_t fV0TimeOffset;          // time offset applied to the times read from the V0 (in ns)
  
     Float_t fFMDLowCut;                    // 
     Float_t fFMDHitCut;                    // 
@@ -104,7 +102,7 @@ class AliTriggerAnalysis : public TObject
     TMap* fTriggerClasses;    // counts the active trigger classes (uses the full string)
     
 
-    ClassDef(AliTriggerAnalysis, 3)
+    ClassDef(AliTriggerAnalysis, 5)
     
   private:
     AliTriggerAnalysis(const AliTriggerAnalysis&);
index a8bcf4440feb58a85fda938d827db92c61a8c7ff..45856bb26aa79bf64bb2da4e354af6f1baf756e3 100644 (file)
@@ -17,6 +17,7 @@
 
 #pragma link C++ class AliPWG0Helper+;
 #pragma link C++ class AliTriggerAnalysis+;
+#pragma link C++ class AliPhysicsSelection+;
 
 #pragma link C++ class AliCorrection+;
 #pragma link C++ class AliUnfolding+;
index 62b7e0a7110a82b9177aa14a9d426230f89f3da7..a114d4cbed4c6e658001fa6da2bcf6550fba2f94 100644 (file)
@@ -12,6 +12,7 @@ SRCS = dNdEta/dNdEtaAnalysis.cxx \
       AliUnfolding.cxx \
       multiplicity/AliMultiplicityCorrection.cxx \
       AliTriggerAnalysis.cxx \
+      AliPhysicsSelection.cxx
 
 HDRS = $(SRCS:.cxx=.h)
 
index 3735a5dab9876e378d04f2338ce3a4b9f6fa801d..09250d85473ae494e9e575eb069dfa7abd25c486 100644 (file)
@@ -236,7 +236,7 @@ void AliTriggerTask::Terminate(Option_t *)
       fStats[i]->Write();
   if (fTrigger)
   {
-    fTrigger->WriteHistograms();
+    fTrigger->SaveHistograms();
     fTrigger->PrintTriggerClasses();
   }
     
index 7505cf32e14599d2909b4858af3a277721b1968a..8b31f1f50d7deab40361690d16cbe836c6914350 100644 (file)
@@ -27,10 +27,10 @@ void GetTimes(UInt_t run, UInt_t* startTime = 0, UInt_t* endTime = 0)
   gSystem->Load("libSTEER");
   
   AliCDBManager * man = AliCDBManager::Instance();
-  man->SetDefaultStorage("raw://");
-  man->SetRun(run);
+  man->SetDefaultStorage("alien://folder=/alice/data/2009/OCDB");
+  //man->SetRun(run);
   AliCDBPath cdb("GRP", "GRP", "Data");
-  obj = man->Get(cdb);
+  obj = man->Get(cdb, run);
   grp = (AliGRPObject*) obj->GetObject();
   
   if (startTime)