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