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