1 /* $Id: AliTriggerAnalysis.cxx 35782 2009-10-22 11:54:31Z jgrosseo $ */
3 /**************************************************************************
4 * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
6 * Author: The ALICE Off-line Project. *
7 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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 //-------------------------------------------------------------------------
25 #include <Riostream.h>
29 #include <TIterator.h>
30 #include "TParameter.h"
33 #include <AliTriggerAnalysis.h>
37 #include <AliESDEvent.h>
39 #include <AliMultiplicity.h>
40 #include <AliESDVZERO.h>
41 #include <AliESDZDC.h>
42 #include <AliESDFMD.h>
44 ClassImp(AliTriggerAnalysis)
46 AliTriggerAnalysis::AliTriggerAnalysis() :
65 AliTriggerAnalysis::~AliTriggerAnalysis()
75 if (fHistFiredBitsSPD)
77 delete fHistFiredBitsSPD;
78 fHistFiredBitsSPD = 0;
113 delete fHistFMDSingle;
125 delete fTriggerClasses;
130 void AliTriggerAnalysis::EnableHistograms()
132 // creates the monitoring histograms
134 // do not add this hists to the directory
135 Bool_t oldStatus = TH1::AddDirectoryStatus();
136 TH1::AddDirectory(kFALSE);
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);
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);
150 fTriggerClasses = new TMap;
151 fTriggerClasses->SetOwner();
153 TH1::AddDirectory(oldStatus);
156 //____________________________________________________________________
157 const char* AliTriggerAnalysis::GetTriggerName(Trigger trigger)
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]
164 UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
166 switch (triggerNoFlags)
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;
187 if (trigger & kOfflineFlag)
193 Bool_t AliTriggerAnalysis::IsTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
195 // checks if an event has been triggered
197 if (trigger & kOfflineFlag)
198 return IsOfflineTriggerFired(aEsd, trigger);
200 return IsTriggerBitFired(aEsd, trigger);
203 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, Trigger trigger) const
205 // checks if an event is fired using the trigger bits
207 return IsTriggerBitFired(aEsd->GetTriggerMask(), trigger);
210 Bool_t AliTriggerAnalysis::IsTriggerBitFired(ULong64_t triggerMask, Trigger trigger) const
212 // checks if an event is fired using the trigger bits
214 // this function needs the branch TriggerMask in the ESD
216 // definitions from p-p.cfg
217 ULong64_t spdFO = (1 << 14);
218 ULong64_t v0left = (1 << 10);
219 ULong64_t v0right = (1 << 11);
230 if (triggerMask & spdFO || ((triggerMask & v0left) || (triggerMask & v0right)))
236 if (triggerMask & spdFO && ((triggerMask & v0left) || (triggerMask & v0right)))
242 if (triggerMask & spdFO && (triggerMask & v0left) && (triggerMask & v0right))
248 if (triggerMask & spdFO)
253 Printf("IsEventTriggered: ERROR: Trigger type %d not implemented in this method", (Int_t) trigger);
260 Bool_t AliTriggerAnalysis::IsTriggerBitFired(const AliESDEvent* aEsd, ULong64_t tclass) const
262 // Checks if corresponding bit in mask is on
264 ULong64_t trigmask = aEsd->GetTriggerMask();
265 return (trigmask & (1ull << (tclass-1)));
268 Bool_t AliTriggerAnalysis::IsOfflineTriggerFired(const AliESDEvent* aEsd, Trigger trigger)
270 // checks if an event has been triggered "offline"
272 UInt_t triggerNoFlags = (UInt_t) trigger % (UInt_t) kStartOfFlags;
274 switch (triggerNoFlags)
283 if (SPDGFOTrigger(aEsd, 0) || V0Trigger(aEsd, kASide) == kV0BB || V0Trigger(aEsd, kCSide) == kV0BB)
289 if (SPDGFOTrigger(aEsd, 0) && (V0Trigger(aEsd, kASide) == kV0BB || V0Trigger(aEsd, kCSide) == kV0BB))
295 if (SPDGFOTrigger(aEsd, 0) && V0Trigger(aEsd, kASide) == kV0BB && V0Trigger(aEsd, kCSide) == kV0BB)
301 if (SPDGFOTrigger(aEsd, 0))
307 if (SPDGFOTrigger(aEsd, 1))
313 if (V0Trigger(aEsd, kASide) == kV0BB)
319 if (V0Trigger(aEsd, kCSide) == kV0BB)
325 if (V0Trigger(aEsd, kASide) == kV0BG)
331 if (V0Trigger(aEsd, kCSide) == kV0BG)
337 if (ZDCTrigger(aEsd, kASide) || ZDCTrigger(aEsd, kCentralBarrel) || ZDCTrigger(aEsd, kCSide))
343 if (ZDCTrigger(aEsd, kASide))
349 if (ZDCTrigger(aEsd, kCSide))
355 if (FMDTrigger(aEsd, kASide))
361 if (FMDTrigger(aEsd, kCSide))
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))
373 AliFatal(Form("Trigger type %d not implemented", triggerNoFlags));
381 Bool_t AliTriggerAnalysis::IsTriggerClassFired(const AliESDEvent* aEsd, const Char_t* tclass) const
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
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);
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;
406 AliError(Form("Syntax error in %s", tclass));
415 Bool_t AliTriggerAnalysis::IsInputFired(const AliESDEvent* aEsd, Char_t level, UInt_t input) const
417 // Checks trigger input of any level
421 case '0': return IsL0InputFired(aEsd,input);
422 case '1': return IsL1InputFired(aEsd,input);
423 case '2': return IsL2InputFired(aEsd,input);
425 AliError(Form("Wrong level %i",level));
430 Bool_t AliTriggerAnalysis::IsL0InputFired(const AliESDEvent* aEsd, UInt_t input) const
432 // Checks if corresponding bit in mask is on
434 UInt_t inpmask = aEsd->GetHeader()->GetL0TriggerInputs();
435 return (inpmask & (1<<(input-1)));
438 Bool_t AliTriggerAnalysis::IsL1InputFired(const AliESDEvent* aEsd, UInt_t input) const
440 // Checks if corresponding bit in mask is on
442 UInt_t inpmask = aEsd->GetHeader()->GetL1TriggerInputs();
443 return (inpmask & (1<<(input-1)));
446 Bool_t AliTriggerAnalysis::IsL2InputFired(const AliESDEvent* aEsd, UInt_t input) const
448 // Checks if corresponding bit in mask is on
450 UInt_t inpmask = aEsd->GetHeader()->GetL2TriggerInputs();
451 return (inpmask & (1<<(input-1)));
454 void AliTriggerAnalysis::FillHistograms(const AliESDEvent* aEsd)
456 // fills the histograms with the info from the ESD
458 fHistBitsSPD->Fill(SPDFiredChips(aEsd, 0), SPDFiredChips(aEsd, 1, kTRUE));
460 V0Trigger(aEsd, kASide, kTRUE);
461 V0Trigger(aEsd, kCSide, kTRUE);
463 AliESDZDC* zdcData = aEsd->GetESDZDC();
466 UInt_t quality = zdcData->GetESDQuality();
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;
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);
486 AliError("AliESDZDC not available");
489 fHistFMDA->Fill(FMDHitCombinations(aEsd, kASide, kTRUE));
490 fHistFMDC->Fill(FMDHitCombinations(aEsd, kCSide, kTRUE));
493 void AliTriggerAnalysis::FillTriggerClasses(const AliESDEvent* aEsd)
495 // fills trigger classes map
497 TParameter<Long64_t>* count = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(aEsd->GetFiredTriggerClasses().Data()));
500 count = new TParameter<Long64_t>(aEsd->GetFiredTriggerClasses(), 0);
501 fTriggerClasses->Add(new TObjString(aEsd->GetFiredTriggerClasses().Data()), count);
503 count->SetVal(count->GetVal() + 1);
505 // TODO add first and last orbit number here
508 Int_t AliTriggerAnalysis::SPDFiredChips(const AliESDEvent* aEsd, Int_t origin, Bool_t fillHists)
510 // returns the number of fired chips in the SPD
512 // origin = 0 --> aEsd->GetMultiplicity()->GetNumberOfFiredChips() (filled from clusters)
513 // origin = 1 --> aEsd->GetMultiplicity()->TestFastOrFiredChips() (from hardware bits)
515 const AliMultiplicity* mult = aEsd->GetMultiplicity();
518 AliError("AliMultiplicity not available");
523 return mult->GetNumberOfFiredChips(0) + mult->GetNumberOfFiredChips(1);
528 for (Int_t i=0; i<1200; i++)
529 if (mult->TestFastOrFiredChips(i) == kTRUE)
533 fHistFiredBitsSPD->Fill(i);
541 Bool_t AliTriggerAnalysis::SPDGFOTrigger(const AliESDEvent* aEsd, Int_t origin)
543 // Returns if the SPD gave a global Fast OR trigger
545 Int_t firedChips = SPDFiredChips(aEsd, origin);
547 if (firedChips >= fSPDGFOThreshold)
552 AliTriggerAnalysis::V0Decision AliTriggerAnalysis::V0Trigger(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists)
554 // Returns the V0 trigger decision in V0A | V0C
556 // Based on algorithm by Cvetan Cheshkov
558 AliESDVZERO* esdV0 = aEsd->GetVZEROData();
561 AliError("AliESDVZERO not available");
573 else if (side == kCSide)
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;
589 weight += timeWeight;
595 time += fV0TimeOffset;
599 if (side == kASide && fHistV0A)
600 fHistV0A->Fill(time);
601 if (side == kCSide && fHistV0C)
602 fHistV0C->Fill(time);
607 if (time > 75 && time < 83)
609 if (time > 54 && time < 57.5)
615 if (time > 75.5 && time < 82)
617 if (time > 69.5 && time < 73)
624 Float_t AliTriggerAnalysis::V0CorrectLeadingTime(Int_t i, Float_t time, Float_t adc) const
626 // Correct for slewing and align the channels
628 // Authors: Cvetan Cheshkov / Raphael Tieulent
630 if (time == 0) return 0;
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};
635 time -= timeShift[i];
637 // Slewing correction
638 if (adc == 0) return time;
640 Float_t p1 = 1.57345e1;
641 Float_t p2 =-4.25603e-1;
643 return (time - p1*TMath::Power(adc,p2));
646 Float_t AliTriggerAnalysis::V0LeadingTimeWeight(Float_t adc) const
648 if (adc < 1e-6) return 0;
651 Float_t p2 =-4.25603e-1;
654 return 1./(p1*p1*TMath::Power(adc,2.*(p2-1.))+p3*p3);
657 Bool_t AliTriggerAnalysis::ZDCTrigger(const AliESDEvent* aEsd, AliceSide side) const
659 // Returns if ZDC triggered
661 AliESDZDC* zdcData = aEsd->GetESDZDC();
664 AliError("AliESDZDC not available");
668 UInt_t quality = zdcData->GetESDQuality();
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;
678 if (side == kASide && ((quality & zpa) || (quality & zna)))
680 if (side == kCentralBarrel && ((quality & zem1) || (quality & zem2)))
682 if (side == kCSide && ((quality & zpc) || (quality & znc)))
688 Int_t AliTriggerAnalysis::FMDHitCombinations(const AliESDEvent* aEsd, AliceSide side, Bool_t fillHists)
690 // returns number of hit combinations agove threshold
692 // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
694 // Workaround for AliESDEvent::GetFMDData is not const!
695 const AliESDFMD* fmdData = (const_cast<AliESDEvent*>(aEsd))->GetFMDData();
698 AliError("AliESDFMD not available");
702 Int_t detFrom = (side == kASide) ? 1 : 3;
703 Int_t detTo = (side == kASide) ? 2 : 3;
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;
719 fHistFMDSingle->Fill(mult);
721 if (mult > fFMDLowCut)
722 totalMult = totalMult + mult;
725 if (totalMult > fFMDHitCut)
729 fHistFMDSum->Fill(totalMult);
741 Bool_t AliTriggerAnalysis::FMDTrigger(const AliESDEvent* aEsd, AliceSide side)
743 // Returns if the FMD triggered
745 // Authors: FMD team, Hans Dalsgaard (code merged from FMD/AliFMDOfflineTrigger)
747 Int_t triggers = FMDHitCombinations(aEsd, side, kFALSE);
755 Long64_t AliTriggerAnalysis::Merge(TCollection* list)
757 // Merge a list of AliMultiplicityCorrection objects with this (needed for
759 // Returns the number of merged objects (including this).
767 TIterator* iter = list->MakeIterator();
770 // collections of all histograms
771 const Int_t nHists = 9;
772 TList collections[nHists];
775 while ((obj = iter->Next())) {
777 AliTriggerAnalysis* entry = dynamic_cast<AliTriggerAnalysis*> (obj);
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);
792 // merge fTriggerClasses
793 TIterator* iter2 = entry->fTriggerClasses->MakeIterator();
795 while ((obj = dynamic_cast<TObjString*> (iter2->Next())))
797 TParameter<Long64_t>* param2 = dynamic_cast<TParameter<Long64_t>*> (entry->fTriggerClasses->GetValue(obj));
799 TParameter<Long64_t>* param1 = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj));
802 param1->SetVal(param1->GetVal() + param2->GetVal());
806 param1 = dynamic_cast<TParameter<Long64_t>*> (param2->Clone());
807 fTriggerClasses->Add(new TObjString(obj->String()), param1);
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++]);
832 void AliTriggerAnalysis::SaveHistograms() const
834 // write histograms to current directory
839 fHistBitsSPD->Write();
840 fHistBitsSPD->ProjectionX();
841 fHistBitsSPD->ProjectionY();
842 fHistFiredBitsSPD->Write();
848 fHistFMDSingle->Write();
849 fHistFMDSum->Write();
851 fTriggerClasses->Write("fTriggerClasses", TObject::kSingleKey);
854 void AliTriggerAnalysis::PrintTriggerClasses() const
856 // print trigger classes
858 Printf("Trigger Classes:");
860 TIterator* iter = fTriggerClasses->MakeIterator();
862 while ((obj = dynamic_cast<TObjString*> (iter->Next())))
864 TParameter<Long64_t>* param = dynamic_cast<TParameter<Long64_t>*> (fTriggerClasses->GetValue(obj));
866 Printf("%s: %ld triggers", obj->String().Data(), param->GetVal());