]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliTriggerAnalysis.cxx
Update the Tag system classes. Reduce memory footprint. Add information from RCT
[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 = aEsd->GetPrimaryVertexSPD();
464     Float_t ptmin, ptmax;
465     fEsdTrackCuts->GetPtRange(ptmin,ptmax);
466     AliDebug(3, Form("ptmin = %f, ptmax = %f\n",ptmin, ptmax));
467
468     if (vertex && vertex->GetNContributors() > 0 && (!vertex->IsFromVertexerZ() || vertex->GetDispersion() < 0.02) && TMath::Abs(vertex->GetZv()) < 10.) {
469       AliDebug(3,Form("Check on the vertex passed\n"));
470       for (Int_t i=0; i<aEsd->GetNumberOfTracks(); ++i){
471         if (fEsdTrackCuts->AcceptTrack(aEsd->GetTrack(i))){
472           AliDebug(2, Form("pt of track = %f --> check passed\n",aEsd->GetTrack(i)->Pt()));
473           decision = kTRUE;
474           break;
475         }
476       }
477     }
478     else{
479       AliDebug(4,Form("Check on the vertex not passed\n"));
480       for (Int_t i=0; i<aEsd->GetNumberOfTracks(); ++i){
481         if (fEsdTrackCuts->AcceptTrack(aEsd->GetTrack(i))){
482           AliDebug(4,Form("pt of track = %f --> check would be passed if the vertex was ok\n",aEsd->GetTrack(i)->Pt()));
483           break;
484         }
485       }
486     }
487     if (!decision) AliDebug(3,("Check for kOneTrack NOT passed\n"));
488   }
489
490   return decision;
491 }
492
493 Bool_t AliTriggerAnalysis::IsTriggerClassFired(const AliESDEvent* aEsd, const Char_t* tclass) const 
494 {
495   // tclass is logical function of inputs, e.g. 01 && 02 || 03 && 11 && 21
496   // = L0 inp 1 && L0 inp 2 || L0 inp 3 && L1 inp 1 && L2 inp 1
497   // NO brackets in logical function !
498   // Spaces between operators and inputs.
499   // Not all logical functions are available in CTP= 
500   // =any function of first 4 inputs; 'AND' of other inputs, check not done
501   // This method will be replaced/complemened by similar one
502   // which works withh class and inputs names as in CTP cfg file
503   
504   TString TClass(tclass);
505   TObjArray* tcltokens = TClass.Tokenize(" ");
506   Char_t level=((TObjString*)tcltokens->At(0))->String()[0];
507   UInt_t input=atoi((((TObjString*)tcltokens->At(0))->String()).Remove(0));
508   Bool_t tcl = IsInputFired(aEsd,level,input);
509  
510   for (Int_t i=1;i<tcltokens->GetEntriesFast();i=i+2) {
511     level=((TObjString*)tcltokens->At(i+1))->String()[0];
512     input=atoi((((TObjString*)tcltokens->At(i+1))->String()).Remove(0));
513     Bool_t inpnext = IsInputFired(aEsd,level,input);
514     Char_t op =((TObjString*)tcltokens->At(i))->String()[0];
515     if (op == '&') tcl=tcl && inpnext;
516     else if (op == '|') tcl =tcl || inpnext;
517     else {
518        AliError(Form("Syntax error in %s", tclass));
519        tcltokens->Delete();
520        return kFALSE;
521     }
522   }
523   tcltokens->Delete();
524   return tcl;
525 }
526
527 Bool_t AliTriggerAnalysis::IsInputFired(const AliESDEvent* aEsd, Char_t level, UInt_t input) const
528 {
529   // Checks trigger input of any level
530   
531   switch (level)
532   {
533     case '0': return IsL0InputFired(aEsd,input);
534     case '1': return IsL1InputFired(aEsd,input);
535     case '2': return IsL2InputFired(aEsd,input);
536     default:
537       AliError(Form("Wrong level %i",level));
538       return kFALSE;
539   }
540 }
541
542 Bool_t AliTriggerAnalysis::IsL0InputFired(const AliESDEvent* aEsd, UInt_t input) const 
543 {
544   // Checks if corresponding bit in mask is on
545   
546   UInt_t inpmask = aEsd->GetHeader()->GetL0TriggerInputs();
547   return (inpmask & (1<<(input-1)));
548 }
549
550 Bool_t AliTriggerAnalysis::IsL1InputFired(const AliESDEvent* aEsd, UInt_t input) const
551 {
552   // Checks if corresponding bit in mask is on
553   
554   UInt_t inpmask = aEsd->GetHeader()->GetL1TriggerInputs();
555   return (inpmask & (1<<(input-1)));
556 }
557
558 Bool_t AliTriggerAnalysis::IsL2InputFired(const AliESDEvent* aEsd, UInt_t input) const 
559 {
560   // Checks if corresponding bit in mask is on
561   
562   UInt_t inpmask = aEsd->GetHeader()->GetL2TriggerInputs();
563   return (inpmask & (1<<(input-1)));
564 }
565
566 void AliTriggerAnalysis::FillHistograms(const AliESDEvent* aEsd) 
567 {
568   // fills the histograms with the info from the ESD
569   
570   fHistBitsSPD->Fill(SPDFiredChips(aEsd, 0), SPDFiredChips(aEsd, 1, kTRUE));
571   
572   V0Trigger(aEsd, kASide, kFALSE, kTRUE);
573   V0Trigger(aEsd, kCSide, kFALSE, kTRUE);
574   
575   AliESDZDC* zdcData = aEsd->GetESDZDC();
576   if (zdcData)
577   {
578     UInt_t quality = zdcData->GetESDQuality();
579     
580     // from Nora's presentation, general first physics meeting 16.10.09
581     static UInt_t zpc  = 0x20;
582     static UInt_t znc  = 0x10;
583     static UInt_t zem1 = 0x08;
584     static UInt_t zem2 = 0x04;
585     static UInt_t zpa  = 0x02;
586     static UInt_t zna  = 0x01;
587    
588     fHistZDC->Fill(1, quality & zna);
589     fHistZDC->Fill(2, quality & zpa);
590     fHistZDC->Fill(3, quality & zem2);
591     fHistZDC->Fill(4, quality & zem1);
592     fHistZDC->Fill(5, quality & znc);
593     fHistZDC->Fill(6, quality & zpc);
594   }
595   else
596   {
597     fHistZDC->Fill(-1);
598     AliError("AliESDZDC not available");
599   }
600   
601   fHistFMDA->Fill(FMDHitCombinations(aEsd, kASide, kTRUE));
602   fHistFMDC->Fill(FMDHitCombinations(aEsd, kCSide, kTRUE));
603 }
604   
605 void AliTriggerAnalysis::FillTriggerClasses(const AliESDEvent* aEsd)
606 {
607   // fills trigger classes map
608   
609   TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(aEsd->GetFiredTriggerClasses().Data()));
610   if (!count)
611   {
612     count = new TParameter<Long64_t>(aEsd->GetFiredTriggerClasses(), 0);
613     fTriggerClasses->Add(new TObjString(aEsd->GetFiredTriggerClasses().Data()), count);
614   }
615   count->SetVal(count->GetVal() + 1);
616   
617   // TODO add first and last orbit number here
618 }
619
620 Int_t AliTriggerAnalysis::SSDClusters(const AliESDEvent* aEsd)
621 {
622   // returns the number of clusters in the SSD
623   const AliMultiplicity* mult = aEsd->GetMultiplicity();
624   Int_t clusters = mult->GetNumberOfITSClusters(4)+mult->GetNumberOfITSClusters(5);
625   return clusters;
626 }
627
628
629 Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHists)
630 {
631   // returns the number of fired chips in the SPD
632   //
633   // origin = 0 --> aEsd->GetMultiplicity()->GetNumberOfFiredChips() (filled from clusters)
634   // origin = 1 --> aEsd->GetMultiplicity()->TestFastOrFiredChips() (from hardware bits)
635   
636   const AliMultiplicity* mult = aEsd->GetMultiplicity();
637   if (!mult)
638   {
639     AliError("AliMultiplicity not available");
640     return -1;
641   }
642   
643   if (origin == 0)
644     return mult->GetNumberOfFiredChips(0) + mult->GetNumberOfFiredChips(1);
645     
646   if (origin == 1)
647   {
648     Int_t nChips = 0;
649     for (Int_t i=0; i<1200; i++)
650     {
651       if (mult->TestFastOrFiredChips(i) == kTRUE)
652       {
653         // efficiency simulation (if enabled)
654         if (fSPDGFOEfficiency)
655         {
656           if (gRandom->Uniform() > fSPDGFOEfficiency->GetBinContent(i+1))
657             continue;
658         }
659         
660         nChips++;
661         if (fillHists)
662           fHistFiredBitsSPD->Fill(i);
663       }
664     }
665     return nChips;
666   }
667   
668   return -1;
669 }
670
671 Bool_t AliTriggerAnalysis::SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin)
672 {
673   // Returns if the SPD gave a global Fast OR trigger
674   
675   Int_t firedChips = SPDFiredChips(aEsd, origin);
676   
677   if (firedChips >= fSPDGFOThreshold)
678     return kTRUE;
679   return kFALSE;
680 }
681
682 AliTriggerAnalysis::V0Decision AliTriggerAnalysis::V0Trigger(const AliESDEvent* aEsd, AliceSide side, Bool_t online, Bool_t fillHists)
683 {
684   // Returns the V0 trigger decision in V0A | V0C
685   //
686   // Returns kV0Fake if the calculated average time is in a window where neither BB nor BG is expected. 
687   // The rate of such triggers can be used to estimate the background. Note that the rate has to be 
688   // rescaled with the size of the windows (numerical values see below in the code)
689   //
690   // argument 'online' is used as a switch between online and offline trigger algorithms
691   //
692   // Based on an algorithm by Cvetan Cheshkov
693   
694   AliESDVZERO* esdV0 = aEsd->GetVZEROData();
695   if (!esdV0)
696   {
697     AliError("AliESDVZERO not available");
698     return kV0Invalid;
699   }
700
701   if (esdV0->TestBit(AliESDVZERO::kDecisionFilled)) {
702     if (online) {
703       AliWarning("V0 online trigger analysis is not yet available!");
704       return kV0BB;
705     }
706     else {
707
708       if (fillHists) {
709         if (side == kASide && fHistV0A)
710           fHistV0A->Fill(esdV0->GetV0ATime());
711         if (side == kCSide && fHistV0C)
712           fHistV0C->Fill(esdV0->GetV0CTime());
713       }
714
715       if (side == kASide) return (V0Decision)esdV0->GetV0ADecision();
716       else if (side == kCSide) return (V0Decision)esdV0->GetV0CDecision();
717       else return kV0Invalid;
718     }
719   }
720
721   Int_t begin = -1;
722   Int_t end = -1;
723   
724   if (side == kASide)
725   {
726     begin = 32;
727     end = 64;
728   } 
729   else if (side == kCSide)
730   {
731     begin = 0;
732     end = 32;
733   }
734   else
735     return kV0Invalid;
736     
737   Float_t time = 0;
738   Float_t weight = 0;
739   if (fMC)
740   {
741     Int_t runRange;
742     if (aEsd->GetRunNumber() <= 104803) runRange = 0;
743     else if (aEsd->GetRunNumber() <= 104876) runRange = 1;
744     else runRange = 2;
745
746     Float_t factors[3][64] = {
747       // runs: 104792-104803
748       {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},
749       // runs: 104841-104876
750       {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},
751       // runs: 104890-92
752       {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}
753     };
754     Float_t dA = 77.4 - 11.0;
755     Float_t dC = 77.4 - 2.9;
756     // Time misalignment
757     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};
758     Float_t dA2 = 2.8, dC2 = 3.3;
759
760     if (online) {
761       for (Int_t i = begin; i < end; ++i) {
762         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
763         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
764         if (esdV0->GetTime(i) >= 1e-6 &&
765             tempTime > fV0HwWinLow && tempTime < fV0HwWinHigh &&
766             tempAdc > fV0HwAdcThr)
767           return kV0BB;
768       }
769       return kV0Empty;
770     }
771     else {
772       for (Int_t i = begin; i < end; ++i) {
773         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
774         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA : esdV0->GetTime(i)+dC;
775         Float_t tempRawTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
776         if (esdV0->GetTime(i) >= 1e-6 &&
777             tempRawTime < 125.0 &&
778             tempAdc > fV0AdcThr) {
779           weight += 1.0;
780           time += tempTime;
781         }
782       }
783     }
784   }
785   else {
786     if (online) {
787       for (Int_t i = begin; i < end; ++i) {
788         if (esdV0->GetTime(i) >= 1e-6 &&
789             esdV0->GetTime(i) > fV0HwWinLow && esdV0->GetTime(i) < fV0HwWinHigh &&
790             esdV0->GetAdc(i) > fV0HwAdcThr)
791           return kV0BB;
792       }
793       return kV0Empty;
794     }
795     else {
796       for (Int_t i = begin; i < end; ++i) {
797         if (esdV0->GetTime(i) > 1e-6 && esdV0->GetAdc(i) > fV0AdcThr) {
798           Float_t correctedTime = V0CorrectLeadingTime(i, esdV0->GetTime(i), esdV0->GetAdc(i),aEsd->GetRunNumber());
799           Float_t timeWeight = V0LeadingTimeWeight(esdV0->GetAdc(i));
800           time += correctedTime*timeWeight;
801             
802           weight += timeWeight;
803         }
804       }
805     }
806   }
807
808   if (weight > 0) 
809     time /= weight;
810   time += fV0TimeOffset;
811
812   if (fillHists)
813   {
814     if (side == kASide && fHistV0A)
815       fHistV0A->Fill(time);
816     if (side == kCSide && fHistV0C)
817       fHistV0C->Fill(time);
818   }
819   
820   if (side == kASide)
821   {
822     if (time > 68 && time < 100)
823       return kV0BB;
824     if (time > 54 && time < 57.5) 
825       return kV0BG;
826     if (time > 57.5 && time < 68)
827       return kV0Fake;
828   }
829   
830   if (side == kCSide)
831   {
832     if (time > 75.5 && time < 100)
833       return kV0BB;
834     if (time > 69.5 && time < 73)
835       return kV0BG; 
836     if (time > 55 && time < 69.5)
837       return kV0Fake;
838   }
839   
840   return kV0Empty;
841 }
842
843 Float_t AliTriggerAnalysis::V0CorrectLeadingTime(Int_t i, Float_t time, Float_t adc, Int_t runNumber) const
844 {
845   // Correct for slewing and align the channels
846   //
847   // Authors: Cvetan Cheshkov / Raphael Tieulent
848
849   if (time == 0) return 0;
850
851   // Time alignment
852   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};
853
854   if(runNumber < 106031)
855     time -= timeShift[i];
856
857   // Slewing correction
858   if (adc == 0) return time;
859
860   Float_t p1 = 1.57345e1;
861   Float_t p2 =-4.25603e-1;
862
863   if(runNumber >= 106031) adc *= (2.5/4.0);
864   return (time - p1*TMath::Power(adc,p2));
865 }
866
867 Float_t AliTriggerAnalysis::V0LeadingTimeWeight(Float_t adc) const
868 {
869   if (adc < 1e-6) return 0;
870
871   Float_t p1 = 40.211;
872   Float_t p2 =-4.25603e-1;
873   Float_t p3 = 0.5646;
874
875   return 1./(p1*p1*TMath::Power(adc,2.*(p2-1.))+p3*p3);
876 }
877
878 Bool_t AliTriggerAnalysis::ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) const
879 {
880   // Returns if ZDC triggered
881   
882   AliESDZDC* zdcData = aEsd->GetESDZDC();
883   if (!zdcData)
884   {
885     AliError("AliESDZDC not available");
886     return kFALSE;
887   }
888   
889   UInt_t quality = zdcData->GetESDQuality();
890   
891   // from Nora's presentation, general first physics meeting 16.10.09
892   static UInt_t zpc  = 0x20;
893   static UInt_t znc  = 0x10;
894   static UInt_t zem1 = 0x08;
895   static UInt_t zem2 = 0x04;
896   static UInt_t zpa  = 0x02;
897   static UInt_t zna  = 0x01;
898   
899   if (side == kASide && ((quality & zpa) || (quality & zna)))
900     return kTRUE;
901   if (side == kCentralBarrel && ((quality & zem1) || (quality & zem2)))
902     return kTRUE;
903   if (side == kCSide && ((quality & zpc) || (quality & znc)))
904     return kTRUE;
905   
906   return kFALSE;
907 }
908
909 Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists)
910 {
911   // returns number of hit combinations agove threshold
912   //
913   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
914
915   // Workaround for AliESDEvent::GetFMDData is not const!
916   const AliESDFMD* fmdData = (const_cast<AliESDEvent*>(aEsd))->GetFMDData();
917   if (!fmdData)
918   {
919     AliError("AliESDFMD not available");
920     return -1;
921   }
922
923   Int_t detFrom = (side == kASide) ? 1 : 3;
924   Int_t detTo   = (side == kASide) ? 2 : 3;
925
926   Int_t triggers = 0;
927   Float_t totalMult = 0;
928   for (UShort_t det=detFrom;det<=detTo;det++) {
929     Int_t nRings = (det == 1 ? 1 : 2);
930     for (UShort_t ir = 0; ir < nRings; ir++) {
931       Char_t   ring = (ir == 0 ? 'I' : 'O');
932       UShort_t nsec = (ir == 0 ? 20  : 40);
933       UShort_t nstr = (ir == 0 ? 512 : 256);
934       for (UShort_t sec =0; sec < nsec;  sec++) {
935         for (UShort_t strip = 0; strip < nstr; strip++) {
936           Float_t mult = fmdData->Multiplicity(det,ring,sec,strip);
937           if (mult == AliESDFMD::kInvalidMult) continue;
938           
939           if (fillHists)
940             fHistFMDSingle->Fill(mult);
941           
942           if (mult > fFMDLowCut)
943             totalMult = totalMult + mult;
944           else
945           {
946             if (totalMult > fFMDHitCut)
947               triggers++;
948               
949             if (fillHists)
950               fHistFMDSum->Fill(totalMult);
951               
952             totalMult = 0;
953           }
954         }
955       }
956     }
957   }
958   
959   return triggers;
960 }
961
962 Bool_t AliTriggerAnalysis::FMDTrigger(const AliESDEvent* aEsd, AliceSide side)
963 {
964   // Returns if the FMD triggered
965   //
966   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
967
968   Int_t triggers = FMDHitCombinations(aEsd, side, kFALSE);
969   
970   if (triggers > 0)
971     return kTRUE;
972     
973   return kFALSE;
974 }
975
976 Long64_t AliTriggerAnalysis::Merge(TCollection* list)
977 {
978   // Merge a list of AliMultiplicityCorrection objects with this (needed for
979   // PROOF).
980   // Returns the number of merged objects (including this).
981
982   if (!list)
983     return 0;
984
985   if (list->IsEmpty())
986     return 1;
987
988   TIterator* iter = list->MakeIterator();
989   TObject* obj;
990
991   // collections of all histograms
992   const Int_t nHists = 9;
993   TList collections[nHists];
994
995   Int_t count = 0;
996   while ((obj = iter->Next())) {
997
998     AliTriggerAnalysis* entry = dynamic_cast<AliTriggerAnalysis*> (obj);
999     if (entry == 0) 
1000       continue;
1001
1002     Int_t n = 0;
1003     collections[n++].Add(entry->fHistV0A);
1004     collections[n++].Add(entry->fHistV0C);
1005     collections[n++].Add(entry->fHistZDC);
1006     collections[n++].Add(entry->fHistFMDA);
1007     collections[n++].Add(entry->fHistFMDC);
1008     collections[n++].Add(entry->fHistFMDSingle);
1009     collections[n++].Add(entry->fHistFMDSum);
1010     collections[n++].Add(entry->fHistBitsSPD);
1011     collections[n++].Add(entry->fHistFiredBitsSPD);
1012
1013     // merge fTriggerClasses
1014     TIterator* iter2 = entry->fTriggerClasses->MakeIterator();
1015     TObjString* obj2 = 0;
1016     while ((obj2 = dynamic_cast<TObjString*> (iter2->Next())))
1017     {
1018       TParameter<Long64_t>* param2 = dynamic_cast<TParameter<Long64_t>*> (entry->fTriggerClasses->GetValue(obj2));
1019       
1020       TParameter<Long64_t>* param1 = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj2));
1021       if (param1)
1022       {
1023         param1->SetVal(param1->GetVal() + param2->GetVal());
1024       }
1025       else
1026       {
1027         param1 = dynamic_cast<TParameter<Long64_t>*> (param2->Clone());
1028         fTriggerClasses->Add(new TObjString(obj2->String()), param1);
1029       }
1030     }
1031     
1032     delete iter2;
1033   
1034     count++;
1035   }
1036
1037   Int_t n = 0;
1038   fHistV0A->Merge(&collections[n++]);
1039   fHistV0C->Merge(&collections[n++]);
1040   fHistZDC->Merge(&collections[n++]);
1041   fHistFMDA->Merge(&collections[n++]);
1042   fHistFMDC->Merge(&collections[n++]);
1043   fHistFMDSingle->Merge(&collections[n++]);
1044   fHistFMDSum->Merge(&collections[n++]);
1045   fHistBitsSPD->Merge(&collections[n++]);
1046   fHistFiredBitsSPD->Merge(&collections[n++]);
1047   
1048   delete iter;
1049
1050   return count+1;
1051 }
1052
1053 void AliTriggerAnalysis::SaveHistograms() const
1054 {
1055   // write histograms to current directory
1056   
1057   if (!fHistBitsSPD)
1058     return;
1059     
1060   fHistBitsSPD->Write();
1061   fHistBitsSPD->ProjectionX();
1062   fHistBitsSPD->ProjectionY();
1063   fHistFiredBitsSPD->Write();
1064   fHistV0A->Write();
1065   fHistV0C->Write();
1066   fHistZDC->Write();
1067   fHistFMDA->Write();
1068   fHistFMDC->Write();
1069   fHistFMDSingle->Write();
1070   fHistFMDSum->Write();
1071   
1072   if (fSPDGFOEfficiency)
1073     fSPDGFOEfficiency->Write("fSPDGFOEfficiency");
1074   
1075   fTriggerClasses->Write("fTriggerClasses", TObject::kSingleKey);
1076 }
1077
1078 void AliTriggerAnalysis::PrintTriggerClasses() const
1079 {
1080   // print trigger classes
1081   
1082   Printf("Trigger Classes:");
1083   
1084   Printf("Event count for trigger combinations:");
1085   
1086   TMap singleTrigger;
1087   singleTrigger.SetOwner();
1088   
1089   TIterator* iter = fTriggerClasses->MakeIterator();
1090   TObjString* obj = 0;
1091   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
1092   {
1093     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj));
1094     
1095     Printf(" %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1096     
1097     TObjArray* tokens = obj->String().Tokenize(" ");
1098     for (Int_t i=0; i<tokens->GetEntries(); i++)
1099     {
1100       TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(((TObjString*) tokens->At(i))->String().Data()));
1101       if (!count)
1102       {
1103         count = new TParameter<Long64_t>(((TObjString*) tokens->At(i))->String().Data(), 0);
1104         singleTrigger.Add(new TObjString(((TObjString*) tokens->At(i))->String().Data()), count);
1105       }
1106       count->SetVal(count->GetVal() + param->GetVal());
1107     }
1108     
1109     delete tokens;
1110   }
1111   delete iter;
1112   
1113   Printf("Event count for single trigger:");
1114   
1115   iter = singleTrigger.MakeIterator();
1116   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
1117   {
1118     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(obj));
1119     
1120     Printf("  %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1121   }
1122   delete iter;
1123   
1124   singleTrigger.DeleteAll();
1125 }