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