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