]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliPhysicsSelection.cxx
adding MC support
[u/mrichter/AliRoot.git] / ANALYSIS / AliPhysicsSelection.cxx
1 /* $Id: AliPhysicsSelection.cxx 35782 2009-10-22 11:54:31Z jgrosseo $ */
2
3 /**************************************************************************
4  * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
5  *                                                                        *
6  * Author: The ALICE Off-line Project.                                    *
7  * Contributors are mentioned in the code where appropriate.              *
8  *                                                                        *
9  * Permission to use, copy, modify and distribute this software and its   *
10  * documentation strictly for non-commercial purposes is hereby granted   *
11  * without fee, provided that the above copyright notice appears in all   *
12  * copies and that both the copyright notice and this permission notice   *
13  * appear in the supporting documentation. The authors make no claims     *
14  * about the suitability of this software for any purpose. It is          *
15  * provided "as is" without express or implied warranty.                  *
16  **************************************************************************/
17
18 //-------------------------------------------------------------------------
19 //                      Implementation of   Class AliPhysicsSelection
20 // This class selects collision candidates from data runs, applying selection cuts on triggers 
21 // and background rejection based on the content of the ESD
22 //
23 // Usage:
24 //
25 // Create the object:
26 //   fPhysicsSelection = new AliPhysicsSelection;
27 //
28 // For MC data, call
29 //   fPhysicsSelection->SetAnalyzeMC()
30 //
31 // To check if an event is a collision candidate, use:
32 //   fPhysicsSelection->IsCollisionCandidate(fESD)
33 //
34 // After processing save the resulting histograms to a file with (a folder physics_selection 
35 //   will be created that contains the histograms):
36 //   fPhysicsSelection->SaveHistograms("physics_selection")
37 //
38 // To print statistics after processing use:
39 //   fPhysicsSelection->Print();
40 //
41 //   Origin: Jan Fiete Grosse-Oetringhaus, CERN
42 //-------------------------------------------------------------------------
43
44 #include <Riostream.h>
45 #include <TH1F.h>
46 #include <TH2F.h>
47 #include <TList.h>
48 #include <TIterator.h>
49 #include <TDirectory.h>
50 #include <TObjArray.h>
51
52 #include <AliPhysicsSelection.h>
53
54 #include <AliTriggerAnalysis.h>
55 #include <AliLog.h>
56
57 #include <AliESDEvent.h>
58
59 ClassImp(AliPhysicsSelection)
60
61 AliPhysicsSelection::AliPhysicsSelection() :
62   AliAnalysisCuts("AliPhysicsSelection", "AliPhysicsSelection"),
63   fCurrentRun(-1),
64   fMC(kFALSE),
65   fCollTrigClasses(),
66   fBGTrigClasses(),
67   fTriggerAnalysis(),
68   fBackgroundIdentification(0),
69   fHistStatistics(0),
70   fHistBunchCrossing(0)
71 {
72   // constructor
73   
74   fCollTrigClasses.SetOwner(1);
75   fBGTrigClasses.SetOwner(1);
76   fTriggerAnalysis.SetOwner(1);
77   
78   AliLog::SetClassDebugLevel("AliPhysicsSelection", AliLog::kWarning);
79 }
80     
81 AliPhysicsSelection::~AliPhysicsSelection()
82 {
83   // destructor
84   
85   fCollTrigClasses.Delete();
86   fBGTrigClasses.Delete();
87   fTriggerAnalysis.Delete();
88
89   if (fHistStatistics)
90   {
91     delete fHistStatistics;
92     fHistStatistics = 0;
93   }
94   
95   if (fHistBunchCrossing)
96   {
97     delete fHistBunchCrossing;
98     fHistBunchCrossing = 0;
99   }
100 }
101
102 Bool_t AliPhysicsSelection::CheckTriggerClass(const AliESDEvent* aEsd, const char* trigger) const
103 {
104   // checks if the given trigger class(es) are found for the current event
105   // format of trigger: +TRIGGER1 -TRIGGER2
106   //   requires TRIGGER1 and rejects TRIGGER2
107   
108   TString str(trigger);
109   TObjArray* tokens = str.Tokenize(" ");
110   
111   for (Int_t i=0; i < tokens->GetEntries(); i++)
112   {
113     TString str2(((TObjString*) tokens->At(i))->String());
114     
115     if (str2[0] != '+' && str2[0] != '-')
116       AliFatal(Form("Invalid trigger syntax: %s", trigger));
117       
118     Bool_t flag = (str2[0] == '+');
119     
120     str2.Remove(0, 1);
121     
122     if (flag && !aEsd->IsTriggerClassFired(str2))
123     {
124       AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is not present", str2.Data()));
125       delete tokens;
126       return kFALSE;
127     }
128     if (!flag && aEsd->IsTriggerClassFired(str2))
129     {
130       AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is present", str2.Data()));
131       delete tokens;
132       return kFALSE;
133     }
134   }
135   
136   delete tokens;
137   return kTRUE;
138 }
139     
140 Bool_t AliPhysicsSelection::IsCollisionCandidate(const AliESDEvent* aEsd)
141 {
142   // checks if the given event is a collision candidate
143   
144   const AliESDHeader* esdHeader = aEsd->GetHeader();
145   if (!esdHeader)
146   {
147     AliError("ESD Header could not be retrieved");
148     return kFALSE;
149   }
150   
151   // check event type; should be PHYSICS = 7 for data and 0 for MC
152   if (!fMC)
153   {
154     if (esdHeader->GetEventType() != 7)
155       return kFALSE;
156   }
157   else
158   {
159     if (esdHeader->GetEventType() != 0)
160       AliFatal(Form("Invalid event type for MC: %d", esdHeader->GetEventType()));
161   }
162     
163   if (fCurrentRun != aEsd->GetRunNumber())
164     if (!Initialize(aEsd->GetRunNumber()))
165       AliFatal(Form("Could not initialize for run %d", aEsd->GetRunNumber()));
166     
167   Bool_t accept = kFALSE;
168     
169   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
170   for (Int_t i=0; i < count; i++)
171   {
172     const char* triggerClass = 0;
173     if (i < fCollTrigClasses.GetEntries())
174       triggerClass = ((TObjString*) fCollTrigClasses.At(i))->String();
175     else
176       triggerClass = ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
177   
178     AliDebug(AliLog::kDebug, Form("Processing trigger class %s", triggerClass));
179   
180     AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
181   
182     triggerAnalysis->FillTriggerClasses(aEsd);
183     
184     if (CheckTriggerClass(aEsd, triggerClass))
185     {
186       triggerAnalysis->FillHistograms(aEsd);
187   
188       fHistStatistics->Fill(1, i);
189     
190       Int_t fastOR = triggerAnalysis->SPDFiredChips(aEsd, 0);
191       Bool_t v0A = triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0A);
192       Bool_t v0C = triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0C);
193       Bool_t v0ABG = triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0ABG);
194       Bool_t v0CBG = triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0CBG);
195       Bool_t v0BG = v0ABG || v0CBG;
196     
197       if (fastOR > 0)
198         fHistStatistics->Fill(2, i);
199       if (fastOR > 1)
200         fHistStatistics->Fill(3, i);
201         
202       if (v0A)
203         fHistStatistics->Fill(4, i);
204       if (v0C)
205         fHistStatistics->Fill(5, i);
206       if (v0ABG)
207         fHistStatistics->Fill(6, i);
208       if (v0CBG)
209         fHistStatistics->Fill(7, i);
210         
211       if (fastOR > 1 && !v0BG)
212         fHistStatistics->Fill(8, i);
213         
214       if ((fastOR > 0 || v0A || v0C) && !v0BG)
215         fHistStatistics->Fill(9, i);
216     
217       if (fastOR > 0 && (v0A || v0C) && !v0BG)
218         fHistStatistics->Fill(10, i);
219   
220       if (v0A && v0C && !v0BG)
221         fHistStatistics->Fill(11, i);
222         
223       if (fastOR > 1 || (fastOR > 0 && (v0A || v0C)) || (v0A && v0C))
224       {
225         if (!v0BG)
226         {
227           fHistStatistics->Fill(12, i);
228       
229           if (fBackgroundIdentification && !fBackgroundIdentification->IsSelected(const_cast<AliESDEvent*> (aEsd)))
230           {
231             AliDebug(AliLog::kDebug, "Rejecting event because of background identification");
232             fHistStatistics->Fill(13, i);
233           }
234           else
235           {
236             AliDebug(AliLog::kDebug, "Accepted event for histograms");
237             
238             fHistStatistics->Fill(14, i);
239             fHistBunchCrossing->Fill(aEsd->GetBunchCrossNumber(), i);
240             if (i < fCollTrigClasses.GetEntries())
241               accept = kTRUE;
242           }
243         }
244         else
245           AliDebug(AliLog::kDebug, "Rejecting event because of V0 BG flag");
246       }
247       else
248         AliDebug(AliLog::kDebug, "Rejecting event because trigger condition is not fulfilled");
249     }
250   }
251  
252   if (accept)
253     AliDebug(AliLog::kDebug, "Accepted event as collision candidate");
254   
255   return accept;
256 }
257
258 Int_t AliPhysicsSelection::GetTriggerScheme(UInt_t /* runNumber */)
259 {
260   // returns the current trigger scheme (classes that are accepted/rejected)
261   
262   if (fMC)
263     return 0;
264     
265   // TODO dependent on run number
266   return 1;
267 }  
268     
269 Bool_t AliPhysicsSelection::Initialize(UInt_t runNumber)
270 {
271   // initializes the object for the given run
272   // TODO having the run number here and parameters hardcoded is clearly temporary, a way needs to be found to have a CDB-like configuration also for analysis
273   
274   Bool_t oldStatus = TH1::AddDirectoryStatus();
275   TH1::AddDirectory(kFALSE);
276   
277   Int_t triggerScheme = GetTriggerScheme(runNumber);
278   
279   if (fCurrentRun != -1 && triggerScheme != GetTriggerScheme(fCurrentRun))
280     AliFatal("Processing several runs with different trigger schemes is not supported");
281   
282   AliInfo(Form("Initializing for run %d", runNumber));
283   
284   Bool_t firstTime = kFALSE;
285   if (fCurrentRun == -1)
286   {
287     // initialize first time
288     fCurrentRun = runNumber;
289     firstTime = kTRUE;
290   }
291   
292   if (firstTime)
293   {
294     switch (triggerScheme)
295     {
296       case 0:
297         fCollTrigClasses.Add(new TObjString(""));
298         break;
299       
300       case 1:
301         fCollTrigClasses.Add(new TObjString("+CINT1B-ABCE-NOPF-ALL"));
302         fBGTrigClasses.Add(new TObjString("+CINT1A-ABCE-NOPF-ALL"));
303         fBGTrigClasses.Add(new TObjString("+CINT1C-ABCE-NOPF-ALL"));
304         fBGTrigClasses.Add(new TObjString("+CINT1-E-NOPF-ALL"));
305         break;
306         
307       default:
308         AliFatal(Form("Unsupported trigger scheme %d", triggerScheme));
309     }
310   
311     Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
312     
313     for (Int_t i=0; i<count; i++)
314     {
315       AliTriggerAnalysis* triggerAnalysis = new AliTriggerAnalysis;
316       triggerAnalysis->SetAnalyzeMC(fMC);
317       triggerAnalysis->EnableHistograms();
318       triggerAnalysis->SetSPDGFOThreshhold(1);
319       fTriggerAnalysis.Add(triggerAnalysis);
320     }
321       
322     if (fHistStatistics)
323       delete fHistStatistics;
324   
325     fHistStatistics = new TH2F("fHistStatistics", "fHistStatistics;;", 14, 0.5, 14.5, count, -0.5, -0.5 + count);
326       
327     Int_t n = 1;
328     fHistStatistics->GetXaxis()->SetBinLabel(n++, "Triggered");
329     fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 1");
330     fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 2");
331     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0A");
332     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0C");
333     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0A BG");
334     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0C BG");
335     fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 2 &!V0 BG");
336     fHistStatistics->GetXaxis()->SetBinLabel(n++, "(FO >= 1 | V0A | V0C) & !V0 BG");
337     fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 1 & (V0A | V0C) & !V0 BG");
338     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0A & V0C & !V0 BG");
339     fHistStatistics->GetXaxis()->SetBinLabel(n++, "(FO >= 2 | (FO >= 1 & (V0A | V0C)) | (V0A & V0C)) & !V0 BG");
340     fHistStatistics->GetXaxis()->SetBinLabel(n++, "Background identification");
341     fHistStatistics->GetXaxis()->SetBinLabel(n++, "Accepted");
342     
343     if (fHistBunchCrossing)
344       delete fHistBunchCrossing;
345   
346     fHistBunchCrossing = new TH2F("fHistBunchCrossing", "fHistBunchCrossing;bunch crossing number;", 4000, -0.5, 3999.5,  count, -0.5, -0.5 + count);
347       
348     n = 1;
349     for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
350     {
351       fHistStatistics->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
352       fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
353       n++;
354     }
355     for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
356     {
357       fHistStatistics->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
358       fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
359       n++;
360     }
361   }
362     
363   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
364   for (Int_t i=0; i<count; i++)
365   {
366     AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
367   
368     switch (runNumber)
369     {
370       case 104315:
371       case 104316:
372       case 104320:
373       case 104321:
374       case 104439:
375         triggerAnalysis->SetV0TimeOffset(7.5);
376         break;
377       default:
378         triggerAnalysis->SetV0TimeOffset(0);
379     }
380   }
381     
382   TH1::AddDirectory(oldStatus);
383   
384   return kTRUE;
385 }
386
387 void AliPhysicsSelection::Print(Option_t* /* option */) const
388 {
389   // print the configuration
390   
391   Printf("Configuration initialized for run %d (MC: %d):", fCurrentRun, fMC);
392   
393   Printf("Collision trigger classes:");
394   for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
395     Printf("%s", ((TObjString*) fCollTrigClasses.At(i))->String().Data());
396   
397   Printf("Background trigger classes:");
398   for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
399     Printf("%s", ((TObjString*) fBGTrigClasses.At(i))->String().Data());
400
401   AliTriggerAnalysis* triggerAnalysis = dynamic_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(0));
402   
403   if (triggerAnalysis)
404   {
405     if (triggerAnalysis->GetV0TimeOffset() > 0)
406       Printf("V0 time offset active: %.2f ns", triggerAnalysis->GetV0TimeOffset());
407     
408     Printf("\nTotal available events:");
409     
410     triggerAnalysis->PrintTriggerClasses();
411   }
412   
413   if (fHistStatistics)
414   {
415     Printf("\nSelection statistics for first collision trigger:");
416     
417     Printf("Total events with correct trigger class: %d", (Int_t) fHistStatistics->GetBinContent(1, 1));
418     Printf("Selected collision candidates: %d", (Int_t) fHistStatistics->GetBinContent(fHistStatistics->GetXaxis()->FindBin("Accepted"), 1));
419   }
420   
421   if (fHistBunchCrossing)
422   {
423     Printf("\nBunch crossing statistics:");
424     
425     for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
426     {
427       TString str;
428       str.Form("Trigger %s has accepted events in the bunch crossings: ", fHistBunchCrossing->GetYaxis()->GetBinLabel(i));
429       
430       for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
431         if (fHistBunchCrossing->GetBinContent(j, i) > 0)
432           str += Form("%d, ", (Int_t) fHistBunchCrossing->GetXaxis()->GetBinCenter(j));
433        
434       Printf("%s", str.Data());
435     }
436     
437     for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
438     {
439       Int_t count = 0;
440       for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
441       {
442         if (fHistBunchCrossing->GetBinContent(j, i) > 0)
443           count++;
444       }
445       if (count > 1)
446         Printf("WARNING: Bunch crossing %d has more than one trigger class active. Check BPTX functioning for this run!", (Int_t) fHistBunchCrossing->GetXaxis()->GetBinCenter(j));
447     }
448   }
449 }
450
451 Long64_t AliPhysicsSelection::Merge(TCollection* list)
452 {
453   // Merge a list of AliMultiplicityCorrection objects with this (needed for
454   // PROOF).
455   // Returns the number of merged objects (including this).
456
457   if (!list)
458     return 0;
459
460   if (list->IsEmpty())
461     return 1;
462
463   TIterator* iter = list->MakeIterator();
464   TObject* obj;
465   
466   // collections of all histograms
467   const Int_t nHists = 9;
468   TList collections[nHists];
469
470   Int_t count = 0;
471   while ((obj = iter->Next())) {
472
473     AliPhysicsSelection* entry = dynamic_cast<AliPhysicsSelection*> (obj);
474     if (entry == 0) 
475       continue;
476       
477     collections[0].Add(&(entry->fTriggerAnalysis));
478     if (entry->fHistStatistics)
479       collections[1].Add(entry->fHistStatistics);
480     if (entry->fHistBunchCrossing)
481       collections[2].Add(entry->fHistBunchCrossing);
482     if (entry->fBackgroundIdentification)
483       collections[3].Add(entry->fBackgroundIdentification);
484
485     count++;
486   }
487
488   fTriggerAnalysis.Merge(&collections[0]);
489   if (fHistStatistics)
490     fHistStatistics->Merge(&collections[1]);
491   if (fHistBunchCrossing)
492     fHistBunchCrossing->Merge(&collections[2]);
493   if (fBackgroundIdentification)
494     fBackgroundIdentification->Merge(&collections[3]);
495   
496   delete iter;
497
498   return count+1;
499 }
500
501 void AliPhysicsSelection::SaveHistograms(const char* folder) const
502 {
503   // write histograms to current directory
504   
505   if (!fHistStatistics)
506     return;
507     
508   if (folder)
509   {
510     gDirectory->mkdir(folder);
511     gDirectory->cd(folder);
512   }
513   
514   fHistStatistics->Write();
515   fHistBunchCrossing->Write();
516   
517   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
518   for (Int_t i=0; i < count; i++)
519   {
520     TString triggerClass = "trigger_histograms_";
521     if (i < fCollTrigClasses.GetEntries())
522       triggerClass += ((TObjString*) fCollTrigClasses.At(i))->String();
523     else
524       triggerClass += ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
525   
526     gDirectory->mkdir(triggerClass);
527     gDirectory->cd(triggerClass);
528   
529     static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i))->SaveHistograms();
530     
531     gDirectory->cd("..");
532   }
533  
534   if (fBackgroundIdentification)
535   {
536     gDirectory->mkdir("background_identification");
537     gDirectory->cd("background_identification");
538       
539     fBackgroundIdentification->GetOutput()->Write();
540       
541     gDirectory->cd("..");
542   }
543   
544   if (folder)
545     gDirectory->cd("..");
546 }