1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
19 // AliAnalysisTriggerScalers : utility class to play with OCDB scalers.
21 // Can use it to e.g. :
23 // - get the integrated luminosity for a given period
24 // (see IntegratedLuminosity method)
26 // - plot the trigger rates (see PlotTriggerEvolution)
28 // Please note that this class is doing an OCDB "scan" (for as many run
29 // as you put in the SetRunList method), so it can be really long
30 // if you're using the raw:// OCDB. If you're planning on using this
31 // class for anything more than a brief test, you'd better think
32 // or making a local copy of the OCDB objects required by this class :
37 // and also (to get the fill and period ranges, mainly for drawing)
43 // Please note also that this class evolved from a set a quick macros to
44 // follow the luminosity and trigger rates during data taking to this stage.
46 // Now (Feb 2013) would probably be a good time to regroup a bit, think about it and
47 // make it more general, more robust and just less "had-hoc"...
49 // Ideas for improvement :
51 // - make sure the "what" parameter mean the same thing in all methods
52 // and so can take the same values in all methods
54 // - get the lumi trigger name and cross-section from e.g. OADB instead of
57 // - organize the class so that the CTP Scalers and Config can be fetched
58 // from elsewhere (e.g. from a run based container available in some
59 // distant future in the ESDs/AODs ?)
61 // - robustify the PAC fraction computation
64 // L. Aphecetche (Subatech)
67 #include "AliAnalysisTriggerScalers.h"
68 #include "AliCDBEntry.h"
69 #include "AliCDBManager.h"
70 #include "AliGRPObject.h"
72 #include "AliTimeStamp.h"
73 #include "AliTriggerBCMask.h"
74 #include "AliTriggerClass.h"
75 #include "AliTriggerCluster.h"
76 #include "AliTriggerConfiguration.h"
77 #include "AliTriggerDescriptor.h"
78 #include "AliTriggerInput.h"
79 #include "AliTriggerRunScalers.h"
80 #include "AliTriggerScalers.h"
81 #include "AliTriggerScalersESD.h"
82 #include "AliTriggerScalersRecord.h"
83 #include "Riostream.h"
87 #include "TGraphErrors.h"
89 #include "TObjArray.h"
90 #include "TObjString.h"
94 #include "AliLHCData.h"
100 using std::make_pair;
101 ClassImp(AliAnalysisTriggerScalers)
105 //______________________________________________________________________________
106 UChar_t GetIndex(ULong64_t mask)
108 for ( Int_t i = 0; i < 64; ++i )
110 if ( mask & ( 1ull << i ) ) return i+1;
115 //______________________________________________________________________________
116 TCanvas* NewCanvas(const char* name)
118 TCanvas* c = new TCanvas(name,name);
119 c->SetLeftMargin(0.15);
123 //______________________________________________________________________________
124 void TimeAxis(TGraph* g)
126 // g->GetXaxis()->SetTimeDisplay(1);
127 // g->GetXaxis()->SetTimeFormat("%d/%m %H:%M%F2010-12-31 24:00:00");
128 // g->GetXaxis()->SetNdivisions(505);
129 g->GetYaxis()->SetTitleOffset(1.5);
131 g->GetXaxis()->SetTimeDisplay(1);
132 // g->GetXaxis()->SetTimeFormat("%d/%m %H:%M%F2010-12-31 24:00:00");
133 g->GetXaxis()->SetTimeFormat("%d/%m %H:%M");
134 g->GetXaxis()->SetTimeOffset(0,"gmt");
135 g->GetXaxis()->SetNdivisions(505);
141 //_____________________________________________________________________________
142 AliAnalysisTriggerScalers::AliAnalysisTriggerScalers(Int_t runNumber, const char* ocdbPath) : fRunList(), fOCDBPath(ocdbPath), fShouldCorrectForPileUp(kTRUE), fCrossSectionUnit("UB")
144 // ctor using a single run number
145 SetRunList(runNumber);
148 //_____________________________________________________________________________
149 AliAnalysisTriggerScalers::AliAnalysisTriggerScalers(const std::vector<int>& runs, const char* ocdbPath) : fRunList(), fOCDBPath(ocdbPath), fShouldCorrectForPileUp(kTRUE), fCrossSectionUnit("UB")
151 // ctor using a vector of run numbers
155 //_____________________________________________________________________________
156 AliAnalysisTriggerScalers::AliAnalysisTriggerScalers(const std::set<int>& runs, const char* ocdbPath) : fRunList(), fOCDBPath(ocdbPath), fShouldCorrectForPileUp(kTRUE), fCrossSectionUnit("UB")
158 // ctor using a set of run numbers
162 //_____________________________________________________________________________
163 AliAnalysisTriggerScalers::AliAnalysisTriggerScalers(const char* runlist, const char* ocdbPath) : fRunList(), fOCDBPath(ocdbPath), fShouldCorrectForPileUp(kTRUE), fCrossSectionUnit("UB")
165 // ctor from a run list (txt file)
169 //_____________________________________________________________________________
170 AliAnalysisTriggerScalers::~AliAnalysisTriggerScalers()
175 //______________________________________________________________________________
176 Bool_t AliAnalysisTriggerScalers::CheckRecord(const AliTriggerScalersRecord& record,
180 UInt_t timelapse) const
182 /// Check if a trigger scaler record should be kept for our
183 /// luminosity computations
185 const AliTriggerScalers* scaler = record.GetTriggerScalersForClass(index);
187 UInt_t a = scaler->GetLOCA() - refa;
188 UInt_t b = scaler->GetLOCB() - refb;
192 if ( b <= 2 || ( a <= 2 ) || timelapse <= 9 ) // empirical cuts
201 //______________________________________________________________________________
202 void AliAnalysisTriggerScalers::DrawFill(Int_t run1, Int_t run2, double ymin, double ymax, const char* label, Int_t color)
204 // Draw a yellow box to indicate a run range
206 AliDebugClass(1,Form("RUN1 %09d RUN2 %09d YMIN %e YMAX %e %s",
207 run1,run2,ymin,ymax,label));
208 TBox* b = new TBox(run1*1.0,ymin,run2*1.0,ymax);
209 b->SetFillColor(color);
211 TText* text = new TText((run1+run2)/2.0,ymin + (ymax-ymin)*0.85,label);
212 text->SetTextSize(0.025);
213 text->SetTextFont(42);
214 text->SetTextAlign(23);
215 text->SetTextAngle(45);
219 //_____________________________________________________________________________
220 void AliAnalysisTriggerScalers::DrawFills(Double_t ymin, Double_t ymax, Int_t color)
222 /// Draw the fill ranges corresponding to the list of runs
223 /// Note that this method will scan the OCDB to get the run -> fill number relationship,
224 /// so it's better in this case to use a local copy of the OCDB if you have one. Otherwise
227 std::map<int, std::pair<int,int> > fills;
229 GetFillBoundaries(fills);
231 std::map<int, std::pair<int,int> >::const_iterator it;
233 for ( it = fills.begin(); it != fills.end(); ++it )
235 const std::pair<int,int>& p = it->second;
237 fillnumber.Form("%d",it->first);
238 DrawFill(p.first,p.second,ymin,ymax,fillnumber.Data(),color);
242 //_____________________________________________________________________________
243 void AliAnalysisTriggerScalers::DrawPeriods(Double_t ymin, Double_t ymax, Int_t color)
245 /// Draw the period ranges corresponding to the list of runs
246 /// Note that this method will scan the OCDB to get the run -> fill number relationship,
247 /// so it's better in this case to use a local copy of the OCDB if you have one. Otherwise
250 std::map<std::string, std::pair<int,int> > periods;
252 GetLHCPeriodBoundaries(periods);
254 for ( std::map<std::string, std::pair<int,int> >::const_iterator it = periods.begin(); it != periods.end(); ++it )
256 const std::pair<int,int>& p = it->second;
257 DrawFill(p.first,p.second,ymin,ymax,it->first.c_str(),color);
262 //_____________________________________________________________________________
263 void AliAnalysisTriggerScalers::GetCTPObjects(Int_t runNumber,
264 AliTriggerConfiguration*& tc,
265 AliTriggerRunScalers*& trs,
266 AliLHCData*& lhc) const
268 /// For a given run, get the CTP objects we need to do our work
270 /// FIXME: this is the method that should probably be virtualized so
271 /// we can get those objects from either OCDB or run based container in AODs/ESDs
273 tc = static_cast<AliTriggerConfiguration*>(GetOCDBObject("GRP/CTP/Config",runNumber));
274 trs = static_cast<AliTriggerRunScalers*>(GetOCDBObject("GRP/CTP/Scalers",runNumber));
275 lhc = static_cast<AliLHCData*>(GetOCDBObject("GRP/GRP/LHCData",runNumber));
278 //_____________________________________________________________________________
279 void AliAnalysisTriggerScalers::GetFillBoundaries(std::map<int, std::pair<int,int> >& fills)
281 /// Given a list of runs get the fills they are in
285 for ( std::vector<int>::const_iterator it = fRunList.begin(); it != fRunList.end(); ++it )
289 int fill = GetFillNumberFromRunNumber(runnumber);
291 if (fills.count(fill))
293 std::pair<int,int>& p = fills[fill];
294 p.first = TMath::Min(runnumber,p.first);
295 p.second = TMath::Max(runnumber,p.second);
299 fills[fill] = make_pair(runnumber,runnumber);
305 //______________________________________________________________________________
306 Int_t AliAnalysisTriggerScalers::GetFillNumberFromRunNumber(Int_t runNumber)
308 /// Get the fill number of a run
310 AliLHCData* lhcData = static_cast<AliLHCData*>(GetOCDBObject("GRP/GRP/LHCData",runNumber));
313 Int_t fillNumber = lhcData->GetFillNumber();
314 if ( fillNumber == 0 && runNumber == 189694)
316 // manual hack because GRP info incorrect for this run ?
325 //______________________________________________________________________________
326 TObject* AliAnalysisTriggerScalers::GetOCDBObject(const char* path, Int_t runNumber) const
328 // Helper method to get an object from the OCDB
330 if ( !AliCDBManager::Instance()->IsDefaultStorageSet() )
332 AliInfo(Form("Setting OCDB default storage to %s",fOCDBPath.Data()));
333 AliCDBManager::Instance()->SetDefaultStorage(fOCDBPath.Data());
336 AliCDBManager::Instance()->SetRun(runNumber);
338 gErrorIgnoreLevel=kWarning; // to avoid all the TAlienFile::Open messages...
340 AliCDBEntry* e = AliCDBManager::Instance()->Get(path);
341 return e ? e->GetObject() : 0x0;
345 //______________________________________________________________________________
346 TString AliAnalysisTriggerScalers::GetLHCPeriodFromRunNumber(Int_t runNumber) const
348 /// Get a list of LHC periods from a list of run numbers
350 AliGRPObject* grp = static_cast<AliGRPObject*>(GetOCDBObject("GRP/GRP/Data",runNumber));
354 return grp->GetLHCPeriod();
357 //_____________________________________________________________________________
359 AliAnalysisTriggerScalers::GetLuminosityTriggerAndCrossSection(Int_t runNumber,
360 TString& lumiTriggerClassName,
361 Double_t& lumiTriggerCrossSection,
362 Double_t& lumiTriggerCrossSectionError) const
364 /// For one given run, get the trigger class name to be used as the luminometer,
365 /// and its corresponding cross-section
367 /// FIXME: all is hard-coded for now... (use OADB for this ?)
369 lumiTriggerClassName="";
370 lumiTriggerCrossSection=0.0;
371 lumiTriggerCrossSectionError=0.0; // FIXME
373 TString period = GetLHCPeriodFromRunNumber(runNumber);
375 if ( period.BeginsWith("LHC11") )
377 AliError("Not implemented yet");
379 else if ( period.BeginsWith("LHC12") )
381 if ( period == "LHC12h" || period == "LHC12i" )
383 // pp 8 TeV main-satellites
384 lumiTriggerClassName = "C0TVX-S-NOPF-ALLNOTRD";
385 lumiTriggerCrossSection = 28.0; // FIXME: not from a vdM yet
389 AliError("Not implemented yet");
392 else if ( period.BeginsWith("LHC13") )
394 if ( period == "LHC13g" )
397 lumiTriggerClassName = "C0TVX-B-NOPF-ALLNOTRD";
398 lumiTriggerCrossSection = 18.0; // FIXME: not from a vdM yet
402 // p-Pb or Pb-p 5.02 TeV
403 lumiTriggerClassName = "C0TVX-B-NOPF-ALLNOTRD";
404 lumiTriggerCrossSection = 0.755*2000; // FIXME: not from a vdM yet
409 AliError("Not implemented yet");
414 if ( fCrossSectionUnit=="PB" )
418 else if (fCrossSectionUnit=="NB")
422 else if ( fCrossSectionUnit=="UB" )
426 else if ( fCrossSectionUnit=="MB" )
431 lumiTriggerCrossSection *= csMult;
432 lumiTriggerCrossSectionError *= csMult;
435 //_____________________________________________________________________________
436 void AliAnalysisTriggerScalers::GetLHCPeriodBoundaries(std::map<std::string, std::pair<int,int> >& periods)
438 /// Given a list of runs get the fills they are in
442 for ( std::vector<int>::const_iterator it = fRunList.begin(); it != fRunList.end(); ++it )
446 std::string period = GetLHCPeriodFromRunNumber(runnumber).Data();
448 if (periods.count(period))
450 std::pair<int,int>& p = periods[period];
451 p.first = TMath::Min(runnumber,p.first);
452 p.second = TMath::Max(runnumber,p.second);
456 periods[period] = make_pair(runnumber,runnumber);
462 //______________________________________________________________________________
463 Int_t AliAnalysisTriggerScalers::GetTriggerInput(Int_t runNumber, const char* inputname)
465 /// Get one trigger input for a given run
466 AliTriggerConfiguration* tc = static_cast<AliTriggerConfiguration*>(GetOCDBObject("GRP/CTP/Config",runNumber));
469 const TObjArray& inputs = tc->GetInputs();
471 AliTriggerInput* ti = static_cast<AliTriggerInput*>(inputs.FindObject(inputname));
475 return ti->GetSignature();
478 //______________________________________________________________________________
480 AliAnalysisTriggerScalers::GetPauseAndConfigCorrection(Int_t runNumber, const char* triggerClassName)
482 /// Tries to estimate the duration of the Pause and Configure phase(s) in a given run
483 /// Probably not really accurate for the moment though.
488 AliTriggerConfiguration* tc = static_cast<AliTriggerConfiguration*>(GetOCDBObject("GRP/CTP/Config",runNumber));
490 AliTriggerRunScalers* trs = static_cast<AliTriggerRunScalers*>(GetOCDBObject("GRP/CTP/Scalers",runNumber));
491 const TObjArray* scalers = trs->GetScalersRecords();
494 AliTriggerScalersRecord* record;
500 Int_t classIndex = tc->GetClassIndexFromName(triggerClassName);
502 while ( ( record = static_cast<AliTriggerScalersRecord*>(next()) ) )
504 const AliTriggerScalers* scaler = record->GetTriggerScalersForClass(classIndex);
506 const AliTimeStamp* ats = record->GetTimeStamp();
508 UInt_t seconds = ats->GetSeconds();// - TTimeStamp::GetZoneOffset();
510 TTimeStamp ts(seconds,ats->GetMicroSecs());
512 UInt_t l0b = scaler->GetLOCB() - refl0b;
513 UInt_t l0a = scaler->GetLOCA() - refl0a;
514 UInt_t timelapse = seconds - reft;
516 if ( l0b <=2 || timelapse < 9 ) continue;
519 refl0b = scaler->GetLOCB();
520 refl0a = scaler->GetLOCA();
536 return total > 0 ? nopac*1.0/total : 0.0;
539 //______________________________________________________________________________
540 void AliAnalysisTriggerScalers::GetPileUpFactor(Int_t runNumber, const char* triggerClassName,
542 Double_t& value, Double_t& error)
544 /// Get the mean pile-up correction factor for the given run
550 AliError(Form("Cannot work with purity=%f for trigger %s in run %d. Should be strictly positive",purity,triggerClassName,runNumber));
554 AliTriggerRunScalers* trs = static_cast<AliTriggerRunScalers*>(GetOCDBObject("GRP/CTP/Scalers",runNumber));
555 const TObjArray* scalers = trs->GetScalersRecords();
556 const AliTriggerScalersRecord* begin = (AliTriggerScalersRecord*)(scalers->First());
557 const AliTriggerScalersRecord* end = (AliTriggerScalersRecord*)(scalers->Last());
559 time_t duration = TMath::Nint((end->GetTimeStamp()->GetBunchCross() - begin->GetTimeStamp()->GetBunchCross())*AliTimeStamp::fNanosecPerBC*1E-9);
563 AliError(Form("Got zero duration for run %d",runNumber));
567 AliAnalysisTriggerScalerItem* item = GetTriggerScaler(runNumber,"L0B",triggerClassName);
570 AliError(Form("Could not get L0B for trigger %s in run %d",triggerClassName,runNumber));
574 AliLHCData* lhc = static_cast<AliLHCData*>(GetOCDBObject("GRP/GRP/LHCData",runNumber));
576 Int_t nbcx = NumberOfInteractingBunches(*lhc,runNumber);
580 AliError(Form("Cannot work with nbcx=%d for trigger %s in run %d. Should be strictly positive",nbcx,triggerClassName,runNumber));
584 Double_t itemValue = purity*item->Value();
588 AliError(Form("Cannot work with value=%f for trigger %s in run %d. Should be strictly positive",itemValue,triggerClassName,runNumber));
592 Double_t mu = Mu(itemValue/duration,nbcx);
594 value = mu/(1.0-TMath::Exp(-mu));
596 error = 0.0; // FIXME
599 //______________________________________________________________________________
600 AliAnalysisTriggerScalerItem*
601 AliAnalysisTriggerScalers::GetTriggerScaler(Int_t runNumber, const char* level, const char* triggerClassName)
603 // Get the scaler for a given trigger class and a given trigger level.
604 // Returned object must be deleted by the user.
606 AliTriggerConfiguration* tc = static_cast<AliTriggerConfiguration*>(GetOCDBObject("GRP/CTP/Config",runNumber));
607 AliTriggerRunScalers* trs = static_cast<AliTriggerRunScalers*>(GetOCDBObject("GRP/CTP/Scalers",runNumber));
608 AliGRPObject* grp = static_cast<AliGRPObject*>(GetOCDBObject("GRP/GRP/Data",runNumber));
610 if (!tc || !trs || !grp) return 0x0;
612 TString diCurrent(Form("L3:%5.0f;DIP:%5.0f [L3:%d;DIP:%d]",
613 grp->GetL3Current((AliGRPObject::Stats)0),
614 grp->GetDipoleCurrent((AliGRPObject::Stats)0),
615 grp->GetL3Polarity(),
616 grp->GetDipolePolarity()));
618 const TObjArray& trClasses = tc->GetClasses();
620 const TObjArray* scalers = trs->GetScalersRecords();
621 const AliTriggerScalersRecord* begin = (AliTriggerScalersRecord*)(scalers->First());
622 const AliTriggerScalersRecord* end = (AliTriggerScalersRecord*)(scalers->Last());
624 time_t duration = TMath::Nint((end->GetTimeStamp()->GetBunchCross() - begin->GetTimeStamp()->GetBunchCross())*AliTimeStamp::fNanosecPerBC*1E-9);
626 AliTriggerClass* triggerClass = static_cast<AliTriggerClass*>(trClasses.FindObject(triggerClassName));
631 UChar_t index = GetIndex(triggerClass->GetMask());
635 const AliTriggerScalers* sbegin = begin->GetTriggerScalersForClass(index);
636 const AliTriggerScalers* send = end->GetTriggerScalersForClass(index);
638 TString swhat(level);
641 if ( swhat.BeginsWith("L0A") )
643 value = send->GetLOCA() - sbegin->GetLOCA();
645 else if ( swhat.BeginsWith("L0B") )
647 value = send->GetLOCB() - sbegin->GetLOCB();
649 else if ( swhat.BeginsWith("L1A") )
651 value = send->GetL1CA() - sbegin->GetL1CA();
653 else if ( swhat.BeginsWith("L1B") )
655 value = send->GetL1CB() - sbegin->GetL1CB();
657 else if ( swhat.BeginsWith("L2A") )
659 value = send->GetL2CA() - sbegin->GetL2CA();
661 else if ( swhat.BeginsWith("L2B") )
663 value = send->GetL2CB() - sbegin->GetL2CB();
667 // FIXME: get the downscaling factor here for all cases (only BC1 supported here so far)
668 if ( TString(triggerClassName).Contains("_B1") )
670 ds = GetTriggerInput(runNumber,"BC1");
674 swhat.ReplaceAll("RATE","");
676 return new AliAnalysisTriggerScalerItem(runNumber,swhat.Data(),diCurrent.Data(),triggerClass->GetName(),value,triggerClass->GetBCMask(),ds,duration);
679 //______________________________________________________________________________
681 AliAnalysisTriggerScalers::IntegratedLuminosityGraph(const char* triggerClassName, const char* triggerClassNameForPACEstimate)
683 // Get the integrated luminosity graph for one trigger
685 if ( fRunList.empty() )
687 AliError("Cannot work without a run list");
691 Double_t offset(0.0);
693 std::vector<double> vx;
694 std::vector<double> vy;
698 for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
700 TGraph* g = IntegratedLuminosityGraph(fRunList[i],triggerClassName,triggerClassNameForPACEstimate);
704 for ( Int_t j = 0 ; j < g->GetN(); ++j )
708 vy.push_back(y+offset);
716 TGraph* g = new TGraph(vx.size(),&vx[0],&vy[0]);
723 //______________________________________________________________________________
725 AliAnalysisTriggerScalers::IntegratedLuminosityGraph(Int_t runNumber, const char* triggerClassName, const char* triggerClassNameForPACEstimate)
727 // Get the integrated luminosity graph for one trigger for the given run
729 // the triggerClassNameForPACEstimate is only used if triggerClassName is,
730 // for the given run, the same as the luminometer trigger class name,
731 // and should be the trigger class with the higher and non-zero L0A rate (except during PACs)
734 AliTriggerConfiguration* tc(0x0);
735 AliTriggerRunScalers* trs(0x0);
736 AliLHCData* lhc(0x0);
738 GetCTPObjects(runNumber,tc,trs,lhc);
740 if (!tc || !trs || !lhc)
745 const TObjArray* scalerRecords= trs->GetScalersRecords();
746 const TObjArray& triggerClasses = tc->GetClasses();
748 TString lumiTriggerClassName;
749 Double_t lumiTriggerCrossSection(0.0);
750 Double_t lumiTriggerCrossSectionError(0.0);
752 GetLuminosityTriggerAndCrossSection(runNumber,lumiTriggerClassName,lumiTriggerCrossSection,lumiTriggerCrossSectionError);
754 AliTriggerClass* lumiTriggerClass = static_cast<AliTriggerClass*>(triggerClasses.FindObject(lumiTriggerClassName));
756 if (!lumiTriggerClass)
758 AliError(Form("Could not find lumiTriggerClass=%s",lumiTriggerClassName.Data()));
762 if (lumiTriggerCrossSection==0.0)
764 AliError(Form("Could not get cross-section for lumiTriggerClass=%s",lumiTriggerClassName.Data()));
768 AliTriggerClass* triggerClass = static_cast<AliTriggerClass*>(triggerClasses.FindObject(triggerClassName));
772 AliError(Form("Could not find triggerClass=%s",triggerClassName));
776 Int_t nbcx = NumberOfInteractingBunches(*lhc,runNumber);
778 if (nbcx <= 0 && ShouldCorrectForPileUp())
780 AliError("Could not get the number of bunches, so won't be able to correct for pile-up");
783 TIter nextScaler(scalerRecords);
784 AliTriggerScalersRecord* record;
787 UInt_t refl0b[] = { 0, 0 };
788 UInt_t refl0a[] = { 0, 0 };
790 UInt_t l0b[] = { 0, 0 };
791 UInt_t l0a[] = { 0, 0 };
795 UChar_t classIndices[2];
797 // class 0 = class for luminomiter
798 // class 1 = class for livetime estimate
799 // or for PAC estimate if triggerClassNameForPACEstimate
800 // is given and triggerClass is the same as the lumi class
802 classIndices[0] = GetIndex(lumiTriggerClass->GetMask());
803 classIndices[1] = GetIndex(triggerClass->GetMask());
805 Bool_t sameClass = ( classIndices[0] == classIndices[1] );
806 Bool_t pacRemoval(kFALSE);
808 if ( sameClass && strlen(triggerClassNameForPACEstimate) > 0 )
810 AliTriggerClass* triggerClassForPACEstimate = static_cast<AliTriggerClass*>(triggerClasses.FindObject(triggerClassNameForPACEstimate));
811 if (!triggerClassForPACEstimate)
813 AliError(Form("Could not find triggerClassForPACEstimate=%s. Will not correct for PAC durations",triggerClassNameForPACEstimate));
816 classIndices[1] = GetIndex(triggerClassForPACEstimate->GetMask());
817 sameClass = ( classIndices[0] == classIndices[1] );
824 ULong64_t counter(0);
826 std::vector<Double_t> vseconds;
827 std::vector<Double_t> vlivetime;
828 std::vector<Double_t> vlumi;
830 while ( ( record = static_cast<AliTriggerScalersRecord*>(nextScaler()) ) )
832 const AliTimeStamp* ats = record->GetTimeStamp();
834 UInt_t seconds = ats->GetSeconds();// - TTimeStamp::GetZoneOffset();
836 TTimeStamp ts(seconds,ats->GetMicroSecs());
838 UInt_t timelapse = seconds - reft;
840 Bool_t ok = sameClass || CheckRecord(*record,classIndices[1],refl0b[1],refl0a[1],timelapse);
842 if (ok) reft = seconds;
844 for ( Int_t i = 0; i < 2; ++i )
846 const AliTriggerScalers* scaler = record->GetTriggerScalersForClass(classIndices[i]);
850 l0b[i] = scaler->GetLOCB() - refl0b[i];
856 l0a[i] = scaler->GetLOCA() - refl0a[i];
858 refl0b[i] = scaler->GetLOCB();
859 refl0a[i] = scaler->GetLOCA();
868 Double_t factor(1.0);
870 if ( ShouldCorrectForPileUp() && nbcx > 0 && l0b[0] > 0 )
872 Double_t mu = Mu(l0b[0]/timelapse,nbcx);
873 factor = mu/(1-TMath::Exp(-mu));
876 vseconds.push_back(seconds);
880 if ( ok && !sameClass && !pacRemoval )
882 lt = (l0a[1]*1.0)/l0b[1];
885 counter += TMath::Nint(factor*l0b[0]*lt);
887 vlumi.push_back( counter / lumiTriggerCrossSection );
891 if ( vseconds.empty() ) return 0x0;
893 TGraph* g = new TGraph(vseconds.size(),&vseconds[0],&vlumi[0]);
900 //______________________________________________________________________________
901 void AliAnalysisTriggerScalers::IntegratedLuminosity(const char* triggerList,
902 const char* lumiTrigger,
903 Double_t lumiCrossSection,
904 const char* csvOutputFile,
908 // Compute the luminosity for a set of triggers
910 // for T0 based lumi (end of pp 2012), use lumiTrigger = C0TVX-S-NOPF-ALLNOTRD and lumiCrossSection = 28 mb (and csUnit="nb")
912 // for T0 based lumi (pPb 2013), use lumiTrigger = C0TVX-B-NOPF-ALLNOTRD
913 // and lumiCrossSection = 0.755*2000 mb (and csUnit="mub")
915 // for T0 based lumi (pp 2.76 TeV 2013), use lumiTrigger = C0TVX-B-NOPF-ALLNOTRD
916 // and lumiCrossSection = 18 mb (and csUnit="nb")
920 std::map<std::string,float> lumiPerTrigger;
921 std::map<int, std::map<std::string,float> > lumiPerFillPerTrigger;
923 std::map<std::string,time_t> durationPerTrigger;
925 TString sTriggerList(triggerList);
927 if ( sTriggerList.Length()==0 )
929 if ( lumiCrossSection < 30 && lumiCrossSection > 20 )
931 // assuming it's pp 2012
932 sTriggerList = "CMUL8-S-NOPF-ALLNOTRD,CMUL7-S-NOPF-ALLNOTRD,CMUL8-S-NOPF-MUON,CMUL7-S-NOPF-MUON";
934 // for C0MUL-SC-NOPF-MUON must use C0TVX-SC as lumiTrigger (with same cross-section as C0TVX=28mb)
936 else if ( lumiCrossSection < 30 )
938 // assuming it's pp 2013
939 sTriggerList = "CMSL7-B-NOPF-MUON,CMSH7-B-NOPF-MUON,CMUL7-B-NOPF-MUON,CMSL7-B-NOPF-ALLNOTRD,CMSH7-B-NOPF-ALLNOTRD,CMUL7-B-NOPF-ALLNOTRD";
943 sTriggerList = "CMSL7-B-NOPF-MUON,CMSH7-B-NOPF-MUON,CMUL7-B-NOPF-MUON,CMSL7-B-NOPF-ALLNOTRD,CMSH7-B-NOPF-ALLNOTRD,CMUL7-B-NOPF-ALLNOTRD";
947 TObjArray* triggerArray = sTriggerList.Tokenize(",");
948 triggerArray->SetOwner(kTRUE);
950 std::ofstream* out(0x0);
952 if ( TString(csvOutputFile).Length() > 0 )
954 out = new std::ofstream(gSystem->ExpandPathName(csvOutputFile));
955 if (!out || out->bad())
962 TIter nextTrigger(triggerArray);
967 (*out) << "Fill number" << sep;
971 std::map<std::string,float>::const_iterator it;
973 while ( ( trigger = static_cast<TObjString*>(nextTrigger())))
975 lumiPerTrigger[trigger->String().Data()] = 0;
978 for ( it = lumiPerTrigger.begin(); it != lumiPerTrigger.end(); ++it )
980 (*out) << "lumi from " << it->first.c_str() << sep;
983 (*out) << "comments" << sep;
985 for ( it = lumiPerTrigger.begin(); it != lumiPerTrigger.end(); ++it )
987 (*out) << "recorded " << it->first.c_str() << " (integrated)" << sep;
990 (*out) << Form("LHC delivered (%s^-1) per fill ",csUnit)
991 << sep << Form("LHC delivered (%s^-1 integrated)",csUnit) << sep;
992 (*out) << "lumi tot muon" << sep << "eff (%)" << sep;
998 Int_t currentFillNumber(-1);
1002 TString scsUnit(csUnit);
1005 if ( scsUnit=="PB" )
1009 else if (scsUnit=="NB")
1013 else if ( scsUnit=="UB" )
1017 else if ( scsUnit=="MB" )
1023 AliError(Form("Don't know how to deal with cross-section unit=%s",csUnit));
1028 for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
1030 Int_t runNumber = fRunList[i];
1031 Bool_t atLeastOneTriggerFound(kFALSE);
1033 // find out which trigger classes to use
1035 AliTriggerConfiguration* tc = static_cast<AliTriggerConfiguration*>(GetOCDBObject("GRP/CTP/Config",runNumber));
1036 const TObjArray& trClasses = tc->GetClasses();
1040 fillNumber = GetFillNumberFromRunNumber(runNumber);
1041 if ( fillNumber == 0 )
1043 AliError(Form("Got fillNumber = 0 for run %09d !",runNumber));
1046 if ( fillNumber != currentFillNumber )
1048 std::map<std::string,float>::const_iterator it;
1050 for ( it = lumiPerTrigger.begin(); it != lumiPerTrigger.end(); ++it )
1052 lumiPerFillPerTrigger[fillNumber][it->first.c_str()] = 0;
1054 currentFillNumber = fillNumber;
1058 AliTriggerRunScalers* trs = static_cast<AliTriggerRunScalers*>(GetOCDBObject("GRP/CTP/Scalers",runNumber));
1059 const TObjArray* scalers = trs->GetScalersRecords();
1060 const AliTriggerScalersRecord* begin = (AliTriggerScalersRecord*)(scalers->First());
1061 const AliTriggerScalersRecord* end = (AliTriggerScalersRecord*)(scalers->Last());
1063 time_t duration = TMath::Nint((end->GetTimeStamp()->GetBunchCross() - begin->GetTimeStamp()->GetBunchCross())*AliTimeStamp::fNanosecPerBC*1E-9);
1065 nextTrigger.Reset();
1067 TString lumiTriggerClassName(lumiTrigger);
1068 float lumiSigma = lumiCrossSection*csMult; //in csUnit
1069 AliAnalysisTriggerScalerItem* lumiB = GetTriggerScaler(runNumber,"L0B",lumiTriggerClassName.Data());
1073 AliError(Form("Did not find lumiTrigger %s for run %09d",lumiTriggerClassName.Data(),runNumber));
1077 Float_t pacCorrection(1.0);
1079 while ( ( trigger = static_cast<TObjString*>(nextTrigger()) ) )
1081 TString muTriggerClassName(trigger->String());
1084 if ( !trClasses.FindObject(muTriggerClassName.Data() ) )
1089 if ( muTriggerClassName.Contains("CMUL8") ) ++n;
1090 if ( muTriggerClassName.Contains("CMUL7") ) ++n;
1094 AliError(Form("More than 1 relevant trigger class found for run %09d ! Check that !",runNumber));
1099 AliAnalysisTriggerScalerItem* muonA = GetTriggerScaler(runNumber,"L0A",muTriggerClassName.Data());
1100 AliAnalysisTriggerScalerItem* muonB = GetTriggerScaler(runNumber,"L0B",muTriggerClassName.Data());
1102 if (!muonA || !muonB) continue;
1104 atLeastOneTriggerFound = kTRUE;
1108 if ( muonB->ValueCorrectedForDownscale() > 0 )
1110 ratio = muonA->ValueCorrectedForDownscale()/muonB->ValueCorrectedForDownscale();
1113 ratio *= lumiB->ValueCorrectedForDownscale()/lumiSigma;
1115 if ( muTriggerClassName.BeginsWith("CMUL"))
1120 if ( muTriggerClassName.Contains("CMSH") )
1122 pacCorrection = GetPauseAndConfigCorrection(runNumber,muTriggerClassName.Data());
1125 lumiPerTrigger[muTriggerClassName.Data()] += ratio;
1126 durationPerTrigger[muTriggerClassName.Data()] += duration;
1127 lumiPerFillPerTrigger[currentFillNumber][muTriggerClassName.Data()] += ratio;
1131 lumiPerTrigger[lumiTriggerClassName.Data()] += lumiB->ValueCorrectedForDownscale()/lumiSigma;
1132 durationPerTrigger[lumiTriggerClassName.Data()] += duration;
1133 lumiPerFillPerTrigger[currentFillNumber][lumiTriggerClassName.Data()] += lumiB->ValueCorrectedForDownscale()/lumiSigma;
1135 TString lumiPACCorrected(Form("%s(-PAC)",lumiTriggerClassName.Data()));
1137 lumiPerTrigger[lumiPACCorrected.Data()] += pacCorrection*lumiB->ValueCorrectedForDownscale()/lumiSigma;
1138 durationPerTrigger[lumiPACCorrected.Data()] += duration;
1139 lumiPerFillPerTrigger[currentFillNumber][lumiPACCorrected.Data()] += pacCorrection*lumiB->ValueCorrectedForDownscale()/lumiSigma;
1141 if (!atLeastOneTriggerFound && sTriggerList.Contains("CMUL") )
1143 AliError(Form("Found no relevant trigger for run %09d",runNumber));
1147 AliInfo(Form("Integrated lumi %7.4f %s^-1",intLumi,csUnit));
1149 std::map<std::string,float>::const_iterator it;
1151 for ( it = lumiPerTrigger.begin(); it != lumiPerTrigger.end(); ++it )
1153 AliInfo(Form("Trigger %30s Lumi %10.4f %s^-1 duration %10ld s",it->first.c_str(),it->second,csUnit,durationPerTrigger[it->first]));
1158 std::map<int, std::map<std::string, float> >::const_iterator fit;
1160 lumiPerTrigger.clear();
1162 for ( fit = lumiPerFillPerTrigger.begin(); fit != lumiPerFillPerTrigger.end(); ++fit )
1164 int fill = fit->first;
1165 std::map<std::string,float>::const_iterator tit;
1166 (*out) << fill << sep;
1168 for ( tit = fit->second.begin(); tit != fit->second.end(); ++tit )
1170 (*out) << Form("%e",tit->second) << sep;
1173 (*out) << sep; // comment (empty)
1175 for ( tit = fit->second.begin(); tit != fit->second.end(); ++tit )
1177 lumiPerTrigger[tit->first] += tit->second;
1179 (*out) << Form("%e",lumiPerTrigger[tit->first]) << sep;
1181 (*out) << sep << "0" << sep << "0" << sep << "0" << sep << "0" << sep << std::endl; // LHC per fill, LHC integrated, lumi tot muon , efficiency
1185 delete triggerArray;
1188 //______________________________________________________________________________
1189 TGraph* AliAnalysisTriggerScalers::MakeGraph(const std::vector<int>& vx,
1190 const std::vector<int>& vex,
1191 const std::vector<double>& vy,
1192 const std::vector<double>& vey)
1194 /// Build a graph from a set of stl vectors
1196 if ( ! ( vx.size() == vex.size() && vx.size() == vy.size() && vx.size() == vey.size() ) )
1198 std::cerr << "incompatible sizes" << std::endl;
1202 double* x = new double[vx.size()];
1203 double* ex = new double[vx.size()];
1204 double* y = new double[vx.size()];
1205 double* ey = new double[vx.size()];
1207 for ( size_t i = 0; i < vx.size(); ++i )
1215 TGraph* g = new TGraphErrors(vx.size(),x,y,ex,ey);
1227 //______________________________________________________________________________
1228 Double_t AliAnalysisTriggerScalers::Mu(Double_t L0B, Double_t Nb)
1230 /// L0B = trigger rate before any veto
1231 /// Nb = number of crossing bunches
1233 Double_t p0 = 1-L0B/(Nb*11245.0); // proba to get *no* collision per bunch crossing
1235 return -TMath::Log(p0);
1238 //______________________________________________________________________________
1239 Int_t AliAnalysisTriggerScalers::NumberOfInteractingBunches(const AliLHCData& lhc, Int_t runNumber) const
1241 /// Extract the number of colliding bunches from the LHC data
1243 Int_t numberOfInteractingBunches(0);
1244 Int_t numberOfInteractingBunchesMeasured(0);
1250 AliLHCDipValI* val = lhc.GetBunchConfigDeclared(beam1,0);
1252 for ( Int_t i = 0; i < val->GetSizeTotal(); ++i )
1254 if ( val->GetValue(i) < 0 ) ++numberOfInteractingBunches;
1257 AliLHCDipValI* valm = lhc.GetBunchConfigMeasured(beam1,0);
1259 for ( Int_t i = 0; i < valm->GetSizeTotal(); ++i )
1261 if ( valm->GetValue(i) < 0 ) ++numberOfInteractingBunchesMeasured;
1264 valm = lhc.GetBunchConfigMeasured(beam2,0);
1266 for ( Int_t i = 0; i < valm->GetSizeTotal(); ++i )
1268 if ( valm->GetValue(i) <= 0 ) ++nIBM2;
1271 if (!numberOfInteractingBunches)
1276 if ( numberOfInteractingBunches != numberOfInteractingBunchesMeasured ||
1277 numberOfInteractingBunches != nIBM2 )
1279 Int_t delta = TMath::Max(numberOfInteractingBunches - numberOfInteractingBunchesMeasured,numberOfInteractingBunches-nIBM2);
1281 if ( 1.0*TMath::Abs(delta)/numberOfInteractingBunches > 0.05 ) // more than 5% difference
1283 AliWarning(Form("Got some different number of interacting bunches for fill %d run %d ! NumberOfInteractingBunches=%3d NumberOfInteractingBunchesMeasured=%3d NIBM2=%3d. Will use %d",
1284 lhc.GetFillNumber(),runNumber,
1285 numberOfInteractingBunches,numberOfInteractingBunchesMeasured,nIBM2,numberOfInteractingBunches));
1289 AliDebug(1,Form("Got some different number of interacting bunches for fill %d run %d ! NumberOfInteractingBunches=%3d NumberOfInteractingBunchesMeasured=%3d NIBM2=%3d. Will use %d",
1290 lhc.GetFillNumber(),runNumber,
1291 numberOfInteractingBunches,numberOfInteractingBunchesMeasured,nIBM2,numberOfInteractingBunches));
1296 return numberOfInteractingBunches;
1300 //______________________________________________________________________________
1301 TGraph* AliAnalysisTriggerScalers::PlotTrigger(const char* triggerClassName,
1304 // Plot one of the scalers (L0A,L0B,L0AOVERB,etc...) of a given triggerClass
1305 // Get one value per run.
1307 std::vector<Float_t> x;
1308 std::vector<Float_t> y;
1313 for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
1315 Int_t runNumber = fRunList[i];
1317 AliAnalysisTriggerScalerItem* s = GetTriggerScaler(runNumber,what,triggerClassName);
1321 x.push_back(runNumber);
1323 Double_t value = s->ValueCorrectedForDownscale();
1325 if ( TString(what).Contains("RATE") )
1336 if ( fRunList.size() ) {
1337 mean /= fRunList.size();
1340 AliInfo(Form("Integral %e mean %e",integral,mean));
1342 return new TGraph(x.size(),&x[0],&y[0]);
1345 //______________________________________________________________________________
1346 TGraph* AliAnalysisTriggerScalers::PlotTriggerEvolution(const char* triggerClassName,
1352 /// Make a graph of "what" for a given trigger class.
1356 /// - L0B (level 0 before any veto applied)
1357 /// - L0A (level 0 after all vetoes)
1358 /// - L0AOVERB (L0A/L0B)
1359 /// - mu ( = -TMath::Log( 1 - P(0) ) where P(0) is the proba to have zero collisions in a bunch crossing)
1360 /// - pileupfactor = mu/(1-exp(-mu)) : the factor to apply to correct the cint1b count rate
1361 /// - vsnb = L0B/(NumberOfInteractingBunches*11245)
1362 /// - NBCX = NumberOfInteractingBunches
1364 TString swhat(what);
1367 if ( swhat.Contains(";"))
1369 swhat.ReplaceAll(";",",");
1370 AliWarningClass("; is not a valid separator, replaced it with ,");
1376 TObjArray* a = swhat.Tokenize(",");
1377 if (a->GetEntries()>1)
1381 TObjString* str(0x0);
1382 Double_t ymin(TMath::Limits<Double_t>::Max());
1384 Double_t xmin(TMath::Limits<Double_t>::Max());
1388 while ( ( str = static_cast<TObjString*>(next())) )
1390 g = PlotTriggerEvolution(triggerClassName,str->String().Data(),false,0x0,removeZeros);
1392 for ( Int_t i = 0; i < g->GetN(); ++i )
1394 ymin = TMath::Min(ymin,g->GetY()[i]);
1395 ymax = TMath::Max(ymax,g->GetY()[i]);
1397 xmin = TMath::Min(xmin,g->GetX()[0]);
1398 xmax = TMath::Max(xmax,g->GetX()[g->GetN()-1]);
1401 gStyle->SetOptTitle(0);
1403 AliInfoClass(Form("x %e ; %e y %e ; %e",xmin,xmax,ymin,ymax));
1404 TH2* h = new TH2F("h",triggerClassName,100,xmin,xmax,100,ymin,ymax);
1405 h->SetStats(kFALSE);
1406 h->GetXaxis()->SetTimeDisplay(1);
1407 h->GetXaxis()->SetTimeFormat("%d/%m %H:%M");
1408 h->GetXaxis()->SetTimeOffset(0,"gmt");
1409 h->GetXaxis()->SetNdivisions(505);
1412 TIter nextGraph(&graphs);
1414 while ( ( g = static_cast<TGraph*>(nextGraph())) )
1416 AliInfoClass(g->GetTitle());
1418 g->SetLineColor(color);
1419 g->SetMarkerColor(color);
1420 g->SetMarkerStyle(marker);
1429 std::vector<int> vx;
1430 std::vector<int> vex;
1431 std::vector<double> vy;
1432 std::vector<double> vey;
1438 for ( std::vector<int>::size_type iRun = 0; iRun < fRunList.size(); ++iRun )
1440 Int_t runNumber = fRunList[iRun];
1442 AliTriggerConfiguration* tc(0x0);
1443 AliTriggerRunScalers* trs(0x0);
1444 AliLHCData* lhc(0x0);
1446 GetCTPObjects(runNumber,tc,trs,lhc);
1448 Int_t numberOfInteractingBunches = NumberOfInteractingBunches(*lhc,runNumber);
1450 const TObjArray* scalers = trs->GetScalersRecords();
1452 const TObjArray& trClasses = tc->GetClasses();
1453 TIter next(&trClasses);
1454 AliTriggerClass* triggerClass;
1456 while ( ( triggerClass = static_cast<AliTriggerClass*>(next()) ) )
1458 UChar_t index = GetIndex(triggerClass->GetMask());
1460 if ( !TString(triggerClass->GetName()).Contains(triggerClassName) ) continue;
1462 TIter nextScaler(scalers);
1463 AliTriggerScalersRecord* record;
1471 Bool_t first(kTRUE);
1473 while ( ( record = static_cast<AliTriggerScalersRecord*>(nextScaler()) ) )
1475 const AliTriggerScalers* scaler = record->GetTriggerScalersForClass(index);
1478 const AliTimeStamp* ats = record->GetTimeStamp();
1480 UInt_t seconds = ats->GetSeconds();// - TTimeStamp::GetZoneOffset();
1482 TTimeStamp ts(seconds,ats->GetMicroSecs());
1484 UInt_t l0b = scaler->GetLOCB() - refl0b;
1485 UInt_t l0a = scaler->GetLOCA() - refl0a;
1486 UInt_t l1b = scaler->GetL1CB() - refl1b;
1487 UInt_t l1a = scaler->GetL1CA() - refl1a;
1488 UInt_t l2b = scaler->GetL2CB() - refl2b;
1489 UInt_t l2a = scaler->GetL2CA() - refl2a;
1490 UInt_t timelapse = seconds - reft;
1492 if ( removeZeros && ( l0b <= 2 || ( l0a <= 2 && l0a != 0 ) || timelapse <= 9 ) )
1494 AliInfo(Form("Skipping point for %s l0b %d l0a %d timelapse %d ts %s",
1495 triggerClassName,l0b,l0a,timelapse,ts.AsString()));
1500 refl0b = scaler->GetLOCB();
1501 refl0a = scaler->GetLOCA();
1502 refl1b = scaler->GetL1CB();
1503 refl1a = scaler->GetL1CA();
1504 refl2b = scaler->GetL2CB();
1505 refl2a = scaler->GetL2CA();
1516 if (swhat.Contains("L0AOVERB") )
1518 value = l0a*1.0/l0b;
1521 error = value*TMath::Sqrt(1.0/l0b+1.0/l0a);
1528 else if ( swhat.Contains("L0B") )
1532 else if (swhat.Contains("L0A") )
1536 else if ( swhat.Contains("L1B") )
1540 else if (swhat.Contains("L1A") )
1544 else if ( swhat.Contains("L2B") )
1548 else if (swhat.Contains("L2A") )
1552 else if ( swhat.Contains("MU") )
1554 if (timelapse==0) continue;
1555 value = Mu(l0b/timelapse,numberOfInteractingBunches);
1556 error = 0.0; // FIXME
1558 else if ( swhat.Contains("PILEUPFACTOR") )
1560 Double_t mu = Mu(l0b/timelapse,numberOfInteractingBunches);
1561 value = mu/(1-TMath::Exp(-mu));
1562 error = -1.0; // FIXME
1564 else if ( swhat.Contains("VSNB") )
1566 value = l0b/(11245.0*numberOfInteractingBunches);
1567 error = -1.0; // FIXME
1569 else if ( swhat.Contains("NBCX"))
1571 value = numberOfInteractingBunches;
1572 error = 1.0; // FIXME
1577 AliInfo(Form("Unknown what %s",what));
1580 if ( ! swhat.Contains("OVER") && ! swhat.Contains("RATIO") &&
1581 ! swhat.Contains("MU") && ! swhat.Contains("PILEUPFACTOR") &&
1582 ! swhat.Contains("RAW") & ! swhat.Contains("NBCX") )
1587 if ( !TMath::Finite(value) ) continue;
1589 if ( !swhat.Contains("L0AOVERB") && error > 0 )
1591 error = value >0 ? 1.0/TMath::Sqrt(value) : 0.0;
1594 if (removeZeros && TMath::Abs(value) < 1E-6 )
1604 vx.push_back(seconds);
1607 vy.push_back(value);
1609 vey.push_back(error);
1617 if (mean && nvalues)
1622 if ( vx.empty() ) return 0;
1624 TGraph* g = MakeGraph(vx,vex,vy,vey);
1625 TString title(Form("TriggerEvolution%s-%s",triggerClassName,swhat.Data()));
1627 g->SetName(title.Data());
1628 g->SetTitle(title.Data());
1630 g->GetYaxis()->SetTitle(title.Data());
1634 TCanvas* c = NewCanvas(g->GetName());
1635 g->SetLineColor(color);
1636 g->SetMarkerColor(color);
1637 g->SetMarkerStyle(marker);
1639 c->SaveAs(Form("%s.png",title.Data()));
1646 //______________________________________________________________________________
1647 TGraph* AliAnalysisTriggerScalers::PlotTriggerRatio(const char* triggerClassName1,
1649 const char* triggerClassName2,
1652 // Plot the ratio of two scalers.
1653 // Get one value per run.
1655 std::vector<Float_t> x;
1656 std::vector<Float_t> y;
1658 for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
1660 Int_t runNumber = fRunList[i];
1662 AliAnalysisTriggerScalerItem* s1 = GetTriggerScaler(runNumber,what1,triggerClassName1);
1663 AliAnalysisTriggerScalerItem* s2 = GetTriggerScaler(runNumber,what2,triggerClassName2);
1665 if (!s1 || !s2) continue;
1667 x.push_back(runNumber);
1671 if ( s2->ValueCorrectedForDownscale() > 0 )
1673 ratio = s1->ValueCorrectedForDownscale()/s2->ValueCorrectedForDownscale();
1678 AliDebug(1,Form("RUN %09d %20s (%s) %12llu (%5d) %20s (%s) %12llu (%5d) R %7.2f",
1680 triggerClassName1,what1,
1681 s1->Value(),s1->DownscalingFactor(),
1682 triggerClassName2,what2,
1683 s2->Value(),s2->DownscalingFactor(),
1687 return new TGraph(x.size(),&x[0],&y[0]);
1690 //______________________________________________________________________________
1691 TGraph* AliAnalysisTriggerScalers::PlotTriggerRatioEvolution(const char* triggerClassName1,
1693 const char* triggerClassName2,
1696 /// Plots the evolution of one trigger ratio
1698 Bool_t draw(kFALSE);
1699 Bool_t removeZeros(kFALSE);
1701 TGraph* g1 = PlotTriggerEvolution(triggerClassName1,what1,draw,0x0,removeZeros);
1702 TGraph* g2 = PlotTriggerEvolution(triggerClassName2,what2,draw,0x0,removeZeros);
1704 if (!g1 || !g2) return 0x0;
1706 if ( g1->GetN() != g2->GetN() )
1708 AliError("Oups. Did not get the same number of points for the 2 graphs ?!");
1714 Double_t x1err,x2err;
1715 Double_t y1err,y2err;
1717 TGraphErrors* g = new TGraphErrors(g1->GetN());
1720 for ( Int_t i = 0; i < g1->GetN(); ++i )
1722 g1->GetPoint(i,x1,y1);
1723 g2->GetPoint(i,x2,y2);
1725 x1err = g1->GetErrorX(i);
1726 x2err = g2->GetErrorX(i);
1728 y1err = g1->GetErrorY(i);
1729 y2err = g2->GetErrorY(i);
1733 AliError(Form("Points at %d don't have the same x ! ",i));
1739 if ( TMath::Abs(y2) > 1E-12 )
1746 if ( TMath::Abs(x1) > 1E-12)
1748 yerr += TMath::Sqrt(x1err*x1err/(x1*x1));
1751 if ( TMath::Abs(x2) > 1E-12)
1753 yerr += TMath::Sqrt(x2err*x2err/(x2*x2));
1758 g->SetPoint(j,x1,y);
1759 g->SetPointError(j,x1err,yerr);
1772 //______________________________________________________________________________
1773 void AliAnalysisTriggerScalers::Print(Option_t* /* opt */) const
1775 /// print our runlist
1776 AliAnalysisTriggerScalers::PrintIntegers(fRunList,',');
1779 //______________________________________________________________________________
1780 void AliAnalysisTriggerScalers::PrintIntegers(const std::vector<int>& integers,
1784 /// print a list of integers
1785 for ( std::vector<int>::size_type i = 0; i < integers.size(); ++i )
1787 out << integers[i] << sep;
1792 //______________________________________________________________________________
1793 void AliAnalysisTriggerScalers::ReadIntegers(const char* filename,
1794 std::vector<int>& integers,
1797 /// Read integers from filename, where integers are either
1798 /// separated by "," or by return carriage
1800 if ( gSystem->AccessPathName(filename)==kTRUE )
1804 std::ifstream in(gSystem->ExpandPathName(filename));
1807 std::set<int> runset;
1811 for ( std::vector<int>::size_type j = 0; j < integers.size(); ++ j )
1813 runset.insert(integers[j]);
1819 in.getline(line,10000,'\n');
1821 TString sline(line);
1823 if (sline.Contains(","))
1825 TObjArray* a = sline.Tokenize(",");
1828 while ( ( s = static_cast<TObjString*>(next()) ) )
1830 runset.insert(s->String().Atoi());
1836 runset.insert(sline.Atoi());
1846 for ( std::set<int>::const_iterator it = runset.begin(); it != runset.end(); ++it )
1848 integers.push_back((*it));
1851 std::sort(integers.begin(),integers.end());
1855 //______________________________________________________________________________
1856 void AliAnalysisTriggerScalers::SetRunList(Int_t runNumber)
1858 // Make the runlist be a single run
1860 fRunList.push_back(runNumber);
1863 //______________________________________________________________________________
1864 void AliAnalysisTriggerScalers::SetRunList(const std::vector<int>& runs)
1866 // Make the runlist be a single run
1870 //______________________________________________________________________________
1871 void AliAnalysisTriggerScalers::SetRunList(const std::set<int>& runs)
1873 // Make the runlist be a single run
1875 for ( std::set<int>::const_iterator it = runs.begin(); it != runs.end(); ++it )
1877 fRunList.push_back(*it);
1881 //______________________________________________________________________________
1882 void AliAnalysisTriggerScalers::SetRunList(const char* runlist)
1884 // Read the runlist from an ASCII file or a comma separated list
1885 // or a space separated list
1889 if ( TString(runlist).Contains(",") || TString(runlist).Contains(" ") )
1891 TObjArray* runs = 0x0;
1892 if ( TString(runlist).Contains(",") )
1894 runs = TString(runlist).Tokenize(",");
1898 runs = TString(runlist).Tokenize(" ");
1902 std::set<int> runset;
1904 while ( ( s = static_cast<TObjString*>(next()) ) )
1906 runset.insert(s->String().Atoi());
1909 for ( std::set<int>::const_iterator it = runset.begin(); it != runset.end(); ++it )
1911 fRunList.push_back((*it));
1914 std::sort(fRunList.begin(),fRunList.end());
1920 ReadIntegers(runlist,fRunList);
1921 if (fRunList.empty())
1923 if ( TString(runlist).IsDigit() )
1925 SetRunList(TString(runlist).Atoi());
1929 AliError("Could not set run list !");
1935 //______________________________________________________________________________
1936 //______________________________________________________________________________
1937 //______________________________________________________________________________
1938 //______________________________________________________________________________
1939 //______________________________________________________________________________
1940 //______________________________________________________________________________
1941 //______________________________________________________________________________
1943 ClassImp(AliAnalysisTriggerScalerItem)
1945 //______________________________________________________________________________
1947 AliAnalysisTriggerScalerItem::BCMaskName() const
1949 // return bc mask name
1950 if ( BCMask() ) return BCMask()->GetName(); else return "";
1953 //______________________________________________________________________________
1954 Int_t AliAnalysisTriggerScalerItem::Compare(const TObject* obj) const
1956 /// compare two scaler items (by means of their run number only)
1957 const AliAnalysisTriggerScalerItem* s = static_cast<const AliAnalysisTriggerScalerItem*>(obj);
1959 if ( s->RunNumber() < RunNumber() )
1963 else if ( s->RunNumber() > RunNumber() )
1970 //______________________________________________________________________________
1971 void AliAnalysisTriggerScalerItem::Print(Option_t* opt) const
1976 if ( RunNumber() > 0 )
1980 if ( sopt.Contains("FULL") )
1985 std::cout << Form("RUN %6d CLASS %24s (%5s %4d BCs) SCALER %s %12llu DS %6d DURATION %ld",
1986 RunNumber(),TriggerClassName(),
1988 BCMask() ? BCMask()->GetNUnmaskedBCs() : 0,
1990 Value(),DownscalingFactor(),Duration()) << std::endl;
1995 std::cout << Form("CLASS %24s SCALER %20llu NRUNS %d",
1996 TriggerClassName(),Value(),NofRuns()) << std::endl;