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