]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliTriggerAnalysis.cxx
Removing the dependence on libESD (removed AliESDCaloClusters (Gustavo) see savannah...
[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
46 ClassImp(AliTriggerAnalysis)
47
48 AliTriggerAnalysis::AliTriggerAnalysis() :
49   fSPDGFOThreshold(2),
50   fSPDGFOEfficiency(0),
51   fV0TimeOffset(0),
52   fV0AdcThr(0),
53   fV0HwAdcThr(2.5),
54   fV0HwWinLow(61.5),
55   fV0HwWinHigh(86.5),
56   fFMDLowCut(0.2),
57   fFMDHitCut(0.5),
58   fHistBitsSPD(0),
59   fHistFiredBitsSPD(0),
60   fHistV0A(0),       
61   fHistV0C(0),
62   fHistZDC(0),    
63   fHistFMDA(0),    
64   fHistFMDC(0),   
65   fHistFMDSingle(0),
66   fHistFMDSum(0),
67   fTriggerClasses(0),
68   fMC(kFALSE)
69 {
70   // constructor
71 }
72
73 AliTriggerAnalysis::~AliTriggerAnalysis()
74 {
75   // destructor
76   
77   if (fHistBitsSPD)
78   {
79     delete fHistBitsSPD;
80     fHistBitsSPD = 0;
81   }
82
83   if (fHistFiredBitsSPD)
84   {
85     delete fHistFiredBitsSPD;
86     fHistFiredBitsSPD = 0;
87   }
88
89   if (fHistV0A)
90   {
91     delete fHistV0A;
92     fHistV0A = 0;
93   }
94
95   if (fHistV0C)
96   {
97     delete fHistV0C;
98     fHistV0C = 0;
99   }
100
101   if (fHistZDC)
102   {
103     delete fHistZDC;
104     fHistZDC = 0;
105   }
106
107   if (fHistFMDA)
108   {
109     delete fHistFMDA;
110     fHistFMDA = 0;
111   }
112
113   if (fHistFMDC)
114   {
115     delete fHistFMDC;
116     fHistFMDC = 0;
117   }
118
119   if (fHistFMDSingle)
120   {
121     delete fHistFMDSingle;
122     fHistFMDSingle = 0;
123   }
124
125   if (fHistFMDSum)
126   {
127     delete fHistFMDSum;
128     fHistFMDSum = 0;
129   }
130
131   if (fTriggerClasses)
132   {
133     fTriggerClasses->DeleteAll();
134     delete fTriggerClasses;
135     fTriggerClasses = 0;
136   }
137 }
138
139 void AliTriggerAnalysis::EnableHistograms()
140 {
141   // creates the monitoring histograms
142   
143   // do not add this hists to the directory
144   Bool_t oldStatus = TH1::AddDirectoryStatus();
145   TH1::AddDirectory(kFALSE);
146   
147   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);
148   fHistFiredBitsSPD = new TH1F("fHistFiredBitsSPD", "SPD GFO Hardware;chip number;events", 1200, -0.5, 1199.5);
149   fHistV0A = new TH1F("fHistV0A", "V0A;leading time (ns);events", 400, -100, 100);
150   fHistV0C = new TH1F("fHistV0C", "V0C;leading time (ns);events", 400, -100, 100);
151   fHistZDC = new TH1F("fHistZDC", "ZDC;trigger bits;events", 8, -1.5, 6.5);
152   
153   // TODO check limits
154   fHistFMDA = new TH1F("fHistFMDA", "FMDA;combinations above threshold;events", 102, -1.5, 100.5);
155   fHistFMDC = new TH1F("fHistFMDC", "FMDC;combinations above threshold;events", 102, -1.5, 100.5);
156   fHistFMDSingle = new TH1F("fHistFMDSingle", "FMD single;multiplicity value;counts", 1000, 0, 10);
157   fHistFMDSum = new TH1F("fHistFMDSum", "FMD sum;multiplicity value;counts", 1000, 0, 10);
158   
159   fTriggerClasses = new TMap;
160   fTriggerClasses->SetOwner();
161   
162   TH1::AddDirectory(oldStatus);
163 }
164
165 //____________________________________________________________________
166 const char* AliTriggerAnalysis::GetTriggerName(Trigger trigger) 
167 {
168   // returns the name of the requested trigger
169   // the returned string will only be valid until the next call to this function [not thread-safe]
170   
171   static TString str;
172   
173   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
174   
175   switch (triggerNoFlags)
176   {
177     case kAcceptAll : str = "ACCEPT ALL (bypass!)"; break;
178     case kMB1 : str = "MB1"; break;
179     case kMB2 : str = "MB2"; break;
180     case kMB3 : str = "MB3"; break;
181     case kSPDGFO : str = "SPD GFO"; break;
182     case kSPDGFOBits : str = "SPD GFO Bits"; break;
183     case kV0A : str = "V0 A BB"; break;
184     case kV0C : str = "V0 C BB"; break;
185     case kV0OR : str = "V0 OR BB"; break;
186     case kV0AND : str = "V0 AND BB"; break;
187     case kV0ABG : str = "V0 A BG"; break;
188     case kV0CBG : str = "V0 C BG"; break;
189     case kZDC : str = "ZDC"; break;
190     case kZDCA : str = "ZDC A"; break;
191     case kZDCC : str = "ZDC C"; break;
192     case kFMDA : str = "FMD A"; break;
193     case kFMDC : str = "FMD C"; break;
194     case kFPANY : str = "SPD GFO | V0 | ZDC | FMD"; break;
195     case kNSD1 : str = "NSD1"; break;
196     case kMB1Prime: str = "MB1prime"; break;
197     default: str = ""; break;
198   }
199    
200   if (trigger & kOfflineFlag)
201     str += " OFFLINE";  
202   
203   if (trigger & kOneParticle)
204     str += " OneParticle";  
205
206   return str;
207 }
208
209 Bool_t AliTriggerAnalysis::IsTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
210 {
211   // checks if an event has been triggered
212
213   if (trigger & kOfflineFlag)
214     return IsOfflineTriggerFired(aEsd, trigger);
215     
216   return IsTriggerBitFired(aEsd, trigger);
217 }
218
219 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, Trigger trigger) const
220 {
221   // checks if an event is fired using the trigger bits
222
223   return IsTriggerBitFired(aEsd->GetTriggerMask(), trigger);
224 }
225
226 Bool_t AliTriggerAnalysis::IsTriggerBitFired(ULong64_t triggerMask, Trigger trigger) const
227 {
228   // checks if an event is fired using the trigger bits
229   //
230   // this function needs the branch TriggerMask in the ESD
231   
232   // definitions from p-p.cfg
233   ULong64_t spdFO = (1 << 14);
234   ULong64_t v0left = (1 << 10);
235   ULong64_t v0right = (1 << 11);
236
237   switch (trigger)
238   {
239     case kAcceptAll:
240     {
241       return kTRUE;
242       break;
243     }
244     case kMB1:
245     {
246       if (triggerMask & spdFO || ((triggerMask & v0left) || (triggerMask & v0right)))
247         return kTRUE;
248       break;
249     }
250     case kMB2:
251     {
252       if (triggerMask & spdFO && ((triggerMask & v0left) || (triggerMask & v0right)))
253         return kTRUE;
254       break;
255     }
256     case kMB3:
257     {
258       if (triggerMask & spdFO && (triggerMask & v0left) && (triggerMask & v0right))
259         return kTRUE;
260       break;
261     }
262     case kSPDGFO:
263     {
264       if (triggerMask & spdFO)
265         return kTRUE;
266       break;
267     }
268     default:
269       Printf("IsEventTriggered: ERROR: Trigger type %d not implemented in this method", (Int_t) trigger);
270       break;
271   }
272
273   return kFALSE;
274 }
275
276 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, ULong64_t tclass) const
277 {
278   // Checks if corresponding bit in mask is on
279   
280   ULong64_t trigmask = aEsd->GetTriggerMask();
281   return (trigmask & (1ull << (tclass-1)));
282 }
283
284 Bool_t AliTriggerAnalysis::IsOfflineTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
285 {
286   // checks if an event has been triggered "offline"
287
288   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
289   
290   Bool_t decision = kFALSE;
291   switch (triggerNoFlags)
292   {
293     case kAcceptAll:
294     {
295       decision = kTRUE;
296       break;
297     }
298     case kMB1:
299     {
300       if (SPDGFOTrigger(aEsd, 0) || V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
301         decision = kTRUE;
302       break;
303     }
304     case kMB2:
305     {
306       if (SPDGFOTrigger(aEsd, 0) && (V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB))
307         decision = kTRUE;
308       break;
309     }
310     case kMB3:
311     {
312       if (SPDGFOTrigger(aEsd, 0) && V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
313         decision = kTRUE;
314       break;
315     }
316     case kSPDGFO:
317     {
318       if (SPDGFOTrigger(aEsd, 0))
319         decision = kTRUE;
320       break;
321     }
322     case kSPDGFOBits:
323     {
324       if (SPDGFOTrigger(aEsd, 1))
325         decision = kTRUE;
326       break;
327     }
328     case kV0A:
329     {
330       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB)
331         decision = kTRUE;
332       break;
333     }
334     case kV0C:
335     {
336       if (V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
337         decision = kTRUE;
338       break;
339     }
340     case kV0OR:
341     {
342       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
343         decision = kTRUE;
344       break;
345     }
346     case kV0AND:
347     {
348       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
349         decision = kTRUE;
350       break;
351     }
352     case kV0ABG:
353     {
354       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BG)
355         decision = kTRUE;
356       break;
357     }
358     case kV0CBG:
359     {
360       if (V0Trigger(aEsd, kCSide, kFALSE) == kV0BG)
361         decision = kTRUE;
362       break;
363     }
364     case kZDC:
365     {
366       if (ZDCTrigger(aEsd, kASide) || ZDCTrigger(aEsd, kCentralBarrel) || ZDCTrigger(aEsd, kCSide))
367         decision = kTRUE;
368       break;
369     }
370     case kZDCA:
371     {
372       if (ZDCTrigger(aEsd, kASide))
373         decision = kTRUE;
374       break;
375     }
376     case kZDCC:
377     {
378       if (ZDCTrigger(aEsd, kCSide))
379         decision = kTRUE;
380       break;
381     }
382     case kFMDA:
383     {
384       if (FMDTrigger(aEsd, kASide))
385         decision = kTRUE;
386       break;
387     }
388     case kFMDC:
389     {
390       if (FMDTrigger(aEsd, kCSide))
391         decision = kTRUE;
392       break;
393     }
394     case kFPANY:
395     {
396       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))
397         decision = kTRUE;
398       break;
399     }
400     case kNSD1:
401     {
402       if (SPDFiredChips(aEsd, 0) >= 5 || (V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB))
403         decision = kTRUE;
404        break;
405     }
406     case kMB1Prime:
407     {
408       Int_t count = 0;
409       if (SPDGFOTrigger(aEsd, 0))
410         count++;
411       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB)
412         count++;
413       if (V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
414         count++;
415       
416       if (count >= 2)
417         decision = kTRUE;
418         
419       break;
420     }
421     default:
422     {
423       AliFatal(Form("Trigger type %d not implemented", triggerNoFlags));
424     }
425   }
426   
427   // hadron-level requirement
428   if (decision && (trigger & kOneParticle))
429   {
430     decision = kFALSE;
431     
432     const AliESDVertex* vertex = aEsd->GetPrimaryVertexSPD();
433     const AliMultiplicity* mult = aEsd->GetMultiplicity();
434
435     if (mult && vertex && vertex->GetNContributors() > 0 && (!vertex->IsFromVertexerZ() || vertex->GetDispersion() < 0.02) && TMath::Abs(vertex->GetZv()) < 5.5) 
436     {
437       for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
438       {
439         if (TMath::Abs(mult->GetEta(i)) < 1)
440         {
441           decision = kTRUE;
442           break;
443         }
444       }
445     }
446   }
447
448   return decision;
449 }
450
451 Bool_t AliTriggerAnalysis::IsTriggerClassFired(const AliESDEvent* aEsd, const Char_t* tclass) const 
452 {
453   // tclass is logical function of inputs, e.g. 01 && 02 || 03 && 11 && 21
454   // = L0 inp 1 && L0 inp 2 || L0 inp 3 && L1 inp 1 && L2 inp 1
455   // NO brackets in logical function !
456   // Spaces between operators and inputs.
457   // Not all logical functions are available in CTP= 
458   // =any function of first 4 inputs; 'AND' of other inputs, check not done
459   // This method will be replaced/complemened by similar one
460   // which works withh class and inputs names as in CTP cfg file
461   
462   TString TClass(tclass);
463   TObjArray* tcltokens = TClass.Tokenize(" ");
464   Char_t level=((TObjString*)tcltokens->At(0))->String()[0];
465   UInt_t input=atoi((((TObjString*)tcltokens->At(0))->String()).Remove(0));
466   Bool_t tcl = IsInputFired(aEsd,level,input);
467  
468   for (Int_t i=1;i<tcltokens->GetEntriesFast();i=i+2) {
469     level=((TObjString*)tcltokens->At(i+1))->String()[0];
470     input=atoi((((TObjString*)tcltokens->At(i+1))->String()).Remove(0));
471     Bool_t inpnext = IsInputFired(aEsd,level,input);
472     Char_t op =((TObjString*)tcltokens->At(i))->String()[0];
473     if (op == '&') tcl=tcl && inpnext;
474     else if (op == '|') tcl =tcl || inpnext;
475     else {
476        AliError(Form("Syntax error in %s", tclass));
477        tcltokens->Delete();
478        return kFALSE;
479     }
480   }
481   tcltokens->Delete();
482   return tcl;
483 }
484
485 Bool_t AliTriggerAnalysis::IsInputFired(const AliESDEvent* aEsd, Char_t level, UInt_t input) const
486 {
487   // Checks trigger input of any level
488   
489   switch (level)
490   {
491     case '0': return IsL0InputFired(aEsd,input);
492     case '1': return IsL1InputFired(aEsd,input);
493     case '2': return IsL2InputFired(aEsd,input);
494     default:
495       AliError(Form("Wrong level %i",level));
496       return kFALSE;
497   }
498 }
499
500 Bool_t AliTriggerAnalysis::IsL0InputFired(const AliESDEvent* aEsd, UInt_t input) const 
501 {
502   // Checks if corresponding bit in mask is on
503   
504   UInt_t inpmask = aEsd->GetHeader()->GetL0TriggerInputs();
505   return (inpmask & (1<<(input-1)));
506 }
507
508 Bool_t AliTriggerAnalysis::IsL1InputFired(const AliESDEvent* aEsd, UInt_t input) const
509 {
510   // Checks if corresponding bit in mask is on
511   
512   UInt_t inpmask = aEsd->GetHeader()->GetL1TriggerInputs();
513   return (inpmask & (1<<(input-1)));
514 }
515
516 Bool_t AliTriggerAnalysis::IsL2InputFired(const AliESDEvent* aEsd, UInt_t input) const 
517 {
518   // Checks if corresponding bit in mask is on
519   
520   UInt_t inpmask = aEsd->GetHeader()->GetL2TriggerInputs();
521   return (inpmask & (1<<(input-1)));
522 }
523
524 void AliTriggerAnalysis::FillHistograms(const AliESDEvent* aEsd) 
525 {
526   // fills the histograms with the info from the ESD
527   
528   fHistBitsSPD->Fill(SPDFiredChips(aEsd, 0), SPDFiredChips(aEsd, 1, kTRUE));
529   
530   V0Trigger(aEsd, kASide, kFALSE, kTRUE);
531   V0Trigger(aEsd, kCSide, kFALSE, kTRUE);
532   
533   AliESDZDC* zdcData = aEsd->GetESDZDC();
534   if (zdcData)
535   {
536     UInt_t quality = zdcData->GetESDQuality();
537     
538     // from Nora's presentation, general first physics meeting 16.10.09
539     static UInt_t zpc  = 0x20;
540     static UInt_t znc  = 0x10;
541     static UInt_t zem1 = 0x08;
542     static UInt_t zem2 = 0x04;
543     static UInt_t zpa  = 0x02;
544     static UInt_t zna  = 0x01;
545    
546     fHistZDC->Fill(1, quality & zna);
547     fHistZDC->Fill(2, quality & zpa);
548     fHistZDC->Fill(3, quality & zem2);
549     fHistZDC->Fill(4, quality & zem1);
550     fHistZDC->Fill(5, quality & znc);
551     fHistZDC->Fill(6, quality & zpc);
552   }
553   else
554   {
555     fHistZDC->Fill(-1);
556     AliError("AliESDZDC not available");
557   }
558   
559   fHistFMDA->Fill(FMDHitCombinations(aEsd, kASide, kTRUE));
560   fHistFMDC->Fill(FMDHitCombinations(aEsd, kCSide, kTRUE));
561 }
562   
563 void AliTriggerAnalysis::FillTriggerClasses(const AliESDEvent* aEsd)
564 {
565   // fills trigger classes map
566   
567   TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(aEsd->GetFiredTriggerClasses().Data()));
568   if (!count)
569   {
570     count = new TParameter<Long64_t>(aEsd->GetFiredTriggerClasses(), 0);
571     fTriggerClasses->Add(new TObjString(aEsd->GetFiredTriggerClasses().Data()), count);
572   }
573   count->SetVal(count->GetVal() + 1);
574   
575   // TODO add first and last orbit number here
576 }
577
578 Int_t AliTriggerAnalysis::SSDClusters(const AliESDEvent* aEsd)
579 {
580   // returns the number of clusters in the SSD
581   const AliMultiplicity* mult = aEsd->GetMultiplicity();
582   Int_t clusters = mult->GetNumberOfITSClusters(4)+mult->GetNumberOfITSClusters(5);
583   return clusters;
584 }
585
586
587 Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHists)
588 {
589   // returns the number of fired chips in the SPD
590   //
591   // origin = 0 --> aEsd->GetMultiplicity()->GetNumberOfFiredChips() (filled from clusters)
592   // origin = 1 --> aEsd->GetMultiplicity()->TestFastOrFiredChips() (from hardware bits)
593   
594   const AliMultiplicity* mult = aEsd->GetMultiplicity();
595   if (!mult)
596   {
597     AliError("AliMultiplicity not available");
598     return -1;
599   }
600   
601   if (origin == 0)
602     return mult->GetNumberOfFiredChips(0) + mult->GetNumberOfFiredChips(1);
603     
604   if (origin == 1)
605   {
606     Int_t nChips = 0;
607     for (Int_t i=0; i<1200; i++)
608     {
609       if (mult->TestFastOrFiredChips(i) == kTRUE)
610       {
611         // efficiency simulation (if enabled)
612         if (fSPDGFOEfficiency)
613         {
614           if (gRandom->Uniform() > fSPDGFOEfficiency->GetBinContent(i+1))
615             continue;
616         }
617         
618         nChips++;
619         if (fillHists)
620           fHistFiredBitsSPD->Fill(i);
621       }
622     }
623     return nChips;
624   }
625   
626   return -1;
627 }
628
629 Bool_t AliTriggerAnalysis::SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin)
630 {
631   // Returns if the SPD gave a global Fast OR trigger
632   
633   Int_t firedChips = SPDFiredChips(aEsd, origin);
634   
635   if (firedChips >= fSPDGFOThreshold)
636     return kTRUE;
637   return kFALSE;
638 }
639
640 AliTriggerAnalysis::V0Decision AliTriggerAnalysis::V0Trigger(const AliESDEvent* aEsd, AliceSide side, Bool_t online, Bool_t fillHists)
641 {
642   // Returns the V0 trigger decision in V0A | V0C
643   //
644   // Returns kV0Fake if the calculated average time is in a window where neither BB nor BG is expected. 
645   // The rate of such triggers can be used to estimate the background. Note that the rate has to be 
646   // rescaled with the size of the windows (numerical values see below in the code)
647   //
648   // argument 'online' is used as a switch between online and offline trigger algorithms
649   //
650   // Based on an algorithm by Cvetan Cheshkov
651   
652   AliESDVZERO* esdV0 = aEsd->GetVZEROData();
653   if (!esdV0)
654   {
655     AliError("AliESDVZERO not available");
656     return kV0Invalid;
657   }
658
659   if (esdV0->TestBit(AliESDVZERO::kDecisionFilled)) {
660     if (online) {
661       AliWarning("V0 online trigger analysis is not yet available!");
662       return kV0BB;
663     }
664     else {
665
666       if (fillHists) {
667         if (side == kASide && fHistV0A)
668           fHistV0A->Fill(esdV0->GetV0ATime());
669         if (side == kCSide && fHistV0C)
670           fHistV0C->Fill(esdV0->GetV0CTime());
671       }
672
673       if (side == kASide) return (V0Decision)esdV0->GetV0ADecision();
674       else if (side == kCSide) return (V0Decision)esdV0->GetV0CDecision();
675       else return kV0Invalid;
676     }
677   }
678
679   Int_t begin = -1;
680   Int_t end = -1;
681   
682   if (side == kASide)
683   {
684     begin = 32;
685     end = 64;
686   } 
687   else if (side == kCSide)
688   {
689     begin = 0;
690     end = 32;
691   }
692   else
693     return kV0Invalid;
694     
695   Float_t time = 0;
696   Float_t weight = 0;
697   if (fMC)
698   {
699     Int_t runRange;
700     if (aEsd->GetRunNumber() <= 104803) runRange = 0;
701     else if (aEsd->GetRunNumber() <= 104876) runRange = 1;
702     else runRange = 2;
703
704     Float_t factors[3][64] = {
705       // runs: 104792-104803
706       {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},
707       // runs: 104841-104876
708       {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},
709       // runs: 104890-92
710       {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}
711     };
712     Float_t dA = 77.4 - 11.0;
713     Float_t dC = 77.4 - 2.9;
714     // Time misalignment
715     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};
716     Float_t dA2 = 2.8, dC2 = 3.3;
717
718     if (online) {
719       for (Int_t i = begin; i < end; ++i) {
720         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
721         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
722         if (esdV0->GetTime(i) >= 1e-6 &&
723             tempTime > fV0HwWinLow && tempTime < fV0HwWinHigh &&
724             tempAdc > fV0HwAdcThr)
725           return kV0BB;
726       }
727       return kV0Empty;
728     }
729     else {
730       for (Int_t i = begin; i < end; ++i) {
731         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
732         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA : esdV0->GetTime(i)+dC;
733         Float_t tempRawTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
734         if (esdV0->GetTime(i) >= 1e-6 &&
735             tempRawTime < 125.0 &&
736             tempAdc > fV0AdcThr) {
737           weight += 1.0;
738           time += tempTime;
739         }
740       }
741     }
742   }
743   else {
744     if (online) {
745       for (Int_t i = begin; i < end; ++i) {
746         if (esdV0->GetTime(i) >= 1e-6 &&
747             esdV0->GetTime(i) > fV0HwWinLow && esdV0->GetTime(i) < fV0HwWinHigh &&
748             esdV0->GetAdc(i) > fV0HwAdcThr)
749           return kV0BB;
750       }
751       return kV0Empty;
752     }
753     else {
754       for (Int_t i = begin; i < end; ++i) {
755         if (esdV0->GetTime(i) > 1e-6 && esdV0->GetAdc(i) > fV0AdcThr) {
756           Float_t correctedTime = V0CorrectLeadingTime(i, esdV0->GetTime(i), esdV0->GetAdc(i),aEsd->GetRunNumber());
757           Float_t timeWeight = V0LeadingTimeWeight(esdV0->GetAdc(i));
758           time += correctedTime*timeWeight;
759             
760           weight += timeWeight;
761         }
762       }
763     }
764   }
765
766   if (weight > 0) 
767     time /= weight;
768   time += fV0TimeOffset;
769
770   if (fillHists)
771   {
772     if (side == kASide && fHistV0A)
773       fHistV0A->Fill(time);
774     if (side == kCSide && fHistV0C)
775       fHistV0C->Fill(time);
776   }
777   
778   if (side == kASide)
779   {
780     if (time > 68 && time < 100)
781       return kV0BB;
782     if (time > 54 && time < 57.5) 
783       return kV0BG;
784     if (time > 57.5 && time < 68)
785       return kV0Fake;
786   }
787   
788   if (side == kCSide)
789   {
790     if (time > 75.5 && time < 100)
791       return kV0BB;
792     if (time > 69.5 && time < 73)
793       return kV0BG; 
794     if (time > 55 && time < 69.5)
795       return kV0Fake;
796   }
797   
798   return kV0Empty;
799 }
800
801 Float_t AliTriggerAnalysis::V0CorrectLeadingTime(Int_t i, Float_t time, Float_t adc, Int_t runNumber) const
802 {
803   // Correct for slewing and align the channels
804   //
805   // Authors: Cvetan Cheshkov / Raphael Tieulent
806
807   if (time == 0) return 0;
808
809   // Time alignment
810   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};
811
812   if(runNumber < 106031)
813     time -= timeShift[i];
814
815   // Slewing correction
816   if (adc == 0) return time;
817
818   Float_t p1 = 1.57345e1;
819   Float_t p2 =-4.25603e-1;
820
821   if(runNumber >= 106031) adc *= (2.5/4.0);
822   return (time - p1*TMath::Power(adc,p2));
823 }
824
825 Float_t AliTriggerAnalysis::V0LeadingTimeWeight(Float_t adc) const
826 {
827   if (adc < 1e-6) return 0;
828
829   Float_t p1 = 40.211;
830   Float_t p2 =-4.25603e-1;
831   Float_t p3 = 0.5646;
832
833   return 1./(p1*p1*TMath::Power(adc,2.*(p2-1.))+p3*p3);
834 }
835
836 Bool_t AliTriggerAnalysis::ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) const
837 {
838   // Returns if ZDC triggered
839   
840   AliESDZDC* zdcData = aEsd->GetESDZDC();
841   if (!zdcData)
842   {
843     AliError("AliESDZDC not available");
844     return kFALSE;
845   }
846   
847   UInt_t quality = zdcData->GetESDQuality();
848   
849   // from Nora's presentation, general first physics meeting 16.10.09
850   static UInt_t zpc  = 0x20;
851   static UInt_t znc  = 0x10;
852   static UInt_t zem1 = 0x08;
853   static UInt_t zem2 = 0x04;
854   static UInt_t zpa  = 0x02;
855   static UInt_t zna  = 0x01;
856   
857   if (side == kASide && ((quality & zpa) || (quality & zna)))
858     return kTRUE;
859   if (side == kCentralBarrel && ((quality & zem1) || (quality & zem2)))
860     return kTRUE;
861   if (side == kCSide && ((quality & zpc) || (quality & znc)))
862     return kTRUE;
863   
864   return kFALSE;
865 }
866
867 Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists)
868 {
869   // returns number of hit combinations agove threshold
870   //
871   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
872
873   // Workaround for AliESDEvent::GetFMDData is not const!
874   const AliESDFMD* fmdData = (const_cast<AliESDEvent*>(aEsd))->GetFMDData();
875   if (!fmdData)
876   {
877     AliError("AliESDFMD not available");
878     return -1;
879   }
880
881   Int_t detFrom = (side == kASide) ? 1 : 3;
882   Int_t detTo   = (side == kASide) ? 2 : 3;
883
884   Int_t triggers = 0;
885   Float_t totalMult = 0;
886   for (UShort_t det=detFrom;det<=detTo;det++) {
887     Int_t nRings = (det == 1 ? 1 : 2);
888     for (UShort_t ir = 0; ir < nRings; ir++) {
889       Char_t   ring = (ir == 0 ? 'I' : 'O');
890       UShort_t nsec = (ir == 0 ? 20  : 40);
891       UShort_t nstr = (ir == 0 ? 512 : 256);
892       for (UShort_t sec =0; sec < nsec;  sec++) {
893         for (UShort_t strip = 0; strip < nstr; strip++) {
894           Float_t mult = fmdData->Multiplicity(det,ring,sec,strip);
895           if (mult == AliESDFMD::kInvalidMult) continue;
896           
897           if (fillHists)
898             fHistFMDSingle->Fill(mult);
899           
900           if (mult > fFMDLowCut)
901             totalMult = totalMult + mult;
902           else
903           {
904             if (totalMult > fFMDHitCut)
905               triggers++;
906               
907             if (fillHists)
908               fHistFMDSum->Fill(totalMult);
909               
910             totalMult = 0;
911           }
912         }
913       }
914     }
915   }
916   
917   return triggers;
918 }
919
920 Bool_t AliTriggerAnalysis::FMDTrigger(const AliESDEvent* aEsd, AliceSide side)
921 {
922   // Returns if the FMD triggered
923   //
924   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
925
926   Int_t triggers = FMDHitCombinations(aEsd, side, kFALSE);
927   
928   if (triggers > 0)
929     return kTRUE;
930     
931   return kFALSE;
932 }
933
934 Long64_t AliTriggerAnalysis::Merge(TCollection* list)
935 {
936   // Merge a list of AliMultiplicityCorrection objects with this (needed for
937   // PROOF).
938   // Returns the number of merged objects (including this).
939
940   if (!list)
941     return 0;
942
943   if (list->IsEmpty())
944     return 1;
945
946   TIterator* iter = list->MakeIterator();
947   TObject* obj;
948
949   // collections of all histograms
950   const Int_t nHists = 9;
951   TList collections[nHists];
952
953   Int_t count = 0;
954   while ((obj = iter->Next())) {
955
956     AliTriggerAnalysis* entry = dynamic_cast<AliTriggerAnalysis*> (obj);
957     if (entry == 0) 
958       continue;
959
960     Int_t n = 0;
961     collections[n++].Add(entry->fHistV0A);
962     collections[n++].Add(entry->fHistV0C);
963     collections[n++].Add(entry->fHistZDC);
964     collections[n++].Add(entry->fHistFMDA);
965     collections[n++].Add(entry->fHistFMDC);
966     collections[n++].Add(entry->fHistFMDSingle);
967     collections[n++].Add(entry->fHistFMDSum);
968     collections[n++].Add(entry->fHistBitsSPD);
969     collections[n++].Add(entry->fHistFiredBitsSPD);
970
971     // merge fTriggerClasses
972     TIterator* iter2 = entry->fTriggerClasses->MakeIterator();
973     TObjString* obj2 = 0;
974     while ((obj2 = dynamic_cast<TObjString*> (iter2->Next())))
975     {
976       TParameter<Long64_t>* param2 = dynamic_cast<TParameter<Long64_t>*> (entry->fTriggerClasses->GetValue(obj2));
977       
978       TParameter<Long64_t>* param1 = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj2));
979       if (param1)
980       {
981         param1->SetVal(param1->GetVal() + param2->GetVal());
982       }
983       else
984       {
985         param1 = dynamic_cast<TParameter<Long64_t>*> (param2->Clone());
986         fTriggerClasses->Add(new TObjString(obj2->String()), param1);
987       }
988     }
989     
990     delete iter2;
991   
992     count++;
993   }
994
995   Int_t n = 0;
996   fHistV0A->Merge(&collections[n++]);
997   fHistV0C->Merge(&collections[n++]);
998   fHistZDC->Merge(&collections[n++]);
999   fHistFMDA->Merge(&collections[n++]);
1000   fHistFMDC->Merge(&collections[n++]);
1001   fHistFMDSingle->Merge(&collections[n++]);
1002   fHistFMDSum->Merge(&collections[n++]);
1003   fHistBitsSPD->Merge(&collections[n++]);
1004   fHistFiredBitsSPD->Merge(&collections[n++]);
1005   
1006   delete iter;
1007
1008   return count+1;
1009 }
1010
1011 void AliTriggerAnalysis::SaveHistograms() const
1012 {
1013   // write histograms to current directory
1014   
1015   if (!fHistBitsSPD)
1016     return;
1017     
1018   fHistBitsSPD->Write();
1019   fHistBitsSPD->ProjectionX();
1020   fHistBitsSPD->ProjectionY();
1021   fHistFiredBitsSPD->Write();
1022   fHistV0A->Write();
1023   fHistV0C->Write();
1024   fHistZDC->Write();
1025   fHistFMDA->Write();
1026   fHistFMDC->Write();
1027   fHistFMDSingle->Write();
1028   fHistFMDSum->Write();
1029   
1030   if (fSPDGFOEfficiency)
1031     fSPDGFOEfficiency->Write("fSPDGFOEfficiency");
1032   
1033   fTriggerClasses->Write("fTriggerClasses", TObject::kSingleKey);
1034 }
1035
1036 void AliTriggerAnalysis::PrintTriggerClasses() const
1037 {
1038   // print trigger classes
1039   
1040   Printf("Trigger Classes:");
1041   
1042   Printf("Event count for trigger combinations:");
1043   
1044   TMap singleTrigger;
1045   singleTrigger.SetOwner();
1046   
1047   TIterator* iter = fTriggerClasses->MakeIterator();
1048   TObjString* obj = 0;
1049   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
1050   {
1051     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj));
1052     
1053     Printf(" %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1054     
1055     TObjArray* tokens = obj->String().Tokenize(" ");
1056     for (Int_t i=0; i<tokens->GetEntries(); i++)
1057     {
1058       TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(((TObjString*) tokens->At(i))->String().Data()));
1059       if (!count)
1060       {
1061         count = new TParameter<Long64_t>(((TObjString*) tokens->At(i))->String().Data(), 0);
1062         singleTrigger.Add(new TObjString(((TObjString*) tokens->At(i))->String().Data()), count);
1063       }
1064       count->SetVal(count->GetVal() + param->GetVal());
1065     }
1066     
1067     delete tokens;
1068   }
1069   delete iter;
1070   
1071   Printf("Event count for single trigger:");
1072   
1073   iter = singleTrigger.MakeIterator();
1074   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
1075   {
1076     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(obj));
1077     
1078     Printf("  %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1079   }
1080   delete iter;
1081   
1082   singleTrigger.DeleteAll();
1083 }