]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliTriggerAnalysis.cxx
moving AliPhysicsSelection and AliTriggerAnalysis to ANALYSIS
[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
33 #include <AliTriggerAnalysis.h>
34
35 #include <AliLog.h>
36
37 #include <AliESDEvent.h>
38
39 #include <AliMultiplicity.h>
40 #include <AliESDVZERO.h>
41 #include <AliESDZDC.h>
42 #include <AliESDFMD.h>
43
44 ClassImp(AliTriggerAnalysis)
45
46 AliTriggerAnalysis::AliTriggerAnalysis() :
47   fSPDGFOThreshold(2),
48   fV0TimeOffset(0),
49   fFMDLowCut(0.2),
50   fFMDHitCut(0.5),
51   fHistBitsSPD(0),
52   fHistFiredBitsSPD(0),
53   fHistV0A(0),       
54   fHistV0C(0),
55   fHistZDC(0),    
56   fHistFMDA(0),    
57   fHistFMDC(0),   
58   fHistFMDSingle(0),
59   fHistFMDSum(0),
60   fTriggerClasses(0)
61 {
62   // constructor
63 }
64
65 AliTriggerAnalysis::~AliTriggerAnalysis()
66 {
67   // destructor
68   
69   if (fHistBitsSPD)
70   {
71     delete fHistBitsSPD;
72     fHistBitsSPD = 0;
73   }
74
75   if (fHistFiredBitsSPD)
76   {
77     delete fHistFiredBitsSPD;
78     fHistFiredBitsSPD = 0;
79   }
80
81   if (fHistV0A)
82   {
83     delete fHistV0A;
84     fHistV0A = 0;
85   }
86
87   if (fHistV0C)
88   {
89     delete fHistV0C;
90     fHistV0C = 0;
91   }
92
93   if (fHistZDC)
94   {
95     delete fHistZDC;
96     fHistZDC = 0;
97   }
98
99   if (fHistFMDA)
100   {
101     delete fHistFMDA;
102     fHistFMDA = 0;
103   }
104
105   if (fHistFMDC)
106   {
107     delete fHistFMDC;
108     fHistFMDC = 0;
109   }
110
111   if (fHistFMDSingle)
112   {
113     delete fHistFMDSingle;
114     fHistFMDSingle = 0;
115   }
116
117   if (fHistFMDSum)
118   {
119     delete fHistFMDSum;
120     fHistFMDSum = 0;
121   }
122
123   if (fTriggerClasses)
124   {
125     delete fTriggerClasses;
126     fTriggerClasses = 0;
127   }
128 }
129
130 void AliTriggerAnalysis::EnableHistograms()
131 {
132   // creates the monitoring histograms
133   
134   // do not add this hists to the directory
135   Bool_t oldStatus = TH1::AddDirectoryStatus();
136   TH1::AddDirectory(kFALSE);
137   
138   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);
139   fHistFiredBitsSPD = new TH1F("fHistFiredBitsSPD", "SPD GFO Hardware;chip number;events", 1200, -0.5, 1199.5);
140   fHistV0A = new TH1F("fHistV0A", "V0A;leading time (ns);events", 200, 0, 100);
141   fHistV0C = new TH1F("fHistV0C", "V0C;leading time (ns);events", 200, 0, 100);
142   fHistZDC = new TH1F("fHistZDC", "ZDC;trigger bits;events", 8, -1.5, 6.5);
143   
144   // TODO check limits
145   fHistFMDA = new TH1F("fHistFMDA", "FMDA;combinations above threshold;events", 102, -1.5, 100.5);
146   fHistFMDC = new TH1F("fHistFMDC", "FMDC;combinations above threshold;events", 102, -1.5, 100.5);
147   fHistFMDSingle = new TH1F("fHistFMDSingle", "FMD single;multiplicity value;counts", 1000, 0, 10);
148   fHistFMDSum = new TH1F("fHistFMDSum", "FMD sum;multiplicity value;counts", 1000, 0, 10);
149   
150   fTriggerClasses = new TMap;
151   fTriggerClasses->SetOwner();
152   
153   TH1::AddDirectory(oldStatus);
154 }
155
156 //____________________________________________________________________
157 const char* AliTriggerAnalysis::GetTriggerName(Trigger trigger) 
158 {
159   // returns the name of the requested trigger
160   // the returned string will only be valid until the next call to this function [not thread-safe]
161   
162   static TString str;
163   
164   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
165   
166   switch (triggerNoFlags)
167   {
168     case kAcceptAll : str = "ACCEPT ALL (bypass!)"; break;
169     case kMB1 : str = "MB1"; break;
170     case kMB2 : str = "MB2"; break;
171     case kMB3 : str = "MB3"; break;
172     case kSPDGFO : str = "SPD GFO"; break;
173     case kSPDGFOBits : str = "SPD GFO Bits"; break;
174     case kV0A : str = "V0 A BB"; break;
175     case kV0C : str = "V0 C BB"; break;
176     case kV0ABG : str = "V0 A BG"; break;
177     case kV0CBG : str = "V0 C BG"; break;
178     case kZDC : str = "ZDC"; break;
179     case kZDCA : str = "ZDC A"; break;
180     case kZDCC : str = "ZDC C"; break;
181     case kFMDA : str = "FMD A"; break;
182     case kFMDC : str = "FMD C"; break;
183     case kFPANY : str = "SPD GFO | V0 | ZDC | FMD"; break;
184     default: str = ""; break;
185   }
186    
187   if (trigger & kOfflineFlag)
188     str += " OFFLINE";  
189   
190   return str;
191 }
192
193 Bool_t AliTriggerAnalysis::IsTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
194 {
195   // checks if an event has been triggered
196
197   if (trigger & kOfflineFlag)
198     return IsOfflineTriggerFired(aEsd, trigger);
199     
200   return IsTriggerBitFired(aEsd, trigger);
201 }
202
203 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, Trigger trigger) const
204 {
205   // checks if an event is fired using the trigger bits
206
207   return IsTriggerBitFired(aEsd->GetTriggerMask(), trigger);
208 }
209
210 Bool_t AliTriggerAnalysis::IsTriggerBitFired(ULong64_t triggerMask, Trigger trigger) const
211 {
212   // checks if an event is fired using the trigger bits
213   //
214   // this function needs the branch TriggerMask in the ESD
215   
216   // definitions from p-p.cfg
217   ULong64_t spdFO = (1 << 14);
218   ULong64_t v0left = (1 << 10);
219   ULong64_t v0right = (1 << 11);
220
221   switch (trigger)
222   {
223     case kAcceptAll:
224     {
225       return kTRUE;
226       break;
227     }
228     case kMB1:
229     {
230       if (triggerMask & spdFO || ((triggerMask & v0left) || (triggerMask & v0right)))
231         return kTRUE;
232       break;
233     }
234     case kMB2:
235     {
236       if (triggerMask & spdFO && ((triggerMask & v0left) || (triggerMask & v0right)))
237         return kTRUE;
238       break;
239     }
240     case kMB3:
241     {
242       if (triggerMask & spdFO && (triggerMask & v0left) && (triggerMask & v0right))
243         return kTRUE;
244       break;
245     }
246     case kSPDGFO:
247     {
248       if (triggerMask & spdFO)
249         return kTRUE;
250       break;
251     }
252     default:
253       Printf("IsEventTriggered: ERROR: Trigger type %d not implemented in this method", (Int_t) trigger);
254       break;
255   }
256
257   return kFALSE;
258 }
259
260 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, ULong64_t tclass) const
261 {
262   // Checks if corresponding bit in mask is on
263   
264   ULong64_t trigmask = aEsd->GetTriggerMask();
265   return (trigmask & (1ull << (tclass-1)));
266 }
267
268 Bool_t AliTriggerAnalysis::IsOfflineTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
269 {
270   // checks if an event has been triggered "offline"
271
272   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
273   
274   switch (triggerNoFlags)
275   {
276     case kAcceptAll:
277     {
278       return kTRUE;
279       break;
280     }
281     case kMB1:
282     {
283       if (SPDGFOTrigger(aEsd, 0) || V0Trigger(aEsd, kASide) == kV0BB || V0Trigger(aEsd, kCSide) == kV0BB)
284         return kTRUE;
285       break;
286     }
287     case kMB2:
288     {
289       if (SPDGFOTrigger(aEsd, 0) && (V0Trigger(aEsd, kASide) == kV0BB || V0Trigger(aEsd, kCSide) == kV0BB))
290         return kTRUE;
291       break;
292     }
293     case kMB3:
294     {
295       if (SPDGFOTrigger(aEsd, 0) && V0Trigger(aEsd, kASide) == kV0BB && V0Trigger(aEsd, kCSide) == kV0BB)
296         return kTRUE;
297       break;
298     }
299     case kSPDGFO:
300     {
301       if (SPDGFOTrigger(aEsd, 0))
302         return kTRUE;
303       break;
304     }
305     case kSPDGFOBits:
306     {
307       if (SPDGFOTrigger(aEsd, 1))
308         return kTRUE;
309       break;
310     }
311     case kV0A:
312     {
313       if (V0Trigger(aEsd, kASide) == kV0BB)
314         return kTRUE;
315       break;
316     }
317     case kV0C:
318     {
319       if (V0Trigger(aEsd, kCSide) == kV0BB)
320         return kTRUE;
321       break;
322     }
323     case kV0ABG:
324     {
325       if (V0Trigger(aEsd, kASide) == kV0BG)
326         return kTRUE;
327       break;
328     }
329     case kV0CBG:
330     {
331       if (V0Trigger(aEsd, kCSide) == kV0BG)
332         return kTRUE;
333       break;
334     }
335     case kZDC:
336     {
337       if (ZDCTrigger(aEsd, kASide) || ZDCTrigger(aEsd, kCentralBarrel) || ZDCTrigger(aEsd, kCSide))
338         return kTRUE;
339       break;
340     }
341     case kZDCA:
342     {
343       if (ZDCTrigger(aEsd, kASide))
344         return kTRUE;
345       break;
346     }
347     case kZDCC:
348     {
349       if (ZDCTrigger(aEsd, kCSide))
350         return kTRUE;
351       break;
352     }
353     case kFMDA:
354     {
355       if (FMDTrigger(aEsd, kASide))
356         return kTRUE;
357       break;
358     }
359     case kFMDC:
360     {
361       if (FMDTrigger(aEsd, kCSide))
362         return kTRUE;
363       break;
364     }
365     case kFPANY:
366     {
367       if (SPDGFOTrigger(aEsd, 0) || V0Trigger(aEsd, kASide) == kV0BB || V0Trigger(aEsd, kCSide) == kV0BB || ZDCTrigger(aEsd, kASide) || ZDCTrigger(aEsd, kCentralBarrel) || ZDCTrigger(aEsd, kCSide) || FMDTrigger(aEsd, kASide) || FMDTrigger(aEsd, kCSide))
368         return kTRUE;
369       break;
370     }
371     default:
372     {
373       AliFatal(Form("Trigger type %d not implemented", triggerNoFlags));
374     }
375   }
376   
377   return kFALSE;
378 }
379
380
381 Bool_t AliTriggerAnalysis::IsTriggerClassFired(const AliESDEvent* aEsd, const Char_t* tclass) const 
382 {
383   // tclass is logical function of inputs, e.g. 01 && 02 || 03 && 11 && 21
384   // = L0 inp 1 && L0 inp 2 || L0 inp 3 && L1 inp 1 && L2 inp 1
385   // NO brackets in logical function !
386   // Spaces between operators and inputs.
387   // Not all logical functions are available in CTP= 
388   // =any function of first 4 inputs; 'AND' of other inputs, check not done
389   // This method will be replaced/complemened by similar one
390   // which works withh class and inputs names as in CTP cfg file
391   
392   TString TClass(tclass);
393   TObjArray* tcltokens = TClass.Tokenize(" ");
394   Char_t level=((TObjString*)tcltokens->At(0))->String()[0];
395   UInt_t input=atoi((((TObjString*)tcltokens->At(0))->String()).Remove(0));
396   Bool_t tcl = IsInputFired(aEsd,level,input);
397  
398   for (Int_t i=1;i<tcltokens->GetEntriesFast();i=i+2) {
399     level=((TObjString*)tcltokens->At(i+1))->String()[0];
400     input=atoi((((TObjString*)tcltokens->At(i+1))->String()).Remove(0));
401     Bool_t inpnext = IsInputFired(aEsd,level,input);
402     Char_t op =((TObjString*)tcltokens->At(i))->String()[0];
403     if (op == '&') tcl=tcl && inpnext;
404     else if (op == '|') tcl =tcl || inpnext;
405     else {
406        AliError(Form("Syntax error in %s", tclass));
407        tcltokens->Delete();
408        return kFALSE;
409     }
410   }
411   tcltokens->Delete();
412   return tcl;
413 }
414
415 Bool_t AliTriggerAnalysis::IsInputFired(const AliESDEvent* aEsd, Char_t level, UInt_t input) const
416 {
417   // Checks trigger input of any level
418   
419   switch (level)
420   {
421     case '0': return IsL0InputFired(aEsd,input);
422     case '1': return IsL1InputFired(aEsd,input);
423     case '2': return IsL2InputFired(aEsd,input);
424     default:
425       AliError(Form("Wrong level %i",level));
426       return kFALSE;
427   }
428 }
429
430 Bool_t AliTriggerAnalysis::IsL0InputFired(const AliESDEvent* aEsd, UInt_t input) const 
431 {
432   // Checks if corresponding bit in mask is on
433   
434   UInt_t inpmask = aEsd->GetHeader()->GetL0TriggerInputs();
435   return (inpmask & (1<<(input-1)));
436 }
437
438 Bool_t AliTriggerAnalysis::IsL1InputFired(const AliESDEvent* aEsd, UInt_t input) const
439 {
440   // Checks if corresponding bit in mask is on
441   
442   UInt_t inpmask = aEsd->GetHeader()->GetL1TriggerInputs();
443   return (inpmask & (1<<(input-1)));
444 }
445
446 Bool_t AliTriggerAnalysis::IsL2InputFired(const AliESDEvent* aEsd, UInt_t input) const 
447 {
448   // Checks if corresponding bit in mask is on
449   
450   UInt_t inpmask = aEsd->GetHeader()->GetL2TriggerInputs();
451   return (inpmask & (1<<(input-1)));
452 }
453
454 void AliTriggerAnalysis::FillHistograms(const AliESDEvent* aEsd) 
455 {
456   // fills the histograms with the info from the ESD
457   
458   fHistBitsSPD->Fill(SPDFiredChips(aEsd, 0), SPDFiredChips(aEsd, 1, kTRUE));
459   
460   V0Trigger(aEsd, kASide, kTRUE);
461   V0Trigger(aEsd, kCSide, kTRUE);
462   
463   AliESDZDC* zdcData = aEsd->GetESDZDC();
464   if (zdcData)
465   {
466     UInt_t quality = zdcData->GetESDQuality();
467     
468     // from Nora's presentation, general first physics meeting 16.10.09
469     static UInt_t zpc  = 0x20;
470     static UInt_t znc  = 0x10;
471     static UInt_t zem1 = 0x08;
472     static UInt_t zem2 = 0x04;
473     static UInt_t zpa  = 0x02;
474     static UInt_t zna  = 0x01;
475    
476     fHistZDC->Fill(1, quality & zna);
477     fHistZDC->Fill(2, quality & zpa);
478     fHistZDC->Fill(3, quality & zem2);
479     fHistZDC->Fill(4, quality & zem1);
480     fHistZDC->Fill(5, quality & znc);
481     fHistZDC->Fill(6, quality & zpc);
482   }
483   else
484   {
485     fHistZDC->Fill(-1);
486     AliError("AliESDZDC not available");
487   }
488   
489   fHistFMDA->Fill(FMDHitCombinations(aEsd, kASide, kTRUE));
490   fHistFMDC->Fill(FMDHitCombinations(aEsd, kCSide, kTRUE));
491 }
492   
493 void AliTriggerAnalysis::FillTriggerClasses(const AliESDEvent* aEsd)
494 {
495   // fills trigger classes map
496   
497   TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(aEsd->GetFiredTriggerClasses().Data()));
498   if (!count)
499   {
500     count = new TParameter<Long64_t>(aEsd->GetFiredTriggerClasses(), 0);
501     fTriggerClasses->Add(new TObjString(aEsd->GetFiredTriggerClasses().Data()), count);
502   }
503   count->SetVal(count->GetVal() + 1);
504   
505   // TODO add first and last orbit number here
506 }
507
508 Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHists)
509 {
510   // returns the number of fired chips in the SPD
511   //
512   // origin = 0 --> aEsd->GetMultiplicity()->GetNumberOfFiredChips() (filled from clusters)
513   // origin = 1 --> aEsd->GetMultiplicity()->TestFastOrFiredChips() (from hardware bits)
514   
515   const AliMultiplicity* mult = aEsd->GetMultiplicity();
516   if (!mult)
517   {
518     AliError("AliMultiplicity not available");
519     return -1;
520   }
521   
522   if (origin == 0)
523     return mult->GetNumberOfFiredChips(0) + mult->GetNumberOfFiredChips(1);
524     
525   if (origin == 1)
526   {
527     Int_t nChips = 0;
528     for (Int_t i=0; i<1200; i++)
529       if (mult->TestFastOrFiredChips(i) == kTRUE)
530       {
531         nChips++;
532         if (fillHists)
533           fHistFiredBitsSPD->Fill(i);
534       }
535     return nChips;
536   }
537   
538   return -1;
539 }
540
541 Bool_t AliTriggerAnalysis::SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin)
542 {
543   // Returns if the SPD gave a global Fast OR trigger
544   
545   Int_t firedChips = SPDFiredChips(aEsd, origin);
546   
547   if (firedChips >= fSPDGFOThreshold)
548     return kTRUE;
549   return kFALSE;
550 }
551
552 AliTriggerAnalysis::V0Decision AliTriggerAnalysis::V0Trigger(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists)
553 {
554   // Returns the V0 trigger decision in V0A | V0C
555   //
556   // Based on algorithm by Cvetan Cheshkov
557   
558   AliESDVZERO* esdV0 = aEsd->GetVZEROData();
559   if (!esdV0)
560   {
561     AliError("AliESDVZERO not available");
562     return kV0Invalid;
563   }
564   
565   Int_t begin = -1;
566   Int_t end = -1;
567   
568   if (side == kASide)
569   {
570     begin = 32;
571     end = 64;
572   } 
573   else if (side == kCSide)
574   {
575     begin = 0;
576     end = 32;
577   }
578   else
579     return kV0Invalid;
580   
581   Float_t time = 0;
582   Float_t weight = 0;
583   for (Int_t i = begin; i < end; ++i) {
584     if (esdV0->GetTime(i) > 1e-6) {
585       Float_t correctedTime = V0CorrectLeadingTime(i, esdV0->GetTime(i), esdV0->GetAdc(i));
586       Float_t timeWeight = V0LeadingTimeWeight(esdV0->GetAdc(i));
587       time += correctedTime*timeWeight;
588       
589       weight += timeWeight;
590     }
591   }
592
593   if (weight > 0) 
594     time /= weight;
595   time += fV0TimeOffset;
596
597   if (fillHists)
598   {
599     if (side == kASide && fHistV0A)
600       fHistV0A->Fill(time);
601     if (side == kCSide && fHistV0C)
602       fHistV0C->Fill(time);
603   }
604   
605   if (side == kASide)
606   {
607     if (time > 75 && time < 83)
608       return kV0BB;
609     if (time > 54 && time < 57.5) 
610       return kV0BG;
611   }
612   
613   if (side == kCSide)
614   {
615     if (time > 75.5 && time < 82)
616       return kV0BB;
617     if (time > 69.5 && time < 73)
618       return kV0BG; 
619   }
620   
621   return kV0Empty;
622 }
623
624 Float_t AliTriggerAnalysis::V0CorrectLeadingTime(Int_t i, Float_t time, Float_t adc) const
625 {
626   // Correct for slewing and align the channels
627   //
628   // Authors: Cvetan Cheshkov / Raphael Tieulent
629
630   if (time == 0) return 0;
631
632   // Time alignment
633   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};
634
635   time -= timeShift[i];
636
637   // Slewing correction
638   if (adc == 0) return time;
639
640   Float_t p1 = 1.57345e1;
641   Float_t p2 =-4.25603e-1;
642
643   return (time - p1*TMath::Power(adc,p2));
644 }
645
646 Float_t AliTriggerAnalysis::V0LeadingTimeWeight(Float_t adc) const
647 {
648   if (adc < 1e-6) return 0;
649
650   Float_t p1 = 40.211;
651   Float_t p2 =-4.25603e-1;
652   Float_t p3 = 0.5646;
653
654   return 1./(p1*p1*TMath::Power(adc,2.*(p2-1.))+p3*p3);
655 }
656
657 Bool_t AliTriggerAnalysis::ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) const
658 {
659   // Returns if ZDC triggered
660   
661   AliESDZDC* zdcData = aEsd->GetESDZDC();
662   if (!zdcData)
663   {
664     AliError("AliESDZDC not available");
665     return kFALSE;
666   }
667   
668   UInt_t quality = zdcData->GetESDQuality();
669   
670   // from Nora's presentation, general first physics meeting 16.10.09
671   static UInt_t zpc  = 0x20;
672   static UInt_t znc  = 0x10;
673   static UInt_t zem1 = 0x08;
674   static UInt_t zem2 = 0x04;
675   static UInt_t zpa  = 0x02;
676   static UInt_t zna  = 0x01;
677   
678   if (side == kASide && ((quality & zpa) || (quality & zna)))
679     return kTRUE;
680   if (side == kCentralBarrel && ((quality & zem1) || (quality & zem2)))
681     return kTRUE;
682   if (side == kCSide && ((quality & zpc) || (quality & znc)))
683     return kTRUE;
684   
685   return kFALSE;
686 }
687
688 Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists)
689 {
690   // returns number of hit combinations agove threshold
691   //
692   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
693
694   // Workaround for AliESDEvent::GetFMDData is not const!
695   const AliESDFMD* fmdData = (const_cast<AliESDEvent*>(aEsd))->GetFMDData();
696   if (!fmdData)
697   {
698     AliError("AliESDFMD not available");
699     return -1;
700   }
701
702   Int_t detFrom = (side == kASide) ? 1 : 3;
703   Int_t detTo   = (side == kASide) ? 2 : 3;
704
705   Int_t triggers = 0;
706   Float_t totalMult = 0;
707   for (UShort_t det=detFrom;det<=detTo;det++) {
708     Int_t nRings = (det == 1 ? 1 : 2);
709     for (UShort_t ir = 0; ir < nRings; ir++) {
710       Char_t   ring = (ir == 0 ? 'I' : 'O');
711       UShort_t nsec = (ir == 0 ? 20  : 40);
712       UShort_t nstr = (ir == 0 ? 512 : 256);
713       for (UShort_t sec =0; sec < nsec;  sec++) {
714         for (UShort_t strip = 0; strip < nstr; strip++) {
715           Float_t mult = fmdData->Multiplicity(det,ring,sec,strip);
716           if (mult == AliESDFMD::kInvalidMult) continue;
717           
718           if (fillHists)
719             fHistFMDSingle->Fill(mult);
720           
721           if (mult > fFMDLowCut)
722             totalMult = totalMult + mult;
723           else
724           {
725             if (totalMult > fFMDHitCut)
726               triggers++;
727               
728             if (fillHists)
729               fHistFMDSum->Fill(totalMult);
730               
731             totalMult = 0;
732           }
733         }
734       }
735     }
736   }
737   
738   return triggers;
739 }
740
741 Bool_t AliTriggerAnalysis::FMDTrigger(const AliESDEvent* aEsd, AliceSide side)
742 {
743   // Returns if the FMD triggered
744   //
745   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
746
747   Int_t triggers = FMDHitCombinations(aEsd, side, kFALSE);
748   
749   if (triggers > 0)
750     return kTRUE;
751     
752   return kFALSE;
753 }
754
755 Long64_t AliTriggerAnalysis::Merge(TCollection* list)
756 {
757   // Merge a list of AliMultiplicityCorrection objects with this (needed for
758   // PROOF).
759   // Returns the number of merged objects (including this).
760
761   if (!list)
762     return 0;
763
764   if (list->IsEmpty())
765     return 1;
766
767   TIterator* iter = list->MakeIterator();
768   TObject* obj;
769
770   // collections of all histograms
771   const Int_t nHists = 9;
772   TList collections[nHists];
773
774   Int_t count = 0;
775   while ((obj = iter->Next())) {
776
777     AliTriggerAnalysis* entry = dynamic_cast<AliTriggerAnalysis*> (obj);
778     if (entry == 0) 
779       continue;
780
781     Int_t n = 0;
782     collections[n++].Add(entry->fHistV0A);
783     collections[n++].Add(entry->fHistV0C);
784     collections[n++].Add(entry->fHistZDC);
785     collections[n++].Add(entry->fHistFMDA);
786     collections[n++].Add(entry->fHistFMDC);
787     collections[n++].Add(entry->fHistFMDSingle);
788     collections[n++].Add(entry->fHistFMDSum);
789     collections[n++].Add(entry->fHistBitsSPD);
790     collections[n++].Add(entry->fHistFiredBitsSPD);
791
792     // merge fTriggerClasses
793     TIterator* iter2 = entry->fTriggerClasses->MakeIterator();
794     TObjString* obj = 0;
795     while ((obj = dynamic_cast<TObjString*> (iter2->Next())))
796     {
797       TParameter<Long64_t>* param2 = dynamic_cast<TParameter<Long64_t>*> (entry->fTriggerClasses->GetValue(obj));
798       
799       TParameter<Long64_t>* param1 = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj));
800       if (param1)
801       {
802         param1->SetVal(param1->GetVal() + param2->GetVal());
803       }
804       else
805       {
806         param1 = dynamic_cast<TParameter<Long64_t>*> (param2->Clone());
807         fTriggerClasses->Add(new TObjString(obj->String()), param1);
808       }
809     }
810     
811     delete iter2;
812   
813     count++;
814   }
815
816   Int_t n = 0;
817   fHistV0A->Merge(&collections[n++]);
818   fHistV0C->Merge(&collections[n++]);
819   fHistZDC->Merge(&collections[n++]);
820   fHistFMDA->Merge(&collections[n++]);
821   fHistFMDC->Merge(&collections[n++]);
822   fHistFMDSingle->Merge(&collections[n++]);
823   fHistFMDSum->Merge(&collections[n++]);
824   fHistBitsSPD->Merge(&collections[n++]);
825   fHistFiredBitsSPD->Merge(&collections[n++]);
826   
827   delete iter;
828
829   return count+1;
830 }
831
832 void AliTriggerAnalysis::SaveHistograms() const
833 {
834   // write histograms to current directory
835   
836   if (!fHistBitsSPD)
837     return;
838     
839   fHistBitsSPD->Write();
840   fHistBitsSPD->ProjectionX();
841   fHistBitsSPD->ProjectionY();
842   fHistFiredBitsSPD->Write();
843   fHistV0A->Write();
844   fHistV0C->Write();
845   fHistZDC->Write();
846   fHistFMDA->Write();
847   fHistFMDC->Write();
848   fHistFMDSingle->Write();
849   fHistFMDSum->Write();
850   
851   fTriggerClasses->Write("fTriggerClasses", TObject::kSingleKey);
852 }
853
854 void AliTriggerAnalysis::PrintTriggerClasses() const
855 {
856   // print trigger classes
857   
858   Printf("Trigger Classes:");
859   
860   TIterator* iter = fTriggerClasses->MakeIterator();
861   TObjString* obj = 0;
862   while ((obj = dynamic_cast<TObjString*> (iter->Next())))
863   {
864     TParameter<Long64_t>* param = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj));
865     
866     Printf("%s: %ld triggers", obj->String().Data(), param->GetVal());
867   }
868 }