]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/dielectron/AliAnalysisTaskMultiDielectron.cxx
framework updates, including Track rotation method (Jens); updated filtering macro
[u/mrichter/AliRoot.git] / PWG3 / dielectron / AliAnalysisTaskMultiDielectron.cxx
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
16 /* $Id$ */
17
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 <AliAODInputHandler.h>
31 #include <AliAnalysisManager.h>
32 #include <AliVEvent.h>
33 #include <AliTriggerAnalysis.h>
34
35 #include "AliDielectron.h"
36 #include "AliDielectronHistos.h"
37 #include "AliDielectronCF.h"
38 #include "AliDielectronMC.h"
39 #include "AliAnalysisTaskMultiDielectron.h"
40
41 ClassImp(AliAnalysisTaskMultiDielectron)
42
43 //_________________________________________________________________________________
44 AliAnalysisTaskMultiDielectron::AliAnalysisTaskMultiDielectron() :
45   AliAnalysisTaskSE(),
46   fListDielectron(),
47   fListHistos(),
48   fListCF(),
49   fSelectPhysics(kFALSE),
50   fTriggerMask(AliVEvent::kMB),
51   fTriggerOnV0AND(kFALSE),
52   fRejectPileup(kFALSE),
53   fTriggerAnalysis(0x0),
54   fEventFilter(0x0),
55   fEventStat(0x0)
56 {
57   //
58   // Constructor
59   //
60 }
61
62 //_________________________________________________________________________________
63 AliAnalysisTaskMultiDielectron::AliAnalysisTaskMultiDielectron(const char *name) :
64   AliAnalysisTaskSE(name),
65   fListDielectron(),
66   fListHistos(),
67   fListCF(),
68   fSelectPhysics(kFALSE),
69   fTriggerMask(AliVEvent::kMB),
70   fTriggerOnV0AND(kFALSE),
71   fRejectPileup(kFALSE),
72   fTriggerAnalysis(0x0),
73   fEventFilter(0x0),
74   fEventStat(0x0)
75 {
76   //
77   // Constructor
78   //
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();
87   fListCF.SetOwner();
88 }
89
90
91 //_________________________________________________________________________________
92 void AliAnalysisTaskMultiDielectron::UserCreateOutputObjects()
93 {
94   //
95   // Add all histogram manager histogram lists to the output TList
96   //
97
98   if (!fListHistos.IsEmpty()||!fListCF.IsEmpty()) return; //already initialised
99
100   AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
101   Bool_t isESD=man->GetInputEventHandler()->IsA()==AliESDInputHandler::Class();
102 //   Bool_t isAOD=man->GetInputEventHandler()->IsA()==AliAODInputHandler::Class();
103   
104   TIter nextDie(&fListDielectron);
105   AliDielectron *die=0;
106   while ( (die=static_cast<AliDielectron*>(nextDie())) ){
107     die->Init();
108     if (die->GetHistogramList()) fListHistos.Add(const_cast<THashList*>(die->GetHistogramList()));
109     if (die->GetCFManagerPair()) fListCF.Add(const_cast<AliCFContainer*>(die->GetCFManagerPair()->GetContainer()));
110   }
111
112   Int_t cuts=fListDielectron.GetEntries();
113   Int_t nbins=kNbinsEvent+2*cuts;
114   if (!fEventStat){
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");
121     
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()));
125     }
126   }
127
128   if (!fTriggerAnalysis) fTriggerAnalysis=new AliTriggerAnalysis;
129   fTriggerAnalysis->EnableHistograms();
130   fTriggerAnalysis->SetAnalyzeMC(AliDielectronMC::Instance()->HasMC());
131   
132   PostData(1, &fListHistos);
133   PostData(2, &fListCF);
134   PostData(3, fEventStat);
135 }
136
137 //_________________________________________________________________________________
138 void AliAnalysisTaskMultiDielectron::UserExec(Option_t *)
139 {
140   //
141   // Main loop. Called for every event
142   //
143
144   if (fListHistos.IsEmpty()&&fListCF.IsEmpty()) return;
145
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());
152   } else {
153     //load esd pid bethe bloch parameters depending on the existance of the MC handler
154     // yes: MC parameters
155     // no:  data parameters
156
157     //ESD case
158     if (isESD){
159       if (!AliDielectronVarManager::GetESDpid()){
160         
161         if (AliDielectronMC::Instance()->HasMC()) {
162           AliDielectronVarManager::InitESDpid();
163         } else {
164           AliDielectronVarManager::InitESDpid(1);
165         }
166       }
167     }
168     //AOD case
169     if (isAOD){
170       if (!AliDielectronVarManager::GetAODpidUtil()){
171         if (AliDielectronMC::Instance()->HasMC()) {
172           AliDielectronVarManager::InitAODpidUtil();
173         } else {
174           AliDielectronVarManager::InitAODpidUtil(1);
175         }
176       }
177     }
178   } 
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;
185   }
186   
187   //Before physics selection
188   fEventStat->Fill(kAllEvents);
189   if (isSelected==0) {
190     PostData(3,fEventStat);
191     return;
192   }
193   //after physics selection
194   fEventStat->Fill(kSelectedEvents);
195
196   //V0and
197   if (fTriggerOnV0AND&&isESD){
198     if (!fTriggerAnalysis->IsOfflineTriggerFired(static_cast<AliESDEvent*>(InputEvent()), AliTriggerAnalysis::kV0AND)) return;
199   }
200   fEventStat->Fill(kV0andEvents);
201   
202   //event filter
203   if (fEventFilter) {
204     if (!fEventFilter->IsSelected(InputEvent())) return;
205   }
206   fEventStat->Fill(kFilteredEvents);
207   
208   //pileup
209   if (fRejectPileup){
210     if (InputEvent()->IsPileupFromSPD(3,0.8,3.,2.,5.)) return;
211   }
212   fEventStat->Fill(kPileupEvents);
213   
214   //bz for AliKF
215   Double_t bz = InputEvent()->GetMagneticField();
216   AliKFParticle::SetField( bz );
217
218   AliDielectronPID::SetCorrVal((Double_t)InputEvent()->GetRunNumber());
219   
220   //Process event in all AliDielectron instances
221   TIter nextDie(&fListDielectron);
222   AliDielectron *die=0;
223   Int_t idie=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);
230     }
231     ++idie;
232   }
233   
234   PostData(1, &fListHistos);
235   PostData(2, &fListCF);
236   PostData(3,fEventStat);
237 }
238
239 //_________________________________________________________________________________
240 void AliAnalysisTaskMultiDielectron::FinishTaskOutput()
241 {
242   //
243   // Write debug tree
244   //
245   TIter nextDie(&fListDielectron);
246   AliDielectron *die=0;
247   while ( (die=static_cast<AliDielectron*>(nextDie())) ){
248     die->SaveDebugTree();
249   }
250 }
251