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