]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG3/dielectron/AliAnalysisTaskDielectronME.cxx
Fix Coverity reports
[u/mrichter/AliRoot.git] / PWG3 / dielectron / AliAnalysisTaskDielectronME.cxx
CommitLineData
ffbede40 1/*************************************************************************
2* Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
3* *
4* Author: The ALICE Off-line Project. *
5* Contributors are mentioned in the code where appropriate. *
6* *
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**************************************************************************/
15
27de2dfb 16/* $Id$ */
17
ffbede40 18///////////////////////////////////////////////////////////////////////////
19// //
20// Basic Analysis Task //
21// //
22///////////////////////////////////////////////////////////////////////////
23
24#include <TChain.h>
25#include <TH1D.h>
26
27#include <AliCFContainer.h>
28#include <AliInputEventHandler.h>
29#include <AliESDInputHandler.h>
30#include <AliAnalysisManager.h>
31#include <AliVEvent.h>
32
33#include "AliDielectron.h"
34#include "AliDielectronHistos.h"
35#include "AliDielectronCF.h"
36#include "AliDielectronMC.h"
37#include "AliAnalysisTaskDielectronME.h"
38
39ClassImp(AliAnalysisTaskDielectronME)
40
41//_________________________________________________________________________________
42AliAnalysisTaskDielectronME::AliAnalysisTaskDielectronME() :
43 AliAnalysisTaskME(),
44 fListDielectron(),
45 fListHistos(),
46 fListCF(),
47 fPoolDepth(2),
48 fSelectPhysics(kFALSE),
49 fTriggerMask(AliVEvent::kMB),
50 fEventStat(0x0)
51{
52 //
53 // Constructor
54 //
55}
56
57//_________________________________________________________________________________
58AliAnalysisTaskDielectronME::AliAnalysisTaskDielectronME(const char *name) :
59 AliAnalysisTaskME(name),
60 fListDielectron(),
61 fListHistos(),
62 fListCF(),
63 fPoolDepth(2),
64 fSelectPhysics(kFALSE),
65 fTriggerMask(AliVEvent::kMB),
66 fEventStat(0x0)
67{
68 //
69 // Constructor
70 //
71 DefineInput(0,TChain::Class());
72 DefineOutput(1, TList::Class());
73 DefineOutput(2, TList::Class());
74 DefineOutput(3, TH1D::Class());
75 fListHistos.SetName("Dielectron_Histos_Multi");
76 fListCF.SetName("Dielectron_CF_Multi");
77}
78
79
80//_________________________________________________________________________________
81void AliAnalysisTaskDielectronME::UserCreateOutputObjects()
82{
83 //
84 // Add all histogram manager histogram lists to the output TList
85 //
86
87 if (!fListHistos.IsEmpty()||!fListCF.IsEmpty()) return; //already initialised
88
89 TIter nextDie(&fListDielectron);
90 AliDielectron *die=0;
91 while ( (die=static_cast<AliDielectron*>(nextDie())) ){
92 die->Init();
93 if (die->GetHistogramList()) fListHistos.Add(const_cast<THashList*>(die->GetHistogramList()));
94 if (die->GetCFManagerPair()) fListCF.Add(const_cast<AliCFContainer*>(die->GetCFManagerPair()->GetContainer()));
95 }
96
97 Int_t cuts=fListDielectron.GetEntries();
98 Int_t nbins=2+2*cuts;
99 if (!fEventStat){
100 fEventStat=new TH1D("hEventStat","Event statistics",nbins,0,nbins);
101 fEventStat->GetXaxis()->SetBinLabel(1,"Before Phys. Sel.");
102 fEventStat->GetXaxis()->SetBinLabel(2,"After Phys. Sel.");
103 for (Int_t i=0; i<cuts; ++i){
104 fEventStat->GetXaxis()->SetBinLabel(3+2*i,Form("#splitline{1 candidate}{%s}",fListDielectron.At(i)->GetName()));
105 fEventStat->GetXaxis()->SetBinLabel(4+2*i,Form("#splitline{With >1 candidate}{%s}",fListDielectron.At(i)->GetName()));
106 }
107 }
108
109 PostData(1, &fListHistos);
110 PostData(2, &fListCF);
111 PostData(3, fEventStat);
112}
113
114//_________________________________________________________________________________
115void AliAnalysisTaskDielectronME::UserExec(Option_t *)
116{
117 //
118 // Main loop. Called for every event
119 //
120
121 if (fListHistos.IsEmpty()&&fListCF.IsEmpty()) return;
122
123 AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
124 AliESDInputHandler *esdHandler=0x0;
125 if ( (esdHandler=dynamic_cast<AliESDInputHandler*>(man->GetInputEventHandler())) && esdHandler->GetESDpid() ){
126 AliDielectronVarManager::SetESDpid(esdHandler->GetESDpid());
127 } else {
128 //load esd pid bethe bloch parameters depending on the existance of the MC handler
129 // yes: MC parameters
130 // no: data parameters
131 if (!AliDielectronVarManager::GetESDpid()){
132 if (AliDielectronMC::Instance()->HasMC()) {
133 AliDielectronVarManager::InitESDpid();
134 } else {
135 AliDielectronVarManager::InitESDpid(1);
136 }
137 }
138 }
139 // Was event selected ?
140 AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
141 UInt_t isSelected = AliVEvent::kAny;
142 if( fSelectPhysics && inputHandler && inputHandler->GetEventSelection() ) {
143 isSelected = inputHandler->IsEventSelected();
144 isSelected&=fTriggerMask;
145 }
146
147 //Before physics selection
148 fEventStat->Fill(0.);
149 if (isSelected==0) {
150 PostData(3,fEventStat);
151 return;
152 }
153 //after physics selection
154 fEventStat->Fill(1.);
155
156 //bz for AliKF
157 Double_t bz = GetEvent(0)->GetMagneticField();
158 AliKFParticle::SetField( bz );
159
160 //Process event in all AliDielectron instances
161 TIter nextDie(&fListDielectron);
162 AliDielectron *die=0;
163 Int_t idie=0;
164 while ( (die=static_cast<AliDielectron*>(nextDie())) ){
165 for (Int_t evt1=0; evt1<fPoolDepth-1; evt1++){
166 for (Int_t evt2=evt1+1; evt2<fPoolDepth; evt2++){
167 die->Process((AliESDEvent*)GetEvent(evt1),(AliESDEvent*)GetEvent(evt2));
168 if (die->HasCandidates()){
169 Int_t ncandidates=die->GetPairArray(1)->GetEntriesFast();
170 if (ncandidates==1) fEventStat->Fill(3+2*idie);
171 else if (ncandidates>1) fEventStat->Fill(4+2*idie);
172 }
173 }
174 }
175 ++idie;
176 }
177
178 PostData(1, &fListHistos);
179 PostData(2, &fListCF);
180 PostData(3,fEventStat);
181}
182
183//_________________________________________________________________________________
184void AliAnalysisTaskDielectronME::FinishTaskOutput()
185{
186 //
187 // Write debug tree
188 //
189 TIter nextDie(&fListDielectron);
190 AliDielectron *die=0;
191 while ( (die=static_cast<AliDielectron*>(nextDie())) ){
192 die->SaveDebugTree();
193 }
194}
195