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),
74 fCollTrigClasses.SetOwner(1);
75 fBGTrigClasses.SetOwner(1);
76 fTriggerAnalysis.SetOwner(1);
78 AliLog::SetClassDebugLevel("AliPhysicsSelection", AliLog::kWarning);
81 AliPhysicsSelection::~AliPhysicsSelection()
85 fCollTrigClasses.Delete();
86 fBGTrigClasses.Delete();
87 fTriggerAnalysis.Delete();
91 delete fHistStatistics;
95 if (fHistBunchCrossing)
97 delete fHistBunchCrossing;
98 fHistBunchCrossing = 0;
102 Bool_t AliPhysicsSelection::CheckTriggerClass(const AliESDEvent* aEsd, const char* trigger) const
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
108 TString str(trigger);
109 TObjArray* tokens = str.Tokenize(" ");
111 for (Int_t i=0; i < tokens->GetEntries(); i++)
113 TString str2(((TObjString*) tokens->At(i))->String());
115 if (str2[0] != '+' && str2[0] != '-')
116 AliFatal(Form("Invalid trigger syntax: %s", trigger));
118 Bool_t flag = (str2[0] == '+');
122 if (flag && !aEsd->IsTriggerClassFired(str2))
124 AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is not present", str2.Data()));
128 if (!flag && aEsd->IsTriggerClassFired(str2))
130 AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is present", str2.Data()));
140 Bool_t AliPhysicsSelection::IsCollisionCandidate(const AliESDEvent* aEsd)
142 // checks if the given event is a collision candidate
144 if (fCurrentRun != aEsd->GetRunNumber())
145 if (!Initialize(aEsd->GetRunNumber()))
146 AliFatal(Form("Could not initialize for run %d", aEsd->GetRunNumber()));
148 const AliESDHeader* esdHeader = aEsd->GetHeader();
151 AliError("ESD Header could not be retrieved");
155 // check event type; should be PHYSICS = 7 for data and 0 for MC
158 if (esdHeader->GetEventType() != 7)
163 if (esdHeader->GetEventType() != 0)
164 AliFatal(Form("Invalid event type for MC: %d", esdHeader->GetEventType()));
167 Bool_t accept = kFALSE;
169 Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
170 for (Int_t i=0; i < count; i++)
172 const char* triggerClass = 0;
173 if (i < fCollTrigClasses.GetEntries())
174 triggerClass = ((TObjString*) fCollTrigClasses.At(i))->String();
176 triggerClass = ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
178 AliDebug(AliLog::kDebug, Form("Processing trigger class %s", triggerClass));
180 AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
182 triggerAnalysis->FillTriggerClasses(aEsd);
184 if (CheckTriggerClass(aEsd, triggerClass))
186 triggerAnalysis->FillHistograms(aEsd);
188 fHistStatistics->Fill(1, i);
190 // hardware trigger (should only remove events for MC)
191 // replay CINT1B hardware trigger
192 // TODO this has to depend on the actual hardware trigger (and that depends on the run...)
193 Int_t fastORHW = triggerAnalysis->SPDFiredChips(aEsd, 1); // SPD number of chips from trigger bits (!)
194 Bool_t v0A = triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0A);
195 Bool_t v0C = triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0C);
197 if (fastORHW == 0 && !v0A && !v0C)
199 AliDebug(AliLog::kDebug, "Rejecting event because hardware trigger is not fired");
203 fHistStatistics->Fill(2, i);
206 Int_t fastOROffline = triggerAnalysis->SPDFiredChips(aEsd, 0); // SPD number of chips from clusters (!)
207 Bool_t v0ABG = triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0ABG);
208 Bool_t v0CBG = triggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0CBG);
209 Bool_t v0BG = v0ABG || v0CBG;
211 if (fastOROffline > 0)
212 fHistStatistics->Fill(3, i);
213 if (fastOROffline > 1)
214 fHistStatistics->Fill(4, i);
217 fHistStatistics->Fill(5, i);
219 fHistStatistics->Fill(6, i);
221 fHistStatistics->Fill(7, i);
223 fHistStatistics->Fill(8, i);
225 if (fastOROffline > 1 && !v0BG)
226 fHistStatistics->Fill(9, i);
228 if ((fastOROffline > 0 || v0A || v0C) && !v0BG)
229 fHistStatistics->Fill(10, i);
231 if (fastOROffline > 0 && (v0A || v0C) && !v0BG)
232 fHistStatistics->Fill(11, i);
234 if (v0A && v0C && !v0BG)
235 fHistStatistics->Fill(12, i);
237 if (fastOROffline > 1 || (fastOROffline > 0 && (v0A || v0C)) || (v0A && v0C))
241 fHistStatistics->Fill(13, i);
243 if (fBackgroundIdentification && !fBackgroundIdentification->IsSelected(const_cast<AliESDEvent*> (aEsd)))
245 AliDebug(AliLog::kDebug, "Rejecting event because of background identification");
246 fHistStatistics->Fill(14, i);
250 AliDebug(AliLog::kDebug, "Accepted event for histograms");
252 fHistStatistics->Fill(15, i);
253 fHistBunchCrossing->Fill(aEsd->GetBunchCrossNumber(), i);
254 if (i < fCollTrigClasses.GetEntries())
259 AliDebug(AliLog::kDebug, "Rejecting event because of V0 BG flag");
262 AliDebug(AliLog::kDebug, "Rejecting event because trigger condition is not fulfilled");
267 AliDebug(AliLog::kDebug, "Accepted event as collision candidate");
272 Int_t AliPhysicsSelection::GetTriggerScheme(UInt_t runNumber)
274 // returns the current trigger scheme (classes that are accepted/rejected)
279 // TODO dependent on run number
290 // default: CINT1 suite
294 Bool_t AliPhysicsSelection::Initialize(UInt_t runNumber)
296 // initializes the object for the given run
297 // 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
299 Bool_t oldStatus = TH1::AddDirectoryStatus();
300 TH1::AddDirectory(kFALSE);
302 Int_t triggerScheme = GetTriggerScheme(runNumber);
304 if (fCurrentRun != -1 && triggerScheme != GetTriggerScheme(fCurrentRun))
305 AliFatal("Processing several runs with different trigger schemes is not supported");
307 AliInfo(Form("Initializing for run %d", runNumber));
309 Bool_t firstTime = kFALSE;
310 if (fCurrentRun == -1)
312 // initialize first time
313 fCurrentRun = runNumber;
319 switch (triggerScheme)
322 fCollTrigClasses.Add(new TObjString(""));
326 fCollTrigClasses.Add(new TObjString("+CINT1B-ABCE-NOPF-ALL"));
327 fBGTrigClasses.Add(new TObjString("+CINT1A-ABCE-NOPF-ALL"));
328 fBGTrigClasses.Add(new TObjString("+CINT1C-ABCE-NOPF-ALL"));
329 fBGTrigClasses.Add(new TObjString("+CINT1-E-NOPF-ALL"));
333 fCollTrigClasses.Add(new TObjString("+CSMBB-ABCE-NOPF-ALL"));
334 fBGTrigClasses.Add(new TObjString("+CSMBA-ABCE-NOPF-ALL -CSMBB-ABCE-NOPF-ALL"));
335 fBGTrigClasses.Add(new TObjString("+CSMBC-ABCE-NOPF-ALL -CSMBB-ABCE-NOPF-ALL"));
339 AliFatal(Form("Unsupported trigger scheme %d", triggerScheme));
342 Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
344 for (Int_t i=0; i<count; i++)
346 AliTriggerAnalysis* triggerAnalysis = new AliTriggerAnalysis;
347 triggerAnalysis->SetAnalyzeMC(fMC);
348 triggerAnalysis->EnableHistograms();
349 triggerAnalysis->SetSPDGFOThreshhold(1);
350 fTriggerAnalysis.Add(triggerAnalysis);
354 delete fHistStatistics;
356 fHistStatistics = new TH2F("fHistStatistics", "fHistStatistics;;", 15, 0.5, 15.5, count, -0.5, -0.5 + count);
359 fHistStatistics->GetXaxis()->SetBinLabel(n++, "Trigger class");
360 fHistStatistics->GetXaxis()->SetBinLabel(n++, "Hardware trigger");
361 fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 1");
362 fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 2");
363 fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0A");
364 fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0C");
365 fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0A BG");
366 fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0C BG");
367 fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 2 &!V0 BG");
368 fHistStatistics->GetXaxis()->SetBinLabel(n++, "(FO >= 1 | V0A | V0C) & !V0 BG");
369 fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO >= 1 & (V0A | V0C) & !V0 BG");
370 fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0A & V0C & !V0 BG");
371 fHistStatistics->GetXaxis()->SetBinLabel(n++, "(FO >= 2 | (FO >= 1 & (V0A | V0C)) | (V0A & V0C)) & !V0 BG");
372 fHistStatistics->GetXaxis()->SetBinLabel(n++, "Background identification");
373 fHistStatistics->GetXaxis()->SetBinLabel(n++, "Accepted");
375 if (fHistBunchCrossing)
376 delete fHistBunchCrossing;
378 fHistBunchCrossing = new TH2F("fHistBunchCrossing", "fHistBunchCrossing;bunch crossing number;", 4000, -0.5, 3999.5, count, -0.5, -0.5 + count);
381 for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
383 fHistStatistics->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
384 fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fCollTrigClasses.At(i))->String());
387 for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
389 fHistStatistics->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
390 fHistBunchCrossing->GetYaxis()->SetBinLabel(n, ((TObjString*) fBGTrigClasses.At(i))->String());
395 Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
396 for (Int_t i=0; i<count; i++)
398 AliTriggerAnalysis* triggerAnalysis = static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i));
407 triggerAnalysis->SetV0TimeOffset(7.5);
410 triggerAnalysis->SetV0TimeOffset(0);
414 TH1::AddDirectory(oldStatus);
419 void AliPhysicsSelection::Print(Option_t* /* option */) const
421 // print the configuration
423 Printf("Configuration initialized for run %d (MC: %d):", fCurrentRun, fMC);
425 Printf("Collision trigger classes:");
426 for (Int_t i=0; i < fCollTrigClasses.GetEntries(); i++)
427 Printf("%s", ((TObjString*) fCollTrigClasses.At(i))->String().Data());
429 Printf("Background trigger classes:");
430 for (Int_t i=0; i < fBGTrigClasses.GetEntries(); i++)
431 Printf("%s", ((TObjString*) fBGTrigClasses.At(i))->String().Data());
433 AliTriggerAnalysis* triggerAnalysis = dynamic_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(0));
437 if (triggerAnalysis->GetV0TimeOffset() > 0)
438 Printf("V0 time offset active: %.2f ns", triggerAnalysis->GetV0TimeOffset());
440 Printf("\nTotal available events:");
442 triggerAnalysis->PrintTriggerClasses();
445 if (fHistStatistics && fCollTrigClasses.GetEntries() > 0)
447 Printf("\nSelection statistics for first collision trigger (%s):", ((TObjString*) fCollTrigClasses.First())->String().Data());
449 Printf("Total events with correct trigger class: %d", (Int_t) fHistStatistics->GetBinContent(1, 1));
450 Printf("Selected collision candidates: %d", (Int_t) fHistStatistics->GetBinContent(fHistStatistics->GetXaxis()->FindBin("Accepted"), 1));
453 if (fHistBunchCrossing)
455 Printf("\nBunch crossing statistics:");
457 for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
460 str.Form("Trigger %s has accepted events in the bunch crossings: ", fHistBunchCrossing->GetYaxis()->GetBinLabel(i));
462 for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
463 if (fHistBunchCrossing->GetBinContent(j, i) > 0)
464 str += Form("%d, ", (Int_t) fHistBunchCrossing->GetXaxis()->GetBinCenter(j));
466 Printf("%s", str.Data());
469 for (Int_t j=1; j<=fHistBunchCrossing->GetNbinsX(); j++)
472 for (Int_t i=1; i<=fHistBunchCrossing->GetNbinsY(); i++)
474 if (fHistBunchCrossing->GetBinContent(j, i) > 0)
478 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 Long64_t AliPhysicsSelection::Merge(TCollection* list)
485 // Merge a list of AliMultiplicityCorrection objects with this (needed for
487 // Returns the number of merged objects (including this).
495 TIterator* iter = list->MakeIterator();
498 // collections of all histograms
499 const Int_t nHists = 9;
500 TList collections[nHists];
503 while ((obj = iter->Next())) {
505 AliPhysicsSelection* entry = dynamic_cast<AliPhysicsSelection*> (obj);
509 collections[0].Add(&(entry->fTriggerAnalysis));
510 if (entry->fHistStatistics)
511 collections[1].Add(entry->fHistStatistics);
512 if (entry->fHistBunchCrossing)
513 collections[2].Add(entry->fHistBunchCrossing);
514 if (entry->fBackgroundIdentification)
515 collections[3].Add(entry->fBackgroundIdentification);
520 fTriggerAnalysis.Merge(&collections[0]);
522 fHistStatistics->Merge(&collections[1]);
523 if (fHistBunchCrossing)
524 fHistBunchCrossing->Merge(&collections[2]);
525 if (fBackgroundIdentification)
526 fBackgroundIdentification->Merge(&collections[3]);
533 void AliPhysicsSelection::SaveHistograms(const char* folder) const
535 // write histograms to current directory
537 if (!fHistStatistics)
542 gDirectory->mkdir(folder);
543 gDirectory->cd(folder);
546 fHistStatistics->Write();
547 fHistBunchCrossing->Write();
549 Int_t count = fCollTrigClasses.GetEntries() + fBGTrigClasses.GetEntries();
550 for (Int_t i=0; i < count; i++)
552 TString triggerClass = "trigger_histograms_";
553 if (i < fCollTrigClasses.GetEntries())
554 triggerClass += ((TObjString*) fCollTrigClasses.At(i))->String();
556 triggerClass += ((TObjString*) fBGTrigClasses.At(i - fCollTrigClasses.GetEntries()))->String();
558 gDirectory->mkdir(triggerClass);
559 gDirectory->cd(triggerClass);
561 static_cast<AliTriggerAnalysis*> (fTriggerAnalysis.At(i))->SaveHistograms();
563 gDirectory->cd("..");
566 if (fBackgroundIdentification)
568 gDirectory->mkdir("background_identification");
569 gDirectory->cd("background_identification");
571 fBackgroundIdentification->GetOutput()->Write();
573 gDirectory->cd("..");
577 gDirectory->cd("..");