adding arguments to AddTaskPhysicsSelection for MC and BG identification (Antonin)
[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   fSkipTriggerClassSelection(0),
72   fUsingCustomClasses(0)
73 {
74   // constructor
75   
76   fCollTrigClasses.SetOwner(1);
77   fBGTrigClasses.SetOwner(1);
78   fTriggerAnalysis.SetOwner(1);
79   
80   AliLog::SetClassDebugLevel("AliPhysicsSelection", AliLog::kWarning);
81 }
82     
83 AliPhysicsSelection::~AliPhysicsSelection()
84 {
85   // destructor
86   
87   fCollTrigClasses.Delete();
88   fBGTrigClasses.Delete();
89   fTriggerAnalysis.Delete();
90
91   if (fHistStatistics)
92   {
93     delete fHistStatistics;
94     fHistStatistics = 0;
95   }
96   
97   if (fHistBunchCrossing)
98   {
99     delete fHistBunchCrossing;
100     fHistBunchCrossing = 0;
101   }
102 }
103
104 Bool_t AliPhysicsSelection::CheckTriggerClass(const AliESDEvent* aEsd, const char* trigger) const
105 {
106   // checks if the given trigger class(es) are found for the current event
107   // format of trigger: +TRIGGER1 -TRIGGER2
108   //   requires TRIGGER1 and rejects TRIGGER2
109   
110   TString str(trigger);
111   TObjArray* tokens = str.Tokenize(" ");
112   
113   for (Int_t i=0; i < tokens->GetEntries(); i++)
114   {
115     TString str2(((TObjString*) tokens->At(i))->String());
116     
117     if (str2[0] != '+' && str2[0] != '-')
118       AliFatal(Form("Invalid trigger syntax: %s", trigger));
119       
120     Bool_t flag = (str2[0] == '+');
121     
122     str2.Remove(0, 1);
123     
124     if (flag && !aEsd->IsTriggerClassFired(str2))
125     {
126       AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is not present", str2.Data()));
127       delete tokens;
128       return kFALSE;
129     }
130     if (!flag && aEsd->IsTriggerClassFired(str2))
131     {
132       AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is present", str2.Data()));
133       delete tokens;
134       return kFALSE;
135     }
136   }
137   
138   delete tokens;
139   return kTRUE;
140 }
141     
142 Bool_t AliPhysicsSelection::IsCollisionCandidate(const AliESDEvent* aEsd)
143 {
144   // checks if the given event is a collision candidate
145   
146   if (fCurrentRun != aEsd->GetRunNumber())
147     if (!Initialize(aEsd->GetRunNumber()))
148       AliFatal(Form("Could not initialize for run %d", aEsd->GetRunNumber()));
149     
150   const AliESDHeader* esdHeader = aEsd->GetHeader();
151   if (!esdHeader)
152   {
153     AliError("ESD Header could not be retrieved");
154     return kFALSE;
155   }
156   
157   // check event type; should be PHYSICS = 7 for data and 0 for MC
158   if (!fMC)
159   {
160     if (esdHeader->GetEventType() != 7)
161       return kFALSE;
162   }
163   else
164   {
165     if (esdHeader->GetEventType() != 0)
166       AliFatal(Form("Invalid event type for MC: %d", esdHeader->GetEventType()));
167   }
168   
169   Bool_t accept = kFALSE;
170     
171   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
172   for (Int_t i=0; i < count; i++)
173   {
174     const char* triggerClass = 0;
175     if (i < fCollTrigClasses.GetEntries())
176       triggerClass = ((TObjString*) fCollTrigClasses.At(i))->String();
177     else
178       triggerClass = ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
179   
180     AliDebug(AliLog::kDebug, Form("Processing trigger class %s", triggerClass));
181   
182     AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
183   
184     triggerAnalysis->FillTriggerClasses(aEsd);
185     
186     if (CheckTriggerClass(aEsd, triggerClass))
187     {
188       triggerAnalysis->FillHistograms(aEsd);
189   
190       fHistStatistics->Fill(1, i);
191       
192       // hardware trigger (should only remove events for MC)
193       // replay CINT1B hardware trigger
194       // TODO this has to depend on the actual hardware trigger (and that depends on the run...)
195       Int_t fastORHW = triggerAnalysis->SPDFiredChips(aEsd, 1); // SPD number of chips from trigger bits (!)
196       Bool_t v0A = triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0A);
197       Bool_t v0C = triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0C);
198         
199       if (fastORHW == 0 && !v0A && !v0C)
200       {
201         AliDebug(AliLog::kDebug, "Rejecting event because hardware trigger is not fired");
202         continue;
203       }
204       
205       fHistStatistics->Fill(2, i);
206     
207       // offline trigger
208       Int_t fastOROffline = triggerAnalysis->SPDFiredChips(aEsd, 0); // SPD number of chips from clusters (!)
209       Bool_t v0ABG = triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0ABG);
210       Bool_t v0CBG = triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0CBG);
211       Bool_t v0BG = v0ABG || v0CBG;
212     
213       if (fastOROffline > 0)
214         fHistStatistics->Fill(3, i);
215       if (fastOROffline > 1)
216         fHistStatistics->Fill(4, i);
217         
218       if (v0A)
219         fHistStatistics->Fill(5, i);
220       if (v0C)
221         fHistStatistics->Fill(6, i);
222       if (v0ABG)
223         fHistStatistics->Fill(7, i);
224       if (v0CBG)
225         fHistStatistics->Fill(8, i);
226         
227       if (fastOROffline > 1 && !v0BG)
228         fHistStatistics->Fill(9, i);
229         
230       if ((fastOROffline > 0 || v0A || v0C) && !v0BG)
231         fHistStatistics->Fill(10, i);
232     
233       if (fastOROffline > 0 && (v0A || v0C) && !v0BG)
234         fHistStatistics->Fill(11, i);
235   
236       if (v0A && v0C && !v0BG)
237         fHistStatistics->Fill(12, i);
238         
239       if (fastOROffline > 1 || (fastOROffline > 0 && (v0A || v0C)) || (v0A && v0C))
240       {
241         if (!v0BG)
242         {
243           fHistStatistics->Fill(13, i);
244       
245           if (fBackgroundIdentification && !fBackgroundIdentification->IsSelected(const_cast<AliESDEvent*> (aEsd)))
246           {
247             AliDebug(AliLog::kDebug, "Rejecting event because of background identification");
248             fHistStatistics->Fill(14, i);
249           }
250           else
251           {
252             AliDebug(AliLog::kDebug, "Accepted event for histograms");
253             
254             fHistStatistics->Fill(15, i);
255             fHistBunchCrossing->Fill(aEsd->GetBunchCrossNumber(), i);
256             if (i < fCollTrigClasses.GetEntries() || fSkipTriggerClassSelection)
257               accept = kTRUE;
258           }
259         }
260         else
261           AliDebug(AliLog::kDebug, "Rejecting event because of V0 BG flag");
262       }
263       else
264         AliDebug(AliLog::kDebug, "Rejecting event because trigger condition is not fulfilled");
265     }
266   }
267  
268   if (accept)
269     AliDebug(AliLog::kDebug, "Accepted event as collision candidate");
270   
271   return accept;
272 }
273
274 Int_t AliPhysicsSelection::GetTriggerScheme(UInt_t runNumber)
275 {
276   // returns the current trigger scheme (classes that are accepted/rejected)
277   
278   if (fMC)
279     return 0;
280     
281   // TODO dependent on run number
282   
283   switch (runNumber)
284   {
285     // CSMBB triggers
286     case 104044:
287     case 105054:
288     case 105057:
289       return 2;
290   }
291   
292   // default: CINT1 suite
293   return 1;
294 }  
295     
296 Bool_t AliPhysicsSelection::Initialize(UInt_t runNumber)
297 {
298   // initializes the object for the given run
299   // 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
300   
301   Bool_t oldStatus = TH1::AddDirectoryStatus();
302   TH1::AddDirectory(kFALSE);
303   
304   Int_t triggerScheme = GetTriggerScheme(runNumber);
305   
306   if (!fUsingCustomClasses && fCurrentRun != -1 && triggerScheme != GetTriggerScheme(fCurrentRun))
307     AliFatal("Processing several runs with different trigger schemes is not supported");
308   
309   AliInfo(Form("Initializing for run %d", runNumber));
310   
311   // initialize first time?
312   if (fCurrentRun == -1)
313   {
314     if (fUsingCustomClasses) {
315       AliInfo("Using user-provided trigger classes");
316     } else {
317       switch (triggerScheme)
318       {
319       case 0:
320         fCollTrigClasses.Add(new TObjString(""));
321         break;
322         
323       case 1:
324         fCollTrigClasses.Add(new TObjString("+CINT1B-ABCE-NOPF-ALL"));
325         fBGTrigClasses.Add(new TObjString("+CINT1A-ABCE-NOPF-ALL"));
326         fBGTrigClasses.Add(new TObjString("+CINT1C-ABCE-NOPF-ALL"));
327         fBGTrigClasses.Add(new TObjString("+CINT1-E-NOPF-ALL"));
328         break;
329         
330       case 2:
331         fCollTrigClasses.Add(new TObjString("+CSMBB-ABCE-NOPF-ALL"));
332         fBGTrigClasses.Add(new TObjString("+CSMBA-ABCE-NOPF-ALL -CSMBB-ABCE-NOPF-ALL"));
333         fBGTrigClasses.Add(new TObjString("+CSMBC-ABCE-NOPF-ALL -CSMBB-ABCE-NOPF-ALL"));
334         break;
335         
336       default:
337         AliFatal(Form("Unsupported trigger scheme %d", triggerScheme));
338       }
339     }
340     
341     Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
342     
343     for (Int_t i=0; i<count; i++)
344     {
345       AliTriggerAnalysis* triggerAnalysis = new AliTriggerAnalysis;
346       triggerAnalysis->SetAnalyzeMC(fMC);
347       triggerAnalysis->EnableHistograms();
348       triggerAnalysis->SetSPDGFOThreshhold(1);
349       fTriggerAnalysis.Add(triggerAnalysis);
350     }
351       
352     if (fHistStatistics)
353       delete fHistStatistics;
354   
355     fHistStatistics = new TH2F("fHistStatistics", "fHistStatistics;;", 15, 0.5, 15.5, count, -0.5, -0.5 + count);
356       
357     Int_t n = 1;
358     fHistStatistics->GetXaxis()->SetBinLabel(n++, "Trigger class");
359     fHistStatistics->GetXaxis()->SetBinLabel(n++, "Hardware trigger");
360     fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 1");
361     fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 2");
362     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0A");
363     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0C");
364     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0A BG");
365     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0C BG");
366     fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 2 &!V0 BG");
367     fHistStatistics->GetXaxis()->SetBinLabel(n++, "(FO >= 1 | V0A | V0C) & !V0 BG");
368     fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 1 & (V0A | V0C) & !V0 BG");
369     fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0A & V0C & !V0 BG");
370     fHistStatistics->GetXaxis()->SetBinLabel(n++, "(FO >= 2 | (FO >= 1 & (V0A | V0C)) | (V0A & V0C)) & !V0 BG");
371     fHistStatistics->GetXaxis()->SetBinLabel(n++, "Background identification");
372     fHistStatistics->GetXaxis()->SetBinLabel(n++, "Accepted");
373     
374     if (fHistBunchCrossing)
375       delete fHistBunchCrossing;
376   
377     fHistBunchCrossing = new TH2F("fHistBunchCrossing", "fHistBunchCrossing;bunch crossing number;", 4000, -0.5, 3999.5,  count, -0.5, -0.5 + count);
378       
379     n = 1;
380     for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
381     {
382       fHistStatistics->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
383       fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
384       n++;
385     }
386     for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
387     {
388       fHistStatistics->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
389       fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
390       n++;
391     }
392   }
393     
394   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
395   for (Int_t i=0; i<count; i++)
396   {
397     AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
398   
399     switch (runNumber)
400     {
401       case 104315:
402       case 104316:
403       case 104320:
404       case 104321:
405       case 104439:
406         triggerAnalysis->SetV0TimeOffset(7.5);
407         break;
408       default:
409         triggerAnalysis->SetV0TimeOffset(0);
410     }
411   }
412     
413   fCurrentRun = runNumber;
414   
415   TH1::AddDirectory(oldStatus);
416   
417   return kTRUE;
418 }
419
420 void AliPhysicsSelection::Print(Option_t* /* option */) const
421 {
422   // print the configuration
423   
424   Printf("Configuration initialized for run %d (MC: %d):", fCurrentRun, fMC);
425   
426   Printf("Collision trigger classes:");
427   for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
428     Printf("%s", ((TObjString*) fCollTrigClasses.At(i))->String().Data());
429   
430   Printf("Background trigger classes:");
431   for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
432     Printf("%s", ((TObjString*) fBGTrigClasses.At(i))->String().Data());
433
434   AliTriggerAnalysis* triggerAnalysis = dynamic_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(0));
435   
436   if (triggerAnalysis)
437   {
438     if (triggerAnalysis->GetV0TimeOffset() > 0)
439       Printf("V0 time offset active: %.2f ns", triggerAnalysis->GetV0TimeOffset());
440     
441     Printf("\nTotal available events:");
442     
443     triggerAnalysis->PrintTriggerClasses();
444   }
445   
446   if (fHistStatistics && fCollTrigClasses.GetEntries() > 0)
447   {
448     Printf("\nSelection statistics for first collision trigger (%s):", ((TObjString*) fCollTrigClasses.First())->String().Data());
449     
450     Printf("Total events with correct trigger class: %d", (Int_t) fHistStatistics->GetBinContent(1, 1));
451     Printf("Selected collision candidates: %d", (Int_t) fHistStatistics->GetBinContent(fHistStatistics->GetXaxis()->FindBin("Accepted"), 1));
452   }
453   
454   if (fHistBunchCrossing)
455   {
456     Printf("\nBunch crossing statistics:");
457     
458     for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
459     {
460       TString str;
461       str.Form("Trigger %s has accepted events in the bunch crossings: ", fHistBunchCrossing->GetYaxis()->GetBinLabel(i));
462       
463       for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
464         if (fHistBunchCrossing->GetBinContent(j, i) > 0)
465           str += Form("%d, ", (Int_t) fHistBunchCrossing->GetXaxis()->GetBinCenter(j));
466        
467       Printf("%s", str.Data());
468     }
469     
470     for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
471     {
472       Int_t count = 0;
473       for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
474       {
475         if (fHistBunchCrossing->GetBinContent(j, i) > 0)
476           count++;
477       }
478       if (count > 1)
479         Printf("WARNING: Bunch crossing %d has more than one trigger class active. Check BPTX functioning for this run!", (Int_t) fHistBunchCrossing->GetXaxis()->GetBinCenter(j));
480     }
481   }
482
483   if (fUsingCustomClasses)        
484     Printf("WARNING: Using custom trigger classes!");
485   if (fSkipTriggerClassSelection) 
486     Printf("WARNING: Skipping trigger class selection!");
487 }
488
489 Long64_t AliPhysicsSelection::Merge(TCollection* list)
490 {
491   // Merge a list of AliMultiplicityCorrection objects with this (needed for
492   // PROOF).
493   // Returns the number of merged objects (including this).
494
495   if (!list)
496     return 0;
497
498   if (list->IsEmpty())
499     return 1;
500
501   TIterator* iter = list->MakeIterator();
502   TObject* obj;
503   
504   // collections of all histograms
505   const Int_t nHists = 9;
506   TList collections[nHists];
507
508   Int_t count = 0;
509   while ((obj = iter->Next())) {
510
511     AliPhysicsSelection* entry = dynamic_cast<AliPhysicsSelection*> (obj);
512     if (entry == 0) 
513       continue;
514       
515     collections[0].Add(&(entry->fTriggerAnalysis));
516     if (entry->fHistStatistics)
517       collections[1].Add(entry->fHistStatistics);
518     if (entry->fHistBunchCrossing)
519       collections[2].Add(entry->fHistBunchCrossing);
520     if (entry->fBackgroundIdentification)
521       collections[3].Add(entry->fBackgroundIdentification);
522
523     count++;
524   }
525
526   fTriggerAnalysis.Merge(&collections[0]);
527   if (fHistStatistics)
528     fHistStatistics->Merge(&collections[1]);
529   if (fHistBunchCrossing)
530     fHistBunchCrossing->Merge(&collections[2]);
531   if (fBackgroundIdentification)
532     fBackgroundIdentification->Merge(&collections[3]);
533   
534   delete iter;
535
536   return count+1;
537 }
538
539 void AliPhysicsSelection::SaveHistograms(const char* folder) const
540 {
541   // write histograms to current directory
542   
543   if (!fHistStatistics)
544     return;
545     
546   if (folder)
547   {
548     gDirectory->mkdir(folder);
549     gDirectory->cd(folder);
550   }
551   
552   fHistStatistics->Write();
553   fHistBunchCrossing->Write();
554   
555   Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
556   for (Int_t i=0; i < count; i++)
557   {
558     TString triggerClass = "trigger_histograms_";
559     if (i < fCollTrigClasses.GetEntries())
560       triggerClass += ((TObjString*) fCollTrigClasses.At(i))->String();
561     else
562       triggerClass += ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
563   
564     gDirectory->mkdir(triggerClass);
565     gDirectory->cd(triggerClass);
566   
567     static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i))->SaveHistograms();
568     
569     gDirectory->cd("..");
570   }
571  
572   if (fBackgroundIdentification)
573   {
574     gDirectory->mkdir("background_identification");
575     gDirectory->cd("background_identification");
576       
577     fBackgroundIdentification->GetOutput()->Write();
578       
579     gDirectory->cd("..");
580   }
581   
582   if (folder)
583     gDirectory->cd("..");
584 }