]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliTriggerAnalysis.cxx
Fix for LEGO train generation: the configuration macro name for a wagon between ...
[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 <Riostream.h>
26 #include <TH1F.h>
27 #include <TH2F.h>
28 #include <TList.h>
29 #include <TIterator.h>
30 #include "TParameter.h"
31 #include <TMap.h>
32 #include <TRandom.h>
33
34 #include <AliTriggerAnalysis.h>
35
36 #include <AliLog.h>
37 #include <AliAnalysisManager.h>
38
39 #include <AliESDEvent.h>
40
41 #include <AliMultiplicity.h>
42 #include <AliESDVZERO.h>
43 #include <AliESDZDC.h>
44 #include <AliESDFMD.h>
45 #include <AliESDVertex.h>
46 #include <AliESDtrackCuts.h>
47
48 ClassImp(AliTriggerAnalysis)
49
50 AliTriggerAnalysis::AliTriggerAnalysis() :
51   fSPDGFOThreshold(2),
52   fSPDGFOEfficiency(0),
53   fV0TimeOffset(0),
54   fV0AdcThr(0),
55   fV0HwAdcThr(2.5),
56   fV0HwWinLow(61.5),
57   fV0HwWinHigh(86.5),
58   fZDCCutRefSum(-568.5),
59   fZDCCutRefDelta(-2.1),
60   fZDCCutSigmaSum(3.25),
61   fZDCCutSigmaDelta(2.25),
62   fZDCCutRefSumCorr(-65.5),
63   fZDCCutRefDeltaCorr(-2.1),
64   fZDCCutSigmaSumCorr(6.0),
65   fZDCCutSigmaDeltaCorr(1.2),
66   fZDCCutZNATimeCorr(2.0),
67   fZDCCutZNCTimeCorr(5.0),
68   fASPDCvsTCut(65),
69   fBSPDCvsTCut(4),
70   fDoFMD(kTRUE),
71   fFMDLowCut(0.2),
72   fFMDHitCut(0.5),
73   fHistBitsSPD(0),
74   fHistFiredBitsSPD(0),
75   fHistSPDClsVsTrk(0),
76   fHistV0A(0),       
77   fHistV0C(0),
78   fHistZDC(0),    
79   fHistTDCZDC(0),    
80   fHistTimeZDC(0),    
81   fHistTimeCorrZDC(0),    
82   fHistFMDA(0),    
83   fHistFMDC(0),   
84   fHistFMDSingle(0),
85   fHistFMDSum(0),
86   fHistT0(0),
87   fTriggerClasses(0),
88   fMC(kFALSE),
89   fEsdTrackCuts(0),
90   fTPCOnly(kFALSE)
91 {
92   // constructor
93 }
94
95 AliTriggerAnalysis::~AliTriggerAnalysis()
96 {
97   // destructor
98   
99   if (fHistBitsSPD)
100   {
101     delete fHistBitsSPD;
102     fHistBitsSPD = 0;
103   }
104
105   if (fHistFiredBitsSPD)
106   {
107     delete fHistFiredBitsSPD;
108     fHistFiredBitsSPD = 0;
109   }
110
111   if (fHistSPDClsVsTrk)
112   {
113     delete fHistSPDClsVsTrk;
114     fHistSPDClsVsTrk = 0;
115   }
116
117   if (fHistV0A)
118   {
119     delete fHistV0A;
120     fHistV0A = 0;
121   }
122
123   if (fHistV0C)
124   {
125     delete fHistV0C;
126     fHistV0C = 0;
127   }
128
129   if (fHistZDC)
130   {
131     delete fHistZDC;
132     fHistZDC = 0;
133   }
134   if (fHistTDCZDC)
135   {
136     delete fHistTDCZDC;
137     fHistTDCZDC = 0;
138   }
139   if (fHistTimeZDC)
140   {
141     delete fHistTimeZDC;
142     fHistTimeZDC = 0;
143   }
144   if (fHistTimeCorrZDC)
145   {
146     delete fHistTimeCorrZDC;
147     fHistTimeCorrZDC = 0;
148   }
149
150   if (fHistFMDA)
151   {
152     delete fHistFMDA;
153     fHistFMDA = 0;
154   }
155
156   if (fHistFMDC)
157   {
158     delete fHistFMDC;
159     fHistFMDC = 0;
160   }
161
162   if (fHistFMDSingle)
163   {
164     delete fHistFMDSingle;
165     fHistFMDSingle = 0;
166   }
167
168   if (fHistFMDSum)
169   {
170     delete fHistFMDSum;
171     fHistFMDSum = 0;
172   }
173
174   if (fHistT0)
175   {
176     delete fHistT0;
177     fHistT0 = 0;
178   }
179
180   if (fTriggerClasses)
181   {
182     fTriggerClasses->DeleteAll();
183     delete fTriggerClasses;
184     fTriggerClasses = 0;
185   }
186
187   if (fEsdTrackCuts){
188     delete fEsdTrackCuts;
189     fEsdTrackCuts =0;
190   }
191 }
192
193 void AliTriggerAnalysis::EnableHistograms()
194 {
195   // creates the monitoring histograms
196   
197   // do not add this hists to the directory
198   Bool_t oldStatus = TH1::AddDirectoryStatus();
199   TH1::AddDirectory(kFALSE);
200   
201   fHistBitsSPD = new TH2F("fHistBitsSPD", "SPD GFO;number of fired chips (offline);number of fired chips (hardware)", 1202, -1.5, 1200.5, 1202, -1.5, 1200.5);
202   fHistFiredBitsSPD = new TH1F("fHistFiredBitsSPD", "SPD GFO Hardware;chip number;events", 1200, -0.5, 1199.5);
203   fHistSPDClsVsTrk = new TH2F("fHistSPDClsVsTrk", "SPD Clusters vs Tracklets", 300, -0.5, 2999.5, 1000, -0.5, 9999.5);
204   fHistV0A = new TH1F("fHistV0A", "V0A;leading time (ns);events", 400, -100, 100);
205   fHistV0C = new TH1F("fHistV0C", "V0C;leading time (ns);events", 400, -100, 100);
206   fHistZDC = new TH1F("fHistZDC", "ZDC;trigger bits;events", 8, -1.5, 6.5);
207   fHistTDCZDC = new TH1F("fHistTDCZDC", "ZDC;TDC bits;events", 32, -0.5, 32-0.5);
208   fHistTimeZDC = new TH2F("fHistTimeZDC", "ZDC;TDC timing A+C vs C-A; events", 120,-30,30,120,-600,-540);
209   fHistTimeCorrZDC = new TH2F("fHistTimeCorrZDC", "ZDC;Corrected TDC timing A+C vs C-A; events", 120,-30,30,260,-100,30);
210   
211   // TODO check limits
212   fHistFMDA = new TH1F("fHistFMDA", "FMDA;combinations above threshold;events", 102, -1.5, 100.5);
213   fHistFMDC = new TH1F("fHistFMDC", "FMDC;combinations above threshold;events", 102, -1.5, 100.5);
214   fHistFMDSingle = new TH1F("fHistFMDSingle", "FMD single;multiplicity value;counts", 1000, 0, 10);
215   fHistFMDSum = new TH1F("fHistFMDSum", "FMD sum;multiplicity value;counts", 1000, 0, 10);
216   fHistT0 = new TH1F("fHistT0", "T0;time (ns);events", 100, -25, 25);
217   
218   fTriggerClasses = new TMap;
219   fTriggerClasses->SetOwner();
220   
221   TH1::AddDirectory(oldStatus);
222 }
223
224 //____________________________________________________________________
225 const char* AliTriggerAnalysis::GetTriggerName(Trigger trigger) 
226 {
227   // returns the name of the requested trigger
228   // the returned string will only be valid until the next call to this function [not thread-safe]
229   
230   static TString str;
231   
232   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
233   
234   switch (triggerNoFlags)
235   {
236     case kAcceptAll : str = "ACCEPT ALL (bypass!)"; break;
237     case kMB1 : str = "MB1"; break;
238     case kMB2 : str = "MB2"; break;
239     case kMB3 : str = "MB3"; break;
240     case kSPDGFO : str = "SPD GFO"; break;
241     case kSPDGFOBits : str = "SPD GFO Bits"; break;
242     case kSPDGFOL0 : str = "SPD GFO L0 (first layer)"; break;
243     case kSPDGFOL1 : str = "SPD GFO L1 (second layer)"; break;
244     case kSPDClsVsTrkBG :  str = "Cluster vs Tracklets"; break;
245     case kV0A : str = "V0 A BB"; break;
246     case kV0C : str = "V0 C BB"; break;
247     case kV0OR : str = "V0 OR BB"; break;
248     case kV0AND : str = "V0 AND BB"; break;
249     case kV0ABG : str = "V0 A BG"; break;
250     case kV0CBG : str = "V0 C BG"; break;
251     case kZDC : str = "ZDC"; break;
252     case kZDCA : str = "ZDC A"; break;
253     case kZDCC : str = "ZDC C"; break;
254     case kZNA : str = "ZN A"; break;
255     case kZNC : str = "ZN C"; break;
256     case kZNABG : str = "ZN A BG"; break;
257     case kZNCBG : str = "ZN C BG"; break;
258     case kFMDA : str = "FMD A"; break;
259     case kFMDC : str = "FMD C"; break;
260     case kFPANY : str = "SPD GFO | V0 | ZDC | FMD"; break;
261     case kNSD1 : str = "NSD1"; break;
262     case kMB1Prime: str = "MB1prime"; break;
263     case kZDCTDCA : str = "ZDC TDC A"; break;
264     case kZDCTDCC : str = "ZDC TDC C"; break;
265     case kZDCTime : str = "ZDC Time Cut"; break;
266     case kCentral : str = "V0 Central"; break;
267     case kSemiCentral : str = "V0 Semi-central"; break;
268     case kEMCAL : str = "EMCAL"; break;
269     default: str = ""; break;
270   }
271    
272   if (trigger & kOfflineFlag)
273     str += " OFFLINE";  
274   
275   if (trigger & kOneParticle)
276     str += " OneParticle";  
277
278   if (trigger & kOneTrack)
279     str += " OneTrack";  
280
281   return str;
282 }
283
284 Bool_t AliTriggerAnalysis::IsTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
285 {
286   // checks if an event has been triggered
287
288   if (trigger & kOfflineFlag)
289     return IsOfflineTriggerFired(aEsd, trigger);
290     
291   return IsTriggerBitFired(aEsd, trigger);
292 }
293
294 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, Trigger trigger) const
295 {
296   // checks if an event is fired using the trigger bits
297
298   return IsTriggerBitFired(aEsd->GetTriggerMask(), trigger);
299 }
300
301 Bool_t AliTriggerAnalysis::IsTriggerBitFired(ULong64_t triggerMask, Trigger trigger) const
302 {
303   // checks if an event is fired using the trigger bits
304   //
305   // this function needs the branch TriggerMask in the ESD
306   
307   // definitions from p-p.cfg
308   ULong64_t spdFO = (1 << 14);
309   ULong64_t v0left = (1 << 10);
310   ULong64_t v0right = (1 << 11);
311
312   switch (trigger)
313   {
314     case kAcceptAll:
315     {
316       return kTRUE;
317       break;
318     }
319     case kMB1:
320     {
321       if (triggerMask & spdFO || ((triggerMask & v0left) || (triggerMask & v0right)))
322         return kTRUE;
323       break;
324     }
325     case kMB2:
326     {
327       if (triggerMask & spdFO && ((triggerMask & v0left) || (triggerMask & v0right)))
328         return kTRUE;
329       break;
330     }
331     case kMB3:
332     {
333       if (triggerMask & spdFO && (triggerMask & v0left) && (triggerMask & v0right))
334         return kTRUE;
335       break;
336     }
337     case kSPDGFO:
338     {
339       if (triggerMask & spdFO)
340         return kTRUE;
341       break;
342     }
343     default:
344       Printf("IsEventTriggered: ERROR: Trigger type %d not implemented in this method", (Int_t) trigger);
345       break;
346   }
347
348   return kFALSE;
349 }
350
351 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, ULong64_t tclass) const
352 {
353   // Checks if corresponding bit in mask is on
354   
355   ULong64_t trigmask = aEsd->GetTriggerMask();
356   return (trigmask & (1ull << (tclass-1)));
357 }
358   
359 Int_t AliTriggerAnalysis::EvaluateTrigger(const AliESDEvent* aEsd, Trigger trigger)
360 {
361   // evaluates a given trigger "offline"
362   // trigger combinations are not supported, for that see IsOfflineTriggerFired
363   
364   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
365   Bool_t offline = kFALSE;
366   if (trigger & kOfflineFlag)
367     offline = kTRUE;
368   
369   Int_t decision = 0;
370   switch (triggerNoFlags)
371   {
372     case kSPDGFO:
373     {
374       decision = SPDFiredChips(aEsd, (offline) ? 0 : 1);
375       break;
376     }
377     case kSPDGFOL0:
378     {
379       decision = SPDFiredChips(aEsd, (offline) ? 0 : 1, kFALSE, 1);
380       break;
381     }
382     case kSPDGFOL1:
383     {
384       decision = SPDFiredChips(aEsd, (offline) ? 0 : 1, kFALSE, 2);
385       break;
386     } 
387     case kSPDClsVsTrkBG:
388     {
389       if(!offline) 
390         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
391       decision = IsSPDClusterVsTrackletBG(aEsd);
392       break;
393     }
394     case kV0A:
395     {
396       if (V0Trigger(aEsd, kASide, !offline) == kV0BB)
397         decision = 1;
398       break;
399     }
400     case kV0C:
401     {
402       if (V0Trigger(aEsd, kCSide, !offline) == kV0BB)
403         decision = 1;
404       break;
405     }
406     case kV0ABG:
407     {
408       if (V0Trigger(aEsd, kASide, !offline) == kV0BG)
409         decision = 1;
410       break;
411     }
412     case kV0CBG:
413     {
414       if (V0Trigger(aEsd, kCSide, !offline) == kV0BG)
415         decision = 1;
416       break;
417     }
418     case kT0:
419     {
420       if (T0Trigger(aEsd, !offline) == kT0BB)
421         decision = 1;
422       break;
423     }
424     case kT0BG:
425     {
426       if (!offline)
427         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
428       if (T0Trigger(aEsd, !offline) == kT0DecBG)
429         decision = 1;
430       break;
431     }
432     case kT0Pileup:
433     {
434       if (!offline)
435         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
436       if (T0Trigger(aEsd, !offline) == kT0DecPileup)
437         decision = 1;
438       break;
439     }
440     case kZDCA:
441     {
442       if (!offline)
443         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
444       if (ZDCTrigger(aEsd, kASide))
445         decision = 1;
446       break;
447     }
448     case kZDCC:
449     {
450       if (!offline)
451         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
452       if (ZDCTrigger(aEsd, kCSide))
453         decision = 1;
454       break;
455     }
456     case kZDCTDCA:
457     {
458       if (!offline)
459         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
460       if (ZDCTDCTrigger(aEsd, kASide))
461         decision = 1;
462       break;
463     }
464     case kZDCTDCC:
465     {
466       if (!offline)
467         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
468       if (ZDCTDCTrigger(aEsd, kCSide))
469         decision = 1;
470       break;
471     }
472     case kZDCTime:
473     {
474       if (!offline)
475         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
476       if (ZDCTimeTrigger(aEsd))
477         decision = 1;
478       break;
479     }
480     case kZNA:
481     {
482       if (!offline)
483         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
484       if (ZDCTDCTrigger(aEsd,kASide,kTRUE,kFALSE,kFALSE))
485         decision = 1;
486       break;
487     }
488     case kZNC:
489     {
490       if (!offline)
491         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
492       if (ZDCTDCTrigger(aEsd,kCSide,kTRUE,kFALSE,kFALSE))
493         decision = 1;
494       break;
495     }
496     case kZNABG:
497     {
498       if (!offline)
499         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
500       if (ZDCTimeBGTrigger(aEsd,kASide))
501         decision = 1;
502       break;
503     }
504     case kZNCBG:
505     {
506       if (!offline)
507         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
508       if (ZDCTimeBGTrigger(aEsd,kCSide))
509         decision = 1;
510       break;
511     }
512     case kFMDA:
513     {
514       if (!offline)
515         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
516       if (FMDTrigger(aEsd, kASide))
517         decision = 1;
518       break;
519     }
520     case kFMDC:
521     {
522       if (!offline)
523         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
524       if (FMDTrigger(aEsd, kCSide))
525         decision = 1;
526       break;
527     }
528     case kCTPV0A:
529     {
530       if (offline)
531         AliFatal(Form("Offline trigger not available for trigger %d", triggerNoFlags));
532       if (IsL0InputFired(aEsd, 2))
533         decision = 1;
534       break;
535     }
536     case kCTPV0C:
537     {
538       if (offline)
539         AliFatal(Form("Offline trigger not available for trigger %d", triggerNoFlags));
540       if (IsL0InputFired(aEsd, 3))
541         decision = 1;
542       break;
543     }
544     case kTPCLaserWarmUp:
545     {
546       if (!offline)
547         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
548       return IsLaserWarmUpTPCEvent(aEsd);
549     }
550     case kCentral:
551     {
552       if (offline)
553         AliFatal(Form("Offline trigger not available for trigger %d - use centrality selection", triggerNoFlags));
554       if (aEsd->GetVZEROData()) {
555         if (aEsd->GetVZEROData()->TestBit(AliESDVZERO::kTriggerChargeBitsFilled)) {
556           if (aEsd->GetVZEROData()->GetTriggerBits() & (1<<AliESDVZERO::kCTA2andCTC2))
557             decision = kTRUE;
558         }
559         else
560           AliWarning("V0 centrality trigger bits were not filled!");
561       }
562       break;
563     }
564     case kSemiCentral:
565     {
566       if (offline)
567         AliFatal(Form("Offline trigger not available for trigger %d - use centrality selection", triggerNoFlags));
568       if (aEsd->GetVZEROData()) {
569         if (aEsd->GetVZEROData()->TestBit(AliESDVZERO::kTriggerChargeBitsFilled)) {
570           if (aEsd->GetVZEROData()->GetTriggerBits() & (1<<AliESDVZERO::kCTA1andCTC1))
571             decision = kTRUE;
572         }
573         else
574           AliWarning("V0 centrality trigger bits were not filled!");
575       }
576       break;
577     }
578   case kEMCAL:
579     {
580       if(!offline) 
581         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
582       if (EMCALCellsTrigger(aEsd)) 
583         decision = kTRUE;
584       break;
585     }
586     default:
587     {
588       AliFatal(Form("Trigger type %d not implemented", triggerNoFlags));
589     }
590   }  
591   
592   return decision;
593 }
594
595 Bool_t AliTriggerAnalysis::IsLaserWarmUpTPCEvent(const AliESDEvent* esd)
596 {
597   //
598   // This function flags noisy TPC events which can happen during laser warm-up.
599   //
600   
601   Int_t trackCounter = 0;
602   for (Int_t i=0; i<esd->GetNumberOfTracks(); ++i) 
603   {
604     AliESDtrack *track = esd->GetTrack(i);
605     if (!track) 
606       continue;
607       
608     if (track->GetTPCNcls() < 30) continue;
609     if (TMath::Abs(track->Eta()) > 0.005) continue;
610     if (track->Pt() < 4) continue;
611     if (track->GetKinkIndex(0) > 0) continue;
612     
613     UInt_t status = track->GetStatus();
614     if ((status&AliESDtrack::kITSrefit)==AliESDtrack::kITSrefit) continue; // explicitly ask for tracks without ITS refit
615     if ((status&AliESDtrack::kTPCrefit)!=AliESDtrack::kTPCrefit) continue;
616     
617     if (track->GetTPCsignal() > 10) continue;          // explicitly ask for tracks without dE/dx
618     
619     trackCounter++;
620   }
621   if (trackCounter > 15) 
622     return kTRUE;
623   return kFALSE;
624 }
625
626 Bool_t AliTriggerAnalysis::IsOfflineTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
627 {
628   // checks if an event has been triggered "offline"
629
630   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
631   
632   Bool_t decision = kFALSE;
633   switch (triggerNoFlags)
634   {
635     case kAcceptAll:
636     {
637       decision = kTRUE;
638       break;
639     }
640     case kMB1:
641     {
642       if (SPDGFOTrigger(aEsd, 0) || V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
643         decision = kTRUE;
644       break;
645     }
646     case kMB2:
647     {
648       if (SPDGFOTrigger(aEsd, 0) && (V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB))
649         decision = kTRUE;
650       break;
651     }
652     case kMB3:
653     {
654       if (SPDGFOTrigger(aEsd, 0) && V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
655         decision = kTRUE;
656       break;
657     }
658     case kSPDGFO:
659     {
660       if (SPDGFOTrigger(aEsd, 0))
661         decision = kTRUE;
662       break;
663     }
664     case kSPDGFOBits:
665     {
666       if (SPDGFOTrigger(aEsd, 1))
667         decision = kTRUE;
668       break;
669     }
670     case kV0A:
671     {
672       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB)
673         decision = kTRUE;
674       break;
675     }
676     case kV0C:
677     {
678       if (V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
679         decision = kTRUE;
680       break;
681     }
682     case kV0OR:
683     {
684       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
685         decision = kTRUE;
686       break;
687     }
688     case kV0AND:
689     {
690       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
691         decision = kTRUE;
692       break;
693     }
694     case kV0ABG:
695     {
696       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BG)
697         decision = kTRUE;
698       break;
699     }
700     case kV0CBG:
701     {
702       if (V0Trigger(aEsd, kCSide, kFALSE) == kV0BG)
703         decision = kTRUE;
704       break;
705     }
706     case kZDC:
707     {
708       if (ZDCTrigger(aEsd, kASide) || ZDCTrigger(aEsd, kCentralBarrel) || ZDCTrigger(aEsd, kCSide))
709         decision = kTRUE;
710       break;
711     }
712     case kZDCA:
713     {
714       if (ZDCTrigger(aEsd, kASide))
715         decision = kTRUE;
716       break;
717     }
718     case kZDCC:
719     {
720       if (ZDCTrigger(aEsd, kCSide))
721         decision = kTRUE;
722       break;
723     }
724     case kZNA:
725     {
726       if (ZDCTDCTrigger(aEsd,kASide,kTRUE,kFALSE,kFALSE))
727         decision = kTRUE;
728       break;
729     }
730     case kZNC:
731     {
732       if (ZDCTDCTrigger(aEsd,kCSide,kTRUE,kFALSE,kFALSE))
733         decision = kTRUE;
734       break;
735     }
736     case kZNABG:
737     {
738       if (ZDCTimeBGTrigger(aEsd,kASide))
739         decision = kTRUE;
740       break;
741     }
742     case kZNCBG:
743     {
744       if (ZDCTimeBGTrigger(aEsd,kCSide))
745         decision = kTRUE;
746       break;
747     }
748     case kFMDA:
749     {
750       if (FMDTrigger(aEsd, kASide))
751         decision = kTRUE;
752       break;
753     }
754     case kFMDC:
755     {
756       if (FMDTrigger(aEsd, kCSide))
757         decision = kTRUE;
758       break;
759     }
760     case kFPANY:
761     {
762       if (SPDGFOTrigger(aEsd, 0) || V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB || ZDCTrigger(aEsd, kASide) || ZDCTrigger(aEsd, kCentralBarrel) || ZDCTrigger(aEsd, kCSide) || FMDTrigger(aEsd, kASide) || FMDTrigger(aEsd, kCSide))
763         decision = kTRUE;
764       break;
765     }
766     case kNSD1:
767     {
768       if (SPDFiredChips(aEsd, 0) >= 5 || (V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB))
769         decision = kTRUE;
770        break;
771     }
772     case kMB1Prime:
773     {
774       Int_t count = 0;
775       if (SPDGFOTrigger(aEsd, 0))
776         count++;
777       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB)
778         count++;
779       if (V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
780         count++;
781       
782       if (count >= 2)
783         decision = kTRUE;
784         
785       break;
786     }
787   case kEMCAL:
788     {
789       if (EMCALCellsTrigger(aEsd))
790         decision = kTRUE;
791       break;
792     }
793   default:
794       {
795       AliFatal(Form("Trigger type %d not implemented", triggerNoFlags));
796     }
797   }
798   
799   // hadron-level requirement
800   if (decision && (trigger & kOneParticle))
801   {
802     decision = kFALSE;
803     
804     const AliESDVertex* vertex = aEsd->GetPrimaryVertexSPD();
805     const AliMultiplicity* mult = aEsd->GetMultiplicity();
806
807     if (mult && vertex && vertex->GetNContributors() > 0 && (!vertex->IsFromVertexerZ() || vertex->GetDispersion() < 0.02) && TMath::Abs(vertex->GetZv()) < 5.5) 
808     {
809       for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
810       {
811         if (TMath::Abs(mult->GetEta(i)) < 1)
812         {
813           decision = kTRUE;
814           break;
815         }
816       }
817     }
818   }
819
820   // hadron level definition for TPC tracks
821
822   if (decision && (trigger & kOneTrack))
823   {
824     decision = kFALSE;
825     const AliESDVertex* vertex =0x0;
826     vertex = aEsd->GetPrimaryVertexTracks();
827     if (!vertex || vertex->GetNContributors() <= 0)
828     {
829       vertex = aEsd->GetPrimaryVertexSPD();
830     }
831     Float_t ptmin, ptmax;
832     fEsdTrackCuts->GetPtRange(ptmin,ptmax);
833     AliDebug(3, Form("ptmin = %f, ptmax = %f\n",ptmin, ptmax));
834
835     if (vertex && vertex->GetNContributors() > 0 && (!vertex->IsFromVertexerZ() || vertex->GetDispersion() < 0.02) && TMath::Abs(vertex->GetZv()) < 10.) {
836       AliDebug(3,Form("Check on the vertex passed\n"));
837       for (Int_t i=0; i<aEsd->GetNumberOfTracks(); ++i){
838         if (fTPCOnly == kFALSE){
839           if (fEsdTrackCuts->AcceptTrack(aEsd->GetTrack(i))){
840             AliDebug(2, Form("pt of track = %f --> check passed\n",aEsd->GetTrack(i)->Pt()));
841             decision = kTRUE;
842              break;
843           }
844         }
845         else {
846           // TPC only tracks
847           AliESDtrack *tpcTrack = fEsdTrackCuts->GetTPCOnlyTrack((AliESDEvent*)aEsd, i);
848           if (!tpcTrack){
849             AliDebug(3,Form("track %d is NOT a TPC track",i));
850             continue;
851           }
852           else{
853             AliDebug(3,Form("track %d IS a TPC track",i));
854             if (!(fEsdTrackCuts->AcceptTrack(tpcTrack))) {
855               AliDebug(2, Form("TPC track %d NOT ACCEPTED, pt = %f, eta = %f",i,tpcTrack->Pt(), tpcTrack->Eta()));
856               delete tpcTrack; tpcTrack = 0x0;
857               continue;
858             } // end if the TPC track is not accepted
859             else{
860               AliDebug(2, Form("TPC track %d ACCEPTED, pt = %f, eta = %f",i,tpcTrack->Pt(), tpcTrack->Eta()));
861               decision = kTRUE;
862               delete tpcTrack; tpcTrack = 0x0;
863               break;
864             } // end if the TPC track is accepted
865           } // end if it is a TPC track
866         } // end if you are looking at TPC only tracks                  
867       } // end loop on tracks
868     } // end check on vertex
869     else{
870       AliDebug(4,Form("Check on the vertex not passed\n"));
871       for (Int_t i=0; i<aEsd->GetNumberOfTracks(); ++i){
872         if (fEsdTrackCuts->AcceptTrack(aEsd->GetTrack(i))){
873           AliDebug(4,Form("pt of track = %f --> check would be passed if the vertex was ok\n",aEsd->GetTrack(i)->Pt()));
874           break;
875         }
876       }
877     }
878     if (!decision) AliDebug(3,("Check for kOneTrack NOT passed\n"));
879   }
880
881   return decision;
882 }
883
884 Bool_t AliTriggerAnalysis::IsTriggerClassFired(const AliESDEvent* aEsd, const Char_t* tclass) const 
885 {
886   // tclass is logical function of inputs, e.g. 01 && 02 || 03 && 11 && 21
887   // = L0 inp 1 && L0 inp 2 || L0 inp 3 && L1 inp 1 && L2 inp 1
888   // NO brackets in logical function !
889   // Spaces between operators and inputs.
890   // Not all logical functions are available in CTP= 
891   // =any function of first 4 inputs; 'AND' of other inputs, check not done
892   // This method will be replaced/complemened by similar one
893   // which works withh class and inputs names as in CTP cfg file
894   
895   TString TClass(tclass);
896   TObjArray* tcltokens = TClass.Tokenize(" ");
897   Char_t level=((TObjString*)tcltokens->At(0))->String()[0];
898   UInt_t input=atoi((((TObjString*)tcltokens->At(0))->String()).Remove(0));
899   Bool_t tcl = IsInputFired(aEsd,level,input);
900  
901   for (Int_t i=1;i<tcltokens->GetEntriesFast();i=i+2) {
902     level=((TObjString*)tcltokens->At(i+1))->String()[0];
903     input=atoi((((TObjString*)tcltokens->At(i+1))->String()).Remove(0));
904     Bool_t inpnext = IsInputFired(aEsd,level,input);
905     Char_t op =((TObjString*)tcltokens->At(i))->String()[0];
906     if (op == '&') tcl=tcl && inpnext;
907     else if (op == '|') tcl =tcl || inpnext;
908     else {
909        AliError(Form("Syntax error in %s", tclass));
910        delete tcltokens;
911        tcltokens = 0;
912        //      tcltokens->Delete();
913        return kFALSE;
914     }
915   }
916   delete tcltokens;
917   tcltokens = 0;
918        //  tcltokens->Delete();
919   return tcl;
920 }
921
922 Bool_t AliTriggerAnalysis::IsInputFired(const AliESDEvent* aEsd, Char_t level, UInt_t input) const
923 {
924   // Checks trigger input of any level
925   
926   switch (level)
927   {
928     case '0': return IsL0InputFired(aEsd,input);
929     case '1': return IsL1InputFired(aEsd,input);
930     case '2': return IsL2InputFired(aEsd,input);
931     default:
932       AliError(Form("Wrong level %i",level));
933       return kFALSE;
934   }
935 }
936
937 Bool_t AliTriggerAnalysis::IsL0InputFired(const AliESDEvent* aEsd, UInt_t input) const 
938 {
939   // Checks if corresponding bit in mask is on
940   
941   UInt_t inpmask = aEsd->GetHeader()->GetL0TriggerInputs();
942   return (inpmask & (1<<(input-1)));
943 }
944
945 Bool_t AliTriggerAnalysis::IsL1InputFired(const AliESDEvent* aEsd, UInt_t input) const
946 {
947   // Checks if corresponding bit in mask is on
948   
949   UInt_t inpmask = aEsd->GetHeader()->GetL1TriggerInputs();
950   return (inpmask & (1<<(input-1)));
951 }
952
953 Bool_t AliTriggerAnalysis::IsL2InputFired(const AliESDEvent* aEsd, UInt_t input) const 
954 {
955   // Checks if corresponding bit in mask is on
956   
957   UInt_t inpmask = aEsd->GetHeader()->GetL2TriggerInputs();
958   return (inpmask & (1<<(input-1)));
959 }
960
961 void AliTriggerAnalysis::FillHistograms(const AliESDEvent* aEsd) 
962 {
963   // fills the histograms with the info from the ESD
964   
965   fHistBitsSPD->Fill(SPDFiredChips(aEsd, 0), SPDFiredChips(aEsd, 1, kTRUE));
966   
967   V0Trigger(aEsd, kASide, kFALSE, kTRUE);
968   V0Trigger(aEsd, kCSide, kFALSE, kTRUE);
969   T0Trigger(aEsd, kFALSE, kTRUE);
970   ZDCTDCTrigger(aEsd,kASide,kFALSE,kFALSE,kTRUE);
971   ZDCTimeTrigger(aEsd,kTRUE);
972   IsSPDClusterVsTrackletBG(aEsd, kTRUE);
973
974   AliESDZDC* zdcData = aEsd->GetESDZDC();
975   if (zdcData)
976   {
977     UInt_t quality = zdcData->GetESDQuality();
978     
979     // from Nora's presentation, general first physics meeting 16.10.09
980     static UInt_t zpc  = 0x20;
981     static UInt_t znc  = 0x10;
982     static UInt_t zem1 = 0x08;
983     static UInt_t zem2 = 0x04;
984     static UInt_t zpa  = 0x02;
985     static UInt_t zna  = 0x01;
986    
987     fHistZDC->Fill(1, (quality & zna)  ? 1 : 0);
988     fHistZDC->Fill(2, (quality & zpa)  ? 1 : 0);
989     fHistZDC->Fill(3, (quality & zem2) ? 1 : 0);
990     fHistZDC->Fill(4, (quality & zem1) ? 1 : 0);
991     fHistZDC->Fill(5, (quality & znc)  ? 1 : 0);
992     fHistZDC->Fill(6, (quality & zpc)  ? 1 : 0);
993   }
994   else
995   {
996     fHistZDC->Fill(-1);
997     AliError("AliESDZDC not available");
998   }
999   
1000   if (fDoFMD) {
1001     fHistFMDA->Fill(FMDHitCombinations(aEsd, kASide, kTRUE));
1002     fHistFMDC->Fill(FMDHitCombinations(aEsd, kCSide, kTRUE));
1003   }
1004 }
1005   
1006 void AliTriggerAnalysis::FillTriggerClasses(const AliESDEvent* aEsd)
1007 {
1008   // fills trigger classes map
1009   
1010   TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(aEsd->GetFiredTriggerClasses().Data()));
1011   if (!count)
1012   {
1013     count = new TParameter<Long64_t>(aEsd->GetFiredTriggerClasses(), 0);
1014     fTriggerClasses->Add(new TObjString(aEsd->GetFiredTriggerClasses().Data()), count);
1015   }
1016   count->SetVal(count->GetVal() + 1);
1017 }
1018
1019 Int_t AliTriggerAnalysis::SSDClusters(const AliESDEvent* aEsd)
1020 {
1021   // returns the number of clusters in the SSD
1022   const AliMultiplicity* mult = aEsd->GetMultiplicity();
1023   Int_t clusters = mult->GetNumberOfITSClusters(4)+mult->GetNumberOfITSClusters(5);
1024   return clusters;
1025 }
1026
1027
1028 Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHists, Int_t layer)
1029 {
1030   // returns the number of fired chips in the SPD
1031   //
1032   // origin = 0 --> aEsd->GetMultiplicity()->GetNumberOfFiredChips() (filled from clusters)
1033   // origin = 1 --> aEsd->GetMultiplicity()->TestFastOrFiredChips() (from hardware bits)
1034   // layer  = 0 --> both layers
1035   // layer  = 1 --> inner
1036   // layer  = 2 --> outer
1037   
1038   const AliMultiplicity* mult = aEsd->GetMultiplicity();
1039   if (!mult)
1040   {
1041     AliError("AliMultiplicity not available");
1042     return -1;
1043   }
1044   
1045   if (origin == 0){
1046     if (layer == 0) 
1047       return mult->GetNumberOfFiredChips(0) + mult->GetNumberOfFiredChips(1);
1048
1049     return mult->GetNumberOfFiredChips(layer-1); 
1050   }
1051     
1052   if (origin == 1)
1053   {
1054     Int_t nChips = 0;
1055     Int_t firstChip = 0;
1056     Int_t lastChip  = 1200;
1057     if(layer == 1)
1058       lastChip  = 400;
1059     if(layer == 2)
1060       firstChip = 400;
1061
1062     for (Int_t i=firstChip; i<lastChip; i++)
1063     {
1064       if (mult->TestFastOrFiredChips(i) == kTRUE)
1065       {
1066         // efficiency simulation (if enabled)
1067         if (fSPDGFOEfficiency)
1068         {
1069           if (gRandom->Uniform() > fSPDGFOEfficiency->GetBinContent(i+1))
1070             continue;
1071         }
1072         
1073         nChips++;
1074         if (fillHists)
1075           fHistFiredBitsSPD->Fill(i);
1076       }
1077     }
1078     return nChips;
1079   }
1080   
1081   return -1;
1082 }
1083
1084 Bool_t AliTriggerAnalysis::SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin)
1085 {
1086   // Returns if the SPD gave a global Fast OR trigger
1087   
1088   Int_t firedChips = SPDFiredChips(aEsd, origin);
1089   
1090   if (firedChips >= fSPDGFOThreshold)
1091     return kTRUE;
1092   return kFALSE;
1093 }
1094
1095 Bool_t AliTriggerAnalysis::IsSPDClusterVsTrackletBG(const AliESDEvent* aEsd, Bool_t fillHists){
1096   //rejects BG based on the cluster vs tracklet correlation
1097   // returns true if the event is BG
1098   const AliMultiplicity* mult = aEsd->GetMultiplicity();
1099   if (!mult){
1100     AliFatal("No multiplicity object"); // TODO: Should this be fatal?
1101   }
1102   Int_t ntracklet = mult->GetNumberOfTracklets();
1103
1104   Int_t spdClusters = 0;
1105   for(Int_t ilayer = 0; ilayer < 2; ilayer++){
1106     spdClusters += mult->GetNumberOfITSClusters(ilayer);
1107   }
1108
1109   if(fillHists) {
1110     fHistSPDClsVsTrk->Fill(ntracklet,spdClusters);
1111   }
1112
1113   Bool_t isCvsTOk = kFALSE;
1114   Float_t limit = Float_t(fASPDCvsTCut) + Float_t(ntracklet) * fBSPDCvsTCut;  
1115   if (spdClusters > limit)        isCvsTOk = kTRUE;
1116   else                            isCvsTOk = kFALSE ;
1117
1118   return isCvsTOk;
1119
1120 }
1121   
1122
1123 AliTriggerAnalysis::V0Decision AliTriggerAnalysis::V0Trigger(const AliESDEvent* aEsd, AliceSide side, Bool_t online, Bool_t fillHists)
1124 {
1125   // Returns the V0 trigger decision in V0A | V0C
1126   //
1127   // Returns kV0Fake if the calculated average time is in a window where neither BB nor BG is expected. 
1128   // The rate of such triggers can be used to estimate the background. Note that the rate has to be 
1129   // rescaled with the size of the windows (numerical values see below in the code)
1130   //
1131   // argument 'online' is used as a switch between online and offline trigger algorithms
1132   //
1133   // Based on an algorithm by Cvetan Cheshkov
1134
1135   AliESDVZERO* esdV0 = aEsd->GetVZEROData();
1136   if (!esdV0)
1137   {
1138     AliError("AliESDVZERO not available");
1139     return kV0Invalid;
1140   }
1141   AliDebug(2,Form("In V0Trigger: %f %f",esdV0->GetV0ATime(),esdV0->GetV0CTime()));
1142
1143   Int_t begin = -1;
1144   Int_t end = -1;
1145   
1146   if (side == kASide)
1147   {
1148     begin = 32;
1149     end = 64;
1150   } 
1151   else if (side == kCSide)
1152   {
1153     begin = 0;
1154     end = 32;
1155   }
1156   else
1157     return kV0Invalid;
1158     
1159    if (esdV0->TestBit(AliESDVZERO::kDecisionFilled)) {
1160     if (online) {
1161       if (esdV0->TestBit(AliESDVZERO::kOnlineBitsFilled)) {
1162         for (Int_t i = begin; i < end; ++i) {
1163           if (esdV0->GetBBFlag(i)) return kV0BB;
1164         }
1165         for (Int_t i = begin; i < end; ++i) {
1166           if (esdV0->GetBGFlag(i)) return kV0BG;
1167         }
1168         return kV0Empty;
1169       }
1170       else {
1171         AliWarning("V0 online trigger analysis is not yet available!");
1172         return kV0BB;
1173       }
1174     }
1175     else {
1176
1177       if (fillHists) {
1178         if (side == kASide && fHistV0A)
1179           fHistV0A->Fill(esdV0->GetV0ATime());
1180         if (side == kCSide && fHistV0C)
1181           fHistV0C->Fill(esdV0->GetV0CTime());
1182       }
1183
1184       if (side == kASide) return (V0Decision)esdV0->GetV0ADecision();
1185       else if (side == kCSide) return (V0Decision)esdV0->GetV0CDecision();
1186       else return kV0Invalid;
1187     }
1188   }
1189
1190   Float_t time = 0;
1191   Float_t weight = 0;
1192   if (fMC)
1193   {
1194     Int_t runRange;
1195     if (aEsd->GetRunNumber() <= 104803) runRange = 0;
1196     else if (aEsd->GetRunNumber() <= 104876) runRange = 1;
1197     else runRange = 2;
1198
1199     Float_t factors[3][64] = {
1200       // runs: 104792-104803
1201       {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},
1202       // runs: 104841-104876
1203       {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},
1204       // runs: 104890-92
1205       {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}
1206     };
1207     Float_t dA = 77.4 - 11.0;
1208     Float_t dC = 77.4 - 2.9;
1209     // Time misalignment
1210     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};
1211     Float_t dA2 = 2.8, dC2 = 3.3;
1212
1213     if (online) {
1214       for (Int_t i = begin; i < end; ++i) {
1215         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
1216         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
1217         if (esdV0->GetTime(i) >= 1e-6 &&
1218             tempTime > fV0HwWinLow && tempTime < fV0HwWinHigh &&
1219             tempAdc > fV0HwAdcThr)
1220           return kV0BB;
1221       }
1222       return kV0Empty;
1223     }
1224     else {
1225       for (Int_t i = begin; i < end; ++i) {
1226         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
1227         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA : esdV0->GetTime(i)+dC;
1228         Float_t tempRawTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
1229         if (esdV0->GetTime(i) >= 1e-6 &&
1230             tempRawTime < 125.0 &&
1231             tempAdc > fV0AdcThr) {
1232           weight += 1.0;
1233           time += tempTime;
1234         }
1235       }
1236     }
1237   }
1238   else {
1239     if (online) {
1240       for (Int_t i = begin; i < end; ++i) {
1241         if (esdV0->GetTime(i) >= 1e-6 &&
1242             esdV0->GetTime(i) > fV0HwWinLow && esdV0->GetTime(i) < fV0HwWinHigh &&
1243             esdV0->GetAdc(i) > fV0HwAdcThr)
1244           return kV0BB;
1245       }
1246       return kV0Empty;
1247     }
1248     else {
1249       for (Int_t i = begin; i < end; ++i) {
1250         if (esdV0->GetTime(i) > 1e-6 && esdV0->GetAdc(i) > fV0AdcThr) {
1251           Float_t correctedTime = V0CorrectLeadingTime(i, esdV0->GetTime(i), esdV0->GetAdc(i),aEsd->GetRunNumber());
1252           Float_t timeWeight = V0LeadingTimeWeight(esdV0->GetAdc(i));
1253           time += correctedTime*timeWeight;
1254             
1255           weight += timeWeight;
1256         }
1257       }
1258     }
1259   }
1260
1261   if (weight > 0) 
1262     time /= weight;
1263   time += fV0TimeOffset;
1264
1265   if (fillHists)
1266   {
1267     if (side == kASide && fHistV0A)
1268       fHistV0A->Fill(time);
1269     if (side == kCSide && fHistV0C)
1270       fHistV0C->Fill(time);
1271   }
1272   
1273   if (side == kASide)
1274   {
1275     if (time > 68 && time < 100)
1276       return kV0BB;
1277     if (time > 54 && time < 57.5) 
1278       return kV0BG;
1279     if (time > 57.5 && time < 68)
1280       return kV0Fake;
1281   }
1282   
1283   if (side == kCSide)
1284   {
1285     if (time > 75.5 && time < 100)
1286       return kV0BB;
1287     if (time > 69.5 && time < 73)
1288       return kV0BG; 
1289     if (time > 55 && time < 69.5)
1290       return kV0Fake;
1291   }
1292   
1293   return kV0Empty;
1294 }
1295
1296 Float_t AliTriggerAnalysis::V0CorrectLeadingTime(Int_t i, Float_t time, Float_t adc, Int_t runNumber) const
1297 {
1298   // Correct for slewing and align the channels
1299   //
1300   // Authors: Cvetan Cheshkov / Raphael Tieulent
1301
1302   if (time == 0) return 0;
1303
1304   // Time alignment
1305   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};
1306
1307   if(runNumber < 106031)
1308     time -= timeShift[i];
1309
1310   // Slewing correction
1311   if (adc == 0) return time;
1312
1313   Float_t p1 = 1.57345e1;
1314   Float_t p2 =-4.25603e-1;
1315
1316   if(runNumber >= 106031) adc *= (2.5/4.0);
1317   return (time - p1*TMath::Power(adc,p2));
1318 }
1319
1320 Float_t AliTriggerAnalysis::V0LeadingTimeWeight(Float_t adc) const
1321 {
1322   if (adc < 1e-6) return 0;
1323
1324   Float_t p1 = 40.211;
1325   Float_t p2 =-4.25603e-1;
1326   Float_t p3 = 0.5646;
1327
1328   return 1./(p1*p1*TMath::Power(adc,2.*(p2-1.))+p3*p3);
1329 }
1330
1331
1332 Bool_t AliTriggerAnalysis::ZDCTDCTrigger(const AliESDEvent* aEsd, AliceSide side, Bool_t useZN, Bool_t useZP, Bool_t fillHists) const
1333 {
1334   // Returns if ZDC triggered, based on TDC information 
1335   
1336   AliESDZDC *esdZDC = aEsd->GetESDZDC();
1337
1338   Bool_t zdcNA = kFALSE;
1339   Bool_t zdcNC = kFALSE;
1340   Bool_t zdcPA = kFALSE;
1341   Bool_t zdcPC = kFALSE;
1342
1343   if (fMC) {
1344         // If it's MC, we use the energy
1345     Double_t minEnergy = 0;
1346     Double_t  zNCEnergy = esdZDC->GetZDCN1Energy();
1347     Double_t  zPCEnergy = esdZDC->GetZDCP1Energy();
1348     Double_t  zNAEnergy = esdZDC->GetZDCN2Energy();
1349     Double_t  zPAEnergy = esdZDC->GetZDCP2Energy();
1350     zdcNA = (zNAEnergy>minEnergy);
1351     zdcNC = (zNCEnergy>minEnergy);
1352     zdcPA = (zPAEnergy>minEnergy);
1353     zdcPC = (zPCEnergy>minEnergy);
1354
1355   }
1356   else {
1357
1358     Bool_t tdc[32] = {kFALSE};
1359     for(Int_t itdc=0; itdc<32; itdc++){
1360       for(Int_t i=0; i<4; i++){
1361         if (esdZDC->GetZDCTDCData(itdc, i) != 0){
1362           tdc[itdc] = kTRUE;
1363         }
1364       }
1365       if(fillHists && tdc[itdc]) {
1366         fHistTDCZDC->Fill(itdc);
1367       }    
1368     }
1369     zdcNA = tdc[12];
1370     zdcNC = tdc[10];
1371     zdcPA = tdc[13];
1372     zdcPC = tdc[11];
1373   }
1374
1375   if (side == kASide) return ((useZP&&zdcPA) || (useZN&&zdcNA)); 
1376   if (side == kCSide) return ((useZP&&zdcPC) || (useZN&&zdcNC)); 
1377   return kFALSE;
1378 }
1379
1380 Bool_t AliTriggerAnalysis::ZDCTimeTrigger(const AliESDEvent *aEsd, Bool_t fillHists) const
1381 {
1382   // This method implements a selection
1383   // based on the timing in both sides of zdcN
1384   // It can be used in order to eliminate
1385   // parasitic collisions
1386   Bool_t zdcAccept = kFALSE;
1387   AliESDZDC *esdZDC = aEsd->GetESDZDC();
1388
1389   if(fMC) {
1390      UInt_t esdFlag =  esdZDC->GetESDQuality();
1391
1392      Bool_t znaFired=kFALSE, zpaFired=kFALSE;
1393      Bool_t zem1Fired=kFALSE, zem2Fired=kFALSE;
1394      Bool_t zncFired=kFALSE, zpcFired=kFALSE;
1395      //
1396      // **** Trigger patterns
1397      if((esdFlag & 0x00000001) == 0x00000001) znaFired=kTRUE;
1398      if((esdFlag & 0x00000002) == 0x00000002) zpaFired=kTRUE;
1399      if((esdFlag & 0x00000004) == 0x00000004) zem1Fired=kTRUE;
1400      if((esdFlag & 0x00000008) == 0x00000008) zem2Fired=kTRUE;
1401      if((esdFlag & 0x00000010) == 0x00000010) zncFired=kTRUE;
1402      if((esdFlag & 0x00000020) == 0x00000020) zpcFired=kTRUE;
1403      zdcAccept = (znaFired | zncFired);
1404   }
1405   else {
1406     for(Int_t i = 0; i < 4; ++i) {
1407       if (esdZDC->GetZDCTDCData(10,i) != 0) {
1408         Float_t tdcC = 0.025*(esdZDC->GetZDCTDCData(10,i)-esdZDC->GetZDCTDCData(14,i)); 
1409         Float_t tdcCcorr = esdZDC->GetZDCTDCCorrected(10,i); 
1410         for(Int_t j = 0; j < 4; ++j) {
1411           if (esdZDC->GetZDCTDCData(12,j) != 0) {
1412             Float_t tdcA = 0.025*(esdZDC->GetZDCTDCData(12,j)-esdZDC->GetZDCTDCData(14,j));
1413
1414             Float_t tdcAcorr = esdZDC->GetZDCTDCCorrected(12,j);
1415             if(fillHists) {
1416               fHistTimeZDC->Fill(tdcC-tdcA,tdcC+tdcA);
1417               fHistTimeCorrZDC->Fill(tdcCcorr-tdcAcorr,tdcCcorr+tdcAcorr);
1418             }
1419             if (esdZDC->TestBit(AliESDZDC::kCorrectedTDCFilled)) {
1420               if (((tdcCcorr-tdcAcorr-fZDCCutRefDeltaCorr)*(tdcCcorr-tdcAcorr-fZDCCutRefDeltaCorr)/(fZDCCutSigmaDeltaCorr*fZDCCutSigmaDeltaCorr) +
1421                    (tdcCcorr+tdcAcorr-fZDCCutRefSumCorr)*(tdcCcorr+tdcAcorr-fZDCCutRefSumCorr)/(fZDCCutSigmaSumCorr*fZDCCutSigmaSumCorr))< 1.0)
1422                 zdcAccept = kTRUE;
1423             }
1424             else {
1425               if (((tdcC-tdcA-fZDCCutRefDelta)*(tdcC-tdcA-fZDCCutRefDelta)/(fZDCCutSigmaDelta*fZDCCutSigmaDelta) +
1426                    (tdcC+tdcA-fZDCCutRefSum)*(tdcC+tdcA-fZDCCutRefSum)/(fZDCCutSigmaSum*fZDCCutSigmaSum))< 1.0)
1427                 zdcAccept = kTRUE;
1428             }
1429           }
1430         }
1431       }
1432     }
1433   }
1434   return zdcAccept;
1435 }
1436
1437 Bool_t AliTriggerAnalysis::ZDCTimeBGTrigger(const AliESDEvent *aEsd, AliceSide side) const
1438 {
1439   // This method implements a selection
1440   // based on the timing in of zdcN
1441   // It can be used in order to flag background
1442   // ** So far only implemented for the 2012 pA run **
1443
1444   if(fMC) return kFALSE;
1445
1446   AliESDZDC *zdcData = aEsd->GetESDZDC();
1447   Bool_t zna = kFALSE;
1448   Bool_t znc = kFALSE;
1449
1450   Float_t tdcC=999, tdcCcorr=999, tdcA=999, tdcAcorr=999;
1451   for(Int_t i = 0; i < 4; ++i) {
1452     if (zdcData->GetZDCTDCData(10,i) != 0) {
1453       znc = kTRUE;
1454       tdcC = 0.025*(zdcData->GetZDCTDCData(10,i)-zdcData->GetZDCTDCData(14,i));
1455       tdcCcorr = zdcData->GetZDCTDCCorrected(10,i);
1456     }
1457   }
1458   for(Int_t j = 0; j < 4; ++j) {
1459     if (zdcData->GetZDCTDCData(12,j) != 0) {
1460       zna = kTRUE;
1461       tdcA = 0.025*(zdcData->GetZDCTDCData(12,j)-zdcData->GetZDCTDCData(14,j));
1462       tdcAcorr = zdcData->GetZDCTDCCorrected(12,j);
1463     }
1464   }
1465
1466   const Int_t runNumber = aEsd->GetRunNumber();
1467   if(runNumber<188124 || (runNumber>188374 && runNumber<194713)){ // FIXME: end of pA-run is not known
1468     AliError(Form(" ZN BG time cut not implemented for run %d",runNumber));
1469     return kFALSE;
1470   }
1471
1472   Bool_t znabg = (zna && (TMath::Abs(tdcAcorr)>fZDCCutZNATimeCorr));
1473   Bool_t zncbg = (znc && (TMath::Abs(tdcCcorr)>fZDCCutZNCTimeCorr));
1474
1475   // Printf("Checking ZN background (time) for run %d, A = %d, time=%2.2f, C = %d, time=%2.2f",runNumber,(Int_t)zna,tdcAcorr,(Int_t)znc,tdcCcorr);
1476   // Printf("Checking ZN background (time) for run %d, A-BG = %d, C-BG = %d",runNumber,(Int_t)znabg,(Int_t)zncbg);
1477
1478   if (side == kASide) return znabg;
1479   if (side == kCSide) return zncbg;
1480   return kFALSE;
1481 }
1482
1483 Bool_t AliTriggerAnalysis::ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) const
1484 {
1485   // Returns if ZDC triggered
1486   
1487   AliESDZDC* zdcData = aEsd->GetESDZDC();
1488   if (!zdcData)
1489   {
1490     AliError("AliESDZDC not available");
1491     return kFALSE;
1492   }
1493   
1494   UInt_t quality = zdcData->GetESDQuality();
1495   
1496   // from Nora's presentation, general first physics meeting 16.10.09
1497   static UInt_t zpc  = 0x20;
1498   static UInt_t znc  = 0x10;
1499   static UInt_t zem1 = 0x08;
1500   static UInt_t zem2 = 0x04;
1501   static UInt_t zpa  = 0x02;
1502   static UInt_t zna  = 0x01;
1503   
1504   if (side == kASide && ((quality & zpa) || (quality & zna)))
1505     return kTRUE;
1506   if (side == kCentralBarrel && ((quality & zem1) || (quality & zem2)))
1507     return kTRUE;
1508   if (side == kCSide && ((quality & zpc) || (quality & znc)))
1509     return kTRUE;
1510   
1511   return kFALSE;
1512 }
1513
1514 Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists)
1515 {
1516   // returns number of hit combinations agove threshold
1517   //
1518   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
1519
1520   if (!fDoFMD)
1521     return -1;
1522
1523   // Workaround for AliESDEvent::GetFMDData is not const!
1524   const AliESDFMD* fmdData = (const_cast<AliESDEvent*>(aEsd))->GetFMDData();
1525   if (!fmdData)
1526   {
1527     AliError("AliESDFMD not available");
1528     return -1;
1529   }
1530
1531   Int_t detFrom = (side == kASide) ? 1 : 3;
1532   Int_t detTo   = (side == kASide) ? 2 : 3;
1533
1534   Int_t triggers = 0;
1535   Float_t totalMult = 0;
1536   for (UShort_t det=detFrom;det<=detTo;det++) {
1537     Int_t nRings = (det == 1 ? 1 : 2);
1538     for (UShort_t ir = 0; ir < nRings; ir++) {
1539       Char_t   ring = (ir == 0 ? 'I' : 'O');
1540       UShort_t nsec = (ir == 0 ? 20  : 40);
1541       UShort_t nstr = (ir == 0 ? 512 : 256);
1542       for (UShort_t sec =0; sec < nsec;  sec++) {
1543         for (UShort_t strip = 0; strip < nstr; strip++) {
1544           Float_t mult = fmdData->Multiplicity(det,ring,sec,strip);
1545           if (mult == AliESDFMD::kInvalidMult) continue;
1546           
1547           if (fillHists)
1548             fHistFMDSingle->Fill(mult);
1549           
1550           if (mult > fFMDLowCut)
1551             totalMult = totalMult + mult;
1552           else
1553           {
1554             if (totalMult > fFMDHitCut)
1555               triggers++;
1556               
1557             if (fillHists)
1558               fHistFMDSum->Fill(totalMult);
1559               
1560             totalMult = 0;
1561           }
1562         }
1563       }
1564     }
1565   }
1566   
1567   return triggers;
1568 }
1569
1570 Bool_t AliTriggerAnalysis::FMDTrigger(const AliESDEvent* aEsd, AliceSide side)
1571 {
1572   // Returns if the FMD triggered
1573   //
1574   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
1575
1576   Int_t triggers = FMDHitCombinations(aEsd, side, kFALSE);
1577   
1578   if (triggers > 0)
1579     return kTRUE;
1580     
1581   return kFALSE;
1582 }
1583
1584 Long64_t AliTriggerAnalysis::Merge(TCollection* list)
1585 {
1586   // Merge a list of AliMultiplicityCorrection objects with this (needed for
1587   // PROOF).
1588   // Returns the number of merged objects (including this).
1589
1590   if (!list)
1591     return 0;
1592
1593   if (list->IsEmpty())
1594     return 1;
1595
1596   TIterator* iter = list->MakeIterator();
1597   TObject* obj;
1598
1599   // collections of all histograms
1600   const Int_t nHists = 14;
1601   TList collections[nHists];
1602
1603   Int_t count = 0;
1604   while ((obj = iter->Next())) {
1605
1606     AliTriggerAnalysis* entry = dynamic_cast<AliTriggerAnalysis*> (obj);
1607     if (entry == 0) 
1608       continue;
1609
1610     Int_t n = 0;
1611     collections[n++].Add(entry->fHistV0A);
1612     collections[n++].Add(entry->fHistV0C);
1613     collections[n++].Add(entry->fHistZDC);
1614     collections[n++].Add(entry->fHistTDCZDC);
1615     collections[n++].Add(entry->fHistTimeZDC);
1616     collections[n++].Add(entry->fHistTimeCorrZDC);
1617     collections[n++].Add(entry->fHistFMDA);
1618     collections[n++].Add(entry->fHistFMDC);
1619     collections[n++].Add(entry->fHistFMDSingle);
1620     collections[n++].Add(entry->fHistFMDSum);
1621     collections[n++].Add(entry->fHistBitsSPD);
1622     collections[n++].Add(entry->fHistFiredBitsSPD);
1623     collections[n++].Add(entry->fHistSPDClsVsTrk);
1624     collections[n++].Add(entry->fHistT0);
1625
1626     // merge fTriggerClasses
1627     TIterator* iter2 = entry->fTriggerClasses->MakeIterator();
1628     TObjString* obj2 = 0;
1629     while ((obj2 = dynamic_cast<TObjString*> (iter2->Next())))
1630     {
1631       TParameter<Long64_t>* param2 = static_cast<TParameter<Long64_t>*> (entry->fTriggerClasses->GetValue(obj2));
1632       
1633       TParameter<Long64_t>* param1 = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj2));
1634       if (param1)
1635       {
1636         param1->SetVal(param1->GetVal() + param2->GetVal());
1637       }
1638       else
1639       {
1640         param1 = dynamic_cast<TParameter<Long64_t>*> (param2->Clone());
1641         fTriggerClasses->Add(new TObjString(obj2->String()), param1);
1642       }
1643     }
1644     
1645     delete iter2;
1646   
1647     count++;
1648   }
1649
1650   Int_t n = 0;
1651   fHistV0A->Merge(&collections[n++]);
1652   fHistV0C->Merge(&collections[n++]);
1653   fHistZDC->Merge(&collections[n++]);
1654   fHistTDCZDC->Merge(&collections[n++]);
1655   if (fHistTimeZDC)
1656     fHistTimeZDC->Merge(&collections[n++]);
1657   else
1658     n++;
1659   if (fHistTimeCorrZDC)
1660     fHistTimeCorrZDC->Merge(&collections[n++]);
1661   else
1662     n++;
1663   fHistFMDA->Merge(&collections[n++]);
1664   fHistFMDC->Merge(&collections[n++]);
1665   fHistFMDSingle->Merge(&collections[n++]);
1666   fHistFMDSum->Merge(&collections[n++]);
1667   fHistBitsSPD->Merge(&collections[n++]);
1668   fHistFiredBitsSPD->Merge(&collections[n++]);
1669   fHistSPDClsVsTrk->Merge(&collections[n++]);
1670   fHistT0->Merge(&collections[n++]);
1671   delete iter;
1672
1673   return count+1;
1674 }
1675
1676 void AliTriggerAnalysis::SaveHistograms() const
1677 {
1678   // write histograms to current directory
1679   
1680   if (!fHistBitsSPD)
1681     return;
1682     
1683   if (fHistBitsSPD) {
1684     fHistBitsSPD->Write();
1685     fHistBitsSPD->ProjectionX();
1686     fHistBitsSPD->ProjectionY();
1687   }
1688   else Printf("Cannot save fHistBitsSPD");
1689   if (fHistFiredBitsSPD) fHistFiredBitsSPD->Write();
1690   else Printf("Cannot save fHistFiredBitsSPD");
1691   if (fHistV0A) fHistV0A->Write();
1692   else Printf("Cannot save fHistV0A");
1693   if (fHistV0C) fHistV0C->Write();
1694   else Printf("Cannot save fHistV0C");
1695   if (fHistZDC) fHistZDC->Write();
1696   else Printf("Cannot save fHistZDC");
1697   if (fHistTDCZDC) fHistTDCZDC->Write();
1698   else Printf("Cannot save fHistTDCZDC");
1699   if (fHistTimeZDC) fHistTimeZDC->Write();
1700   else Printf("Cannot save fHistTimeZDC");
1701   if (fHistTimeCorrZDC) fHistTimeCorrZDC->Write();
1702   else Printf("Cannot save fHistTimeCorrZDC");
1703   if (fHistFMDA) fHistFMDA->Write();
1704   else Printf("Cannot save fHistFMDA");
1705   if (fHistFMDC) fHistFMDC->Write();
1706   else Printf("Cannot save fHistFMDC");
1707   if (fHistFMDSingle) fHistFMDSingle->Write();
1708   else Printf("Cannot save fHistFMDSingle");
1709   if (fHistFMDSum) fHistFMDSum->Write();
1710   else Printf("Cannot save fHistFMDSum");
1711   if (fSPDGFOEfficiency) fSPDGFOEfficiency->Write("fSPDGFOEfficiency");
1712   if (fHistSPDClsVsTrk) fHistSPDClsVsTrk->Write("fHistSPDClsVsTrk");
1713   if (fHistT0) fHistT0->Write("fHistT0");
1714
1715   //  else Printf("Cannot save fSPDGFOEfficiency");
1716   
1717   fTriggerClasses->Write("fTriggerClasses", TObject::kSingleKey);
1718 }
1719
1720 void AliTriggerAnalysis::PrintTriggerClasses() const
1721 {
1722   // print trigger classes
1723   
1724   Printf("Trigger Classes:");
1725   
1726   Printf("Event count for trigger combinations:");
1727   
1728   TMap singleTrigger;
1729   singleTrigger.SetOwner();
1730   
1731   TIterator* iter = fTriggerClasses->MakeIterator();
1732   TObjString* obj = 0;
1733   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
1734   {
1735     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj));
1736     
1737     Printf(" %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1738     
1739     TObjArray* tokens = obj->String().Tokenize(" ");
1740     for (Int_t i=0; i<tokens->GetEntries(); i++)
1741     {
1742       TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(((TObjString*) tokens->At(i))->String().Data()));
1743       if (!count)
1744       {
1745         count = new TParameter<Long64_t>(((TObjString*) tokens->At(i))->String().Data(), 0);
1746         singleTrigger.Add(new TObjString(((TObjString*) tokens->At(i))->String().Data()), count);
1747       }
1748       count->SetVal(count->GetVal() + param->GetVal());
1749     }
1750     
1751     delete tokens;
1752   }
1753   delete iter;
1754   
1755   Printf("Event count for single trigger:");
1756   
1757   iter = singleTrigger.MakeIterator();
1758   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
1759   {
1760     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(obj));
1761     
1762     Printf("  %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1763   }
1764   delete iter;
1765   
1766   singleTrigger.DeleteAll();
1767 }
1768
1769
1770 //----------------------------------------------------------------------------------------------------
1771 AliTriggerAnalysis::T0Decision AliTriggerAnalysis::T0Trigger(const AliESDEvent* aEsd, Bool_t online, Bool_t fillHists) 
1772 {
1773   // Returns the T0 TVDC trigger decision
1774   //  
1775   // argument 'online' is used as a switch between online and offline trigger algorithms
1776   // in online mode return 0TVX 
1777   // in offline mode in addition check pile-up and background :
1778   // pile-up readed from ESD: check if TVDC (0TVX module name) has more 1 hit;
1779   // backgroud flag readed from ESD : check in given time interval OrA and OrC were correct but TVDC not  
1780   // 
1781   // Based on an algorithm by Alla Maevskaya 
1782
1783   const AliESDTZERO* esdT0 = aEsd->GetESDTZERO();
1784   if (!esdT0)
1785   {
1786     AliError("AliESDTZERO not available");
1787     return kT0Invalid;
1788   }
1789   //????  AliDebug(2,Form("In T0Trigger: %f %f",esdV0->GetV0ATime(),esdV0->GetV0CTime()));
1790   Float_t  tvdc[5] ;
1791   for (Int_t ii=0; ii<5; ii++)
1792     tvdc[ii] = esdT0->GetTVDC(ii);
1793   //  Int_t trig=esdT0->GetT0Trig();
1794   //  cout<<" T0 trig "<<trig<<endl;    
1795
1796   if(fillHists) fHistT0->Fill(tvdc[0]);
1797     
1798   if (online) {
1799     if(aEsd->GetHeader()->GetFiredTriggerInputs().Contains("0TVX") ) return kT0BB;
1800   }
1801   else {
1802   
1803     if (esdT0->GetPileupFlag()) return kT0DecPileup;
1804     if (esdT0->GetBackgroundFlag()) return kT0DecBG;
1805     if (tvdc[0]>-5 && tvdc[0]<5 && tvdc[0] != 0) return kT0BB; 
1806   }
1807
1808   if (fMC)
1809     if( esdT0->GetT0zVertex()>-12.3 && esdT0->GetT0zVertex() < 10.3) return kT0BB; 
1810  
1811   return kT0Empty;
1812 }
1813
1814 //----------------------------------------------------------------------------------------------------
1815 Bool_t AliTriggerAnalysis::EMCALCellsTrigger(const AliESDEvent *aEsd)
1816 {
1817   //
1818   // Returns the EMCAL trigger decision
1819   // so far only implemented for LHC11a data
1820   // see http://alisoft.cern.ch/viewvc/trunk/PWGGA/EMCALTasks/AliEmcalPhysicsSelection.cxx?view=markup&root=AliRoot Revision 56136
1821   //
1822
1823   Bool_t isFired = kTRUE;
1824   const Int_t runNumber = aEsd->GetRunNumber();
1825
1826   /*
1827   // Load EMCAL branches from the manager
1828   AliAnalysisManager *am = AliAnalysisManager::GetAnalysisManager();
1829   am->LoadBranch("EMCALCells.");
1830   am->LoadBranch("CaloClusters");
1831   */
1832
1833   // Get EMCAL cells
1834   AliVCaloCells *cells = aEsd->GetEMCALCells();
1835   const Short_t nCells = cells->GetNumberOfCells();
1836
1837   // count cells above threshold per sm
1838   Int_t nCellCount[10] = {0,0,0,0,0,0,0,0,0,0};
1839   for(Int_t iCell=0; iCell<nCells; ++iCell) {
1840     Short_t cellId = cells->GetCellNumber(iCell);
1841     Double_t cellE = cells->GetCellAmplitude(cellId);
1842     Int_t sm       = cellId / (24*48);
1843     if (cellE>0.1)
1844       ++nCellCount[sm];
1845   }
1846
1847   // Trigger decision for LHC11a
1848   Bool_t isLedEvent = kFALSE;
1849   if ((runNumber>=144871) && (runNumber<=146860)) {
1850     if (nCellCount[4] > 100)
1851       isLedEvent = kTRUE;
1852     else {
1853       if ((runNumber>=146858) && (runNumber<=146860)) {
1854         if (nCellCount[3]>=35)
1855           isLedEvent = kTRUE;
1856       }
1857     }
1858   }
1859   
1860   if (isLedEvent) {
1861     isFired = kFALSE;
1862   }
1863
1864   return isFired;
1865 }