]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/AliTriggerAnalysis.cxx
first paper analysis code:
[u/mrichter/AliRoot.git] / PWG0 / 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
31 #include <AliTriggerAnalysis.h>
32
33 #include <AliLog.h>
34
35 #include <AliESDEvent.h>
36
37 #include <AliMultiplicity.h>
38 #include <AliESDVZERO.h>
39 #include <AliESDZDC.h>
40 #include <AliESDFMD.h>
41
42 ClassImp(AliTriggerAnalysis)
43
44 AliTriggerAnalysis::AliTriggerAnalysis() :
45   fSPDGFOThreshold(2),
46   fV0AThreshold(1),
47   fV0CThreshold(1),
48   fFMDLowCut(0.2),
49   fFMDHitCut(0.5),
50   fHistSPD(0),
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 {
61 }
62
63 void AliTriggerAnalysis::EnableHistograms()
64 {
65   // creates the monitoring histograms
66   
67   fHistSPD = new TH1F("fHistSPD", "SPD GFO;number of fired chips;events", 1202, -1.5, 1200.5);
68   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);
69   fHistFiredBitsSPD = new TH1F("fHistFiredBitsSPD", "SPD GFO HARDWARE;chip number;events", 1200, -0.5, 1199.5);
70   fHistV0A = new TH1F("fHistV0A", "V0A;number of BB triggers;events", 34, -1.5, 32.5);
71   fHistV0C = new TH1F("fHistV0C", "V0C;number of BB triggers;events", 34, -1.5, 32.5);
72   fHistZDC = new TH1F("fHistZDC", "ZDC;trigger bits;events", 8, -1.5, 6.5);
73   
74   // TODO check limits
75   fHistFMDA = new TH1F("fHistFMDA", "FMDA;combinations above threshold;events", 102, -1.5, 100.5);
76   fHistFMDC = new TH1F("fHistFMDC", "FMDC;combinations above threshold;events", 102, -1.5, 100.5);
77   fHistFMDSingle = new TH1F("fHistFMDSingle", "FMD single;multiplicity value;counts", 1000, 0, 10);
78   fHistFMDSum = new TH1F("fHistFMDSum", "FMD sum;multiplicity value;counts", 1000, 0, 10);
79 }
80
81 //____________________________________________________________________
82 const char* AliTriggerAnalysis::GetTriggerName(Trigger trigger) 
83 {
84   // returns the name of the requested trigger
85   // the returned string will only be valid until the next call to this function [not thread-safe]
86   
87   static TString str;
88   
89   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
90   
91   switch (triggerNoFlags)
92   {
93     case kAcceptAll : str = "ACCEPT ALL (bypass!)"; break;
94     case kMB1 : str = "MB1"; break;
95     case kMB2 : str = "MB2"; break;
96     case kMB3 : str = "MB3"; break;
97     case kSPDGFO : str = "SPD GFO"; break;
98     case kSPDGFOBits : str = "SPD GFO Bits"; break;
99     case kV0A : str = "V0 A"; break;
100     case kV0C : str = "V0 C"; break;
101     case kZDC : str = "ZDC"; break;
102     case kZDCA : str = "ZDC A"; break;
103     case kZDCC : str = "ZDC C"; break;
104     case kFMDA : str = "FMD A"; break;
105     case kFMDC : str = "FMD C"; break;
106     case kFPANY : str = "SPD GFO | V0 | ZDC | FMD"; break;
107     default: str = ""; break;
108   }
109    
110   if (trigger & kOfflineFlag)
111     str += " OFFLINE";  
112   
113   return str;
114 }
115
116 Bool_t AliTriggerAnalysis::IsTriggerFired(const AliESDEvent* aEsd, Trigger trigger) const
117 {
118   // checks if an event has been triggered
119
120   if (trigger & kOfflineFlag)
121     return IsOfflineTriggerFired(aEsd, trigger);
122     
123   return IsTriggerBitFired(aEsd, trigger);
124 }
125
126 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, Trigger trigger) const
127 {
128   // checks if an event is fired using the trigger bits
129
130   return IsTriggerBitFired(aEsd->GetTriggerMask(), trigger);
131 }
132
133 Bool_t AliTriggerAnalysis::IsTriggerBitFired(ULong64_t triggerMask, Trigger trigger) const
134 {
135   // checks if an event is fired using the trigger bits
136   //
137   // this function needs the branch TriggerMask in the ESD
138   
139   // definitions from p-p.cfg
140   ULong64_t spdFO = (1 << 14);
141   ULong64_t v0left = (1 << 10);
142   ULong64_t v0right = (1 << 11);
143
144   switch (trigger)
145   {
146     case kAcceptAll:
147     {
148       return kTRUE;
149       break;
150     }
151     case kMB1:
152     {
153       if (triggerMask & spdFO || ((triggerMask & v0left) || (triggerMask & v0right)))
154         return kTRUE;
155       break;
156     }
157     case kMB2:
158     {
159       if (triggerMask & spdFO && ((triggerMask & v0left) || (triggerMask & v0right)))
160         return kTRUE;
161       break;
162     }
163     case kMB3:
164     {
165       if (triggerMask & spdFO && (triggerMask & v0left) && (triggerMask & v0right))
166         return kTRUE;
167       break;
168     }
169     case kSPDGFO:
170     {
171       if (triggerMask & spdFO)
172         return kTRUE;
173       break;
174     }
175     default:
176       Printf("IsEventTriggered: ERROR: Trigger type %d not implemented in this method", (Int_t) trigger);
177       break;
178   }
179
180   return kFALSE;
181 }
182
183 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, ULong64_t tclass) const
184 {
185   // Checks if corresponding bit in mask is on
186   
187   ULong64_t trigmask = aEsd->GetTriggerMask();
188   return (trigmask & (1ull << (tclass-1)));
189 }
190
191 Bool_t AliTriggerAnalysis::IsOfflineTriggerFired(const AliESDEvent* aEsd, Trigger trigger) const
192 {
193   // checks if an event has been triggered "offline"
194
195   UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
196   
197   switch (triggerNoFlags)
198   {
199     case kAcceptAll:
200     {
201       return kTRUE;
202       break;
203     }
204     case kMB1:
205     {
206       if (SPDGFOTrigger(aEsd, 0) || V0Trigger(aEsd, kASide) || V0Trigger(aEsd, kCSide))
207         return kTRUE;
208       break;
209     }
210     case kMB2:
211     {
212       if (SPDGFOTrigger(aEsd, 0) && (V0Trigger(aEsd, kASide) || V0Trigger(aEsd, kCSide)))
213         return kTRUE;
214       break;
215     }
216     case kMB3:
217     {
218       if (SPDGFOTrigger(aEsd, 0) && V0Trigger(aEsd, kASide) && V0Trigger(aEsd, kCSide))
219         return kTRUE;
220       break;
221     }
222     case kSPDGFO:
223     {
224       if (SPDGFOTrigger(aEsd, 0))
225         return kTRUE;
226       break;
227     }
228     case kSPDGFOBits:
229     {
230       if (SPDGFOTrigger(aEsd, 1))
231         return kTRUE;
232       break;
233     }
234     case kV0A:
235     {
236       if (V0Trigger(aEsd, kASide))
237         return kTRUE;
238       break;
239     }
240     case kV0C:
241     {
242       if (V0Trigger(aEsd, kCSide))
243         return kTRUE;
244       break;
245     }
246     case kZDC:
247     {
248       if (ZDCTrigger(aEsd, kASide) || ZDCTrigger(aEsd, kCentralBarrel) || ZDCTrigger(aEsd, kCSide))
249         return kTRUE;
250       break;
251     }
252     case kZDCA:
253     {
254       if (ZDCTrigger(aEsd, kASide))
255         return kTRUE;
256       break;
257     }
258     case kZDCC:
259     {
260       if (ZDCTrigger(aEsd, kCSide))
261         return kTRUE;
262       break;
263     }
264     case kFMDA:
265     {
266       if (FMDTrigger(aEsd, kASide))
267         return kTRUE;
268       break;
269     }
270     case kFMDC:
271     {
272       if (FMDTrigger(aEsd, kCSide))
273         return kTRUE;
274       break;
275     }
276     case kFPANY:
277     {
278       if (SPDGFOTrigger(aEsd, 0) || V0Trigger(aEsd, kASide) || V0Trigger(aEsd, kCSide) || ZDCTrigger(aEsd, kASide) || ZDCTrigger(aEsd, kCentralBarrel) || ZDCTrigger(aEsd, kCSide) || FMDTrigger(aEsd, kASide) || FMDTrigger(aEsd, kCSide))
279         return kTRUE;
280       break;
281     }
282     default:
283     {
284       AliFatal(Form("Trigger type %d not implemented", triggerNoFlags));
285     }
286   }
287   
288   return kFALSE;
289 }
290
291
292 Bool_t AliTriggerAnalysis::IsTriggerClassFired(const AliESDEvent* aEsd, const Char_t* tclass) const 
293 {
294   // tclass is logical function of inputs, e.g. 01 && 02 || 03 && 11 && 21
295   // = L0 inp 1 && L0 inp 2 || L0 inp 3 && L1 inp 1 && L2 inp 1
296   // NO brackets in logical function !
297   // Spaces between operators and inputs.
298   // Not all logical functions are available in CTP= 
299   // =any function of first 4 inputs; 'AND' of other inputs, check not done
300   // This method will be replaced/complemened by similar one
301   // which works withh class and inputs names as in CTP cfg file
302   
303   TString TClass(tclass);
304   TObjArray* tcltokens = TClass.Tokenize(" ");
305   Char_t level=((TObjString*)tcltokens->At(0))->String()[0];
306   UInt_t input=atoi((((TObjString*)tcltokens->At(0))->String()).Remove(0));
307   Bool_t tcl = IsInputFired(aEsd,level,input);
308  
309   for (Int_t i=1;i<tcltokens->GetEntriesFast();i=i+2) {
310     level=((TObjString*)tcltokens->At(i+1))->String()[0];
311     input=atoi((((TObjString*)tcltokens->At(i+1))->String()).Remove(0));
312     Bool_t inpnext = IsInputFired(aEsd,level,input);
313     Char_t op =((TObjString*)tcltokens->At(i))->String()[0];
314     if (op == '&') tcl=tcl && inpnext;
315     else if (op == '|') tcl =tcl || inpnext;
316     else {
317        AliError(Form("Syntax error in %s", tclass));
318        tcltokens->Delete();
319        return kFALSE;
320     }
321   }
322   tcltokens->Delete();
323   return tcl;
324 }
325
326 Bool_t AliTriggerAnalysis::IsInputFired(const AliESDEvent* aEsd, Char_t level, UInt_t input) const
327 {
328   // Checks trigger input of any level
329   
330   switch (level)
331   {
332     case '0': return IsL0InputFired(aEsd,input);
333     case '1': return IsL1InputFired(aEsd,input);
334     case '2': return IsL2InputFired(aEsd,input);
335     default:
336       AliError(Form("Wrong level %i",level));
337       return kFALSE;
338   }
339 }
340
341 Bool_t AliTriggerAnalysis::IsL0InputFired(const AliESDEvent* aEsd, UInt_t input) const 
342 {
343   // Checks if corresponding bit in mask is on
344   
345   UInt_t inpmask = aEsd->GetHeader()->GetL0TriggerInputs();
346   return (inpmask & (1<<(input-1)));
347 }
348
349 Bool_t AliTriggerAnalysis::IsL1InputFired(const AliESDEvent* aEsd, UInt_t input) const
350 {
351   // Checks if corresponding bit in mask is on
352   
353   UInt_t inpmask = aEsd->GetHeader()->GetL1TriggerInputs();
354   return (inpmask & (1<<(input-1)));
355 }
356
357 Bool_t AliTriggerAnalysis::IsL2InputFired(const AliESDEvent* aEsd, UInt_t input) const 
358 {
359   // Checks if corresponding bit in mask is on
360   
361   UInt_t inpmask = aEsd->GetHeader()->GetL2TriggerInputs();
362   return (inpmask & (1<<(input-1)));
363 }
364
365 void AliTriggerAnalysis::FillHistograms(const AliESDEvent* aEsd) 
366 {
367   // fills the histograms with the info from the ESD
368   
369   fHistSPD->Fill(SPDFiredChips(aEsd, 0));
370   fHistBitsSPD->Fill(SPDFiredChips(aEsd, 0), SPDFiredChips(aEsd, 1, kTRUE));
371   
372   fHistV0A->Fill(V0BBTriggers(aEsd, kASide));  
373   fHistV0C->Fill(V0BBTriggers(aEsd, kCSide));
374   
375   AliESDZDC* zdcData = aEsd->GetESDZDC();
376   if (zdcData)
377   {
378     UInt_t quality = zdcData->GetESDQuality();
379     
380     // from Nora's presentation, general first physics meeting 16.10.09
381     static UInt_t zpc  = 0x20;
382     static UInt_t znc  = 0x10;
383     static UInt_t zem1 = 0x08;
384     static UInt_t zem2 = 0x04;
385     static UInt_t zpa  = 0x02;
386     static UInt_t zna  = 0x01;
387    
388     fHistZDC->Fill(1, quality & zna);
389     fHistZDC->Fill(2, quality & zpa);
390     fHistZDC->Fill(3, quality & zem2);
391     fHistZDC->Fill(4, quality & zem1);
392     fHistZDC->Fill(5, quality & znc);
393     fHistZDC->Fill(6, quality & zpc);
394   }
395   else
396   {
397     fHistZDC->Fill(-1);
398     AliError("AliESDZDC not available");
399   }
400   
401   fHistFMDA->Fill(FMDHitCombinations(aEsd, kASide, kTRUE));
402   fHistFMDC->Fill(FMDHitCombinations(aEsd, kCSide, kTRUE));
403 }
404
405 Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHist) const
406 {
407   // returns the number of fired chips in the SPD
408   //
409   // origin = 0 --> aEsd->GetMultiplicity()->GetNumberOfFiredChips() (filled from clusters)
410   // origin = 1 --> aEsd->GetMultiplicity()->TestFastOrFiredChips() (from hardware bits)
411   
412   const AliMultiplicity* mult = aEsd->GetMultiplicity();
413   if (!mult)
414   {
415     AliError("AliMultiplicity not available");
416     return -1;
417   }
418   
419   if (origin == 0)
420     return mult->GetNumberOfFiredChips(0) + mult->GetNumberOfFiredChips(1);
421     
422   if (origin == 1)
423   {
424     Int_t nChips = 0;
425     for (Int_t i=0; i<1200; i++)
426       if (mult->TestFastOrFiredChips(i) == kTRUE)
427       {
428         nChips++;
429         if (fillHist)
430           fHistFiredBitsSPD->Fill(i);
431       }
432     return nChips;
433   }
434   
435   return -1;
436 }
437
438 Bool_t AliTriggerAnalysis::SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin) const
439 {
440   // Returns if the SPD gave a global Fast OR trigger
441   
442   Int_t firedChips = SPDFiredChips(aEsd, origin);
443   
444   if (firedChips >= fSPDGFOThreshold)
445     return kTRUE;
446   return kFALSE;
447 }
448
449 Int_t AliTriggerAnalysis::V0BBTriggers(const AliESDEvent* aEsd, AliceSide side) const
450 {
451   // returns the number of BB triggers in V0A | V0C
452   
453   AliESDVZERO* v0Data = aEsd->GetVZEROData();
454   if (!v0Data)
455   {
456     AliError("AliESDVZERO not available");
457     return -1;
458   }
459   
460   Int_t count = 0;
461   for (Int_t i=0; i<32; i++)
462   {
463     if (side == kASide && v0Data->BBTriggerV0A(i))
464       count++;
465     if (side == kCSide && v0Data->BBTriggerV0C(i))
466       count++;
467   }
468   
469   return count;
470 }
471
472 Bool_t AliTriggerAnalysis::V0Trigger(const AliESDEvent* aEsd, AliceSide side) const
473 {
474   // Returns if the V0 triggered
475   
476   Int_t count = V0BBTriggers(aEsd, side);
477   
478   if (side == kASide && count >= fV0AThreshold)
479     return kTRUE;
480   if (side == kCSide && count >= fV0CThreshold)
481     return kTRUE;
482   return kFALSE;
483 }
484
485 Bool_t AliTriggerAnalysis::ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) const
486 {
487   // Returns if ZDC triggered
488   
489   AliESDZDC* zdcData = aEsd->GetESDZDC();
490   if (!zdcData)
491   {
492     AliError("AliESDZDC not available");
493     return kFALSE;
494   }
495   
496   UInt_t quality = zdcData->GetESDQuality();
497   
498   // from Nora's presentation, general first physics meeting 16.10.09
499   static UInt_t zpc  = 0x20;
500   static UInt_t znc  = 0x10;
501   static UInt_t zem1 = 0x08;
502   static UInt_t zem2 = 0x04;
503   static UInt_t zpa  = 0x02;
504   static UInt_t zna  = 0x01;
505   
506   if (side == kASide && ((quality & zpa) || (quality & zna)))
507     return kTRUE;
508   if (side == kCentralBarrel && ((quality & zem1) || (quality & zem2)))
509     return kTRUE;
510   if (side == kCSide && ((quality & zpc) || (quality & znc)))
511     return kTRUE;
512   
513   return kFALSE;
514 }
515
516 Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHistograms) const
517 {
518   // returns number of hit combinations agove threshold
519   //
520   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
521
522   // Workaround for AliESDEvent::GetFMDData is not const!
523   const AliESDFMD* fmdData = (const_cast<AliESDEvent*>(aEsd))->GetFMDData();
524   if (!fmdData)
525   {
526     AliError("AliESDFMD not available");
527     return -1;
528   }
529
530   Int_t detFrom = (side == kASide) ? 1 : 3;
531   Int_t detTo   = (side == kASide) ? 2 : 3;
532
533   Int_t triggers = 0;
534   Float_t totalMult = 0;
535   for (UShort_t det=detFrom;det<=detTo;det++) {
536     Int_t nRings = (det == 1 ? 1 : 2);
537     for (UShort_t ir = 0; ir < nRings; ir++) {
538       Char_t   ring = (ir == 0 ? 'I' : 'O');
539       UShort_t nsec = (ir == 0 ? 20  : 40);
540       UShort_t nstr = (ir == 0 ? 512 : 256);
541       for (UShort_t sec =0; sec < nsec;  sec++) {
542         for (UShort_t strip = 0; strip < nstr; strip++) {
543           Float_t mult = fmdData->Multiplicity(det,ring,sec,strip);
544           if (mult == AliESDFMD::kInvalidMult) continue;
545           
546           if (fillHistograms)
547             fHistFMDSingle->Fill(mult);
548           
549           if (mult > fFMDLowCut)
550             totalMult = totalMult + mult;
551           else
552           {
553             if (totalMult > fFMDHitCut)
554               triggers++;
555               
556             if (fillHistograms)
557               fHistFMDSum->Fill(totalMult);
558               
559             totalMult = 0;
560           }
561         }
562       }
563     }
564   }
565   
566   return triggers;
567 }
568
569 Bool_t AliTriggerAnalysis::FMDTrigger(const AliESDEvent* aEsd, AliceSide side) const
570 {
571   // Returns if the FMD triggered
572   //
573   // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
574
575   Int_t triggers = FMDHitCombinations(aEsd, side, kFALSE);
576   
577   if (triggers > 0)
578     return kTRUE;
579     
580   return kFALSE;
581 }
582
583 Long64_t AliTriggerAnalysis::Merge(TCollection* list)
584 {
585   // Merge a list of AliMultiplicityCorrection objects with this (needed for
586   // PROOF).
587   // Returns the number of merged objects (including this).
588
589   if (!list)
590     return 0;
591
592   if (list->IsEmpty())
593     return 1;
594
595   TIterator* iter = list->MakeIterator();
596   TObject* obj;
597
598   // collections of all histograms
599   const Int_t nHists = 10;
600   TList collections[nHists];
601
602   Int_t count = 0;
603   while ((obj = iter->Next())) {
604
605     AliTriggerAnalysis* entry = dynamic_cast<AliTriggerAnalysis*> (obj);
606     if (entry == 0) 
607       continue;
608
609     collections[0].Add(entry->fHistSPD);
610     collections[1].Add(entry->fHistV0A);
611     collections[2].Add(entry->fHistV0C);
612     collections[3].Add(entry->fHistZDC);
613     collections[4].Add(entry->fHistFMDA);
614     collections[5].Add(entry->fHistFMDC);
615     collections[6].Add(entry->fHistFMDSingle);
616     collections[7].Add(entry->fHistFMDSum);
617     collections[8].Add(entry->fHistBitsSPD);
618     collections[9].Add(entry->fHistFiredBitsSPD);
619
620     count++;
621   }
622
623   fHistSPD->Merge(&collections[0]);
624   fHistV0A->Merge(&collections[1]);
625   fHistV0C->Merge(&collections[2]);
626   fHistZDC->Merge(&collections[3]);
627   fHistFMDA->Merge(&collections[4]);
628   fHistFMDC->Merge(&collections[5]);
629   fHistFMDSingle->Merge(&collections[6]);
630   fHistFMDSum->Merge(&collections[7]);
631   fHistBitsSPD->Merge(&collections[8]);
632   fHistFiredBitsSPD->Merge(&collections[9]);
633
634   delete iter;
635
636   return count+1;
637 }
638
639 void AliTriggerAnalysis::WriteHistograms() const
640 {
641   // write histograms to current directory
642   
643   if (!fHistSPD)
644     return;
645     
646   fHistSPD->Write();
647   fHistBitsSPD->Write();
648   fHistFiredBitsSPD->Write();
649   fHistV0A->Write();
650   fHistV0C->Write();
651   fHistZDC->Write();
652   fHistFMDA->Write();
653   fHistFMDC->Write();
654   fHistFMDSingle->Write();
655   fHistFMDSum->Write();
656 }