]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliTriggerAnalysis.cxx
fixes for uploading/compiling packages in plugin proof mode
[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   // TODO add first and last orbit number here
626 }
627
628 Int_t AliTriggerAnalysis::SSDClusters(const AliESDEvent* aEsd)
629 {
630   // returns the number of clusters in the SSD
631   const AliMultiplicity* mult = aEsd->GetMultiplicity();
632   Int_t clusters = mult->GetNumberOfITSClusters(4)+mult->GetNumberOfITSClusters(5);
633   return clusters;
634 }
635
636
637 Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHists, Int_t layer)
638 {
639   // returns the number of fired chips in the SPD
640   //
641   // origin = 0 --> aEsd->GetMultiplicity()->GetNumberOfFiredChips() (filled from clusters)
642   // origin = 1 --> aEsd->GetMultiplicity()->TestFastOrFiredChips() (from hardware bits)
643   // layer  = 0 --> both layers
644   // layer  = 1 --> inner
645   // layer  = 2 --> outer
646   
647   const AliMultiplicity* mult = aEsd->GetMultiplicity();
648   if (!mult)
649   {
650     AliError("AliMultiplicity not available");
651     return -1;
652   }
653   
654   if (origin == 0){
655     if (layer == 0) 
656       return mult->GetNumberOfFiredChips(0) + mult->GetNumberOfFiredChips(1);
657
658     return mult->GetNumberOfFiredChips(layer-1); 
659   }
660     
661   if (origin == 1)
662   {
663     Int_t nChips = 0;
664     Int_t firstChip = 0;
665     Int_t lastChip  = 1200;
666     if(layer == 1)
667       lastChip  = 400;
668     if(layer == 2)
669       firstChip = 400;
670
671     for (Int_t i=0; i<1200; i++)
672     {
673       if (mult->TestFastOrFiredChips(i) == kTRUE)
674       {
675         // efficiency simulation (if enabled)
676         if (fSPDGFOEfficiency)
677         {
678           if (gRandom->Uniform() > fSPDGFOEfficiency->GetBinContent(i+1))
679             continue;
680         }
681         
682         nChips++;
683         if (fillHists)
684           fHistFiredBitsSPD->Fill(i);
685       }
686     }
687     return nChips;
688   }
689   
690   return -1;
691 }
692
693 Bool_t AliTriggerAnalysis::SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin)
694 {
695   // Returns if the SPD gave a global Fast OR trigger
696   
697   Int_t firedChips = SPDFiredChips(aEsd, origin);
698   
699   if (firedChips >= fSPDGFOThreshold)
700     return kTRUE;
701   return kFALSE;
702 }
703
704 AliTriggerAnalysis::V0Decision AliTriggerAnalysis::V0Trigger(const AliESDEvent* aEsd, AliceSide side, Bool_t online, Bool_t fillHists)
705 {
706   // Returns the V0 trigger decision in V0A | V0C
707   //
708   // Returns kV0Fake if the calculated average time is in a window where neither BB nor BG is expected. 
709   // The rate of such triggers can be used to estimate the background. Note that the rate has to be 
710   // rescaled with the size of the windows (numerical values see below in the code)
711   //
712   // argument 'online' is used as a switch between online and offline trigger algorithms
713   //
714   // Based on an algorithm by Cvetan Cheshkov
715
716   AliESDVZERO* esdV0 = aEsd->GetVZEROData();
717   if (!esdV0)
718   {
719     AliError("AliESDVZERO not available");
720     return kV0Invalid;
721   }
722   AliDebug(2,Form("In V0Trigger: %f %f",esdV0->GetV0ATime(),esdV0->GetV0CTime()));
723
724   Int_t begin = -1;
725   Int_t end = -1;
726   
727   if (side == kASide)
728   {
729     begin = 32;
730     end = 64;
731   } 
732   else if (side == kCSide)
733   {
734     begin = 0;
735     end = 32;
736   }
737   else
738     return kV0Invalid;
739     
740    if (esdV0->TestBit(AliESDVZERO::kDecisionFilled)) {
741     if (online) {
742       if (esdV0->TestBit(AliESDVZERO::kOnlineBitsFilled)) {
743         for (Int_t i = begin; i < end; ++i) {
744           if (esdV0->GetBBFlag(i)) return kV0BB;
745         }
746         for (Int_t i = begin; i < end; ++i) {
747           if (esdV0->GetBGFlag(i)) return kV0BG;
748         }
749         return kV0Empty;
750       }
751       else {
752         AliWarning("V0 online trigger analysis is not yet available!");
753         return kV0BB;
754       }
755     }
756     else {
757
758       if (fillHists) {
759         if (side == kASide && fHistV0A)
760           fHistV0A->Fill(esdV0->GetV0ATime());
761         if (side == kCSide && fHistV0C)
762           fHistV0C->Fill(esdV0->GetV0CTime());
763       }
764
765       if (side == kASide) return (V0Decision)esdV0->GetV0ADecision();
766       else if (side == kCSide) return (V0Decision)esdV0->GetV0CDecision();
767       else return kV0Invalid;
768     }
769   }
770
771   Float_t time = 0;
772   Float_t weight = 0;
773   if (fMC)
774   {
775     Int_t runRange;
776     if (aEsd->GetRunNumber() <= 104803) runRange = 0;
777     else if (aEsd->GetRunNumber() <= 104876) runRange = 1;
778     else runRange = 2;
779
780     Float_t factors[3][64] = {
781       // runs: 104792-104803
782       {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},
783       // runs: 104841-104876
784       {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},
785       // runs: 104890-92
786       {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}
787     };
788     Float_t dA = 77.4 - 11.0;
789     Float_t dC = 77.4 - 2.9;
790     // Time misalignment
791     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};
792     Float_t dA2 = 2.8, dC2 = 3.3;
793
794     if (online) {
795       for (Int_t i = begin; i < end; ++i) {
796         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
797         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
798         if (esdV0->GetTime(i) >= 1e-6 &&
799             tempTime > fV0HwWinLow && tempTime < fV0HwWinHigh &&
800             tempAdc > fV0HwAdcThr)
801           return kV0BB;
802       }
803       return kV0Empty;
804     }
805     else {
806       for (Int_t i = begin; i < end; ++i) {
807         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
808         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA : esdV0->GetTime(i)+dC;
809         Float_t tempRawTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
810         if (esdV0->GetTime(i) >= 1e-6 &&
811             tempRawTime < 125.0 &&
812             tempAdc > fV0AdcThr) {
813           weight += 1.0;
814           time += tempTime;
815         }
816       }
817     }
818   }
819   else {
820     if (online) {
821       for (Int_t i = begin; i < end; ++i) {
822         if (esdV0->GetTime(i) >= 1e-6 &&
823             esdV0->GetTime(i) > fV0HwWinLow && esdV0->GetTime(i) < fV0HwWinHigh &&
824             esdV0->GetAdc(i) > fV0HwAdcThr)
825           return kV0BB;
826       }
827       return kV0Empty;
828     }
829     else {
830       for (Int_t i = begin; i < end; ++i) {
831         if (esdV0->GetTime(i) > 1e-6 && esdV0->GetAdc(i) > fV0AdcThr) {
832           Float_t correctedTime = V0CorrectLeadingTime(i, esdV0->GetTime(i), esdV0->GetAdc(i),aEsd->GetRunNumber());
833           Float_t timeWeight = V0LeadingTimeWeight(esdV0->GetAdc(i));
834           time += correctedTime*timeWeight;
835             
836           weight += timeWeight;
837         }
838       }
839     }
840   }
841
842   if (weight > 0) 
843     time /= weight;
844   time += fV0TimeOffset;
845
846   if (fillHists)
847   {
848     if (side == kASide && fHistV0A)
849       fHistV0A->Fill(time);
850     if (side == kCSide && fHistV0C)
851       fHistV0C->Fill(time);
852   }
853   
854   if (side == kASide)
855   {
856     if (time > 68 && time < 100)
857       return kV0BB;
858     if (time > 54 && time < 57.5) 
859       return kV0BG;
860     if (time > 57.5 && time < 68)
861       return kV0Fake;
862   }
863   
864   if (side == kCSide)
865   {
866     if (time > 75.5 && time < 100)
867       return kV0BB;
868     if (time > 69.5 && time < 73)
869       return kV0BG; 
870     if (time > 55 && time < 69.5)
871       return kV0Fake;
872   }
873   
874   return kV0Empty;
875 }
876
877 Float_t AliTriggerAnalysis::V0CorrectLeadingTime(Int_t i, Float_t time, Float_t adc, Int_t runNumber) const
878 {
879   // Correct for slewing and align the channels
880   //
881   // Authors: Cvetan Cheshkov / Raphael Tieulent
882
883   if (time == 0) return 0;
884
885   // Time alignment
886   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};
887
888   if(runNumber < 106031)
889     time -= timeShift[i];
890
891   // Slewing correction
892   if (adc == 0) return time;
893
894   Float_t p1 = 1.57345e1;
895   Float_t p2 =-4.25603e-1;
896
897   if(runNumber >= 106031) adc *= (2.5/4.0);
898   return (time - p1*TMath::Power(adc,p2));
899 }
900
901 Float_t AliTriggerAnalysis::V0LeadingTimeWeight(Float_t adc) const
902 {
903   if (adc < 1e-6) return 0;
904
905   Float_t p1 = 40.211;
906   Float_t p2 =-4.25603e-1;
907   Float_t p3 = 0.5646;
908
909   return 1./(p1*p1*TMath::Power(adc,2.*(p2-1.))+p3*p3);
910 }
911
912 Bool_t AliTriggerAnalysis::ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) const
913 {
914   // Returns if ZDC triggered
915   
916   AliESDZDC* zdcData = aEsd->GetESDZDC();
917   if (!zdcData)
918   {
919     AliError("AliESDZDC not available");
920     return kFALSE;
921   }
922   
923   UInt_t quality = zdcData->GetESDQuality();
924   
925   // from Nora's presentation, general first physics meeting 16.10.09
926   static UInt_t zpc  = 0x20;
927   static UInt_t znc  = 0x10;
928   static UInt_t zem1 = 0x08;
929   static UInt_t zem2 = 0x04;
930   static UInt_t zpa  = 0x02;
931   static UInt_t zna  = 0x01;
932   
933   if (side == kASide && ((quality & zpa) || (quality & zna)))
934     return kTRUE;
935   if (side == kCentralBarrel && ((quality & zem1) || (quality & zem2)))
936     return kTRUE;
937   if (side == kCSide && ((quality & zpc) || (quality & znc)))
938     return kTRUE;
939   
940   return kFALSE;
941 }
942
943 Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists)
944 {
945   // returns number of hit combinations agove threshold
946   //
947   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
948
949   if (!fDoFMD)
950     return -1;
951
952   // Workaround for AliESDEvent::GetFMDData is not const!
953   const AliESDFMD* fmdData = (const_cast<AliESDEvent*>(aEsd))->GetFMDData();
954   if (!fmdData)
955   {
956     AliError("AliESDFMD not available");
957     return -1;
958   }
959
960   Int_t detFrom = (side == kASide) ? 1 : 3;
961   Int_t detTo   = (side == kASide) ? 2 : 3;
962
963   Int_t triggers = 0;
964   Float_t totalMult = 0;
965   for (UShort_t det=detFrom;det<=detTo;det++) {
966     Int_t nRings = (det == 1 ? 1 : 2);
967     for (UShort_t ir = 0; ir < nRings; ir++) {
968       Char_t   ring = (ir == 0 ? 'I' : 'O');
969       UShort_t nsec = (ir == 0 ? 20  : 40);
970       UShort_t nstr = (ir == 0 ? 512 : 256);
971       for (UShort_t sec =0; sec < nsec;  sec++) {
972         for (UShort_t strip = 0; strip < nstr; strip++) {
973           Float_t mult = fmdData->Multiplicity(det,ring,sec,strip);
974           if (mult == AliESDFMD::kInvalidMult) continue;
975           
976           if (fillHists)
977             fHistFMDSingle->Fill(mult);
978           
979           if (mult > fFMDLowCut)
980             totalMult = totalMult + mult;
981           else
982           {
983             if (totalMult > fFMDHitCut)
984               triggers++;
985               
986             if (fillHists)
987               fHistFMDSum->Fill(totalMult);
988               
989             totalMult = 0;
990           }
991         }
992       }
993     }
994   }
995   
996   return triggers;
997 }
998
999 Bool_t AliTriggerAnalysis::FMDTrigger(const AliESDEvent* aEsd, AliceSide side)
1000 {
1001   // Returns if the FMD triggered
1002   //
1003   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
1004
1005   Int_t triggers = FMDHitCombinations(aEsd, side, kFALSE);
1006   
1007   if (triggers > 0)
1008     return kTRUE;
1009     
1010   return kFALSE;
1011 }
1012
1013 Long64_t AliTriggerAnalysis::Merge(TCollection* list)
1014 {
1015   // Merge a list of AliMultiplicityCorrection objects with this (needed for
1016   // PROOF).
1017   // Returns the number of merged objects (including this).
1018
1019   if (!list)
1020     return 0;
1021
1022   if (list->IsEmpty())
1023     return 1;
1024
1025   TIterator* iter = list->MakeIterator();
1026   TObject* obj;
1027
1028   // collections of all histograms
1029   const Int_t nHists = 9;
1030   TList collections[nHists];
1031
1032   Int_t count = 0;
1033   while ((obj = iter->Next())) {
1034
1035     AliTriggerAnalysis* entry = dynamic_cast<AliTriggerAnalysis*> (obj);
1036     if (entry == 0) 
1037       continue;
1038
1039     Int_t n = 0;
1040     collections[n++].Add(entry->fHistV0A);
1041     collections[n++].Add(entry->fHistV0C);
1042     collections[n++].Add(entry->fHistZDC);
1043     collections[n++].Add(entry->fHistFMDA);
1044     collections[n++].Add(entry->fHistFMDC);
1045     collections[n++].Add(entry->fHistFMDSingle);
1046     collections[n++].Add(entry->fHistFMDSum);
1047     collections[n++].Add(entry->fHistBitsSPD);
1048     collections[n++].Add(entry->fHistFiredBitsSPD);
1049
1050     // merge fTriggerClasses
1051     TIterator* iter2 = entry->fTriggerClasses->MakeIterator();
1052     TObjString* obj2 = 0;
1053     while ((obj2 = dynamic_cast<TObjString*> (iter2->Next())))
1054     {
1055       TParameter<Long64_t>* param2 = dynamic_cast<TParameter<Long64_t>*> (entry->fTriggerClasses->GetValue(obj2));
1056       
1057       TParameter<Long64_t>* param1 = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj2));
1058       if (param1)
1059       {
1060         param1->SetVal(param1->GetVal() + param2->GetVal());
1061       }
1062       else
1063       {
1064         param1 = dynamic_cast<TParameter<Long64_t>*> (param2->Clone());
1065         fTriggerClasses->Add(new TObjString(obj2->String()), param1);
1066       }
1067     }
1068     
1069     delete iter2;
1070   
1071     count++;
1072   }
1073
1074   Int_t n = 0;
1075   fHistV0A->Merge(&collections[n++]);
1076   fHistV0C->Merge(&collections[n++]);
1077   fHistZDC->Merge(&collections[n++]);
1078   fHistFMDA->Merge(&collections[n++]);
1079   fHistFMDC->Merge(&collections[n++]);
1080   fHistFMDSingle->Merge(&collections[n++]);
1081   fHistFMDSum->Merge(&collections[n++]);
1082   fHistBitsSPD->Merge(&collections[n++]);
1083   fHistFiredBitsSPD->Merge(&collections[n++]);
1084   
1085   delete iter;
1086
1087   return count+1;
1088 }
1089
1090 void AliTriggerAnalysis::SaveHistograms() const
1091 {
1092   // write histograms to current directory
1093   
1094   if (!fHistBitsSPD)
1095     return;
1096     
1097   fHistBitsSPD->Write();
1098   fHistBitsSPD->ProjectionX();
1099   fHistBitsSPD->ProjectionY();
1100   fHistFiredBitsSPD->Write();
1101   fHistV0A->Write();
1102   fHistV0C->Write();
1103   fHistZDC->Write();
1104   fHistFMDA->Write();
1105   fHistFMDC->Write();
1106   fHistFMDSingle->Write();
1107   fHistFMDSum->Write();
1108   
1109   if (fSPDGFOEfficiency)
1110     fSPDGFOEfficiency->Write("fSPDGFOEfficiency");
1111   
1112   fTriggerClasses->Write("fTriggerClasses", TObject::kSingleKey);
1113 }
1114
1115 void AliTriggerAnalysis::PrintTriggerClasses() const
1116 {
1117   // print trigger classes
1118   
1119   Printf("Trigger Classes:");
1120   
1121   Printf("Event count for trigger combinations:");
1122   
1123   TMap singleTrigger;
1124   singleTrigger.SetOwner();
1125   
1126   TIterator* iter = fTriggerClasses->MakeIterator();
1127   TObjString* obj = 0;
1128   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
1129   {
1130     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj));
1131     
1132     Printf(" %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1133     
1134     TObjArray* tokens = obj->String().Tokenize(" ");
1135     for (Int_t i=0; i<tokens->GetEntries(); i++)
1136     {
1137       TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(((TObjString*) tokens->At(i))->String().Data()));
1138       if (!count)
1139       {
1140         count = new TParameter<Long64_t>(((TObjString*) tokens->At(i))->String().Data(), 0);
1141         singleTrigger.Add(new TObjString(((TObjString*) tokens->At(i))->String().Data()), count);
1142       }
1143       count->SetVal(count->GetVal() + param->GetVal());
1144     }
1145     
1146     delete tokens;
1147   }
1148   delete iter;
1149   
1150   Printf("Event count for single trigger:");
1151   
1152   iter = singleTrigger.MakeIterator();
1153   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
1154   {
1155     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(obj));
1156     
1157     Printf("  %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1158   }
1159   delete iter;
1160   
1161   singleTrigger.DeleteAll();
1162 }