]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/muondep/AliAnalysisTriggerScalers.cxx
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWG / muondep / AliAnalysisTriggerScalers.cxx
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
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 **************************************************************************/
15
16 // $Id$
17
18 //
19 // AliAnalysisTriggerScalers : utility class to play with OCDB scalers.
20 //
21 // Can use it to e.g. :
22 //
23 // - get the integrated luminosity for a given period
24 //   (see IntegratedLuminosity method)
25 //
26 // - plot the trigger rates (see PlotTriggerEvolution)
27 //
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 :
33 //
34 // GRP/CTP/Config
35 // GRP/CTP/Scalers
36 //
37 // and also (to get the fill and period ranges, mainly for drawing)
38 //
39 // GRP/GRP/Data
40 // GRP/GRP/LHCData
41 //
42 //
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.
45 //
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"...
48 //
49 // Ideas for improvement :
50 //
51 // - make sure the "what" parameter mean the same thing in all methods
52 //   and so can take the same values in all methods
53 //
54 // - get the lumi trigger name and cross-section from e.g. OADB instead of
55 //   hard-coding them ?
56 //
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 ?)
60 //
61 // - robustify the PAC fraction computation
62 //   
63 //
64 // L. Aphecetche (Subatech)
65 //
66
67 #include "AliAnalysisTriggerScalers.h"
68 #include "AliCDBEntry.h"
69 #include "AliCDBManager.h"
70 #include "AliGRPObject.h"
71 #include "AliLog.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"
84 #include "TDatime.h"
85 #include "TError.h"
86 #include "TGraph.h"
87 #include "TGraphErrors.h"
88 #include "TH2.h"
89 #include "TObjArray.h"
90 #include "TObjString.h"
91 #include <set>
92 #include <map>
93 #include "TStyle.h"
94 #include "AliLHCData.h"
95 #include "TMath.h"
96 #include "TAxis.h"
97 #include "TBox.h"
98 #include "TCanvas.h"
99 #include "TText.h"
100 #include <sstream>
101
102 using std::make_pair;
103
104 ClassImp(AliAnalysisTriggerScalers)
105
106 namespace {
107   
108   //______________________________________________________________________________
109   UChar_t GetIndex(ULong64_t mask)
110   {
111     for ( Int_t i = 0; i < 64; ++i )
112     {
113       if ( mask & ( 1ull << i ) ) return i+1;
114     }
115     return 0;
116   }
117   
118   //______________________________________________________________________________
119   TCanvas* NewCanvas(const char* name)
120   {
121     TCanvas* c = new TCanvas(name,name);
122     c->SetLeftMargin(0.15);
123     return c;
124   }
125   
126   //______________________________________________________________________________
127   void TimeAxis(TGraph* g)
128   {
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);
133     
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);
139
140   }
141
142   //______________________________________________________________________________
143   Bool_t IsMainSatelliteCollision(const char* triggerClassName)
144   {
145     TString tcn(triggerClassName);
146     tcn.ToUpper();
147     
148     return (tcn.Contains("-S-") ||  tcn.Contains("-SC-") || tcn.Contains("-SA-"));
149   }
150 }
151
152 //_____________________________________________________________________________
153 AliAnalysisTriggerScalers::AliAnalysisTriggerScalers(Int_t runNumber, const char* ocdbPath) : fRunList(), fOCDBPath(ocdbPath), fShouldCorrectForPileUp(kTRUE), fCrossSectionUnit("UB")
154 {
155   // ctor using a single run number
156   SetRunList(runNumber);
157 }
158
159 //_____________________________________________________________________________
160 AliAnalysisTriggerScalers::AliAnalysisTriggerScalers(const std::vector<int>& runs, const char* ocdbPath) : fRunList(), fOCDBPath(ocdbPath), fShouldCorrectForPileUp(kTRUE), fCrossSectionUnit("UB")
161 {
162   // ctor using a vector of run numbers
163   SetRunList(runs);
164 }
165
166 //_____________________________________________________________________________
167 AliAnalysisTriggerScalers::AliAnalysisTriggerScalers(const std::set<int>& runs, const char* ocdbPath) : fRunList(), fOCDBPath(ocdbPath), fShouldCorrectForPileUp(kTRUE), fCrossSectionUnit("UB")
168 {
169   // ctor using a set of run numbers
170   SetRunList(runs);
171 }
172
173 //_____________________________________________________________________________
174 AliAnalysisTriggerScalers::AliAnalysisTriggerScalers(const char* runlist, const char* ocdbPath) : fRunList(), fOCDBPath(ocdbPath), fShouldCorrectForPileUp(kTRUE), fCrossSectionUnit("UB")
175 {
176   // ctor from a run list (txt file)
177   SetRunList(runlist);
178 }
179
180 //_____________________________________________________________________________
181 AliAnalysisTriggerScalers::~AliAnalysisTriggerScalers()
182 {
183   // dtor
184 }
185
186 //______________________________________________________________________________
187 Bool_t AliAnalysisTriggerScalers::CheckRecord(const AliTriggerScalersRecord& record,
188                                               UInt_t index,
189                                               UInt_t refb,
190                                               UInt_t refa,
191                                               UInt_t timelapse) const
192 {
193   /// Check if a trigger scaler record should be kept for our
194   /// luminosity computations
195   
196   const AliTriggerScalers* scaler = record.GetTriggerScalersForClass(index);
197   
198   UInt_t a = scaler->GetLOCA() - refa;
199   UInt_t b = scaler->GetLOCB() - refb;
200   
201   Bool_t ok(kTRUE);
202   
203   if ( b <= 2 || ( a <= 2 ) || timelapse <= 9 ) // empirical cuts
204   {
205     ok = kFALSE;
206   }
207   
208   return ok;
209 }
210
211
212 //______________________________________________________________________________
213 void AliAnalysisTriggerScalers::DrawFill(Int_t run1, Int_t run2, double ymin, double ymax, const char* label, Int_t color)
214 {
215   // Draw a yellow box to indicate a run range
216   
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);
221   b->Draw();
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);
227   text->Draw();
228 }
229
230 //_____________________________________________________________________________
231 void AliAnalysisTriggerScalers::DrawFills(Double_t ymin, Double_t ymax, Int_t color)
232 {
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
236   /// it will be long.
237   ///
238   std::map<int, std::pair<int,int> > fills;
239   
240   GetFillBoundaries(fills);
241   
242   std::map<int, std::pair<int,int> >::const_iterator it;
243   
244   for ( it = fills.begin(); it != fills.end(); ++it )
245   {
246     const std::pair<int,int>& p = it->second;
247     TString fillnumber;
248     fillnumber.Form("%d",it->first);
249     DrawFill(p.first,p.second,ymin,ymax,fillnumber.Data(),color);
250   }
251 }
252
253 //_____________________________________________________________________________
254 void AliAnalysisTriggerScalers::DrawPeriods(Double_t ymin, Double_t ymax, Int_t color)
255 {
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
259   /// it will be long.
260
261   std::map<std::string, std::pair<int,int> > periods;
262   
263   GetLHCPeriodBoundaries(periods);
264   
265   for ( std::map<std::string, std::pair<int,int> >::const_iterator it = periods.begin(); it != periods.end(); ++it )
266   {
267     const std::pair<int,int>& p = it->second;
268     DrawFill(p.first,p.second,ymin,ymax,it->first.c_str(),color);
269   }
270
271 }
272
273 //_____________________________________________________________________________
274 void AliAnalysisTriggerScalers::GetCTPObjects(Int_t runNumber,
275                                               AliTriggerConfiguration*& tc,
276                                               AliTriggerRunScalers*& trs,
277                                               AliLHCData*& lhc) const
278 {
279   /// For a given run, get the CTP objects we need to do our work
280   ///
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
283   
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));
287 }
288
289 //_____________________________________________________________________________
290 void AliAnalysisTriggerScalers::GetFillBoundaries(std::map<int, std::pair<int,int> >& fills)
291 {
292   /// Given a list of runs get the fills they are in
293   
294   fills.clear();
295   
296   for ( std::vector<int>::const_iterator it = fRunList.begin(); it != fRunList.end(); ++it )
297   {
298     int runnumber = *it;
299     
300     int fill = GetFillNumberFromRunNumber(runnumber);
301     
302     if (fills.count(fill))
303     {
304       std::pair<int,int>& p = fills[fill];
305       p.first = TMath::Min(runnumber,p.first);
306       p.second = TMath::Max(runnumber,p.second);
307     }
308     else
309     {
310       fills[fill] = make_pair(runnumber,runnumber);
311     }
312   }
313 }
314
315
316 //______________________________________________________________________________
317 Int_t AliAnalysisTriggerScalers::GetFillNumberFromRunNumber(Int_t runNumber)
318 {
319   /// Get the fill number of a run
320   
321   AliLHCData* lhcData = static_cast<AliLHCData*>(GetOCDBObject("GRP/GRP/LHCData",runNumber));
322   if (lhcData)
323   {
324     Int_t fillNumber = lhcData->GetFillNumber();
325     if ( fillNumber == 0 && runNumber == 189694)
326     {
327       // manual hack because GRP info incorrect for this run ?
328       fillNumber = 3135;
329     }
330
331     return fillNumber;
332   }
333   return -1;
334 }
335
336 //______________________________________________________________________________
337 TObject* AliAnalysisTriggerScalers::GetOCDBObject(const char* path, Int_t runNumber) const
338 {
339   // Helper method to get an object from the OCDB
340   
341   if ( !AliCDBManager::Instance()->IsDefaultStorageSet() )
342   {
343     AliInfo(Form("Setting OCDB default storage to %s",fOCDBPath.Data()));
344     AliCDBManager::Instance()->SetDefaultStorage(fOCDBPath.Data());
345   }
346   
347   AliCDBManager::Instance()->SetRun(runNumber);
348   
349   gErrorIgnoreLevel=kWarning; // to avoid all the TAlienFile::Open messages...
350   
351   AliCDBEntry* e = AliCDBManager::Instance()->Get(path);
352   return e ? e->GetObject() : 0x0;
353 }
354
355
356 //______________________________________________________________________________
357 TString AliAnalysisTriggerScalers::GetLHCPeriodFromRunNumber(Int_t runNumber) const
358 {
359   /// Get a list of LHC periods from a list of run numbers
360   
361   AliGRPObject* grp = static_cast<AliGRPObject*>(GetOCDBObject("GRP/GRP/Data",runNumber));
362
363   if (!grp) return "";
364   
365   return grp->GetLHCPeriod();
366 }
367
368 //_____________________________________________________________________________
369 void
370 AliAnalysisTriggerScalers::GetLuminosityTriggerAndCrossSection(Int_t runNumber,
371                                                                TString& lumiTriggerClassName,
372                                                                Double_t& lumiTriggerCrossSection,
373                                                                Double_t& lumiTriggerCrossSectionError) const
374 {
375   /// For one given run, get the trigger class name to be used as the luminometer,
376   /// and its corresponding cross-section
377   ///
378   /// FIXME: all is hard-coded for now... (use OADB for this ?)
379   
380   lumiTriggerClassName="";
381   lumiTriggerCrossSection=0.0;
382   lumiTriggerCrossSectionError=0.0; // FIXME
383   
384   TString period = GetLHCPeriodFromRunNumber(runNumber);
385   
386   if ( period.BeginsWith("LHC11") )
387   {
388     AliError("Not implemented yet");
389   }
390   else if ( period.BeginsWith("LHC12") )
391   {
392      if ( period == "LHC12h" || period == "LHC12i" )
393      {
394        // pp 8 TeV main-satellites
395        lumiTriggerClassName = "C0TVX-S-NOPF-ALLNOTRD";
396        lumiTriggerCrossSection = 28.0; // FIXME: not from a vdM yet       
397      }
398     else
399     {
400       AliError("Not implemented yet");
401     }
402   }
403   else if ( period.BeginsWith("LHC13") )
404   {
405     if ( period == "LHC13g" )
406     {
407       // pp 2.76 TeV
408       lumiTriggerClassName = "C0TVX-B-NOPF-ALLNOTRD";
409       lumiTriggerCrossSection = 18.0; // FIXME: not from a vdM yet
410     }
411     else
412     {
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
416     }
417   }
418   else
419   {
420     AliError("Not implemented yet");
421   }
422   
423   float csMult(1.0);
424   
425   if ( fCrossSectionUnit=="PB" )
426   {
427     csMult=1E9;
428   }
429   else if (fCrossSectionUnit=="NB")
430   {
431     csMult=1E6;
432   }
433   else if ( fCrossSectionUnit=="UB" )
434   {
435     csMult=1E3;
436   }
437   else if ( fCrossSectionUnit=="MB" )
438   {
439     csMult=1.0;
440   }
441
442   lumiTriggerCrossSection *= csMult;
443   lumiTriggerCrossSectionError *= csMult;
444 }
445
446 //_____________________________________________________________________________
447 void AliAnalysisTriggerScalers::GetLHCPeriodBoundaries(std::map<std::string, std::pair<int,int> >& periods)
448 {
449   /// Given a list of runs get the fills they are in
450   
451   periods.clear();
452   
453   for ( std::vector<int>::const_iterator it = fRunList.begin(); it != fRunList.end(); ++it )
454   {
455     int runnumber = *it;
456     
457     std::string period = GetLHCPeriodFromRunNumber(runnumber).Data();
458     
459     if (periods.count(period))
460     {
461       std::pair<int,int>& p = periods[period];
462       p.first = TMath::Min(runnumber,p.first);
463       p.second = TMath::Max(runnumber,p.second);
464     }
465     else
466     {
467       periods[period] = make_pair(runnumber,runnumber);
468     }
469   }
470 }
471
472
473 //______________________________________________________________________________
474 Int_t AliAnalysisTriggerScalers::GetTriggerInput(Int_t runNumber, const char* inputname)
475 {
476   /// Get one trigger input for a given run
477   AliTriggerConfiguration* tc = static_cast<AliTriggerConfiguration*>(GetOCDBObject("GRP/CTP/Config",runNumber));
478   if (!tc) return -1;
479   
480   const TObjArray& inputs = tc->GetInputs();
481   
482   AliTriggerInput* ti = static_cast<AliTriggerInput*>(inputs.FindObject(inputname));
483   
484   if (!ti) return -1;
485   
486   return ti->GetSignature();
487 }
488
489 //______________________________________________________________________________
490 Float_t
491 AliAnalysisTriggerScalers::GetPauseAndConfigCorrection(Int_t runNumber, const char* triggerClassName)
492 {
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.
495   
496   ULong64_t total(0);
497   ULong64_t nopac(0);
498   
499   AliTriggerConfiguration* tc = static_cast<AliTriggerConfiguration*>(GetOCDBObject("GRP/CTP/Config",runNumber));
500
501   AliTriggerRunScalers* trs = static_cast<AliTriggerRunScalers*>(GetOCDBObject("GRP/CTP/Scalers",runNumber));
502   const TObjArray* scalers = trs->GetScalersRecords();
503   
504   TIter next(scalers);
505   AliTriggerScalersRecord* record;
506   UInt_t reft(0);
507   UInt_t refl0b(0);
508   UInt_t refl0a(0);
509   Bool_t first(kTRUE);
510
511   Int_t classIndex = tc->GetClassIndexFromName(triggerClassName);
512   
513   while ( ( record = static_cast<AliTriggerScalersRecord*>(next()) ) )
514   {
515     const AliTriggerScalers* scaler = record->GetTriggerScalersForClass(classIndex);
516         
517     const AliTimeStamp* ats = record->GetTimeStamp();
518     
519     UInt_t seconds = ats->GetSeconds();// - TTimeStamp::GetZoneOffset();
520     
521     TTimeStamp ts(seconds,ats->GetMicroSecs());
522     
523     UInt_t l0b = scaler->GetLOCB() - refl0b;
524     UInt_t l0a = scaler->GetLOCA() - refl0a;
525     UInt_t timelapse = seconds - reft;
526
527     if ( l0b <=2 || timelapse < 9 ) continue;
528
529     reft = seconds;
530     refl0b = scaler->GetLOCB();
531     refl0a = scaler->GetLOCA();
532
533     if ( first )
534     {
535       first = kFALSE;
536       continue;
537     }
538     
539     total += l0b;
540     
541     if ( l0a > 0 )
542     {
543       nopac += l0b;
544     }
545   }
546
547   return total > 0 ? nopac*1.0/total : 0.0;
548 }
549
550 //______________________________________________________________________________
551 void AliAnalysisTriggerScalers::ShowPileUpFactors(const char* triggerClassName, Double_t purity)
552 {
553   /// Printout the pile-up factors for the runlist, including some purity if needed
554   
555   if (purity<=0.0)
556   {
557     AliError(Form("Cannot work with purity=%f for trigger %s. Should be strictly positive",purity,triggerClassName));
558     return;
559   }
560
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
564   
565   for ( std::vector<int>::const_iterator it = rl.begin(); it != rl.end(); ++it )
566   {
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());
572   }
573
574   for ( std::vector<std::string>::size_type i = 0; i  < lines.size(); ++i )
575   {
576     std::cout << lines[i].c_str() << std::endl;
577   }
578 }
579
580 //______________________________________________________________________________
581 void AliAnalysisTriggerScalers::GetPileUpFactor(Int_t runNumber, const char* triggerClassName,
582                                                 Double_t purity,
583                                                 Double_t& value, Double_t& error)
584 {
585   /// Get the mean pile-up correction factor for the given run
586
587   value = error = 0.0;
588
589   if (purity<=0.0)
590   {
591     AliError(Form("Cannot work with purity=%f for trigger %s in run %d. Should be strictly positive",purity,triggerClassName,runNumber));
592     return;             
593   }
594   
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());
599   
600   time_t duration = TMath::Nint((end->GetTimeStamp()->GetBunchCross() - begin->GetTimeStamp()->GetBunchCross())*AliTimeStamp::fNanosecPerBC*1E-9);
601   
602   if (!duration)
603   {
604     AliError(Form("Got zero duration for run %d",runNumber));
605     return;
606   }
607   
608   AliAnalysisTriggerScalerItem* item = GetTriggerScaler(runNumber,"L0B",triggerClassName);
609   if (!item)
610   {
611     AliError(Form("Could not get L0B for trigger %s in run %d",triggerClassName,runNumber));
612     return;
613   }
614   
615   AliLHCData* lhc = static_cast<AliLHCData*>(GetOCDBObject("GRP/GRP/LHCData",runNumber));
616   
617   Bool_t mainSat = IsMainSatelliteCollision(triggerClassName);
618   
619   Int_t nbcx = NumberOfInteractingBunches(*lhc,runNumber,mainSat);
620   
621   if ( nbcx<=0.0 )
622   {
623     AliError(Form("Cannot work with nbcx=%d for trigger %s in run %d. Should be strictly positive",nbcx,triggerClassName,runNumber));
624     return;
625   }
626   
627   Double_t itemValue = purity*item->Value();
628   
629   if (itemValue<=0.0)
630   {
631     AliError(Form("Cannot work with value=%f for trigger %s in run %d. Should be strictly positive",itemValue,triggerClassName,runNumber));
632     return;
633   }
634
635   Double_t mu = Mu(itemValue/duration,nbcx);
636   
637   value = mu/(1.0-TMath::Exp(-mu));
638   
639   error = 0.0; // FIXME
640 }
641
642 //______________________________________________________________________________
643 AliAnalysisTriggerScalerItem* 
644 AliAnalysisTriggerScalers::GetTriggerScaler(Int_t runNumber, const char* level, const char* triggerClassName)
645 {
646   // Get the scaler for a given trigger class and a given trigger level.
647   // Returned object must be deleted by the user.
648   
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));
652
653   if (!tc || !trs || !grp) return 0x0;
654   
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()));
660   
661   const TObjArray& trClasses = tc->GetClasses();
662   
663   const TObjArray* scalers = trs->GetScalersRecords();  
664   const AliTriggerScalersRecord* begin = (AliTriggerScalersRecord*)(scalers->First());
665   const AliTriggerScalersRecord* end = (AliTriggerScalersRecord*)(scalers->Last());
666
667   time_t duration = TMath::Nint((end->GetTimeStamp()->GetBunchCross() - begin->GetTimeStamp()->GetBunchCross())*AliTimeStamp::fNanosecPerBC*1E-9);
668   
669   AliTriggerClass* triggerClass = static_cast<AliTriggerClass*>(trClasses.FindObject(triggerClassName));
670   if (!triggerClass)
671   {
672     return 0x0;
673   }
674   UChar_t index = GetIndex(triggerClass->GetMask());
675   
676   ULong64_t value(0);
677   
678   const AliTriggerScalers* sbegin = begin->GetTriggerScalersForClass(index);
679   const AliTriggerScalers* send = end->GetTriggerScalersForClass(index);
680   
681   TString swhat(level);
682   swhat.ToUpper();
683   
684   if ( swhat.BeginsWith("L0A") )
685   {
686     value = send->GetLOCA() - sbegin->GetLOCA();
687   }
688   else if ( swhat.BeginsWith("L0B") )
689   {
690     value = send->GetLOCB() - sbegin->GetLOCB();
691   }
692   else if ( swhat.BeginsWith("L1A") )
693   {
694     value = send->GetL1CA() - sbegin->GetL1CA();
695   }
696   else if ( swhat.BeginsWith("L1B") )
697   {
698     value = send->GetL1CB() - sbegin->GetL1CB();
699   }
700   else if ( swhat.BeginsWith("L2A") )
701   {
702     value = send->GetL2CA() - sbegin->GetL2CA();
703   }
704   else if ( swhat.BeginsWith("L2B") )
705   {
706     value = send->GetL2CB() - sbegin->GetL2CB();
707   }
708   Int_t ds(1);
709   
710   // FIXME: get the downscaling factor here for all cases (only BC1 supported here so far)
711   if ( TString(triggerClassName).Contains("_B1") )
712   {
713     ds = GetTriggerInput(runNumber,"BC1");
714     if (!ds) ds=1;
715   }
716   
717   swhat.ReplaceAll("RATE","");
718   
719   return new AliAnalysisTriggerScalerItem(runNumber,swhat.Data(),diCurrent.Data(),triggerClass->GetName(),value,triggerClass->GetBCMask(),ds,duration);
720 }
721
722 //______________________________________________________________________________
723 TGraph*
724 AliAnalysisTriggerScalers::IntegratedLuminosityGraph(const char* triggerClassName, const char* triggerClassNameForPACEstimate)
725 {
726   // Get the integrated luminosity graph for one trigger
727   
728   if ( fRunList.empty() )
729   {
730     AliError("Cannot work without a run list");
731     return 0x0;
732   }
733   
734   Double_t offset(0.0);
735   
736   std::vector<double> vx;
737   std::vector<double> vy;
738   
739   Double_t x,y;
740   
741   for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
742   {
743     TGraph* g = IntegratedLuminosityGraph(fRunList[i],triggerClassName,triggerClassNameForPACEstimate);
744     
745     if (!g) continue;
746     
747     for ( Int_t j = 0 ; j < g->GetN(); ++j )
748     {
749       g->GetPoint(j,x,y);
750       vx.push_back(x);
751       vy.push_back(y+offset);
752     }
753     
754     offset += y;
755     
756     delete g;
757   }
758   
759   TGraph* g = new TGraph(vx.size(),&vx[0],&vy[0]);
760   
761   TimeAxis(g);
762   
763   return g;
764 }
765
766 //______________________________________________________________________________
767 TGraph*
768 AliAnalysisTriggerScalers::IntegratedLuminosityGraph(Int_t runNumber, const char* triggerClassName, const char* triggerClassNameForPACEstimate)
769 {
770   // Get the integrated luminosity graph for one trigger for the given run
771   //
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)
775   //
776   
777   AliTriggerConfiguration* tc(0x0);
778   AliTriggerRunScalers* trs(0x0);
779   AliLHCData* lhc(0x0);
780   
781   GetCTPObjects(runNumber,tc,trs,lhc);
782   
783   if (!tc || !trs || !lhc)
784   {
785     return 0x0;
786   }
787   
788   const TObjArray* scalerRecords= trs->GetScalersRecords();
789   const TObjArray& triggerClasses = tc->GetClasses();
790   
791   TString lumiTriggerClassName;
792   Double_t lumiTriggerCrossSection(0.0);
793   Double_t lumiTriggerCrossSectionError(0.0);
794   
795   GetLuminosityTriggerAndCrossSection(runNumber,lumiTriggerClassName,lumiTriggerCrossSection,lumiTriggerCrossSectionError);
796
797   AliTriggerClass* lumiTriggerClass = static_cast<AliTriggerClass*>(triggerClasses.FindObject(lumiTriggerClassName));
798
799   if (!lumiTriggerClass)
800   {
801     AliError(Form("Could not find lumiTriggerClass=%s",lumiTriggerClassName.Data()));
802     return 0x0;
803   }
804   
805   if (lumiTriggerCrossSection==0.0)
806   {
807     AliError(Form("Could not get cross-section for lumiTriggerClass=%s",lumiTriggerClassName.Data()));
808     return 0x0;
809   }
810   
811   AliTriggerClass* triggerClass = static_cast<AliTriggerClass*>(triggerClasses.FindObject(triggerClassName));
812
813   if (!triggerClass)
814   {
815     AliError(Form("Could not find triggerClass=%s",triggerClassName));
816     return 0x0;
817   }
818
819   Int_t nbcx = NumberOfInteractingBunches(*lhc,runNumber);
820   
821   if (nbcx <= 0 && ShouldCorrectForPileUp())
822   {
823     AliError("Could not get the number of bunches, so won't be able to correct for pile-up");
824   }
825   
826   TIter nextScaler(scalerRecords);
827   AliTriggerScalersRecord* record;
828   UInt_t reft(0);
829   
830   UInt_t refl0b[] = { 0, 0 };
831   UInt_t refl0a[] = { 0, 0 };
832
833   UInt_t l0b[] = { 0, 0 };
834   UInt_t l0a[] = { 0, 0 };
835   
836   Bool_t first(kTRUE);
837   
838   UChar_t classIndices[2];
839   
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
844   
845   classIndices[0] = GetIndex(lumiTriggerClass->GetMask());
846   classIndices[1] = GetIndex(triggerClass->GetMask());
847
848   Bool_t sameClass = ( classIndices[0] == classIndices[1] );
849   Bool_t pacRemoval(kFALSE);
850   
851   if ( sameClass && strlen(triggerClassNameForPACEstimate) > 0 )
852   {
853     AliTriggerClass* triggerClassForPACEstimate = static_cast<AliTriggerClass*>(triggerClasses.FindObject(triggerClassNameForPACEstimate));
854     if (!triggerClassForPACEstimate)
855     {
856       AliError(Form("Could not find triggerClassForPACEstimate=%s. Will not correct for PAC durations",triggerClassNameForPACEstimate));
857       return 0x0;
858     }
859     classIndices[1] = GetIndex(triggerClassForPACEstimate->GetMask());
860     sameClass = ( classIndices[0] == classIndices[1] );
861     if (!sameClass)
862     {
863       pacRemoval=kTRUE;
864     }
865   }
866   
867   ULong64_t counter(0);
868   
869   std::vector<Double_t> vseconds;
870   std::vector<Double_t> vlivetime;
871   std::vector<Double_t> vlumi;
872   
873   while ( ( record = static_cast<AliTriggerScalersRecord*>(nextScaler()) ) )
874   {
875     const AliTimeStamp* ats = record->GetTimeStamp();
876
877     UInt_t seconds = ats->GetSeconds();// - TTimeStamp::GetZoneOffset();
878     
879     TTimeStamp ts(seconds,ats->GetMicroSecs());
880
881     UInt_t timelapse = seconds - reft;
882
883     Bool_t ok = sameClass || CheckRecord(*record,classIndices[1],refl0b[1],refl0a[1],timelapse);
884     
885     if (ok) reft = seconds;
886     
887     for ( Int_t i = 0; i < 2; ++i )
888     {
889       const AliTriggerScalers* scaler = record->GetTriggerScalersForClass(classIndices[i]);
890       
891       if (ok)
892       {
893         l0b[i] = scaler->GetLOCB() - refl0b[i];
894       }
895       else
896       {
897         l0b[i] = 0;
898       }
899       l0a[i] = scaler->GetLOCA() - refl0a[i];
900
901       refl0b[i] = scaler->GetLOCB();
902       refl0a[i] = scaler->GetLOCA();
903     }
904     
905     if ( first )
906     {
907       first = kFALSE;
908       continue;
909     }
910     
911     Double_t factor(1.0);
912     
913     if ( ShouldCorrectForPileUp() && nbcx > 0 && l0b[0] > 0 )
914     {
915       Double_t mu = Mu(l0b[0]/timelapse,nbcx);
916       factor = mu/(1-TMath::Exp(-mu));
917     }
918     
919     vseconds.push_back(seconds);
920     
921     Double_t lt(1.0);
922     
923     if ( ok && !sameClass && !pacRemoval )
924     {
925       lt = (l0a[1]*1.0)/l0b[1];
926     }
927
928     counter += TMath::Nint(factor*l0b[0]*lt);
929
930     vlumi.push_back( counter / lumiTriggerCrossSection );
931     
932   }
933   
934   if ( vseconds.empty() ) return 0x0;
935   
936   TGraph* g = new TGraph(vseconds.size(),&vseconds[0],&vlumi[0]);
937   
938   TimeAxis(g);
939
940   return g;
941 }
942
943 //______________________________________________________________________________
944 void AliAnalysisTriggerScalers::IntegratedLuminosity(const char* triggerList,
945                                                      const char* lumiTrigger,
946                                                      Double_t lumiCrossSection,
947                                                      const char* csvOutputFile,
948                                                      const char sep,
949                                                      const char* csUnit)
950 {
951   // Compute the luminosity for a set of triggers
952   
953   // for T0 based lumi (end of pp 2012), use lumiTrigger = C0TVX-S-NOPF-ALLNOTRD and lumiCrossSection = 28 mb (and csUnit="nb")
954   //
955   // for T0 based lumi (pPb 2013), use lumiTrigger = C0TVX-B-NOPF-ALLNOTRD
956   // and lumiCrossSection = 0.755*2000 mb (and csUnit="mub")
957   
958   // for T0 based lumi (pp 2.76 TeV 2013), use lumiTrigger = C0TVX-B-NOPF-ALLNOTRD
959   // and lumiCrossSection = 18 mb (and csUnit="nb")
960   
961   double intLumi(0);
962   
963   std::map<std::string,float> lumiPerTrigger;
964   std::map<int, std::map<std::string,float> > lumiPerFillPerTrigger;
965
966   std::map<std::string,time_t> durationPerTrigger;
967   
968   TString sTriggerList(triggerList);
969   
970   if ( sTriggerList.Length()==0 )
971   {
972     if ( lumiCrossSection < 30 && lumiCrossSection  > 20 )
973     {
974       // assuming it's pp 2012
975       sTriggerList = "CMUL8-S-NOPF-ALLNOTRD,CMUL7-S-NOPF-ALLNOTRD,CMUL8-S-NOPF-MUON,CMUL7-S-NOPF-MUON";
976     
977     // for C0MUL-SC-NOPF-MUON must use C0TVX-SC as lumiTrigger (with same cross-section as C0TVX=28mb)
978     }
979     else if ( lumiCrossSection < 30 )
980     {
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";
983     }
984     else
985     {
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";
987     }
988   }
989   
990   TObjArray* triggerArray = sTriggerList.Tokenize(",");
991   triggerArray->SetOwner(kTRUE);
992
993   std::ofstream* out(0x0);
994
995   if ( TString(csvOutputFile).Length() > 0 )
996   {
997     out = new std::ofstream(gSystem->ExpandPathName(csvOutputFile));
998     if (!out || out->bad())
999     {
1000       delete out;
1001       out = 0x0;
1002     }
1003   }
1004   
1005   TIter nextTrigger(triggerArray);
1006   TObjString* trigger;
1007   
1008   if (out)
1009   {
1010     (*out) << "Fill number" << sep;
1011     
1012     nextTrigger.Reset();
1013     
1014     std::map<std::string,float>::const_iterator it;
1015
1016     while ( ( trigger = static_cast<TObjString*>(nextTrigger())))
1017     {
1018       lumiPerTrigger[trigger->String().Data()] = 0;
1019     }
1020     
1021     for ( it = lumiPerTrigger.begin(); it != lumiPerTrigger.end(); ++it )
1022     {
1023       (*out) << "lumi from " << it->first.c_str() << sep;
1024     }
1025
1026     (*out) << "comments" << sep;
1027     
1028     for ( it = lumiPerTrigger.begin(); it != lumiPerTrigger.end(); ++it )
1029     {
1030       (*out) << "recorded " <<  it->first.c_str() << " (integrated)" << sep;
1031     }
1032
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;
1037     
1038     nextTrigger.Reset();
1039   }
1040   
1041   Int_t currentFillNumber(-1);
1042   Int_t fillNumber(0);
1043   
1044   float csMult(1.0);
1045   TString scsUnit(csUnit);
1046   scsUnit.ToUpper();
1047   
1048   if ( scsUnit=="PB" )
1049   {
1050     csMult=1E9;
1051   }
1052   else if (scsUnit=="NB")
1053   {
1054     csMult=1E6;
1055   }
1056   else if ( scsUnit=="UB" )
1057   {
1058     csMult=1E3;
1059   }
1060   else if ( scsUnit=="MB" )
1061   {
1062     csMult=1.0;
1063   }
1064   else
1065   {
1066     AliError(Form("Don't know how to deal with cross-section unit=%s",csUnit));
1067     return;
1068   }
1069   
1070
1071   for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
1072   {
1073     Int_t runNumber = fRunList[i];
1074     Bool_t atLeastOneTriggerFound(kFALSE);
1075     
1076     // find out which trigger classes to use
1077     
1078     AliTriggerConfiguration* tc = static_cast<AliTriggerConfiguration*>(GetOCDBObject("GRP/CTP/Config",runNumber));
1079     const TObjArray& trClasses = tc->GetClasses();
1080         
1081     if (out)
1082     {
1083       fillNumber = GetFillNumberFromRunNumber(runNumber);
1084       if ( fillNumber == 0 )
1085       {
1086         AliError(Form("Got fillNumber = 0 for run %09d !",runNumber));
1087       }
1088       
1089       if ( fillNumber != currentFillNumber )
1090       {
1091         std::map<std::string,float>::const_iterator it;
1092         
1093         for ( it = lumiPerTrigger.begin(); it != lumiPerTrigger.end(); ++it )
1094         {
1095           lumiPerFillPerTrigger[fillNumber][it->first.c_str()] = 0;
1096         }
1097         currentFillNumber = fillNumber;
1098       }
1099     }
1100
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());
1105     
1106     time_t duration = TMath::Nint((end->GetTimeStamp()->GetBunchCross() - begin->GetTimeStamp()->GetBunchCross())*AliTimeStamp::fNanosecPerBC*1E-9);
1107
1108     nextTrigger.Reset();
1109
1110     TString lumiTriggerClassName(lumiTrigger);
1111     float lumiSigma = lumiCrossSection*csMult; //in csUnit
1112     AliAnalysisTriggerScalerItem* lumiB = GetTriggerScaler(runNumber,"L0B",lumiTriggerClassName.Data());
1113
1114     if (!lumiB)
1115     {
1116       AliError(Form("Did not find lumiTrigger %s for run %09d",lumiTriggerClassName.Data(),runNumber));
1117       continue;
1118     }
1119         
1120     Float_t pacCorrection(1.0);
1121     
1122     while ( ( trigger = static_cast<TObjString*>(nextTrigger()) ) )
1123     {
1124       TString muTriggerClassName(trigger->String());
1125       Int_t n(0);
1126       
1127       if ( !trClasses.FindObject(muTriggerClassName.Data() ) )
1128       {
1129         continue;
1130       }
1131       
1132       if ( muTriggerClassName.Contains("CMUL8") ) ++n;
1133       if ( muTriggerClassName.Contains("CMUL7") ) ++n;
1134       
1135       if ( n>1 )
1136       {
1137         AliError(Form("More than 1 relevant trigger class found for run %09d ! Check that !",runNumber));
1138         trClasses.Print();
1139         continue;
1140       }
1141       
1142       AliAnalysisTriggerScalerItem* muonA = GetTriggerScaler(runNumber,"L0A",muTriggerClassName.Data());
1143       AliAnalysisTriggerScalerItem* muonB = GetTriggerScaler(runNumber,"L0B",muTriggerClassName.Data());
1144       
1145       if (!muonA || !muonB) continue;
1146       
1147       atLeastOneTriggerFound = kTRUE;
1148       
1149       Float_t ratio(0);
1150       
1151       if ( muonB->ValueCorrectedForDownscale() > 0 )
1152       {
1153         ratio = muonA->ValueCorrectedForDownscale()/muonB->ValueCorrectedForDownscale();
1154       }
1155       
1156       ratio *= lumiB->ValueCorrectedForDownscale()/lumiSigma;
1157       
1158       if ( muTriggerClassName.BeginsWith("CMUL"))
1159       {
1160         intLumi += ratio;
1161       }
1162       
1163       if ( muTriggerClassName.Contains("CMSH") )
1164       {
1165         pacCorrection = GetPauseAndConfigCorrection(runNumber,muTriggerClassName.Data());
1166       }
1167       
1168       lumiPerTrigger[muTriggerClassName.Data()] += ratio;
1169       durationPerTrigger[muTriggerClassName.Data()] += duration;
1170       lumiPerFillPerTrigger[currentFillNumber][muTriggerClassName.Data()] += ratio;
1171       
1172     }
1173     
1174     lumiPerTrigger[lumiTriggerClassName.Data()] += lumiB->ValueCorrectedForDownscale()/lumiSigma;
1175     durationPerTrigger[lumiTriggerClassName.Data()] += duration;
1176     lumiPerFillPerTrigger[currentFillNumber][lumiTriggerClassName.Data()] += lumiB->ValueCorrectedForDownscale()/lumiSigma;
1177
1178     TString lumiPACCorrected(Form("%s(-PAC)",lumiTriggerClassName.Data()));
1179     
1180     lumiPerTrigger[lumiPACCorrected.Data()] += pacCorrection*lumiB->ValueCorrectedForDownscale()/lumiSigma;
1181     durationPerTrigger[lumiPACCorrected.Data()] += duration;
1182     lumiPerFillPerTrigger[currentFillNumber][lumiPACCorrected.Data()] += pacCorrection*lumiB->ValueCorrectedForDownscale()/lumiSigma;
1183
1184     if (!atLeastOneTriggerFound && sTriggerList.Contains("CMUL") )
1185     {
1186       AliError(Form("Found no relevant trigger for run %09d",runNumber));
1187     }
1188   }
1189   
1190   AliInfo(Form("Integrated lumi %7.4f %s^-1",intLumi,csUnit));
1191   
1192   std::map<std::string,float>::const_iterator it;
1193   
1194   for ( it = lumiPerTrigger.begin(); it != lumiPerTrigger.end(); ++it )
1195   {
1196     AliInfo(Form("Trigger %30s Lumi %10.4f %s^-1 duration %10ld s",it->first.c_str(),it->second,csUnit,durationPerTrigger[it->first]));
1197   }
1198   
1199   if (out)
1200   {
1201     std::map<int, std::map<std::string, float> >::const_iterator fit;
1202     
1203     lumiPerTrigger.clear();
1204     
1205     for ( fit = lumiPerFillPerTrigger.begin(); fit != lumiPerFillPerTrigger.end(); ++fit )
1206     {
1207       int fill = fit->first;
1208       std::map<std::string,float>::const_iterator tit;
1209       (*out) << fill << sep;
1210       
1211       for ( tit = fit->second.begin(); tit != fit->second.end(); ++tit )
1212       {
1213         (*out) << Form("%e",tit->second) << sep;
1214       }
1215       
1216       (*out) << sep; // comment (empty)
1217       
1218       for ( tit = fit->second.begin(); tit != fit->second.end(); ++tit )
1219       {
1220         lumiPerTrigger[tit->first] += tit->second;        
1221         
1222         (*out) <<  Form("%e",lumiPerTrigger[tit->first]) << sep;
1223       }
1224       (*out) << sep << "0" << sep << "0" << sep << "0" << sep << "0" << sep << std::endl; // LHC per fill, LHC integrated, lumi tot muon , efficiency
1225     }
1226   }
1227   //
1228   delete triggerArray;
1229 }
1230
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)
1236 {
1237   /// Build a graph from a set of stl vectors
1238   
1239   if ( ! ( vx.size() == vex.size() && vx.size() == vy.size() && vx.size() == vey.size() ) )
1240   {
1241     std::cerr << "incompatible sizes" << std::endl;
1242     return 0x0;
1243   }
1244   
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()];
1249   
1250   for ( size_t i = 0; i < vx.size(); ++i )
1251   {
1252     x[i] = vx[i];
1253     ex[i] = vex[i];
1254     y[i] = vy[i];
1255     ey[i] = vey[i];
1256   }
1257   
1258   TGraph* g = new TGraphErrors(vx.size(),x,y,ex,ey);
1259   
1260   TimeAxis(g);
1261   
1262   delete[] x;
1263   delete[] y;
1264   delete[] ex;
1265   delete[] ey;
1266   
1267   return g;
1268 }
1269
1270 //______________________________________________________________________________
1271 Double_t AliAnalysisTriggerScalers::Mu(Double_t L0B, Double_t Nb)
1272 {
1273   /// L0B = trigger rate before any veto
1274   /// Nb = number of crossing bunches
1275   
1276   Double_t p0 = 1-L0B/(Nb*11245.0); // proba to get *no* collision per bunch crossing
1277   
1278   return -TMath::Log(p0);
1279 }
1280
1281 //______________________________________________________________________________
1282 Int_t AliAnalysisTriggerScalers::NumberOfInteractingBunches(const AliLHCData& lhc, Int_t runNumber, Bool_t mainSat) const
1283 {
1284   /// Extract the number of colliding bunches from the LHC data
1285   /// Use mainSat = true if the collisons were main - satellite ones
1286   
1287   Int_t numberOfInteractingBunches(0);
1288   Int_t numberOfInteractingBunchesMeasured(0);
1289   Int_t nIBM2(0);
1290
1291   int beam1(0);
1292   int beam2(1);
1293
1294   AliLHCDipValI* val = lhc.GetBunchConfigDeclared(beam1,0);
1295   AliLHCDipValI* valm = lhc.GetBunchConfigMeasured(beam1,0);
1296
1297   
1298   if ( mainSat )
1299   {
1300     numberOfInteractingBunches = val->GetSizeTotal();
1301     numberOfInteractingBunchesMeasured = valm->GetSizeTotal();
1302     nIBM2 = numberOfInteractingBunches;
1303   }
1304   else
1305   {
1306     for ( Int_t i = 0; i < val->GetSizeTotal(); ++i )
1307     {
1308       if ( val->GetValue(i) < 0 ) ++numberOfInteractingBunches;
1309     }
1310
1311     for ( Int_t i = 0; i < valm->GetSizeTotal(); ++i )
1312     {
1313       if ( valm->GetValue(i) < 0 ) ++numberOfInteractingBunchesMeasured;
1314     }
1315
1316     valm = lhc.GetBunchConfigMeasured(beam2,0);
1317
1318     for ( Int_t i = 0; i < valm->GetSizeTotal(); ++i )
1319     {
1320       if ( valm->GetValue(i) <= 0 ) ++nIBM2;
1321     }
1322   }
1323   
1324   if (!numberOfInteractingBunches)
1325   {
1326     return 0;
1327   }
1328
1329   if ( numberOfInteractingBunches != numberOfInteractingBunchesMeasured ||
1330        numberOfInteractingBunches != nIBM2 )
1331   {
1332     Int_t delta = TMath::Max(numberOfInteractingBunches - numberOfInteractingBunchesMeasured,numberOfInteractingBunches-nIBM2);
1333     
1334     if ( 1.0*TMath::Abs(delta)/numberOfInteractingBunches > 0.05 ) // more than 5% difference
1335     {
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));
1339     }
1340     else
1341     {
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));
1345       
1346     }
1347   }
1348   
1349   return numberOfInteractingBunches;
1350 }
1351
1352
1353 //______________________________________________________________________________
1354 TGraph* AliAnalysisTriggerScalers::PlotTrigger(const char* triggerClassName,
1355                                                const char* what)
1356 {
1357   // Plot one of the scalers (L0A,L0B,L0AOVERB,etc...) of a given triggerClass
1358   // Get one value per run.
1359   
1360   std::vector<Float_t> x;
1361   std::vector<Float_t> y;
1362
1363   double integral(0);
1364   double mean(0);
1365   
1366   for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
1367   {
1368     Int_t runNumber = fRunList[i];
1369     
1370     AliAnalysisTriggerScalerItem* s = GetTriggerScaler(runNumber,what,triggerClassName);
1371     
1372     if (!s) continue;
1373     
1374     x.push_back(runNumber);
1375
1376     Double_t value = s->ValueCorrectedForDownscale();
1377     
1378     if ( TString(what).Contains("RATE") )
1379     {
1380       value = s->Rate();
1381     }
1382   
1383     integral += value;
1384     mean += value;
1385     
1386     y.push_back(value);
1387   }
1388   
1389   if ( fRunList.size() ) {
1390     mean /= fRunList.size();
1391   }
1392   
1393   AliInfo(Form("Integral %e mean %e",integral,mean));
1394   
1395   return new TGraph(x.size(),&x[0],&y[0]);
1396 }
1397
1398 //______________________________________________________________________________
1399 TGraph* AliAnalysisTriggerScalers::PlotTriggerEvolution(const char* triggerClassName,
1400                                                         const char* what,
1401                                                         bool draw,
1402                                                         double* mean,
1403                                                         bool removeZeros)
1404 {
1405   /// Make a graph of "what" for a given trigger class.
1406   ///
1407   /// What can be :
1408   ///
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
1416   
1417   TString swhat(what);
1418   swhat.ToUpper();
1419   
1420   if ( swhat.Contains(";"))
1421   {
1422     swhat.ReplaceAll(";",",");
1423     AliWarningClass("; is not a valid separator, replaced it with ,");
1424   }
1425   
1426   int color(1);
1427   int marker(20);
1428   
1429   TObjArray* a = swhat.Tokenize(",");
1430   if (a->GetEntries()>1)
1431   {
1432     TObjArray graphs;
1433     TIter next(a);
1434     TObjString* str(0x0);
1435     Double_t ymin(TMath::Limits<Double_t>::Max());
1436     Double_t ymax(0);
1437     Double_t xmin(TMath::Limits<Double_t>::Max());
1438     Double_t xmax(0);
1439     TGraph* g(0x0);
1440     
1441     while ( ( str = static_cast<TObjString*>(next())) )
1442     {
1443       g = PlotTriggerEvolution(triggerClassName,str->String().Data(),false,0x0,removeZeros);
1444       graphs.Add(g);
1445       for ( Int_t i = 0; i < g->GetN(); ++i )
1446       {
1447         ymin = TMath::Min(ymin,g->GetY()[i]);
1448         ymax = TMath::Max(ymax,g->GetY()[i]);
1449       }
1450       xmin = TMath::Min(xmin,g->GetX()[0]);
1451       xmax = TMath::Max(xmax,g->GetX()[g->GetN()-1]);
1452     }
1453     
1454     gStyle->SetOptTitle(0);
1455     
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);
1463     h->Draw();
1464     
1465     TIter nextGraph(&graphs);
1466     
1467     while ( ( g = static_cast<TGraph*>(nextGraph())) )
1468     {
1469       AliInfoClass(g->GetTitle());
1470       g->Draw("lp");
1471       g->SetLineColor(color);
1472       g->SetMarkerColor(color);
1473       g->SetMarkerStyle(marker);
1474       ++color;
1475       ++marker;
1476     }
1477     delete a;
1478     return 0x0;
1479   }
1480   delete a;
1481   
1482   std::vector<int> vx;
1483   std::vector<int> vex;
1484   std::vector<double> vy;
1485   std::vector<double> vey;
1486   
1487   
1488   if (mean) *mean=0;
1489   double nvalues(0);
1490   
1491   for ( std::vector<int>::size_type iRun = 0; iRun < fRunList.size(); ++iRun )
1492   {
1493     Int_t runNumber = fRunList[iRun];
1494     
1495     AliTriggerConfiguration* tc(0x0);
1496     AliTriggerRunScalers* trs(0x0);
1497     AliLHCData* lhc(0x0);
1498     
1499     GetCTPObjects(runNumber,tc,trs,lhc);
1500     
1501     const TObjArray* scalers = trs->GetScalersRecords();
1502     
1503     const TObjArray& trClasses = tc->GetClasses();
1504     TIter next(&trClasses);
1505     AliTriggerClass* triggerClass;
1506     
1507     while ( ( triggerClass = static_cast<AliTriggerClass*>(next()) ) )
1508     {
1509       UChar_t index = GetIndex(triggerClass->GetMask());
1510       
1511       if ( !TString(triggerClass->GetName()).Contains(triggerClassName) ) continue;
1512       
1513       Bool_t mainSat = IsMainSatelliteCollision(triggerClass->GetName());
1514       
1515       Int_t numberOfInteractingBunches = NumberOfInteractingBunches(*lhc,runNumber,mainSat);
1516       
1517       TIter nextScaler(scalers);
1518       AliTriggerScalersRecord* record;
1519       UInt_t reft(0);
1520       UInt_t refl0b(0);
1521       UInt_t refl1b(0);
1522       UInt_t refl2b(0);
1523       UInt_t refl0a(0);
1524       UInt_t refl1a(0);
1525       UInt_t refl2a(0);
1526       Bool_t first(kTRUE);
1527       
1528       while ( ( record = static_cast<AliTriggerScalersRecord*>(nextScaler()) ) )
1529       {
1530         const AliTriggerScalers* scaler = record->GetTriggerScalersForClass(index);
1531         
1532         
1533         const AliTimeStamp* ats = record->GetTimeStamp();
1534         
1535         UInt_t seconds = ats->GetSeconds();// - TTimeStamp::GetZoneOffset();
1536         
1537         TTimeStamp ts(seconds,ats->GetMicroSecs());
1538         
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;
1546         
1547         if ( removeZeros && ( l0b <= 2 || ( l0a <= 2 && l0a != 0 ) || timelapse <= 9 ) )
1548         {
1549           AliInfo(Form("Skipping point for %s l0b %d l0a %d timelapse %d ts %s",
1550                        triggerClassName,l0b,l0a,timelapse,ts.AsString()));
1551           continue;
1552         }
1553         
1554         reft = seconds;
1555         refl0b = scaler->GetLOCB();
1556         refl0a = scaler->GetLOCA();
1557         refl1b = scaler->GetL1CB();
1558         refl1a = scaler->GetL1CA();
1559         refl2b = scaler->GetL2CB();
1560         refl2a = scaler->GetL2CA();
1561         
1562         if ( first )
1563         {
1564           first = kFALSE;
1565           continue;
1566         }
1567         
1568         double value(1.0);
1569         double error(0.0);
1570         
1571         if (swhat.Contains("L0AOVERB") )
1572         {
1573           value = l0a*1.0/l0b;
1574           if ( l0a > 0 )
1575           {
1576             error = value*TMath::Sqrt(1.0/l0b+1.0/l0a);
1577           }
1578           else
1579           {
1580             error = 0.0;
1581           }
1582         }
1583         else if ( swhat.Contains("L0B") )
1584         {
1585           value = l0b;
1586         }
1587         else if (swhat.Contains("L0A") )
1588         {
1589           value = l0a;
1590         }
1591         else if ( swhat.Contains("L1B") )
1592         {
1593           value = l1b;
1594         }
1595         else if (swhat.Contains("L1A") )
1596         {
1597           value = l1a;
1598         }
1599         else if ( swhat.Contains("L2B") )
1600         {
1601           value = l2b;
1602         }
1603         else if (swhat.Contains("L2A") )
1604         {
1605           value = l2a;
1606         }
1607         else if ( swhat.Contains("MU") )
1608         {
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
1613         }
1614         else if ( swhat.Contains("PILEUPFACTOR") )
1615         {
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
1620         }
1621         else if ( swhat.Contains("VSNB") )
1622         {
1623           value = l0b/(11245.0*numberOfInteractingBunches);
1624           error = -1.0; // FIXME
1625         }
1626         else if ( swhat.Contains("NBCX"))
1627         {
1628           value = numberOfInteractingBunches;
1629           error = 1.0; // FIXME          
1630         }
1631         else
1632         {
1633           value = timelapse;
1634           AliInfo(Form("Unknown what %s",what));
1635         }
1636         
1637         if ( ! swhat.Contains("OVER") && ! swhat.Contains("RATIO") &&
1638             ! swhat.Contains("MU") && ! swhat.Contains("PILEUPFACTOR") &&
1639             ! swhat.Contains("RAW") & ! swhat.Contains("NBCX") )
1640         {
1641           value /= timelapse;
1642         }
1643         
1644         if ( !TMath::Finite(value) ) continue;
1645         
1646         if ( !swhat.Contains("L0AOVERB") && error > 0 )
1647         {
1648           error = value >0 ? 1.0/TMath::Sqrt(value) : 0.0;
1649         }
1650         
1651         if (removeZeros && TMath::Abs(value) < 1E-6 )
1652         {
1653           continue;
1654         }
1655         if (mean)
1656         {
1657           (*mean) += value;
1658           nvalues += 1.0;
1659         }
1660         
1661         vx.push_back(seconds);
1662         vex.push_back(1);
1663         
1664         vy.push_back(value);
1665         
1666         vey.push_back(error);
1667         
1668       }
1669       
1670     }
1671     
1672   }
1673   
1674   if (mean && nvalues)
1675   {
1676     (*mean) /= nvalues;
1677   }
1678   
1679   if ( vx.empty() ) return 0;
1680   
1681   TGraph* g = MakeGraph(vx,vex,vy,vey);
1682   TString title(Form("TriggerEvolution%s-%s",triggerClassName,swhat.Data()));
1683   
1684   g->SetName(title.Data());
1685   g->SetTitle(title.Data());
1686   
1687   g->GetYaxis()->SetTitle(title.Data());
1688   
1689   if (draw)
1690   {
1691     TCanvas* c = NewCanvas(g->GetName());
1692     g->SetLineColor(color);
1693     g->SetMarkerColor(color);
1694     g->SetMarkerStyle(marker);
1695     g->Draw("ALP");
1696     c->SaveAs(Form("%s.png",title.Data()));
1697   }
1698   
1699   
1700   return g;
1701 }
1702
1703 //______________________________________________________________________________
1704 TGraph* AliAnalysisTriggerScalers::PlotTriggerRatio(const char* triggerClassName1,
1705                                                     const char* what1,
1706                                                     const char* triggerClassName2,
1707                                                     const char* what2)
1708 {
1709   // Plot the ratio of two scalers.
1710   // Get one value per run.
1711
1712   std::vector<Float_t> x;
1713   std::vector<Float_t> y;
1714   
1715   for ( std::vector<int>::size_type i = 0; i < fRunList.size(); ++i )
1716   {
1717     Int_t runNumber = fRunList[i];
1718     
1719     AliAnalysisTriggerScalerItem* s1 = GetTriggerScaler(runNumber,what1,triggerClassName1);
1720     AliAnalysisTriggerScalerItem* s2 = GetTriggerScaler(runNumber,what2,triggerClassName2);
1721     
1722     if (!s1 || !s2) continue;
1723     
1724     x.push_back(runNumber);
1725     
1726     Float_t ratio(0);
1727     
1728     if ( s2->ValueCorrectedForDownscale() > 0 )
1729     {
1730       ratio = s1->ValueCorrectedForDownscale()/s2->ValueCorrectedForDownscale();
1731     }
1732     
1733     y.push_back(ratio);
1734     
1735       AliDebug(1,Form("RUN %09d %20s (%s) %12llu (%5d) %20s (%s) %12llu (%5d) R %7.2f",
1736                    runNumber,
1737                    triggerClassName1,what1,
1738                    s1->Value(),s1->DownscalingFactor(),
1739                    triggerClassName2,what2,
1740                    s2->Value(),s2->DownscalingFactor(),
1741                    ratio));
1742   }
1743   
1744   return new TGraph(x.size(),&x[0],&y[0]);
1745 }
1746
1747 //______________________________________________________________________________
1748 TGraph* AliAnalysisTriggerScalers::PlotTriggerRatioEvolution(const char* triggerClassName1,
1749                                                              const char* what1,
1750                                                              const char* triggerClassName2,
1751                                                              const char* what2)
1752 {
1753   /// Plots the evolution of one trigger ratio
1754   
1755   Bool_t draw(kFALSE);
1756   Bool_t removeZeros(kFALSE);
1757   
1758   TGraph* g1 = PlotTriggerEvolution(triggerClassName1,what1,draw,0x0,removeZeros);
1759   TGraph* g2 = PlotTriggerEvolution(triggerClassName2,what2,draw,0x0,removeZeros);
1760   
1761   if (!g1 || !g2) return 0x0;
1762   
1763   if ( g1->GetN() != g2->GetN() )
1764   {
1765     AliError("Oups. Did not get the same number of points for the 2 graphs ?!");
1766     return 0x0;
1767   }
1768
1769   Double_t x1,x2;
1770   Double_t y1,y2;
1771   Double_t x1err,x2err;
1772   Double_t y1err,y2err;
1773   
1774   TGraphErrors* g = new TGraphErrors(g1->GetN());
1775   Int_t j(0);
1776   
1777   for ( Int_t i = 0; i < g1->GetN(); ++i )
1778   {
1779     g1->GetPoint(i,x1,y1);
1780     g2->GetPoint(i,x2,y2);
1781     
1782     x1err = g1->GetErrorX(i);
1783     x2err = g2->GetErrorX(i);
1784
1785     y1err = g1->GetErrorY(i);
1786     y2err = g2->GetErrorY(i);
1787     
1788     if  (x1!=x2)
1789     {
1790       AliError(Form("Points at %d don't have the same x ! ",i));
1791       continue;
1792     }
1793     
1794     Double_t y(0.0);
1795     
1796     if ( TMath::Abs(y2) > 1E-12 )
1797     {
1798       y = y1/y2;
1799     }
1800     
1801     Double_t yerr(0.0);
1802     
1803     if ( TMath::Abs(x1) > 1E-12)
1804     {
1805       yerr += TMath::Sqrt(x1err*x1err/(x1*x1));
1806     }
1807
1808     if ( TMath::Abs(x2) > 1E-12)
1809     {
1810       yerr += TMath::Sqrt(x2err*x2err/(x2*x2));
1811     }
1812
1813     yerr *= y;
1814       
1815     g->SetPoint(j,x1,y);
1816     g->SetPointError(j,x1err,yerr);
1817     
1818     ++j;
1819   }
1820   
1821   delete g1;
1822   delete g2;
1823   
1824   TimeAxis(g);
1825   
1826   return g;
1827 }
1828
1829 //______________________________________________________________________________
1830 void AliAnalysisTriggerScalers::Print(Option_t* /* opt */) const
1831 {
1832   /// print our runlist
1833   AliAnalysisTriggerScalers::PrintIntegers(fRunList,',');
1834 }
1835
1836 //______________________________________________________________________________
1837 void AliAnalysisTriggerScalers::PrintIntegers(const std::vector<int>& integers,
1838                                               char sep,
1839                                               std::ostream& out)
1840 {
1841   /// print a list of integers
1842   for ( std::vector<int>::size_type i = 0; i < integers.size(); ++i )
1843   {
1844     out << integers[i] << sep;
1845   }
1846   out << std::endl;
1847 }
1848
1849 //______________________________________________________________________________
1850 void AliAnalysisTriggerScalers::ReadIntegers(const char* filename,
1851                                              std::vector<int>& integers,
1852                                              Bool_t resetVector)
1853 {
1854   /// Read integers from filename, where integers are either
1855   /// separated by "," or by return carriage
1856   
1857   if ( gSystem->AccessPathName(filename)==kTRUE ) 
1858   {
1859     return;
1860   }
1861   std::ifstream in(gSystem->ExpandPathName(filename));
1862   int i;
1863   
1864   std::set<int> runset;
1865   
1866   if (!resetVector)
1867   {
1868     for ( std::vector<int>::size_type j = 0; j < integers.size(); ++ j )
1869     {
1870       runset.insert(integers[j]);
1871     }
1872   }
1873   
1874   char line[10000];
1875   
1876   in.getline(line,10000,'\n');
1877   
1878   TString sline(line);
1879   
1880   if (sline.Contains(","))
1881   {
1882     TObjArray* a = sline.Tokenize(",");
1883     TIter next(a);
1884     TObjString* s;
1885     while ( ( s = static_cast<TObjString*>(next()) ) )
1886     {
1887       runset.insert(s->String().Atoi());
1888     }
1889     delete a;
1890   }
1891   else
1892   {
1893     runset.insert(sline.Atoi());
1894     
1895     while ( in >> i )
1896     {
1897       runset.insert(i);
1898     }
1899   }
1900   
1901   integers.clear();
1902   
1903   for ( std::set<int>::const_iterator it = runset.begin(); it != runset.end(); ++it ) 
1904   {
1905     integers.push_back((*it)); 
1906   }
1907   
1908   std::sort(integers.begin(),integers.end());
1909 }
1910
1911
1912 //______________________________________________________________________________
1913 void AliAnalysisTriggerScalers::SetRunList(Int_t runNumber)
1914 {
1915   // Make the runlist be a single run
1916   fRunList.clear();
1917   fRunList.push_back(runNumber);
1918 }
1919
1920 //______________________________________________________________________________
1921 void AliAnalysisTriggerScalers::SetRunList(const std::vector<int>& runs)
1922 {
1923   // Make the runlist be a single run
1924   fRunList = runs;
1925 }
1926
1927 //______________________________________________________________________________
1928 void AliAnalysisTriggerScalers::SetRunList(const std::set<int>& runs)
1929 {
1930   // Make the runlist be a single run
1931   fRunList.clear();
1932   for ( std::set<int>::const_iterator it = runs.begin(); it != runs.end(); ++it )
1933   {
1934     fRunList.push_back(*it);
1935   }
1936 }
1937
1938 //______________________________________________________________________________
1939 void AliAnalysisTriggerScalers::SetRunList(const char* runlist)
1940 {
1941   // Read the runlist from an ASCII file or a comma separated list
1942   // or a space separated list
1943   
1944   fRunList.clear();
1945   
1946   if ( TString(runlist).Contains(",") || TString(runlist).Contains(" ") )
1947   {
1948     TObjArray* runs = 0x0;
1949     if ( TString(runlist).Contains(",") ) 
1950     {
1951       runs = TString(runlist).Tokenize(",");
1952     }
1953     else
1954     {
1955       runs = TString(runlist).Tokenize(" ");
1956     }
1957     TIter next(runs);
1958     TObjString* s;
1959     std::set<int> runset;
1960     
1961     while ( ( s = static_cast<TObjString*>(next()) ) )
1962     {
1963       runset.insert(s->String().Atoi());
1964     }
1965
1966     for ( std::set<int>::const_iterator it = runset.begin(); it != runset.end(); ++it ) 
1967     {
1968       fRunList.push_back((*it)); 
1969     }
1970     
1971     std::sort(fRunList.begin(),fRunList.end());
1972
1973     delete runs;
1974   }
1975   else
1976   {
1977     ReadIntegers(runlist,fRunList);
1978     if (fRunList.empty())
1979     {
1980       if ( TString(runlist).IsDigit() )
1981       {
1982         SetRunList(TString(runlist).Atoi());
1983       }
1984       else
1985       {
1986         AliError("Could not set run list !");
1987       }
1988     }
1989   }
1990 }
1991
1992 //______________________________________________________________________________
1993 //______________________________________________________________________________
1994 //______________________________________________________________________________
1995 //______________________________________________________________________________
1996 //______________________________________________________________________________
1997 //______________________________________________________________________________
1998 //______________________________________________________________________________
1999
2000 ClassImp(AliAnalysisTriggerScalerItem)
2001
2002 //______________________________________________________________________________
2003 const char* 
2004 AliAnalysisTriggerScalerItem::BCMaskName() const 
2005 {
2006   // return bc mask name
2007   if ( BCMask() ) return BCMask()->GetName(); else return ""; 
2008 }
2009
2010 //______________________________________________________________________________
2011 Int_t AliAnalysisTriggerScalerItem::Compare(const TObject* obj) const
2012 {
2013   /// compare two scaler items (by means of their run number only)
2014   const AliAnalysisTriggerScalerItem* s = static_cast<const AliAnalysisTriggerScalerItem*>(obj);
2015   
2016   if ( s->RunNumber() < RunNumber() ) 
2017   {
2018     return -1;
2019   }
2020   else if ( s->RunNumber() > RunNumber() )
2021   {
2022     return 1;
2023   }
2024   return 0;
2025 }
2026
2027 //______________________________________________________________________________
2028 void AliAnalysisTriggerScalerItem::Print(Option_t* opt) const
2029 {
2030   /// Printout
2031   TString sopt(opt);
2032   
2033   if ( RunNumber() > 0 )
2034   {
2035     sopt.ToUpper();
2036     
2037     if ( sopt.Contains("FULL") )
2038     {
2039     }
2040     else
2041     {
2042       std::cout << Form("RUN %6d CLASS %24s (%5s %4d BCs) SCALER %s %12llu DS %6d DURATION %ld",
2043                    RunNumber(),TriggerClassName(),
2044                    BCMaskName(),
2045                    BCMask() ? BCMask()->GetNUnmaskedBCs() : 0,
2046                    Level(),
2047                    Value(),DownscalingFactor(),Duration()) << std::endl;
2048     }
2049   }
2050   else
2051   {
2052     std::cout << Form("CLASS %24s SCALER %20llu NRUNS %d",
2053                  TriggerClassName(),Value(),NofRuns()) << std::endl;    
2054   }
2055 }