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