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