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"
102 using std::make_pair;
104 ClassImp(AliAnalysisTriggerScalers)
108 //______________________________________________________________________________
109 UChar_t GetIndex(ULong64_t mask)
111 for ( Int_t i = 0; i < 64; ++i )
113 if ( mask & ( 1ull << i ) ) return i+1;
118 //______________________________________________________________________________
119 TCanvas* NewCanvas(const char* name)
121 TCanvas* c = new TCanvas(name,name);
122 c->SetLeftMargin(0.15);
126 //______________________________________________________________________________
127 void TimeAxis(TGraph* g)
129 // g->GetXaxis()->SetTimeDisplay(1);
130 // g->GetXaxis()->SetTimeFormat("%d/%m %H:%M%F2010-12-31 24:00:00");
131 // g->GetXaxis()->SetNdivisions(505);
132 g->GetYaxis()->SetTitleOffset(1.5);
134 g->GetXaxis()->SetTimeDisplay(1);
135 // g->GetXaxis()->SetTimeFormat("%d/%m %H:%M%F2010-12-31 24:00:00");
136 g->GetXaxis()->SetTimeFormat("%d/%m %H:%M");
137 g->GetXaxis()->SetTimeOffset(0,"gmt");
138 g->GetXaxis()->SetNdivisions(505);
142 //______________________________________________________________________________
143 Bool_t IsMainSatelliteCollision(const char* triggerClassName)
145 TString tcn(triggerClassName);
148 return (tcn.Contains("-S-") || tcn.Contains("-SC-") || tcn.Contains("-SA-"));
152 //_____________________________________________________________________________
153 AliAnalysisTriggerScalers::AliAnalysisTriggerScalers(Int_t runNumber, const char* ocdbPath) : fRunList(), fOCDBPath(ocdbPath), fShouldCorrectForPileUp(kTRUE), fCrossSectionUnit("UB")
155 // ctor using a single run number
156 SetRunList(runNumber);
159 //_____________________________________________________________________________
160 AliAnalysisTriggerScalers::AliAnalysisTriggerScalers(const std::vector<int>& runs, const char* ocdbPath) : fRunList(), fOCDBPath(ocdbPath), fShouldCorrectForPileUp(kTRUE), fCrossSectionUnit("UB")
162 // ctor using a vector of run numbers
166 //_____________________________________________________________________________
167 AliAnalysisTriggerScalers::AliAnalysisTriggerScalers(const std::set<int>& runs, const char* ocdbPath) : fRunList(), fOCDBPath(ocdbPath), fShouldCorrectForPileUp(kTRUE), fCrossSectionUnit("UB")
169 // ctor using a set of run numbers
173 //_____________________________________________________________________________
174 AliAnalysisTriggerScalers::AliAnalysisTriggerScalers(const char* runlist, const char* ocdbPath) : fRunList(), fOCDBPath(ocdbPath), fShouldCorrectForPileUp(kTRUE), fCrossSectionUnit("UB")
176 // ctor from a run list (txt file)
180 //_____________________________________________________________________________
181 AliAnalysisTriggerScalers::~AliAnalysisTriggerScalers()
186 //______________________________________________________________________________
187 Bool_t AliAnalysisTriggerScalers::CheckRecord(const AliTriggerScalersRecord& record,
191 UInt_t timelapse) const
193 /// Check if a trigger scaler record should be kept for our
194 /// luminosity computations
196 const AliTriggerScalers* scaler = record.GetTriggerScalersForClass(index);
198 UInt_t a = scaler->GetLOCA() - refa;
199 UInt_t b = scaler->GetLOCB() - refb;
203 if ( b <= 2 || ( a <= 2 ) || timelapse <= 9 ) // empirical cuts
212 //______________________________________________________________________________
213 void AliAnalysisTriggerScalers::DrawFill(Int_t run1, Int_t run2, double ymin, double ymax, const char* label, Int_t color)
215 // Draw a yellow box to indicate a run range
217 AliDebugClass(1,Form("RUN1 %09d RUN2 %09d YMIN %e YMAX %e %s",
218 run1,run2,ymin,ymax,label));
219 TBox* b = new TBox(run1*1.0,ymin,run2*1.0,ymax);
220 b->SetFillColor(color);
222 TText* text = new TText((run1+run2)/2.0,ymin + (ymax-ymin)*0.85,label);
223 text->SetTextSize(0.025);
224 text->SetTextFont(42);
225 text->SetTextAlign(23);
226 text->SetTextAngle(45);
230 //_____________________________________________________________________________
231 void AliAnalysisTriggerScalers::DrawFills(Double_t ymin, Double_t ymax, Int_t color)
233 /// Draw the fill ranges corresponding to the list of runs
234 /// Note that this method will scan the OCDB to get the run -> fill number relationship,
235 /// so it's better in this case to use a local copy of the OCDB if you have one. Otherwise
238 std::map<int, std::pair<int,int> > fills;
240 GetFillBoundaries(fills);
242 std::map<int, std::pair<int,int> >::const_iterator it;
244 for ( it = fills.begin(); it != fills.end(); ++it )
246 const std::pair<int,int>& p = it->second;
248 fillnumber.Form("%d",it->first);
249 DrawFill(p.first,p.second,ymin,ymax,fillnumber.Data(),color);
253 //_____________________________________________________________________________
254 void AliAnalysisTriggerScalers::DrawPeriods(Double_t ymin, Double_t ymax, Int_t color)
256 /// Draw the period ranges corresponding to the list of runs
257 /// Note that this method will scan the OCDB to get the run -> fill number relationship,
258 /// so it's better in this case to use a local copy of the OCDB if you have one. Otherwise
261 std::map<std::string, std::pair<int,int> > periods;
263 GetLHCPeriodBoundaries(periods);
265 for ( std::map<std::string, std::pair<int,int> >::const_iterator it = periods.begin(); it != periods.end(); ++it )
267 const std::pair<int,int>& p = it->second;
268 DrawFill(p.first,p.second,ymin,ymax,it->first.c_str(),color);
273 //_____________________________________________________________________________
274 void AliAnalysisTriggerScalers::GetCTPObjects(Int_t runNumber,
275 AliTriggerConfiguration*& tc,
276 AliTriggerRunScalers*& trs,
277 AliLHCData*& lhc) const
279 /// For a given run, get the CTP objects we need to do our work
281 /// FIXME: this is the method that should probably be virtualized so
282 /// we can get those objects from either OCDB or run based container in AODs/ESDs
284 tc = static_cast<AliTriggerConfiguration*>(GetOCDBObject("GRP/CTP/Config",runNumber));
285 trs = static_cast<AliTriggerRunScalers*>(GetOCDBObject("GRP/CTP/Scalers",runNumber));
286 lhc = static_cast<AliLHCData*>(GetOCDBObject("GRP/GRP/LHCData",runNumber));
289 //_____________________________________________________________________________
290 void AliAnalysisTriggerScalers::GetFillBoundaries(std::map<int, std::pair<int,int> >& fills)
292 /// Given a list of runs get the fills they are in
296 for ( std::vector<int>::const_iterator it = fRunList.begin(); it != fRunList.end(); ++it )
300 int fill = GetFillNumberFromRunNumber(runnumber);
302 if (fills.count(fill))
304 std::pair<int,int>& p = fills[fill];
305 p.first = TMath::Min(runnumber,p.first);
306 p.second = TMath::Max(runnumber,p.second);
310 fills[fill] = make_pair(runnumber,runnumber);
316 //______________________________________________________________________________
317 Int_t AliAnalysisTriggerScalers::GetFillNumberFromRunNumber(Int_t runNumber)
319 /// Get the fill number of a run
321 AliLHCData* lhcData = static_cast<AliLHCData*>(GetOCDBObject("GRP/GRP/LHCData",runNumber));
324 Int_t fillNumber = lhcData->GetFillNumber();
325 if ( fillNumber == 0 && runNumber == 189694)
327 // manual hack because GRP info incorrect for this run ?
336 //______________________________________________________________________________
337 TObject* AliAnalysisTriggerScalers::GetOCDBObject(const char* path, Int_t runNumber) const
339 // Helper method to get an object from the OCDB
341 if ( !AliCDBManager::Instance()->IsDefaultStorageSet() )
343 AliInfo(Form("Setting OCDB default storage to %s",fOCDBPath.Data()));
344 AliCDBManager::Instance()->SetDefaultStorage(fOCDBPath.Data());
347 AliCDBManager::Instance()->SetRun(runNumber);
349 gErrorIgnoreLevel=kWarning; // to avoid all the TAlienFile::Open messages...
351 AliCDBEntry* e = AliCDBManager::Instance()->Get(path);
352 return e ? e->GetObject() : 0x0;
356 //______________________________________________________________________________
357 TString AliAnalysisTriggerScalers::GetLHCPeriodFromRunNumber(Int_t runNumber) const
359 /// Get a list of LHC periods from a list of run numbers
361 AliGRPObject* grp = static_cast<AliGRPObject*>(GetOCDBObject("GRP/GRP/Data",runNumber));
365 return grp->GetLHCPeriod();
368 //_____________________________________________________________________________
370 AliAnalysisTriggerScalers::GetLuminosityTriggerAndCrossSection(Int_t runNumber,
371 TString& lumiTriggerClassName,
372 Double_t& lumiTriggerCrossSection,
373 Double_t& lumiTriggerCrossSectionError) const
375 /// For one given run, get the trigger class name to be used as the luminometer,
376 /// and its corresponding cross-section
378 /// FIXME: all is hard-coded for now... (use OADB for this ?)
380 lumiTriggerClassName="";
381 lumiTriggerCrossSection=0.0;
382 lumiTriggerCrossSectionError=0.0; // FIXME
384 TString period = GetLHCPeriodFromRunNumber(runNumber);
386 if ( period.BeginsWith("LHC11") )
388 AliError("Not implemented yet");
390 else if ( period.BeginsWith("LHC12") )
392 if ( period == "LHC12h" || period == "LHC12i" )
394 // pp 8 TeV main-satellites
395 lumiTriggerClassName = "C0TVX-S-NOPF-ALLNOTRD";
396 lumiTriggerCrossSection = 28.0; // FIXME: not from a vdM yet
400 AliError("Not implemented yet");
403 else if ( period.BeginsWith("LHC13") )
405 if ( period == "LHC13g" )
408 lumiTriggerClassName = "C0TVX-B-NOPF-ALLNOTRD";
409 lumiTriggerCrossSection = 18.0; // FIXME: not from a vdM yet
413 // p-Pb or Pb-p 5.02 TeV
414 lumiTriggerClassName = "C0TVX-B-NOPF-ALLNOTRD";
415 lumiTriggerCrossSection = 0.755*2000; // FIXME: not from a vdM yet
420 AliError("Not implemented yet");
425 if ( fCrossSectionUnit=="PB" )
429 else if (fCrossSectionUnit=="NB")
433 else if ( fCrossSectionUnit=="UB" )
437 else if ( fCrossSectionUnit=="MB" )
442 lumiTriggerCrossSection *= csMult;
443 lumiTriggerCrossSectionError *= csMult;
446 //_____________________________________________________________________________
447 void AliAnalysisTriggerScalers::GetLHCPeriodBoundaries(std::map<std::string, std::pair<int,int> >& periods)
449 /// Given a list of runs get the fills they are in
453 for ( std::vector<int>::const_iterator it = fRunList.begin(); it != fRunList.end(); ++it )
457 std::string period = GetLHCPeriodFromRunNumber(runnumber).Data();
459 if (periods.count(period))
461 std::pair<int,int>& p = periods[period];
462 p.first = TMath::Min(runnumber,p.first);
463 p.second = TMath::Max(runnumber,p.second);
467 periods[period] = make_pair(runnumber,runnumber);
473 //______________________________________________________________________________
474 Int_t AliAnalysisTriggerScalers::GetTriggerInput(Int_t runNumber, const char* inputname)
476 /// Get one trigger input for a given run
477 AliTriggerConfiguration* tc = static_cast<AliTriggerConfiguration*>(GetOCDBObject("GRP/CTP/Config",runNumber));
480 const TObjArray& inputs = tc->GetInputs();
482 AliTriggerInput* ti = static_cast<AliTriggerInput*>(inputs.FindObject(inputname));
486 return ti->GetSignature();
489 //______________________________________________________________________________
491 AliAnalysisTriggerScalers::GetPauseAndConfigCorrection(Int_t runNumber, const char* triggerClassName)
493 /// Tries to estimate the duration of the Pause and Configure phase(s) in a given run
494 /// Probably not really accurate for the moment though.
499 AliTriggerConfiguration* tc = static_cast<AliTriggerConfiguration*>(GetOCDBObject("GRP/CTP/Config",runNumber));
501 AliTriggerRunScalers* trs = static_cast<AliTriggerRunScalers*>(GetOCDBObject("GRP/CTP/Scalers",runNumber));
502 const TObjArray* scalers = trs->GetScalersRecords();
505 AliTriggerScalersRecord* record;
511 Int_t classIndex = tc->GetClassIndexFromName(triggerClassName);
513 while ( ( record = static_cast<AliTriggerScalersRecord*>(next()) ) )
515 const AliTriggerScalers* scaler = record->GetTriggerScalersForClass(classIndex);
517 const AliTimeStamp* ats = record->GetTimeStamp();
519 UInt_t seconds = ats->GetSeconds();// - TTimeStamp::GetZoneOffset();
521 TTimeStamp ts(seconds,ats->GetMicroSecs());
523 UInt_t l0b = scaler->GetLOCB() - refl0b;
524 UInt_t l0a = scaler->GetLOCA() - refl0a;
525 UInt_t timelapse = seconds - reft;
527 if ( l0b <=2 || timelapse < 9 ) continue;
530 refl0b = scaler->GetLOCB();
531 refl0a = scaler->GetLOCA();
547 return total > 0 ? nopac*1.0/total : 0.0;
550 //______________________________________________________________________________
551 void AliAnalysisTriggerScalers::ShowPileUpFactors(const char* triggerClassName, Double_t purity)
553 /// Printout the pile-up factors for the runlist, including some purity if needed
557 AliError(Form("Cannot work with purity=%f for trigger %s. Should be strictly positive",purity,triggerClassName));
561 const std::vector<int>& rl = GetRunList();
562 Double_t value(0.0),error(0.0);
563 std::vector<std::string> lines; // put the data into lines so we can get the OCDB access printout grouped before the important stuff
565 for ( std::vector<int>::const_iterator it = rl.begin(); it != rl.end(); ++it )
567 GetPileUpFactor(*it,triggerClassName,purity,value,error);
568 TString period = GetLHCPeriodFromRunNumber(*it);
569 std::ostringstream str;
570 str << Form("RUN %6d PERIOD %6s PILE-UP CORRECTION FACTOR (mu/(1-exp(-mu)) = %7.4f",*it,period.Data(),value);
571 lines.push_back(str.str());
574 for ( std::vector<std::string>::size_type i = 0; i < lines.size(); ++i )
576 std::cout << lines[i].c_str() << std::endl;
580 //______________________________________________________________________________
581 void AliAnalysisTriggerScalers::GetPileUpFactor(Int_t runNumber, const char* triggerClassName,
583 Double_t& value, Double_t& error)
585 /// Get the mean pile-up correction factor for the given run
591 AliError(Form("Cannot work with purity=%f for trigger %s in run %d. Should be strictly positive",purity,triggerClassName,runNumber));
595 AliTriggerRunScalers* trs = static_cast<AliTriggerRunScalers*>(GetOCDBObject("GRP/CTP/Scalers",runNumber));
596 const TObjArray* scalers = trs->GetScalersRecords();
597 const AliTriggerScalersRecord* begin = (AliTriggerScalersRecord*)(scalers->First());
598 const AliTriggerScalersRecord* end = (AliTriggerScalersRecord*)(scalers->Last());
600 time_t duration = TMath::Nint((end->GetTimeStamp()->GetBunchCross() - begin->GetTimeStamp()->GetBunchCross())*AliTimeStamp::fNanosecPerBC*1E-9);
604 AliError(Form("Got zero duration for run %d",runNumber));
608 AliAnalysisTriggerScalerItem* item = GetTriggerScaler(runNumber,"L0B",triggerClassName);
611 AliError(Form("Could not get L0B for trigger %s in run %d",triggerClassName,runNumber));
615 AliLHCData* lhc = static_cast<AliLHCData*>(GetOCDBObject("GRP/GRP/LHCData",runNumber));
617 Bool_t mainSat = IsMainSatelliteCollision(triggerClassName);
619 Int_t nbcx = NumberOfInteractingBunches(*lhc,runNumber,mainSat);
623 AliError(Form("Cannot work with nbcx=%d for trigger %s in run %d. Should be strictly positive",nbcx,triggerClassName,runNumber));
627 Double_t itemValue = purity*item->Value();
631 AliError(Form("Cannot work with value=%f for trigger %s in run %d. Should be strictly positive",itemValue,triggerClassName,runNumber));
635 Double_t mu = Mu(itemValue/duration,nbcx);
637 value = mu/(1.0-TMath::Exp(-mu));
639 error = 0.0; // FIXME
642 //______________________________________________________________________________
643 AliAnalysisTriggerScalerItem*
644 AliAnalysisTriggerScalers::GetTriggerScaler(Int_t runNumber, const char* level, const char* triggerClassName)
646 // Get the scaler for a given trigger class and a given trigger level.
647 // Returned object must be deleted by the user.
649 AliTriggerConfiguration* tc = static_cast<AliTriggerConfiguration*>(GetOCDBObject("GRP/CTP/Config",runNumber));
650 AliTriggerRunScalers* trs = static_cast<AliTriggerRunScalers*>(GetOCDBObject("GRP/CTP/Scalers",runNumber));
651 AliGRPObject* grp = static_cast<AliGRPObject*>(GetOCDBObject("GRP/GRP/Data",runNumber));
653 if (!tc || !trs || !grp) return 0x0;
655 TString diCurrent(Form("L3:%5.0f;DIP:%5.0f [L3:%d;DIP:%d]",
656 grp->GetL3Current((AliGRPObject::Stats)0),
657 grp->GetDipoleCurrent((AliGRPObject::Stats)0),
658 grp->GetL3Polarity(),
659 grp->GetDipolePolarity()));
661 const TObjArray& trClasses = tc->GetClasses();
663 const TObjArray* scalers = trs->GetScalersRecords();
664 const AliTriggerScalersRecord* begin = (AliTriggerScalersRecord*)(scalers->First());
665 const AliTriggerScalersRecord* end = (AliTriggerScalersRecord*)(scalers->Last());
667 time_t duration = TMath::Nint((end->GetTimeStamp()->GetBunchCross() - begin->GetTimeStamp()->GetBunchCross())*AliTimeStamp::fNanosecPerBC*1E-9);
669 AliTriggerClass* triggerClass = static_cast<AliTriggerClass*>(trClasses.FindObject(triggerClassName));
674 UChar_t index = GetIndex(triggerClass->GetMask());
678 const AliTriggerScalers* sbegin = begin->GetTriggerScalersForClass(index);
679 const AliTriggerScalers* send = end->GetTriggerScalersForClass(index);
681 TString swhat(level);
684 if ( swhat.BeginsWith("L0A") )
686 value = send->GetLOCA() - sbegin->GetLOCA();
688 else if ( swhat.BeginsWith("L0B") )
690 value = send->GetLOCB() - sbegin->GetLOCB();
692 else if ( swhat.BeginsWith("L1A") )
694 value = send->GetL1CA() - sbegin->GetL1CA();
696 else if ( swhat.BeginsWith("L1B") )
698 value = send->GetL1CB() - sbegin->GetL1CB();
700 else if ( swhat.BeginsWith("L2A") )
702 value = send->GetL2CA() - sbegin->GetL2CA();
704 else if ( swhat.BeginsWith("L2B") )
706 value = send->GetL2CB() - sbegin->GetL2CB();
710 // FIXME: get the downscaling factor here for all cases (only BC1 supported here so far)
711 if ( TString(triggerClassName).Contains("_B1") )
713 ds = GetTriggerInput(runNumber,"BC1");
717 swhat.ReplaceAll("RATE","");
719 return new AliAnalysisTriggerScalerItem(runNumber,swhat.Data(),diCurrent.Data(),triggerClass->GetName(),value,triggerClass->GetBCMask(),ds,duration);
722 //______________________________________________________________________________
724 AliAnalysisTriggerScalers::IntegratedLuminosityGraph(const char* triggerClassName, const char* triggerClassNameForPACEstimate)
726 // Get the integrated luminosity graph for one trigger
728 if ( fRunList.empty() )
730 AliError("Cannot work without a run list");
734 Double_t offset(0.0);
736 std::vector<double> vx;
737 std::vector<double> vy;
741 for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
743 TGraph* g = IntegratedLuminosityGraph(fRunList[i],triggerClassName,triggerClassNameForPACEstimate);
747 for ( Int_t j = 0 ; j < g->GetN(); ++j )
751 vy.push_back(y+offset);
759 TGraph* g = new TGraph(vx.size(),&vx[0],&vy[0]);
766 //______________________________________________________________________________
768 AliAnalysisTriggerScalers::IntegratedLuminosityGraph(Int_t runNumber, const char* triggerClassName, const char* triggerClassNameForPACEstimate)
770 // Get the integrated luminosity graph for one trigger for the given run
772 // the triggerClassNameForPACEstimate is only used if triggerClassName is,
773 // for the given run, the same as the luminometer trigger class name,
774 // and should be the trigger class with the higher and non-zero L0A rate (except during PACs)
777 AliTriggerConfiguration* tc(0x0);
778 AliTriggerRunScalers* trs(0x0);
779 AliLHCData* lhc(0x0);
781 GetCTPObjects(runNumber,tc,trs,lhc);
783 if (!tc || !trs || !lhc)
788 const TObjArray* scalerRecords= trs->GetScalersRecords();
789 const TObjArray& triggerClasses = tc->GetClasses();
791 TString lumiTriggerClassName;
792 Double_t lumiTriggerCrossSection(0.0);
793 Double_t lumiTriggerCrossSectionError(0.0);
795 GetLuminosityTriggerAndCrossSection(runNumber,lumiTriggerClassName,lumiTriggerCrossSection,lumiTriggerCrossSectionError);
797 AliTriggerClass* lumiTriggerClass = static_cast<AliTriggerClass*>(triggerClasses.FindObject(lumiTriggerClassName));
799 if (!lumiTriggerClass)
801 AliError(Form("Could not find lumiTriggerClass=%s",lumiTriggerClassName.Data()));
805 if (lumiTriggerCrossSection==0.0)
807 AliError(Form("Could not get cross-section for lumiTriggerClass=%s",lumiTriggerClassName.Data()));
811 AliTriggerClass* triggerClass = static_cast<AliTriggerClass*>(triggerClasses.FindObject(triggerClassName));
815 AliError(Form("Could not find triggerClass=%s",triggerClassName));
819 Int_t nbcx = NumberOfInteractingBunches(*lhc,runNumber);
821 if (nbcx <= 0 && ShouldCorrectForPileUp())
823 AliError("Could not get the number of bunches, so won't be able to correct for pile-up");
826 TIter nextScaler(scalerRecords);
827 AliTriggerScalersRecord* record;
830 UInt_t refl0b[] = { 0, 0 };
831 UInt_t refl0a[] = { 0, 0 };
833 UInt_t l0b[] = { 0, 0 };
834 UInt_t l0a[] = { 0, 0 };
838 UChar_t classIndices[2];
840 // class 0 = class for luminomiter
841 // class 1 = class for livetime estimate
842 // or for PAC estimate if triggerClassNameForPACEstimate
843 // is given and triggerClass is the same as the lumi class
845 classIndices[0] = GetIndex(lumiTriggerClass->GetMask());
846 classIndices[1] = GetIndex(triggerClass->GetMask());
848 Bool_t sameClass = ( classIndices[0] == classIndices[1] );
849 Bool_t pacRemoval(kFALSE);
851 if ( sameClass && strlen(triggerClassNameForPACEstimate) > 0 )
853 AliTriggerClass* triggerClassForPACEstimate = static_cast<AliTriggerClass*>(triggerClasses.FindObject(triggerClassNameForPACEstimate));
854 if (!triggerClassForPACEstimate)
856 AliError(Form("Could not find triggerClassForPACEstimate=%s. Will not correct for PAC durations",triggerClassNameForPACEstimate));
859 classIndices[1] = GetIndex(triggerClassForPACEstimate->GetMask());
860 sameClass = ( classIndices[0] == classIndices[1] );
867 ULong64_t counter(0);
869 std::vector<Double_t> vseconds;
870 std::vector<Double_t> vlivetime;
871 std::vector<Double_t> vlumi;
873 while ( ( record = static_cast<AliTriggerScalersRecord*>(nextScaler()) ) )
875 const AliTimeStamp* ats = record->GetTimeStamp();
877 UInt_t seconds = ats->GetSeconds();// - TTimeStamp::GetZoneOffset();
879 TTimeStamp ts(seconds,ats->GetMicroSecs());
881 UInt_t timelapse = seconds - reft;
883 Bool_t ok = sameClass || CheckRecord(*record,classIndices[1],refl0b[1],refl0a[1],timelapse);
885 if (ok) reft = seconds;
887 for ( Int_t i = 0; i < 2; ++i )
889 const AliTriggerScalers* scaler = record->GetTriggerScalersForClass(classIndices[i]);
893 l0b[i] = scaler->GetLOCB() - refl0b[i];
899 l0a[i] = scaler->GetLOCA() - refl0a[i];
901 refl0b[i] = scaler->GetLOCB();
902 refl0a[i] = scaler->GetLOCA();
911 Double_t factor(1.0);
913 if ( ShouldCorrectForPileUp() && nbcx > 0 && l0b[0] > 0 )
915 Double_t mu = Mu(l0b[0]/timelapse,nbcx);
916 factor = mu/(1-TMath::Exp(-mu));
919 vseconds.push_back(seconds);
923 if ( ok && !sameClass && !pacRemoval )
925 lt = (l0a[1]*1.0)/l0b[1];
928 counter += TMath::Nint(factor*l0b[0]*lt);
930 vlumi.push_back( counter / lumiTriggerCrossSection );
934 if ( vseconds.empty() ) return 0x0;
936 TGraph* g = new TGraph(vseconds.size(),&vseconds[0],&vlumi[0]);
943 //______________________________________________________________________________
944 void AliAnalysisTriggerScalers::IntegratedLuminosity(const char* triggerList,
945 const char* lumiTrigger,
946 Double_t lumiCrossSection,
947 const char* csvOutputFile,
951 // Compute the luminosity for a set of triggers
953 // for T0 based lumi (end of pp 2012), use lumiTrigger = C0TVX-S-NOPF-ALLNOTRD and lumiCrossSection = 28 mb (and csUnit="nb")
955 // for T0 based lumi (pPb 2013), use lumiTrigger = C0TVX-B-NOPF-ALLNOTRD
956 // and lumiCrossSection = 0.755*2000 mb (and csUnit="mub")
958 // for T0 based lumi (pp 2.76 TeV 2013), use lumiTrigger = C0TVX-B-NOPF-ALLNOTRD
959 // and lumiCrossSection = 18 mb (and csUnit="nb")
963 std::map<std::string,float> lumiPerTrigger;
964 std::map<int, std::map<std::string,float> > lumiPerFillPerTrigger;
966 std::map<std::string,time_t> durationPerTrigger;
968 TString sTriggerList(triggerList);
970 if ( sTriggerList.Length()==0 )
972 if ( lumiCrossSection < 30 && lumiCrossSection > 20 )
974 // assuming it's pp 2012
975 sTriggerList = "CMUL8-S-NOPF-ALLNOTRD,CMUL7-S-NOPF-ALLNOTRD,CMUL8-S-NOPF-MUON,CMUL7-S-NOPF-MUON";
977 // for C0MUL-SC-NOPF-MUON must use C0TVX-SC as lumiTrigger (with same cross-section as C0TVX=28mb)
979 else if ( lumiCrossSection < 30 )
981 // assuming it's pp 2013
982 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";
986 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";
990 TObjArray* triggerArray = sTriggerList.Tokenize(",");
991 triggerArray->SetOwner(kTRUE);
993 std::ofstream* out(0x0);
995 if ( TString(csvOutputFile).Length() > 0 )
997 out = new std::ofstream(gSystem->ExpandPathName(csvOutputFile));
998 if (!out || out->bad())
1005 TIter nextTrigger(triggerArray);
1006 TObjString* trigger;
1010 (*out) << "Fill number" << sep;
1012 nextTrigger.Reset();
1014 std::map<std::string,float>::const_iterator it;
1016 while ( ( trigger = static_cast<TObjString*>(nextTrigger())))
1018 lumiPerTrigger[trigger->String().Data()] = 0;
1021 for ( it = lumiPerTrigger.begin(); it != lumiPerTrigger.end(); ++it )
1023 (*out) << "lumi from " << it->first.c_str() << sep;
1026 (*out) << "comments" << sep;
1028 for ( it = lumiPerTrigger.begin(); it != lumiPerTrigger.end(); ++it )
1030 (*out) << "recorded " << it->first.c_str() << " (integrated)" << sep;
1033 (*out) << Form("LHC delivered (%s^-1) per fill ",csUnit)
1034 << sep << Form("LHC delivered (%s^-1 integrated)",csUnit) << sep;
1035 (*out) << "lumi tot muon" << sep << "eff (%)" << sep;
1036 (*out) << std::endl;
1038 nextTrigger.Reset();
1041 Int_t currentFillNumber(-1);
1042 Int_t fillNumber(0);
1045 TString scsUnit(csUnit);
1048 if ( scsUnit=="PB" )
1052 else if (scsUnit=="NB")
1056 else if ( scsUnit=="UB" )
1060 else if ( scsUnit=="MB" )
1066 AliError(Form("Don't know how to deal with cross-section unit=%s",csUnit));
1071 for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
1073 Int_t runNumber = fRunList[i];
1074 Bool_t atLeastOneTriggerFound(kFALSE);
1076 // find out which trigger classes to use
1078 AliTriggerConfiguration* tc = static_cast<AliTriggerConfiguration*>(GetOCDBObject("GRP/CTP/Config",runNumber));
1079 const TObjArray& trClasses = tc->GetClasses();
1083 fillNumber = GetFillNumberFromRunNumber(runNumber);
1084 if ( fillNumber == 0 )
1086 AliError(Form("Got fillNumber = 0 for run %09d !",runNumber));
1089 if ( fillNumber != currentFillNumber )
1091 std::map<std::string,float>::const_iterator it;
1093 for ( it = lumiPerTrigger.begin(); it != lumiPerTrigger.end(); ++it )
1095 lumiPerFillPerTrigger[fillNumber][it->first.c_str()] = 0;
1097 currentFillNumber = fillNumber;
1101 AliTriggerRunScalers* trs = static_cast<AliTriggerRunScalers*>(GetOCDBObject("GRP/CTP/Scalers",runNumber));
1102 const TObjArray* scalers = trs->GetScalersRecords();
1103 const AliTriggerScalersRecord* begin = (AliTriggerScalersRecord*)(scalers->First());
1104 const AliTriggerScalersRecord* end = (AliTriggerScalersRecord*)(scalers->Last());
1106 time_t duration = TMath::Nint((end->GetTimeStamp()->GetBunchCross() - begin->GetTimeStamp()->GetBunchCross())*AliTimeStamp::fNanosecPerBC*1E-9);
1108 nextTrigger.Reset();
1110 TString lumiTriggerClassName(lumiTrigger);
1111 float lumiSigma = lumiCrossSection*csMult; //in csUnit
1112 AliAnalysisTriggerScalerItem* lumiB = GetTriggerScaler(runNumber,"L0B",lumiTriggerClassName.Data());
1116 AliError(Form("Did not find lumiTrigger %s for run %09d",lumiTriggerClassName.Data(),runNumber));
1120 Float_t pacCorrection(1.0);
1122 while ( ( trigger = static_cast<TObjString*>(nextTrigger()) ) )
1124 TString muTriggerClassName(trigger->String());
1127 if ( !trClasses.FindObject(muTriggerClassName.Data() ) )
1132 if ( muTriggerClassName.Contains("CMUL8") ) ++n;
1133 if ( muTriggerClassName.Contains("CMUL7") ) ++n;
1137 AliError(Form("More than 1 relevant trigger class found for run %09d ! Check that !",runNumber));
1142 AliAnalysisTriggerScalerItem* muonA = GetTriggerScaler(runNumber,"L0A",muTriggerClassName.Data());
1143 AliAnalysisTriggerScalerItem* muonB = GetTriggerScaler(runNumber,"L0B",muTriggerClassName.Data());
1145 if (!muonA || !muonB) continue;
1147 atLeastOneTriggerFound = kTRUE;
1151 if ( muonB->ValueCorrectedForDownscale() > 0 )
1153 ratio = muonA->ValueCorrectedForDownscale()/muonB->ValueCorrectedForDownscale();
1156 ratio *= lumiB->ValueCorrectedForDownscale()/lumiSigma;
1158 if ( muTriggerClassName.BeginsWith("CMUL"))
1163 if ( muTriggerClassName.Contains("CMSH") )
1165 pacCorrection = GetPauseAndConfigCorrection(runNumber,muTriggerClassName.Data());
1168 lumiPerTrigger[muTriggerClassName.Data()] += ratio;
1169 durationPerTrigger[muTriggerClassName.Data()] += duration;
1170 lumiPerFillPerTrigger[currentFillNumber][muTriggerClassName.Data()] += ratio;
1174 lumiPerTrigger[lumiTriggerClassName.Data()] += lumiB->ValueCorrectedForDownscale()/lumiSigma;
1175 durationPerTrigger[lumiTriggerClassName.Data()] += duration;
1176 lumiPerFillPerTrigger[currentFillNumber][lumiTriggerClassName.Data()] += lumiB->ValueCorrectedForDownscale()/lumiSigma;
1178 TString lumiPACCorrected(Form("%s(-PAC)",lumiTriggerClassName.Data()));
1180 lumiPerTrigger[lumiPACCorrected.Data()] += pacCorrection*lumiB->ValueCorrectedForDownscale()/lumiSigma;
1181 durationPerTrigger[lumiPACCorrected.Data()] += duration;
1182 lumiPerFillPerTrigger[currentFillNumber][lumiPACCorrected.Data()] += pacCorrection*lumiB->ValueCorrectedForDownscale()/lumiSigma;
1184 if (!atLeastOneTriggerFound && sTriggerList.Contains("CMUL") )
1186 AliError(Form("Found no relevant trigger for run %09d",runNumber));
1190 AliInfo(Form("Integrated lumi %7.4f %s^-1",intLumi,csUnit));
1192 std::map<std::string,float>::const_iterator it;
1194 for ( it = lumiPerTrigger.begin(); it != lumiPerTrigger.end(); ++it )
1196 AliInfo(Form("Trigger %30s Lumi %10.4f %s^-1 duration %10ld s",it->first.c_str(),it->second,csUnit,durationPerTrigger[it->first]));
1201 std::map<int, std::map<std::string, float> >::const_iterator fit;
1203 lumiPerTrigger.clear();
1205 for ( fit = lumiPerFillPerTrigger.begin(); fit != lumiPerFillPerTrigger.end(); ++fit )
1207 int fill = fit->first;
1208 std::map<std::string,float>::const_iterator tit;
1209 (*out) << fill << sep;
1211 for ( tit = fit->second.begin(); tit != fit->second.end(); ++tit )
1213 (*out) << Form("%e",tit->second) << sep;
1216 (*out) << sep; // comment (empty)
1218 for ( tit = fit->second.begin(); tit != fit->second.end(); ++tit )
1220 lumiPerTrigger[tit->first] += tit->second;
1222 (*out) << Form("%e",lumiPerTrigger[tit->first]) << sep;
1224 (*out) << sep << "0" << sep << "0" << sep << "0" << sep << "0" << sep << std::endl; // LHC per fill, LHC integrated, lumi tot muon , efficiency
1228 delete triggerArray;
1231 //______________________________________________________________________________
1232 TGraph* AliAnalysisTriggerScalers::MakeGraph(const std::vector<int>& vx,
1233 const std::vector<int>& vex,
1234 const std::vector<double>& vy,
1235 const std::vector<double>& vey)
1237 /// Build a graph from a set of stl vectors
1239 if ( ! ( vx.size() == vex.size() && vx.size() == vy.size() && vx.size() == vey.size() ) )
1241 std::cerr << "incompatible sizes" << std::endl;
1245 double* x = new double[vx.size()];
1246 double* ex = new double[vx.size()];
1247 double* y = new double[vx.size()];
1248 double* ey = new double[vx.size()];
1250 for ( size_t i = 0; i < vx.size(); ++i )
1258 TGraph* g = new TGraphErrors(vx.size(),x,y,ex,ey);
1270 //______________________________________________________________________________
1271 Double_t AliAnalysisTriggerScalers::Mu(Double_t L0B, Double_t Nb)
1273 /// L0B = trigger rate before any veto
1274 /// Nb = number of crossing bunches
1276 Double_t p0 = 1-L0B/(Nb*11245.0); // proba to get *no* collision per bunch crossing
1278 return -TMath::Log(p0);
1281 //______________________________________________________________________________
1282 Int_t AliAnalysisTriggerScalers::NumberOfInteractingBunches(const AliLHCData& lhc, Int_t runNumber, Bool_t mainSat) const
1284 /// Extract the number of colliding bunches from the LHC data
1285 /// Use mainSat = true if the collisons were main - satellite ones
1287 Int_t numberOfInteractingBunches(0);
1288 Int_t numberOfInteractingBunchesMeasured(0);
1294 AliLHCDipValI* val = lhc.GetBunchConfigDeclared(beam1,0);
1295 AliLHCDipValI* valm = lhc.GetBunchConfigMeasured(beam1,0);
1300 numberOfInteractingBunches = val->GetSizeTotal();
1301 numberOfInteractingBunchesMeasured = valm->GetSizeTotal();
1302 nIBM2 = numberOfInteractingBunches;
1306 for ( Int_t i = 0; i < val->GetSizeTotal(); ++i )
1308 if ( val->GetValue(i) < 0 ) ++numberOfInteractingBunches;
1311 for ( Int_t i = 0; i < valm->GetSizeTotal(); ++i )
1313 if ( valm->GetValue(i) < 0 ) ++numberOfInteractingBunchesMeasured;
1316 valm = lhc.GetBunchConfigMeasured(beam2,0);
1318 for ( Int_t i = 0; i < valm->GetSizeTotal(); ++i )
1320 if ( valm->GetValue(i) <= 0 ) ++nIBM2;
1324 if (!numberOfInteractingBunches)
1329 if ( numberOfInteractingBunches != numberOfInteractingBunchesMeasured ||
1330 numberOfInteractingBunches != nIBM2 )
1332 Int_t delta = TMath::Max(numberOfInteractingBunches - numberOfInteractingBunchesMeasured,numberOfInteractingBunches-nIBM2);
1334 if ( 1.0*TMath::Abs(delta)/numberOfInteractingBunches > 0.05 ) // more than 5% difference
1336 AliWarning(Form("Got some different number of interacting bunches for fill %d run %d ! NumberOfInteractingBunches=%3d NumberOfInteractingBunchesMeasured=%3d NIBM2=%3d. Will use %d",
1337 lhc.GetFillNumber(),runNumber,
1338 numberOfInteractingBunches,numberOfInteractingBunchesMeasured,nIBM2,numberOfInteractingBunches));
1342 AliDebug(1,Form("Got some different number of interacting bunches for fill %d run %d ! NumberOfInteractingBunches=%3d NumberOfInteractingBunchesMeasured=%3d NIBM2=%3d. Will use %d",
1343 lhc.GetFillNumber(),runNumber,
1344 numberOfInteractingBunches,numberOfInteractingBunchesMeasured,nIBM2,numberOfInteractingBunches));
1349 return numberOfInteractingBunches;
1353 //______________________________________________________________________________
1354 TGraph* AliAnalysisTriggerScalers::PlotTrigger(const char* triggerClassName,
1357 // Plot one of the scalers (L0A,L0B,L0AOVERB,etc...) of a given triggerClass
1358 // Get one value per run.
1360 std::vector<Float_t> x;
1361 std::vector<Float_t> y;
1366 for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
1368 Int_t runNumber = fRunList[i];
1370 AliAnalysisTriggerScalerItem* s = GetTriggerScaler(runNumber,what,triggerClassName);
1374 x.push_back(runNumber);
1376 Double_t value = s->ValueCorrectedForDownscale();
1378 if ( TString(what).Contains("RATE") )
1389 if ( fRunList.size() ) {
1390 mean /= fRunList.size();
1393 AliInfo(Form("Integral %e mean %e",integral,mean));
1395 return new TGraph(x.size(),&x[0],&y[0]);
1398 //______________________________________________________________________________
1399 TGraph* AliAnalysisTriggerScalers::PlotTriggerEvolution(const char* triggerClassName,
1405 /// Make a graph of "what" for a given trigger class.
1409 /// - L0B (level 0 before any veto applied)
1410 /// - L0A (level 0 after all vetoes)
1411 /// - L0AOVERB (L0A/L0B)
1412 /// - mu ( = -TMath::Log( 1 - P(0) ) where P(0) is the proba to have zero collisions in a bunch crossing)
1413 /// - pileupfactor = mu/(1-exp(-mu)) : the factor to apply to correct the cint1b count rate
1414 /// - vsnb = L0B/(NumberOfInteractingBunches*11245)
1415 /// - NBCX = NumberOfInteractingBunches
1417 TString swhat(what);
1420 if ( swhat.Contains(";"))
1422 swhat.ReplaceAll(";",",");
1423 AliWarningClass("; is not a valid separator, replaced it with ,");
1429 TObjArray* a = swhat.Tokenize(",");
1430 if (a->GetEntries()>1)
1434 TObjString* str(0x0);
1435 Double_t ymin(TMath::Limits<Double_t>::Max());
1437 Double_t xmin(TMath::Limits<Double_t>::Max());
1441 while ( ( str = static_cast<TObjString*>(next())) )
1443 g = PlotTriggerEvolution(triggerClassName,str->String().Data(),false,0x0,removeZeros);
1445 for ( Int_t i = 0; i < g->GetN(); ++i )
1447 ymin = TMath::Min(ymin,g->GetY()[i]);
1448 ymax = TMath::Max(ymax,g->GetY()[i]);
1450 xmin = TMath::Min(xmin,g->GetX()[0]);
1451 xmax = TMath::Max(xmax,g->GetX()[g->GetN()-1]);
1454 gStyle->SetOptTitle(0);
1456 AliInfoClass(Form("x %e ; %e y %e ; %e",xmin,xmax,ymin,ymax));
1457 TH2* h = new TH2F("h",triggerClassName,100,xmin,xmax,100,ymin,ymax);
1458 h->SetStats(kFALSE);
1459 h->GetXaxis()->SetTimeDisplay(1);
1460 h->GetXaxis()->SetTimeFormat("%d/%m %H:%M");
1461 h->GetXaxis()->SetTimeOffset(0,"gmt");
1462 h->GetXaxis()->SetNdivisions(505);
1465 TIter nextGraph(&graphs);
1467 while ( ( g = static_cast<TGraph*>(nextGraph())) )
1469 AliInfoClass(g->GetTitle());
1471 g->SetLineColor(color);
1472 g->SetMarkerColor(color);
1473 g->SetMarkerStyle(marker);
1482 std::vector<int> vx;
1483 std::vector<int> vex;
1484 std::vector<double> vy;
1485 std::vector<double> vey;
1491 for ( std::vector<int>::size_type iRun = 0; iRun < fRunList.size(); ++iRun )
1493 Int_t runNumber = fRunList[iRun];
1495 AliTriggerConfiguration* tc(0x0);
1496 AliTriggerRunScalers* trs(0x0);
1497 AliLHCData* lhc(0x0);
1499 GetCTPObjects(runNumber,tc,trs,lhc);
1501 const TObjArray* scalers = trs->GetScalersRecords();
1503 const TObjArray& trClasses = tc->GetClasses();
1504 TIter next(&trClasses);
1505 AliTriggerClass* triggerClass;
1507 while ( ( triggerClass = static_cast<AliTriggerClass*>(next()) ) )
1509 UChar_t index = GetIndex(triggerClass->GetMask());
1511 if ( !TString(triggerClass->GetName()).Contains(triggerClassName) ) continue;
1513 Bool_t mainSat = IsMainSatelliteCollision(triggerClass->GetName());
1515 Int_t numberOfInteractingBunches = NumberOfInteractingBunches(*lhc,runNumber,mainSat);
1517 TIter nextScaler(scalers);
1518 AliTriggerScalersRecord* record;
1526 Bool_t first(kTRUE);
1528 while ( ( record = static_cast<AliTriggerScalersRecord*>(nextScaler()) ) )
1530 const AliTriggerScalers* scaler = record->GetTriggerScalersForClass(index);
1533 const AliTimeStamp* ats = record->GetTimeStamp();
1535 UInt_t seconds = ats->GetSeconds();// - TTimeStamp::GetZoneOffset();
1537 TTimeStamp ts(seconds,ats->GetMicroSecs());
1539 UInt_t l0b = scaler->GetLOCB() - refl0b;
1540 UInt_t l0a = scaler->GetLOCA() - refl0a;
1541 UInt_t l1b = scaler->GetL1CB() - refl1b;
1542 UInt_t l1a = scaler->GetL1CA() - refl1a;
1543 UInt_t l2b = scaler->GetL2CB() - refl2b;
1544 UInt_t l2a = scaler->GetL2CA() - refl2a;
1545 UInt_t timelapse = seconds - reft;
1547 if ( removeZeros && ( l0b <= 2 || ( l0a <= 2 && l0a != 0 ) || timelapse <= 9 ) )
1549 AliInfo(Form("Skipping point for %s l0b %d l0a %d timelapse %d ts %s",
1550 triggerClassName,l0b,l0a,timelapse,ts.AsString()));
1555 refl0b = scaler->GetLOCB();
1556 refl0a = scaler->GetLOCA();
1557 refl1b = scaler->GetL1CB();
1558 refl1a = scaler->GetL1CA();
1559 refl2b = scaler->GetL2CB();
1560 refl2a = scaler->GetL2CA();
1571 if (swhat.Contains("L0AOVERB") )
1573 value = l0a*1.0/l0b;
1576 error = value*TMath::Sqrt(1.0/l0b+1.0/l0a);
1583 else if ( swhat.Contains("L0B") )
1587 else if (swhat.Contains("L0A") )
1591 else if ( swhat.Contains("L1B") )
1595 else if (swhat.Contains("L1A") )
1599 else if ( swhat.Contains("L2B") )
1603 else if (swhat.Contains("L2A") )
1607 else if ( swhat.Contains("MU") )
1609 AliInfo(Form("run %d seconds %d L0b %d NB %d",runNumber,timelapse,l0b,numberOfInteractingBunches));
1610 if (timelapse==0) continue;
1611 value = Mu(l0b/timelapse,numberOfInteractingBunches);
1612 error = 0.0; // FIXME
1614 else if ( swhat.Contains("PILEUPFACTOR") )
1616 if (timelapse==0) continue;
1617 Double_t mu = Mu(l0b/timelapse,numberOfInteractingBunches);
1618 value = mu/(1-TMath::Exp(-mu));
1619 error = -1.0; // FIXME
1621 else if ( swhat.Contains("VSNB") )
1623 value = l0b/(11245.0*numberOfInteractingBunches);
1624 error = -1.0; // FIXME
1626 else if ( swhat.Contains("NBCX"))
1628 value = numberOfInteractingBunches;
1629 error = 1.0; // FIXME
1634 AliInfo(Form("Unknown what %s",what));
1637 if ( ! swhat.Contains("OVER") && ! swhat.Contains("RATIO") &&
1638 ! swhat.Contains("MU") && ! swhat.Contains("PILEUPFACTOR") &&
1639 ! swhat.Contains("RAW") & ! swhat.Contains("NBCX") )
1644 if ( !TMath::Finite(value) ) continue;
1646 if ( !swhat.Contains("L0AOVERB") && error > 0 )
1648 error = value >0 ? 1.0/TMath::Sqrt(value) : 0.0;
1651 if (removeZeros && TMath::Abs(value) < 1E-6 )
1661 vx.push_back(seconds);
1664 vy.push_back(value);
1666 vey.push_back(error);
1674 if (mean && nvalues)
1679 if ( vx.empty() ) return 0;
1681 TGraph* g = MakeGraph(vx,vex,vy,vey);
1682 TString title(Form("TriggerEvolution%s-%s",triggerClassName,swhat.Data()));
1684 g->SetName(title.Data());
1685 g->SetTitle(title.Data());
1687 g->GetYaxis()->SetTitle(title.Data());
1691 TCanvas* c = NewCanvas(g->GetName());
1692 g->SetLineColor(color);
1693 g->SetMarkerColor(color);
1694 g->SetMarkerStyle(marker);
1696 c->SaveAs(Form("%s.png",title.Data()));
1703 //______________________________________________________________________________
1704 TGraph* AliAnalysisTriggerScalers::PlotTriggerRatio(const char* triggerClassName1,
1706 const char* triggerClassName2,
1709 // Plot the ratio of two scalers.
1710 // Get one value per run.
1712 std::vector<Float_t> x;
1713 std::vector<Float_t> y;
1715 for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
1717 Int_t runNumber = fRunList[i];
1719 AliAnalysisTriggerScalerItem* s1 = GetTriggerScaler(runNumber,what1,triggerClassName1);
1720 AliAnalysisTriggerScalerItem* s2 = GetTriggerScaler(runNumber,what2,triggerClassName2);
1722 if (!s1 || !s2) continue;
1724 x.push_back(runNumber);
1728 if ( s2->ValueCorrectedForDownscale() > 0 )
1730 ratio = s1->ValueCorrectedForDownscale()/s2->ValueCorrectedForDownscale();
1735 AliDebug(1,Form("RUN %09d %20s (%s) %12llu (%5d) %20s (%s) %12llu (%5d) R %7.2f",
1737 triggerClassName1,what1,
1738 s1->Value(),s1->DownscalingFactor(),
1739 triggerClassName2,what2,
1740 s2->Value(),s2->DownscalingFactor(),
1744 return new TGraph(x.size(),&x[0],&y[0]);
1747 //______________________________________________________________________________
1748 TGraph* AliAnalysisTriggerScalers::PlotTriggerRatioEvolution(const char* triggerClassName1,
1750 const char* triggerClassName2,
1753 /// Plots the evolution of one trigger ratio
1755 Bool_t draw(kFALSE);
1756 Bool_t removeZeros(kFALSE);
1758 TGraph* g1 = PlotTriggerEvolution(triggerClassName1,what1,draw,0x0,removeZeros);
1759 TGraph* g2 = PlotTriggerEvolution(triggerClassName2,what2,draw,0x0,removeZeros);
1761 if (!g1 || !g2) return 0x0;
1763 if ( g1->GetN() != g2->GetN() )
1765 AliError("Oups. Did not get the same number of points for the 2 graphs ?!");
1771 Double_t x1err,x2err;
1772 Double_t y1err,y2err;
1774 TGraphErrors* g = new TGraphErrors(g1->GetN());
1777 for ( Int_t i = 0; i < g1->GetN(); ++i )
1779 g1->GetPoint(i,x1,y1);
1780 g2->GetPoint(i,x2,y2);
1782 x1err = g1->GetErrorX(i);
1783 x2err = g2->GetErrorX(i);
1785 y1err = g1->GetErrorY(i);
1786 y2err = g2->GetErrorY(i);
1790 AliError(Form("Points at %d don't have the same x ! ",i));
1796 if ( TMath::Abs(y2) > 1E-12 )
1803 if ( TMath::Abs(x1) > 1E-12)
1805 yerr += TMath::Sqrt(x1err*x1err/(x1*x1));
1808 if ( TMath::Abs(x2) > 1E-12)
1810 yerr += TMath::Sqrt(x2err*x2err/(x2*x2));
1815 g->SetPoint(j,x1,y);
1816 g->SetPointError(j,x1err,yerr);
1829 //______________________________________________________________________________
1830 void AliAnalysisTriggerScalers::Print(Option_t* /* opt */) const
1832 /// print our runlist
1833 AliAnalysisTriggerScalers::PrintIntegers(fRunList,',');
1836 //______________________________________________________________________________
1837 void AliAnalysisTriggerScalers::PrintIntegers(const std::vector<int>& integers,
1841 /// print a list of integers
1842 for ( std::vector<int>::size_type i = 0; i < integers.size(); ++i )
1844 out << integers[i] << sep;
1849 //______________________________________________________________________________
1850 void AliAnalysisTriggerScalers::ReadIntegers(const char* filename,
1851 std::vector<int>& integers,
1854 /// Read integers from filename, where integers are either
1855 /// separated by "," or by return carriage
1857 if ( gSystem->AccessPathName(filename)==kTRUE )
1861 std::ifstream in(gSystem->ExpandPathName(filename));
1864 std::set<int> runset;
1868 for ( std::vector<int>::size_type j = 0; j < integers.size(); ++ j )
1870 runset.insert(integers[j]);
1876 in.getline(line,10000,'\n');
1878 TString sline(line);
1880 if (sline.Contains(","))
1882 TObjArray* a = sline.Tokenize(",");
1885 while ( ( s = static_cast<TObjString*>(next()) ) )
1887 runset.insert(s->String().Atoi());
1893 runset.insert(sline.Atoi());
1903 for ( std::set<int>::const_iterator it = runset.begin(); it != runset.end(); ++it )
1905 integers.push_back((*it));
1908 std::sort(integers.begin(),integers.end());
1912 //______________________________________________________________________________
1913 void AliAnalysisTriggerScalers::SetRunList(Int_t runNumber)
1915 // Make the runlist be a single run
1917 fRunList.push_back(runNumber);
1920 //______________________________________________________________________________
1921 void AliAnalysisTriggerScalers::SetRunList(const std::vector<int>& runs)
1923 // Make the runlist be a single run
1927 //______________________________________________________________________________
1928 void AliAnalysisTriggerScalers::SetRunList(const std::set<int>& runs)
1930 // Make the runlist be a single run
1932 for ( std::set<int>::const_iterator it = runs.begin(); it != runs.end(); ++it )
1934 fRunList.push_back(*it);
1938 //______________________________________________________________________________
1939 void AliAnalysisTriggerScalers::SetRunList(const char* runlist)
1941 // Read the runlist from an ASCII file or a comma separated list
1942 // or a space separated list
1946 if ( TString(runlist).Contains(",") || TString(runlist).Contains(" ") )
1948 TObjArray* runs = 0x0;
1949 if ( TString(runlist).Contains(",") )
1951 runs = TString(runlist).Tokenize(",");
1955 runs = TString(runlist).Tokenize(" ");
1959 std::set<int> runset;
1961 while ( ( s = static_cast<TObjString*>(next()) ) )
1963 runset.insert(s->String().Atoi());
1966 for ( std::set<int>::const_iterator it = runset.begin(); it != runset.end(); ++it )
1968 fRunList.push_back((*it));
1971 std::sort(fRunList.begin(),fRunList.end());
1977 ReadIntegers(runlist,fRunList);
1978 if (fRunList.empty())
1980 if ( TString(runlist).IsDigit() )
1982 SetRunList(TString(runlist).Atoi());
1986 AliError("Could not set run list !");
1992 //______________________________________________________________________________
1993 //______________________________________________________________________________
1994 //______________________________________________________________________________
1995 //______________________________________________________________________________
1996 //______________________________________________________________________________
1997 //______________________________________________________________________________
1998 //______________________________________________________________________________
2000 ClassImp(AliAnalysisTriggerScalerItem)
2002 //______________________________________________________________________________
2004 AliAnalysisTriggerScalerItem::BCMaskName() const
2006 // return bc mask name
2007 if ( BCMask() ) return BCMask()->GetName(); else return "";
2010 //______________________________________________________________________________
2011 Int_t AliAnalysisTriggerScalerItem::Compare(const TObject* obj) const
2013 /// compare two scaler items (by means of their run number only)
2014 const AliAnalysisTriggerScalerItem* s = static_cast<const AliAnalysisTriggerScalerItem*>(obj);
2016 if ( s->RunNumber() < RunNumber() )
2020 else if ( s->RunNumber() > RunNumber() )
2027 //______________________________________________________________________________
2028 void AliAnalysisTriggerScalerItem::Print(Option_t* opt) const
2033 if ( RunNumber() > 0 )
2037 if ( sopt.Contains("FULL") )
2042 std::cout << Form("RUN %6d CLASS %24s (%5s %4d BCs) SCALER %s %12llu DS %6d DURATION %ld",
2043 RunNumber(),TriggerClassName(),
2045 BCMask() ? BCMask()->GetNUnmaskedBCs() : 0,
2047 Value(),DownscalingFactor(),Duration()) << std::endl;
2052 std::cout << Form("CLASS %24s SCALER %20llu NRUNS %d",
2053 TriggerClassName(),Value(),NofRuns()) << std::endl;