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