]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliTriggerAnalysis.cxx
7e5b17c55ac29fd889e6a8703d71976274772ed0
[u/mrichter/AliRoot.git] / ANALYSIS / AliTriggerAnalysis.cxx
1 /* $Id: AliTriggerAnalysis.cxx 35782 2009-10-22 11:54:31Z jgrosseo $ */
2
3 /**************************************************************************
4  * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
5  *                                                                        *
6  * Author: The ALICE Off-line Project.                                    *
7  * Contributors are mentioned in the code where appropriate.              *
8  *                                                                        *
9  * Permission to use, copy, modify and distribute this software and its   *
10  * documentation strictly for non-commercial purposes is hereby granted   *
11  * without fee, provided that the above copyright notice appears in all   *
12  * copies and that both the copyright notice and this permission notice   *
13  * appear in the supporting documentation. The authors make no claims     *
14  * about the suitability of this software for any purpose. It is          *
15  * provided "as is" without express or implied warranty.                  *
16  **************************************************************************/
17
18 //-------------------------------------------------------------------------
19 //                      Implementation of   Class AliTriggerAnalysis
20 //   This class provides function to check if events have been triggered based on the data in the ESD
21 //   The trigger bits, trigger class inputs and only the data (offline trigger) can be used
22 //   Origin: Jan Fiete Grosse-Oetringhaus, CERN
23 //-------------------------------------------------------------------------
24
25 #include <Riostream.h>
26 #include <TH1F.h>
27 #include <TH2F.h>
28 #include <TList.h>
29 #include <TIterator.h>
30 #include "TParameter.h"
31 #include <TMap.h>
32 #include <TRandom.h>
33
34 #include <AliTriggerAnalysis.h>
35
36 #include <AliLog.h>
37
38 #include <AliESDEvent.h>
39
40 #include <AliMultiplicity.h>
41 #include <AliESDVZERO.h>
42 #include <AliESDZDC.h>
43 #include <AliESDFMD.h>
44 #include <AliESDVertex.h>
45 #include <AliESDtrackCuts.h>
46
47 ClassImp(AliTriggerAnalysis)
48
49 AliTriggerAnalysis::AliTriggerAnalysis() :
50   fSPDGFOThreshold(2),
51   fSPDGFOEfficiency(0),
52   fV0TimeOffset(0),
53   fV0AdcThr(0),
54   fV0HwAdcThr(2.5),
55   fV0HwWinLow(61.5),
56   fV0HwWinHigh(86.5),
57   fDoFMD(kTRUE),
58   fFMDLowCut(0.2),
59   fFMDHitCut(0.5),
60   fHistBitsSPD(0),
61   fHistFiredBitsSPD(0),
62   fHistV0A(0),       
63   fHistV0C(0),
64   fHistZDC(0),    
65   fHistTDCZDC(0),    
66   fHistTimeZDC(0),    
67   fHistFMDA(0),    
68   fHistFMDC(0),   
69   fHistFMDSingle(0),
70   fHistFMDSum(0),
71   fTriggerClasses(0),
72   fMC(kFALSE),
73   fEsdTrackCuts(0)
74 {
75   // constructor
76 }
77
78 AliTriggerAnalysis::~AliTriggerAnalysis()
79 {
80   // destructor
81   
82   if (fHistBitsSPD)
83   {
84     delete fHistBitsSPD;
85     fHistBitsSPD = 0;
86   }
87
88   if (fHistFiredBitsSPD)
89   {
90     delete fHistFiredBitsSPD;
91     fHistFiredBitsSPD = 0;
92   }
93
94   if (fHistV0A)
95   {
96     delete fHistV0A;
97     fHistV0A = 0;
98   }
99
100   if (fHistV0C)
101   {
102     delete fHistV0C;
103     fHistV0C = 0;
104   }
105
106   if (fHistZDC)
107   {
108     delete fHistZDC;
109     fHistZDC = 0;
110   }
111   if (fHistTDCZDC)
112   {
113     delete fHistTDCZDC;
114     fHistTDCZDC = 0;
115   }
116   if (fHistTimeZDC)
117   {
118     delete fHistTimeZDC;
119     fHistTimeZDC = 0;
120   }
121
122   if (fHistFMDA)
123   {
124     delete fHistFMDA;
125     fHistFMDA = 0;
126   }
127
128   if (fHistFMDC)
129   {
130     delete fHistFMDC;
131     fHistFMDC = 0;
132   }
133
134   if (fHistFMDSingle)
135   {
136     delete fHistFMDSingle;
137     fHistFMDSingle = 0;
138   }
139
140   if (fHistFMDSum)
141   {
142     delete fHistFMDSum;
143     fHistFMDSum = 0;
144   }
145
146   if (fTriggerClasses)
147   {
148     fTriggerClasses->DeleteAll();
149     delete fTriggerClasses;
150     fTriggerClasses = 0;
151   }
152
153   if (fEsdTrackCuts){
154     delete fEsdTrackCuts;
155     fEsdTrackCuts =0;
156   }
157 }
158
159 void AliTriggerAnalysis::EnableHistograms()
160 {
161   // creates the monitoring histograms
162   
163   // do not add this hists to the directory
164   Bool_t oldStatus = TH1::AddDirectoryStatus();
165   TH1::AddDirectory(kFALSE);
166   
167   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);
168   fHistFiredBitsSPD = new TH1F("fHistFiredBitsSPD", "SPD GFO Hardware;chip number;events", 1200, -0.5, 1199.5);
169   fHistV0A = new TH1F("fHistV0A", "V0A;leading time (ns);events", 400, -100, 100);
170   fHistV0C = new TH1F("fHistV0C", "V0C;leading time (ns);events", 400, -100, 100);
171   fHistZDC = new TH1F("fHistZDC", "ZDC;trigger bits;events", 8, -1.5, 6.5);
172   fHistTDCZDC = new TH1F("fHistTDCZDC", "ZDC;TDC bits;events", 32, -0.5, 32-0.5);
173   fHistTimeZDC = new TH2F("fHistTimeZDC", "ZDC;TDC timing A+C vs C-A; events", 120,-30,30,120,-600,-540);
174   
175   // TODO check limits
176   fHistFMDA = new TH1F("fHistFMDA", "FMDA;combinations above threshold;events", 102, -1.5, 100.5);
177   fHistFMDC = new TH1F("fHistFMDC", "FMDC;combinations above threshold;events", 102, -1.5, 100.5);
178   fHistFMDSingle = new TH1F("fHistFMDSingle", "FMD single;multiplicity value;counts", 1000, 0, 10);
179   fHistFMDSum = new TH1F("fHistFMDSum", "FMD sum;multiplicity value;counts", 1000, 0, 10);
180   
181   fTriggerClasses = new TMap;
182   fTriggerClasses->SetOwner();
183   
184   TH1::AddDirectory(oldStatus);
185 }
186
187 //____________________________________________________________________
188 const char* AliTriggerAnalysis::GetTriggerName(Trigger trigger) 
189 {
190   // returns the name of the requested trigger
191   // the returned string will only be valid until the next call to this function [not thread-safe]
192   
193   static TString str;
194   
195   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
196   
197   switch (triggerNoFlags)
198   {
199     case kAcceptAll : str = "ACCEPT ALL (bypass!)"; break;
200     case kMB1 : str = "MB1"; break;
201     case kMB2 : str = "MB2"; break;
202     case kMB3 : str = "MB3"; break;
203     case kSPDGFO : str = "SPD GFO"; break;
204     case kSPDGFOBits : str = "SPD GFO Bits"; break;
205     case kV0A : str = "V0 A BB"; break;
206     case kV0C : str = "V0 C BB"; break;
207     case kV0OR : str = "V0 OR BB"; break;
208     case kV0AND : str = "V0 AND BB"; break;
209     case kV0ABG : str = "V0 A BG"; break;
210     case kV0CBG : str = "V0 C BG"; break;
211     case kZDC : str = "ZDC"; break;
212     case kZDCA : str = "ZDC A"; break;
213     case kZDCC : str = "ZDC C"; break;
214     case kFMDA : str = "FMD A"; break;
215     case kFMDC : str = "FMD C"; break;
216     case kFPANY : str = "SPD GFO | V0 | ZDC | FMD"; break;
217     case kNSD1 : str = "NSD1"; break;
218     case kMB1Prime: str = "MB1prime"; break;
219     default: str = ""; break;
220   }
221    
222   if (trigger & kOfflineFlag)
223     str += " OFFLINE";  
224   
225   if (trigger & kOneParticle)
226     str += " OneParticle";  
227
228   if (trigger & kOneTrack)
229     str += " OneTrack";  
230
231   return str;
232 }
233
234 Bool_t AliTriggerAnalysis::IsTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
235 {
236   // checks if an event has been triggered
237
238   if (trigger & kOfflineFlag)
239     return IsOfflineTriggerFired(aEsd, trigger);
240     
241   return IsTriggerBitFired(aEsd, trigger);
242 }
243
244 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, Trigger trigger) const
245 {
246   // checks if an event is fired using the trigger bits
247
248   return IsTriggerBitFired(aEsd->GetTriggerMask(), trigger);
249 }
250
251 Bool_t AliTriggerAnalysis::IsTriggerBitFired(ULong64_t triggerMask, Trigger trigger) const
252 {
253   // checks if an event is fired using the trigger bits
254   //
255   // this function needs the branch TriggerMask in the ESD
256   
257   // definitions from p-p.cfg
258   ULong64_t spdFO = (1 << 14);
259   ULong64_t v0left = (1 << 10);
260   ULong64_t v0right = (1 << 11);
261
262   switch (trigger)
263   {
264     case kAcceptAll:
265     {
266       return kTRUE;
267       break;
268     }
269     case kMB1:
270     {
271       if (triggerMask & spdFO || ((triggerMask & v0left) || (triggerMask & v0right)))
272         return kTRUE;
273       break;
274     }
275     case kMB2:
276     {
277       if (triggerMask & spdFO && ((triggerMask & v0left) || (triggerMask & v0right)))
278         return kTRUE;
279       break;
280     }
281     case kMB3:
282     {
283       if (triggerMask & spdFO && (triggerMask & v0left) && (triggerMask & v0right))
284         return kTRUE;
285       break;
286     }
287     case kSPDGFO:
288     {
289       if (triggerMask & spdFO)
290         return kTRUE;
291       break;
292     }
293     default:
294       Printf("IsEventTriggered: ERROR: Trigger type %d not implemented in this method", (Int_t) trigger);
295       break;
296   }
297
298   return kFALSE;
299 }
300
301 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, ULong64_t tclass) const
302 {
303   // Checks if corresponding bit in mask is on
304   
305   ULong64_t trigmask = aEsd->GetTriggerMask();
306   return (trigmask & (1ull << (tclass-1)));
307 }
308
309 Bool_t AliTriggerAnalysis::IsOfflineTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
310 {
311   // checks if an event has been triggered "offline"
312
313   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
314   
315   Bool_t decision = kFALSE;
316   switch (triggerNoFlags)
317   {
318     case kAcceptAll:
319     {
320       decision = kTRUE;
321       break;
322     }
323     case kMB1:
324     {
325       if (SPDGFOTrigger(aEsd, 0) || V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
326         decision = kTRUE;
327       break;
328     }
329     case kMB2:
330     {
331       if (SPDGFOTrigger(aEsd, 0) && (V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB))
332         decision = kTRUE;
333       break;
334     }
335     case kMB3:
336     {
337       if (SPDGFOTrigger(aEsd, 0) && V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
338         decision = kTRUE;
339       break;
340     }
341     case kSPDGFO:
342     {
343       if (SPDGFOTrigger(aEsd, 0))
344         decision = kTRUE;
345       break;
346     }
347     case kSPDGFOBits:
348     {
349       if (SPDGFOTrigger(aEsd, 1))
350         decision = kTRUE;
351       break;
352     }
353     case kV0A:
354     {
355       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB)
356         decision = kTRUE;
357       break;
358     }
359     case kV0C:
360     {
361       if (V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
362         decision = kTRUE;
363       break;
364     }
365     case kV0OR:
366     {
367       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
368         decision = kTRUE;
369       break;
370     }
371     case kV0AND:
372     {
373       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
374         decision = kTRUE;
375       break;
376     }
377     case kV0ABG:
378     {
379       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BG)
380         decision = kTRUE;
381       break;
382     }
383     case kV0CBG:
384     {
385       if (V0Trigger(aEsd, kCSide, kFALSE) == kV0BG)
386         decision = kTRUE;
387       break;
388     }
389     case kZDC:
390     {
391       if (ZDCTrigger(aEsd, kASide) || ZDCTrigger(aEsd, kCentralBarrel) || ZDCTrigger(aEsd, kCSide))
392         decision = kTRUE;
393       break;
394     }
395     case kZDCA:
396     {
397       if (ZDCTrigger(aEsd, kASide))
398         decision = kTRUE;
399       break;
400     }
401     case kZDCC:
402     {
403       if (ZDCTrigger(aEsd, kCSide))
404         decision = kTRUE;
405       break;
406     }
407     case kFMDA:
408     {
409       if (FMDTrigger(aEsd, kASide))
410         decision = kTRUE;
411       break;
412     }
413     case kFMDC:
414     {
415       if (FMDTrigger(aEsd, kCSide))
416         decision = kTRUE;
417       break;
418     }
419     case kFPANY:
420     {
421       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))
422         decision = kTRUE;
423       break;
424     }
425     case kNSD1:
426     {
427       if (SPDFiredChips(aEsd, 0) >= 5 || (V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB))
428         decision = kTRUE;
429        break;
430     }
431     case kMB1Prime:
432     {
433       Int_t count = 0;
434       if (SPDGFOTrigger(aEsd, 0))
435         count++;
436       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB)
437         count++;
438       if (V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
439         count++;
440       
441       if (count >= 2)
442         decision = kTRUE;
443         
444       break;
445     }
446     default:
447     {
448       AliFatal(Form("Trigger type %d not implemented", triggerNoFlags));
449     }
450   }
451   
452   // hadron-level requirement
453   if (decision && (trigger & kOneParticle))
454   {
455     decision = kFALSE;
456     
457     const AliESDVertex* vertex = aEsd->GetPrimaryVertexSPD();
458     const AliMultiplicity* mult = aEsd->GetMultiplicity();
459
460     if (mult && vertex && vertex->GetNContributors() > 0 && (!vertex->IsFromVertexerZ() || vertex->GetDispersion() < 0.02) && TMath::Abs(vertex->GetZv()) < 5.5) 
461     {
462       for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
463       {
464         if (TMath::Abs(mult->GetEta(i)) < 1)
465         {
466           decision = kTRUE;
467           break;
468         }
469       }
470     }
471   }
472
473   // hadron level definition for TPC tracks
474
475   if (decision && (trigger & kOneTrack))
476   {
477     decision = kFALSE;
478     const AliESDVertex* vertex =0x0;
479     vertex = aEsd->GetPrimaryVertexTracks();
480     if (!vertex || vertex->GetNContributors() <= 0)
481     {
482       vertex = aEsd->GetPrimaryVertexSPD();
483     }
484     Float_t ptmin, ptmax;
485     fEsdTrackCuts->GetPtRange(ptmin,ptmax);
486     AliDebug(3, Form("ptmin = %f, ptmax = %f\n",ptmin, ptmax));
487
488     if (vertex && vertex->GetNContributors() > 0 && (!vertex->IsFromVertexerZ() || vertex->GetDispersion() < 0.02) && TMath::Abs(vertex->GetZv()) < 10.) {
489       AliDebug(3,Form("Check on the vertex passed\n"));
490       for (Int_t i=0; i<aEsd->GetNumberOfTracks(); ++i){
491         if (fEsdTrackCuts->AcceptTrack(aEsd->GetTrack(i))){
492           AliDebug(2, Form("pt of track = %f --> check passed\n",aEsd->GetTrack(i)->Pt()));
493           decision = kTRUE;
494           break;
495         }
496       }
497     }
498     else{
499       AliDebug(4,Form("Check on the vertex not passed\n"));
500       for (Int_t i=0; i<aEsd->GetNumberOfTracks(); ++i){
501         if (fEsdTrackCuts->AcceptTrack(aEsd->GetTrack(i))){
502           AliDebug(4,Form("pt of track = %f --> check would be passed if the vertex was ok\n",aEsd->GetTrack(i)->Pt()));
503           break;
504         }
505       }
506     }
507     if (!decision) AliDebug(3,("Check for kOneTrack NOT passed\n"));
508   }
509
510   return decision;
511 }
512
513 Bool_t AliTriggerAnalysis::IsTriggerClassFired(const AliESDEvent* aEsd, const Char_t* tclass) const 
514 {
515   // tclass is logical function of inputs, e.g. 01 && 02 || 03 && 11 && 21
516   // = L0 inp 1 && L0 inp 2 || L0 inp 3 && L1 inp 1 && L2 inp 1
517   // NO brackets in logical function !
518   // Spaces between operators and inputs.
519   // Not all logical functions are available in CTP= 
520   // =any function of first 4 inputs; 'AND' of other inputs, check not done
521   // This method will be replaced/complemened by similar one
522   // which works withh class and inputs names as in CTP cfg file
523   
524   TString TClass(tclass);
525   TObjArray* tcltokens = TClass.Tokenize(" ");
526   Char_t level=((TObjString*)tcltokens->At(0))->String()[0];
527   UInt_t input=atoi((((TObjString*)tcltokens->At(0))->String()).Remove(0));
528   Bool_t tcl = IsInputFired(aEsd,level,input);
529  
530   for (Int_t i=1;i<tcltokens->GetEntriesFast();i=i+2) {
531     level=((TObjString*)tcltokens->At(i+1))->String()[0];
532     input=atoi((((TObjString*)tcltokens->At(i+1))->String()).Remove(0));
533     Bool_t inpnext = IsInputFired(aEsd,level,input);
534     Char_t op =((TObjString*)tcltokens->At(i))->String()[0];
535     if (op == '&') tcl=tcl && inpnext;
536     else if (op == '|') tcl =tcl || inpnext;
537     else {
538        AliError(Form("Syntax error in %s", tclass));
539        tcltokens->Delete();
540        return kFALSE;
541     }
542   }
543   tcltokens->Delete();
544   return tcl;
545 }
546
547 Bool_t AliTriggerAnalysis::IsInputFired(const AliESDEvent* aEsd, Char_t level, UInt_t input) const
548 {
549   // Checks trigger input of any level
550   
551   switch (level)
552   {
553     case '0': return IsL0InputFired(aEsd,input);
554     case '1': return IsL1InputFired(aEsd,input);
555     case '2': return IsL2InputFired(aEsd,input);
556     default:
557       AliError(Form("Wrong level %i",level));
558       return kFALSE;
559   }
560 }
561
562 Bool_t AliTriggerAnalysis::IsL0InputFired(const AliESDEvent* aEsd, UInt_t input) const 
563 {
564   // Checks if corresponding bit in mask is on
565   
566   UInt_t inpmask = aEsd->GetHeader()->GetL0TriggerInputs();
567   return (inpmask & (1<<(input-1)));
568 }
569
570 Bool_t AliTriggerAnalysis::IsL1InputFired(const AliESDEvent* aEsd, UInt_t input) const
571 {
572   // Checks if corresponding bit in mask is on
573   
574   UInt_t inpmask = aEsd->GetHeader()->GetL1TriggerInputs();
575   return (inpmask & (1<<(input-1)));
576 }
577
578 Bool_t AliTriggerAnalysis::IsL2InputFired(const AliESDEvent* aEsd, UInt_t input) const 
579 {
580   // Checks if corresponding bit in mask is on
581   
582   UInt_t inpmask = aEsd->GetHeader()->GetL2TriggerInputs();
583   return (inpmask & (1<<(input-1)));
584 }
585
586 void AliTriggerAnalysis::FillHistograms(const AliESDEvent* aEsd) 
587 {
588   // fills the histograms with the info from the ESD
589   
590   fHistBitsSPD->Fill(SPDFiredChips(aEsd, 0), SPDFiredChips(aEsd, 1, kTRUE));
591   
592   V0Trigger(aEsd, kASide, kFALSE, kTRUE);
593   V0Trigger(aEsd, kCSide, kFALSE, kTRUE);
594   ZDCTDCTrigger(aEsd,kASide,kFALSE,kFALSE,kTRUE);
595   ZDCTimeTrigger(aEsd,kTRUE);
596
597   AliESDZDC* zdcData = aEsd->GetESDZDC();
598   if (zdcData)
599   {
600     UInt_t quality = zdcData->GetESDQuality();
601     
602     // from Nora's presentation, general first physics meeting 16.10.09
603     static UInt_t zpc  = 0x20;
604     static UInt_t znc  = 0x10;
605     static UInt_t zem1 = 0x08;
606     static UInt_t zem2 = 0x04;
607     static UInt_t zpa  = 0x02;
608     static UInt_t zna  = 0x01;
609    
610     fHistZDC->Fill(1, quality & zna);
611     fHistZDC->Fill(2, quality & zpa);
612     fHistZDC->Fill(3, quality & zem2);
613     fHistZDC->Fill(4, quality & zem1);
614     fHistZDC->Fill(5, quality & znc);
615     fHistZDC->Fill(6, quality & zpc);
616   }
617   else
618   {
619     fHistZDC->Fill(-1);
620     AliError("AliESDZDC not available");
621   }
622   
623   if (fDoFMD) {
624     fHistFMDA->Fill(FMDHitCombinations(aEsd, kASide, kTRUE));
625     fHistFMDC->Fill(FMDHitCombinations(aEsd, kCSide, kTRUE));
626   }
627 }
628   
629 void AliTriggerAnalysis::FillTriggerClasses(const AliESDEvent* aEsd)
630 {
631   // fills trigger classes map
632   
633   TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(aEsd->GetFiredTriggerClasses().Data()));
634   if (!count)
635   {
636     count = new TParameter<Long64_t>(aEsd->GetFiredTriggerClasses(), 0);
637     fTriggerClasses->Add(new TObjString(aEsd->GetFiredTriggerClasses().Data()), count);
638   }
639   count->SetVal(count->GetVal() + 1);
640 }
641
642 Int_t AliTriggerAnalysis::SSDClusters(const AliESDEvent* aEsd)
643 {
644   // returns the number of clusters in the SSD
645   const AliMultiplicity* mult = aEsd->GetMultiplicity();
646   Int_t clusters = mult->GetNumberOfITSClusters(4)+mult->GetNumberOfITSClusters(5);
647   return clusters;
648 }
649
650
651 Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHists, Int_t layer)
652 {
653   // returns the number of fired chips in the SPD
654   //
655   // origin = 0 --> aEsd->GetMultiplicity()->GetNumberOfFiredChips() (filled from clusters)
656   // origin = 1 --> aEsd->GetMultiplicity()->TestFastOrFiredChips() (from hardware bits)
657   // layer  = 0 --> both layers
658   // layer  = 1 --> inner
659   // layer  = 2 --> outer
660   
661   const AliMultiplicity* mult = aEsd->GetMultiplicity();
662   if (!mult)
663   {
664     AliError("AliMultiplicity not available");
665     return -1;
666   }
667   
668   if (origin == 0){
669     if (layer == 0) 
670       return mult->GetNumberOfFiredChips(0) + mult->GetNumberOfFiredChips(1);
671
672     return mult->GetNumberOfFiredChips(layer-1); 
673   }
674     
675   if (origin == 1)
676   {
677     Int_t nChips = 0;
678     Int_t firstChip = 0;
679     Int_t lastChip  = 1200;
680     if(layer == 1)
681       lastChip  = 400;
682     if(layer == 2)
683       firstChip = 400;
684
685     for (Int_t i=0; i<1200; i++)
686     {
687       if (mult->TestFastOrFiredChips(i) == kTRUE)
688       {
689         // efficiency simulation (if enabled)
690         if (fSPDGFOEfficiency)
691         {
692           if (gRandom->Uniform() > fSPDGFOEfficiency->GetBinContent(i+1))
693             continue;
694         }
695         
696         nChips++;
697         if (fillHists)
698           fHistFiredBitsSPD->Fill(i);
699       }
700     }
701     return nChips;
702   }
703   
704   return -1;
705 }
706
707 Bool_t AliTriggerAnalysis::SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin)
708 {
709   // Returns if the SPD gave a global Fast OR trigger
710   
711   Int_t firedChips = SPDFiredChips(aEsd, origin);
712   
713   if (firedChips >= fSPDGFOThreshold)
714     return kTRUE;
715   return kFALSE;
716 }
717
718 AliTriggerAnalysis::V0Decision AliTriggerAnalysis::V0Trigger(const AliESDEvent* aEsd, AliceSide side, Bool_t online, Bool_t fillHists)
719 {
720   // Returns the V0 trigger decision in V0A | V0C
721   //
722   // Returns kV0Fake if the calculated average time is in a window where neither BB nor BG is expected. 
723   // The rate of such triggers can be used to estimate the background. Note that the rate has to be 
724   // rescaled with the size of the windows (numerical values see below in the code)
725   //
726   // argument 'online' is used as a switch between online and offline trigger algorithms
727   //
728   // Based on an algorithm by Cvetan Cheshkov
729
730   AliESDVZERO* esdV0 = aEsd->GetVZEROData();
731   if (!esdV0)
732   {
733     AliError("AliESDVZERO not available");
734     return kV0Invalid;
735   }
736   AliDebug(2,Form("In V0Trigger: %f %f",esdV0->GetV0ATime(),esdV0->GetV0CTime()));
737
738   Int_t begin = -1;
739   Int_t end = -1;
740   
741   if (side == kASide)
742   {
743     begin = 32;
744     end = 64;
745   } 
746   else if (side == kCSide)
747   {
748     begin = 0;
749     end = 32;
750   }
751   else
752     return kV0Invalid;
753     
754    if (esdV0->TestBit(AliESDVZERO::kDecisionFilled)) {
755     if (online) {
756       if (esdV0->TestBit(AliESDVZERO::kOnlineBitsFilled)) {
757         for (Int_t i = begin; i < end; ++i) {
758           if (esdV0->GetBBFlag(i)) return kV0BB;
759         }
760         for (Int_t i = begin; i < end; ++i) {
761           if (esdV0->GetBGFlag(i)) return kV0BG;
762         }
763         return kV0Empty;
764       }
765       else {
766         AliWarning("V0 online trigger analysis is not yet available!");
767         return kV0BB;
768       }
769     }
770     else {
771
772       if (fillHists) {
773         if (side == kASide && fHistV0A)
774           fHistV0A->Fill(esdV0->GetV0ATime());
775         if (side == kCSide && fHistV0C)
776           fHistV0C->Fill(esdV0->GetV0CTime());
777       }
778
779       if (side == kASide) return (V0Decision)esdV0->GetV0ADecision();
780       else if (side == kCSide) return (V0Decision)esdV0->GetV0CDecision();
781       else return kV0Invalid;
782     }
783   }
784
785   Float_t time = 0;
786   Float_t weight = 0;
787   if (fMC)
788   {
789     Int_t runRange;
790     if (aEsd->GetRunNumber() <= 104803) runRange = 0;
791     else if (aEsd->GetRunNumber() <= 104876) runRange = 1;
792     else runRange = 2;
793
794     Float_t factors[3][64] = {
795       // runs: 104792-104803
796       {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},
797       // runs: 104841-104876
798       {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},
799       // runs: 104890-92
800       {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}
801     };
802     Float_t dA = 77.4 - 11.0;
803     Float_t dC = 77.4 - 2.9;
804     // Time misalignment
805     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};
806     Float_t dA2 = 2.8, dC2 = 3.3;
807
808     if (online) {
809       for (Int_t i = begin; i < end; ++i) {
810         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
811         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
812         if (esdV0->GetTime(i) >= 1e-6 &&
813             tempTime > fV0HwWinLow && tempTime < fV0HwWinHigh &&
814             tempAdc > fV0HwAdcThr)
815           return kV0BB;
816       }
817       return kV0Empty;
818     }
819     else {
820       for (Int_t i = begin; i < end; ++i) {
821         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
822         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA : esdV0->GetTime(i)+dC;
823         Float_t tempRawTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
824         if (esdV0->GetTime(i) >= 1e-6 &&
825             tempRawTime < 125.0 &&
826             tempAdc > fV0AdcThr) {
827           weight += 1.0;
828           time += tempTime;
829         }
830       }
831     }
832   }
833   else {
834     if (online) {
835       for (Int_t i = begin; i < end; ++i) {
836         if (esdV0->GetTime(i) >= 1e-6 &&
837             esdV0->GetTime(i) > fV0HwWinLow && esdV0->GetTime(i) < fV0HwWinHigh &&
838             esdV0->GetAdc(i) > fV0HwAdcThr)
839           return kV0BB;
840       }
841       return kV0Empty;
842     }
843     else {
844       for (Int_t i = begin; i < end; ++i) {
845         if (esdV0->GetTime(i) > 1e-6 && esdV0->GetAdc(i) > fV0AdcThr) {
846           Float_t correctedTime = V0CorrectLeadingTime(i, esdV0->GetTime(i), esdV0->GetAdc(i),aEsd->GetRunNumber());
847           Float_t timeWeight = V0LeadingTimeWeight(esdV0->GetAdc(i));
848           time += correctedTime*timeWeight;
849             
850           weight += timeWeight;
851         }
852       }
853     }
854   }
855
856   if (weight > 0) 
857     time /= weight;
858   time += fV0TimeOffset;
859
860   if (fillHists)
861   {
862     if (side == kASide && fHistV0A)
863       fHistV0A->Fill(time);
864     if (side == kCSide && fHistV0C)
865       fHistV0C->Fill(time);
866   }
867   
868   if (side == kASide)
869   {
870     if (time > 68 && time < 100)
871       return kV0BB;
872     if (time > 54 && time < 57.5) 
873       return kV0BG;
874     if (time > 57.5 && time < 68)
875       return kV0Fake;
876   }
877   
878   if (side == kCSide)
879   {
880     if (time > 75.5 && time < 100)
881       return kV0BB;
882     if (time > 69.5 && time < 73)
883       return kV0BG; 
884     if (time > 55 && time < 69.5)
885       return kV0Fake;
886   }
887   
888   return kV0Empty;
889 }
890
891 Float_t AliTriggerAnalysis::V0CorrectLeadingTime(Int_t i, Float_t time, Float_t adc, Int_t runNumber) const
892 {
893   // Correct for slewing and align the channels
894   //
895   // Authors: Cvetan Cheshkov / Raphael Tieulent
896
897   if (time == 0) return 0;
898
899   // Time alignment
900   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};
901
902   if(runNumber < 106031)
903     time -= timeShift[i];
904
905   // Slewing correction
906   if (adc == 0) return time;
907
908   Float_t p1 = 1.57345e1;
909   Float_t p2 =-4.25603e-1;
910
911   if(runNumber >= 106031) adc *= (2.5/4.0);
912   return (time - p1*TMath::Power(adc,p2));
913 }
914
915 Float_t AliTriggerAnalysis::V0LeadingTimeWeight(Float_t adc) const
916 {
917   if (adc < 1e-6) return 0;
918
919   Float_t p1 = 40.211;
920   Float_t p2 =-4.25603e-1;
921   Float_t p3 = 0.5646;
922
923   return 1./(p1*p1*TMath::Power(adc,2.*(p2-1.))+p3*p3);
924 }
925
926
927 Bool_t AliTriggerAnalysis::ZDCTDCTrigger(const AliESDEvent* aEsd, AliceSide side, Bool_t useZN, Bool_t useZP, Bool_t fillHists) const
928 {
929   // Returns if ZDC triggered, based on TDC information 
930   
931   AliESDZDC *esdZDC = aEsd->GetESDZDC();
932
933   Bool_t zdcNA = kFALSE;
934   Bool_t zdcNC = kFALSE;
935   Bool_t zdcPA = kFALSE;
936   Bool_t zdcPC = kFALSE;
937
938   if (fMC) {
939         // If it's MC, we use the energy
940     Double_t minEnergy = 0;
941     Double_t  zNCEnergy = esdZDC->GetZDCN1Energy();
942     Double_t  zPCEnergy = esdZDC->GetZDCP1Energy();
943     Double_t  zNAEnergy = esdZDC->GetZDCN2Energy();
944     Double_t  zPAEnergy = esdZDC->GetZDCP2Energy();
945     zdcNA = (zNAEnergy>minEnergy);
946     zdcNC = (zNCEnergy>minEnergy);
947     zdcPA = (zPAEnergy>minEnergy);
948     zdcPC = (zPCEnergy>minEnergy);
949
950   }
951   else {
952
953     Bool_t tdc[32] = {kFALSE};
954     for(Int_t itdc=0; itdc<32; itdc++){
955       for(Int_t i=0; i<4; i++){
956         if (esdZDC->GetZDCTDCData(itdc, i) != 0){
957           tdc[itdc] = kTRUE;
958         }
959       }
960       if(fillHists && tdc[itdc]) {
961         fHistTDCZDC->Fill(itdc);
962       }    
963     }
964     zdcNA = tdc[12];
965     zdcNC = tdc[10];
966     zdcPA = tdc[13];
967     zdcPC = tdc[11];
968   }
969
970   if (side == kASide) return ((useZP&&zdcPA) || (useZN&&zdcNA)); 
971   if (side == kCSide) return ((useZP&&zdcPC) || (useZN&&zdcNC)); 
972   return kFALSE;
973 }
974
975 Bool_t AliTriggerAnalysis::ZDCTimeTrigger(const AliESDEvent *aEsd, Bool_t fillHists) const
976 {
977   // This method implements a selection
978   // based on the timing in both sides of zdcN
979   // It can be used in order to eliminate
980   // parasitic collisions
981   Bool_t zdcAccept = kFALSE;
982   AliESDZDC *esdZDC = aEsd->GetESDZDC();
983
984   if(fMC) {
985      UInt_t esdFlag =  esdZDC->GetESDQuality();
986
987      Bool_t znaFired=kFALSE, zpaFired=kFALSE;
988      Bool_t zem1Fired=kFALSE, zem2Fired=kFALSE;
989      Bool_t zncFired=kFALSE, zpcFired=kFALSE;
990      //
991      // **** Trigger patterns
992      if((esdFlag & 0x00000001) == 0x00000001) znaFired=kTRUE;
993      if((esdFlag & 0x00000002) == 0x00000002) zpaFired=kTRUE;
994      if((esdFlag & 0x00000004) == 0x00000004) zem1Fired=kTRUE;
995      if((esdFlag & 0x00000008) == 0x00000008) zem2Fired=kTRUE;
996      if((esdFlag & 0x00000010) == 0x00000010) zncFired=kTRUE;
997      if((esdFlag & 0x00000020) == 0x00000020) zpcFired=kTRUE;
998      zdcAccept = (znaFired | zncFired);
999   }
1000   else {
1001     const Float_t refSum = -568.5;
1002     const Float_t refDelta = -2.1;
1003     const Float_t sigmaSum = 3.25;
1004     const Float_t sigmaDelta = 2.25;
1005     for(Int_t i = 0; i < 4; ++i) {
1006       if (esdZDC->GetZDCTDCData(10,i) != 0) {
1007         Float_t tdcC = 0.025*(esdZDC->GetZDCTDCData(10,i)-esdZDC->GetZDCTDCData(14,i)); 
1008         for(Int_t j = 0; j < 4; ++j) {
1009           if (esdZDC->GetZDCTDCData(12,j) != 0) {
1010             Float_t tdcA = 0.025*(esdZDC->GetZDCTDCData(12,j)-esdZDC->GetZDCTDCData(14,j));
1011             if(fillHists) fHistTimeZDC->Fill(tdcC-tdcA,tdcC+tdcA);
1012             if (((tdcC-tdcA-refDelta)*(tdcC-tdcA-refDelta)/(sigmaDelta*sigmaDelta) +
1013                  (tdcC+tdcA-refSum)*(tdcC+tdcA-refSum)/(sigmaSum*sigmaSum))< 1.0)
1014               zdcAccept = kTRUE;
1015           }
1016         }
1017       }
1018     }
1019   }
1020   return zdcAccept;
1021 }
1022
1023 Bool_t AliTriggerAnalysis::ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) const
1024 {
1025   // Returns if ZDC triggered
1026   
1027   AliESDZDC* zdcData = aEsd->GetESDZDC();
1028   if (!zdcData)
1029   {
1030     AliError("AliESDZDC not available");
1031     return kFALSE;
1032   }
1033   
1034   UInt_t quality = zdcData->GetESDQuality();
1035   
1036   // from Nora's presentation, general first physics meeting 16.10.09
1037   static UInt_t zpc  = 0x20;
1038   static UInt_t znc  = 0x10;
1039   static UInt_t zem1 = 0x08;
1040   static UInt_t zem2 = 0x04;
1041   static UInt_t zpa  = 0x02;
1042   static UInt_t zna  = 0x01;
1043   
1044   if (side == kASide && ((quality & zpa) || (quality & zna)))
1045     return kTRUE;
1046   if (side == kCentralBarrel && ((quality & zem1) || (quality & zem2)))
1047     return kTRUE;
1048   if (side == kCSide && ((quality & zpc) || (quality & znc)))
1049     return kTRUE;
1050   
1051   return kFALSE;
1052 }
1053
1054 Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists)
1055 {
1056   // returns number of hit combinations agove threshold
1057   //
1058   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
1059
1060   if (!fDoFMD)
1061     return -1;
1062
1063   // Workaround for AliESDEvent::GetFMDData is not const!
1064   const AliESDFMD* fmdData = (const_cast<AliESDEvent*>(aEsd))->GetFMDData();
1065   if (!fmdData)
1066   {
1067     AliError("AliESDFMD not available");
1068     return -1;
1069   }
1070
1071   Int_t detFrom = (side == kASide) ? 1 : 3;
1072   Int_t detTo   = (side == kASide) ? 2 : 3;
1073
1074   Int_t triggers = 0;
1075   Float_t totalMult = 0;
1076   for (UShort_t det=detFrom;det<=detTo;det++) {
1077     Int_t nRings = (det == 1 ? 1 : 2);
1078     for (UShort_t ir = 0; ir < nRings; ir++) {
1079       Char_t   ring = (ir == 0 ? 'I' : 'O');
1080       UShort_t nsec = (ir == 0 ? 20  : 40);
1081       UShort_t nstr = (ir == 0 ? 512 : 256);
1082       for (UShort_t sec =0; sec < nsec;  sec++) {
1083         for (UShort_t strip = 0; strip < nstr; strip++) {
1084           Float_t mult = fmdData->Multiplicity(det,ring,sec,strip);
1085           if (mult == AliESDFMD::kInvalidMult) continue;
1086           
1087           if (fillHists)
1088             fHistFMDSingle->Fill(mult);
1089           
1090           if (mult > fFMDLowCut)
1091             totalMult = totalMult + mult;
1092           else
1093           {
1094             if (totalMult > fFMDHitCut)
1095               triggers++;
1096               
1097             if (fillHists)
1098               fHistFMDSum->Fill(totalMult);
1099               
1100             totalMult = 0;
1101           }
1102         }
1103       }
1104     }
1105   }
1106   
1107   return triggers;
1108 }
1109
1110 Bool_t AliTriggerAnalysis::FMDTrigger(const AliESDEvent* aEsd, AliceSide side)
1111 {
1112   // Returns if the FMD triggered
1113   //
1114   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
1115
1116   Int_t triggers = FMDHitCombinations(aEsd, side, kFALSE);
1117   
1118   if (triggers > 0)
1119     return kTRUE;
1120     
1121   return kFALSE;
1122 }
1123
1124 Long64_t AliTriggerAnalysis::Merge(TCollection* list)
1125 {
1126   // Merge a list of AliMultiplicityCorrection objects with this (needed for
1127   // PROOF).
1128   // Returns the number of merged objects (including this).
1129
1130   if (!list)
1131     return 0;
1132
1133   if (list->IsEmpty())
1134     return 1;
1135
1136   TIterator* iter = list->MakeIterator();
1137   TObject* obj;
1138
1139   // collections of all histograms
1140   const Int_t nHists = 11;
1141   TList collections[nHists];
1142
1143   Int_t count = 0;
1144   while ((obj = iter->Next())) {
1145
1146     AliTriggerAnalysis* entry = dynamic_cast<AliTriggerAnalysis*> (obj);
1147     if (entry == 0) 
1148       continue;
1149
1150     Int_t n = 0;
1151     collections[n++].Add(entry->fHistV0A);
1152     collections[n++].Add(entry->fHistV0C);
1153     collections[n++].Add(entry->fHistZDC);
1154     collections[n++].Add(entry->fHistTDCZDC);
1155     collections[n++].Add(entry->fHistTimeZDC);
1156     collections[n++].Add(entry->fHistFMDA);
1157     collections[n++].Add(entry->fHistFMDC);
1158     collections[n++].Add(entry->fHistFMDSingle);
1159     collections[n++].Add(entry->fHistFMDSum);
1160     collections[n++].Add(entry->fHistBitsSPD);
1161     collections[n++].Add(entry->fHistFiredBitsSPD);
1162
1163     // merge fTriggerClasses
1164     TIterator* iter2 = entry->fTriggerClasses->MakeIterator();
1165     TObjString* obj2 = 0;
1166     while ((obj2 = dynamic_cast<TObjString*> (iter2->Next())))
1167     {
1168       TParameter<Long64_t>* param2 = dynamic_cast<TParameter<Long64_t>*> (entry->fTriggerClasses->GetValue(obj2));
1169       
1170       TParameter<Long64_t>* param1 = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj2));
1171       if (param1)
1172       {
1173         param1->SetVal(param1->GetVal() + param2->GetVal());
1174       }
1175       else
1176       {
1177         param1 = dynamic_cast<TParameter<Long64_t>*> (param2->Clone());
1178         fTriggerClasses->Add(new TObjString(obj2->String()), param1);
1179       }
1180     }
1181     
1182     delete iter2;
1183   
1184     count++;
1185   }
1186
1187   Int_t n = 0;
1188   fHistV0A->Merge(&collections[n++]);
1189   fHistV0C->Merge(&collections[n++]);
1190   fHistZDC->Merge(&collections[n++]);
1191   fHistTDCZDC->Merge(&collections[n++]);
1192   fHistTimeZDC->Merge(&collections[n++]);
1193   fHistFMDA->Merge(&collections[n++]);
1194   fHistFMDC->Merge(&collections[n++]);
1195   fHistFMDSingle->Merge(&collections[n++]);
1196   fHistFMDSum->Merge(&collections[n++]);
1197   fHistBitsSPD->Merge(&collections[n++]);
1198   fHistFiredBitsSPD->Merge(&collections[n++]);
1199   
1200   delete iter;
1201
1202   return count+1;
1203 }
1204
1205 void AliTriggerAnalysis::SaveHistograms() const
1206 {
1207   // write histograms to current directory
1208   
1209   if (!fHistBitsSPD)
1210     return;
1211     
1212   if (fHistBitsSPD) {
1213     fHistBitsSPD->Write();
1214     fHistBitsSPD->ProjectionX();
1215     fHistBitsSPD->ProjectionY();
1216   }
1217   else Printf("Cannot save fHistBitsSPD");
1218   if (fHistFiredBitsSPD) fHistFiredBitsSPD->Write();
1219   else Printf("Cannot save fHistFiredBitsSPD");
1220   if (fHistV0A) fHistV0A->Write();
1221   else Printf("Cannot save fHistV0A");
1222   if (fHistV0C) fHistV0C->Write();
1223   else Printf("Cannot save fHistV0C");
1224   if (fHistZDC) fHistZDC->Write();
1225   else Printf("Cannot save fHistZDC");
1226   if (fHistTDCZDC) fHistTDCZDC->Write();
1227   else Printf("Cannot save fHistTDCZDC");
1228   if (fHistTimeZDC) fHistTimeZDC->Write();
1229   else Printf("Cannot save fHistTimeZDC");
1230   if (fHistFMDA) fHistFMDA->Write();
1231   else Printf("Cannot save fHistFMDA");
1232   if (fHistFMDC) fHistFMDC->Write();
1233   else Printf("Cannot save fHistFMDC");
1234   if (fHistFMDSingle) fHistFMDSingle->Write();
1235   else Printf("Cannot save fHistFMDSingle");
1236   if (fHistFMDSum) fHistFMDSum->Write();
1237   else Printf("Cannot save fHistFMDSum");
1238   if (fSPDGFOEfficiency) fSPDGFOEfficiency->Write("fSPDGFOEfficiency");
1239   else Printf("Cannot save fSPDGFOEfficiency");
1240   
1241   fTriggerClasses->Write("fTriggerClasses", TObject::kSingleKey);
1242 }
1243
1244 void AliTriggerAnalysis::PrintTriggerClasses() const
1245 {
1246   // print trigger classes
1247   
1248   Printf("Trigger Classes:");
1249   
1250   Printf("Event count for trigger combinations:");
1251   
1252   TMap singleTrigger;
1253   singleTrigger.SetOwner();
1254   
1255   TIterator* iter = fTriggerClasses->MakeIterator();
1256   TObjString* obj = 0;
1257   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
1258   {
1259     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj));
1260     
1261     Printf(" %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1262     
1263     TObjArray* tokens = obj->String().Tokenize(" ");
1264     for (Int_t i=0; i<tokens->GetEntries(); i++)
1265     {
1266       TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(((TObjString*) tokens->At(i))->String().Data()));
1267       if (!count)
1268       {
1269         count = new TParameter<Long64_t>(((TObjString*) tokens->At(i))->String().Data(), 0);
1270         singleTrigger.Add(new TObjString(((TObjString*) tokens->At(i))->String().Data()), count);
1271       }
1272       count->SetVal(count->GetVal() + param->GetVal());
1273     }
1274     
1275     delete tokens;
1276   }
1277   delete iter;
1278   
1279   Printf("Event count for single trigger:");
1280   
1281   iter = singleTrigger.MakeIterator();
1282   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
1283   {
1284     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(obj));
1285     
1286     Printf("  %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1287   }
1288   delete iter;
1289   
1290   singleTrigger.DeleteAll();
1291 }