1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
17 // Filling an AliCFContainer with the quantities pt, eta and phi
18 // for tracks which survivied the particle cuts (MC resp. ESD tracks)
19 // Track selection is done using the AliHFE package
22 // Raphaelle Bailhache <R.Bailhache@gsi.de>
23 // Markus Fasel <M.Fasel@gsi.de>
24 // Matus Kalisky <matus.kalisky@cern.ch>
25 // MinJung Kweon <minjung@physi.uni-heidelberg.de>
31 #include <TDirectory.h>
34 #include <TIterator.h>
38 #include <TObjArray.h>
39 #include <TObjString.h>
40 #include <TParticle.h>
46 #include "AliAnalysisManager.h"
47 #include "AliAODInputHandler.h"
48 #include "AliAODMCParticle.h"
49 #include "AliAODTrack.h"
50 #include "AliAODVertex.h"
51 #include "AliCentrality.h"
52 #include "AliCFContainer.h"
53 #include "AliCFManager.h"
54 #include "AliESDEvent.h"
55 #include "AliESDInputHandler.h"
56 #include "AliESDtrack.h"
58 #include "AliMCEvent.h"
59 #include "AliMCEventHandler.h"
60 #include "AliMCParticle.h"
61 #include "AliMultiplicity.h"
63 #include "AliPIDResponse.h"
64 #include "AliOADBContainer.h"
66 #include "AliTriggerAnalysis.h"
67 #include "AliVVertex.h"
69 #include "AliHFEcollection.h"
70 #include "AliHFEcontainer.h"
71 #include "AliHFEcuts.h"
72 #include "AliHFEelecbackground.h"
73 #include "AliHFEmcQA.h"
74 #include "AliHFEpairs.h"
75 #include "AliHFEpid.h"
76 #include "AliHFEpidQAmanager.h"
77 #include "AliHFEpostAnalysis.h"
78 #include "AliHFEsecVtxs.h"
79 #include "AliHFEsecVtx.h"
80 #include "AliHFEsignalCuts.h"
81 #include "AliHFEtaggedTrackAnalysis.h"
82 #include "AliHFEtools.h"
83 #include "AliHFEvarManager.h"
84 #include "AliAnalysisTaskHFE.h"
86 ClassImp(AliAnalysisTaskHFE)
88 //____________________________________________________________
89 AliAnalysisTaskHFE::AliAnalysisTaskHFE():
90 AliAnalysisTaskSE("PID efficiency Analysis")
93 , fFillSignalOnly(kTRUE)
96 , fApplyCutAOD(kFALSE)
98 , fBackGroundFactorApply(kFALSE)
99 , fRemovePileUp(kFALSE)
100 , fIdentifiedAsPileUp(kFALSE)
101 , fIdentifiedAsOutInz(kFALSE)
102 , fPassTheEventCut(kFALSE)
103 , fRejectKinkMother(kTRUE)
104 , fisppMultiBin(kFALSE)
105 , fPbPbUserCentralityBinning(kFALSE)
106 , fisNonHFEsystematics(kFALSE)
107 , fSpecialTrigger(NULL)
109 , fCentralityPercent(-1)
111 , fWeightBackGround(0.)
113 , fHadronBackgroundOADB(NULL)
118 , fTriggerAnalysis(NULL)
121 , fPIDpreselect(NULL)
123 , fTaggedTrackCuts(NULL)
124 , fCleanTaggedTrack(kFALSE)
125 , fVariablesTRDTaggedTrack(kFALSE)
126 , fCutspreselect(NULL)
128 , fElecBackGround(NULL)
130 , fTaggedTrackAnalysis(NULL)
136 , fHistELECBACKGROUND(NULL)
137 , fQACollection(NULL)
142 memset(fElecBackgroundFactor, 0, sizeof(Double_t) * kElecBgSpecies * kBgPtBins * kCentBins * kBgLevels);
143 memset(fkBackGroundFactorArray, 0, sizeof(TF1 *) * 12);
144 memset(fBinLimit, 0, sizeof(Double_t) * (kBgPtBins+1));
145 memset(&fisppMultiBin, kFALSE, sizeof(fisppMultiBin));
146 memset(fCentralityLimits, 0, sizeof(Float_t) * 12);
151 //____________________________________________________________
152 AliAnalysisTaskHFE::AliAnalysisTaskHFE(const char * name):
153 AliAnalysisTaskSE(name)
156 , fFillSignalOnly(kTRUE)
157 , fFillNoCuts(kFALSE)
159 , fApplyCutAOD(kFALSE)
161 , fBackGroundFactorApply(kFALSE)
162 , fRemovePileUp(kFALSE)
163 , fIdentifiedAsPileUp(kFALSE)
164 , fIdentifiedAsOutInz(kFALSE)
165 , fPassTheEventCut(kFALSE)
166 , fRejectKinkMother(kTRUE)
167 , fisppMultiBin(kFALSE)
168 , fPbPbUserCentralityBinning(kFALSE)
169 , fisNonHFEsystematics(kFALSE)
170 , fSpecialTrigger(NULL)
172 , fCentralityPercent(-1)
174 , fWeightBackGround(0.)
176 , fHadronBackgroundOADB(NULL)
181 , fTriggerAnalysis(NULL)
184 , fPIDpreselect(NULL)
186 , fTaggedTrackCuts(NULL)
187 , fCleanTaggedTrack(kFALSE)
188 , fVariablesTRDTaggedTrack(kFALSE)
189 , fCutspreselect(NULL)
191 , fElecBackGround(NULL)
193 , fTaggedTrackAnalysis(NULL)
199 , fHistELECBACKGROUND(NULL)
203 // Default constructor
205 DefineOutput(1, TList::Class());
206 DefineOutput(2, TList::Class());
208 fPID = new AliHFEpid("hfePid");
209 fPIDqa = new AliHFEpidQAmanager;
210 fVarManager = new AliHFEvarManager("hfeVarManager");
212 memset(fElecBackgroundFactor, 0, sizeof(Double_t) * kElecBgSpecies * kBgPtBins * kCentBins * kBgLevels);
213 memset(fkBackGroundFactorArray, 0, sizeof(TF1 *) * 12);
214 memset(fBinLimit, 0, sizeof(Double_t) * (kBgPtBins+1));
215 memset(&fisppMultiBin, kFALSE, sizeof(fisppMultiBin));
216 memset(fCentralityLimits, 0, sizeof(Float_t) * 12);
220 //____________________________________________________________
221 AliAnalysisTaskHFE::AliAnalysisTaskHFE(const AliAnalysisTaskHFE &ref):
222 AliAnalysisTaskSE(ref)
225 , fFillSignalOnly(ref.fFillSignalOnly)
226 , fFillNoCuts(ref.fFillNoCuts)
227 , fUseFlagAOD(ref.fUseFlagAOD)
228 , fApplyCutAOD(ref.fApplyCutAOD)
230 , fBackGroundFactorApply(ref.fBackGroundFactorApply)
231 , fRemovePileUp(ref.fRemovePileUp)
232 , fIdentifiedAsPileUp(ref.fIdentifiedAsPileUp)
233 , fIdentifiedAsOutInz(ref.fIdentifiedAsOutInz)
234 , fPassTheEventCut(ref.fPassTheEventCut)
235 , fRejectKinkMother(ref.fRejectKinkMother)
236 , fisppMultiBin(ref.fisppMultiBin)
237 , fPbPbUserCentralityBinning(ref.fPbPbUserCentralityBinning)
238 , fisNonHFEsystematics(ref.fisNonHFEsystematics)
239 , fSpecialTrigger(ref.fSpecialTrigger)
240 , fCentralityF(ref.fCentralityF)
241 , fCentralityPercent(ref.fCentralityPercent)
242 , fContributors(ref.fContributors)
243 , fWeightBackGround(ref.fWeightBackGround)
245 , fHadronBackgroundOADB(ref.fHadronBackgroundOADB)
250 , fTriggerAnalysis(NULL)
253 , fPIDpreselect(NULL)
255 , fTaggedTrackCuts(NULL)
256 , fCleanTaggedTrack(ref.fCleanTaggedTrack)
257 , fVariablesTRDTaggedTrack(ref.fVariablesTRDTaggedTrack)
258 , fCutspreselect(NULL)
260 , fElecBackGround(NULL)
262 , fTaggedTrackAnalysis(NULL)
268 , fHistELECBACKGROUND(NULL)
269 , fQACollection(NULL)
277 //____________________________________________________________
278 AliAnalysisTaskHFE &AliAnalysisTaskHFE::operator=(const AliAnalysisTaskHFE &ref){
280 // Assignment operator
287 //____________________________________________________________
288 void AliAnalysisTaskHFE::Copy(TObject &o) const {
290 // Copy into object o
292 AliAnalysisTaskHFE &target = dynamic_cast<AliAnalysisTaskHFE &>(o);
293 target.fQAlevel = fQAlevel;
294 target.fPlugins = fPlugins;
295 target.fFillSignalOnly = fFillSignalOnly;
296 target.fFillNoCuts = fFillNoCuts;
297 target.fUseFlagAOD = fUseFlagAOD;
298 target.fApplyCutAOD = fApplyCutAOD;
299 target.fFlags = fFlags;
300 target.fBackGroundFactorApply = fBackGroundFactorApply;
301 target.fRemovePileUp = fRemovePileUp;
302 target.fIdentifiedAsPileUp = fIdentifiedAsPileUp;
303 target.fIdentifiedAsOutInz = fIdentifiedAsOutInz;
304 target.fPassTheEventCut = fPassTheEventCut;
305 target.fRejectKinkMother = fRejectKinkMother;
306 target.fisppMultiBin = fisppMultiBin;
307 target.fPbPbUserCentralityBinning = fPbPbUserCentralityBinning;
308 target.fisNonHFEsystematics = fisNonHFEsystematics;
309 target.fSpecialTrigger = fSpecialTrigger;
310 target.fCentralityF = fCentralityF;
311 target.fCentralityPercent = fCentralityPercent;
312 target.fContributors = fContributors;
313 target.fWeightBackGround = fWeightBackGround;
315 target.fHadronBackgroundOADB = fHadronBackgroundOADB;
316 target.fContainer = fContainer;
317 target.fVarManager = fVarManager;
318 target.fSignalCuts = fSignalCuts;
320 target.fTriggerAnalysis = fTriggerAnalysis;
322 target.fPIDqa = fPIDqa;
323 target.fPIDpreselect = fPIDpreselect;
324 target.fCuts = fCuts;
325 target.fTaggedTrackCuts = fTaggedTrackCuts;
326 target.fCleanTaggedTrack = fCleanTaggedTrack;
327 target.fVariablesTRDTaggedTrack = fVariablesTRDTaggedTrack;
328 target.fCutspreselect = fCutspreselect;
329 target.fSecVtx = fSecVtx;
330 target.fElecBackGround = fElecBackGround;
331 target.fMCQA = fMCQA;
332 target.fTaggedTrackAnalysis = fTaggedTrackAnalysis;
333 target.fExtraCuts = fExtraCuts;
335 target.fOutput = fOutput;
336 target.fHistMCQA = fHistMCQA;
337 target.fHistSECVTX = fHistSECVTX;
338 target.fHistELECBACKGROUND = fHistELECBACKGROUND;
339 target.fQACollection = fQACollection;
342 //____________________________________________________________
343 AliAnalysisTaskHFE::~AliAnalysisTaskHFE(){
347 if(fPID) delete fPID;
348 if(fPIDpreselect) delete fPIDpreselect;
349 if(fVarManager) delete fVarManager;
350 if(fCFM) delete fCFM;
351 if(fTriggerAnalysis) delete fTriggerAnalysis;
352 if(fSignalCuts) delete fSignalCuts;
353 if(fSecVtx) delete fSecVtx;
354 if(fMCQA) delete fMCQA;
355 if(fElecBackGround) delete fElecBackGround;
356 if(fSpecialTrigger) delete fSpecialTrigger;
357 // Delete output objects only if we are not running in PROOF mode because otherwise this produces a crash during merging
358 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
359 if(mgr && mgr->GetAnalysisType() != AliAnalysisManager::kProofAnalysis){
360 if(fPIDqa) delete fPIDqa;
361 if(fOutput) delete fOutput;
366 //____________________________________________________________
367 void AliAnalysisTaskHFE::UserCreateOutputObjects(){
369 // Creating output container and output objects
370 // Here we also Initialize the correction framework container and
375 // QA histograms are created if requested
376 // Called once per worker
378 AliDebug(3, "Creating Output Objects");
379 // Automatic determination of the analysis mode
380 AliVEventHandler *inputHandler = dynamic_cast<AliVEventHandler *>(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
381 if(!TString(inputHandler->IsA()->GetName()).CompareTo("AliAODInputHandler")){
385 if(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())
388 printf("Analysis Mode: %s Analysis\n", IsAODanalysis() ? "AOD" : "ESD");
389 printf("MC Data available %s\n", HasMCData() ? "Yes" : "No");
391 // Enable Trigger Analysis
392 fTriggerAnalysis = new AliTriggerAnalysis;
393 fTriggerAnalysis->EnableHistograms();
394 fTriggerAnalysis->SetAnalyzeMC(HasMCData());
397 // Make lists for Output
398 if(!fQA) fQA = new TList;
400 if(!fOutput) fOutput = new TList;
403 // First Part: Make QA histograms
404 fQACollection = new AliHFEcollection("TaskQA", "QA histos from the Electron Task");
405 fQACollection->CreateTH1F("nElectronTracksEvent", "Number of Electron Candidates", 100, 0, 100);
406 fQACollection->CreateTH1F("nElectron", "Number of electrons", 100, 0, 100);
407 fQACollection->CreateTH2F("radius", "Production Vertex", 100, 0.0, 5.0, 100, 0.0, 5.0);
408 fQACollection->CreateTH2F("TPCdEdxBeforePID", "TPC dE/dx; p (GeV/c); TPC dE/dx (a.u.)", 1000, 0., 10., 200, 0., 200.);
409 fQACollection->CreateTH2F("TPCnSigmaBeforePID", "TPC dE/dx; p (GeV/c); TPC dE/dx - <TPC dE/dx>|_{el} (#sigma)", 1000, 0., 10., 1000, -10., 10.);
411 InitPIDperformanceQA();
412 InitContaminationQA();
413 InitHistoITScluster();
414 fQA->Add(fQACollection);
417 fPID->SetHasMCData(HasMCData());
418 if(!fPID->GetNumberOfPIDdetectors()) fPID->AddDetector("TPC", 0);
420 AliInfo("PID QA switched on");
421 fPIDqa->Initialize(fPID);
422 fQA->Add(fPIDqa->MakeList("HFEpidQA"));
424 fPID->SortDetectors();
426 // Initialize correction Framework and Cuts
427 const Int_t kNcutSteps = AliHFEcuts::kNcutStepsMCTrack + AliHFEcuts::kNcutStepsRecTrack + AliHFEcuts::kNcutStepsDETrack + AliHFEcuts::kNcutStepsSecvtxTrack;
428 fCFM = new AliCFManager;
429 fCFM->SetNStepParticle(kNcutSteps);
430 MakeParticleContainer();
431 MakeEventContainer();
432 // Temporary fix: Initialize particle cuts with NULL
433 for(Int_t istep = 0; istep < kNcutSteps; istep++)
434 fCFM->SetParticleCutsList(istep, NULL);
436 AliWarning("Cuts not available. Default cuts will be used");
437 fCuts = new AliHFEcuts;
438 fCuts->CreateStandardCuts();
440 if(IsAODanalysis()) fCuts->SetAOD();
441 // Make clone for V0 tagging step
442 fCuts->Initialize(fCFM);
443 if(fCuts->IsQAOn()) fQA->Add(fCuts->GetQAhistograms());
444 fSignalCuts = new AliHFEsignalCuts("HFEsignalCuts", "HFE MC Signal definition");
445 fVarManager->SetSignalCuts(fSignalCuts);
447 // add output objects to the List
448 fOutput->AddAt(fContainer, 0);
449 fOutput->AddAt(fCFM->GetEventContainer(), 1);
451 // mcQA----------------------------------
452 if (HasMCData() && IsQAOn(kMCqa)) {
454 if(!fMCQA) fMCQA = new AliHFEmcQA;
455 if(!fHistMCQA) fHistMCQA = new TList();
456 fHistMCQA->SetOwner();
457 if(IsPbPb()) fMCQA->SetPbPb();
458 if(fisppMultiBin) fMCQA->SetPPMultiBin();
459 if(TestBit(kTreeStream)){
460 fMCQA->EnableDebugStreamer();
462 fMCQA->CreatDefaultHistograms(fHistMCQA);
463 fMCQA->SetBackgroundWeightFactor(fElecBackgroundFactor[0][0][0],fBinLimit);
467 // secvtx----------------------------------
468 if (GetPlugin(kSecVtx)) {
469 AliInfo("Secondary Vertex Analysis on");
470 if(!fSecVtx) fSecVtx = new AliHFEsecVtx;
471 fSecVtx->SetHasMCData(HasMCData());
473 if(!fHistSECVTX) fHistSECVTX = new TList();
474 fHistSECVTX->SetOwner();
475 fSecVtx->CreateHistograms(fHistSECVTX);
476 fOutput->Add(fHistSECVTX);
479 // background----------------------------------
480 if (GetPlugin(kIsElecBackGround)) {
481 AliInfo("Electron BackGround Analysis on");
482 if(!fElecBackGround){
483 AliWarning("ElecBackGround not available. Default elecbackground will be used");
484 fElecBackGround = new AliHFEelecbackground;
486 fElecBackGround->SetHasMCData(HasMCData());
488 if(!fHistELECBACKGROUND) fHistELECBACKGROUND = new TList();
489 fHistELECBACKGROUND->SetOwner();
490 fElecBackGround->CreateHistograms(fHistELECBACKGROUND);
491 fOutput->Add(fHistELECBACKGROUND);
495 if(GetPlugin(kTaggedTrackAnalysis)){
496 AliInfo("Analysis on V0-tagged tracks enabled");
497 fTaggedTrackAnalysis = new AliHFEtaggedTrackAnalysis(Form("taggedTrackAnalysis%s", GetName()));
498 fTaggedTrackAnalysis->SetCuts(fTaggedTrackCuts);
499 fTaggedTrackAnalysis->SetClean(fCleanTaggedTrack);
500 AliHFEvarManager *varManager = fTaggedTrackAnalysis->GetVarManager();
501 TObjArray *array = fVarManager->GetVariables();
502 Int_t nvars = array->GetEntriesFast();
504 for(Int_t v = 0; v < nvars; v++) {
505 AliHFEvarManager::AliHFEvariable *variable = (AliHFEvarManager::AliHFEvariable *) array->At(v);
506 if(!variable) continue;
507 TString name(((AliHFEvarManager::AliHFEvariable *)variable)->GetName());
508 if(!name.CompareTo("source")) namee = TString("species");
509 else namee = TString(name);
510 Int_t nbins = variable->GetNumberOfBins();
511 if(variable->HasUserDefinedBinning()){
512 varManager->AddVariable(namee, nbins, variable->GetBinning());
514 varManager->AddVariable(namee, nbins, variable->GetMinimum(), variable->GetMaximum());
516 //printf("For AliTaggedTrackAnalysis, had the variable %s and the one used %s\n",(const char*)variable->GetName(),(const char*) namee);
518 if(fPIDqa->HasHighResolutionHistos())
519 fTaggedTrackAnalysis->GetPIDqa()->SetHighResolutionHistos();
520 fTaggedTrackAnalysis->SetPID(fPID);
521 fTaggedTrackAnalysis->SetVariablesTRD(fVariablesTRDTaggedTrack);
522 fTaggedTrackAnalysis->InitContainer();
523 fOutput->Add(fTaggedTrackAnalysis->GetContainer());
524 fQA->Add(fTaggedTrackAnalysis->GetPIDQA());
525 fQA->Add(fTaggedTrackAnalysis->GetCutQA());
526 fQA->Add(fTaggedTrackAnalysis->GetQAcollection());
531 PostData(1, fOutput);
535 //____________________________________________________________
536 void AliAnalysisTaskHFE::UserExec(Option_t *){
540 AliDebug(3, "Starting Single Event Analysis");
542 AliError("Reconstructed Event not available");
546 AliDebug(4, Form("MC Event: %p", fMCEvent));
548 AliError("No MC Event, but MC Data required");
553 AliError("HFE cuts not available");
556 if(!fPID->IsInitialized()){
557 // Initialize PID with the given run number
558 fPID->InitializePID(fInputEvent->GetRunNumber());
561 // Initialize hadronic background from OADB Container
562 AliDebug(2, Form("Apply background factors: %s, OADB Container %p", fBackGroundFactorApply ? "Yes" : "No", fHadronBackgroundOADB));
563 if(fBackGroundFactorApply && !TestBit(kBackgroundInitialized)){
564 AliDebug(2, "Initializing Background from OADB");
565 if(!InitializeHadronBackground(fInputEvent->GetRunNumber())) AliError("Failed initializing hadronic background parameterization from OADB");
566 else AliDebug(2, "Successfully loaded Background from OADB");
567 SetBit(kBackgroundInitialized);
570 if(IsESDanalysis() && HasMCData()){
571 // Protect against missing MC trees
572 AliMCEventHandler *mcH = dynamic_cast<AliMCEventHandler *>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
574 AliError("No MC Event Handler available");
577 if(!mcH->InitOk()) return;
578 if(!mcH->TreeK()) return;
579 if(!mcH->TreeTR()) return;
582 // need the centrality for everything (MC also)
584 if(!ReadCentrality()) fCentralityF = -1;
585 //printf("pass centrality\n");
586 //printf("Reading fCentralityF %f\n",fCentralityF);
588 // See if pile up and z in the range
589 RejectionPileUpVertexRangeEventCut();
591 // Protect agains missing
593 //printf("Has MC data\n");
594 fSignalCuts->SetMCEvent(fMCEvent);
595 ProcessMC(); // Run the MC loop + MC QA in case MC Data are available
598 AliPIDResponse *pidResponse = fInputHandler->GetPIDResponse();
600 AliDebug(1, "Using default PID Response");
601 pidResponse = AliHFEtools::GetDefaultPID(HasMCData(), fInputEvent->IsA() == AliAODEvent::Class());
603 fPID->SetPIDResponse(pidResponse);
604 if(fPIDpreselect) fPIDpreselect->SetPIDResponse(pidResponse);
611 const char *specialTrigger = GetSpecialTrigger(fInputEvent->GetRunNumber());
612 // Check Trigger selection
614 AliDebug(2, Form("Special Trigger requested: %s", specialTrigger));
615 AliESDEvent *ev = dynamic_cast<AliESDEvent *>(fInputEvent);
616 if(!(ev && ev->IsTriggerClassFired(specialTrigger))){
617 AliDebug(2, "Event not selected");
619 } else AliDebug(2, "Event Selected");
620 } else AliDebug(2, "No Special Trigger requested");
625 PostData(1, fOutput);
629 //____________________________________________________________
630 void AliAnalysisTaskHFE::Terminate(Option_t *){
632 // Terminate not implemented at the moment
634 if(GetPlugin(kPostProcess)){
635 fOutput = dynamic_cast<TList *>(GetOutputData(1));
636 fQA = dynamic_cast<TList *>(GetOutputData(2));
638 AliError("Results not available");
642 AliError("QA output not available");
645 fContainer = dynamic_cast<AliHFEcontainer *>(fOutput->FindObject("trackContainer"));
647 AliError("Track container not found");
650 AliHFEpostAnalysis postanalysis;
651 postanalysis.SetTaskResults(fContainer);
652 TList *qalist = dynamic_cast<TList *>(fQA->FindObject("list_TaskQA"));
654 AliError("QA List not found");
657 postanalysis.SetTaskQA(qalist);
658 printf("Running post analysis\n");
660 postanalysis.DrawMCSignal2Background();
661 postanalysis.DrawEfficiency();
662 postanalysis.DrawPIDperformance();
663 postanalysis.DrawCutEfficiency();
665 if (GetPlugin(kIsElecBackGround)) {
666 AliHFEelecbackground elecBackGround;
668 if(!(oe = (TList*)dynamic_cast<TList *>(fOutput->FindObject("HFEelecbackground")))){
671 elecBackGround.Load(oe);
672 elecBackGround.Plot();
673 elecBackGround.PostProcess();
677 //_______________________________________________________________
678 Bool_t AliAnalysisTaskHFE::IsEventInBinZero() {
683 //printf("test in IsEventInBinZero\n");
685 AliError("Reconstructed Event not available");
690 const AliVVertex *vertex = fInputEvent->GetPrimaryVertex();
691 if(!vertex) return kTRUE;
692 //if(vertex) return kTRUE;
695 if(fInputEvent->GetNumberOfTracks()<=0) return kTRUE;
696 //if(fInputEvent->GetNumberOfTracks()>0) return kTRUE;
702 //____________________________________________________________
703 void AliAnalysisTaskHFE::ProcessMC(){
705 // Runs the MC Loop (filling the container for the MC Cut Steps with the observables pt, eta and phi)
706 // In case MC QA is on also MC QA loop is done
708 AliDebug(3, "Processing MC Information");
709 Double_t eventContainer [4];
710 eventContainer[0] = fMCEvent->GetPrimaryVertex()->GetZ();
711 eventContainer[2] = fCentralityF;
712 eventContainer[3] = fContributors;
713 fVz = eventContainer[0];
714 //printf("z position is %f\n",eventContainer[0]);
715 //if(fCFM->CheckEventCuts(AliHFEcuts::kEventStepGenerated, fMCEvent))
716 fCFM->GetEventContainer()->Fill(eventContainer,AliHFEcuts::kEventStepGenerated);
717 Int_t nElectrons = 0;
719 if(!((fIdentifiedAsPileUp) || (TMath::Abs(fVz) > fCuts->GetVertexRange()) || (fCentralityF < 0))){ //kStepMCGeneratedZOutNoPileUpCentralityFine
720 if (HasMCData() && IsQAOn(kMCqa)) {
721 AliDebug(2, "Running MC QA");
723 if(fMCEvent->Stack()){
724 fMCQA->SetMCEvent(fMCEvent);
725 fMCQA->SetGenEventHeader(fMCEvent->GenEventHeader());
726 fMCQA->SetCentrality(fCentralityF);
727 fMCQA->SetPercentrality(fCentralityPercent);
728 if(IsPbPb()) { fMCQA->SetPbPb();}
731 if(fisppMultiBin) fMCQA->SetPPMultiBin();
736 fMCQA->GetMesonKine();
738 // loop over all tracks for decayed electrons
739 for (Int_t igen = 0; igen < fMCEvent->GetNumberOfTracks(); igen++){
740 TParticle* mcpart = fMCEvent->Stack()->Particle(igen);
741 if(!mcpart) continue;
742 fMCQA->GetQuarkKine(mcpart, igen, AliHFEmcQA::kCharm);
743 fMCQA->GetQuarkKine(mcpart, igen, AliHFEmcQA::kBeauty);
744 fMCQA->GetHadronKine(mcpart, AliHFEmcQA::kCharm);
745 fMCQA->GetHadronKine(mcpart, AliHFEmcQA::kBeauty);
746 fMCQA->GetDecayedKine(mcpart, AliHFEmcQA::kCharm, AliHFEmcQA::kElectronPDG); // no accept cut
747 fMCQA->GetDecayedKine(mcpart, AliHFEmcQA::kBeauty, AliHFEmcQA::kElectronPDG); // no accept cut
748 fMCQA->GetDecayedKine(mcpart, AliHFEmcQA::kOthers, AliHFEmcQA::kElectronPDG); // no accept cut
750 //fMCQA->EndOfEventAna(AliHFEmcQA::kCharm);
751 //fMCQA->EndOfEventAna(AliHFEmcQA::kBeauty);
754 } // end of MC QA loop
756 // -----------------------------------------------------------------
757 fCFM->SetMCEventInfo(fMCEvent);
758 // fCFM->CheckEventCuts(AliCFManager::kEvtRecCuts, fESD);
760 fCFM->SetMCEventInfo(fInputEvent);
763 AliVParticle *mctrack = NULL;
764 AliDebug(3, Form("Number of Tracks: %d", fMCEvent->GetNumberOfTracks()));
765 for(Int_t imc = 0; imc <fMCEvent->GetNumberOfTracks(); imc++){
766 if(!(mctrack = fMCEvent->GetTrack(imc))) continue;
767 AliDebug(4, "Next MC Track");
768 if(ProcessMCtrack(mctrack)) nElectrons++;
771 // fCFM->CheckEventCuts(AliCFManager::kEvtRecCuts, fESD);
772 fQACollection->Fill("nElectron", nElectrons);
775 //____________________________________________________________
776 void AliAnalysisTaskHFE::ProcessESD(){
778 // Run Analysis of reconstructed event in ESD Mode
779 // Loop over Tracks, filter according cut steps defined in AliHFEcuts
781 AliDebug(3, "Processing ESD Event");
782 AliESDEvent *fESD = dynamic_cast<AliESDEvent *>(fInputEvent);
784 AliError("ESD Event required for ESD Analysis");
788 // Set magnetic field if V0 task on
789 if(fTaggedTrackAnalysis) {
790 fTaggedTrackAnalysis->SetMagneticField(fESD->GetMagneticField());
791 fTaggedTrackAnalysis->SetCentrality(fCentralityF);
792 if(IsPbPb()) fTaggedTrackAnalysis->SetPbPb();
793 else fTaggedTrackAnalysis->SetPP();
796 // Do event Normalization
797 Double_t eventContainer[4];
798 eventContainer[0] = 0.;
799 if(HasMCData()) eventContainer[0] = fVz;
801 if(fESD->GetPrimaryVertexSPD()) eventContainer[0] = fESD->GetPrimaryVertexSPD()->GetZ();
803 eventContainer[1] = 0.;
804 eventContainer[2] = fCentralityF;
805 eventContainer[3] = fContributors;
806 if(fTriggerAnalysis->IsOfflineTriggerFired(fESD, AliTriggerAnalysis::kV0AND))
807 eventContainer[1] = 1.;
810 fCFM->GetEventContainer()->Fill(eventContainer, AliHFEcuts::kEventStepRecNoCut);
813 if(fIdentifiedAsPileUp) return;
814 fCFM->GetEventContainer()->Fill(eventContainer, AliHFEcuts::kEventStepRecNoPileUp);
817 if(TMath::Abs(fCentralityF) < 0) return;
818 fCFM->GetEventContainer()->Fill(eventContainer, AliHFEcuts::kEventStepRecCentralityOk);
819 //printf("In ProcessESD %f\n",fCentralityF);
822 if(fIdentifiedAsOutInz) return;
823 fCFM->GetEventContainer()->Fill(eventContainer, AliHFEcuts::kEventStepZRange);
826 if(!fPassTheEventCut) return;
827 fCFM->GetEventContainer()->Fill(eventContainer, AliHFEcuts::kEventStepReconstructed);
831 fContainer->NewEvent();
833 if (GetPlugin(kIsElecBackGround)) {
834 fElecBackGround->SetEvent(fESD);
836 if (GetPlugin(kSecVtx)) {
837 fSecVtx->SetEvent(fESD);
838 fSecVtx->GetPrimaryCondition();
842 if (GetPlugin(kSecVtx)) {
843 fSecVtx->SetMCEvent(fMCEvent);
844 fSecVtx->SetMCQA(fMCQA);
846 if (GetPlugin(kIsElecBackGround)) {
847 fElecBackGround->SetMCEvent(fMCEvent);
851 Double_t container[10];
852 memset(container, 0, sizeof(Double_t) * 10);
853 // container for the output THnSparse
854 Double_t dataE[6]; // [pT, eta, Phi, type, 'C' or 'B']
855 Double_t dataDca[6]; // [source, pT, dca, centrality]
856 Int_t nElectronCandidates = 0;
857 AliESDtrack *track = NULL, *htrack = NULL;
858 AliMCParticle *mctrack = NULL;
859 AliMCParticle *mctrackmother = NULL;
862 Bool_t signal = kTRUE;
864 fCFM->SetRecEventInfo(fESD);
866 // minjung for IP QA(temporary ~ 2weeks)
868 fExtraCuts = new AliHFEextraCuts("hfetmpCuts","HFE tmp Cuts");
870 fExtraCuts->SetRecEventInfo(fESD);
872 // Electron background analysis
873 if (GetPlugin(kIsElecBackGround)) {
875 AliDebug(2, "Running BackGround Analysis");
877 fElecBackGround->Reset();
879 } // end of electron background analysis
883 AliDebug(3, Form("Number of Tracks: %d", fESD->GetNumberOfTracks()));
884 for(Int_t itrack = 0; itrack < fESD->GetNumberOfTracks(); itrack++){
885 AliDebug(4, "New ESD track");
886 track = fESD->GetTrack(itrack);
887 track->SetESDEvent(fESD);
889 // fill counts of v0-identified particles
891 if(track->TestBit(BIT(14))) v0pid = AliPID::kElectron;
892 else if(track->TestBit(BIT(15))) v0pid = AliPID::kPion;
893 else if(track->TestBit(BIT(16))) v0pid = AliPID::kProton;
894 // here the tagged track analysis will run
895 if(fTaggedTrackAnalysis && v0pid > -1){
896 AliDebug(1, Form("Track identified as %s", AliPID::ParticleName(v0pid)));
897 fTaggedTrackAnalysis->ProcessTrack(track, v0pid);
900 AliDebug(3, Form("Doing track %d, %p", itrack, track));
902 //////////////////////////////////////
904 /////////////////////////////////////
905 if(fPIDpreselect && fCutspreselect) {
906 if(!PreSelectTrack(track)) continue;
911 // Fill step without any cut
914 // Check if it is electrons near the vertex
915 if(!(mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(track->GetLabel()))))) continue;
917 if(fFillSignalOnly && !fCFM->CheckParticleCuts(AliHFEcuts::kStepMCGenerated, mctrack)) signal = kFALSE;
918 else AliDebug(3, "Signal Electron");
920 // Fill K pt for Ke3 contributions
921 if(mctrack && (TMath::Abs(mctrack->Particle()->GetPdgCode())==321)) fQACollection->Fill("Kptspectra",mctrack->Pt());
922 else if(mctrack && (TMath::Abs(mctrack->Particle()->GetPdgCode())==130)) fQACollection->Fill("K0Lptspectra",mctrack->Pt());
924 // Cache new Track information inside the var manager
925 fVarManager->NewTrack(track, mctrack, fCentralityF, -1, signal);
928 if(signal || !fFillSignalOnly){
929 fVarManager->FillContainer(fContainer, "recTrackContReco", AliHFEcuts::kStepRecNoCut, kFALSE);
930 fVarManager->FillContainer(fContainer, "recTrackContMC", AliHFEcuts::kStepRecNoCut, kTRUE);
934 // RecKine: ITSTPC cuts
935 if(!ProcessCutStep(AliHFEcuts::kStepRecKineITSTPC, track)) continue;
939 if(fRejectKinkMother) {
940 if(track->GetKinkIndex(0) != 0) continue; } // Quick and dirty fix to reject both kink mothers and daughters
941 if(!ProcessCutStep(AliHFEcuts::kStepRecPrim, track)) continue;
943 // HFEcuts: ITS layers cuts
944 if(!ProcessCutStep(AliHFEcuts::kStepHFEcutsITS, track)) continue;
946 // HFE cuts: TOF PID and mismatch flag
947 if(!ProcessCutStep(AliHFEcuts::kStepHFEcutsTOF, track)) continue;
949 // HFE cuts: TPC PID cleanup
950 if(!ProcessCutStep(AliHFEcuts::kStepHFEcutsTPC, track)) continue;
952 // HFEcuts: Nb of tracklets TRD0
953 if(!ProcessCutStep(AliHFEcuts::kStepHFEcutsTRD, track)) continue;
955 // Fill correlation maps before PID
956 if(signal && fContainer->GetCorrelationMatrix("correlationstepbeforePID")) {
957 //printf("Fill correlation maps before PID\n");
958 fVarManager->FillCorrelationMatrix(fContainer->GetCorrelationMatrix("correlationstepbeforePID"));
962 FillProductionVertex(track);
965 fMCQA->SetCentrality(fCentralityF);
966 if(mctrack && (TMath::Abs(mctrack->Particle()->GetPdgCode()) == 11)){
967 Double_t weightElecBgV0[kBgLevels] = {0.,0.,0.};
968 Double_t hfeimpactRtmp=0., hfeimpactnsigmaRtmp=0.;
969 fExtraCuts->GetHFEImpactParameters(track, hfeimpactRtmp, hfeimpactnsigmaRtmp);
970 UChar_t itsPixel = track->GetITSClusterMap();
971 Double_t ilyrhit=0, ilyrstat=0;
972 for(Int_t ilyr=0; ilyr<6; ilyr++){
973 if(TESTBIT(itsPixel, ilyr)) ilyrhit += TMath::Power(2,ilyr);
974 if(fExtraCuts->CheckITSstatus(fExtraCuts->GetITSstatus(track,ilyr))) ilyrstat += TMath::Power(2,ilyr);
976 fMCQA->SetITSInfo(ilyrhit,ilyrstat);
977 fMCQA->SetHFEImpactParameters(hfeimpactRtmp, hfeimpactnsigmaRtmp);
978 fMCQA->SetTrkKine(track->Pt(),track->Eta(), track->Phi());
979 fMCQA->SetContainerStep(3);
980 for(Int_t iLevel = 0; iLevel < kBgLevels; iLevel++){
981 weightElecBgV0[iLevel] = fMCQA->GetWeightFactor(mctrack, iLevel); // positive:conversion e, negative: nonHFE
982 if(!fisNonHFEsystematics)break;
985 if(fisNonHFEsystematics){
986 //Fill additional containers for electron source distinction
987 Int_t elecSource = 0;
988 elecSource = fMCQA->GetElecSource(mctrack->Particle());
989 const Char_t *sourceName[kElecBgSpecies]={"Pion","Eta","Omega","Phi","EtaPrime","Rho"};
990 const Char_t *levelName[kBgLevels]={"Best","Lower","Upper"};
992 for(Int_t iSource = AliHFEmcQA::kPi0; iSource <= AliHFEmcQA::kGammaRho0; iSource++){
993 if((iSource == AliHFEmcQA::kElse)||(iSource == AliHFEmcQA::kMisID)) continue;
994 if(elecSource == iSource){
995 for(Int_t iLevel = 0; iLevel < kBgLevels; iLevel++){
996 if(weightElecBgV0[iLevel]>0){
997 fVarManager->FillContainer(fContainer, Form("conversionElecs%s%s",sourceName[iName], levelName[iLevel]), 3, kFALSE, weightElecBgV0[iLevel]);
999 else if(weightElecBgV0[iLevel]<0){
1000 fVarManager->FillContainer(fContainer, Form("mesonElecs%s%s",sourceName[iName], levelName[iLevel]), 3, kFALSE, -1*weightElecBgV0[iLevel]);
1006 if(iName == kElecBgSpecies)iName = 0;
1010 if(weightElecBgV0[0]>0) {
1011 fVarManager->FillContainer(fContainer, "conversionElecs", 3, kFALSE, weightElecBgV0[0]);
1012 fVarManager->FillContainer(fContainer, "conversionElecs", 4, kTRUE, weightElecBgV0[0]);
1014 else if(weightElecBgV0[0]<0) {
1015 fVarManager->FillContainer(fContainer, "mesonElecs", 3, kFALSE, -1*weightElecBgV0[0]);
1016 fVarManager->FillContainer(fContainer, "mesonElecs", 4, kTRUE, -1*weightElecBgV0[0]);
1022 Double_t hfeimpactR4all=0., hfeimpactnsigmaR4all=0.;
1023 Int_t sourceDca =-1;
1024 if(mctrack && (TMath::Abs(mctrack->Particle()->GetPdgCode()) == 211)){
1026 fExtraCuts->GetHFEImpactParameters(track, hfeimpactR4all, hfeimpactnsigmaR4all);
1027 dataDca[0]=0; //pion
1028 dataDca[1]=track->Pt();
1029 dataDca[2]=hfeimpactR4all;
1030 dataDca[3]=fCentralityF;
1032 dataDca[5] = double(track->Charge());
1033 fQACollection->Fill("Dca", dataDca);
1036 else if(mctrack && (TMath::Abs(mctrack->Particle()->GetPdgCode()) == 11)){ // to increas statistics for Martin
1038 fExtraCuts->GetHFEImpactParameters(track, hfeimpactR4all, hfeimpactnsigmaR4all);
1039 if(fSignalCuts->IsCharmElectron(track)){
1042 else if(fSignalCuts->IsBeautyElectron(track)){
1045 else if(fSignalCuts->IsGammaElectron(track)){
1048 else if(fSignalCuts->IsNonHFElectron(track)){
1051 else if(fSignalCuts->IsJpsiElectron(track)){
1057 dataDca[0]=sourceDca;
1058 dataDca[1]=track->Pt();
1059 dataDca[2]=hfeimpactR4all;
1060 dataDca[3]=fCentralityF;
1062 dataDca[5] = double(track->Charge());
1063 if(signal) fQACollection->Fill("Dca", dataDca);
1068 if(TMath::Abs(track->Eta()) < 0.5){
1069 if(track->GetInnerParam())
1070 fQACollection->Fill("TPCdEdxBeforePID", track->GetInnerParam()->P(), track->GetTPCsignal());
1071 fQACollection->Fill("TPCnSigmaBeforePID", track->P(), fInputHandler->GetPIDResponse()->NumberOfSigmasTPC(track, AliPID::kElectron));
1074 AliHFEpidObject hfetrack;
1075 hfetrack.SetAnalysisType(AliHFEpidObject::kESDanalysis);
1076 hfetrack.SetRecTrack(track);
1077 if(HasMCData()) hfetrack.SetMCTrack(mctrack);
1078 hfetrack.SetCentrality(fCentralityF);
1079 if(IsPbPb()) hfetrack.SetPbPb();
1080 else hfetrack.SetPP();
1081 fPID->SetVarManager(fVarManager);
1082 if(!fPID->IsSelected(&hfetrack, fContainer, "recTrackCont", fPIDqa)) continue;
1083 nElectronCandidates++;
1085 // Temporary histogram for chi2/ITS cluster
1087 TBits shared = track->GetTPCSharedMap();
1089 if(shared.CountBits() >= 2) sharebit=1;
1091 Double_t itschi2percluster = 0.0;
1092 Double_t itsnbcls = static_cast<Double_t>(track->GetNcls(0));
1093 if(itsnbcls > 0) itschi2percluster = track->GetITSchi2()/itsnbcls;
1095 Double_t itsChi2[7] = {track->Pt(),track->Eta(), track->Phi(),
1096 fCentralityF,track->GetTPCsignalN(), sharebit,itschi2percluster};
1097 fQACollection->Fill("fChi2perITScluster", itsChi2);
1101 Double_t itschi2percluster = 0.0;
1102 Double_t itsnbcls = static_cast<Double_t>(track->GetNcls(0));
1103 if(itsnbcls > 0) itschi2percluster = track->GetITSchi2()/itsnbcls;
1105 Double_t itsChi2[3] = {track->Pt(), fCentralityF, itschi2percluster};
1106 fQACollection->Fill("fChi2perITScluster", itsChi2);
1109 // Fill Histogram for Hadronic Background
1111 if(mctrack && (TMath::Abs(mctrack->Particle()->GetPdgCode()) != 11))
1112 fVarManager->FillContainer(fContainer, "hadronicBackground", UInt_t(0), kFALSE);
1114 // Fill Ke3 contributions
1115 Int_t glabel=TMath::Abs(mctrack->GetMother());
1116 if((mctrackmother = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(glabel)))){
1117 if(TMath::Abs(mctrackmother->Particle()->GetPdgCode())==321)
1118 fQACollection->Fill("Ke3Kecorr",mctrack->Pt(),mctrackmother->Pt());
1119 else if(TMath::Abs(mctrackmother->Particle()->GetPdgCode())==130)
1120 fQACollection->Fill("Ke3K0Lecorr",mctrack->Pt(),mctrackmother->Pt());
1127 // Apply weight for background contamination
1128 if(fBackGroundFactorApply) {
1129 if(IsPbPb() && fCentralityF >= 0) fWeightBackGround = fkBackGroundFactorArray[fCentralityF >= 0 ? fCentralityF : 0]->Eval(TMath::Abs(track->P()));
1130 else fWeightBackGround = fkBackGroundFactorArray[0]->Eval(TMath::Abs(track->P())); // pp case
1132 if(fWeightBackGround < 0.0) fWeightBackGround = 0.0;
1133 else if(fWeightBackGround > 1.0) fWeightBackGround = 1.0;
1134 // weightBackGround as special weight
1135 fVarManager->FillContainer(fContainer, "hadronicBackground", 1, kFALSE, fWeightBackGround);
1137 fVarManager->FillCorrelationMatrix(fContainer->GetCorrelationMatrix("correlationstepafterPID"));
1140 Bool_t bTagged=kFALSE;
1141 if(GetPlugin(kSecVtx)) {
1142 AliDebug(2, "Running Secondary Vertex Analysis");
1143 if(fSecVtx->Process(track) && signal) {
1144 fVarManager->FillContainer(fContainer, "recTrackContSecvtxReco", AliHFEcuts::kStepHFEcutsSecvtx, kFALSE);
1145 fVarManager->FillContainer(fContainer, "recTrackContSecvtxMC", AliHFEcuts::kStepHFEcutsSecvtx, kTRUE);
1151 dataE[0] = track->Pt();
1152 dataE[1] = track->Eta();
1153 dataE[2] = track->Phi();
1154 dataE[3] = track->Charge();
1158 // Track selected: distinguish between true and fake
1159 AliDebug(1, Form("Candidate Selected, filling THnSparse, PID: %d\n", mctrack->Particle()->GetPdgCode()));
1160 if((pid = TMath::Abs(mctrack->Particle()->GetPdgCode())) == 11){
1162 if(fSignalCuts->IsCharmElectron(track))
1164 else if(fSignalCuts->IsBeautyElectron(track))
1166 AliDebug(1, Form("Type: %d\n", type));
1168 dataE[5] = type; // beauty[1] or charm[2]
1169 dataE[4] = 2; // signal electron
1172 dataE[4] = 1; // not a signal electron
1177 // Fill THnSparse with the information for Fake Electrons
1181 // fill the performance THnSparse, if the mc origin could be defined
1183 AliDebug(1, Form("Entries: [%.3f|%.3f|%.3f|%f|%f|%f]\n", dataE[0],dataE[1],dataE[2],dataE[3],dataE[4],dataE[5]));
1184 fQACollection->Fill("PIDperformance", dataE);
1188 // Electron background analysis
1189 if (GetPlugin(kIsElecBackGround)) {
1191 AliDebug(2, "Running BackGround Analysis");
1193 for(Int_t jtrack = 0; jtrack < fESD->GetNumberOfTracks(); jtrack++){
1194 htrack = fESD->GetTrack(jtrack);
1195 if ( itrack == jtrack ) continue;
1196 fElecBackGround->PairAnalysis(track, htrack);
1198 } // end of electron background analysis
1200 if (GetPlugin(kDEstep)) {
1201 Double_t weightElecBgV0[kBgLevels] = {0.,0.,0.,};
1202 Int_t elecSource = 0;
1203 Double_t hfeimpactR=0., hfeimpactnsigmaR=0.;
1204 fExtraCuts->GetHFEImpactParameters(track, hfeimpactR, hfeimpactnsigmaR);
1207 if(mctrack && (TMath::Abs(mctrack->Particle()->GetPdgCode()) != 11)){
1208 fQACollection->Fill("hadronsBeforeIPcut",track->Pt());
1210 if(fMCQA && signal) {
1212 fMCQA->SetContainerStep(0);
1213 for(Int_t iLevel = 0; iLevel < kBgLevels; iLevel++){
1214 weightElecBgV0[iLevel] = fMCQA->GetWeightFactor(mctrack, iLevel); // positive:conversion e, negative: nonHFE
1215 if(!fisNonHFEsystematics)break;
1218 if(fisNonHFEsystematics){
1219 //Fill additional containers for electron source distinction
1220 elecSource = fMCQA->GetElecSource(mctrack->Particle());
1221 const Char_t *sourceName[kElecBgSpecies]={"Pion","Eta","Omega","Phi","EtaPrime","Rho"};
1222 const Char_t *levelName[kBgLevels]={"Best","Lower","Upper"};
1224 for(Int_t iSource = AliHFEmcQA::kPi0; iSource <= AliHFEmcQA::kGammaRho0; iSource++){
1225 if((iSource == AliHFEmcQA::kElse)||(iSource == AliHFEmcQA::kMisID)) continue;
1226 if(elecSource == iSource){
1227 for(Int_t iLevel = 0; iLevel < kBgLevels; iLevel++){
1228 if(weightElecBgV0[iLevel]>0) fVarManager->FillContainer(fContainer, Form("conversionElecs%s%s",sourceName[iName], levelName[iLevel]), 0, kFALSE, weightElecBgV0[iLevel]);
1229 else if(weightElecBgV0[iLevel]<0) fVarManager->FillContainer(fContainer, Form("mesonElecs%s%s",sourceName[iName], levelName[iLevel]), 0, kFALSE, -1*weightElecBgV0[iLevel]);
1234 if(iName == kElecBgSpecies)iName = 0;
1238 if(weightElecBgV0[0]>0) {
1239 fVarManager->FillContainer(fContainer, "conversionElecs", 0, kFALSE, weightElecBgV0[0]);
1240 fVarManager->FillContainer(fContainer, "conversionElecs", 5, kTRUE, weightElecBgV0[0]);
1242 else if(weightElecBgV0[0]<0) {
1243 fVarManager->FillContainer(fContainer, "mesonElecs", 0, kFALSE, -1*weightElecBgV0[0]);
1244 fVarManager->FillContainer(fContainer, "mesonElecs", 5, kTRUE, -1*weightElecBgV0[0]);
1247 if(bTagged){ // bg estimation for the secondary vertex tagged signals
1248 if(weightElecBgV0[0]>0) fVarManager->FillContainer(fContainer, "conversionElecs", 2, kFALSE, weightElecBgV0[0]);
1249 else if(weightElecBgV0[0]<0) fVarManager->FillContainer(fContainer, "mesonElecs", 2, kFALSE, -1*weightElecBgV0[0]);
1254 dataDca[0]=-1; //for data, don't know the origin
1255 dataDca[1]=track->Pt();
1256 dataDca[2]=hfeimpactR;
1257 dataDca[3]=fCentralityF;
1259 dataDca[5] = double(track->Charge());
1260 if (!HasMCData()) fQACollection->Fill("Dca", dataDca);
1262 // Fill Containers for impact parameter analysis
1263 if(!fCFM->CheckParticleCuts(AliHFEcuts::kStepHFEcutsDca + AliHFEcuts::kNcutStepsMCTrack + AliHFEcuts::kNcutStepsRecTrack,track)) continue;
1265 if(fMCQA && signal) {
1266 fMCQA->SetContainerStep(1);
1267 for(Int_t iLevel = 0; iLevel < kBgLevels; iLevel++){
1268 weightElecBgV0[iLevel] = fMCQA->GetWeightFactor(mctrack, iLevel); // positive:conversion e, negative: nonHFE
1269 if(!fisNonHFEsystematics)break;
1271 if(fisNonHFEsystematics){
1272 //Fill additional containers for electron source distinction
1273 elecSource = fMCQA->GetElecSource(mctrack->Particle());
1274 const Char_t *sourceName[kElecBgSpecies]={"Pion","Eta","Omega","Phi","EtaPrime","Rho"};
1275 const Char_t *levelName[kBgLevels]={"Best","Lower","Upper"};
1277 for(Int_t iSource = AliHFEmcQA::kPi0; iSource <= AliHFEmcQA::kGammaRho0; iSource++){
1278 if((iSource == AliHFEmcQA::kElse)||(iSource == AliHFEmcQA::kMisID)) continue;
1279 if(elecSource == iSource){
1280 for(Int_t iLevel = 0; iLevel < kBgLevels; iLevel++){
1281 if(weightElecBgV0[iLevel]>0) fVarManager->FillContainer(fContainer, Form("conversionElecs%s%s",sourceName[iName], levelName[iLevel]), 1, kFALSE, weightElecBgV0[iLevel]);
1282 else if(weightElecBgV0[iLevel]<0) fVarManager->FillContainer(fContainer, Form("mesonElecs%s%s",sourceName[iName], levelName[iLevel]), 1, kFALSE, -1*weightElecBgV0[iLevel]);
1287 if(iName == kElecBgSpecies)iName = 0;
1291 if(weightElecBgV0[0]>0) {
1292 fVarManager->FillContainer(fContainer, "conversionElecs", 1, kFALSE, weightElecBgV0[0]);
1293 fVarManager->FillContainer(fContainer, "conversionElecs", 6, kTRUE, weightElecBgV0[0]);
1295 else if(weightElecBgV0[0]<0) {
1296 fVarManager->FillContainer(fContainer, "mesonElecs", 1, kFALSE, -1*weightElecBgV0[0]);
1297 fVarManager->FillContainer(fContainer, "mesonElecs", 6, kTRUE, -1*weightElecBgV0[0]);
1303 fVarManager->FillContainer(fContainer, "recTrackContDEReco", AliHFEcuts::kStepHFEcutsDca, kFALSE);
1304 fVarManager->FillContainer(fContainer, "recTrackContDEMC", AliHFEcuts::kStepHFEcutsDca, kTRUE);
1305 fVarManager->FillCorrelationMatrix(fContainer->GetCorrelationMatrix("correlationstepafterDE"));
1308 if(mctrack && (TMath::Abs(mctrack->Particle()->GetPdgCode()) != 11)){
1309 fQACollection->Fill("hadronsAfterIPcut",track->Pt());
1315 fQACollection->Fill("nElectronTracksEvent", nElectronCandidates);
1318 //____________________________________________________________
1319 void AliAnalysisTaskHFE::ProcessAOD(){
1321 // Run Analysis in AOD Mode
1322 // Function is still in development
1324 //printf("Process AOD\n");
1325 AliDebug(3, "Processing AOD Event");
1326 Double_t eventContainer[4];
1327 if(HasMCData()) eventContainer[0] = fVz;
1329 eventContainer[0] = fInputEvent->GetPrimaryVertex()->GetZ();
1331 eventContainer[1] = 1.; // No Information available in AOD analysis, assume all events have V0AND
1332 eventContainer[2] = fCentralityF;
1333 eventContainer[3] = fContributors;
1335 //printf("value event container %f, %f, %f, %f\n",eventContainer[0],eventContainer[1],eventContainer[2],eventContainer[3]);
1337 AliAODEvent *fAOD = dynamic_cast<AliAODEvent *>(fInputEvent);
1339 AliError("AOD Event required for AOD Analysis");
1343 //printf("Will fill\n");
1345 fCFM->GetEventContainer()->Fill(eventContainer, AliHFEcuts::kEventStepRecNoCut);
1348 if(fIdentifiedAsPileUp) return;
1349 fCFM->GetEventContainer()->Fill(eventContainer, AliHFEcuts::kEventStepRecNoPileUp);
1352 if(fIdentifiedAsOutInz) return;
1353 fCFM->GetEventContainer()->Fill(eventContainer, AliHFEcuts::kEventStepZRange);
1356 if(!fPassTheEventCut) return;
1357 fCFM->GetEventContainer()->Fill(eventContainer, AliHFEcuts::kEventStepReconstructed);
1360 fContainer->NewEvent();
1362 fCFM->SetRecEventInfo(fAOD);
1364 // Look for kink mother
1365 Int_t numberofvertices = fAOD->GetNumberOfVertices();
1366 Double_t listofmotherkink[numberofvertices];
1367 Int_t numberofmotherkink = 0;
1368 for(Int_t ivertex=0; ivertex < numberofvertices; ivertex++) {
1369 AliAODVertex *aodvertex = fAOD->GetVertex(ivertex);
1370 if(!aodvertex) continue;
1371 if(aodvertex->GetType()==AliAODVertex::kKink) {
1372 AliAODTrack *mother = (AliAODTrack *) aodvertex->GetParent();
1373 if(!mother) continue;
1374 Int_t idmother = mother->GetID();
1375 listofmotherkink[numberofmotherkink] = idmother;
1376 //printf("ID %d\n",idmother);
1377 numberofmotherkink++;
1380 //printf("Number of kink mother in the events %d\n",numberofmotherkink);
1384 AliAODTrack *track = NULL;
1385 AliAODMCParticle *mctrack = NULL;
1386 Double_t dataE[6]; // [pT, eta, Phi, Charge, type, 'C' or 'B']
1387 Int_t nElectronCandidates = 0;
1390 //printf("Number of track %d\n",(Int_t) fAOD->GetNumberOfTracks());
1391 for(Int_t itrack = 0; itrack < fAOD->GetNumberOfTracks(); itrack++){
1392 track = fAOD->GetTrack(itrack); mctrack = NULL;
1393 if(!track) continue;
1395 if(track->GetFlags() != fFlags) continue; // Only process AOD tracks where the HFE is set
1397 //printf("Pass the flag\n");
1402 Int_t label = TMath::Abs(track->GetLabel());
1404 mctrack = dynamic_cast<AliAODMCParticle *>(fMCEvent->GetTrack(label));
1405 if(fFillSignalOnly && !fCFM->CheckParticleCuts(AliHFEcuts::kStepMCGenerated, mctrack)) signal = kFALSE;
1407 fVarManager->NewTrack(track, mctrack, fCentralityF, -1, kTRUE);
1410 if(signal || !fFillSignalOnly){
1411 fVarManager->FillContainer(fContainer, "recTrackContReco", AliHFEcuts::kStepRecNoCut, kFALSE);
1412 fVarManager->FillContainer(fContainer, "recTrackContMC", AliHFEcuts::kStepRecNoCut, kTRUE);
1417 // RecKine: ITSTPC cuts
1418 if(!ProcessCutStep(AliHFEcuts::kStepRecKineITSTPC, track)) continue;
1420 // Reject kink mother
1421 Bool_t kinkmotherpass = kTRUE;
1422 for(Int_t kinkmother = 0; kinkmother < numberofmotherkink; kinkmother++) {
1423 if(track->GetID() == listofmotherkink[kinkmother]) {
1424 kinkmotherpass = kFALSE;
1428 if(!kinkmotherpass) continue;
1431 if(!ProcessCutStep(AliHFEcuts::kStepRecPrim, track)) continue;
1433 // HFEcuts: ITS layers cuts
1434 if(!ProcessCutStep(AliHFEcuts::kStepHFEcutsITS, track)) continue;
1436 // HFE cuts: TOF PID and mismatch flag
1437 if(!ProcessCutStep(AliHFEcuts::kStepHFEcutsTOF, track)) continue;
1439 // HFE cuts: TPC PID cleanup
1440 if(!ProcessCutStep(AliHFEcuts::kStepHFEcutsTPC, track)) continue;
1442 // HFEcuts: Nb of tracklets TRD0
1443 if(!ProcessCutStep(AliHFEcuts::kStepHFEcutsTRD, track)) continue;
1446 // Fill correlation maps before PID
1447 if(signal && fContainer->GetCorrelationMatrix("correlationstepbeforePID")) {
1448 //printf("Fill correlation maps before PID\n");
1449 fVarManager->FillCorrelationMatrix(fContainer->GetCorrelationMatrix("correlationstepbeforePID"));
1452 //printf("Will process to PID\n");
1454 // track accepted, do PID
1455 AliHFEpidObject hfetrack;
1456 hfetrack.SetAnalysisType(AliHFEpidObject::kAODanalysis);
1457 hfetrack.SetRecTrack(track);
1458 if(HasMCData()) hfetrack.SetMCTrack(mctrack);
1459 hfetrack.SetCentrality(fCentralityF);
1460 if(IsPbPb()) hfetrack.SetPbPb();
1461 else hfetrack.SetPP();
1462 fPID->SetVarManager(fVarManager);
1463 if(!fPID->IsSelected(&hfetrack, fContainer, "recTrackCont", fPIDqa)) continue; // we will do PID here as soon as possible
1466 // Apply weight for background contamination
1467 //Double_t weightBackGround = 1.0;
1469 // Apply weight for background contamination
1470 if(fBackGroundFactorApply) {
1471 if(IsPbPb() && fCentralityF >= 0) fWeightBackGround = fkBackGroundFactorArray[fCentralityF >= 0 ? fCentralityF : 0]->Eval(TMath::Abs(track->P()));
1472 else fWeightBackGround = fkBackGroundFactorArray[0]->Eval(TMath::Abs(track->P())); // pp case
1474 if(fWeightBackGround < 0.0) fWeightBackGround = 0.0;
1475 else if(fWeightBackGround > 1.0) fWeightBackGround = 1.0;
1476 // weightBackGround as special weight
1477 fVarManager->FillContainer(fContainer, "hadronicBackground", 1, kFALSE, fWeightBackGround);
1479 fVarManager->FillCorrelationMatrix(fContainer->GetCorrelationMatrix("correlationstepafterPID"));
1482 nElectronCandidates++;
1484 dataE[0] = track->Pt();
1485 dataE[1] = track->Eta();
1486 dataE[2] = track->Phi();
1487 dataE[3] = track->Charge();
1490 // Track selected: distinguish between true and fake
1491 AliDebug(1, Form("Candidate Selected, filling THnSparse, PID: %d\n", mctrack ? mctrack->GetPdgCode(): -1));
1492 if(mctrack && ((pid = TMath::Abs(mctrack->GetPdgCode())) == 11)){
1495 if(fSignalCuts->IsCharmElectron(track))
1497 else if(fSignalCuts->IsBeautyElectron(track))
1499 AliDebug(1, Form("Type: %d\n", type));
1501 dataE[5] = type; // beauty[1] or charm[2]
1502 dataE[4] = 2; // signal electron
1505 dataE[4] = 1; // not a signal electron
1510 // Fill THnSparse with the information for Fake Electrons
1514 // fill the performance THnSparse, if the mc origin could be defined
1516 AliDebug(1, Form("Entries: [%.3f|%.3f|%.3f|%f|%f|%f]\n", dataE[0],dataE[1],dataE[2],dataE[3],dataE[4],dataE[5]));
1517 fQACollection->Fill("PIDperformance", dataE);
1521 fQACollection->Fill("nElectronTracksEvent", nElectronCandidates);
1524 //____________________________________________________________
1525 Bool_t AliAnalysisTaskHFE::ProcessMCtrack(AliVParticle *track){
1527 // Filter the Monte Carlo Track
1528 // Additionally Fill a THnSparse for Signal To Background Studies
1529 // Works for AOD and MC analysis Type
1531 fVarManager->NewTrack(track, NULL, fCentralityF, -1, kTRUE);
1532 Double_t signalContainer[6];
1534 signalContainer[0] = track->Pt();
1535 signalContainer[1] = track->Eta();
1536 signalContainer[2] = track->Phi();
1537 signalContainer[3] = track->Charge()/3;
1539 Double_t vertex[3] = {0.,0.,0.}; // Production vertex cut to mask gammas which are NOT supposed to have hits in the first ITS layer(s)
1540 if(IsESDanalysis()){
1541 AliMCParticle *mctrack = dynamic_cast<AliMCParticle *>(track);
1543 vertex[0] = mctrack->Particle()->Vx();
1544 vertex[1] = mctrack->Particle()->Vy();
1547 AliAODMCParticle *aodmctrack = dynamic_cast<AliAODMCParticle *>(track);
1548 if(aodmctrack) aodmctrack->XvYvZv(vertex);
1551 if(!fCFM->CheckParticleCuts(AliHFEcuts::kStepMCGenerated, track)) return kFALSE;
1552 fVarManager->FillContainer(fContainer, "MCTrackCont", AliHFEcuts::kStepMCGenerated, kFALSE);
1553 signalContainer[4] = 0;
1554 if(fSignalCuts->IsSelected(track)){
1555 //fVarManager->FillContainer(fContainer, "MCTrackCont", AliHFEcuts::kStepMCsignal, kFALSE);
1556 // Filling of the Signal/Background histogram using the
1557 // definition of the codes for charm and beauty as below in
1558 // th crearion of the histogram
1559 if(fSignalCuts->IsCharmElectron(track))
1560 signalContainer[4] = 1;
1562 signalContainer[4] = 2;
1564 signalContainer[4] = 0; // (and other background)
1566 signalContainer[5] = 0;
1567 // apply cut on the sqrt of the production vertex
1568 Double_t radVertex = TMath::Sqrt(vertex[0]*vertex[0] + vertex[1] * vertex[1]);
1569 if(radVertex < 3.5){
1570 // Within first ITS layer(2) -> Background we cannot reject by ITS cut, let it pass
1571 signalContainer[5] = 1;
1572 } else if (radVertex < 7.5){
1573 signalContainer[5] = 2;
1575 fQACollection->Fill("SignalToBackgroundMC", signalContainer);
1577 // Step GeneratedZOutNoPileUp
1578 if((fIdentifiedAsPileUp) || (TMath::Abs(fVz) > fCuts->GetVertexRange()) || (fCentralityF < 0)) return kFALSE;
1579 fVarManager->FillContainer(fContainer, "MCTrackCont", AliHFEcuts::kStepMCGeneratedZOutNoPileUpCentralityFine, kFALSE);
1580 //printf("In ProcessMCtrack %f\n",fCentralityF);
1582 // Step Generated Event Cut
1583 if(!fPassTheEventCut) return kFALSE;
1584 fVarManager->FillContainer(fContainer, "MCTrackCont", AliHFEcuts::kStepMCGeneratedEventCut, kFALSE);
1586 if(!fCFM->CheckParticleCuts(AliHFEcuts::kStepMCInAcceptance, track)) return kFALSE;
1587 fVarManager->FillContainer(fContainer, "MCTrackCont", AliHFEcuts::kStepMCInAcceptance, kFALSE);
1591 //____________________________________________________________
1592 Bool_t AliAnalysisTaskHFE::PreSelectTrack(AliESDtrack *track) const {
1598 Bool_t survived = kTRUE;
1600 if(!fCutspreselect->CheckParticleCuts(AliHFEcuts::kStepRecKineITSTPC, track)) {
1602 //printf("Did not pass AliHFEcuts::kStepRecKineITSTPC\n");
1604 //else printf("Pass AliHFEcuts::kStepRecKineITSTPC\n");
1605 if(!fCutspreselect->CheckParticleCuts(AliHFEcuts::kStepRecPrim, track)) {
1607 //printf("Did not pass AliHFEcuts::kStepRecPrim\n");
1609 //else printf("Pass AliHFEcuts::kStepRecPrim\n");
1610 if(!fCutspreselect->CheckParticleCuts(AliHFEcuts::kStepHFEcutsITS, track)) {
1612 //printf("Did not pass AliHFEcuts::kStepHFEcutsITS\n");
1614 //else printf("Pass AliHFEcuts::kStepHFEcutsITS\n");
1615 if(!fCutspreselect->CheckParticleCuts(AliHFEcuts::kStepHFEcutsTOF, track)) {
1617 //printf("Did not pass AliHFEcuts::kStepHFEcutsTOF\n");
1619 //else printf("Pass AliHFEcuts::kStepHFEcutsTOF\n");
1620 if(!fCutspreselect->CheckParticleCuts(AliHFEcuts::kStepHFEcutsTRD, track)) {
1622 //printf("Did not pass AliHFEcuts::kStepHFEcutsTRD\n");
1624 //else printf("Pass AliHFEcuts::kStepHFEcutsTRD\n");
1628 AliHFEpidObject hfetrack;
1629 hfetrack.SetAnalysisType(AliHFEpidObject::kESDanalysis);
1630 hfetrack.SetRecTrack(track);
1631 if(!fPIDpreselect->IsSelected(&hfetrack)) {
1632 //printf("Did not pass AliHFEcuts::kPID\n");
1635 //else printf("Pass AliHFEcuts::kPID\n");
1641 //____________________________________________________________
1642 void AliAnalysisTaskHFE::MakeEventContainer(){
1644 // Create the event container for the correction framework and link it
1645 // 1st bin: Vertex z-position
1646 // 2nd bin: V0AND decision (normalization to sigma_inel)
1647 // 3rd bin: Centrality class (for pp defined as number of contributors in vertex.)
1648 // 4th bin: Number of contributors > 0
1651 const Int_t kNvar = 4; // number of variables on the grid:
1652 Int_t nBins[kNvar] = {120, 2, 11, 2};
1653 Double_t binMin[kNvar] = {-30. , 0., 0.0, 0.};
1654 Double_t binMax[kNvar] = {30., 2., 11.0, 2.};
1656 AliCFContainer *evCont = new AliCFContainer("eventContainer", "Container for events", AliHFEcuts::kNcutStepsEvent, kNvar, nBins);
1658 Double_t *vertexBins = AliHFEtools::MakeLinearBinning(nBins[0], binMin[0], binMax[0]);
1659 Double_t *v0andBins = AliHFEtools::MakeLinearBinning(nBins[1], binMin[1], binMax[1]);
1660 Double_t *centralityBins = AliHFEtools::MakeLinearBinning(nBins[2], binMin[2], binMax[2]);
1661 Double_t *contributorsBins = AliHFEtools::MakeLinearBinning(nBins[3], binMin[3], binMax[3]);
1662 evCont->SetBinLimits(0, vertexBins);
1663 evCont->SetBinLimits(1, v0andBins);
1664 evCont->SetBinLimits(2, centralityBins);
1665 evCont->SetBinLimits(3, contributorsBins);
1666 delete[] vertexBins; delete[] v0andBins; delete[] centralityBins; delete[] contributorsBins;
1668 fCFM->SetEventContainer(evCont);
1671 //____________________________________________________________
1672 void AliAnalysisTaskHFE::MakeParticleContainer(){
1674 // Create the particle container for the correction framework manager and
1677 if(!fContainer) fContainer = new AliHFEcontainer("trackContainer");
1678 fVarManager->DefineVariables(fContainer);
1680 // Create Correction Framework containers
1681 fContainer->CreateContainer("MCTrackCont", "Track Container filled with MC information", AliHFEcuts::kNcutStepsMCTrack);
1682 fContainer->CreateContainer("recTrackContReco", "Track Container filled with MC information", AliHFEcuts::kNcutStepsRecTrack + fPID->GetNumberOfPIDdetectors());
1683 fContainer->CreateContainer("recTrackContMC", "Track Container filled with MC information", AliHFEcuts::kNcutStepsRecTrack + fPID->GetNumberOfPIDdetectors());
1685 fContainer->CreateContainer("hadronicBackground", "Container for Hadronic Background", 2);
1686 fContainer->CreateContainer("recTrackContDEReco", "Container for displaced electron analysis with Reco information", 1);
1687 fContainer->CreateContainer("recTrackContDEMC", "Container for displaced electron analysis with MC information", 1);
1688 fContainer->CreateContainer("recTrackContSecvtxReco", "Container for secondary vertexing analysis with Reco information", 1);
1689 fContainer->CreateContainer("recTrackContSecvtxMC", "Container for secondary vertexing analysis with MC information", 1);
1692 fContainer->CreateContainer("conversionElecs", "Container for weighted conversion electrons",7);
1693 fContainer->CreateContainer("mesonElecs", "Container for weighted electrons from meson decays",7);
1694 fContainer->Sumw2("conversionElecs");
1695 fContainer->Sumw2("mesonElecs");
1697 if(fisNonHFEsystematics){
1698 const Char_t *sourceName[kElecBgSpecies]={"Pion","Eta","Omega","Phi","EtaPrime","Rho"};
1699 const Char_t *levelName[kBgLevels]={"Best","Lower","Upper"};
1700 for(Int_t iSource = 0; iSource < kElecBgSpecies; iSource++){
1701 for(Int_t iLevel = 0; iLevel < kBgLevels; iLevel++){
1702 fContainer->CreateContainer(Form("conversionElecs%s%s",sourceName[iSource],levelName[iLevel]), Form("Container for weighted conversion electrons from %s grandm., %s level",sourceName[iSource],levelName[iLevel]),4);
1703 fContainer->CreateContainer(Form("mesonElecs%s%s",sourceName[iSource],levelName[iLevel]), Form("Container for weighted electrons from %s decays, %s level",sourceName[iSource],levelName[iLevel]),4);
1704 fContainer->Sumw2(Form("conversionElecs%s%s",sourceName[iSource],levelName[iLevel]));
1705 fContainer->Sumw2(Form("mesonElecs%s%s",sourceName[iSource],levelName[iLevel]));
1709 //fContainer->CreateContainer("charmElecs", "Container for weighted charm electrons",2);
1712 fContainer->CreateCorrelationMatrix("correlationstepafterPID","THnSparse with correlations");
1713 fContainer->CreateCorrelationMatrix("correlationstepafterDE","THnSparse with correlations");
1714 if(!fVarManager->IsVariableDefined("centrality")) {
1715 //printf("Create the two other correlation maps\n");
1716 fContainer->CreateCorrelationMatrix("correlationstepbeforePID","THnSparse with correlations");
1717 fContainer->CreateCorrelationMatrix("correlationstepafterTOF","THnSparse with correlations");
1720 // Define the step names
1721 for(UInt_t istep = 0; istep < AliHFEcuts::kNcutStepsMCTrack; istep++){
1722 fContainer->SetStepTitle("MCTrackCont", AliHFEcuts::MCCutName(istep), istep);
1724 for(UInt_t istep = 0; istep < AliHFEcuts::kNcutStepsRecTrack; istep++){
1725 fContainer->SetStepTitle("recTrackContReco", AliHFEcuts::RecoCutName(istep), istep);
1726 fContainer->SetStepTitle("recTrackContMC", AliHFEcuts::RecoCutName(istep), istep);
1728 for(UInt_t ipid = 0; ipid < fPID->GetNumberOfPIDdetectors(); ipid++){
1729 fContainer->SetStepTitle("recTrackContReco", fPID->SortedDetectorName(ipid), AliHFEcuts::kNcutStepsRecTrack + ipid);
1730 fContainer->SetStepTitle("recTrackContMC", fPID->SortedDetectorName(ipid), AliHFEcuts::kNcutStepsRecTrack + ipid);
1734 //____________________________________________________________
1735 void AliAnalysisTaskHFE::InitPIDperformanceQA(){
1736 // Add a histogram for Fake electrons
1738 Int_t nBin[nDim] = {40, 8, 18, 2, 3, 3};
1739 //number of variables on the grid:pt,eta,phi,charge,
1740 const Double_t kPtbound[2] = {0.1, 20.};
1741 const Double_t kEtabound[2] = {-0.8, 0.8};
1742 const Double_t kPhibound[2] = {0., 2. * TMath::Pi()};
1743 const Double_t kChargebound[2] = {-1.1, 1.1};
1744 const Double_t kAddInf1bound[2] = {0., 3.};
1745 const Double_t kAddInf2bound[2] = {0., 3.};
1746 Double_t minima[nDim] = {kPtbound[0], kEtabound[0], kPhibound[0], kChargebound[0], kAddInf1bound[0], kAddInf2bound[0]};
1747 Double_t maxima[nDim] = {kPtbound[1], kEtabound[1], kPhibound[1], kChargebound[1], kAddInf1bound[1], kAddInf2bound[1]};
1749 fQACollection->CreateTHnSparse("PIDperformance", "PID performance; pT [GeV/c]; theta [rad]; phi [rad]; charge; type (0 - not el, 1 - other el, 2 - HF el; flavor (0 - no, 1 - charm, 2 - bottom)", nDim, nBin, minima, maxima);
1750 fQACollection->CreateTHnSparse("SignalToBackgroundMC", "PID performance; pT [GeV/c]; theta [rad]; phi [rad]; charge; flavor (0 - no, 1 - charm, 2 - bottom); ITS Cluster (0 - no, 1 - first (and maybe second), 2 - second)", nDim, nBin, minima, maxima);
1752 fQACollection->BinLogAxis("PIDperformance", 0);
1753 fQACollection->BinLogAxis("SignalToBackgroundMC", 0);
1754 fQACollection->Sumw2("PIDperformance");
1755 fQACollection->Sumw2("SignalToBackgroundMC");
1758 //____________________________________________________________
1759 void AliAnalysisTaskHFE::InitContaminationQA(){
1761 // Add QA for Impact Parameter cut
1764 TObjArray *array = fVarManager->GetVariables();
1765 Int_t nvars = array->GetEntriesFast();
1766 for(Int_t v = 0; v < nvars; v++) {
1767 AliHFEvarManager::AliHFEvariable *variable = (AliHFEvarManager::AliHFEvariable *) array->At(v);
1768 if(!variable) continue;
1769 TString name(((AliHFEvarManager::AliHFEvariable *)variable)->GetName());
1770 if(!name.CompareTo("pt")) {
1771 const Int_t nBinPt = variable->GetNumberOfBins();
1772 const Double_t *kPtRange = variable->GetBinning();
1774 fQACollection->CreateTH1Farray("hadronsBeforeIPcut", "Hadrons before IP cut", nBinPt, kPtRange);
1775 fQACollection->CreateTH1Farray("hadronsAfterIPcut", "Hadrons after IP cut", nBinPt, kPtRange);
1777 fQACollection->CreateTH2Farray("Ke3Kecorr", "Ke3 decay e and K correlation; Ke3K p_{t}; Ke3e p_{t}; ", nBinPt, kPtRange, 20,0.,20.);
1778 fQACollection->CreateTH2Farray("Ke3K0Lecorr", "Ke3 decay e and K0L correlation; Ke3K0L p_{t}; Ke3e p_{t}; ", nBinPt, kPtRange, 20,0.,20.);
1779 fQACollection->CreateTH1Farray("Kptspectra", "Charged Kaons: MC p_{t} ", nBinPt, kPtRange);
1780 fQACollection->CreateTH1Farray("K0Lptspectra", "K0L: MC p_{t} ", nBinPt, kPtRange);
1782 const Double_t kDCAbound[2] = {-0.2, 0.2};
1784 const Int_t nDimDca=6;
1785 const Int_t nBinDca[nDimDca] = { 8, nBinPt, 800, 12, 6, 2};
1786 Double_t minimaDca[nDimDca] = { -1., 0., kDCAbound[0], -1., -1, -1.1};
1787 Double_t maximaDca[nDimDca] = { 7., 20., kDCAbound[1], 11., 5, 1.1};
1789 Double_t *sourceBins = AliHFEtools::MakeLinearBinning(nBinDca[0], minimaDca[0], maximaDca[0]);
1790 Double_t *dcaBins = AliHFEtools::MakeLinearBinning(nBinDca[2], minimaDca[2], maximaDca[2]);
1791 Double_t *centralityBins = AliHFEtools::MakeLinearBinning(nBinDca[3], minimaDca[3], maximaDca[3]);
1792 Double_t *v0PIDBins = AliHFEtools::MakeLinearBinning(nBinDca[4], minimaDca[4], maximaDca[4]);
1793 Double_t *chargeBins = AliHFEtools::MakeLinearBinning(nBinDca[5], minimaDca[5], maximaDca[5]);
1795 fQACollection->CreateTHnSparseNoLimits("Dca", "Dca; source (0-all, 1-charm,etc); pT [GeV/c]; dca; centrality bin; v0pid; charge", nDimDca, nBinDca);
1796 ((THnSparse*)(fQACollection->Get("Dca")))->SetBinEdges(0, sourceBins);
1797 ((THnSparse*)(fQACollection->Get("Dca")))->SetBinEdges(1, kPtRange);
1798 ((THnSparse*)(fQACollection->Get("Dca")))->SetBinEdges(2, dcaBins);
1799 ((THnSparse*)(fQACollection->Get("Dca")))->SetBinEdges(3, centralityBins);
1800 ((THnSparse*)(fQACollection->Get("Dca")))->SetBinEdges(4, v0PIDBins);
1801 ((THnSparse*)(fQACollection->Get("Dca")))->SetBinEdges(5, chargeBins);
1809 //____________________________________________________________
1810 void AliAnalysisTaskHFE::InitHistoITScluster(){
1812 // Initialize a temporary histogram to monitor the chi2/ITS cluster
1814 const Int_t kNDim = 7;
1815 const Int_t kNBins[kNDim] = {88, 20,90,11, 160, 2, 1000};
1816 const Double_t kMin[kNDim] = {0.1, -1,0, 0.,0., 0, 0.};
1817 const Double_t kMax[kNDim] = {20., 1, 2.*TMath::Pi(), 11.,160, 2, 100.};
1818 fQACollection->CreateTHnSparse("fChi2perITScluster", "chi2/ITS cluster; p_{T} (GeV/c);eta;phi; centrality class;nclus;sharebit; #chi^{2}/ITS cluster", kNDim, kNBins, kMin, kMax);
1819 fQACollection->BinLogAxis("fChi2perITScluster", 0);
1823 const Int_t kNDim = 3;
1824 const Int_t kNBins[kNDim] = {44, 11, 1000};
1825 const Double_t kMin[kNDim] = {0.1, 0., 0.};
1826 const Double_t kMax[kNDim] = {20., 11., 100.};
1827 fQACollection->CreateTHnSparse("fChi2perITScluster", "chi2/ITS cluster; p_{T} (GeV/c); centrality class; #chi^{2}/ITS cluster", kNDim, kNBins, kMin, kMax);
1828 fQACollection->BinLogAxis("fChi2perITScluster", 0);
1832 //____________________________________________________________
1833 void AliAnalysisTaskHFE::SelectSpecialTrigger(const Char_t *trgclust, Int_t runMin, Int_t runMax){
1835 // Select only events triggered by a special trigeer cluster
1837 if(!fSpecialTrigger) fSpecialTrigger = new AliOADBContainer("SpecialTrigger");
1838 fSpecialTrigger->AppendObject(new TObjString(trgclust), runMin, runMax);
1841 //____________________________________________________________
1842 const Char_t * AliAnalysisTaskHFE::GetSpecialTrigger(Int_t run){
1844 // Derive selected trigger string for given run
1846 if(!fSpecialTrigger) return NULL;
1847 TObjString *trg = dynamic_cast<TObjString *>(fSpecialTrigger->GetObject(run));
1848 if(!trg) return NULL;
1849 return trg->String().Data();
1852 //____________________________________________________________
1853 void AliAnalysisTaskHFE::PrintStatus() const {
1855 // Print Analysis status
1857 printf("\n\tAnalysis Settings\n\t========================================\n\n");
1858 printf("\tSecondary Vertex finding: %s\n", GetPlugin(kSecVtx) ? "YES" : "NO");
1859 printf("\tPrimary Vertex resolution: %s\n", GetPlugin(kPriVtx) ? "YES" : "NO");
1860 printf("\tDisplaced electron analysis step: %s\n", GetPlugin(kDEstep) ? "YES" : "NO");
1861 printf("\tTagged Track Analysis: %s\n", GetPlugin(kTaggedTrackAnalysis) ? "YES" : "NO");
1863 printf("\tParticle Identification Detectors:\n");
1864 fPID->PrintStatus();
1867 printf("\t\tPID: %s\n", IsQAOn(kPIDqa) ? "YES" : "NO");
1868 printf("\t\tCUTS: %s\n", (fCuts != NULL && fCuts->IsQAOn()) ? "YES" : "NO");
1869 printf("\t\tMC: %s\n", IsQAOn(kMCqa) ? "YES" : "NO");
1873 //____________________________________________________________
1874 Bool_t AliAnalysisTaskHFE::FillProductionVertex(const AliVParticle * const track) const{
1876 // Find the production vertex of the associated MC track
1878 if(!fMCEvent) return kFALSE;
1879 const AliVParticle *mctrack = NULL;
1880 TString objectType = track->IsA()->GetName();
1881 if(objectType.CompareTo("AliESDtrack") == 0 || objectType.CompareTo("AliAODTrack") == 0){
1882 // Reconstructed track
1883 mctrack = fMCEvent->GetTrack(TMath::Abs(track->GetLabel()));
1889 if(!mctrack) return kFALSE;
1894 if(TString(mctrack->IsA()->GetName()).CompareTo("AliMCParticle") == 0){
1896 const AliMCParticle *mcpart = dynamic_cast<const AliMCParticle *>(mctrack);
1902 // case AODMCParticle
1903 const AliAODMCParticle *mcpart = dynamic_cast<const AliAODMCParticle *>(mctrack);
1910 //printf("xv %f, yv %f\n",xv,yv);
1911 fQACollection->Fill("radius", TMath::Abs(xv),TMath::Abs(yv));
1916 //__________________________________________
1917 void AliAnalysisTaskHFE::SwitchOnPlugin(Int_t plug){
1921 // - Primary vertex studies
1922 // - Secondary vertex Studies
1923 // - Post Processing
1926 case kPriVtx: SETBIT(fPlugins, plug); break;
1927 case kSecVtx: SETBIT(fPlugins, plug); break;
1928 case kIsElecBackGround: SETBIT(fPlugins, plug); break;
1929 case kPostProcess: SETBIT(fPlugins, plug); break;
1930 case kDEstep: SETBIT(fPlugins, plug); break;
1931 case kTaggedTrackAnalysis: SETBIT(fPlugins, plug); break;
1932 default: AliError("Unknown Plugin");
1935 //__________________________________________
1936 Bool_t AliAnalysisTaskHFE::ProcessCutStep(Int_t cutStep, AliVParticle *track){
1938 // Check single track cuts for a given cut step
1939 // Fill the particle container
1941 const Int_t kMCOffset = AliHFEcuts::kNcutStepsMCTrack;
1942 if(!fCFM->CheckParticleCuts(cutStep + kMCOffset, track)) return kFALSE;
1943 if(fVarManager->IsSignalTrack()) {
1944 fVarManager->FillContainer(fContainer, "recTrackContReco", cutStep, kFALSE);
1945 fVarManager->FillContainer(fContainer, "recTrackContMC", cutStep, kTRUE);
1949 //___________________________________________________
1950 Bool_t AliAnalysisTaskHFE::ReadCentrality() {
1952 // Recover the centrality of the event from ESD or AOD
1955 Float_t fCentralityLimitstemp[12];
1956 Float_t fCentralityLimitsdefault[12]= {0.,5.,10., 20., 30., 40., 50., 60.,70.,80., 90., 100.};
1957 if(!fPbPbUserCentralityBinning) memcpy(fCentralityLimitstemp,fCentralityLimitsdefault,sizeof(fCentralityLimitsdefault));
1958 else memcpy(fCentralityLimitstemp,fCentralityLimits,sizeof(fCentralityLimitsdefault));
1964 AliCentrality *centrality = fInputEvent->GetCentrality();
1965 fCentralityPercent = centrality->GetCentralityPercentile("V0M");
1966 //printf("centrality %f\n",fCentralityPercent);
1968 for(Int_t ibin = 0; ibin < 11; ibin++){
1969 if(fCentralityPercent >= fCentralityLimitstemp[ibin] && fCentralityPercent < fCentralityLimitstemp[ibin+1]){
1971 //printf("test bin %f, low %f, high %f, %d\n",fCentralityPercent,fCentralityLimitstemp[ibin],fCentralityLimitstemp[ibin+1],ibin);
1976 if(bin == -1) bin = 11; // Overflow
1978 // PP: Tracklet multiplicity, use common definition
1979 Int_t itsMultiplicity = GetITSMultiplicity(fInputEvent);
1980 Int_t multiplicityLimits[8] = {0, 1, 9, 17, 25, 36, 60, 500};
1981 for(Int_t ibin = 0; ibin < 7; ibin++){
1982 if(itsMultiplicity >= multiplicityLimits[ibin] && itsMultiplicity < multiplicityLimits[ibin + 1]){
1987 if(bin == -1) bin = 7; // Overflow
1990 AliDebug(2, Form("Centrality class %d\n", fCentralityF));
1993 // contributors, to be outsourced
1994 const AliVVertex *vtx;
1995 if(IsAODanalysis()){
1996 AliAODEvent *fAOD = dynamic_cast<AliAODEvent *>(fInputEvent);
1998 AliError("AOD Event required for AOD Analysis");
2001 vtx = fAOD->GetPrimaryVertex();
2003 AliESDEvent *fESD = dynamic_cast<AliESDEvent *>(fInputEvent);
2005 AliError("ESD Event required for ESD Analysis");
2008 vtx = fESD->GetPrimaryVertexSPD();
2011 fContributors = 0.5;
2015 Int_t contributorstemp = vtx->GetNContributors();
2016 if( contributorstemp <= 0) fContributors = 0.5;
2017 else fContributors = 1.5;
2018 //printf("Number of contributors %d\n",contributorstemp);
2023 //___________________________________________________
2024 Int_t AliAnalysisTaskHFE::GetITSMultiplicity(AliVEvent *ev){
2026 // Definition of the Multiplicity according to the JPSI group (F. Kramer)
2028 Int_t nTracklets = 0;
2030 Double_t etaRange = 1.6;
2032 if (ev->IsA() == AliAODEvent::Class()) {
2033 AliAODTracklets *tracklets = ((AliAODEvent*)ev)->GetTracklets();
2034 nTracklets = tracklets->GetNumberOfTracklets();
2035 for (Int_t nn = 0; nn < nTracklets; nn++) {
2036 Double_t theta = tracklets->GetTheta(nn);
2037 Double_t eta = -TMath::Log(TMath::Tan(theta/2.0));
2038 if (TMath::Abs(eta) < etaRange) nAcc++;
2040 } else if (ev->IsA() == AliESDEvent::Class()) {
2041 nTracklets = ((AliESDEvent*)ev)->GetMultiplicity()->GetNumberOfTracklets();
2042 for (Int_t nn = 0; nn < nTracklets; nn++) {
2043 Double_t eta = ((AliESDEvent*)ev)->GetMultiplicity()->GetEta(nn);
2044 if (TMath::Abs(eta) < etaRange) nAcc++;
2051 //___________________________________________________
2052 Bool_t AliAnalysisTaskHFE::InitializeHadronBackground(Int_t run){
2054 // Initialize background factors array from the AliOADBContainer
2055 // The container is expected to provide a TObjArray with the name
2056 // "hadronBackground" and the size 12 for the given run number
2058 AliDebug(1, "Deriving hadronic background parameterization from OADB container");
2059 TObjArray *functions = dynamic_cast<TObjArray *>(fHadronBackgroundOADB->GetObject(run, "hadronBackground"));
2061 AliDebug(1, "Content in the OADB Container is not a TObjArray");
2062 fBackGroundFactorApply = kFALSE;
2065 if(functions->GetSize() < 12){
2066 AliDebug(1, Form("Size not matching: 12 expected, %d provided", functions->GetSize()));
2067 fBackGroundFactorApply = kFALSE;
2070 for(Int_t icent = 0; icent < 12; icent++) fkBackGroundFactorArray[icent] = dynamic_cast<const TF1 *>(functions->UncheckedAt(icent));
2074 //___________________________________________________
2075 void AliAnalysisTaskHFE::RejectionPileUpVertexRangeEventCut() {
2077 // Recover the centrality of the event from ESD or AOD
2079 if(IsAODanalysis()){
2081 AliAODEvent *fAOD = dynamic_cast<AliAODEvent *>(fInputEvent);
2083 AliError("AOD Event required for AOD Analysis");
2087 fIdentifiedAsPileUp = kFALSE;
2088 if(fRemovePileUp && fAOD->IsPileupFromSPD()) fIdentifiedAsPileUp = kTRUE;
2090 fIdentifiedAsOutInz = kFALSE;
2091 //printf("Z vertex %f and out %f\n",fAOD->GetPrimaryVertex()->GetZ(),fCuts->GetVertexRange());
2092 if(TMath::Abs(fAOD->GetPrimaryVertex()->GetZ()) > fCuts->GetVertexRange()) fIdentifiedAsOutInz = kTRUE;
2094 fPassTheEventCut = kTRUE;
2095 if(!fCFM->CheckEventCuts(AliHFEcuts::kEventStepReconstructed, fAOD)) fPassTheEventCut = kFALSE;
2100 AliDebug(3, "Processing ESD Centrality");
2101 AliESDEvent *fESD = dynamic_cast<AliESDEvent *>(fInputEvent);
2103 AliError("ESD Event required for ESD Analysis");
2107 fIdentifiedAsPileUp = kFALSE;
2108 if(fRemovePileUp && fESD->IsPileupFromSPD()) fIdentifiedAsPileUp = kTRUE;
2110 fIdentifiedAsOutInz = kFALSE;
2113 if(fESD->GetPrimaryVertex()){
2114 if(TMath::Abs(fESD->GetPrimaryVertex()->GetZ()) > fCuts->GetVertexRange()) fIdentifiedAsOutInz = kTRUE;
2119 if(fESD->GetPrimaryVertexTracks()){
2120 if(TMath::Abs(fESD->GetPrimaryVertexTracks()->GetZ()) > fCuts->GetVertexRange()) fIdentifiedAsOutInz = kTRUE;
2124 fPassTheEventCut = kTRUE;
2125 if(!fCFM->CheckEventCuts(AliHFEcuts::kEventStepReconstructed, fESD)) fPassTheEventCut = kFALSE;