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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////
20 // Basic Analysis Task //
22 ///////////////////////////////////////////////////////////////////////////
27 #include <AliCFContainer.h>
28 #include <AliInputEventHandler.h>
29 #include <AliESDInputHandler.h>
30 #include <AliAODInputHandler.h>
31 #include <AliAnalysisManager.h>
32 #include <AliVEvent.h>
33 #include <AliTriggerAnalysis.h>
35 #include "AliDielectron.h"
36 #include "AliDielectronHistos.h"
37 #include "AliDielectronCF.h"
38 #include "AliDielectronMC.h"
39 #include "AliAnalysisTaskMultiDielectron.h"
41 ClassImp(AliAnalysisTaskMultiDielectron)
43 //_________________________________________________________________________________
44 AliAnalysisTaskMultiDielectron::AliAnalysisTaskMultiDielectron() :
49 fSelectPhysics(kFALSE),
50 fTriggerMask(AliVEvent::kMB),
51 fTriggerOnV0AND(kFALSE),
52 fRejectPileup(kFALSE),
53 fTriggerAnalysis(0x0),
62 //_________________________________________________________________________________
63 AliAnalysisTaskMultiDielectron::AliAnalysisTaskMultiDielectron(const char *name) :
64 AliAnalysisTaskSE(name),
68 fSelectPhysics(kFALSE),
69 fTriggerMask(AliVEvent::kMB),
70 fTriggerOnV0AND(kFALSE),
71 fRejectPileup(kFALSE),
72 fTriggerAnalysis(0x0),
79 DefineInput(0,TChain::Class());
80 DefineOutput(1, TList::Class());
81 DefineOutput(2, TList::Class());
82 DefineOutput(3, TH1D::Class());
83 fListHistos.SetName("Dielectron_Histos_Multi");
84 fListCF.SetName("Dielectron_CF_Multi");
85 fListDielectron.SetOwner();
86 fListHistos.SetOwner();
91 //_________________________________________________________________________________
92 void AliAnalysisTaskMultiDielectron::UserCreateOutputObjects()
95 // Add all histogram manager histogram lists to the output TList
98 if (!fListHistos.IsEmpty()||!fListCF.IsEmpty()) return; //already initialised
100 AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
101 Bool_t isESD=man->GetInputEventHandler()->IsA()==AliESDInputHandler::Class();
102 // Bool_t isAOD=man->GetInputEventHandler()->IsA()==AliAODInputHandler::Class();
104 TIter nextDie(&fListDielectron);
105 AliDielectron *die=0;
106 while ( (die=static_cast<AliDielectron*>(nextDie())) ){
108 if (die->GetHistogramList()) fListHistos.Add(const_cast<THashList*>(die->GetHistogramList()));
109 if (die->GetCFManagerPair()) fListCF.Add(const_cast<AliCFContainer*>(die->GetCFManagerPair()->GetContainer()));
112 Int_t cuts=fListDielectron.GetEntries();
113 Int_t nbins=kNbinsEvent+2*cuts;
115 fEventStat=new TH1D("hEventStat","Event statistics",nbins,0,nbins);
116 fEventStat->GetXaxis()->SetBinLabel(1,"Before Phys. Sel.");
117 fEventStat->GetXaxis()->SetBinLabel(2,"After Phys. Sel.");
118 if (fTriggerOnV0AND&&isESD) fEventStat->GetXaxis()->SetBinLabel(3,"V0and triggers");
119 if (fEventFilter) fEventStat->GetXaxis()->SetBinLabel(4,"After Event Filter");
120 if (fRejectPileup) fEventStat->GetXaxis()->SetBinLabel(5,"After Pileup rejection");
122 for (Int_t i=0; i<cuts; ++i){
123 fEventStat->GetXaxis()->SetBinLabel((kNbinsEvent+1)+2*i,Form("#splitline{1 candidate}{%s}",fListDielectron.At(i)->GetName()));
124 fEventStat->GetXaxis()->SetBinLabel((kNbinsEvent+2)+2*i,Form("#splitline{With >1 candidate}{%s}",fListDielectron.At(i)->GetName()));
128 if (!fTriggerAnalysis) fTriggerAnalysis=new AliTriggerAnalysis;
129 fTriggerAnalysis->EnableHistograms();
130 fTriggerAnalysis->SetAnalyzeMC(AliDielectronMC::Instance()->HasMC());
132 PostData(1, &fListHistos);
133 PostData(2, &fListCF);
134 PostData(3, fEventStat);
137 //_________________________________________________________________________________
138 void AliAnalysisTaskMultiDielectron::UserExec(Option_t *)
141 // Main loop. Called for every event
144 if (fListHistos.IsEmpty()&&fListCF.IsEmpty()) return;
146 AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
147 AliESDInputHandler *esdHandler=0x0;
148 Bool_t isESD=man->GetInputEventHandler()->IsA()==AliESDInputHandler::Class();
149 Bool_t isAOD=man->GetInputEventHandler()->IsA()==AliAODInputHandler::Class();
150 if ( (esdHandler=dynamic_cast<AliESDInputHandler*>(man->GetInputEventHandler())) && esdHandler->GetESDpid() ){
151 AliDielectronVarManager::SetESDpid(esdHandler->GetESDpid());
153 //load esd pid bethe bloch parameters depending on the existance of the MC handler
154 // yes: MC parameters
155 // no: data parameters
159 if (!AliDielectronVarManager::GetESDpid()){
161 if (AliDielectronMC::Instance()->HasMC()) {
162 AliDielectronVarManager::InitESDpid();
164 AliDielectronVarManager::InitESDpid(1);
170 if (!AliDielectronVarManager::GetAODpidUtil()){
171 if (AliDielectronMC::Instance()->HasMC()) {
172 AliDielectronVarManager::InitAODpidUtil();
174 AliDielectronVarManager::InitAODpidUtil(1);
179 // Was event selected ?
180 AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
181 UInt_t isSelected = AliVEvent::kAny;
182 if( fSelectPhysics && inputHandler && inputHandler->GetEventSelection() ) {
183 isSelected = inputHandler->IsEventSelected();
184 isSelected&=fTriggerMask;
187 //Before physics selection
188 fEventStat->Fill(kAllEvents);
190 PostData(3,fEventStat);
193 //after physics selection
194 fEventStat->Fill(kSelectedEvents);
197 if (fTriggerOnV0AND&&isESD){
198 if (!fTriggerAnalysis->IsOfflineTriggerFired(static_cast<AliESDEvent*>(InputEvent()), AliTriggerAnalysis::kV0AND)) return;
200 fEventStat->Fill(kV0andEvents);
204 if (!fEventFilter->IsSelected(InputEvent())) return;
206 fEventStat->Fill(kFilteredEvents);
210 if (InputEvent()->IsPileupFromSPD(3,0.8,3.,2.,5.)) return;
212 fEventStat->Fill(kPileupEvents);
215 Double_t bz = InputEvent()->GetMagneticField();
216 AliKFParticle::SetField( bz );
218 AliDielectronPID::SetCorrVal((Double_t)InputEvent()->GetRunNumber());
220 //Process event in all AliDielectron instances
221 TIter nextDie(&fListDielectron);
222 AliDielectron *die=0;
224 while ( (die=static_cast<AliDielectron*>(nextDie())) ){
225 die->Process(InputEvent());
226 if (die->HasCandidates()){
227 Int_t ncandidates=die->GetPairArray(1)->GetEntriesFast();
228 if (ncandidates==1) fEventStat->Fill((kNbinsEvent)+2*idie);
229 else if (ncandidates>1) fEventStat->Fill((kNbinsEvent+1)+2*idie);
234 PostData(1, &fListHistos);
235 PostData(2, &fListCF);
236 PostData(3,fEventStat);
239 //_________________________________________________________________________________
240 void AliAnalysisTaskMultiDielectron::FinishTaskOutput()
245 TIter nextDie(&fListDielectron);
246 AliDielectron *die=0;
247 while ( (die=static_cast<AliDielectron*>(nextDie())) ){
248 die->SaveDebugTree();