1 /* $Id: AliPhysicsSelection.cxx 35782 2009-10-22 11:54:31Z jgrosseo $ */
3 /**************************************************************************
4 * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
6 * Author: The ALICE Off-line Project. *
7 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
26 // fPhysicsSelection = new AliPhysicsSelection;
29 // fPhysicsSelection->SetAnalyzeMC()
31 // To check if an event is a collision candidate, use:
32 // fPhysicsSelection->IsCollisionCandidate(fESD)
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")
38 // To print statistics after processing use:
39 // fPhysicsSelection->Print();
41 // Origin: Jan Fiete Grosse-Oetringhaus, CERN
42 //-------------------------------------------------------------------------
44 #include <Riostream.h>
48 #include <TIterator.h>
49 #include <TDirectory.h>
50 #include <TObjArray.h>
52 #include <AliPhysicsSelection.h>
54 #include <AliTriggerAnalysis.h>
57 #include <AliESDEvent.h>
59 ClassImp(AliPhysicsSelection)
61 AliPhysicsSelection::AliPhysicsSelection() :
62 AliAnalysisCuts("AliPhysicsSelection", "AliPhysicsSelection"),
68 fBackgroundIdentification(0),
70 fHistBunchCrossing(0),
71 fSkipTriggerClassSelection(0),
72 fUsingCustomClasses(0)
76 fCollTrigClasses.SetOwner(1);
77 fBGTrigClasses.SetOwner(1);
78 fTriggerAnalysis.SetOwner(1);
80 AliLog::SetClassDebugLevel("AliPhysicsSelection", AliLog::kWarning);
83 AliPhysicsSelection::~AliPhysicsSelection()
87 fCollTrigClasses.Delete();
88 fBGTrigClasses.Delete();
89 fTriggerAnalysis.Delete();
93 delete fHistStatistics;
97 if (fHistBunchCrossing)
99 delete fHistBunchCrossing;
100 fHistBunchCrossing = 0;
104 Bool_t AliPhysicsSelection::CheckTriggerClass(const AliESDEvent* aEsd, const char* trigger) const
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
110 TString str(trigger);
111 TObjArray* tokens = str.Tokenize(" ");
113 for (Int_t i=0; i < tokens->GetEntries(); i++)
115 TString str2(((TObjString*) tokens->At(i))->String());
117 if (str2[0] != '+' && str2[0] != '-')
118 AliFatal(Form("Invalid trigger syntax: %s", trigger));
120 Bool_t flag = (str2[0] == '+');
124 if (flag && !aEsd->IsTriggerClassFired(str2))
126 AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is not present", str2.Data()));
130 if (!flag && aEsd->IsTriggerClassFired(str2))
132 AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is present", str2.Data()));
142 Bool_t AliPhysicsSelection::IsCollisionCandidate(const AliESDEvent* aEsd)
144 // checks if the given event is a collision candidate
146 if (fCurrentRun != aEsd->GetRunNumber())
147 if (!Initialize(aEsd->GetRunNumber()))
148 AliFatal(Form("Could not initialize for run %d", aEsd->GetRunNumber()));
150 const AliESDHeader* esdHeader = aEsd->GetHeader();
153 AliError("ESD Header could not be retrieved");
157 // check event type; should be PHYSICS = 7 for data and 0 for MC
160 if (esdHeader->GetEventType() != 7)
165 if (esdHeader->GetEventType() != 0)
166 AliFatal(Form("Invalid event type for MC: %d", esdHeader->GetEventType()));
169 Bool_t accept = kFALSE;
171 Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
172 for (Int_t i=0; i < count; i++)
174 const char* triggerClass = 0;
175 if (i < fCollTrigClasses.GetEntries())
176 triggerClass = ((TObjString*) fCollTrigClasses.At(i))->String();
178 triggerClass = ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
180 AliDebug(AliLog::kDebug, Form("Processing trigger class %s", triggerClass));
182 AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
184 triggerAnalysis->FillTriggerClasses(aEsd);
186 if (CheckTriggerClass(aEsd, triggerClass))
188 triggerAnalysis->FillHistograms(aEsd);
190 fHistStatistics->Fill(1, i);
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);
199 if (fastORHW == 0 && !v0A && !v0C)
201 AliDebug(AliLog::kDebug, "Rejecting event because hardware trigger is not fired");
205 fHistStatistics->Fill(2, i);
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;
213 if (fastOROffline > 0)
214 fHistStatistics->Fill(3, i);
215 if (fastOROffline > 1)
216 fHistStatistics->Fill(4, i);
219 fHistStatistics->Fill(5, i);
221 fHistStatistics->Fill(6, i);
223 fHistStatistics->Fill(7, i);
225 fHistStatistics->Fill(8, i);
227 if (fastOROffline > 1 && !v0BG)
228 fHistStatistics->Fill(9, i);
230 if ((fastOROffline > 0 || v0A || v0C) && !v0BG)
231 fHistStatistics->Fill(10, i);
233 if (fastOROffline > 0 && (v0A || v0C) && !v0BG)
234 fHistStatistics->Fill(11, i);
236 if (v0A && v0C && !v0BG)
237 fHistStatistics->Fill(12, i);
239 if (fastOROffline > 1 || (fastOROffline > 0 && (v0A || v0C)) || (v0A && v0C))
243 fHistStatistics->Fill(13, i);
245 if (fBackgroundIdentification && !fBackgroundIdentification->IsSelected(const_cast<AliESDEvent*> (aEsd)))
247 AliDebug(AliLog::kDebug, "Rejecting event because of background identification");
248 fHistStatistics->Fill(14, i);
252 AliDebug(AliLog::kDebug, "Accepted event for histograms");
254 fHistStatistics->Fill(15, i);
255 fHistBunchCrossing->Fill(aEsd->GetBunchCrossNumber(), i);
256 if (i < fCollTrigClasses.GetEntries() || fSkipTriggerClassSelection)
261 AliDebug(AliLog::kDebug, "Rejecting event because of V0 BG flag");
264 AliDebug(AliLog::kDebug, "Rejecting event because trigger condition is not fulfilled");
269 AliDebug(AliLog::kDebug, "Accepted event as collision candidate");
274 Int_t AliPhysicsSelection::GetTriggerScheme(UInt_t runNumber)
276 // returns the current trigger scheme (classes that are accepted/rejected)
281 // TODO dependent on run number
292 // default: CINT1 suite
296 Bool_t AliPhysicsSelection::Initialize(UInt_t runNumber)
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
301 Bool_t oldStatus = TH1::AddDirectoryStatus();
302 TH1::AddDirectory(kFALSE);
304 Int_t triggerScheme = GetTriggerScheme(runNumber);
306 if (!fUsingCustomClasses && fCurrentRun != -1 && triggerScheme != GetTriggerScheme(fCurrentRun))
307 AliFatal("Processing several runs with different trigger schemes is not supported");
309 AliInfo(Form("Initializing for run %d", runNumber));
311 // initialize first time?
312 if (fCurrentRun == -1)
314 if (fUsingCustomClasses) {
315 AliInfo("Using user-provided trigger classes");
317 switch (triggerScheme)
320 fCollTrigClasses.Add(new TObjString(""));
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"));
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"));
337 AliFatal(Form("Unsupported trigger scheme %d", triggerScheme));
341 Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
343 for (Int_t i=0; i<count; i++)
345 AliTriggerAnalysis* triggerAnalysis = new AliTriggerAnalysis;
346 triggerAnalysis->SetAnalyzeMC(fMC);
347 triggerAnalysis->EnableHistograms();
348 triggerAnalysis->SetSPDGFOThreshhold(1);
349 fTriggerAnalysis.Add(triggerAnalysis);
353 delete fHistStatistics;
355 fHistStatistics = new TH2F("fHistStatistics", "fHistStatistics;;", 15, 0.5, 15.5, count, -0.5, -0.5 + count);
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");
374 if (fHistBunchCrossing)
375 delete fHistBunchCrossing;
377 fHistBunchCrossing = new TH2F("fHistBunchCrossing", "fHistBunchCrossing;bunch crossing number;", 4000, -0.5, 3999.5, count, -0.5, -0.5 + count);
380 for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
382 fHistStatistics->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
383 fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
386 for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
388 fHistStatistics->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
389 fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
394 Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
395 for (Int_t i=0; i<count; i++)
397 AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
406 triggerAnalysis->SetV0TimeOffset(7.5);
409 triggerAnalysis->SetV0TimeOffset(0);
413 fCurrentRun = runNumber;
415 TH1::AddDirectory(oldStatus);
420 void AliPhysicsSelection::Print(Option_t* /* option */) const
422 // print the configuration
424 Printf("Configuration initialized for run %d (MC: %d):", fCurrentRun, fMC);
426 Printf("Collision trigger classes:");
427 for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
428 Printf("%s", ((TObjString*) fCollTrigClasses.At(i))->String().Data());
430 Printf("Background trigger classes:");
431 for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
432 Printf("%s", ((TObjString*) fBGTrigClasses.At(i))->String().Data());
434 AliTriggerAnalysis* triggerAnalysis = dynamic_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(0));
438 if (triggerAnalysis->GetV0TimeOffset() > 0)
439 Printf("V0 time offset active: %.2f ns", triggerAnalysis->GetV0TimeOffset());
441 Printf("\nTotal available events:");
443 triggerAnalysis->PrintTriggerClasses();
446 if (fHistStatistics && fCollTrigClasses.GetEntries() > 0)
448 Printf("\nSelection statistics for first collision trigger (%s):", ((TObjString*) fCollTrigClasses.First())->String().Data());
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));
454 if (fHistBunchCrossing)
456 Printf("\nBunch crossing statistics:");
458 for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
461 str.Form("Trigger %s has accepted events in the bunch crossings: ", fHistBunchCrossing->GetYaxis()->GetBinLabel(i));
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));
467 Printf("%s", str.Data());
470 for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
473 for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
475 if (fHistBunchCrossing->GetBinContent(j, i) > 0)
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));
483 if (fUsingCustomClasses)
484 Printf("WARNING: Using custom trigger classes!");
485 if (fSkipTriggerClassSelection)
486 Printf("WARNING: Skipping trigger class selection!");
489 Long64_t AliPhysicsSelection::Merge(TCollection* list)
491 // Merge a list of AliMultiplicityCorrection objects with this (needed for
493 // Returns the number of merged objects (including this).
501 TIterator* iter = list->MakeIterator();
504 // collections of all histograms
505 const Int_t nHists = 9;
506 TList collections[nHists];
509 while ((obj = iter->Next())) {
511 AliPhysicsSelection* entry = dynamic_cast<AliPhysicsSelection*> (obj);
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);
526 fTriggerAnalysis.Merge(&collections[0]);
528 fHistStatistics->Merge(&collections[1]);
529 if (fHistBunchCrossing)
530 fHistBunchCrossing->Merge(&collections[2]);
531 if (fBackgroundIdentification)
532 fBackgroundIdentification->Merge(&collections[3]);
539 void AliPhysicsSelection::SaveHistograms(const char* folder) const
541 // write histograms to current directory
543 if (!fHistStatistics)
548 gDirectory->mkdir(folder);
549 gDirectory->cd(folder);
552 fHistStatistics->Write();
553 fHistBunchCrossing->Write();
555 Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
556 for (Int_t i=0; i < count; i++)
558 TString triggerClass = "trigger_histograms_";
559 if (i < fCollTrigClasses.GetEntries())
560 triggerClass += ((TObjString*) fCollTrigClasses.At(i))->String();
562 triggerClass += ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
564 gDirectory->mkdir(triggerClass);
565 gDirectory->cd(triggerClass);
567 static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i))->SaveHistograms();
569 gDirectory->cd("..");
572 if (fBackgroundIdentification)
574 gDirectory->mkdir("background_identification");
575 gDirectory->cd("background_identification");
577 fBackgroundIdentification->GetOutput()->Write();
579 gDirectory->cd("..");
583 gDirectory->cd("..");