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