]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliTriggerAnalysis.cxx
move printf to debug
[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        delete tcltokens;
540        tcltokens = 0;
541        //      tcltokens->Delete();
542        return kFALSE;
543     }
544   }
545   delete tcltokens;
546   tcltokens = 0;
547        //  tcltokens->Delete();
548   return tcl;
549 }
550
551 Bool_t AliTriggerAnalysis::IsInputFired(const AliESDEvent* aEsd, Char_t level, UInt_t input) const
552 {
553   // Checks trigger input of any level
554   
555   switch (level)
556   {
557     case '0': return IsL0InputFired(aEsd,input);
558     case '1': return IsL1InputFired(aEsd,input);
559     case '2': return IsL2InputFired(aEsd,input);
560     default:
561       AliError(Form("Wrong level %i",level));
562       return kFALSE;
563   }
564 }
565
566 Bool_t AliTriggerAnalysis::IsL0InputFired(const AliESDEvent* aEsd, UInt_t input) const 
567 {
568   // Checks if corresponding bit in mask is on
569   
570   UInt_t inpmask = aEsd->GetHeader()->GetL0TriggerInputs();
571   return (inpmask & (1<<(input-1)));
572 }
573
574 Bool_t AliTriggerAnalysis::IsL1InputFired(const AliESDEvent* aEsd, UInt_t input) const
575 {
576   // Checks if corresponding bit in mask is on
577   
578   UInt_t inpmask = aEsd->GetHeader()->GetL1TriggerInputs();
579   return (inpmask & (1<<(input-1)));
580 }
581
582 Bool_t AliTriggerAnalysis::IsL2InputFired(const AliESDEvent* aEsd, UInt_t input) const 
583 {
584   // Checks if corresponding bit in mask is on
585   
586   UInt_t inpmask = aEsd->GetHeader()->GetL2TriggerInputs();
587   return (inpmask & (1<<(input-1)));
588 }
589
590 void AliTriggerAnalysis::FillHistograms(const AliESDEvent* aEsd) 
591 {
592   // fills the histograms with the info from the ESD
593   
594   fHistBitsSPD->Fill(SPDFiredChips(aEsd, 0), SPDFiredChips(aEsd, 1, kTRUE));
595   
596   V0Trigger(aEsd, kASide, kFALSE, kTRUE);
597   V0Trigger(aEsd, kCSide, kFALSE, kTRUE);
598   ZDCTDCTrigger(aEsd,kASide,kFALSE,kFALSE,kTRUE);
599   ZDCTimeTrigger(aEsd,kTRUE);
600
601   AliESDZDC* zdcData = aEsd->GetESDZDC();
602   if (zdcData)
603   {
604     UInt_t quality = zdcData->GetESDQuality();
605     
606     // from Nora's presentation, general first physics meeting 16.10.09
607     static UInt_t zpc  = 0x20;
608     static UInt_t znc  = 0x10;
609     static UInt_t zem1 = 0x08;
610     static UInt_t zem2 = 0x04;
611     static UInt_t zpa  = 0x02;
612     static UInt_t zna  = 0x01;
613    
614     fHistZDC->Fill(1, quality & zna);
615     fHistZDC->Fill(2, quality & zpa);
616     fHistZDC->Fill(3, quality & zem2);
617     fHistZDC->Fill(4, quality & zem1);
618     fHistZDC->Fill(5, quality & znc);
619     fHistZDC->Fill(6, quality & zpc);
620   }
621   else
622   {
623     fHistZDC->Fill(-1);
624     AliError("AliESDZDC not available");
625   }
626   
627   if (fDoFMD) {
628     fHistFMDA->Fill(FMDHitCombinations(aEsd, kASide, kTRUE));
629     fHistFMDC->Fill(FMDHitCombinations(aEsd, kCSide, kTRUE));
630   }
631 }
632   
633 void AliTriggerAnalysis::FillTriggerClasses(const AliESDEvent* aEsd)
634 {
635   // fills trigger classes map
636   
637   TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(aEsd->GetFiredTriggerClasses().Data()));
638   if (!count)
639   {
640     count = new TParameter<Long64_t>(aEsd->GetFiredTriggerClasses(), 0);
641     fTriggerClasses->Add(new TObjString(aEsd->GetFiredTriggerClasses().Data()), count);
642   }
643   count->SetVal(count->GetVal() + 1);
644 }
645
646 Int_t AliTriggerAnalysis::SSDClusters(const AliESDEvent* aEsd)
647 {
648   // returns the number of clusters in the SSD
649   const AliMultiplicity* mult = aEsd->GetMultiplicity();
650   Int_t clusters = mult->GetNumberOfITSClusters(4)+mult->GetNumberOfITSClusters(5);
651   return clusters;
652 }
653
654
655 Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHists, Int_t layer)
656 {
657   // returns the number of fired chips in the SPD
658   //
659   // origin = 0 --> aEsd->GetMultiplicity()->GetNumberOfFiredChips() (filled from clusters)
660   // origin = 1 --> aEsd->GetMultiplicity()->TestFastOrFiredChips() (from hardware bits)
661   // layer  = 0 --> both layers
662   // layer  = 1 --> inner
663   // layer  = 2 --> outer
664   
665   const AliMultiplicity* mult = aEsd->GetMultiplicity();
666   if (!mult)
667   {
668     AliError("AliMultiplicity not available");
669     return -1;
670   }
671   
672   if (origin == 0){
673     if (layer == 0) 
674       return mult->GetNumberOfFiredChips(0) + mult->GetNumberOfFiredChips(1);
675
676     return mult->GetNumberOfFiredChips(layer-1); 
677   }
678     
679   if (origin == 1)
680   {
681     Int_t nChips = 0;
682     Int_t firstChip = 0;
683     Int_t lastChip  = 1200;
684     if(layer == 1)
685       lastChip  = 400;
686     if(layer == 2)
687       firstChip = 400;
688
689     for (Int_t i=firstChip; i<lastChip; i++)
690     {
691       if (mult->TestFastOrFiredChips(i) == kTRUE)
692       {
693         // efficiency simulation (if enabled)
694         if (fSPDGFOEfficiency)
695         {
696           if (gRandom->Uniform() > fSPDGFOEfficiency->GetBinContent(i+1))
697             continue;
698         }
699         
700         nChips++;
701         if (fillHists)
702           fHistFiredBitsSPD->Fill(i);
703       }
704     }
705     return nChips;
706   }
707   
708   return -1;
709 }
710
711 Bool_t AliTriggerAnalysis::SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin)
712 {
713   // Returns if the SPD gave a global Fast OR trigger
714   
715   Int_t firedChips = SPDFiredChips(aEsd, origin);
716   
717   if (firedChips >= fSPDGFOThreshold)
718     return kTRUE;
719   return kFALSE;
720 }
721
722 AliTriggerAnalysis::V0Decision AliTriggerAnalysis::V0Trigger(const AliESDEvent* aEsd, AliceSide side, Bool_t online, Bool_t fillHists)
723 {
724   // Returns the V0 trigger decision in V0A | V0C
725   //
726   // Returns kV0Fake if the calculated average time is in a window where neither BB nor BG is expected. 
727   // The rate of such triggers can be used to estimate the background. Note that the rate has to be 
728   // rescaled with the size of the windows (numerical values see below in the code)
729   //
730   // argument 'online' is used as a switch between online and offline trigger algorithms
731   //
732   // Based on an algorithm by Cvetan Cheshkov
733
734   AliESDVZERO* esdV0 = aEsd->GetVZEROData();
735   if (!esdV0)
736   {
737     AliError("AliESDVZERO not available");
738     return kV0Invalid;
739   }
740   AliDebug(2,Form("In V0Trigger: %f %f",esdV0->GetV0ATime(),esdV0->GetV0CTime()));
741
742   Int_t begin = -1;
743   Int_t end = -1;
744   
745   if (side == kASide)
746   {
747     begin = 32;
748     end = 64;
749   } 
750   else if (side == kCSide)
751   {
752     begin = 0;
753     end = 32;
754   }
755   else
756     return kV0Invalid;
757     
758    if (esdV0->TestBit(AliESDVZERO::kDecisionFilled)) {
759     if (online) {
760       if (esdV0->TestBit(AliESDVZERO::kOnlineBitsFilled)) {
761         for (Int_t i = begin; i < end; ++i) {
762           if (esdV0->GetBBFlag(i)) return kV0BB;
763         }
764         for (Int_t i = begin; i < end; ++i) {
765           if (esdV0->GetBGFlag(i)) return kV0BG;
766         }
767         return kV0Empty;
768       }
769       else {
770         AliWarning("V0 online trigger analysis is not yet available!");
771         return kV0BB;
772       }
773     }
774     else {
775
776       if (fillHists) {
777         if (side == kASide && fHistV0A)
778           fHistV0A->Fill(esdV0->GetV0ATime());
779         if (side == kCSide && fHistV0C)
780           fHistV0C->Fill(esdV0->GetV0CTime());
781       }
782
783       if (side == kASide) return (V0Decision)esdV0->GetV0ADecision();
784       else if (side == kCSide) return (V0Decision)esdV0->GetV0CDecision();
785       else return kV0Invalid;
786     }
787   }
788
789   Float_t time = 0;
790   Float_t weight = 0;
791   if (fMC)
792   {
793     Int_t runRange;
794     if (aEsd->GetRunNumber() <= 104803) runRange = 0;
795     else if (aEsd->GetRunNumber() <= 104876) runRange = 1;
796     else runRange = 2;
797
798     Float_t factors[3][64] = {
799       // runs: 104792-104803
800       {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},
801       // runs: 104841-104876
802       {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},
803       // runs: 104890-92
804       {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}
805     };
806     Float_t dA = 77.4 - 11.0;
807     Float_t dC = 77.4 - 2.9;
808     // Time misalignment
809     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};
810     Float_t dA2 = 2.8, dC2 = 3.3;
811
812     if (online) {
813       for (Int_t i = begin; i < end; ++i) {
814         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
815         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
816         if (esdV0->GetTime(i) >= 1e-6 &&
817             tempTime > fV0HwWinLow && tempTime < fV0HwWinHigh &&
818             tempAdc > fV0HwAdcThr)
819           return kV0BB;
820       }
821       return kV0Empty;
822     }
823     else {
824       for (Int_t i = begin; i < end; ++i) {
825         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
826         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA : esdV0->GetTime(i)+dC;
827         Float_t tempRawTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
828         if (esdV0->GetTime(i) >= 1e-6 &&
829             tempRawTime < 125.0 &&
830             tempAdc > fV0AdcThr) {
831           weight += 1.0;
832           time += tempTime;
833         }
834       }
835     }
836   }
837   else {
838     if (online) {
839       for (Int_t i = begin; i < end; ++i) {
840         if (esdV0->GetTime(i) >= 1e-6 &&
841             esdV0->GetTime(i) > fV0HwWinLow && esdV0->GetTime(i) < fV0HwWinHigh &&
842             esdV0->GetAdc(i) > fV0HwAdcThr)
843           return kV0BB;
844       }
845       return kV0Empty;
846     }
847     else {
848       for (Int_t i = begin; i < end; ++i) {
849         if (esdV0->GetTime(i) > 1e-6 && esdV0->GetAdc(i) > fV0AdcThr) {
850           Float_t correctedTime = V0CorrectLeadingTime(i, esdV0->GetTime(i), esdV0->GetAdc(i),aEsd->GetRunNumber());
851           Float_t timeWeight = V0LeadingTimeWeight(esdV0->GetAdc(i));
852           time += correctedTime*timeWeight;
853             
854           weight += timeWeight;
855         }
856       }
857     }
858   }
859
860   if (weight > 0) 
861     time /= weight;
862   time += fV0TimeOffset;
863
864   if (fillHists)
865   {
866     if (side == kASide && fHistV0A)
867       fHistV0A->Fill(time);
868     if (side == kCSide && fHistV0C)
869       fHistV0C->Fill(time);
870   }
871   
872   if (side == kASide)
873   {
874     if (time > 68 && time < 100)
875       return kV0BB;
876     if (time > 54 && time < 57.5) 
877       return kV0BG;
878     if (time > 57.5 && time < 68)
879       return kV0Fake;
880   }
881   
882   if (side == kCSide)
883   {
884     if (time > 75.5 && time < 100)
885       return kV0BB;
886     if (time > 69.5 && time < 73)
887       return kV0BG; 
888     if (time > 55 && time < 69.5)
889       return kV0Fake;
890   }
891   
892   return kV0Empty;
893 }
894
895 Float_t AliTriggerAnalysis::V0CorrectLeadingTime(Int_t i, Float_t time, Float_t adc, Int_t runNumber) const
896 {
897   // Correct for slewing and align the channels
898   //
899   // Authors: Cvetan Cheshkov / Raphael Tieulent
900
901   if (time == 0) return 0;
902
903   // Time alignment
904   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};
905
906   if(runNumber < 106031)
907     time -= timeShift[i];
908
909   // Slewing correction
910   if (adc == 0) return time;
911
912   Float_t p1 = 1.57345e1;
913   Float_t p2 =-4.25603e-1;
914
915   if(runNumber >= 106031) adc *= (2.5/4.0);
916   return (time - p1*TMath::Power(adc,p2));
917 }
918
919 Float_t AliTriggerAnalysis::V0LeadingTimeWeight(Float_t adc) const
920 {
921   if (adc < 1e-6) return 0;
922
923   Float_t p1 = 40.211;
924   Float_t p2 =-4.25603e-1;
925   Float_t p3 = 0.5646;
926
927   return 1./(p1*p1*TMath::Power(adc,2.*(p2-1.))+p3*p3);
928 }
929
930
931 Bool_t AliTriggerAnalysis::ZDCTDCTrigger(const AliESDEvent* aEsd, AliceSide side, Bool_t useZN, Bool_t useZP, Bool_t fillHists) const
932 {
933   // Returns if ZDC triggered, based on TDC information 
934   
935   AliESDZDC *esdZDC = aEsd->GetESDZDC();
936
937   Bool_t zdcNA = kFALSE;
938   Bool_t zdcNC = kFALSE;
939   Bool_t zdcPA = kFALSE;
940   Bool_t zdcPC = kFALSE;
941
942   if (fMC) {
943         // If it's MC, we use the energy
944     Double_t minEnergy = 0;
945     Double_t  zNCEnergy = esdZDC->GetZDCN1Energy();
946     Double_t  zPCEnergy = esdZDC->GetZDCP1Energy();
947     Double_t  zNAEnergy = esdZDC->GetZDCN2Energy();
948     Double_t  zPAEnergy = esdZDC->GetZDCP2Energy();
949     zdcNA = (zNAEnergy>minEnergy);
950     zdcNC = (zNCEnergy>minEnergy);
951     zdcPA = (zPAEnergy>minEnergy);
952     zdcPC = (zPCEnergy>minEnergy);
953
954   }
955   else {
956
957     Bool_t tdc[32] = {kFALSE};
958     for(Int_t itdc=0; itdc<32; itdc++){
959       for(Int_t i=0; i<4; i++){
960         if (esdZDC->GetZDCTDCData(itdc, i) != 0){
961           tdc[itdc] = kTRUE;
962         }
963       }
964       if(fillHists && tdc[itdc]) {
965         fHistTDCZDC->Fill(itdc);
966       }    
967     }
968     zdcNA = tdc[12];
969     zdcNC = tdc[10];
970     zdcPA = tdc[13];
971     zdcPC = tdc[11];
972   }
973
974   if (side == kASide) return ((useZP&&zdcPA) || (useZN&&zdcNA)); 
975   if (side == kCSide) return ((useZP&&zdcPC) || (useZN&&zdcNC)); 
976   return kFALSE;
977 }
978
979 Bool_t AliTriggerAnalysis::ZDCTimeTrigger(const AliESDEvent *aEsd, Bool_t fillHists) const
980 {
981   // This method implements a selection
982   // based on the timing in both sides of zdcN
983   // It can be used in order to eliminate
984   // parasitic collisions
985   Bool_t zdcAccept = kFALSE;
986   AliESDZDC *esdZDC = aEsd->GetESDZDC();
987
988   if(fMC) {
989      UInt_t esdFlag =  esdZDC->GetESDQuality();
990
991      Bool_t znaFired=kFALSE, zpaFired=kFALSE;
992      Bool_t zem1Fired=kFALSE, zem2Fired=kFALSE;
993      Bool_t zncFired=kFALSE, zpcFired=kFALSE;
994      //
995      // **** Trigger patterns
996      if((esdFlag & 0x00000001) == 0x00000001) znaFired=kTRUE;
997      if((esdFlag & 0x00000002) == 0x00000002) zpaFired=kTRUE;
998      if((esdFlag & 0x00000004) == 0x00000004) zem1Fired=kTRUE;
999      if((esdFlag & 0x00000008) == 0x00000008) zem2Fired=kTRUE;
1000      if((esdFlag & 0x00000010) == 0x00000010) zncFired=kTRUE;
1001      if((esdFlag & 0x00000020) == 0x00000020) zpcFired=kTRUE;
1002      zdcAccept = (znaFired | zncFired);
1003   }
1004   else {
1005     const Float_t refSum = -568.5;
1006     const Float_t refDelta = -2.1;
1007     const Float_t sigmaSum = 3.25;
1008     const Float_t sigmaDelta = 2.25;
1009     for(Int_t i = 0; i < 4; ++i) {
1010       if (esdZDC->GetZDCTDCData(10,i) != 0) {
1011         Float_t tdcC = 0.025*(esdZDC->GetZDCTDCData(10,i)-esdZDC->GetZDCTDCData(14,i)); 
1012         for(Int_t j = 0; j < 4; ++j) {
1013           if (esdZDC->GetZDCTDCData(12,j) != 0) {
1014             Float_t tdcA = 0.025*(esdZDC->GetZDCTDCData(12,j)-esdZDC->GetZDCTDCData(14,j));
1015             if(fillHists) fHistTimeZDC->Fill(tdcC-tdcA,tdcC+tdcA);
1016             if (((tdcC-tdcA-refDelta)*(tdcC-tdcA-refDelta)/(sigmaDelta*sigmaDelta) +
1017                  (tdcC+tdcA-refSum)*(tdcC+tdcA-refSum)/(sigmaSum*sigmaSum))< 1.0)
1018               zdcAccept = kTRUE;
1019           }
1020         }
1021       }
1022     }
1023   }
1024   return zdcAccept;
1025 }
1026
1027 Bool_t AliTriggerAnalysis::ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) const
1028 {
1029   // Returns if ZDC triggered
1030   
1031   AliESDZDC* zdcData = aEsd->GetESDZDC();
1032   if (!zdcData)
1033   {
1034     AliError("AliESDZDC not available");
1035     return kFALSE;
1036   }
1037   
1038   UInt_t quality = zdcData->GetESDQuality();
1039   
1040   // from Nora's presentation, general first physics meeting 16.10.09
1041   static UInt_t zpc  = 0x20;
1042   static UInt_t znc  = 0x10;
1043   static UInt_t zem1 = 0x08;
1044   static UInt_t zem2 = 0x04;
1045   static UInt_t zpa  = 0x02;
1046   static UInt_t zna  = 0x01;
1047   
1048   if (side == kASide && ((quality & zpa) || (quality & zna)))
1049     return kTRUE;
1050   if (side == kCentralBarrel && ((quality & zem1) || (quality & zem2)))
1051     return kTRUE;
1052   if (side == kCSide && ((quality & zpc) || (quality & znc)))
1053     return kTRUE;
1054   
1055   return kFALSE;
1056 }
1057
1058 Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists)
1059 {
1060   // returns number of hit combinations agove threshold
1061   //
1062   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
1063
1064   if (!fDoFMD)
1065     return -1;
1066
1067   // Workaround for AliESDEvent::GetFMDData is not const!
1068   const AliESDFMD* fmdData = (const_cast<AliESDEvent*>(aEsd))->GetFMDData();
1069   if (!fmdData)
1070   {
1071     AliError("AliESDFMD not available");
1072     return -1;
1073   }
1074
1075   Int_t detFrom = (side == kASide) ? 1 : 3;
1076   Int_t detTo   = (side == kASide) ? 2 : 3;
1077
1078   Int_t triggers = 0;
1079   Float_t totalMult = 0;
1080   for (UShort_t det=detFrom;det<=detTo;det++) {
1081     Int_t nRings = (det == 1 ? 1 : 2);
1082     for (UShort_t ir = 0; ir < nRings; ir++) {
1083       Char_t   ring = (ir == 0 ? 'I' : 'O');
1084       UShort_t nsec = (ir == 0 ? 20  : 40);
1085       UShort_t nstr = (ir == 0 ? 512 : 256);
1086       for (UShort_t sec =0; sec < nsec;  sec++) {
1087         for (UShort_t strip = 0; strip < nstr; strip++) {
1088           Float_t mult = fmdData->Multiplicity(det,ring,sec,strip);
1089           if (mult == AliESDFMD::kInvalidMult) continue;
1090           
1091           if (fillHists)
1092             fHistFMDSingle->Fill(mult);
1093           
1094           if (mult > fFMDLowCut)
1095             totalMult = totalMult + mult;
1096           else
1097           {
1098             if (totalMult > fFMDHitCut)
1099               triggers++;
1100               
1101             if (fillHists)
1102               fHistFMDSum->Fill(totalMult);
1103               
1104             totalMult = 0;
1105           }
1106         }
1107       }
1108     }
1109   }
1110   
1111   return triggers;
1112 }
1113
1114 Bool_t AliTriggerAnalysis::FMDTrigger(const AliESDEvent* aEsd, AliceSide side)
1115 {
1116   // Returns if the FMD triggered
1117   //
1118   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
1119
1120   Int_t triggers = FMDHitCombinations(aEsd, side, kFALSE);
1121   
1122   if (triggers > 0)
1123     return kTRUE;
1124     
1125   return kFALSE;
1126 }
1127
1128 Long64_t AliTriggerAnalysis::Merge(TCollection* list)
1129 {
1130   // Merge a list of AliMultiplicityCorrection objects with this (needed for
1131   // PROOF).
1132   // Returns the number of merged objects (including this).
1133
1134   if (!list)
1135     return 0;
1136
1137   if (list->IsEmpty())
1138     return 1;
1139
1140   TIterator* iter = list->MakeIterator();
1141   TObject* obj;
1142
1143   // collections of all histograms
1144   const Int_t nHists = 11;
1145   TList collections[nHists];
1146
1147   Int_t count = 0;
1148   while ((obj = iter->Next())) {
1149
1150     AliTriggerAnalysis* entry = dynamic_cast<AliTriggerAnalysis*> (obj);
1151     if (entry == 0) 
1152       continue;
1153
1154     Int_t n = 0;
1155     collections[n++].Add(entry->fHistV0A);
1156     collections[n++].Add(entry->fHistV0C);
1157     collections[n++].Add(entry->fHistZDC);
1158     collections[n++].Add(entry->fHistTDCZDC);
1159     collections[n++].Add(entry->fHistTimeZDC);
1160     collections[n++].Add(entry->fHistFMDA);
1161     collections[n++].Add(entry->fHistFMDC);
1162     collections[n++].Add(entry->fHistFMDSingle);
1163     collections[n++].Add(entry->fHistFMDSum);
1164     collections[n++].Add(entry->fHistBitsSPD);
1165     collections[n++].Add(entry->fHistFiredBitsSPD);
1166
1167     // merge fTriggerClasses
1168     TIterator* iter2 = entry->fTriggerClasses->MakeIterator();
1169     TObjString* obj2 = 0;
1170     while ((obj2 = dynamic_cast<TObjString*> (iter2->Next())))
1171     {
1172       TParameter<Long64_t>* param2 = static_cast<TParameter<Long64_t>*> (entry->fTriggerClasses->GetValue(obj2));
1173       
1174       TParameter<Long64_t>* param1 = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj2));
1175       if (param1)
1176       {
1177         param1->SetVal(param1->GetVal() + param2->GetVal());
1178       }
1179       else
1180       {
1181         param1 = dynamic_cast<TParameter<Long64_t>*> (param2->Clone());
1182         fTriggerClasses->Add(new TObjString(obj2->String()), param1);
1183       }
1184     }
1185     
1186     delete iter2;
1187   
1188     count++;
1189   }
1190
1191   Int_t n = 0;
1192   fHistV0A->Merge(&collections[n++]);
1193   fHistV0C->Merge(&collections[n++]);
1194   fHistZDC->Merge(&collections[n++]);
1195   fHistTDCZDC->Merge(&collections[n++]);
1196   if (fHistTimeZDC)
1197     fHistTimeZDC->Merge(&collections[n++]);
1198   else
1199     n++;
1200   fHistFMDA->Merge(&collections[n++]);
1201   fHistFMDC->Merge(&collections[n++]);
1202   fHistFMDSingle->Merge(&collections[n++]);
1203   fHistFMDSum->Merge(&collections[n++]);
1204   fHistBitsSPD->Merge(&collections[n++]);
1205   fHistFiredBitsSPD->Merge(&collections[n++]);
1206   
1207   delete iter;
1208
1209   return count+1;
1210 }
1211
1212 void AliTriggerAnalysis::SaveHistograms() const
1213 {
1214   // write histograms to current directory
1215   
1216   if (!fHistBitsSPD)
1217     return;
1218     
1219   if (fHistBitsSPD) {
1220     fHistBitsSPD->Write();
1221     fHistBitsSPD->ProjectionX();
1222     fHistBitsSPD->ProjectionY();
1223   }
1224   else Printf("Cannot save fHistBitsSPD");
1225   if (fHistFiredBitsSPD) fHistFiredBitsSPD->Write();
1226   else Printf("Cannot save fHistFiredBitsSPD");
1227   if (fHistV0A) fHistV0A->Write();
1228   else Printf("Cannot save fHistV0A");
1229   if (fHistV0C) fHistV0C->Write();
1230   else Printf("Cannot save fHistV0C");
1231   if (fHistZDC) fHistZDC->Write();
1232   else Printf("Cannot save fHistZDC");
1233   if (fHistTDCZDC) fHistTDCZDC->Write();
1234   else Printf("Cannot save fHistTDCZDC");
1235   if (fHistTimeZDC) fHistTimeZDC->Write();
1236   else Printf("Cannot save fHistTimeZDC");
1237   if (fHistFMDA) fHistFMDA->Write();
1238   else Printf("Cannot save fHistFMDA");
1239   if (fHistFMDC) fHistFMDC->Write();
1240   else Printf("Cannot save fHistFMDC");
1241   if (fHistFMDSingle) fHistFMDSingle->Write();
1242   else Printf("Cannot save fHistFMDSingle");
1243   if (fHistFMDSum) fHistFMDSum->Write();
1244   else Printf("Cannot save fHistFMDSum");
1245   if (fSPDGFOEfficiency) fSPDGFOEfficiency->Write("fSPDGFOEfficiency");
1246   //  else Printf("Cannot save fSPDGFOEfficiency");
1247   
1248   fTriggerClasses->Write("fTriggerClasses", TObject::kSingleKey);
1249 }
1250
1251 void AliTriggerAnalysis::PrintTriggerClasses() const
1252 {
1253   // print trigger classes
1254   
1255   Printf("Trigger Classes:");
1256   
1257   Printf("Event count for trigger combinations:");
1258   
1259   TMap singleTrigger;
1260   singleTrigger.SetOwner();
1261   
1262   TIterator* iter = fTriggerClasses->MakeIterator();
1263   TObjString* obj = 0;
1264   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
1265   {
1266     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj));
1267     
1268     Printf(" %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1269     
1270     TObjArray* tokens = obj->String().Tokenize(" ");
1271     for (Int_t i=0; i<tokens->GetEntries(); i++)
1272     {
1273       TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(((TObjString*) tokens->At(i))->String().Data()));
1274       if (!count)
1275       {
1276         count = new TParameter<Long64_t>(((TObjString*) tokens->At(i))->String().Data(), 0);
1277         singleTrigger.Add(new TObjString(((TObjString*) tokens->At(i))->String().Data()), count);
1278       }
1279       count->SetVal(count->GetVal() + param->GetVal());
1280     }
1281     
1282     delete tokens;
1283   }
1284   delete iter;
1285   
1286   Printf("Event count for single trigger:");
1287   
1288   iter = singleTrigger.MakeIterator();
1289   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
1290   {
1291     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(obj));
1292     
1293     Printf("  %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1294   }
1295   delete iter;
1296   
1297   singleTrigger.DeleteAll();
1298 }