]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FORWARD/analysis2/AliForwardMultiplicity.cxx
Added some more scripts
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / AliForwardMultiplicity.cxx
CommitLineData
7e4038b5 1#include "AliForwardMultiplicity.h"
2#include "AliTriggerAnalysis.h"
3#include "AliPhysicsSelection.h"
4#include "AliLog.h"
5#include "AliFMDAnaParameters.h"
6#include "AliESDEvent.h"
7#include "AliAODHandler.h"
8#include "AliMultiplicity.h"
9#include "AliInputEventHandler.h"
10#include <TH1.h>
11#include <TDirectory.h>
12#include <TTree.h>
13
14//====================================================================
15AliForwardMultiplicity::AliForwardMultiplicity()
16 : AliAnalysisTaskSE(),
7e4038b5 17 fHData(0),
18 fFirstEvent(true),
7e4038b5 19 fESDFMD(),
20 fHistos(),
21 fAODFMD(),
fffea31d 22 fEventInspector(),
23 fEnergyFitter(),
7e4038b5 24 fSharingFilter(),
25 fDensityCalculator(),
26 fCorrections(),
27 fHistCollector(),
fffea31d 28 fList(0)
7e4038b5 29{
30}
31
32//____________________________________________________________________
33AliForwardMultiplicity::AliForwardMultiplicity(const char* name)
34 : AliAnalysisTaskSE(name),
7e4038b5 35 fHData(0),
36 fFirstEvent(true),
7e4038b5 37 fESDFMD(),
38 fHistos(),
39 fAODFMD(kTRUE),
fffea31d 40 fEventInspector("event"),
41 fEnergyFitter("energy"),
7e4038b5 42 fSharingFilter("sharing"),
43 fDensityCalculator("density"),
44 fCorrections("corrections"),
45 fHistCollector("collector"),
fffea31d 46 fList(0)
7e4038b5 47{
48 DefineOutput(1, TList::Class());
7e4038b5 49}
50
51//____________________________________________________________________
52AliForwardMultiplicity::AliForwardMultiplicity(const AliForwardMultiplicity& o)
53 : AliAnalysisTaskSE(o),
7e4038b5 54 fHData(o.fHData),
55 fFirstEvent(true),
7e4038b5 56 fESDFMD(o.fESDFMD),
57 fHistos(o.fHistos),
58 fAODFMD(o.fAODFMD),
fffea31d 59 fEventInspector(o.fEventInspector),
60 fEnergyFitter(o.fEnergyFitter),
7e4038b5 61 fSharingFilter(o.fSharingFilter),
62 fDensityCalculator(o.fDensityCalculator),
63 fCorrections(o.fCorrections),
64 fHistCollector(o.fHistCollector),
fffea31d 65 fList(o.fList)
7e4038b5 66{
67}
68
69//____________________________________________________________________
70AliForwardMultiplicity&
71AliForwardMultiplicity::operator=(const AliForwardMultiplicity& o)
72{
fffea31d 73 AliAnalysisTaskSE::operator=(o);
74
7e4038b5 75 fHData = o.fHData;
76 fFirstEvent = o.fFirstEvent;
fffea31d 77 fEventInspector = o.fEventInspector;
78 fEnergyFitter = o.fEnergyFitter;
7e4038b5 79 fSharingFilter = o.fSharingFilter;
80 fDensityCalculator = o.fDensityCalculator;
81 fCorrections = o.fCorrections;
82 fHistCollector = o.fHistCollector;
83 fHistos = o.fHistos;
84 fAODFMD = o.fAODFMD;
85 fList = o.fList;
7e4038b5 86
87 return *this;
88}
89
fffea31d 90//____________________________________________________________________
91void
92AliForwardMultiplicity::SetDebug(Int_t dbg)
93{
94 fEventInspector.SetDebug(dbg);
95 fEnergyFitter.SetDebug(dbg);
96 fSharingFilter.SetDebug(dbg);
97 fDensityCalculator.SetDebug(dbg);
98 fCorrections.SetDebug(dbg);
99 fHistCollector.SetDebug(dbg);
100}
101
7e4038b5 102//____________________________________________________________________
103void
104AliForwardMultiplicity::Init()
105{
106 fFirstEvent = true;
107}
108
109//____________________________________________________________________
110void
111AliForwardMultiplicity::InitializeSubs()
112{
113 AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
114 pars->Init(kTRUE);
115
fffea31d 116
117 TAxis e(pars->GetNetaBins(), pars->GetEtaMin(), pars->GetEtaMax());
118 TAxis v(pars->GetNvtxBins(), -pars->GetVtxCutZ(), pars->GetVtxCutZ());
119
7e4038b5 120 fHistos.Init(e);
121 fAODFMD.Init(e);
122
123 fHData = static_cast<TH2D*>(fAODFMD.GetHistogram().Clone("d2Ndetadphi"));
124 fHData->SetStats(0);
125 fHData->SetDirectory(0);
9d99b0dd 126 fList->Add(fHData);
127
fffea31d 128 fEnergyFitter.Init(e);
129 fEventInspector.Init(v);
130 fHistCollector.Init(v);
7e4038b5 131}
132
133//____________________________________________________________________
134void
135AliForwardMultiplicity::UserCreateOutputObjects()
136{
137 fList = new TList;
138
139 AliAnalysisManager* am = AliAnalysisManager::GetAnalysisManager();
fffea31d 140 AliAODHandler* ah = dynamic_cast<AliAODHandler*>(am->GetOutputEventHandler());
141 if (!ah) AliFatal("No AOD output handler set in analysis manager");
7e4038b5 142
143
144 TObject* obj = &fAODFMD;
145 ah->AddBranch("AliAODForwardMult", &obj);
146
fffea31d 147 fEventInspector.DefineOutput(fList);
148 fEnergyFitter.DefineOutput(fList);
9d99b0dd 149 fSharingFilter.DefineOutput(fList);
150 fDensityCalculator.DefineOutput(fList);
151 fCorrections.DefineOutput(fList);
7e4038b5 152}
153//____________________________________________________________________
154void
155AliForwardMultiplicity::UserExec(Option_t*)
156{
157 // Get the input data
158 AliESDEvent* esd = dynamic_cast<AliESDEvent*>(InputEvent());
159 if (!esd) {
160 AliWarning("No ESD event found for input event");
161 return;
162 }
163
7e4038b5 164 // On the first event, initialize the parameters
fffea31d 165 if (fFirstEvent && esd->GetESDRun()) {
7e4038b5 166 AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
fffea31d 167 AliInfo(Form("Initializing with parameters from the ESD:\n"
168 " AliESDEvent::GetBeamEnergy() ->%f\n"
169 " AliESDEvent::GetBeamType() ->%s\n"
170 " AliESDEvent::GetCurrentL3() ->%f\n"
171 " AliESDEvent::GetMagneticField()->%f\n"
172 " AliESDEvent::GetRunNumber() ->%d\n",
173 esd->GetBeamEnergy(),
174 esd->GetBeamType(),
175 esd->GetCurrentL3(),
176 esd->GetMagneticField(),
177 esd->GetRunNumber()));
7e4038b5 178 pars->SetParametersFromESD(esd);
179 pars->PrintStatus();
180 fFirstEvent = false;
181
182 InitializeSubs();
183 }
184 // Clear stuff
185 fHistos.Clear();
186 fESDFMD.Clear();
187 fAODFMD.Clear();
188
fffea31d 189 Bool_t lowFlux = kFALSE;
190 UInt_t triggers = 0;
191 Int_t ivz = -1;
192 Double_t vz = 0;
193 UInt_t found = fEventInspector.Process(esd, triggers, lowFlux, ivz, vz);
194 if (found & AliFMDEventInspector::kNoEvent) return;
195 if (found & AliFMDEventInspector::kNoTriggers) return;
196
197 // Set trigger bits, and mark this event for storage
9d99b0dd 198 fAODFMD.SetTriggerBits(triggers);
7e4038b5 199 MarkEventForStore();
200
fffea31d 201 if (found & AliFMDEventInspector::kNoSPD) return;
202 if (found & AliFMDEventInspector::kNoFMD) return;
203 if (found & AliFMDEventInspector::kNoVertex) return;
204 fAODFMD.SetIpZ(vz);
7e4038b5 205
fffea31d 206 if (found & AliFMDEventInspector::kBadVertex) return;
7e4038b5 207
fffea31d 208 // Get FMD data
209 AliESDFMD* esdFMD = esd->GetFMDData();
210 // Apply the sharing filter (or hit merging or clustering if you like)
211 if (!fSharingFilter.Filter(*esdFMD, lowFlux, fESDFMD)) {
212 AliWarning("Sharing filter failed!");
7e4038b5 213 return;
214 }
215
fffea31d 216 // Do the energy stuff
217 if (!fEnergyFitter.Accumulate(*esdFMD, triggers & AliAODForwardMult::kEmpty)){
218 AliWarning("Energy fitter failed");
7e4038b5 219 return;
220 }
221
222 // Calculate the inclusive charged particle density
223 if (!fDensityCalculator.Calculate(fESDFMD, fHistos, ivz, lowFlux)) {
224 AliWarning("Density calculator failed!");
225 return;
226 }
227
228 // Do the secondary and other corrections.
229 if (!fCorrections.Correct(fHistos, ivz)) {
230 AliWarning("Corrections failed");
231 return;
232 }
233
234 if (!fHistCollector.Collect(fHistos, ivz, fAODFMD.GetHistogram())) {
235 AliWarning("Histogram collector failed");
236 return;
237 }
238
239 if (fAODFMD.IsTriggerBits(AliAODForwardMult::kInel))
240 fHData->Add(&(fAODFMD.GetHistogram()));
2d68d438 241
242 PostData(1, fList);
7e4038b5 243}
244
245//____________________________________________________________________
246void
247AliForwardMultiplicity::Terminate(Option_t*)
248{
249 TList* list = dynamic_cast<TList*>(GetOutputData(1));
250 if (!list) {
2d68d438 251 AliError(Form("No output list defined (%p)", GetOutputData(1)));
252 if (GetOutputData(1)) GetOutputData(1)->Print();
7e4038b5 253 return;
254 }
9d99b0dd 255
256 // Get our histograms from the container
fffea31d 257 TH1I* hEventsTr = 0;//static_cast<TH1I*>(list->FindObject("nEventsTr"));
258 TH1I* hEventsTrVtx = 0;//static_cast<TH1I*>(list->FindObject("nEventsTrVtx"));
259 TH1I* hTriggers = 0;
260 if (!fEventInspector.FetchHistograms(list, hEventsTr,
261 hEventsTrVtx, hTriggers)) {
262 AliError(Form("Didn't get histograms from event selector "
263 "(hEventsTr=%p,hEventsTrVtx=%p)",
264 hEventsTr, hEventsTrVtx));
265 list->ls();
266 return;
267 }
268
9d99b0dd 269 TH2D* hData = static_cast<TH2D*>(list->FindObject("d2Ndetadphi"));
fffea31d 270 if (!hData) {
271 AliError(Form("Couldn't get our summed histogram from output "
272 "list %s (d2Ndetadphi=%p)", list->GetName(), hData));
2d68d438 273 list->ls();
274 return;
275 }
9d99b0dd 276
7e4038b5 277 // TH1D* dNdeta = fHData->ProjectionX("dNdeta", 0, -1, "e");
9d99b0dd 278 TH1D* dNdeta = hData->ProjectionX("dNdeta", 1, -1, "e");
279 TH1D* norm = hData->ProjectionX("dNdeta", 0, 1, "");
7e4038b5 280 dNdeta->SetTitle("dN_{ch}/d#eta in the forward regions");
281 dNdeta->SetYTitle("#frac{1}{N}#frac{dN_{ch}}{d#eta}");
282 dNdeta->Divide(norm);
283 dNdeta->SetStats(0);
9d99b0dd 284 dNdeta->Scale(Double_t(hEventsTrVtx->GetEntries())/hEventsTr->GetEntries(),
7e4038b5 285 "width");
7e4038b5 286 list->Add(dNdeta);
fffea31d 287
288 fEnergyFitter.Fit(list);
9d99b0dd 289 fSharingFilter.ScaleHistograms(list,hEventsTr->Integral());
290 fDensityCalculator.ScaleHistograms(list,hEventsTrVtx->Integral());
291 fCorrections.ScaleHistograms(list,hEventsTrVtx->Integral());
7e4038b5 292}
293
294//____________________________________________________________________
295void
296AliForwardMultiplicity::MarkEventForStore() const
297{
298 // Make sure the AOD tree is filled
299 AliAnalysisManager* am = AliAnalysisManager::GetAnalysisManager();
300 AliAODHandler* ah =
301 dynamic_cast<AliAODHandler*>(am->GetOutputEventHandler());
302 if (!ah)
303 AliFatal("No AOD output handler set in analysis manager");
304
305 ah->SetFillAOD(kTRUE);
306}
7e4038b5 307
308//____________________________________________________________________
309void
310AliForwardMultiplicity::Print(Option_t*) const
311{
312}
313
314//
315// EOF
316//