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