]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/muondep/AliAnalysisTaskPileup.cxx
o small fix
[u/mrichter/AliRoot.git] / PWG / muondep / AliAnalysisTaskPileup.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 //             Implementation of class AliAnalysisTaskPileup
20 //
21 // The class re-weights the number of analysed events with the correction
22 // factor to account for the event pileup.
23 //
24 // The final idea will be to include all of the methods for the pileup
25 // calculation.
26 //
27 // For the moment, the implemented method is based on a work 
28 // by L. Aphecetche.
29 // The probability of having at least 1 collision is given by:
30 // MB L0 trigger / max collision rate
31 // The max collision rate is given by:
32 //    number of bunches X revolution frequency ~ CBEAMB L0 trigger
33 // Assuming a poissonian distribution of the collision, with:
34 //    mu = mean number of collisions
35 // the probability of having 0 collisions is:
36 //    exp{-mu} = 1 - L0b_MB / max collision rate
37 // Since the MB trigger is measuring the probability of having at least
38 // 1 collisions, the number should be corrected by the factor:
39 //    CF = mu / ( 1 - exp{-mu} )
40 // The class weights the number of events with this correction factor
41 //
42 // CAVEATS:
43 // - so far the class needs to access the OCDB
44 //   (hopefully it will be changed in the future)
45 //   Hence the following libraries need to be included (in addition to
46 //   the standard analysis ones):
47 //   libRAWDatabase.so libCDB.so libSTEER.so libPWGHFbase.so
48 //-------------------------------------------------------------------------
49
50 #include <Riostream.h>
51
52 // ROOT includes
53 #include "TH1.h"
54 #include "TH2.h"
55 #include "TCanvas.h"
56 #include "TROOT.h"
57 #include "TString.h"
58 #include "TObjArray.h"
59 #include "TObjString.h"
60 #include "TMath.h"
61
62 // STEER includes
63 #include "AliESDEvent.h"
64 #include "AliESDInputHandler.h"
65 #include "AliAODEvent.h"
66 #include "AliAODInputHandler.h"
67 #include "AliTriggerScalersESD.h"
68 #include "AliCentrality.h"
69
70 // ANALYSIS includes
71 #include "AliAnalysisTaskSE.h"
72 #include "AliAnalysisManager.h"
73 #include "AliAnalysisDataSlot.h"
74 #include "AliAnalysisDataContainer.h"
75
76 #include "AliCounterCollection.h"
77
78 #include "AliAnalysisTaskPileup.h"
79
80 #if defined(READOCDB)
81 #include "AliTriggerRunScalers.h"
82 #include "AliTriggerConfiguration.h"
83 #include "AliTriggerClass.h"
84 #include "AliCDBManager.h"
85 #include "AliCDBEntry.h"
86 #include "AliCDBPath.h"
87 #endif
88
89
90 ClassImp(AliAnalysisTaskPileup)
91
92
93 //________________________________________________________________________
94 AliAnalysisTaskPileup::AliAnalysisTaskPileup() :
95 AliAnalysisTaskSE(),
96   fEventCounters(0x0),
97   fHistoEventsList(0x0),
98   fTriggerClasses(0x0),
99   fTriggerClassIndex(0x0),
100   fIsInitCDB(0),
101   fCentralityClasses(0x0),
102   fStorageList("")
103 #if defined(READOCDB)
104   , fTriggerRunScalers(0x0)
105 #endif
106
107 {
108   /// Default Constructor
109 }
110
111
112 //________________________________________________________________________
113 AliAnalysisTaskPileup::AliAnalysisTaskPileup(const char *name) :
114   AliAnalysisTaskSE(name),
115   fEventCounters(0x0),
116   fHistoEventsList(0x0),
117   fTriggerClasses(0x0),
118   fTriggerClassIndex(0x0),
119   fIsInitCDB(0),
120   fCentralityClasses(0x0),
121   fStorageList("")
122 #if defined(READOCDB)
123   , fTriggerRunScalers(0x0)
124 #endif
125 {
126   /// Constructor
127
128   DefineOutput(1,AliCounterCollection::Class());
129   DefineOutput(2,TObjArray::Class());
130 }
131
132 //________________________________________________________________________
133 AliAnalysisTaskPileup::~AliAnalysisTaskPileup()
134 {
135   /// Destructor
136
137   // For proof: do not delete output containers
138   if ( ! AliAnalysisManager::GetAnalysisManager() || ! AliAnalysisManager::GetAnalysisManager()->IsProofMode() ) {
139     delete fEventCounters;
140   }
141
142   delete fHistoEventsList;
143   delete fTriggerClasses;
144   delete fTriggerClassIndex;
145
146   //#if defined(READOCDB)
147   //delete fTriggerRunScalers; // Not owner -> Owned by OCDB
148   //#endif
149
150 }
151
152
153 //___________________________________________________________________________
154 void AliAnalysisTaskPileup::NotifyRun()
155 {
156   /// Notify run
157
158 #if defined(READOCDB)
159   if ( fStorageList.IsNull() ) {
160     AliError("Default storage not set! Pileup corrections won't be calculated!");
161     return;
162   }
163   if ( ! AliCDBManager::Instance()->GetDefaultStorage() ) {
164     TObjArray* storageList = fStorageList.Tokenize(" ");
165     storageList->SetOwner();
166     for ( Int_t ientry=0; ientry<storageList->GetEntries(); ientry++ ) {
167       TObjString* calibStr = (TObjString*)storageList->At(ientry);
168       ientry++;
169       TObjString* dbStr = (TObjString*)storageList->At(ientry);
170       TString calibName = calibStr->GetString();
171       if ( ! calibName.CompareTo("default") ) {
172         AliCDBManager::Instance()->SetDefaultStorage(dbStr->GetName());
173       }
174       else {
175         AliCDBManager::Instance()->SetSpecificStorage(calibStr->GetName(), dbStr->GetName());
176       }
177     }
178     delete storageList;
179   }
180
181   // Default storage was not correclty set: nothing done
182   if ( ! AliCDBManager::Instance()->GetDefaultStorage() ) return;
183
184   AliCDBManager::Instance()->SetRun(InputEvent()->GetRunNumber());
185
186   AliCDBEntry *entry = 0x0;
187   entry = AliCDBManager::Instance()->Get("GRP/CTP/Config");
188   if ( ! entry ) return;
189   AliTriggerConfiguration* trigConf = (AliTriggerConfiguration*)entry->GetObject();
190   const TObjArray& classesArray = trigConf->GetClasses();
191   Int_t nclasses = classesArray.GetEntriesFast();
192
193   if ( fTriggerClasses ) delete fTriggerClasses;
194
195   fTriggerClasses = new TObjArray(nclasses);
196   fTriggerClasses->SetOwner();
197
198   Int_t currActive = -1;
199   for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
200     AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
201     if (trclass && trclass->GetMask()>0) {
202       currActive++;
203       Int_t currPos = currActive;
204
205       // Store the CBEAMB class at the first position
206       TString trigName = trclass->GetName();
207       if ( ( trigName.Contains("CBEAMB") ||  trigName.Contains("CTRUE") ) && ! trigName.Contains("WU") )
208         currPos = 0;
209       else if ( ! fTriggerClasses->At(0) )
210         currPos++;
211
212       Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
213       TObjString* objString = new TObjString(trigName.Data());
214       fTriggerClasses->AddAtAndExpand(objString, currPos);
215       (*fTriggerClassIndex)[currPos] = trindex;
216       if ( fDebug >= 3 ) printf("AliAnalysisTaskPileup: Current class %s  index %i  position %i\n", trigName.Data(), trindex, currPos);
217     } // is class active
218   } // loop on trigger classes
219
220   entry = AliCDBManager::Instance()->Get("GRP/CTP/Scalers");
221   if ( ! entry ) return;
222   AliInfo("Found an AliTriggerRunScalers in GRP/CTP/Scalers, reading it");
223   fTriggerRunScalers = static_cast<AliTriggerRunScalers*> (entry->GetObject());
224   if (fTriggerRunScalers->CorrectScalersOverflow() == 0) AliInfo("32bit Trigger counters corrected for overflow");
225
226   fIsInitCDB = kTRUE;
227
228 #endif
229
230 }
231
232
233 //___________________________________________________________________________
234 void AliAnalysisTaskPileup::UserCreateOutputObjects()
235 {
236   /// Create histograms and counters
237
238   fTriggerClassIndex = new TArrayI(50);
239   fTriggerClassIndex->Reset(-1);
240
241   fCentralityClasses = new TAxis(20, 0., 100.);
242
243   TString centralityClassesStr = "", currClass = "";
244   for ( Int_t ibin=1; ibin<=fCentralityClasses->GetNbins(); ibin++ ){
245     if ( ! centralityClassesStr.IsNull() )
246       centralityClassesStr.Append("/");
247     currClass = Form("%.0f-%.0f",fCentralityClasses->GetBinLowEdge(ibin),fCentralityClasses->GetBinUpEdge(ibin));
248     centralityClassesStr += currClass;
249     fCentralityClasses->SetBinLabel(ibin, currClass.Data());
250   }
251   
252   // The framework has problems if the name of the object
253   // and the one of container differ
254   // To be complaint, get the name from container and set it
255   TString containerName = GetOutputSlot(1)->GetContainer()->GetName();
256
257   // initialize event counters
258   fEventCounters = new AliCounterCollection(containerName.Data());
259   fEventCounters->AddRubric("event", "any/correctedL0");
260   fEventCounters->AddRubric("trigger", 1000000);
261   fEventCounters->AddRubric("vtxSelection", "any/hasVtxContrib/nonPileupSPD");
262   fEventCounters->AddRubric("selected", "yes/no");
263   fEventCounters->AddRubric("centrality", centralityClassesStr.Data());
264   fEventCounters->Init(kTRUE);
265
266   // Post data at least once per task to ensure data synchronisation (required for merging)
267   PostData(1, fEventCounters);
268 }
269
270 //________________________________________________________________________
271 void AliAnalysisTaskPileup::UserExec(Option_t *)
272 {
273   /// Called for each event
274
275   AliAODEvent* aodEvent = 0x0;
276
277   AliESDEvent* esdEvent = dynamic_cast<AliESDEvent*> (InputEvent());
278   if ( ! esdEvent ) {
279     aodEvent = dynamic_cast<AliAODEvent*> (InputEvent());
280   }
281
282   if ( ! aodEvent && ! esdEvent ) {
283     AliError ("AOD or ESD event not found. Nothing done!");
284     return;
285   }
286  
287   // check physics selection
288   Bool_t isPhysicsSelected = (fInputHandler && fInputHandler->IsEventSelected());
289   TString selected = ( isPhysicsSelected ) ? "yes" : "no";
290
291   TString firedTrigClasses = ( esdEvent ) ? esdEvent->GetFiredTriggerClasses() : aodEvent->GetFiredTriggerClasses();
292
293   if ( ! fIsInitCDB ) {
294     delete fTriggerClasses;
295     fTriggerClasses = firedTrigClasses.Tokenize(" ");
296     fTriggerClasses->SetOwner();
297   }
298
299 #if defined(READOCDB)
300
301   //if ( firedTrigClasses.IsNull() ) return; // Reject non-physics events
302
303   if ( fDebug >= 2 ) printf("\nAliAnalysisTaskPileup: Event %lli\n", Entry());
304
305   // If scalers is not correctly loaded, the task becomes a mere event counter
306   Int_t nPoints = ( fTriggerRunScalers ) ? fTriggerRunScalers->GetScalersRecordsESD()->GetEntriesFast() : 0;
307
308   AliTriggerScalersRecordESD* trigScalerRecords1 = 0x0;
309   AliTriggerScalersRecordESD* trigScalerRecords2 = 0x0;
310   if ( nPoints > 1 ) {
311     // Add protection for MC (no scalers there!)
312
313     AliTimeStamp timeStamp(InputEvent()->GetOrbitNumber(), InputEvent()->GetPeriodNumber(), InputEvent()->GetBunchCrossNumber());
314     Int_t position = fTriggerRunScalers->FindNearestScalersRecord(&timeStamp);
315     if ( position < 0 ) {
316       AliWarning("Position out of range: put to 1");
317       position = 1;
318     } 
319     if ( position == 0 ) position++; // Don't trust the first one
320     else if ( position + 1 >= nPoints ) position--;
321     if ( fDebug >= 2 ) printf("AliAnalysisTaskPileup: position %i\n", position);
322     trigScalerRecords1 = (AliTriggerScalersRecordESD*)fTriggerRunScalers->GetScalersRecordsESD()->At(position);
323
324     // Sometimes scalers are filled very close to each others
325     // in this case skip to the next entry
326     for ( Int_t ipos=position+1; ipos<nPoints; ipos++ ) {
327       trigScalerRecords2 = (AliTriggerScalersRecordESD*)fTriggerRunScalers->GetScalersRecordsESD()->At(ipos);
328       Double_t deltaTime = (Double_t)( trigScalerRecords2->GetTimeStamp()->GetSeconds() - trigScalerRecords1->GetTimeStamp()->GetSeconds() );
329       if ( fDebug >= 2 ) printf("AliAnalysisTaskPileup: Pos %i  TimeStamp %u - %u = %.0f\n", ipos, trigScalerRecords2->GetTimeStamp()->GetSeconds(), trigScalerRecords1->GetTimeStamp()->GetSeconds(), deltaTime);
330
331       if ( deltaTime > 1 )
332         break;
333     } // loop on position
334   } // nPoins > 0
335
336 #endif
337
338   ULong64_t trigMask = 0;
339   Double_t correctFactor = 1.;
340
341   Int_t nVtxContrib = ( esdEvent ) ? esdEvent->GetPrimaryVertex()->GetNContributors() : aodEvent->GetPrimaryVertex()->GetNContributors();
342   Bool_t isPileupSPD = ( esdEvent ) ? esdEvent->IsPileupFromSPD(3,0.8) : aodEvent->IsPileupFromSPD(3,0.8);
343
344   Double_t centralityClass = InputEvent()->GetCentrality()->GetCentralityPercentile("V0M");
345   // If percentile is exactly 100. the bin chosen is in overflow
346   // fix it by setting 99.999 instead of 100.
347   if ( centralityClass >= 100. ) centralityClass = 99.999;
348   Int_t centralityBin = fCentralityClasses->FindBin(centralityClass);
349
350   TString vtxSelKey[3] = {"any","hasVtxContrib","nonPileupSPD"};
351   Bool_t fillSel[3] = {kTRUE, ( nVtxContrib > 0 ), ( ( nVtxContrib > 0 ) && ( ! isPileupSPD ) )};
352
353   //const AliTriggerScalersRecordESD* trigScalerRecords = esdEvent->GetHeader()->GetTriggerScalersRecord(); // REMEMBER TO CUT
354
355   TString trigName = "";
356   TString eventType = "";
357
358   // Get entries counts the number of entries != 0
359   // However the position 0 can be empty in MC (it is reserved to CBEAMB)
360   // In this case with GetEntries the last entry is lost
361   // Use GetEntriesFast instead
362   Int_t nTriggerClasses = fTriggerClasses->GetEntriesFast();
363   Int_t classIndex = -1;
364 #if defined(READOCDB)
365   Double_t deltaScalersBeam = 0., deltaScalers = 0.;
366 #endif
367   //Bool_t isFiredOnce = kFALSE;
368   for (Int_t itrig=0; itrig<nTriggerClasses+1; itrig++) {
369
370     Double_t correctFactorL0 = 1.;
371
372     Bool_t isClassFired = kFALSE;
373
374     if ( itrig < nTriggerClasses ) {
375
376       if ( fIsInitCDB ) {
377         // Check if current mask contains trigger
378         classIndex = (*fTriggerClassIndex)[itrig];
379         if ( classIndex < 0 ) continue; // Protection for MC (where BEAMB not present
380         trigMask = ( 1ull << classIndex );
381         isClassFired = ( trigMask & InputEvent()->GetTriggerMask() );
382       }
383       else
384         isClassFired = kTRUE;
385
386       trigName = ((TObjString*)fTriggerClasses->At(itrig))->GetString();
387
388 #if defined(READOCDB)
389       if ( trigScalerRecords2 && ( isClassFired || itrig == 0 ) ) {
390         // Get scalers
391         const AliTriggerScalersESD* scaler1 = trigScalerRecords1->GetTriggerScalersForClass(classIndex+1);
392         const AliTriggerScalersESD* scaler2 = trigScalerRecords2->GetTriggerScalersForClass(classIndex+1);
393         deltaScalers = scaler2->GetLOCB() - scaler1->GetLOCB();
394  
395         if ( itrig == 0 )
396           deltaScalersBeam = deltaScalers;
397         else if ( isClassFired ) {
398           correctFactorL0 = GetL0Correction(deltaScalers, deltaScalersBeam);
399           if ( fDebug >= 2 ) printf("AliAnalysisTaskPileup: Scalers: %s %.0f  %s %.0f -> CF %f\n", fTriggerClasses->At(itrig)->GetName(), deltaScalers, fTriggerClasses->At(0)->GetName(), deltaScalersBeam, correctFactorL0);
400         }
401       }
402 #endif
403     } // if ( itrig < nTriggerClasses )
404     else {
405       classIndex = -1;
406       trigName = "any";
407       isClassFired = kTRUE; // isFiredOnce;
408     }
409
410     if ( ! isClassFired ) continue;
411     //isFiredOnce = kTRUE;
412
413     //const AliTriggerScalersESD* trigScaler = trigScalerRecords->GetTriggerScalersForClass(classIndex+1); // REMEMBER TO CUT
414     //if ( classIndex > 1 ) printf("Index: trigger %i  scaler %i\n", classIndex+1, trigScaler->GetClassIndex()); // REMEMBER TO CUT
415
416     if ( fDebug >= 2 ) printf("AliAnalysisTaskPileup: Fired trig %s\n", trigName.Data());
417
418     for ( Int_t ievType=0; ievType<2; ievType++ ){
419       switch ( ievType ) {
420       case kHeventsCorrectL0:
421         correctFactor = correctFactorL0;
422         eventType = "correctedL0";
423         break;
424       default:
425         correctFactor = 1.;
426         eventType = "any";
427         break;
428       }
429
430       for ( Int_t isel=0; isel<3; isel++ ) {
431         if ( ! fillSel[isel] ) continue;
432         fEventCounters->Count(Form("event:%s/trigger:%s/vtxSelection:%s/selected:%s/centrality:%s",eventType.Data(),trigName.Data(), vtxSelKey[isel].Data(), selected.Data(),fCentralityClasses->GetBinLabel(centralityBin)),correctFactor);
433       } // loop on vertex selection
434     } // loop on event type
435   } // loop on trigger classes
436
437   // Post final data. It will be written to a file with option "RECREATE"
438   PostData(1, fEventCounters);
439 }
440
441
442 //________________________________________________________________________
443 void AliAnalysisTaskPileup::Terminate(Option_t *)
444 {
445   //
446   /// Save final histograms
447   /// and draw result to the screen
448   //
449
450   // Fill the container only at the very last step
451   // i.e. in local, when done interactively
452   if ( gROOT->IsBatch() )
453     return;
454
455   fEventCounters = dynamic_cast<AliCounterCollection*>(GetOutputData(1));
456   if ( ! fEventCounters ) return;
457
458   if ( ! fHistoEventsList ) fHistoEventsList = new TObjArray(0);
459   fHistoEventsList->SetOwner();
460
461   TH2D* histo = 0x0;
462   const Int_t kNevtTypes = 2;
463   TString evtSel[kNevtTypes] = {"event:any", "event:correctedL0"};
464   TString evtName[kNevtTypes] = {"", "CorrectL0"};
465   TString evtTitle[kNevtTypes] = {"Events", "L0 corrected events"};
466   const Int_t kNphysSel = 2;
467   TString physSel[kNphysSel] = {"selected:any","selected:yes"};
468   TString physName[kNphysSel] = {"", "PhysSel"};
469   TString physTitle[kNphysSel] = {"", "w/ physics selection"};
470
471   TString typeName[2] = {"", "Centrality"};
472
473   Int_t currHisto = -1;
474   TString currName = "";
475   for ( Int_t itype=0; itype<2; itype++ ) {
476     for ( Int_t isel=0; isel<kNphysSel; isel++ ) {
477       for ( Int_t iev=0; iev<kNevtTypes; iev++ ) {
478         currName = Form("%s/%s", evtSel[iev].Data(), physSel[isel].Data());
479         if ( itype == 0 )
480           histo = fEventCounters->Get("trigger","vtxSelection",currName.Data());
481         else {
482           currName.Append("/vtxSelection:any");
483           histo = fEventCounters->Get("trigger","centrality",currName.Data());
484         }
485         if ( ! histo ) continue;
486         currHisto++;
487         currName = Form("hEvents%s%s%s", typeName[itype].Data(), evtName[iev].Data(), physName[isel].Data());
488         histo->SetName(currName.Data());
489         currName = Form("%s %s", evtTitle[iev].Data(), physTitle[isel].Data());
490         histo->SetTitle(currName.Data());
491         fHistoEventsList->AddAtAndExpand(histo, currHisto);
492       } // loop on event type
493       TH2D* num = (TH2D*)fHistoEventsList->At(currHisto);
494       TH2D* den = (TH2D*)fHistoEventsList->At(currHisto-1);
495       if ( ! num || ! den ) continue;
496       currName = Form("hPileupL0%s%s_%s", typeName[itype].Data(), physName[isel].Data(), GetName());
497       TH2D* histoPileupL0 = (TH2D*)num->Clone(currName.Data());
498       histoPileupL0->Divide(den);
499       currName.ReplaceAll("hPileupL0","canPileupL0");
500       TCanvas *can = new TCanvas(currName.Data(),"Pileup",10,10,310,310);
501       can->SetFillColor(10); can->SetHighLightColor(10);
502       can->SetLeftMargin(0.15); can->SetBottomMargin(0.15);
503       histoPileupL0->DrawCopy("text");
504       delete histoPileupL0;
505     } // loop on physics selection
506   } // loop on histo type
507
508   PostData(2, fHistoEventsList);
509 }
510
511
512 //________________________________________________________________________
513 Double_t AliAnalysisTaskPileup::GetL0Correction(Double_t nCINT1B, Double_t nCBEAMB)
514 {
515   /// Get the correction factor for L0 calculation
516
517   if ( nCBEAMB == 0. )
518     return 1.;
519
520   Double_t ratio = nCINT1B / nCBEAMB;
521
522   if ( ratio >= 1. || ratio == 0. )
523     return 1.;
524
525   Double_t mu = -TMath::Log(1-ratio);
526
527   return mu / ( 1. - TMath::Exp(-mu) );
528
529 }
530
531 //________________________________________________________________________
532 void AliAnalysisTaskPileup::SetDefaultStorage(TString dbString)
533 {
534   /// Set default storage
535   SetSpecificStorage("default", dbString);
536 }
537
538
539 //________________________________________________________________________
540 void AliAnalysisTaskPileup::SetSpecificStorage(TString calibType, TString dbString)
541 {
542   /// Set specific storage
543 #if defined(READOCDB)
544   if ( ! fStorageList.IsNull() ) fStorageList += " ";
545   fStorageList += Form("%s %s", calibType.Data(), dbString.Data());
546   if ( fDebug >= 3 ) printf("AliAnalysisTaskPileup: Setting storage %s %s", calibType.Data(), dbString.Data());
547 #else
548   calibType = "";
549   dbString  = "";
550   AliWarning("Class was not compiled to run on OCDB. Command will not have effect");
551 #endif
552 }
553