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