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