]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/dielectron/AliAnalysisTaskDielectronSE.cxx
9fe22fa0ddfce28b1c1f2c5643851586acd688df
[u/mrichter/AliRoot.git] / PWG3 / dielectron / AliAnalysisTaskDielectronSE.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 //                      for Dielectron Analysis                          //
20 //                                                                       //
21 ///////////////////////////////////////////////////////////////////////////
22
23 #include <TChain.h>
24 #include <TH1D.h>
25
26 #include <AliCFContainer.h>
27 #include <AliVEvent.h>
28 #include <AliInputEventHandler.h>
29 #include <AliESDInputHandler.h>
30 #include <AliAnalysisManager.h>
31
32 #include "AliDielectron.h"
33 #include "AliDielectronHistos.h"
34 #include "AliDielectronCF.h"
35 #include "AliAnalysisTaskDielectronSE.h"
36
37 ClassImp(AliAnalysisTaskDielectronSE)
38
39 //_________________________________________________________________________________
40 AliAnalysisTaskDielectronSE::AliAnalysisTaskDielectronSE() :
41   AliAnalysisTaskSE(),
42   fDielectron(0),
43   fSelectPhysics(kFALSE),
44   fTriggerMask(AliVEvent::kMB),
45   fEventStat(0x0)
46 {
47   //
48   // Constructor
49   //
50 }
51
52 //_________________________________________________________________________________
53 AliAnalysisTaskDielectronSE::AliAnalysisTaskDielectronSE(const char *name) :
54   AliAnalysisTaskSE(name),
55   fDielectron(0),
56   fSelectPhysics(kFALSE),
57   fTriggerMask(AliVEvent::kMB),
58   fEventStat(0x0)
59 {
60   //
61   // Constructor
62   //
63   DefineInput(0,TChain::Class());
64   DefineOutput(1, THashList::Class());
65   DefineOutput(2, AliCFContainer::Class());
66   DefineOutput(3, TH1D::Class());
67 }
68
69 //_________________________________________________________________________________
70 void AliAnalysisTaskDielectronSE::UserCreateOutputObjects()
71 {
72   //
73   // Initialise the framework objects
74   //
75   if (!fDielectron){
76     AliError("No Dielectron framework object set !!!");
77     return;
78   }
79   fDielectron->Init();
80   if (fDielectron->GetHistogramList()){
81     PostData(1, const_cast<THashList*>(fDielectron->GetHistogramList()));
82   }
83   if (fDielectron->GetCFManagerPair()){
84     PostData(2, const_cast<AliCFContainer*>(fDielectron->GetCFManagerPair()->GetContainer()));
85   }
86   
87   if (!fEventStat){
88     fEventStat=new TH1D("hEventStat","Event statistics",5,0,5);
89     fEventStat->GetXaxis()->SetBinLabel(1,"Before Phys. Sel.");
90     fEventStat->GetXaxis()->SetBinLabel(2,"After Phys. Sel.");
91   }
92   
93   PostData(3,fEventStat);
94   
95 }
96
97 //_________________________________________________________________________________
98 void AliAnalysisTaskDielectronSE::UserExec(Option_t *)
99 {
100   //
101   // Main loop. Called for every event
102   //
103
104   if (!fDielectron) return;
105   
106   AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
107   AliESDInputHandler *esdHandler=0x0;
108   if ( (esdHandler=dynamic_cast<AliESDInputHandler*>(man->GetInputEventHandler())) && esdHandler->GetESDpid() ){
109     AliDielectronVarManager::SetESDpid(esdHandler->GetESDpid());
110   } else {
111     //load esd pid bethe bloch parameters depending on the existance of the MC handler
112     // yes: MC parameters
113     // no:  data parameters
114     if (!AliDielectronVarManager::GetESDpid()){
115       if (AliDielectronMC::Instance()->HasMC()) {
116         AliDielectronVarManager::InitESDpid();
117       } else {
118         AliDielectronVarManager::InitESDpid(1);
119       }
120     }
121   }
122   // Was event selected ?
123   AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
124   UInt_t isSelected = AliVEvent::kAny;
125   if( fSelectPhysics && inputHandler && inputHandler->GetEventSelection() ) {
126     isSelected = inputHandler->IsEventSelected();
127     isSelected&=fTriggerMask;
128   }
129   
130   //Before physics selection
131   fEventStat->Fill(0.);
132   if (isSelected==0) {
133     PostData(3,fEventStat);
134     return;
135   }
136   //after physics selection
137   fEventStat->Fill(1.);
138   
139   //bz for AliKF
140   Double_t bz = InputEvent()->GetMagneticField();
141   AliKFParticle::SetField( bz );
142   
143   fDielectron->Process(InputEvent());
144
145   if (fDielectron->GetHistogramList()){
146     PostData(1, const_cast<THashList*>(fDielectron->GetHistogramList()));
147   }
148   if (fDielectron->GetCFManagerPair()){
149     PostData(2, const_cast<AliCFContainer*>(fDielectron->GetCFManagerPair()->GetContainer()));
150   }
151   PostData(3,fEventStat);
152 }
153