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