1 /*************************************************************************
2 * Copyright(c) 1998-2009, 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 **************************************************************************/
16 ///////////////////////////////////////////////////////////////////////////
18 // Basic Analysis Task //
20 ///////////////////////////////////////////////////////////////////////////
25 #include <AliCFContainer.h>
26 #include <AliInputEventHandler.h>
27 #include <AliESDInputHandler.h>
28 #include <AliAODHandler.h>
29 #include <AliAnalysisManager.h>
30 #include <AliVEvent.h>
31 #include <AliTriggerAnalysis.h>
33 #include "AliDielectron.h"
34 #include "AliDielectronHistos.h"
35 #include "AliDielectronCF.h"
36 #include "AliDielectronMC.h"
37 #include "AliAnalysisTaskMultiDielectron.h"
39 ClassImp(AliAnalysisTaskMultiDielectron)
41 //_________________________________________________________________________________
42 AliAnalysisTaskMultiDielectron::AliAnalysisTaskMultiDielectron() :
47 fSelectPhysics(kFALSE),
48 fTriggerMask(AliVEvent::kMB),
49 fTriggerOnV0AND(kFALSE),
50 fRejectPileup(kFALSE),
51 fTriggerAnalysis(0x0),
60 //_________________________________________________________________________________
61 AliAnalysisTaskMultiDielectron::AliAnalysisTaskMultiDielectron(const char *name) :
62 AliAnalysisTaskSE(name),
66 fSelectPhysics(kFALSE),
67 fTriggerMask(AliVEvent::kMB),
68 fTriggerOnV0AND(kFALSE),
69 fRejectPileup(kFALSE),
70 fTriggerAnalysis(0x0),
77 DefineInput(0,TChain::Class());
78 DefineOutput(1, TList::Class());
79 DefineOutput(2, TList::Class());
80 DefineOutput(3, TH1D::Class());
81 fListHistos.SetName("Dielectron_Histos_Multi");
82 fListCF.SetName("Dielectron_CF_Multi");
83 fListDielectron.SetOwner();
84 fListHistos.SetOwner();
89 //_________________________________________________________________________________
90 void AliAnalysisTaskMultiDielectron::UserCreateOutputObjects()
93 // Add all histogram manager histogram lists to the output TList
96 if (!fListHistos.IsEmpty()||!fListCF.IsEmpty()) return; //already initialised
98 AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
99 Bool_t isESD=man->GetInputEventHandler()->IsA()==AliESDInputHandler::Class();
100 // Bool_t isAOD=man->GetInputEventHandler()->IsA()==AliAODHandler::Class();
102 TIter nextDie(&fListDielectron);
103 AliDielectron *die=0;
104 while ( (die=static_cast<AliDielectron*>(nextDie())) ){
106 if (die->GetHistogramList()) fListHistos.Add(const_cast<THashList*>(die->GetHistogramList()));
107 if (die->GetCFManagerPair()) fListCF.Add(const_cast<AliCFContainer*>(die->GetCFManagerPair()->GetContainer()));
110 Int_t cuts=fListDielectron.GetEntries();
111 Int_t nbins=kNbinsEvent+2*cuts;
113 fEventStat=new TH1D("hEventStat","Event statistics",nbins,0,nbins);
114 fEventStat->GetXaxis()->SetBinLabel(1,"Before Phys. Sel.");
115 fEventStat->GetXaxis()->SetBinLabel(2,"After Phys. Sel.");
116 if (fTriggerOnV0AND&&isESD) fEventStat->GetXaxis()->SetBinLabel(3,"V0and triggers");
117 if (fEventFilter) fEventStat->GetXaxis()->SetBinLabel(4,"After Event Filter");
118 if (fRejectPileup) fEventStat->GetXaxis()->SetBinLabel(5,"After Pileup rejection");
120 for (Int_t i=0; i<cuts; ++i){
121 fEventStat->GetXaxis()->SetBinLabel((kNbinsEvent+1)+2*i,Form("#splitline{1 candidate}{%s}",fListDielectron.At(i)->GetName()));
122 fEventStat->GetXaxis()->SetBinLabel((kNbinsEvent+2)+2*i,Form("#splitline{With >1 candidate}{%s}",fListDielectron.At(i)->GetName()));
126 if (!fTriggerAnalysis) fTriggerAnalysis=new AliTriggerAnalysis;
127 fTriggerAnalysis->EnableHistograms();
128 fTriggerAnalysis->SetAnalyzeMC(AliDielectronMC::Instance()->HasMC());
130 PostData(1, &fListHistos);
131 PostData(2, &fListCF);
132 PostData(3, fEventStat);
135 //_________________________________________________________________________________
136 void AliAnalysisTaskMultiDielectron::UserExec(Option_t *)
139 // Main loop. Called for every event
142 if (fListHistos.IsEmpty()&&fListCF.IsEmpty()) return;
144 AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
145 AliESDInputHandler *esdHandler=0x0;
146 Bool_t isESD=man->GetInputEventHandler()->IsA()==AliESDInputHandler::Class();
147 Bool_t isAOD=man->GetInputEventHandler()->IsA()==AliAODHandler::Class();
148 if ( (esdHandler=dynamic_cast<AliESDInputHandler*>(man->GetInputEventHandler())) && esdHandler->GetESDpid() ){
149 AliDielectronVarManager::SetESDpid(esdHandler->GetESDpid());
151 //load esd pid bethe bloch parameters depending on the existance of the MC handler
152 // yes: MC parameters
153 // no: data parameters
157 if (!AliDielectronVarManager::GetESDpid()){
159 if (AliDielectronMC::Instance()->HasMC()) {
160 AliDielectronVarManager::InitESDpid();
162 AliDielectronVarManager::InitESDpid(1);
168 if (!AliDielectronVarManager::GetAODpidUtil()){
169 if (AliDielectronMC::Instance()->HasMC()) {
170 AliDielectronVarManager::InitAODpidUtil();
172 AliDielectronVarManager::InitAODpidUtil(1);
177 // Was event selected ?
178 AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
179 UInt_t isSelected = AliVEvent::kAny;
180 if( fSelectPhysics && inputHandler && inputHandler->GetEventSelection() ) {
181 isSelected = inputHandler->IsEventSelected();
182 isSelected&=fTriggerMask;
185 //Before physics selection
186 fEventStat->Fill(kAllEvents);
188 PostData(3,fEventStat);
191 //after physics selection
192 fEventStat->Fill(kSelectedEvents);
195 if (fTriggerOnV0AND&&isESD){
196 if (!fTriggerAnalysis->IsOfflineTriggerFired(static_cast<AliESDEvent*>(InputEvent()), AliTriggerAnalysis::kV0AND)) return;
198 fEventStat->Fill(kV0andEvents);
202 if (!fEventFilter->IsSelected(InputEvent())) return;
204 fEventStat->Fill(kFilteredEvents);
208 if (InputEvent()->IsPileupFromSPD(3,0.8,3.,2.,5.)) return;
210 fEventStat->Fill(kPileupEvents);
213 Double_t bz = InputEvent()->GetMagneticField();
214 AliKFParticle::SetField( bz );
216 AliDielectronPID::SetCorrVal((Double_t)InputEvent()->GetRunNumber());
218 //Process event in all AliDielectron instances
219 TIter nextDie(&fListDielectron);
220 AliDielectron *die=0;
222 while ( (die=static_cast<AliDielectron*>(nextDie())) ){
223 die->Process(InputEvent());
224 if (die->HasCandidates()){
225 Int_t ncandidates=die->GetPairArray(1)->GetEntriesFast();
226 if (ncandidates==1) fEventStat->Fill((kNbinsEvent+1)+2*idie);
227 else if (ncandidates>1) fEventStat->Fill((kNbinsEvent+2)+2*idie);
232 PostData(1, &fListHistos);
233 PostData(2, &fListCF);
234 PostData(3,fEventStat);
237 //_________________________________________________________________________________
238 void AliAnalysisTaskMultiDielectron::FinishTaskOutput()
243 TIter nextDie(&fListDielectron);
244 AliDielectron *die=0;
245 while ( (die=static_cast<AliDielectron*>(nextDie())) ){
246 die->SaveDebugTree();