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