]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliCentralMCMultiplicityTask.cxx
Some changes to make it possible to run the DA's
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliCentralMCMultiplicityTask.cxx
1 //====================================================================
2 // 
3 // Base class for classes that calculate the multiplicity in the
4 // central region event-by-event
5 // 
6 // Inputs: 
7 //   - AliESDEvent 
8 //
9 // Outputs: 
10 //   - AliAODCentralMult 
11 // 
12 // Histograms 
13 //   
14 // Corrections used 
15 #include "AliCentralMCMultiplicityTask.h"
16 #include "AliForwardCorrectionManager.h"
17 #include "AliForwardUtil.h"
18 #include "AliLog.h"
19 #include "AliAODHandler.h"
20 #include "AliInputEventHandler.h"
21 #include "AliESDInputHandler.h"
22 #include "AliAnalysisManager.h"
23 #include "AliESDEvent.h"
24 #include "AliMultiplicity.h"
25 #include "AliFMDEventInspector.h"
26 #include <AliMCEvent.h>
27 #include <AliTrackReference.h>
28 #include <AliStack.h>
29 #include <TROOT.h>
30 #include <TH1D.h>
31 #include <TH2D.h>
32 #include <TH3D.h>
33 #include <TFile.h>
34 #include <TError.h>
35 #include <iostream>
36 #include <iomanip>
37
38 //====================================================================
39 AliCentralMCMultiplicityTask::AliCentralMCMultiplicityTask(const char* name) 
40   : AliCentralMultiplicityTask(name),
41     fTrackDensity(name),
42     fAODMCCentral(kTRUE)
43 {
44   // 
45   // Constructor 
46   //   
47   DGUARD(fDebug,3,"Named CTOR of AliCentralMCMultiplicityTask: %s", 
48          name);
49   fBranchNames = 
50     "ESD:AliESDRun.,AliESDHeader.,AliMultiplicity.,"
51     "SPDVertex.,PrimaryVertex.";
52 }
53 //____________________________________________________________________
54 AliCentralMCMultiplicityTask::AliCentralMCMultiplicityTask() 
55   : AliCentralMultiplicityTask(),
56     fTrackDensity(),
57     fAODMCCentral(kTRUE)
58 {
59   // 
60   // Constructor 
61   // 
62   DGUARD(fDebug, 3,"Default CTOR of AliCentralMCMultiplicityTask");
63 }
64 //____________________________________________________________________
65 AliCentralMCMultiplicityTask::AliCentralMCMultiplicityTask(const AliCentralMCMultiplicityTask& o)
66   : AliCentralMultiplicityTask(o),
67     fTrackDensity(o.fTrackDensity),
68     fAODMCCentral(o.fAODMCCentral)
69 {
70   //
71   // Copy constructor 
72   // 
73   DGUARD(fDebug, 3,"COPY CTOR of AliCentralMCMultiplicityTask");
74 }
75 //____________________________________________________________________
76 AliCentralMCMultiplicityTask&
77 AliCentralMCMultiplicityTask::operator=(const AliCentralMCMultiplicityTask& o)
78 {
79   // 
80   // Assignment operator 
81   //
82   DGUARD(fDebug,3,"Assignment of AliCentralMCMultiplicityTask");
83   if (&o == this) return *this; 
84   AliCentralMultiplicityTask::operator=(o);
85   fAODMCCentral     = o.fAODMCCentral;
86   fTrackDensity     = o.fTrackDensity;
87   return *this;
88 }
89 //____________________________________________________________________
90 void AliCentralMCMultiplicityTask::UserCreateOutputObjects() 
91 {
92   // 
93   // Create output objects 
94   // 
95   //
96   DGUARD(fDebug,1,"Create user output in AliCentralMCMultiplicityTask");
97   AliCentralMultiplicityTask::UserCreateOutputObjects();
98
99   AliAnalysisManager* am = AliAnalysisManager::GetAnalysisManager();
100   AliAODHandler*      ah = 
101     dynamic_cast<AliAODHandler*>(am->GetOutputEventHandler());
102   if (ah)
103  {       
104         AliFatal("No AOD output handler set in analysis manager");
105   
106   
107         TObject* obj = &fAODMCCentral;
108         ah->AddBranch("AliAODCentralMult", &obj);
109
110   }
111   fTrackDensity.CreateOutputObjects(fList);
112
113 }
114 //____________________________________________________________________
115 void AliCentralMCMultiplicityTask::UserExec(Option_t* option) 
116 {
117   // 
118   // Process each event 
119   // 
120   // Parameters:
121   //    option Not used
122   //  
123   DGUARD(fDebug,1,"Process event in AliCentralMCMultiplicityTask");
124   fAODMCCentral.Clear("");
125   // Call base class 
126   AliCentralMultiplicityTask::UserExec(option);
127 fAODMCCentral.Init(*(fAODCentral.GetHistogram().GetXaxis()));
128   // check if we need this event 
129   AliAnalysisManager* am = AliAnalysisManager::GetAnalysisManager();
130   AliAODHandler*      ah = 
131     dynamic_cast<AliAODHandler*>(am->GetOutputEventHandler());
132   if (ah)
133 {  
134   //  AliFatal("No AOD output handler set in analysis manager");
135
136   // if base class did not want this event, then neither to we 
137   if (!ah->GetFillAOD() || fIvz <= 0) return;
138  } 
139   const AliMCEvent*  mcEvent = MCEvent();
140   if (!mcEvent) return;
141   TH2D&              hist    = fAODMCCentral.GetHistogram();
142
143   Double_t vz = GetManager().GetSecMap()->GetVertexAxis().GetBinCenter(fIvz);
144
145   fTrackDensity.Calculate(*mcEvent, vz, hist, NULL);
146
147   CorrectData(hist, fIvz);
148
149 }
150
151 //____________________________________________________________________
152 void AliCentralMCMultiplicityTask::Terminate(Option_t* option) 
153 {
154   // 
155   // End of job
156   // 
157   // Parameters:
158   //    option Not used 
159   //
160   DGUARD(fDebug,1,"Final analysis of merge in AliCentralMCMultiplicityTask");
161   AliCentralMultiplicityTask::Terminate(option);
162 }
163 //____________________________________________________________________
164 void
165 AliCentralMCMultiplicityTask::Print(Option_t* option) const
166 {
167   // 
168   // Print information 
169   // 
170   // Parameters:
171   //    option Not used
172   //
173   AliCentralMultiplicityTask::Print(option);
174   gROOT->IncreaseDirLevel();
175   fTrackDensity.Print(option);
176   gROOT->DecreaseDirLevel();
177 }
178 //
179 // EOF
180 //