]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliTriggerAnalysis.cxx
one track trigger using tpc only tracks (Chiara)
[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   fZDCCutRefSum(-568.5),
58   fZDCCutRefDelta(-2.1),
59   fZDCCutSigmaSum(3.25),
60   fZDCCutSigmaDelta(2.25),
61   fZDCCutRefSumCorr(-65.5),
62   fZDCCutRefDeltaCorr(-2.1),
63   fZDCCutSigmaSumCorr(6.0),
64   fZDCCutSigmaDeltaCorr(1.2),
65   fDoFMD(kTRUE),
66   fFMDLowCut(0.2),
67   fFMDHitCut(0.5),
68   fHistBitsSPD(0),
69   fHistFiredBitsSPD(0),
70   fHistV0A(0),       
71   fHistV0C(0),
72   fHistZDC(0),    
73   fHistTDCZDC(0),    
74   fHistTimeZDC(0),    
75   fHistTimeCorrZDC(0),    
76   fHistFMDA(0),    
77   fHistFMDC(0),   
78   fHistFMDSingle(0),
79   fHistFMDSum(0),
80   fTriggerClasses(0),
81   fMC(kFALSE),
82   fEsdTrackCuts(0),
83   fTPCOnly(kFALSE)
84 {
85   // constructor
86 }
87
88 AliTriggerAnalysis::~AliTriggerAnalysis()
89 {
90   // destructor
91   
92   if (fHistBitsSPD)
93   {
94     delete fHistBitsSPD;
95     fHistBitsSPD = 0;
96   }
97
98   if (fHistFiredBitsSPD)
99   {
100     delete fHistFiredBitsSPD;
101     fHistFiredBitsSPD = 0;
102   }
103
104   if (fHistV0A)
105   {
106     delete fHistV0A;
107     fHistV0A = 0;
108   }
109
110   if (fHistV0C)
111   {
112     delete fHistV0C;
113     fHistV0C = 0;
114   }
115
116   if (fHistZDC)
117   {
118     delete fHistZDC;
119     fHistZDC = 0;
120   }
121   if (fHistTDCZDC)
122   {
123     delete fHistTDCZDC;
124     fHistTDCZDC = 0;
125   }
126   if (fHistTimeZDC)
127   {
128     delete fHistTimeZDC;
129     fHistTimeZDC = 0;
130   }
131   if (fHistTimeCorrZDC)
132   {
133     delete fHistTimeCorrZDC;
134     fHistTimeCorrZDC = 0;
135   }
136
137   if (fHistFMDA)
138   {
139     delete fHistFMDA;
140     fHistFMDA = 0;
141   }
142
143   if (fHistFMDC)
144   {
145     delete fHistFMDC;
146     fHistFMDC = 0;
147   }
148
149   if (fHistFMDSingle)
150   {
151     delete fHistFMDSingle;
152     fHistFMDSingle = 0;
153   }
154
155   if (fHistFMDSum)
156   {
157     delete fHistFMDSum;
158     fHistFMDSum = 0;
159   }
160
161   if (fTriggerClasses)
162   {
163     fTriggerClasses->DeleteAll();
164     delete fTriggerClasses;
165     fTriggerClasses = 0;
166   }
167
168   if (fEsdTrackCuts){
169     delete fEsdTrackCuts;
170     fEsdTrackCuts =0;
171   }
172 }
173
174 void AliTriggerAnalysis::EnableHistograms()
175 {
176   // creates the monitoring histograms
177   
178   // do not add this hists to the directory
179   Bool_t oldStatus = TH1::AddDirectoryStatus();
180   TH1::AddDirectory(kFALSE);
181   
182   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);
183   fHistFiredBitsSPD = new TH1F("fHistFiredBitsSPD", "SPD GFO Hardware;chip number;events", 1200, -0.5, 1199.5);
184   fHistV0A = new TH1F("fHistV0A", "V0A;leading time (ns);events", 400, -100, 100);
185   fHistV0C = new TH1F("fHistV0C", "V0C;leading time (ns);events", 400, -100, 100);
186   fHistZDC = new TH1F("fHistZDC", "ZDC;trigger bits;events", 8, -1.5, 6.5);
187   fHistTDCZDC = new TH1F("fHistTDCZDC", "ZDC;TDC bits;events", 32, -0.5, 32-0.5);
188   fHistTimeZDC = new TH2F("fHistTimeZDC", "ZDC;TDC timing A+C vs C-A; events", 120,-30,30,120,-600,-540);
189   fHistTimeCorrZDC = new TH2F("fHistTimeCorrZDC", "ZDC;Corrected TDC timing A+C vs C-A; events", 120,-30,30,120,-100,-40);
190   
191   // TODO check limits
192   fHistFMDA = new TH1F("fHistFMDA", "FMDA;combinations above threshold;events", 102, -1.5, 100.5);
193   fHistFMDC = new TH1F("fHistFMDC", "FMDC;combinations above threshold;events", 102, -1.5, 100.5);
194   fHistFMDSingle = new TH1F("fHistFMDSingle", "FMD single;multiplicity value;counts", 1000, 0, 10);
195   fHistFMDSum = new TH1F("fHistFMDSum", "FMD sum;multiplicity value;counts", 1000, 0, 10);
196   
197   fTriggerClasses = new TMap;
198   fTriggerClasses->SetOwner();
199   
200   TH1::AddDirectory(oldStatus);
201 }
202
203 //____________________________________________________________________
204 const char* AliTriggerAnalysis::GetTriggerName(Trigger trigger) 
205 {
206   // returns the name of the requested trigger
207   // the returned string will only be valid until the next call to this function [not thread-safe]
208   
209   static TString str;
210   
211   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
212   
213   switch (triggerNoFlags)
214   {
215     case kAcceptAll : str = "ACCEPT ALL (bypass!)"; break;
216     case kMB1 : str = "MB1"; break;
217     case kMB2 : str = "MB2"; break;
218     case kMB3 : str = "MB3"; break;
219     case kSPDGFO : str = "SPD GFO"; break;
220     case kSPDGFOBits : str = "SPD GFO Bits"; break;
221     case kSPDGFOL0 : str = "SPD GFO L0 (first layer)"; break;
222     case kSPDGFOL1 : str = "SPD GFO L1 (second layer)"; break;
223     case kV0A : str = "V0 A BB"; break;
224     case kV0C : str = "V0 C BB"; break;
225     case kV0OR : str = "V0 OR BB"; break;
226     case kV0AND : str = "V0 AND BB"; break;
227     case kV0ABG : str = "V0 A BG"; break;
228     case kV0CBG : str = "V0 C BG"; break;
229     case kZDC : str = "ZDC"; break;
230     case kZDCA : str = "ZDC A"; break;
231     case kZDCC : str = "ZDC C"; break;
232     case kFMDA : str = "FMD A"; break;
233     case kFMDC : str = "FMD C"; break;
234     case kFPANY : str = "SPD GFO | V0 | ZDC | FMD"; break;
235     case kNSD1 : str = "NSD1"; break;
236     case kMB1Prime: str = "MB1prime"; break;
237     case kZDCTDCA : str = "ZDC TDC A"; break;
238     case kZDCTDCC : str = "ZDC TDC C"; break;
239     case kZDCTime : str = "ZDC Time Cut"; break;
240     default: str = ""; break;
241   }
242    
243   if (trigger & kOfflineFlag)
244     str += " OFFLINE";  
245   
246   if (trigger & kOneParticle)
247     str += " OneParticle";  
248
249   if (trigger & kOneTrack)
250     str += " OneTrack";  
251
252   return str;
253 }
254
255 Bool_t AliTriggerAnalysis::IsTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
256 {
257   // checks if an event has been triggered
258
259   if (trigger & kOfflineFlag)
260     return IsOfflineTriggerFired(aEsd, trigger);
261     
262   return IsTriggerBitFired(aEsd, trigger);
263 }
264
265 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, Trigger trigger) const
266 {
267   // checks if an event is fired using the trigger bits
268
269   return IsTriggerBitFired(aEsd->GetTriggerMask(), trigger);
270 }
271
272 Bool_t AliTriggerAnalysis::IsTriggerBitFired(ULong64_t triggerMask, Trigger trigger) const
273 {
274   // checks if an event is fired using the trigger bits
275   //
276   // this function needs the branch TriggerMask in the ESD
277   
278   // definitions from p-p.cfg
279   ULong64_t spdFO = (1 << 14);
280   ULong64_t v0left = (1 << 10);
281   ULong64_t v0right = (1 << 11);
282
283   switch (trigger)
284   {
285     case kAcceptAll:
286     {
287       return kTRUE;
288       break;
289     }
290     case kMB1:
291     {
292       if (triggerMask & spdFO || ((triggerMask & v0left) || (triggerMask & v0right)))
293         return kTRUE;
294       break;
295     }
296     case kMB2:
297     {
298       if (triggerMask & spdFO && ((triggerMask & v0left) || (triggerMask & v0right)))
299         return kTRUE;
300       break;
301     }
302     case kMB3:
303     {
304       if (triggerMask & spdFO && (triggerMask & v0left) && (triggerMask & v0right))
305         return kTRUE;
306       break;
307     }
308     case kSPDGFO:
309     {
310       if (triggerMask & spdFO)
311         return kTRUE;
312       break;
313     }
314     default:
315       Printf("IsEventTriggered: ERROR: Trigger type %d not implemented in this method", (Int_t) trigger);
316       break;
317   }
318
319   return kFALSE;
320 }
321
322 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, ULong64_t tclass) const
323 {
324   // Checks if corresponding bit in mask is on
325   
326   ULong64_t trigmask = aEsd->GetTriggerMask();
327   return (trigmask & (1ull << (tclass-1)));
328 }
329   
330 Int_t AliTriggerAnalysis::EvaluateTrigger(const AliESDEvent* aEsd, Trigger trigger)
331 {
332   // evaluates a given trigger "offline"
333   // trigger combinations are not supported, for that see IsOfflineTriggerFired
334   
335   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
336   Bool_t offline = kFALSE;
337   if (trigger & kOfflineFlag)
338     offline = kTRUE;
339   
340   Int_t decision = 0;
341   switch (triggerNoFlags)
342   {
343     case kSPDGFO:
344     {
345       decision = SPDFiredChips(aEsd, (offline) ? 0 : 1);
346       break;
347     }
348     case kSPDGFOL0:
349     {
350       decision = SPDFiredChips(aEsd, (offline) ? 0 : 1, kFALSE, 1);
351       break;
352     }
353     case kSPDGFOL1:
354     {
355       decision = SPDFiredChips(aEsd, (offline) ? 0 : 1, kFALSE, 2);
356       break;
357     }
358     case kV0A:
359     {
360       if (V0Trigger(aEsd, kASide, !offline) == kV0BB)
361         decision = 1;
362       break;
363     }
364     case kV0C:
365     {
366       if (V0Trigger(aEsd, kCSide, !offline) == kV0BB)
367         decision = 1;
368       break;
369     }
370     case kV0ABG:
371     {
372       if (V0Trigger(aEsd, kASide, !offline) == kV0BG)
373         decision = 1;
374       break;
375     }
376     case kV0CBG:
377     {
378       if (V0Trigger(aEsd, kCSide, !offline) == kV0BG)
379         decision = 1;
380       break;
381     }
382     case kZDCA:
383     {
384       if (!offline)
385         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
386       if (ZDCTrigger(aEsd, kASide))
387         decision = 1;
388       break;
389     }
390     case kZDCC:
391     {
392       if (!offline)
393         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
394       if (ZDCTrigger(aEsd, kCSide))
395         decision = 1;
396       break;
397     }
398     case kZDCTDCA:
399     {
400       if (!offline)
401         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
402       if (ZDCTDCTrigger(aEsd, kASide))
403         decision = 1;
404       break;
405     }
406     case kZDCTDCC:
407     {
408       if (!offline)
409         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
410       if (ZDCTDCTrigger(aEsd, kCSide))
411         decision = 1;
412       break;
413     }
414     case kZDCTime:
415     {
416       if (!offline)
417         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
418       if (ZDCTimeTrigger(aEsd))
419         decision = 1;
420       break;
421     }
422     case kFMDA:
423     {
424       if (!offline)
425         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
426       if (FMDTrigger(aEsd, kASide))
427         decision = 1;
428       break;
429     }
430     case kFMDC:
431     {
432       if (!offline)
433         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
434       if (FMDTrigger(aEsd, kCSide))
435         decision = 1;
436       break;
437     }
438     case kCTPV0A:
439     {
440       if (offline)
441         AliFatal(Form("Offline trigger not available for trigger %d", triggerNoFlags));
442       if (IsL0InputFired(aEsd, 2))
443         decision = 1;
444       break;
445     }
446     case kCTPV0C:
447     {
448       if (offline)
449         AliFatal(Form("Offline trigger not available for trigger %d", triggerNoFlags));
450       if (IsL0InputFired(aEsd, 3))
451         decision = 1;
452       break;
453     }
454     case kTPCLaserWarmUp:
455     {
456       if (!offline)
457         AliFatal(Form("Online trigger not available for trigger %d", triggerNoFlags));
458       return IsLaserWarmUpTPCEvent(aEsd);
459     }
460     default:
461     {
462       AliFatal(Form("Trigger type %d not implemented", triggerNoFlags));
463     }
464   }  
465   
466   return decision;
467 }
468
469 Bool_t AliTriggerAnalysis::IsLaserWarmUpTPCEvent(const AliESDEvent* esd)
470 {
471   //
472   // This function flags noisy TPC events which can happen during laser warm-up.
473   //
474   
475   Int_t trackCounter = 0;
476   for (Int_t i=0; i<esd->GetNumberOfTracks(); ++i) 
477   {
478     AliESDtrack *track = esd->GetTrack(i);
479     if (!track) 
480       continue;
481       
482     if (track->GetTPCNcls() < 30) continue;
483     if (TMath::Abs(track->Eta()) > 0.005) continue;
484     if (track->Pt() < 4) continue;
485     if (track->GetKinkIndex(0) > 0) continue;
486     
487     UInt_t status = track->GetStatus();
488     if ((status&AliESDtrack::kITSrefit)==AliESDtrack::kITSrefit) continue; // explicitly ask for tracks without ITS refit
489     if ((status&AliESDtrack::kTPCrefit)!=AliESDtrack::kTPCrefit) continue;
490     
491     if (track->GetTPCsignal() > 10) continue;          // explicitly ask for tracks without dE/dx
492     
493     trackCounter++;
494   }
495   if (trackCounter > 15) 
496     return kTRUE;
497   return kFALSE;
498 }
499
500 Bool_t AliTriggerAnalysis::IsOfflineTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
501 {
502   // checks if an event has been triggered "offline"
503
504   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
505   
506   Bool_t decision = kFALSE;
507   switch (triggerNoFlags)
508   {
509     case kAcceptAll:
510     {
511       decision = kTRUE;
512       break;
513     }
514     case kMB1:
515     {
516       if (SPDGFOTrigger(aEsd, 0) || V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
517         decision = kTRUE;
518       break;
519     }
520     case kMB2:
521     {
522       if (SPDGFOTrigger(aEsd, 0) && (V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB))
523         decision = kTRUE;
524       break;
525     }
526     case kMB3:
527     {
528       if (SPDGFOTrigger(aEsd, 0) && V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
529         decision = kTRUE;
530       break;
531     }
532     case kSPDGFO:
533     {
534       if (SPDGFOTrigger(aEsd, 0))
535         decision = kTRUE;
536       break;
537     }
538     case kSPDGFOBits:
539     {
540       if (SPDGFOTrigger(aEsd, 1))
541         decision = kTRUE;
542       break;
543     }
544     case kV0A:
545     {
546       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB)
547         decision = kTRUE;
548       break;
549     }
550     case kV0C:
551     {
552       if (V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
553         decision = kTRUE;
554       break;
555     }
556     case kV0OR:
557     {
558       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB || V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
559         decision = kTRUE;
560       break;
561     }
562     case kV0AND:
563     {
564       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
565         decision = kTRUE;
566       break;
567     }
568     case kV0ABG:
569     {
570       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BG)
571         decision = kTRUE;
572       break;
573     }
574     case kV0CBG:
575     {
576       if (V0Trigger(aEsd, kCSide, kFALSE) == kV0BG)
577         decision = kTRUE;
578       break;
579     }
580     case kZDC:
581     {
582       if (ZDCTrigger(aEsd, kASide) || ZDCTrigger(aEsd, kCentralBarrel) || ZDCTrigger(aEsd, kCSide))
583         decision = kTRUE;
584       break;
585     }
586     case kZDCA:
587     {
588       if (ZDCTrigger(aEsd, kASide))
589         decision = kTRUE;
590       break;
591     }
592     case kZDCC:
593     {
594       if (ZDCTrigger(aEsd, kCSide))
595         decision = kTRUE;
596       break;
597     }
598     case kFMDA:
599     {
600       if (FMDTrigger(aEsd, kASide))
601         decision = kTRUE;
602       break;
603     }
604     case kFMDC:
605     {
606       if (FMDTrigger(aEsd, kCSide))
607         decision = kTRUE;
608       break;
609     }
610     case kFPANY:
611     {
612       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))
613         decision = kTRUE;
614       break;
615     }
616     case kNSD1:
617     {
618       if (SPDFiredChips(aEsd, 0) >= 5 || (V0Trigger(aEsd, kASide, kFALSE) == kV0BB && V0Trigger(aEsd, kCSide, kFALSE) == kV0BB))
619         decision = kTRUE;
620        break;
621     }
622     case kMB1Prime:
623     {
624       Int_t count = 0;
625       if (SPDGFOTrigger(aEsd, 0))
626         count++;
627       if (V0Trigger(aEsd, kASide, kFALSE) == kV0BB)
628         count++;
629       if (V0Trigger(aEsd, kCSide, kFALSE) == kV0BB)
630         count++;
631       
632       if (count >= 2)
633         decision = kTRUE;
634         
635       break;
636     }
637     default:
638     {
639       AliFatal(Form("Trigger type %d not implemented", triggerNoFlags));
640     }
641   }
642   
643   // hadron-level requirement
644   if (decision && (trigger & kOneParticle))
645   {
646     decision = kFALSE;
647     
648     const AliESDVertex* vertex = aEsd->GetPrimaryVertexSPD();
649     const AliMultiplicity* mult = aEsd->GetMultiplicity();
650
651     if (mult && vertex && vertex->GetNContributors() > 0 && (!vertex->IsFromVertexerZ() || vertex->GetDispersion() < 0.02) && TMath::Abs(vertex->GetZv()) < 5.5) 
652     {
653       for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
654       {
655         if (TMath::Abs(mult->GetEta(i)) < 1)
656         {
657           decision = kTRUE;
658           break;
659         }
660       }
661     }
662   }
663
664   // hadron level definition for TPC tracks
665
666   if (decision && (trigger & kOneTrack))
667   {
668     decision = kFALSE;
669     const AliESDVertex* vertex =0x0;
670     vertex = aEsd->GetPrimaryVertexTracks();
671     if (!vertex || vertex->GetNContributors() <= 0)
672     {
673       vertex = aEsd->GetPrimaryVertexSPD();
674     }
675     Float_t ptmin, ptmax;
676     fEsdTrackCuts->GetPtRange(ptmin,ptmax);
677     AliDebug(3, Form("ptmin = %f, ptmax = %f\n",ptmin, ptmax));
678
679     if (vertex && vertex->GetNContributors() > 0 && (!vertex->IsFromVertexerZ() || vertex->GetDispersion() < 0.02) && TMath::Abs(vertex->GetZv()) < 10.) {
680       AliDebug(3,Form("Check on the vertex passed\n"));
681       for (Int_t i=0; i<aEsd->GetNumberOfTracks(); ++i){
682         if (fTPCOnly == kFALSE){
683           if (fEsdTrackCuts->AcceptTrack(aEsd->GetTrack(i))){
684             AliDebug(2, Form("pt of track = %f --> check passed\n",aEsd->GetTrack(i)->Pt()));
685             decision = kTRUE;
686              break;
687           }
688         }
689         else {
690           // TPC only tracks
691           AliESDtrack *tpcTrack = fEsdTrackCuts->GetTPCOnlyTrack((AliESDEvent*)aEsd, i);
692           if (!tpcTrack){
693             AliDebug(3,Form("track %d is NOT a TPC track",i));
694             continue;
695           }
696           else{
697             AliDebug(3,Form("track %d IS a TPC track",i));
698             if (!(fEsdTrackCuts->AcceptTrack(tpcTrack))) {
699               AliDebug(2, Form("TPC track %d NOT ACCEPTED, pt = %f, eta = %f",i,tpcTrack->Pt(), tpcTrack->Eta()));
700               delete tpcTrack; tpcTrack = 0x0;
701               continue;
702             } // end if the TPC track is not accepted
703             else{
704               AliDebug(2, Form("TPC track %d ACCEPTED, pt = %f, eta = %f",i,tpcTrack->Pt(), tpcTrack->Eta()));
705               decision = kTRUE;
706               delete tpcTrack; tpcTrack = 0x0;
707               break;
708             } // end if the TPC track is accepted
709           } // end if it is a TPC track
710         } // end if you are looking at TPC only tracks                  
711       } // end loop on tracks
712     } // end check on vertex
713     else{
714       AliDebug(4,Form("Check on the vertex not passed\n"));
715       for (Int_t i=0; i<aEsd->GetNumberOfTracks(); ++i){
716         if (fEsdTrackCuts->AcceptTrack(aEsd->GetTrack(i))){
717           AliDebug(4,Form("pt of track = %f --> check would be passed if the vertex was ok\n",aEsd->GetTrack(i)->Pt()));
718           break;
719         }
720       }
721     }
722     if (!decision) AliDebug(3,("Check for kOneTrack NOT passed\n"));
723   }
724
725   return decision;
726 }
727
728 Bool_t AliTriggerAnalysis::IsTriggerClassFired(const AliESDEvent* aEsd, const Char_t* tclass) const 
729 {
730   // tclass is logical function of inputs, e.g. 01 && 02 || 03 && 11 && 21
731   // = L0 inp 1 && L0 inp 2 || L0 inp 3 && L1 inp 1 && L2 inp 1
732   // NO brackets in logical function !
733   // Spaces between operators and inputs.
734   // Not all logical functions are available in CTP= 
735   // =any function of first 4 inputs; 'AND' of other inputs, check not done
736   // This method will be replaced/complemened by similar one
737   // which works withh class and inputs names as in CTP cfg file
738   
739   TString TClass(tclass);
740   TObjArray* tcltokens = TClass.Tokenize(" ");
741   Char_t level=((TObjString*)tcltokens->At(0))->String()[0];
742   UInt_t input=atoi((((TObjString*)tcltokens->At(0))->String()).Remove(0));
743   Bool_t tcl = IsInputFired(aEsd,level,input);
744  
745   for (Int_t i=1;i<tcltokens->GetEntriesFast();i=i+2) {
746     level=((TObjString*)tcltokens->At(i+1))->String()[0];
747     input=atoi((((TObjString*)tcltokens->At(i+1))->String()).Remove(0));
748     Bool_t inpnext = IsInputFired(aEsd,level,input);
749     Char_t op =((TObjString*)tcltokens->At(i))->String()[0];
750     if (op == '&') tcl=tcl && inpnext;
751     else if (op == '|') tcl =tcl || inpnext;
752     else {
753        AliError(Form("Syntax error in %s", tclass));
754        delete tcltokens;
755        tcltokens = 0;
756        //      tcltokens->Delete();
757        return kFALSE;
758     }
759   }
760   delete tcltokens;
761   tcltokens = 0;
762        //  tcltokens->Delete();
763   return tcl;
764 }
765
766 Bool_t AliTriggerAnalysis::IsInputFired(const AliESDEvent* aEsd, Char_t level, UInt_t input) const
767 {
768   // Checks trigger input of any level
769   
770   switch (level)
771   {
772     case '0': return IsL0InputFired(aEsd,input);
773     case '1': return IsL1InputFired(aEsd,input);
774     case '2': return IsL2InputFired(aEsd,input);
775     default:
776       AliError(Form("Wrong level %i",level));
777       return kFALSE;
778   }
779 }
780
781 Bool_t AliTriggerAnalysis::IsL0InputFired(const AliESDEvent* aEsd, UInt_t input) const 
782 {
783   // Checks if corresponding bit in mask is on
784   
785   UInt_t inpmask = aEsd->GetHeader()->GetL0TriggerInputs();
786   return (inpmask & (1<<(input-1)));
787 }
788
789 Bool_t AliTriggerAnalysis::IsL1InputFired(const AliESDEvent* aEsd, UInt_t input) const
790 {
791   // Checks if corresponding bit in mask is on
792   
793   UInt_t inpmask = aEsd->GetHeader()->GetL1TriggerInputs();
794   return (inpmask & (1<<(input-1)));
795 }
796
797 Bool_t AliTriggerAnalysis::IsL2InputFired(const AliESDEvent* aEsd, UInt_t input) const 
798 {
799   // Checks if corresponding bit in mask is on
800   
801   UInt_t inpmask = aEsd->GetHeader()->GetL2TriggerInputs();
802   return (inpmask & (1<<(input-1)));
803 }
804
805 void AliTriggerAnalysis::FillHistograms(const AliESDEvent* aEsd) 
806 {
807   // fills the histograms with the info from the ESD
808   
809   fHistBitsSPD->Fill(SPDFiredChips(aEsd, 0), SPDFiredChips(aEsd, 1, kTRUE));
810   
811   V0Trigger(aEsd, kASide, kFALSE, kTRUE);
812   V0Trigger(aEsd, kCSide, kFALSE, kTRUE);
813   ZDCTDCTrigger(aEsd,kASide,kFALSE,kFALSE,kTRUE);
814   ZDCTimeTrigger(aEsd,kTRUE);
815
816   AliESDZDC* zdcData = aEsd->GetESDZDC();
817   if (zdcData)
818   {
819     UInt_t quality = zdcData->GetESDQuality();
820     
821     // from Nora's presentation, general first physics meeting 16.10.09
822     static UInt_t zpc  = 0x20;
823     static UInt_t znc  = 0x10;
824     static UInt_t zem1 = 0x08;
825     static UInt_t zem2 = 0x04;
826     static UInt_t zpa  = 0x02;
827     static UInt_t zna  = 0x01;
828    
829     fHistZDC->Fill(1, quality & zna);
830     fHistZDC->Fill(2, quality & zpa);
831     fHistZDC->Fill(3, quality & zem2);
832     fHistZDC->Fill(4, quality & zem1);
833     fHistZDC->Fill(5, quality & znc);
834     fHistZDC->Fill(6, quality & zpc);
835   }
836   else
837   {
838     fHistZDC->Fill(-1);
839     AliError("AliESDZDC not available");
840   }
841   
842   if (fDoFMD) {
843     fHistFMDA->Fill(FMDHitCombinations(aEsd, kASide, kTRUE));
844     fHistFMDC->Fill(FMDHitCombinations(aEsd, kCSide, kTRUE));
845   }
846 }
847   
848 void AliTriggerAnalysis::FillTriggerClasses(const AliESDEvent* aEsd)
849 {
850   // fills trigger classes map
851   
852   TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(aEsd->GetFiredTriggerClasses().Data()));
853   if (!count)
854   {
855     count = new TParameter<Long64_t>(aEsd->GetFiredTriggerClasses(), 0);
856     fTriggerClasses->Add(new TObjString(aEsd->GetFiredTriggerClasses().Data()), count);
857   }
858   count->SetVal(count->GetVal() + 1);
859 }
860
861 Int_t AliTriggerAnalysis::SSDClusters(const AliESDEvent* aEsd)
862 {
863   // returns the number of clusters in the SSD
864   const AliMultiplicity* mult = aEsd->GetMultiplicity();
865   Int_t clusters = mult->GetNumberOfITSClusters(4)+mult->GetNumberOfITSClusters(5);
866   return clusters;
867 }
868
869
870 Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHists, Int_t layer)
871 {
872   // returns the number of fired chips in the SPD
873   //
874   // origin = 0 --> aEsd->GetMultiplicity()->GetNumberOfFiredChips() (filled from clusters)
875   // origin = 1 --> aEsd->GetMultiplicity()->TestFastOrFiredChips() (from hardware bits)
876   // layer  = 0 --> both layers
877   // layer  = 1 --> inner
878   // layer  = 2 --> outer
879   
880   const AliMultiplicity* mult = aEsd->GetMultiplicity();
881   if (!mult)
882   {
883     AliError("AliMultiplicity not available");
884     return -1;
885   }
886   
887   if (origin == 0){
888     if (layer == 0) 
889       return mult->GetNumberOfFiredChips(0) + mult->GetNumberOfFiredChips(1);
890
891     return mult->GetNumberOfFiredChips(layer-1); 
892   }
893     
894   if (origin == 1)
895   {
896     Int_t nChips = 0;
897     Int_t firstChip = 0;
898     Int_t lastChip  = 1200;
899     if(layer == 1)
900       lastChip  = 400;
901     if(layer == 2)
902       firstChip = 400;
903
904     for (Int_t i=firstChip; i<lastChip; i++)
905     {
906       if (mult->TestFastOrFiredChips(i) == kTRUE)
907       {
908         // efficiency simulation (if enabled)
909         if (fSPDGFOEfficiency)
910         {
911           if (gRandom->Uniform() > fSPDGFOEfficiency->GetBinContent(i+1))
912             continue;
913         }
914         
915         nChips++;
916         if (fillHists)
917           fHistFiredBitsSPD->Fill(i);
918       }
919     }
920     return nChips;
921   }
922   
923   return -1;
924 }
925
926 Bool_t AliTriggerAnalysis::SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin)
927 {
928   // Returns if the SPD gave a global Fast OR trigger
929   
930   Int_t firedChips = SPDFiredChips(aEsd, origin);
931   
932   if (firedChips >= fSPDGFOThreshold)
933     return kTRUE;
934   return kFALSE;
935 }
936
937 AliTriggerAnalysis::V0Decision AliTriggerAnalysis::V0Trigger(const AliESDEvent* aEsd, AliceSide side, Bool_t online, Bool_t fillHists)
938 {
939   // Returns the V0 trigger decision in V0A | V0C
940   //
941   // Returns kV0Fake if the calculated average time is in a window where neither BB nor BG is expected. 
942   // The rate of such triggers can be used to estimate the background. Note that the rate has to be 
943   // rescaled with the size of the windows (numerical values see below in the code)
944   //
945   // argument 'online' is used as a switch between online and offline trigger algorithms
946   //
947   // Based on an algorithm by Cvetan Cheshkov
948
949   AliESDVZERO* esdV0 = aEsd->GetVZEROData();
950   if (!esdV0)
951   {
952     AliError("AliESDVZERO not available");
953     return kV0Invalid;
954   }
955   AliDebug(2,Form("In V0Trigger: %f %f",esdV0->GetV0ATime(),esdV0->GetV0CTime()));
956
957   Int_t begin = -1;
958   Int_t end = -1;
959   
960   if (side == kASide)
961   {
962     begin = 32;
963     end = 64;
964   } 
965   else if (side == kCSide)
966   {
967     begin = 0;
968     end = 32;
969   }
970   else
971     return kV0Invalid;
972     
973    if (esdV0->TestBit(AliESDVZERO::kDecisionFilled)) {
974     if (online) {
975       if (esdV0->TestBit(AliESDVZERO::kOnlineBitsFilled)) {
976         for (Int_t i = begin; i < end; ++i) {
977           if (esdV0->GetBBFlag(i)) return kV0BB;
978         }
979         for (Int_t i = begin; i < end; ++i) {
980           if (esdV0->GetBGFlag(i)) return kV0BG;
981         }
982         return kV0Empty;
983       }
984       else {
985         AliWarning("V0 online trigger analysis is not yet available!");
986         return kV0BB;
987       }
988     }
989     else {
990
991       if (fillHists) {
992         if (side == kASide && fHistV0A)
993           fHistV0A->Fill(esdV0->GetV0ATime());
994         if (side == kCSide && fHistV0C)
995           fHistV0C->Fill(esdV0->GetV0CTime());
996       }
997
998       if (side == kASide) return (V0Decision)esdV0->GetV0ADecision();
999       else if (side == kCSide) return (V0Decision)esdV0->GetV0CDecision();
1000       else return kV0Invalid;
1001     }
1002   }
1003
1004   Float_t time = 0;
1005   Float_t weight = 0;
1006   if (fMC)
1007   {
1008     Int_t runRange;
1009     if (aEsd->GetRunNumber() <= 104803) runRange = 0;
1010     else if (aEsd->GetRunNumber() <= 104876) runRange = 1;
1011     else runRange = 2;
1012
1013     Float_t factors[3][64] = {
1014       // runs: 104792-104803
1015       {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},
1016       // runs: 104841-104876
1017       {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},
1018       // runs: 104890-92
1019       {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}
1020     };
1021     Float_t dA = 77.4 - 11.0;
1022     Float_t dC = 77.4 - 2.9;
1023     // Time misalignment
1024     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};
1025     Float_t dA2 = 2.8, dC2 = 3.3;
1026
1027     if (online) {
1028       for (Int_t i = begin; i < end; ++i) {
1029         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
1030         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
1031         if (esdV0->GetTime(i) >= 1e-6 &&
1032             tempTime > fV0HwWinLow && tempTime < fV0HwWinHigh &&
1033             tempAdc > fV0HwAdcThr)
1034           return kV0BB;
1035       }
1036       return kV0Empty;
1037     }
1038     else {
1039       for (Int_t i = begin; i < end; ++i) {
1040         Float_t tempAdc = esdV0->GetAdc(i)/factors[runRange][i];
1041         Float_t tempTime = (i >= 32) ? esdV0->GetTime(i)+dA : esdV0->GetTime(i)+dC;
1042         Float_t tempRawTime = (i >= 32) ? esdV0->GetTime(i)+dA+timeShift[i]+dA2 : esdV0->GetTime(i)+dC+timeShift[i]+dC2;
1043         if (esdV0->GetTime(i) >= 1e-6 &&
1044             tempRawTime < 125.0 &&
1045             tempAdc > fV0AdcThr) {
1046           weight += 1.0;
1047           time += tempTime;
1048         }
1049       }
1050     }
1051   }
1052   else {
1053     if (online) {
1054       for (Int_t i = begin; i < end; ++i) {
1055         if (esdV0->GetTime(i) >= 1e-6 &&
1056             esdV0->GetTime(i) > fV0HwWinLow && esdV0->GetTime(i) < fV0HwWinHigh &&
1057             esdV0->GetAdc(i) > fV0HwAdcThr)
1058           return kV0BB;
1059       }
1060       return kV0Empty;
1061     }
1062     else {
1063       for (Int_t i = begin; i < end; ++i) {
1064         if (esdV0->GetTime(i) > 1e-6 && esdV0->GetAdc(i) > fV0AdcThr) {
1065           Float_t correctedTime = V0CorrectLeadingTime(i, esdV0->GetTime(i), esdV0->GetAdc(i),aEsd->GetRunNumber());
1066           Float_t timeWeight = V0LeadingTimeWeight(esdV0->GetAdc(i));
1067           time += correctedTime*timeWeight;
1068             
1069           weight += timeWeight;
1070         }
1071       }
1072     }
1073   }
1074
1075   if (weight > 0) 
1076     time /= weight;
1077   time += fV0TimeOffset;
1078
1079   if (fillHists)
1080   {
1081     if (side == kASide && fHistV0A)
1082       fHistV0A->Fill(time);
1083     if (side == kCSide && fHistV0C)
1084       fHistV0C->Fill(time);
1085   }
1086   
1087   if (side == kASide)
1088   {
1089     if (time > 68 && time < 100)
1090       return kV0BB;
1091     if (time > 54 && time < 57.5) 
1092       return kV0BG;
1093     if (time > 57.5 && time < 68)
1094       return kV0Fake;
1095   }
1096   
1097   if (side == kCSide)
1098   {
1099     if (time > 75.5 && time < 100)
1100       return kV0BB;
1101     if (time > 69.5 && time < 73)
1102       return kV0BG; 
1103     if (time > 55 && time < 69.5)
1104       return kV0Fake;
1105   }
1106   
1107   return kV0Empty;
1108 }
1109
1110 Float_t AliTriggerAnalysis::V0CorrectLeadingTime(Int_t i, Float_t time, Float_t adc, Int_t runNumber) const
1111 {
1112   // Correct for slewing and align the channels
1113   //
1114   // Authors: Cvetan Cheshkov / Raphael Tieulent
1115
1116   if (time == 0) return 0;
1117
1118   // Time alignment
1119   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};
1120
1121   if(runNumber < 106031)
1122     time -= timeShift[i];
1123
1124   // Slewing correction
1125   if (adc == 0) return time;
1126
1127   Float_t p1 = 1.57345e1;
1128   Float_t p2 =-4.25603e-1;
1129
1130   if(runNumber >= 106031) adc *= (2.5/4.0);
1131   return (time - p1*TMath::Power(adc,p2));
1132 }
1133
1134 Float_t AliTriggerAnalysis::V0LeadingTimeWeight(Float_t adc) const
1135 {
1136   if (adc < 1e-6) return 0;
1137
1138   Float_t p1 = 40.211;
1139   Float_t p2 =-4.25603e-1;
1140   Float_t p3 = 0.5646;
1141
1142   return 1./(p1*p1*TMath::Power(adc,2.*(p2-1.))+p3*p3);
1143 }
1144
1145
1146 Bool_t AliTriggerAnalysis::ZDCTDCTrigger(const AliESDEvent* aEsd, AliceSide side, Bool_t useZN, Bool_t useZP, Bool_t fillHists) const
1147 {
1148   // Returns if ZDC triggered, based on TDC information 
1149   
1150   AliESDZDC *esdZDC = aEsd->GetESDZDC();
1151
1152   Bool_t zdcNA = kFALSE;
1153   Bool_t zdcNC = kFALSE;
1154   Bool_t zdcPA = kFALSE;
1155   Bool_t zdcPC = kFALSE;
1156
1157   if (fMC) {
1158         // If it's MC, we use the energy
1159     Double_t minEnergy = 0;
1160     Double_t  zNCEnergy = esdZDC->GetZDCN1Energy();
1161     Double_t  zPCEnergy = esdZDC->GetZDCP1Energy();
1162     Double_t  zNAEnergy = esdZDC->GetZDCN2Energy();
1163     Double_t  zPAEnergy = esdZDC->GetZDCP2Energy();
1164     zdcNA = (zNAEnergy>minEnergy);
1165     zdcNC = (zNCEnergy>minEnergy);
1166     zdcPA = (zPAEnergy>minEnergy);
1167     zdcPC = (zPCEnergy>minEnergy);
1168
1169   }
1170   else {
1171
1172     Bool_t tdc[32] = {kFALSE};
1173     for(Int_t itdc=0; itdc<32; itdc++){
1174       for(Int_t i=0; i<4; i++){
1175         if (esdZDC->GetZDCTDCData(itdc, i) != 0){
1176           tdc[itdc] = kTRUE;
1177         }
1178       }
1179       if(fillHists && tdc[itdc]) {
1180         fHistTDCZDC->Fill(itdc);
1181       }    
1182     }
1183     zdcNA = tdc[12];
1184     zdcNC = tdc[10];
1185     zdcPA = tdc[13];
1186     zdcPC = tdc[11];
1187   }
1188
1189   if (side == kASide) return ((useZP&&zdcPA) || (useZN&&zdcNA)); 
1190   if (side == kCSide) return ((useZP&&zdcPC) || (useZN&&zdcNC)); 
1191   return kFALSE;
1192 }
1193
1194 Bool_t AliTriggerAnalysis::ZDCTimeTrigger(const AliESDEvent *aEsd, Bool_t fillHists) const
1195 {
1196   // This method implements a selection
1197   // based on the timing in both sides of zdcN
1198   // It can be used in order to eliminate
1199   // parasitic collisions
1200   Bool_t zdcAccept = kFALSE;
1201   AliESDZDC *esdZDC = aEsd->GetESDZDC();
1202
1203   if(fMC) {
1204      UInt_t esdFlag =  esdZDC->GetESDQuality();
1205
1206      Bool_t znaFired=kFALSE, zpaFired=kFALSE;
1207      Bool_t zem1Fired=kFALSE, zem2Fired=kFALSE;
1208      Bool_t zncFired=kFALSE, zpcFired=kFALSE;
1209      //
1210      // **** Trigger patterns
1211      if((esdFlag & 0x00000001) == 0x00000001) znaFired=kTRUE;
1212      if((esdFlag & 0x00000002) == 0x00000002) zpaFired=kTRUE;
1213      if((esdFlag & 0x00000004) == 0x00000004) zem1Fired=kTRUE;
1214      if((esdFlag & 0x00000008) == 0x00000008) zem2Fired=kTRUE;
1215      if((esdFlag & 0x00000010) == 0x00000010) zncFired=kTRUE;
1216      if((esdFlag & 0x00000020) == 0x00000020) zpcFired=kTRUE;
1217      zdcAccept = (znaFired | zncFired);
1218   }
1219   else {
1220     for(Int_t i = 0; i < 4; ++i) {
1221       if (esdZDC->GetZDCTDCData(10,i) != 0) {
1222         Float_t tdcC = 0.025*(esdZDC->GetZDCTDCData(10,i)-esdZDC->GetZDCTDCData(14,i)); 
1223         Float_t tdcCcorr = esdZDC->GetZDCTDCCorrected(10,i); 
1224         for(Int_t j = 0; j < 4; ++j) {
1225           if (esdZDC->GetZDCTDCData(12,j) != 0) {
1226             Float_t tdcA = 0.025*(esdZDC->GetZDCTDCData(12,j)-esdZDC->GetZDCTDCData(14,j));
1227
1228             Float_t tdcAcorr = esdZDC->GetZDCTDCCorrected(12,j);
1229             if(fillHists) {
1230               fHistTimeZDC->Fill(tdcC-tdcA,tdcC+tdcA);
1231               fHistTimeCorrZDC->Fill(tdcCcorr-tdcAcorr,tdcCcorr+tdcAcorr);
1232             }
1233             if (esdZDC->TestBit(AliESDZDC::kCorrectedTDCFilled)) {
1234               if (((tdcCcorr-tdcAcorr-fZDCCutRefDeltaCorr)*(tdcCcorr-tdcAcorr-fZDCCutRefDeltaCorr)/(fZDCCutSigmaDeltaCorr*fZDCCutSigmaDeltaCorr) +
1235                    (tdcCcorr+tdcAcorr-fZDCCutRefSumCorr)*(tdcCcorr+tdcAcorr-fZDCCutRefSumCorr)/(fZDCCutSigmaSumCorr*fZDCCutSigmaSumCorr))< 1.0)
1236                 zdcAccept = kTRUE;
1237             }
1238             else {
1239               if (((tdcC-tdcA-fZDCCutRefDelta)*(tdcC-tdcA-fZDCCutRefDelta)/(fZDCCutSigmaDelta*fZDCCutSigmaDelta) +
1240                    (tdcC+tdcA-fZDCCutRefSum)*(tdcC+tdcA-fZDCCutRefSum)/(fZDCCutSigmaSum*fZDCCutSigmaSum))< 1.0)
1241                 zdcAccept = kTRUE;
1242             }
1243           }
1244         }
1245       }
1246     }
1247   }
1248   return zdcAccept;
1249 }
1250
1251 Bool_t AliTriggerAnalysis::ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) const
1252 {
1253   // Returns if ZDC triggered
1254   
1255   AliESDZDC* zdcData = aEsd->GetESDZDC();
1256   if (!zdcData)
1257   {
1258     AliError("AliESDZDC not available");
1259     return kFALSE;
1260   }
1261   
1262   UInt_t quality = zdcData->GetESDQuality();
1263   
1264   // from Nora's presentation, general first physics meeting 16.10.09
1265   static UInt_t zpc  = 0x20;
1266   static UInt_t znc  = 0x10;
1267   static UInt_t zem1 = 0x08;
1268   static UInt_t zem2 = 0x04;
1269   static UInt_t zpa  = 0x02;
1270   static UInt_t zna  = 0x01;
1271   
1272   if (side == kASide && ((quality & zpa) || (quality & zna)))
1273     return kTRUE;
1274   if (side == kCentralBarrel && ((quality & zem1) || (quality & zem2)))
1275     return kTRUE;
1276   if (side == kCSide && ((quality & zpc) || (quality & znc)))
1277     return kTRUE;
1278   
1279   return kFALSE;
1280 }
1281
1282 Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists)
1283 {
1284   // returns number of hit combinations agove threshold
1285   //
1286   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
1287
1288   if (!fDoFMD)
1289     return -1;
1290
1291   // Workaround for AliESDEvent::GetFMDData is not const!
1292   const AliESDFMD* fmdData = (const_cast<AliESDEvent*>(aEsd))->GetFMDData();
1293   if (!fmdData)
1294   {
1295     AliError("AliESDFMD not available");
1296     return -1;
1297   }
1298
1299   Int_t detFrom = (side == kASide) ? 1 : 3;
1300   Int_t detTo   = (side == kASide) ? 2 : 3;
1301
1302   Int_t triggers = 0;
1303   Float_t totalMult = 0;
1304   for (UShort_t det=detFrom;det<=detTo;det++) {
1305     Int_t nRings = (det == 1 ? 1 : 2);
1306     for (UShort_t ir = 0; ir < nRings; ir++) {
1307       Char_t   ring = (ir == 0 ? 'I' : 'O');
1308       UShort_t nsec = (ir == 0 ? 20  : 40);
1309       UShort_t nstr = (ir == 0 ? 512 : 256);
1310       for (UShort_t sec =0; sec < nsec;  sec++) {
1311         for (UShort_t strip = 0; strip < nstr; strip++) {
1312           Float_t mult = fmdData->Multiplicity(det,ring,sec,strip);
1313           if (mult == AliESDFMD::kInvalidMult) continue;
1314           
1315           if (fillHists)
1316             fHistFMDSingle->Fill(mult);
1317           
1318           if (mult > fFMDLowCut)
1319             totalMult = totalMult + mult;
1320           else
1321           {
1322             if (totalMult > fFMDHitCut)
1323               triggers++;
1324               
1325             if (fillHists)
1326               fHistFMDSum->Fill(totalMult);
1327               
1328             totalMult = 0;
1329           }
1330         }
1331       }
1332     }
1333   }
1334   
1335   return triggers;
1336 }
1337
1338 Bool_t AliTriggerAnalysis::FMDTrigger(const AliESDEvent* aEsd, AliceSide side)
1339 {
1340   // Returns if the FMD triggered
1341   //
1342   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
1343
1344   Int_t triggers = FMDHitCombinations(aEsd, side, kFALSE);
1345   
1346   if (triggers > 0)
1347     return kTRUE;
1348     
1349   return kFALSE;
1350 }
1351
1352 Long64_t AliTriggerAnalysis::Merge(TCollection* list)
1353 {
1354   // Merge a list of AliMultiplicityCorrection objects with this (needed for
1355   // PROOF).
1356   // Returns the number of merged objects (including this).
1357
1358   if (!list)
1359     return 0;
1360
1361   if (list->IsEmpty())
1362     return 1;
1363
1364   TIterator* iter = list->MakeIterator();
1365   TObject* obj;
1366
1367   // collections of all histograms
1368   const Int_t nHists = 12;
1369   TList collections[nHists];
1370
1371   Int_t count = 0;
1372   while ((obj = iter->Next())) {
1373
1374     AliTriggerAnalysis* entry = dynamic_cast<AliTriggerAnalysis*> (obj);
1375     if (entry == 0) 
1376       continue;
1377
1378     Int_t n = 0;
1379     collections[n++].Add(entry->fHistV0A);
1380     collections[n++].Add(entry->fHistV0C);
1381     collections[n++].Add(entry->fHistZDC);
1382     collections[n++].Add(entry->fHistTDCZDC);
1383     collections[n++].Add(entry->fHistTimeZDC);
1384     collections[n++].Add(entry->fHistTimeCorrZDC);
1385     collections[n++].Add(entry->fHistFMDA);
1386     collections[n++].Add(entry->fHistFMDC);
1387     collections[n++].Add(entry->fHistFMDSingle);
1388     collections[n++].Add(entry->fHistFMDSum);
1389     collections[n++].Add(entry->fHistBitsSPD);
1390     collections[n++].Add(entry->fHistFiredBitsSPD);
1391
1392     // merge fTriggerClasses
1393     TIterator* iter2 = entry->fTriggerClasses->MakeIterator();
1394     TObjString* obj2 = 0;
1395     while ((obj2 = dynamic_cast<TObjString*> (iter2->Next())))
1396     {
1397       TParameter<Long64_t>* param2 = static_cast<TParameter<Long64_t>*> (entry->fTriggerClasses->GetValue(obj2));
1398       
1399       TParameter<Long64_t>* param1 = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj2));
1400       if (param1)
1401       {
1402         param1->SetVal(param1->GetVal() + param2->GetVal());
1403       }
1404       else
1405       {
1406         param1 = dynamic_cast<TParameter<Long64_t>*> (param2->Clone());
1407         fTriggerClasses->Add(new TObjString(obj2->String()), param1);
1408       }
1409     }
1410     
1411     delete iter2;
1412   
1413     count++;
1414   }
1415
1416   Int_t n = 0;
1417   fHistV0A->Merge(&collections[n++]);
1418   fHistV0C->Merge(&collections[n++]);
1419   fHistZDC->Merge(&collections[n++]);
1420   fHistTDCZDC->Merge(&collections[n++]);
1421   if (fHistTimeZDC)
1422     fHistTimeZDC->Merge(&collections[n++]);
1423   else
1424     n++;
1425   if (fHistTimeCorrZDC)
1426     fHistTimeCorrZDC->Merge(&collections[n++]);
1427   else
1428     n++;
1429   fHistFMDA->Merge(&collections[n++]);
1430   fHistFMDC->Merge(&collections[n++]);
1431   fHistFMDSingle->Merge(&collections[n++]);
1432   fHistFMDSum->Merge(&collections[n++]);
1433   fHistBitsSPD->Merge(&collections[n++]);
1434   fHistFiredBitsSPD->Merge(&collections[n++]);
1435   
1436   delete iter;
1437
1438   return count+1;
1439 }
1440
1441 void AliTriggerAnalysis::SaveHistograms() const
1442 {
1443   // write histograms to current directory
1444   
1445   if (!fHistBitsSPD)
1446     return;
1447     
1448   if (fHistBitsSPD) {
1449     fHistBitsSPD->Write();
1450     fHistBitsSPD->ProjectionX();
1451     fHistBitsSPD->ProjectionY();
1452   }
1453   else Printf("Cannot save fHistBitsSPD");
1454   if (fHistFiredBitsSPD) fHistFiredBitsSPD->Write();
1455   else Printf("Cannot save fHistFiredBitsSPD");
1456   if (fHistV0A) fHistV0A->Write();
1457   else Printf("Cannot save fHistV0A");
1458   if (fHistV0C) fHistV0C->Write();
1459   else Printf("Cannot save fHistV0C");
1460   if (fHistZDC) fHistZDC->Write();
1461   else Printf("Cannot save fHistZDC");
1462   if (fHistTDCZDC) fHistTDCZDC->Write();
1463   else Printf("Cannot save fHistTDCZDC");
1464   if (fHistTimeZDC) fHistTimeZDC->Write();
1465   else Printf("Cannot save fHistTimeZDC");
1466   if (fHistTimeCorrZDC) fHistTimeCorrZDC->Write();
1467   else Printf("Cannot save fHistTimeCorrZDC");
1468   if (fHistFMDA) fHistFMDA->Write();
1469   else Printf("Cannot save fHistFMDA");
1470   if (fHistFMDC) fHistFMDC->Write();
1471   else Printf("Cannot save fHistFMDC");
1472   if (fHistFMDSingle) fHistFMDSingle->Write();
1473   else Printf("Cannot save fHistFMDSingle");
1474   if (fHistFMDSum) fHistFMDSum->Write();
1475   else Printf("Cannot save fHistFMDSum");
1476   if (fSPDGFOEfficiency) fSPDGFOEfficiency->Write("fSPDGFOEfficiency");
1477   //  else Printf("Cannot save fSPDGFOEfficiency");
1478   
1479   fTriggerClasses->Write("fTriggerClasses", TObject::kSingleKey);
1480 }
1481
1482 void AliTriggerAnalysis::PrintTriggerClasses() const
1483 {
1484   // print trigger classes
1485   
1486   Printf("Trigger Classes:");
1487   
1488   Printf("Event count for trigger combinations:");
1489   
1490   TMap singleTrigger;
1491   singleTrigger.SetOwner();
1492   
1493   TIterator* iter = fTriggerClasses->MakeIterator();
1494   TObjString* obj = 0;
1495   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
1496   {
1497     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj));
1498     
1499     Printf(" %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1500     
1501     TObjArray* tokens = obj->String().Tokenize(" ");
1502     for (Int_t i=0; i<tokens->GetEntries(); i++)
1503     {
1504       TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(((TObjString*) tokens->At(i))->String().Data()));
1505       if (!count)
1506       {
1507         count = new TParameter<Long64_t>(((TObjString*) tokens->At(i))->String().Data(), 0);
1508         singleTrigger.Add(new TObjString(((TObjString*) tokens->At(i))->String().Data()), count);
1509       }
1510       count->SetVal(count->GetVal() + param->GetVal());
1511     }
1512     
1513     delete tokens;
1514   }
1515   delete iter;
1516   
1517   Printf("Event count for single trigger:");
1518   
1519   iter = singleTrigger.MakeIterator();
1520   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
1521   {
1522     TParameter<Long64_t>* param = static_cast<TParameter<Long64_t>*> (singleTrigger.GetValue(obj));
1523     
1524     Printf("  %s: %ld triggers", obj->String().Data(), (Long_t)param->GetVal());
1525   }
1526   delete iter;
1527   
1528   singleTrigger.DeleteAll();
1529 }
1530
1531