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