]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliTriggerAnalysis.cxx
a3afe3ce5fe9e2c6b99c3b56d5167fcdcc6e9c28
[u/mrichter/AliRoot.git] / ANALYSIS / AliTriggerAnalysis.cxx
1 /* $Id: AliTriggerAnalysis.cxx 35782 2009-10-22 11:54:31Z jgrosseo $ */
2
3 /**************************************************************************
4  * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
5  *                                                                        *
6  * Author: The ALICE Off-line Project.                                    *
7  * Contributors are mentioned in the code where appropriate.              *
8  *                                                                        *
9  * Permission to use, copy, modify and distribute this software and its   *
10  * documentation strictly for non-commercial purposes is hereby granted   *
11  * without fee, provided that the above copyright notice appears in all   *
12  * copies and that both the copyright notice and this permission notice   *
13  * appear in the supporting documentation. The authors make no claims     *
14  * about the suitability of this software for any purpose. It is          *
15  * provided "as is" without express or implied warranty.                  *
16  **************************************************************************/
17
18 //-------------------------------------------------------------------------
19 //                      Implementation of   Class AliTriggerAnalysis
20 //   This class provides function to check if events have been triggered based on the data in the ESD
21 //   The trigger bits, trigger class inputs and only the data (offline trigger) can be used
22 //   Origin: Jan Fiete Grosse-Oetringhaus, CERN
23 //-------------------------------------------------------------------------
24
25 #include "TH1F.h"
26 #include "TH2F.h"
27 #include "TList.h"
28 #include "TIterator.h"
29 #include "TParameter.h"
30 #include "TMap.h"
31 #include "TRandom.h"
32 #include "AliTriggerAnalysis.h"
33 #include "AliLog.h"
34 #include "AliVEvent.h"
35 #include "AliESDEvent.h"
36 #include "AliMultiplicity.h"
37 #include "AliESDVZERO.h"
38 #include "AliESDZDC.h"
39 #include "AliESDFMD.h"
40 #include "AliESDVertex.h"
41 #include "AliESDtrackCuts.h"
42 #include "AliDAQ.h"
43 #include "AliESDTrdTrack.h"
44 #include "AliVCaloTrigger.h"
45
46 ClassImp(AliTriggerAnalysis)
47
48 AliTriggerAnalysis::AliTriggerAnalysis() :
49 fSPDGFOThreshold(2),
50 fSPDGFOEfficiency(0),
51 fV0TimeOffset(0),
52 fV0AdcThr(0),
53 fV0HwAdcThr(2.5),
54 fV0HwWinLow(61.5),
55 fV0HwWinHigh(86.5),
56 fZDCCutRefSum(-568.5),
57 fZDCCutRefDelta(-2.1),
58 fZDCCutSigmaSum(3.25),
59 fZDCCutSigmaDelta(2.25),
60 fZDCCutRefSumCorr(-65.5),
61 fZDCCutRefDeltaCorr(-2.1),
62 fZDCCutSigmaSumCorr(6.0),
63 fZDCCutSigmaDeltaCorr(1.2),
64 fZDCCutZNATimeCorrMin(0.0),
65 fZDCCutZNATimeCorrMax(2.0),
66 fZDCCutZNCTimeCorrMin(0.0),
67 fZDCCutZNCTimeCorrMax(5.0),
68 fASPDCvsTCut(65),
69 fBSPDCvsTCut(4),
70 fTRDptHSE(3.),
71 fTRDpidHSE(144),
72 fTRDptHQU(2.),
73 fTRDpidHQU(164.),
74 fTRDptHEE(3.),
75 fTRDpidHEE(144),
76 fTRDminSectorHEE(6),
77 fTRDmaxSectorHEE(8),
78 fTRDptHJT(3.),
79 fTRDnHJT(3),
80 fDoFMD(kTRUE),
81 fFMDLowCut(0.2),
82 fFMDHitCut(0.5),
83 fHistBitsSPD(0),
84 fHistFiredBitsSPD(0),
85 fHistSPDClsVsTrk(0),
86 fHistV0A(0),
87 fHistV0C(0),
88 fHistZDC(0),
89 fHistTDCZDC(0),
90 fHistTimeZDC(0),
91 fHistTimeCorrZDC(0),
92 fHistFMDA(0),
93 fHistFMDC(0),
94 fHistFMDSingle(0),
95 fHistFMDSum(0),
96 fHistT0(0),
97 fTriggerClasses(0),
98 fMC(kFALSE),
99 fEsdTrackCuts(0),
100 fTPCOnly(kFALSE)
101 {
102   // constructor
103 }
104
105 //-------------------------------------------------------------------------------------------------
106 AliTriggerAnalysis::~AliTriggerAnalysis(){
107   // destructor
108   if (fHistBitsSPD)      { delete fHistBitsSPD;      fHistBitsSPD = 0;      }
109   if (fHistFiredBitsSPD) { delete fHistFiredBitsSPD; fHistFiredBitsSPD = 0; }
110   if (fHistSPDClsVsTrk)  { delete fHistSPDClsVsTrk;  fHistSPDClsVsTrk = 0;  }
111   if (fHistV0A)          { delete fHistV0A;          fHistV0A = 0;          }
112   if (fHistV0C)          { delete fHistV0C;          fHistV0C = 0;          }
113   if (fHistZDC)          { delete fHistZDC;          fHistZDC = 0;          }
114   if (fHistTDCZDC)       { delete fHistTDCZDC;       fHistTDCZDC = 0;       }
115   if (fHistTimeZDC)      { delete fHistTimeZDC;      fHistTimeZDC = 0;      }
116   if (fHistTimeCorrZDC)  { delete fHistTimeCorrZDC;  fHistTimeCorrZDC = 0;  }
117   if (fHistFMDA)         { delete fHistFMDA;         fHistFMDA = 0;         }
118   if (fHistFMDC)         { delete fHistFMDC;         fHistFMDC = 0;         }
119   if (fHistFMDSingle)    { delete fHistFMDSingle;    fHistFMDSingle = 0;    }
120   if (fHistFMDSum)       { delete fHistFMDSum;       fHistFMDSum = 0;       }
121   if (fHistT0)           { delete fHistT0;           fHistT0 = 0;           }
122   if (fEsdTrackCuts)     { delete fEsdTrackCuts;     fEsdTrackCuts =0;      }
123   if (fTriggerClasses)   { fTriggerClasses->DeleteAll(); delete fTriggerClasses; fTriggerClasses = 0; }
124 }
125
126
127 //-------------------------------------------------------------------------------------------------
128 void AliTriggerAnalysis::EnableHistograms(Bool_t isLowFlux){
129   // creates the monitoring histograms 
130   // dynamical range of histograms can be adapted for pp and pPb via isLowFlux flag)
131   // TODO check limits for FMD
132   
133   Int_t nBins  = isLowFlux ?  600 : 1202;
134   Int_t nBinsX = isLowFlux ?  100 :  300;
135   Int_t nBinsY = isLowFlux ?  500 : 1000;
136   Float_t xMax = isLowFlux ?  400 : 2999.5;
137   Float_t yMax = isLowFlux ? 4000 : 9999.5;
138   
139   // do not add these hists to the directory
140   Bool_t oldStatus = TH1::AddDirectoryStatus();
141   TH1::AddDirectory(kFALSE); 
142   fHistBitsSPD      = new TH2F("fHistBitsSPD", "SPD GFO;number of fired chips (offline);number of fired chips (hardware)", nBins, -1.5, -1.5 + nBins, nBins, -1.5, -1.5+nBins);
143   fHistFiredBitsSPD = new TH1F("fHistFiredBitsSPD", "SPD GFO Hardware;chip number;events", 1200, -0.5, 1199.5);
144   fHistSPDClsVsTrk  = new TH2F("fHistSPDClsVsTrk", "SPD Clusters vs Tracklets; n tracklets; n clusters", nBinsX, -0.5, xMax, nBinsY, -0.5, yMax);
145   fHistV0A          = new TH1F("fHistV0A", "V0A;leading time (ns);events", 400, -100, 100);
146   fHistV0C          = new TH1F("fHistV0C", "V0C;leading time (ns);events", 400, -100, 100);
147   fHistZDC          = new TH1F("fHistZDC", "ZDC;trigger bits;events", 8, -1.5, 6.5);
148   fHistTDCZDC       = new TH1F("fHistTDCZDC", "ZDC;TDC bits;events", 32, -0.5, 32-0.5);
149   fHistTimeZDC      = new TH2F("fHistTimeZDC", "ZDC;TDC timing C-A;TDC timing C+A", 120,-30,30,120,-600,-540);
150   fHistTimeCorrZDC  = new TH2F("fHistTimeCorrZDC", "ZDC;Corrected TDC timing C-A; Corrected TDC timing C+A", 120,-30,30,260,-100,30);
151   fHistFMDA         = new TH1F("fHistFMDA", "FMDA;combinations above threshold;events", 102, -1.5, 100.5);
152   fHistFMDC         = new TH1F("fHistFMDC", "FMDC;combinations above threshold;events", 102, -1.5, 100.5);
153   fHistFMDSingle    = new TH1F("fHistFMDSingle", "FMD single;multiplicity value;counts", 1000, 0, 10);
154   fHistFMDSum       = new TH1F("fHistFMDSum", "FMD sum;multiplicity value;counts", 1000, 0, 10);
155   fHistT0           = new TH1F("fHistT0", "T0;time (ns);events", 100, -25, 25);
156   fTriggerClasses = new TMap;
157   fTriggerClasses->SetOwner();
158   TH1::AddDirectory(oldStatus);
159 }
160
161
162 //-------------------------------------------------------------------------------------------------
163 const char* AliTriggerAnalysis::GetTriggerName(Trigger trigger){
164   // returns the name of the requested trigger
165   // the returned string will only be valid until the next call to this function [not thread-safe]
166   
167   static TString str;
168   
169   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
170   
171   switch (triggerNoFlags)  {
172     case kAcceptAll :      str = "ACCEPT ALL (bypass!)";      break;
173     case kMB1 :            str = "MB1";                       break;
174     case kMB2 :            str = "MB2";                       break;
175     case kMB3 :            str = "MB3";                       break;
176     case kSPDGFO :         str = "SPD GFO";                   break;
177     case kSPDGFOBits :     str = "SPD GFO Bits";              break;
178     case kSPDGFOL0 :       str = "SPD GFO L0 (first layer)";  break;
179     case kSPDGFOL1 :       str = "SPD GFO L1 (second layer)"; break;
180     case kSPDClsVsTrkBG :  str = "Cluster vs Tracklets";      break;
181     case kV0A :            str = "V0 A BB";                   break;
182     case kV0C :            str = "V0 C BB";                   break;
183     case kV0OR :           str = "V0 OR BB";                  break;
184     case kV0AND :          str = "V0 AND BB";                 break;
185     case kV0ABG :          str = "V0 A BG";                   break;
186     case kV0CBG :          str = "V0 C BG";                   break;
187     case kZDC :            str = "ZDC";                       break;
188     case kZDCA :           str = "ZDC A";                     break;
189     case kZDCC :           str = "ZDC C";                     break;
190     case kZNA :            str = "ZN A";                      break;
191     case kZNC :            str = "ZN C";                      break;
192     case kZNABG :          str = "ZN A BG";                   break;
193     case kZNCBG :          str = "ZN C BG";                   break;
194     case kFMDA :           str = "FMD A";                     break;
195     case kFMDC :           str = "FMD C";                     break;
196     case kFPANY :          str = "SPD GFO | V0 | ZDC | FMD";  break;
197     case kNSD1 :           str = "NSD1";                      break;
198     case kMB1Prime:        str = "MB1prime";                  break;
199     case kZDCTDCA :        str = "ZDC TDC A";                 break;
200     case kZDCTDCC :        str = "ZDC TDC C";                 break;
201     case kZDCTime :        str = "ZDC Time Cut";              break;
202     case kCentral :        str = "V0 Central";                break;
203     case kSemiCentral :    str = "V0 Semi-central";           break;
204     case kEmcalL0 :        str = "EMCAL";                     break;
205     case kTRDHCO :         str = "TRD cosmics";               break;
206     case kTRDHJT :         str = "TRD Jet";                   break;
207     case kTRDHSE :         str = "TRD Single Electron";       break;
208     case kTRDHQU :         str = "TRD Quarkonia";             break;
209     case kTRDHEE :         str = "TRD Dielectron";            break;
210     default:               str = "";                          break;
211   }
212   if (trigger & kOfflineFlag) str += " OFFLINE";  
213   if (trigger & kOneParticle) str += " OneParticle";  
214   if (trigger & kOneTrack)    str += " OneTrack";  
215   return str;
216 }
217
218
219 //-------------------------------------------------------------------------------------------------
220 Bool_t AliTriggerAnalysis::IsTriggerFired(const AliESDEvent* aEsd, Trigger trigger){
221   // checks if an event has been triggered
222   if (trigger & kOfflineFlag) return IsOfflineTriggerFired(aEsd, trigger);
223   return IsTriggerBitFired(aEsd, trigger);
224 }
225
226
227 //-------------------------------------------------------------------------------------------------
228 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* /*aEsd*/, Trigger /*trigger*/) const { 
229   AliFatal("This IsTriggerBitFired function is obsolete.\n"); 
230   return 0; 
231 }
232
233
234 //-------------------------------------------------------------------------------------------------
235 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, ULong64_t tclass) const {
236   // Checks if corresponding bit in mask is on
237   ULong64_t trigmask = aEsd->GetTriggerMask();
238   return (trigmask & (1ull << (tclass-1)));
239 }
240
241
242 //-------------------------------------------------------------------------------------------------
243 Int_t AliTriggerAnalysis::EvaluateTrigger(const AliESDEvent* aEsd, Trigger trigger){
244   // evaluates a given trigger
245   // trigger combinations are not supported, for that see IsOfflineTriggerFired
246
247   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
248   Bool_t offline = trigger & kOfflineFlag;
249   
250   if (!offline) {
251     if ( triggerNoFlags==kT0BG
252       || triggerNoFlags==kT0Pileup
253       || triggerNoFlags==kSPDClsVsTrkBG
254       || triggerNoFlags==kZDCA
255       || triggerNoFlags==kZDCC
256       || triggerNoFlags==kZDCTDCA
257       || triggerNoFlags==kZDCTDCC
258       || triggerNoFlags==kZDCTime
259       || triggerNoFlags==kZNA
260       || triggerNoFlags==kZNC
261       || triggerNoFlags==kZNABG
262       || triggerNoFlags==kZNCBG
263       || triggerNoFlags==kFMDA
264       || triggerNoFlags==kFMDC
265       || triggerNoFlags==kTPCLaserWarmUp
266       || triggerNoFlags==kTPCHVdip
267       || triggerNoFlags==kIncompleteEvent
268       || triggerNoFlags==kEMCAL
269       || triggerNoFlags==kEmcalL0
270       || triggerNoFlags==kEmcalL1GammaHigh
271       || triggerNoFlags==kEmcalL1GammaLow
272       || triggerNoFlags==kEmcalL1JetHigh
273       || triggerNoFlags==kEmcalL1JetLow
274       || triggerNoFlags==kTRDHCO
275       || triggerNoFlags==kTRDHJT
276       || triggerNoFlags==kTRDHSE
277       || triggerNoFlags==kTRDHQU
278       || triggerNoFlags==kTRDHEE
279       ) AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
280   } else {
281     if (  triggerNoFlags==kCTPV0A 
282         ||triggerNoFlags==kCTPV0C
283         ||triggerNoFlags==kCentral
284         ||triggerNoFlags==kSemiCentral
285       ) AliFatal(Form("Offline trigger not available for trigger %d", triggerNoFlags));
286   }
287   
288   switch (triggerNoFlags) {
289     case kCTPV0A:          return aEsd->GetHeader()->IsTriggerInputFired("V0A");
290     case kCTPV0C:          return aEsd->GetHeader()->IsTriggerInputFired("V0C");
291     case kSPDGFO:          return SPDFiredChips(aEsd, !offline, kFALSE, 0); 
292     case kSPDGFOL0:        return SPDFiredChips(aEsd, !offline, kFALSE, 1);
293     case kSPDGFOL1:        return SPDFiredChips(aEsd, !offline, kFALSE, 2);
294     case kSPDClsVsTrkBG:   return IsSPDClusterVsTrackletBG(aEsd);
295     case kV0A:             return V0Trigger(aEsd, kASide, !offline) == kV0BB; 
296     case kV0C:             return V0Trigger(aEsd, kCSide, !offline) == kV0BB;
297     case kV0ABG:           return V0Trigger(aEsd, kASide, !offline) == kV0BG;
298     case kV0CBG:           return V0Trigger(aEsd, kCSide, !offline) == kV0BG;
299     case kT0:              return T0Trigger(aEsd, !offline) == kT0BB;
300     case kT0BG:            return T0Trigger(aEsd, !offline) == kT0DecBG;
301     case kT0Pileup:        return T0Trigger(aEsd, !offline) == kT0DecPileup;
302     case kZDCA:            return ZDCTrigger(aEsd, kASide);
303     case kZDCC:            return ZDCTrigger(aEsd, kCSide);
304     case kZDCTDCA:         return ZDCTDCTrigger(aEsd, kASide);
305     case kZDCTDCC:         return ZDCTDCTrigger(aEsd, kCSide);
306     case kZDCTime:         return ZDCTimeTrigger(aEsd);
307     case kZNA:             return ZDCTDCTrigger(aEsd,kASide,kTRUE,kFALSE,kFALSE);
308     case kZNC:             return ZDCTDCTrigger(aEsd,kCSide,kTRUE,kFALSE,kFALSE);
309     case kZNABG:           return ZDCTimeBGTrigger(aEsd,kASide);
310     case kZNCBG:           return ZDCTimeBGTrigger(aEsd,kCSide);
311     case kFMDA:            return FMDTrigger(aEsd, kASide);
312     case kFMDC:            return FMDTrigger(aEsd, kCSide);
313     case kTPCLaserWarmUp:  return IsLaserWarmUpTPCEvent(aEsd);
314     case kTPCHVdip:        return IsHVdipTPCEvent(aEsd);
315     case kIncompleteEvent: return IsIncompleteEvent(aEsd);
316     case kEMCAL:           return EMCALCellsTrigger(aEsd);
317     case kEmcalL0:         return EMCALTrigger(aEsd,kEmcalL0);
318     case kEmcalL1GammaHigh:return EMCALTrigger(aEsd,kEmcalL1GammaHigh);
319     case kEmcalL1GammaLow: return EMCALTrigger(aEsd,kEmcalL1GammaLow);
320     case kEmcalL1JetHigh:  return EMCALTrigger(aEsd,kEmcalL1JetHigh);
321     case kEmcalL1JetLow:   return EMCALTrigger(aEsd,kEmcalL1JetLow);
322     case kTRDHCO:          return TRDTrigger(aEsd,kTRDHCO);
323     case kTRDHJT:          return TRDTrigger(aEsd,kTRDHJT);
324     case kTRDHSE:          return TRDTrigger(aEsd,kTRDHSE);
325     case kTRDHQU:          return TRDTrigger(aEsd,kTRDHQU);
326     case kTRDHEE:          return TRDTrigger(aEsd,kTRDHEE);
327     case kCentral: {
328       if (!aEsd->GetVZEROData()) { AliWarning("V0 centrality trigger bits were not filled!"); return kFALSE; }
329       if (!aEsd->GetVZEROData()->TestBit(AliESDVZERO::kTriggerChargeBitsFilled)) return kFALSE;
330       return aEsd->GetVZEROData()->GetTriggerBits() & (1<<AliESDVZERO::kCTA2andCTC2);
331     }
332     case kSemiCentral: {
333       if (!aEsd->GetVZEROData()) { AliWarning("V0 centrality trigger bits were not filled!"); return kFALSE; }
334       if (!aEsd->GetVZEROData()->TestBit(AliESDVZERO::kTriggerChargeBitsFilled)) return kFALSE;
335       return aEsd->GetVZEROData()->GetTriggerBits() & (1<<AliESDVZERO::kCTA1andCTC1);
336     }
337     default: AliFatal(Form("Trigger type %d not implemented", triggerNoFlags));
338   }
339   return kFALSE;
340 }
341
342
343 //-------------------------------------------------------------------------------------------------
344 Bool_t AliTriggerAnalysis::IsOfflineTriggerFired(const AliESDEvent* aEsd, Trigger trigger){
345   // checks if an event has been triggered "offline"
346   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
347   
348   Bool_t decision = kFALSE;
349   switch (triggerNoFlags) {
350     case kAcceptAll:        decision = kTRUE; break;
351     case kMB1:              decision = SPDGFOTrigger(aEsd, 0) ||  V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB;  break;
352     case kMB2:              decision = SPDGFOTrigger(aEsd, 0) && (V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB); break;
353     case kMB3:              decision = SPDGFOTrigger(aEsd, 0) &&  V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB;  break;
354     case kSPDGFO:           decision = SPDGFOTrigger(aEsd, 0); break;
355     case kSPDGFOBits:       decision = SPDGFOTrigger(aEsd, 1); break;
356     case kV0A:              decision = V0Trigger(aEsd, kASide, kFALSE) == kV0BB; break;
357     case kV0C:              decision = V0Trigger(aEsd, kCSide, kFALSE) == kV0BB; break;
358     case kV0OR:             decision = V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB; break;
359     case kV0AND:            decision = V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB; break;
360     case kV0ABG:            decision = V0Trigger(aEsd, kASide, kFALSE) == kV0BG; break;
361     case kV0CBG:            decision = V0Trigger(aEsd, kCSide, kFALSE) == kV0BG; break;
362     case kZDC:              decision = ZDCTrigger(aEsd, kASide) || ZDCTrigger(aEsd, kCentralBarrel) || ZDCTrigger(aEsd, kCSide); break;
363     case kZDCA:             decision = ZDCTrigger(aEsd, kASide); break;
364     case kZDCC:             decision = ZDCTrigger(aEsd, kCSide); break;
365     case kZNA:              decision = ZDCTDCTrigger(aEsd,kASide,kTRUE,kFALSE,kFALSE); break;
366     case kZNC:              decision = ZDCTDCTrigger(aEsd,kCSide,kTRUE,kFALSE,kFALSE); break;
367     case kZNABG:            decision = ZDCTimeBGTrigger(aEsd,kASide); break;
368     case kZNCBG:            decision = ZDCTimeBGTrigger(aEsd,kCSide); break;
369     case kFMDA:             decision = FMDTrigger(aEsd, kASide); break;
370     case kFMDC:             decision = FMDTrigger(aEsd, kCSide); break;
371     case kEMCAL:            decision = EMCALCellsTrigger(aEsd); break;
372     case kEmcalL0:          decision = EMCALTrigger(aEsd,kEmcalL0);          break;
373     case kEmcalL1GammaHigh: decision = EMCALTrigger(aEsd,kEmcalL1GammaHigh); break;
374     case kEmcalL1GammaLow:  decision = EMCALTrigger(aEsd,kEmcalL1GammaLow);  break;
375     case kEmcalL1JetHigh:   decision = EMCALTrigger(aEsd,kEmcalL1JetHigh);   break;
376     case kEmcalL1JetLow:    decision = EMCALTrigger(aEsd,kEmcalL1JetLow);    break;
377     case kTRDHCO:           decision = TRDTrigger(aEsd,kTRDHCO); break;
378     case kTRDHJT:           decision = TRDTrigger(aEsd,kTRDHJT); break;
379     case kTRDHSE:           decision = TRDTrigger(aEsd,kTRDHSE); break;
380     case kTRDHQU:           decision = TRDTrigger(aEsd,kTRDHQU); break;
381     case kTRDHEE:           decision = TRDTrigger(aEsd,kTRDHEE); break;
382     case kNSD1:             decision = SPDFiredChips(aEsd, 0) >= 5 || (V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB); break;
383     case kFPANY:            decision |= SPDGFOTrigger(aEsd, 0); 
384                             decision |= V0Trigger(aEsd, kASide, kFALSE) == kV0BB;
385                             decision |= V0Trigger(aEsd, kCSide, kFALSE) == kV0BB;
386                             decision |= ZDCTrigger(aEsd, kASide);
387                             decision |= ZDCTrigger(aEsd, kCentralBarrel);
388                             decision |= ZDCTrigger(aEsd, kCSide);
389                             decision |= FMDTrigger(aEsd, kASide);
390                             decision |= FMDTrigger(aEsd, kCSide);
391                             break; 
392     case kMB1Prime:         decision |= SPDGFOTrigger(aEsd, 0) && V0Trigger(aEsd, kASide, kFALSE) == kV0BB;
393                             decision |= SPDGFOTrigger(aEsd, 0) && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB;
394                             decision |= V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB;
395                             break; 
396     default:                AliFatal(Form("Trigger type %d not implemented", triggerNoFlags));
397   }
398   
399   // hadron-level requirement: SPD tracklets
400   if (decision && (trigger & kOneParticle))  {
401     decision = kFALSE;
402     
403     const AliESDVertex* vertex = aEsd->GetPrimaryVertexSPD();
404     const AliMultiplicity* mult = aEsd->GetMultiplicity();
405     
406     if (mult && vertex && vertex->GetNContributors() > 0 && (!vertex->IsFromVertexerZ() || vertex->GetDispersion() < 0.02) && TMath::Abs(vertex->GetZv()) < 5.5)    {
407       for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i) {
408         if (TMath::Abs(mult->GetEta(i)) < 1) {
409           decision = kTRUE;
410           break;
411         }
412       }
413     }
414   }
415   
416   // hadron level requirement: TPC tracks
417   if (decision && (trigger & kOneTrack)) {
418     decision = kFALSE;
419     const AliESDVertex* vertex =0x0;
420     vertex = aEsd->GetPrimaryVertexTracks();
421     if (!vertex || vertex->GetNContributors() <= 0) vertex = aEsd->GetPrimaryVertexSPD();
422     Float_t ptmin, ptmax;
423     fEsdTrackCuts->GetPtRange(ptmin,ptmax);
424     AliDebug(3, Form("ptmin = %f, ptmax = %f\n",ptmin, ptmax));
425     
426     if (vertex && vertex->GetNContributors() > 0 && (!vertex->IsFromVertexerZ() || vertex->GetDispersion() < 0.02) && TMath::Abs(vertex->GetZv()) < 10.) {
427       AliDebug(3,Form("Check on the vertex passed\n"));
428       for (Int_t i=0; i<aEsd->GetNumberOfTracks(); ++i){
429         if (fTPCOnly == kFALSE){
430           if (fEsdTrackCuts->AcceptTrack(aEsd->GetTrack(i))){
431             AliDebug(2, Form("pt of track = %f --> check passed\n",aEsd->GetTrack(i)->Pt()));
432             decision = kTRUE;
433             break;
434           }
435         }
436         else {
437           // TPC only tracks
438           AliESDtrack *tpcTrack = fEsdTrackCuts->GetTPCOnlyTrack((AliESDEvent*)aEsd, i);
439           if (!tpcTrack){
440             AliDebug(3,Form("track %d is NOT a TPC track",i));
441             continue;
442           }
443           else{
444             AliDebug(3,Form("track %d IS a TPC track",i));
445             if (!(fEsdTrackCuts->AcceptTrack(tpcTrack))) {
446               AliDebug(2, Form("TPC track %d NOT ACCEPTED, pt = %f, eta = %f",i,tpcTrack->Pt(), tpcTrack->Eta()));
447               delete tpcTrack; tpcTrack = 0x0;
448               continue;
449             } // end if the TPC track is not accepted
450             else{
451               AliDebug(2, Form("TPC track %d ACCEPTED, pt = %f, eta = %f",i,tpcTrack->Pt(), tpcTrack->Eta()));
452               decision = kTRUE;
453               delete tpcTrack; tpcTrack = 0x0;
454               break;
455             } // end if the TPC track is accepted
456           } // end if it is a TPC track
457         } // end if you are looking at TPC only tracks                  
458       } // end loop on tracks
459     } // end check on vertex
460     else{
461       AliDebug(4,Form("Check on the vertex not passed\n"));
462       for (Int_t i=0; i<aEsd->GetNumberOfTracks(); ++i){
463         if (fEsdTrackCuts->AcceptTrack(aEsd->GetTrack(i))){
464           AliDebug(4,Form("pt of track = %f --> check would be passed if the vertex was ok\n",aEsd->GetTrack(i)->Pt()));
465           break;
466         }
467       }
468     }
469     if (!decision) AliDebug(3,("Check for kOneTrack NOT passed\n"));
470   }
471   
472   return decision;
473 }
474
475
476
477 //-------------------------------------------------------------------------------------------------
478 Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHists, Int_t layer){
479   // returns the number of fired chips in the SPD
480   //
481   // origin = 0 --> aEsd->GetMultiplicity()->GetNumberOfFiredChips() (filled from clusters)
482   // origin = 1 --> aEsd->GetMultiplicity()->TestFastOrFiredChips() (from hardware bits)
483   // layer  = 0 --> both layers
484   // layer  = 1 --> inner
485   // layer  = 2 --> outer
486   
487   const AliMultiplicity* mult = aEsd->GetMultiplicity();
488   if (!mult) { AliFatal("AliMultiplicity not available"); }
489   
490   if (origin == 0) {
491     if (layer == 0) return mult->GetNumberOfFiredChips(0) + mult->GetNumberOfFiredChips(1);
492     return mult->GetNumberOfFiredChips(layer-1); 
493   }
494   
495   if (origin == 1) {
496     Int_t nChips = 0;
497     Int_t firstChip = 0;
498     Int_t lastChip  = 1200;
499     if(layer == 1) lastChip  = 400;
500     if(layer == 2) firstChip = 400;
501
502     for (Int_t i=firstChip; i<lastChip; i++) {
503       if (mult->TestFastOrFiredChips(i)) {
504         // efficiency simulation (if enabled)
505         if (fSPDGFOEfficiency) if (gRandom->Uniform() > fSPDGFOEfficiency->GetBinContent(i+1)) continue;
506         nChips++;
507         if (fillHists) fHistFiredBitsSPD->Fill(i);
508       }
509     }
510     return nChips;
511   }
512   
513   return -1;
514 }
515
516
517 //-------------------------------------------------------------------------------------------------
518 Int_t AliTriggerAnalysis::SSDClusters(const AliVEvent* event){ 
519   return event->GetNumberOfITSClusters(4)+event->GetNumberOfITSClusters(5); 
520 }
521
522
523 //-------------------------------------------------------------------------------------------------
524 AliTriggerAnalysis::V0Decision AliTriggerAnalysis::V0Trigger(const AliESDEvent* aEsd, AliceSide side, Bool_t online, Bool_t fillHists){
525   // Returns the V0 trigger decision in V0A | V0C
526   //
527   // Returns kV0Fake if the calculated average time is in a window where neither BB nor BG is expected. 
528   // The rate of such triggers can be used to estimate the background. Note that the rate has to be 
529   // rescaled with the size of the windows (numerical values see below in the code)
530   //
531   // argument 'online' is used as a switch between online and offline trigger algorithms
532   //
533   // Based on an algorithm by Cvetan Cheshkov
534   
535   AliESDVZERO* esdV0 = aEsd->GetVZEROData();
536   if (!esdV0) { AliError("AliESDVZERO not available");  return kV0Invalid; }
537   if (side != kASide && side != kCSide) return kV0Invalid;
538   
539   AliDebug(2,Form("In V0Trigger: %f %f",esdV0->GetV0ATime(),esdV0->GetV0CTime()));
540   
541   Int_t begin = (side == kASide) ? 32 :  0;
542   Int_t end   = (side == kASide) ? 64 : 32;
543   
544   if (esdV0->TestBit(AliESDVZERO::kDecisionFilled)) {
545     if (online) {
546       if (esdV0->TestBit(AliESDVZERO::kOnlineBitsFilled)) {
547         for (Int_t i = begin; i < end; ++i) if (esdV0->GetBBFlag(i)) return kV0BB;
548         for (Int_t i = begin; i < end; ++i) if (esdV0->GetBGFlag(i)) return kV0BG;
549         return kV0Empty;
550       }
551       else {
552         AliWarning("V0 online trigger analysis is not yet available!");
553         return kV0BB;
554       }
555     }
556     else {
557       if (fillHists) {
558         if (side == kASide && fHistV0A) fHistV0A->Fill(esdV0->GetV0ATime());
559         if (side == kCSide && fHistV0C) fHistV0C->Fill(esdV0->GetV0CTime());
560       }
561       if      (side == kASide) return (V0Decision) esdV0->GetV0ADecision();
562       else if (side == kCSide) return (V0Decision) esdV0->GetV0CDecision();
563     }
564   }
565   
566   Float_t time = 0;
567   Float_t weight = 0;
568   if (fMC) {
569     Int_t runRange;
570     if      (aEsd->GetRunNumber() <= 104803) runRange = 0;
571     else if (aEsd->GetRunNumber() <= 104876) runRange = 1;
572     else runRange = 2;
573     
574     Float_t factors[3][64] = {
575       /*104792-104803*/ {4.6,5.9,6.3,6.0,4.7,5.9,4.9,5.4,4.8,4.1,4.9,4.6,4.5,5.5,5.1,5.8,4.3,4.0,4.0,3.3,3.1,2.9,3.0,5.6,3.3,4.9,3.9,5.3,4.1,4.4,3.9,5.5,5.7,9.5,5.1,5.3,6.6,7.1,8.9,4.4,4.1,5.9,9.0,4.5,4.1,6.0,4.7,7.1,4.2,4.7,3.9,6.3,5.9,4.8,4.7,4.5,4.7,5.4,5.8,5.0,5.1,5.9,5.3,3.6},
576       /*104841-104876*/ {4.6,4.8,4.9,4.8,4.3,4.9,4.4,4.5,4.6,5.0,4.7,4.6,4.7,4.6,4.6,5.5,4.7,4.5,4.7,5.0,6.5,7.6,5.3,4.9,5.5,4.8,4.6,4.9,4.5,4.5,4.6,4.9,5.7,9.8,4.9,5.2,7.1,7.1,8.1,4.4,4.0,6.0,8.3,4.6,4.2,5.6,4.6,6.4,4.4,4.7,4.5,6.5,6.0,4.7,4.5,4.4,4.8,5.5,5.9,5.3,5.0,5.7,5.1,3.6},
577       /*104890-104892*/ {4.7,5.2,4.8,5.0,4.4,5.0,4.4,4.6,4.6,4.5,4.4,4.6,4.5,4.6,4.8,5.5,4.8,4.5,4.4,4.3,5.4,7.7,5.6,5.0,5.4,4.3,4.5,4.8,4.5,4.5,4.6,5.3,5.7,9.6,4.9,5.4,6.1,7.2,8.6,4.4,4.0,5.4,8.8,4.4,4.2,5.8,4.7,6.7,4.3,4.7,4.0,6.1,6.0,4.9,4.8,4.6,4.7,5.2,5.7,5.0,5.0,5.8,5.3,3.6}
578     };
579     Float_t dA = 77.4 - 11.0;
580     Float_t dC = 77.4 - 2.9;
581     // Time misalignment
582     Float_t timeShift[64] = {0.477957 , 0.0889999 , 0.757669 , 0.205439 , 0.239666 , -0.183705 , 0.442873 , -0.281366 , 0.260976 , 0.788995 , 0.974758 , 0.548532 , 0.495023 , 0.868472 , 0.661167 , 0.358307 , 0.221243 , 0.530179 , 1.26696 , 1.33082 , 1.27086 , 1.77133 , 1.10253 , 0.634806 , 2.14838 , 1.50212 , 1.59253 , 1.66122 , 1.16957 , 1.52056 , 1.47791 , 1.81905 , -1.94123 , -1.29124 , -2.16045 , -1.78939 , -3.11111 , -1.87178 , -1.57671 , -1.70311 , -1.81208 , -1.94475 , -2.53058 , -1.7042 , -2.08109 , -1.84416 , -0.61073 , -1.77145 , 0.16999 , -0.0585339 , 0.00401133 , 0.397726 , 0.851111 , 0.264187 , 0.59573 , -0.158263 , 0.584362 , 1.20835 , 0.927573 , 1.13895 , 0.64648 , 2.18747 , 1.68909 , 0.451194};
583     Float_t dA2 = 2.8, dC2 = 3.3;
584     
585     if (online) {
586       for (Int_t i = begin; i < end; ++i) {
587         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
588         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
589         if (esdV0->GetTime(i) >= 1e-6 && tempTime > fV0HwWinLow && tempTime < fV0HwWinHigh && tempAdc > fV0HwAdcThr) return kV0BB;
590       }
591       return kV0Empty;
592     }
593     else {
594       for (Int_t i = begin; i < end; ++i) {
595         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
596         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA : esdV0->GetTime(i)+dC;
597         Float_t tempRawTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
598         if (esdV0->GetTime(i) >= 1e-6 && tempRawTime < 125.0 && tempAdc > fV0AdcThr) {
599           weight += 1.0;
600           time += tempTime;
601         }
602       }
603     }
604   }
605   else {
606     if (online) {
607       for (Int_t i = begin; i < end; ++i) {
608         if (esdV0->GetTime(i) >= 1e-6 && esdV0->GetTime(i) > fV0HwWinLow && esdV0->GetTime(i) < fV0HwWinHigh && esdV0->GetAdc(i) > fV0HwAdcThr) return kV0BB;
609       }
610       return kV0Empty;
611     }
612     else {
613       for (Int_t i = begin; i < end; ++i) {
614         if (esdV0->GetTime(i) > 1e-6 && esdV0->GetAdc(i) > fV0AdcThr) {
615           Float_t correctedTime = V0CorrectLeadingTime(i, esdV0->GetTime(i), esdV0->GetAdc(i),aEsd->GetRunNumber());
616           Float_t timeWeight = V0LeadingTimeWeight(esdV0->GetAdc(i));
617           time += correctedTime*timeWeight;
618           weight += timeWeight;
619         }
620       }
621     }
622   }
623   
624   if (weight > 0) time /= weight;
625   time += fV0TimeOffset;
626   
627   if (fillHists) {
628     if (side == kASide && fHistV0A) fHistV0A->Fill(time);
629     if (side == kCSide && fHistV0C) fHistV0C->Fill(time);
630   }
631   
632   if (side == kASide) {
633     if (time > 68 && time < 100)  return kV0BB;
634     if (time > 54 && time < 57.5) return kV0BG;
635     if (time > 57.5 && time < 68) return kV0Fake;
636   }
637   if (side == kCSide) {
638     if (time > 75.5 && time < 100) return kV0BB;
639     if (time > 69.5 && time < 73)  return kV0BG; 
640     if (time > 55 && time < 69.5)  return kV0Fake;
641   }
642   
643   return kV0Empty;
644 }
645
646
647 //-------------------------------------------------------------------------------------------------
648 Float_t AliTriggerAnalysis::V0CorrectLeadingTime(Int_t i, Float_t time, Float_t adc, Int_t runNumber) const {
649   // Correct for slewing and align the channels
650   // Authors: Cvetan Cheshkov / Raphael Tieulent
651   if (time == 0) return 0;
652
653   // Time alignment
654   Float_t timeShift[64] = {0.477957 , 0.0889999 , 0.757669 , 0.205439 , 0.239666 , -0.183705 , 0.442873 , -0.281366 , 0.260976 , 0.788995 , 0.974758 , 0.548532 , 0.495023 , 0.868472 , 0.661167 , 0.358307 , 0.221243 , 0.530179 , 1.26696 , 1.33082 , 1.27086 , 1.77133 , 1.10253 , 0.634806 , 2.14838 , 1.50212 , 1.59253 , 1.66122 , 1.16957 , 1.52056 , 1.47791 , 1.81905 , -1.94123 , -1.29124 , -2.16045 , -1.78939 , -3.11111 , -1.87178 , -1.57671 , -1.70311 , -1.81208 , -1.94475 , -2.53058 , -1.7042 , -2.08109 , -1.84416 , -0.61073 , -1.77145 , 0.16999 , -0.0585339 , 0.00401133 , 0.397726 , 0.851111 , 0.264187 , 0.59573 , -0.158263 , 0.584362 , 1.20835 , 0.927573 , 1.13895 , 0.64648 , 2.18747 , 1.68909 , 0.451194};
655   
656   if(runNumber < 106031) time -= timeShift[i];
657   
658   // Slewing correction
659   if (adc == 0) return time;
660   
661   Double_t p1 = 1.57345e+1;
662   Double_t p2 =-4.25603e-1;
663   
664   if(runNumber >= 106031) adc *= (2.5/4.0);
665   return (time - p1*TMath::Power(adc,p2));
666 }
667
668
669 //-------------------------------------------------------------------------------------------------
670 Float_t AliTriggerAnalysis::V0LeadingTimeWeight(Float_t adc) const {
671   if (adc < 1e-6) return 0;
672   
673   Float_t p1 = 40.211;
674   Float_t p2 =-4.25603e-1;
675   Float_t p3 = 0.5646;
676   
677   return 1./(p1*p1*TMath::Power(adc,2.*(p2-1.))+p3*p3);
678 }
679
680
681 //-------------------------------------------------------------------------------------------------
682 Bool_t AliTriggerAnalysis::ZDCTDCTrigger(const AliESDEvent* aEsd, AliceSide side, Bool_t useZN, Bool_t useZP, Bool_t fillHists) const{
683   // Returns if ZDC triggered, based on TDC information 
684   
685   AliESDZDC *esdZDC = aEsd->GetESDZDC();
686   
687   Bool_t zdcNA = kFALSE;
688   Bool_t zdcNC = kFALSE;
689   Bool_t zdcPA = kFALSE;
690   Bool_t zdcPC = kFALSE;
691   
692   if (fMC) { // If it's MC, we use the energy
693     Double_t minEnergy = 0;
694     zdcNA = esdZDC->GetZDCN2Energy()>minEnergy;
695     zdcNC = esdZDC->GetZDCN1Energy()>minEnergy;
696     zdcPA = esdZDC->GetZDCP2Energy()>minEnergy;
697     zdcPC = esdZDC->GetZDCP1Energy()>minEnergy;
698   }
699   else {
700     Bool_t tdc[32] = {kFALSE};
701     for(Int_t itdc=0; itdc<32; itdc++){
702       for(Int_t i=0; i<4; i++) tdc[itdc] |= esdZDC->GetZDCTDCData(itdc, i)!=0;
703       if(fillHists && tdc[itdc]) fHistTDCZDC->Fill(itdc);
704     }
705     zdcNA = tdc[12];
706     zdcNC = tdc[10];
707     zdcPA = tdc[13];
708     zdcPC = tdc[11];
709   }
710   
711   if (side == kASide) return ((useZP && zdcPA) || (useZN && zdcNA)); 
712   if (side == kCSide) return ((useZP && zdcPC) || (useZN && zdcNC)); 
713   return kFALSE;
714 }
715
716
717 //-------------------------------------------------------------------------------------------------
718 Bool_t AliTriggerAnalysis::ZDCTimeTrigger(const AliESDEvent *aEsd, Bool_t fillHists) const {
719   // This method implements a selection based on the timing in both sides of zdcN
720   // It can be used in order to eliminate parasitic collisions
721   AliESDZDC *esdZDC = aEsd->GetESDZDC();
722   
723   if(fMC) {
724     UInt_t esdFlag =  esdZDC->GetESDQuality();
725     Bool_t znaFired  = (esdFlag & 0x01) == 0x01;
726     Bool_t zncFired  = (esdFlag & 0x10) == 0x10;
727     return znaFired | zncFired;
728   }
729   else {
730     for(Int_t i=0;i<4;++i) {
731       if (esdZDC->GetZDCTDCData(10,i)==0) continue;
732       Float_t tdcC = 0.025*(esdZDC->GetZDCTDCData(10,i)-esdZDC->GetZDCTDCData(14,i)); 
733       Float_t tdcCcorr = esdZDC->GetZDCTDCCorrected(10,i); 
734       for(Int_t j=0;j<4;++j) {
735         if (esdZDC->GetZDCTDCData(12,j)==0) continue;
736         Float_t tdcA = 0.025*(esdZDC->GetZDCTDCData(12,j)-esdZDC->GetZDCTDCData(14,j));
737         Float_t tdcAcorr = esdZDC->GetZDCTDCCorrected(12,j);
738         if(fillHists) {
739           fHistTimeZDC->Fill(tdcC-tdcA,tdcC+tdcA);
740           fHistTimeCorrZDC->Fill(tdcCcorr-tdcAcorr,tdcCcorr+tdcAcorr);
741         }
742         if (esdZDC->TestBit(AliESDZDC::kCorrectedTDCFilled)) {
743           if (TMath::Power((tdcCcorr-tdcAcorr-fZDCCutRefDeltaCorr)/fZDCCutSigmaDeltaCorr,2.)+
744               TMath::Power((tdcCcorr+tdcAcorr-fZDCCutRefSumCorr  )/fZDCCutSigmaSumCorr,2.) < 1.) return kTRUE;
745         }
746         else {
747           if (TMath::Power((tdcC-tdcA-fZDCCutRefDelta)/fZDCCutSigmaDelta,2.)+
748               TMath::Power((tdcC+tdcA-fZDCCutRefSum  )/fZDCCutSigmaSum,2.  )<1.0) return kTRUE;
749         }
750       }
751     }
752   }
753   return kFALSE;
754 }
755
756
757 //-------------------------------------------------------------------------------------------------
758 Bool_t AliTriggerAnalysis::ZDCTimeBGTrigger(const AliESDEvent *aEsd, AliceSide side) const{
759   // This method implements a selection based on the timing in zdcN
760   // It can be used in order to flag background
761   if(fMC) return kFALSE;
762   
763   AliESDZDC* zdcData = aEsd->GetESDZDC();
764   Bool_t znabadhit = kFALSE;
765   Bool_t zncbadhit = kFALSE;
766   
767   Float_t tdcCcorr=999, tdcAcorr=999;
768   for(Int_t i = 0; i < 4; ++i) {
769     if (zdcData->GetZDCTDCData(10,i)==0) continue;
770     tdcCcorr = TMath::Abs(zdcData->GetZDCTDCCorrected(10,i));
771     if(tdcCcorr<fZDCCutZNCTimeCorrMax && tdcCcorr>fZDCCutZNCTimeCorrMin) zncbadhit = kTRUE;
772   }
773   for(Int_t i = 0; i < 4; ++i) {
774     if (zdcData->GetZDCTDCData(12,i)==0) continue;
775     tdcAcorr = TMath::Abs(zdcData->GetZDCTDCCorrected(12,i));
776     if(tdcAcorr<fZDCCutZNATimeCorrMax && tdcAcorr>fZDCCutZNATimeCorrMin) znabadhit = kTRUE;
777   }
778   
779   if (side == kASide) return znabadhit;
780   if (side == kCSide) return zncbadhit;
781
782   return kFALSE;
783 }
784
785
786 //-------------------------------------------------------------------------------------------------
787 Bool_t AliTriggerAnalysis::ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) const {
788   // Returns if ZDC triggered
789   
790   AliESDZDC* zdcData = aEsd->GetESDZDC();
791   UInt_t quality = zdcData->GetESDQuality();
792
793   // from Nora's presentation, general first physics meeting 16.10.09
794   static UInt_t zpc  = 0x20;
795   static UInt_t znc  = 0x10;
796   static UInt_t zem1 = 0x08;
797   static UInt_t zem2 = 0x04;
798   static UInt_t zpa  = 0x02;
799   static UInt_t zna  = 0x01;
800   
801   if (side == kASide         && ((quality & zpa)  || (quality & zna ))) return kTRUE;
802   if (side == kCentralBarrel && ((quality & zem1) || (quality & zem2))) return kTRUE;
803   if (side == kCSide         && ((quality & zpc)  || (quality & znc ))) return kTRUE;
804   
805   return kFALSE;
806 }
807
808
809 //-------------------------------------------------------------------------------------------------
810 Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists){
811   // returns number of hit combinations above threshold
812   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
813   if (!fDoFMD) return -1;
814   
815   // Workaround for AliESDEvent::GetFMDData is not const!
816   const AliESDFMD* fmdData = (const_cast<AliESDEvent*>(aEsd))->GetFMDData();
817   if (!fmdData) {
818     AliError("AliESDFMD not available");
819     return -1;
820   }
821   
822   Int_t detFrom = (side == kASide) ? 1 : 3;
823   Int_t detTo   = (side == kASide) ? 2 : 3;
824   
825   Int_t triggers = 0;
826   Float_t totalMult = 0;
827   for (UShort_t det=detFrom;det<=detTo;det++) {
828     Int_t nRings = (det == 1 ? 1 : 2);
829     for (UShort_t ir = 0; ir < nRings; ir++) {
830       Char_t   ring = (ir == 0 ? 'I' : 'O');
831       UShort_t nsec = (ir == 0 ? 20  : 40);
832       UShort_t nstr = (ir == 0 ? 512 : 256);
833       for (UShort_t sec =0; sec < nsec;  sec++) {
834         for (UShort_t strip = 0; strip < nstr; strip++) {
835           Float_t mult = fmdData->Multiplicity(det,ring,sec,strip);
836           if (mult == AliESDFMD::kInvalidMult) continue;
837           if (fillHists) fHistFMDSingle->Fill(mult);
838           if (mult > fFMDLowCut)
839             totalMult = totalMult + mult;
840           else {
841             if (totalMult > fFMDHitCut) triggers++;
842             if (fillHists) fHistFMDSum->Fill(totalMult);
843             totalMult = 0;
844           }
845         }
846       }
847     }
848   }
849   return triggers;
850 }
851
852
853 //-------------------------------------------------------------------------------------------------
854 Bool_t AliTriggerAnalysis::FMDTrigger(const AliESDEvent* aEsd, AliceSide side){
855   // Returns if the FMD triggered
856   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
857   return FMDHitCombinations(aEsd, side, kFALSE);
858 }
859
860
861 //-------------------------------------------------------------------------------------------------
862 AliTriggerAnalysis::T0Decision AliTriggerAnalysis::T0Trigger(const AliESDEvent* aEsd, Bool_t online, Bool_t fillHists){
863   // Returns the T0 TVDC trigger decision
864   //  
865   // argument 'online' is used as a switch between online and offline trigger algorithms
866   // in online mode return 0TVX 
867   // in offline mode in addition check pile-up and background :
868   // pile-up readed from ESD: check if TVDC (0TVX module name) has more 1 hit;
869   // backgroud flag readed from ESD : check in given time interval OrA and OrC were correct but TVDC not  
870   // 
871   // Based on an algorithm by Alla Maevskaya 
872   
873   const AliESDTZERO* esdT0 = aEsd->GetESDTZERO();
874   if (!esdT0) {
875     AliError("AliESDTZERO not available");
876     return kT0Invalid;
877   }
878
879   Float_t  tvdc0 = esdT0->GetTVDC(0);
880   if(fillHists) fHistT0->Fill(tvdc0);
881   
882   if (online) {
883     if( aEsd->GetHeader()->IsTriggerInputFired("0TVX")) return kT0BB;
884   }
885   else {
886     if (esdT0->GetPileupFlag()) return kT0DecPileup;
887     if (esdT0->GetBackgroundFlag()) return kT0DecBG;
888     if (tvdc0>-5 && tvdc0<5 && tvdc0!=0) return kT0BB; 
889   }
890   
891   if (fMC) if(esdT0->GetT0zVertex()>-12.3 && esdT0->GetT0zVertex() < 10.3) return kT0BB; 
892   return kT0Empty;
893 }
894
895 //-------------------------------------------------------------------------------------------------
896 Bool_t AliTriggerAnalysis::EMCALCellsTrigger(const AliESDEvent *aEsd){
897   //
898   // Returns the EMCAL trigger decision
899   // so far only implemented for LHC11a data
900   // see http://alisoft.cern.ch/viewvc/trunk/PWGGA/EMCALTasks/AliEmcalPhysicsSelection.cxx?view=markup&root=AliRoot Revision 56136
901   //
902   
903   Bool_t isFired = kTRUE;
904   const Int_t runNumber = aEsd->GetRunNumber();
905   
906   /*
907   // Load EMCAL branches from the manager
908   AliAnalysisManager *am = AliAnalysisManager::GetAnalysisManager();
909   am->LoadBranch("EMCALCells.");
910   am->LoadBranch("CaloClusters");
911    */
912   
913   // Get EMCAL cells
914   AliVCaloCells *cells = aEsd->GetEMCALCells();
915   const Short_t nCells = cells->GetNumberOfCells();
916   
917   // count cells above threshold per sm
918   Int_t nCellCount[10] = {0,0,0,0,0,0,0,0,0,0};
919   for(Int_t iCell=0; iCell<nCells; ++iCell) {
920     Short_t cellId = cells->GetCellNumber(iCell);
921     Double_t cellE = cells->GetCellAmplitude(cellId);
922     Int_t sm       = cellId / (24*48);
923     if (cellE>0.1)
924       ++nCellCount[sm];
925   }
926   
927   // Trigger decision for LHC11a
928   Bool_t isLedEvent = kFALSE;
929   if ((runNumber>=144871) && (runNumber<=146860)) {
930     if (nCellCount[4] > 100)
931       isLedEvent = kTRUE;
932     else {
933       if ((runNumber>=146858) && (runNumber<=146860)) {
934         if (nCellCount[3]>=35)
935           isLedEvent = kTRUE;
936       }
937     }
938   }
939   
940   if (isLedEvent) {
941     isFired = kFALSE;
942   }
943   
944   return isFired;
945 }
946
947
948 //----------------------------------------------------------------------------------------------------
949 Bool_t AliTriggerAnalysis::TRDTrigger(const AliESDEvent *esd, Trigger trigger){
950   // evaluate the TRD trigger conditions,
951   // so far HCO, HSE, HQU, HJT, HEE
952   if(trigger!=kTRDHCO && trigger!=kTRDHJT && trigger!=kTRDHSE && trigger!=kTRDHQU && trigger!=kTRDHEE) {
953     AliWarning("Beware you are erroneously trying to use this function (wrong trigger)");
954     return kFALSE;
955   }
956   
957   Int_t nTrdTracks = esd->GetNumberOfTrdTracks();
958   if (nTrdTracks<=0) return kFALSE;
959   if      (trigger==kTRDHCO) return kTRUE;
960   else if (trigger!=kTRDHJT) {
961     for (Int_t iTrack = 0; iTrack < nTrdTracks; ++iTrack) {
962       AliESDTrdTrack *trdTrack = esd->GetTrdTrack(iTrack);
963       if (!trdTrack) continue;
964       // for the electron triggers we only consider matched tracks
965       if(trigger==kTRDHQU) if (TMath::Abs(trdTrack->Pt())>fTRDptHQU && trdTrack->GetPID()>fTRDpidHQU) return kTRUE;
966       if(trigger==kTRDHSE) if (TMath::Abs(trdTrack->Pt())>fTRDptHSE && trdTrack->GetPID()>fTRDpidHSE) return kTRUE; 
967       if(trigger==kTRDHEE) if (TMath::Abs(trdTrack->Pt())>fTRDptHSE && trdTrack->GetPID()>fTRDpidHSE && trdTrack->GetSector()>=fTRDminSectorHEE && trdTrack->GetSector()<=fTRDmaxSectorHEE) return kTRUE;
968     }
969   } 
970   else if (trigger==kTRDHJT) {
971     Int_t nTracks[90] = { 0 }; // stack-wise counted number of tracks above pt threshold
972     for (Int_t iTrack = 0; iTrack < nTrdTracks; ++iTrack) {
973       AliESDTrdTrack *trdTrack = esd->GetTrdTrack(iTrack);    
974       if (!trdTrack) continue;
975       Int_t globalStack = 5*trdTrack->GetSector() + trdTrack->GetStack();
976       // stack-wise counting of tracks above pt threshold for jet trigger
977       if (TMath::Abs(trdTrack->GetPt()) >= fTRDptHJT) ++nTracks[globalStack];
978     }
979     // check if HJT condition is fulfilled in any stack
980     for (Int_t iStack = 0; iStack < 90; iStack++) if (nTracks[iStack] >= fTRDnHJT) return kTRUE;
981   }
982   return kFALSE;
983 }
984
985
986 //-------------------------------------------------------------------------------------------------
987 Bool_t AliTriggerAnalysis::EMCALTrigger(const AliVEvent* event, Trigger trigger){
988   AliVCaloTrigger* emcalTrigger = event->GetCaloTrigger("EMCAL");
989   if (!emcalTrigger) return kFALSE;
990   Int_t emcalTriggerBits = 0;
991   emcalTrigger->GetTriggerBits(emcalTriggerBits);
992   if      (trigger==kEmcalL0         ) { return emcalTriggerBits & 1<<0; }
993   else if (trigger==kEmcalL1GammaHigh) { return emcalTriggerBits & 1<<1; }
994   else if (trigger==kEmcalL1GammaLow ) { return emcalTriggerBits & 1<<2; }
995   else if (trigger==kEmcalL1JetHigh  ) { return emcalTriggerBits & 1<<3; }
996   else if (trigger==kEmcalL1JetLow   ) { return emcalTriggerBits & 1<<4; }
997   else {
998     AliWarning("Beware you are erroneously trying to use this function (wrong trigger)");
999     return kFALSE;
1000   }
1001   return kFALSE;
1002 }
1003
1004
1005 //-------------------------------------------------------------------------------------------------
1006 Bool_t AliTriggerAnalysis::IsSPDClusterVsTrackletBG(const AliVEvent* event, Bool_t fillHists){
1007   // rejects BG based on the cluster vs tracklet correlation
1008   // returns true if the event is BG
1009   const AliVMultiplicity* mult = event->GetMultiplicity();
1010   if (!mult) { AliFatal("No multiplicity object"); }
1011   Int_t ntracklet   = mult->GetNumberOfTracklets();
1012   Int_t spdClusters = event->GetNumberOfITSClusters(0) + event->GetNumberOfITSClusters(1);
1013   if(fillHists) fHistSPDClsVsTrk->Fill(ntracklet,spdClusters);
1014   return spdClusters > Float_t(fASPDCvsTCut) + Float_t(ntracklet)*fBSPDCvsTCut;
1015 }
1016
1017
1018 //-------------------------------------------------------------------------------------------------
1019 Bool_t AliTriggerAnalysis::IsLaserWarmUpTPCEvent(const AliESDEvent* esd){
1020   // This function flags noisy TPC events which can happen during laser warm-up.
1021   Int_t trackCounter = 0;
1022   for (Int_t i=0; i<esd->GetNumberOfTracks(); ++i) {
1023     AliESDtrack *track = esd->GetTrack(i);
1024     if (!track) continue;
1025     if (track->GetTPCNcls() < 30) continue;
1026     if (TMath::Abs(track->Eta()) > 0.005) continue;
1027     if (track->Pt() < 4) continue;
1028     if (track->GetKinkIndex(0) > 0) continue;
1029     UInt_t status = track->GetStatus();
1030     if ((status&AliESDtrack::kITSrefit)==AliESDtrack::kITSrefit) continue; // explicitly ask for tracks without ITS refit
1031     if ((status&AliESDtrack::kTPCrefit)!=AliESDtrack::kTPCrefit) continue;
1032     if (track->GetTPCsignal() > 10) continue;          // explicitly ask for tracks without dE/dx
1033     trackCounter++;
1034   }
1035   if (trackCounter > 15) return kTRUE;
1036   return kFALSE;
1037 }
1038
1039
1040 //-------------------------------------------------------------------------------------------------
1041 Bool_t AliTriggerAnalysis::IsHVdipTPCEvent(const AliESDEvent* esd) {
1042   // This function flags events in which the TPC chamber HV is not at its nominal value
1043   if (fMC) return kFALSE; // there are no dip events in MC
1044   if (!esd->IsDetectorOn(AliDAQ::kTPC)) return kTRUE;
1045   return kFALSE;
1046 }
1047
1048
1049 //-------------------------------------------------------------------------------------------------
1050 Bool_t AliTriggerAnalysis::IsIncompleteEvent(const AliESDEvent* esd){
1051   // Check whether the event is incomplete 
1052   // (due to DAQ-HLT issues, it could be only part of the event was saved)
1053   if (fMC) return kFALSE; // there are no incomplete events on MC
1054   if ((esd->GetEventType() == 7) &&
1055       (esd->GetDAQDetectorPattern() & (1<<4)) &&
1056       !(esd->GetDAQAttributes() & (1<<7))) return kTRUE;
1057   return kFALSE;
1058 }
1059
1060
1061 //-------------------------------------------------------------------------------------------------
1062 Long64_t AliTriggerAnalysis::Merge(TCollection* list){
1063   // Merge a list of objects with this (needed for PROOF).
1064   // Returns the number of merged objects (including this).
1065   if (!list) return 0;
1066   if (list->IsEmpty()) return 1;
1067   TIterator* iter = list->MakeIterator();
1068   TObject* obj;
1069   
1070   // collections of all histograms
1071   const Int_t nHists = 14;
1072   TList collections[nHists];
1073   
1074   Int_t count = 0;
1075   while ((obj = iter->Next())) {
1076     AliTriggerAnalysis* entry = dynamic_cast<AliTriggerAnalysis*> (obj);
1077     if (entry == 0) continue;
1078     Int_t n = 0;
1079     collections[n++].Add(entry->fHistV0A);
1080     collections[n++].Add(entry->fHistV0C);
1081     collections[n++].Add(entry->fHistZDC);
1082     collections[n++].Add(entry->fHistTDCZDC);
1083     collections[n++].Add(entry->fHistTimeZDC);
1084     collections[n++].Add(entry->fHistTimeCorrZDC);
1085     collections[n++].Add(entry->fHistFMDA);
1086     collections[n++].Add(entry->fHistFMDC);
1087     collections[n++].Add(entry->fHistFMDSingle);
1088     collections[n++].Add(entry->fHistFMDSum);
1089     collections[n++].Add(entry->fHistBitsSPD);
1090     collections[n++].Add(entry->fHistFiredBitsSPD);
1091     collections[n++].Add(entry->fHistSPDClsVsTrk);
1092     collections[n++].Add(entry->fHistT0);
1093     
1094     // merge fTriggerClasses
1095     TIterator* iter2 = entry->fTriggerClasses->MakeIterator();
1096     TObjString* obj2 = 0;
1097     while ((obj2 = dynamic_cast<TObjString*> (iter2->Next()))) {
1098       TParameter<Long64_t>* param2 = static_cast<TParameter<Long64_t>*> (entry->fTriggerClasses->GetValue(obj2));
1099       TParameter<Long64_t>* param1 = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj2));
1100       if (param1) { param1->SetVal(param1->GetVal() + param2->GetVal()); }
1101       else {        
1102         param1 = dynamic_cast<TParameter<Long64_t>*> (param2->Clone());
1103         fTriggerClasses->Add(new TObjString(obj2->String()), param1);
1104       }
1105     }
1106     delete iter2;
1107     count++;
1108   }
1109   
1110   Int_t n = 0;
1111   fHistV0A->Merge(&collections[n++]);
1112   fHistV0C->Merge(&collections[n++]);
1113   fHistZDC->Merge(&collections[n++]);
1114   fHistTDCZDC->Merge(&collections[n++]);
1115   if (fHistTimeZDC)     fHistTimeZDC->Merge(&collections[n++]);     else n++;
1116   if (fHistTimeCorrZDC) fHistTimeCorrZDC->Merge(&collections[n++]); else n++;
1117   fHistFMDA->Merge(&collections[n++]);
1118   fHistFMDC->Merge(&collections[n++]);
1119   fHistFMDSingle->Merge(&collections[n++]);
1120   fHistFMDSum->Merge(&collections[n++]);
1121   fHistBitsSPD->Merge(&collections[n++]);
1122   fHistFiredBitsSPD->Merge(&collections[n++]);
1123   fHistSPDClsVsTrk->Merge(&collections[n++]);
1124   fHistT0->Merge(&collections[n++]);
1125   delete iter;
1126   return count+1;
1127 }
1128
1129
1130 //-------------------------------------------------------------------------------------------------
1131 void AliTriggerAnalysis::FillHistograms(const AliESDEvent* aEsd){
1132   // fills the histograms with the info from the ESD
1133   
1134   fHistBitsSPD->Fill(SPDFiredChips(aEsd, 0), SPDFiredChips(aEsd, 1, kTRUE));
1135   
1136   V0Trigger(aEsd, kASide, kFALSE, kTRUE);
1137   V0Trigger(aEsd, kCSide, kFALSE, kTRUE);
1138   T0Trigger(aEsd, kFALSE, kTRUE);
1139   ZDCTDCTrigger(aEsd,kASide,kFALSE,kFALSE,kTRUE);
1140   ZDCTimeTrigger(aEsd,kTRUE);
1141   IsSPDClusterVsTrackletBG(aEsd, kTRUE);
1142   
1143   AliESDZDC* zdcData = aEsd->GetESDZDC();
1144   if (zdcData)  {
1145     UInt_t quality = zdcData->GetESDQuality();
1146     
1147     // from Nora's presentation, general first physics meeting 16.10.09
1148     static UInt_t zpc  = 0x20;
1149     static UInt_t znc  = 0x10;
1150     static UInt_t zem1 = 0x08;
1151     static UInt_t zem2 = 0x04;
1152     static UInt_t zpa  = 0x02;
1153     static UInt_t zna  = 0x01;
1154     
1155     fHistZDC->Fill(1, (quality & zna)  ? 1 : 0);
1156     fHistZDC->Fill(2, (quality & zpa)  ? 1 : 0);
1157     fHistZDC->Fill(3, (quality & zem2) ? 1 : 0);
1158     fHistZDC->Fill(4, (quality & zem1) ? 1 : 0);
1159     fHistZDC->Fill(5, (quality & znc)  ? 1 : 0);
1160     fHistZDC->Fill(6, (quality & zpc)  ? 1 : 0);
1161   }
1162   else {
1163     fHistZDC->Fill(-1);
1164     AliError("AliESDZDC not available");
1165   }
1166   
1167   if (fDoFMD) {
1168     fHistFMDA->Fill(FMDHitCombinations(aEsd, kASide, kTRUE));
1169     fHistFMDC->Fill(FMDHitCombinations(aEsd, kCSide, kTRUE));
1170   }
1171 }
1172
1173
1174 //-------------------------------------------------------------------------------------------------
1175 void AliTriggerAnalysis::SaveHistograms() const {
1176   // write histograms to current directory
1177   if (fHistBitsSPD)      fHistBitsSPD->Write();
1178   if (fHistFiredBitsSPD) fHistFiredBitsSPD->Write();
1179   if (fHistV0A)          fHistV0A->Write();
1180   if (fHistV0C)          fHistV0C->Write();
1181   if (fHistZDC)          fHistZDC->Write();
1182   if (fHistTDCZDC)       fHistTDCZDC->Write();
1183   if (fHistTimeZDC)      fHistTimeZDC->Write();
1184   if (fHistTimeCorrZDC)  fHistTimeCorrZDC->Write();
1185   if (fHistFMDA)         fHistFMDA->Write();
1186   if (fHistFMDC)         fHistFMDC->Write();
1187   if (fHistFMDSingle)    fHistFMDSingle->Write();
1188   if (fHistFMDSum)       fHistFMDSum->Write();
1189   if (fSPDGFOEfficiency) fSPDGFOEfficiency->Write("fSPDGFOEfficiency");
1190   if (fHistSPDClsVsTrk)  fHistSPDClsVsTrk->Write("fHistSPDClsVsTrk");
1191   if (fHistT0)           fHistT0->Write("fHistT0");
1192   fTriggerClasses->Write("fTriggerClasses", TObject::kSingleKey);
1193 }
1194
1195
1196 //-------------------------------------------------------------------------------------------------
1197 void AliTriggerAnalysis::FillTriggerClasses(const AliESDEvent* aEsd){
1198   // fills trigger classes map
1199   TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(aEsd->GetFiredTriggerClasses().Data()));
1200   if (!count) {
1201     count = new TParameter<Long64_t>(aEsd->GetFiredTriggerClasses(), 0);
1202     fTriggerClasses->Add(new TObjString(aEsd->GetFiredTriggerClasses().Data()), count);
1203   }
1204   count->SetVal(count->GetVal() + 1);
1205 }
1206
1207
1208 //-------------------------------------------------------------------------------------------------
1209 void AliTriggerAnalysis::PrintTriggerClasses() const {
1210   // print trigger classes
1211   
1212   Printf("Trigger Classes:");
1213   Printf("Event count for trigger combinations:");
1214   TMap singleTrigger;
1215   singleTrigger.SetOwner();
1216   TIterator* iter = fTriggerClasses->MakeIterator();
1217   TObjString* obj = 0;
1218   while ((obj = dynamic_cast<TObjString*> (iter->Next()))) {
1219     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj));
1220     Printf(" %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1221     TObjArray* tokens = obj->String().Tokenize(" ");
1222     for (Int_t i=0; i<tokens->GetEntries(); i++) {
1223       TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(((TObjString*) tokens->At(i))->String().Data()));
1224       if (!count) {
1225         count = new TParameter<Long64_t>(((TObjString*) tokens->At(i))->String().Data(), 0);
1226         singleTrigger.Add(new TObjString(((TObjString*) tokens->At(i))->String().Data()), count);
1227       }
1228       count->SetVal(count->GetVal() + param->GetVal());
1229     }
1230     delete tokens;
1231   }
1232   delete iter;
1233   
1234   Printf("Event count for single trigger:");
1235   iter = singleTrigger.MakeIterator();
1236   while ((obj = dynamic_cast<TObjString*> (iter->Next()))) {
1237     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(obj));
1238     Printf("  %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1239   }
1240   delete iter;
1241   singleTrigger.DeleteAll();
1242 }