]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGLF/FORWARD/analysis2/AliForwardMCCorrectionsTask.cxx
Updates
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliForwardMCCorrectionsTask.cxx
CommitLineData
7984e5f7 1//
2// Calculate the corrections in the forward regions
3//
4// Inputs:
5// - AliESDEvent
6//
7// Outputs:
8// - AliAODForwardMult
9//
10// Histograms
11//
12// Corrections used
13//
0bd4b00f 14#include "AliForwardMCCorrectionsTask.h"
c8b1a7db 15#include "AliForwardCorrectionManager.h"
563a673f 16#include "AliTriggerAnalysis.h"
17#include "AliPhysicsSelection.h"
18#include "AliLog.h"
0bd4b00f 19#include "AliHeader.h"
20#include "AliGenEventHeader.h"
563a673f 21#include "AliESDEvent.h"
22#include "AliAODHandler.h"
23#include "AliMultiplicity.h"
24#include "AliInputEventHandler.h"
0bd4b00f 25#include "AliStack.h"
26#include "AliMCEvent.h"
5bb5d1f6 27#include "AliAODForwardMult.h"
28#include "AliFMDStripIndex.h"
f91358e1 29#include "AliFMDCorrSecondaryMap.h"
563a673f 30#include <TH1.h>
fb3430ac 31#include <TH2D.h>
563a673f 32#include <TDirectory.h>
33#include <TTree.h>
fb3430ac 34#include <TList.h>
f91358e1 35#include <TROOT.h>
36#include <iostream>
563a673f 37
563a673f 38//====================================================================
0bd4b00f 39AliForwardMCCorrectionsTask::AliForwardMCCorrectionsTask()
c8b1a7db 40 : AliBaseMCCorrectionsTask(),
f91358e1 41 fTrackDensity(),
42 fESDFMD(),
c8b1a7db 43 fSecCorr(0)
563a673f 44{
7984e5f7 45 //
46 // Constructor
47 //
48 // Parameters:
49 // name Name of task
50 //
563a673f 51}
52
53//____________________________________________________________________
0bd4b00f 54AliForwardMCCorrectionsTask::AliForwardMCCorrectionsTask(const char* name)
c8b1a7db 55 : AliBaseMCCorrectionsTask(name, &(AliForwardCorrectionManager::Instance())),
f91358e1 56 fTrackDensity("trackDensity"),
57 fESDFMD(),
c8b1a7db 58 fSecCorr(0)
563a673f 59{
7984e5f7 60 //
61 // Constructor
62 //
63 // Parameters:
64 // name Name of task
563a673f 65}
66
563a673f 67
68//____________________________________________________________________
c8b1a7db 69AliBaseMCCorrectionsTask::VtxBin*
70AliForwardMCCorrectionsTask::CreateVtxBin(Double_t low, Double_t high)
563a673f 71{
c8b1a7db 72 return new AliForwardMCCorrectionsTask::VtxBin(low,high, fEtaAxis);
563a673f 73}
74
75//____________________________________________________________________
c8b1a7db 76Bool_t
77AliForwardMCCorrectionsTask::PreEvent()
563a673f 78{
c8b1a7db 79 // Clear our ESD object
80 fESDFMD.Clear();
81 return true;
563a673f 82}
83
84//____________________________________________________________________
c8b1a7db 85Bool_t
86AliForwardMCCorrectionsTask::ProcessESD(const AliESDEvent& esd,
87 const AliMCEvent& mc,
88 AliBaseMCCorrectionsTask::VtxBin& bin,
89 Double_t vz)
563a673f 90{
c8b1a7db 91 AliESDFMD* esdFMD = esd.GetFMDData();
f91358e1 92
c8b1a7db 93 fTrackDensity.Calculate(*esdFMD, mc, vz, fESDFMD, bin.fPrimary);
94 bin.fCounts->Fill(0.5);
7984e5f7 95
c8b1a7db 96 AliForwardMCCorrectionsTask::VtxBin& vb =
97 static_cast<AliForwardMCCorrectionsTask::VtxBin&>(bin);
f91358e1 98
99 // And then bin the data in our vtxbin
100 for (UShort_t d=1; d<=3; d++) {
101 UShort_t nr = (d == 1 ? 1 : 2);
102 for (UShort_t q=0; q<nr; q++) {
103 Char_t r = (q == 0 ? 'I' : 'O');
104 UShort_t ns= (q == 0 ? 20 : 40);
105 UShort_t nt= (q == 0 ? 512 : 256);
c8b1a7db 106 TH2D* h = vb.fHists.Get(d,r);
f91358e1 107
108 for (UShort_t s=0; s<ns; s++) {
109 for (UShort_t t=0; t<nt; t++) {
110 Float_t mult = fESDFMD.Multiplicity(d,r,s,t);
111
112 if (mult == 0 || mult > 20) continue;
113
114 Float_t phi = fESDFMD.Phi(d,r,s,t) / 180 * TMath::Pi();
115 Float_t eta = fESDFMD.Eta(d,r,s,t);
116 h->Fill(eta,phi,mult);
117 } // for t
118 } // for s
119 } // for q
120 } // for d
c8b1a7db 121 return true;
bcd522ff 122}
563a673f 123//____________________________________________________________________
124void
c8b1a7db 125AliForwardMCCorrectionsTask::CreateCorrections(TList* results)
563a673f 126{
c8b1a7db 127 fSecCorr = new AliFMDCorrSecondaryMap;
128 fSecCorr->SetVertexAxis(fVtxAxis);
129 fSecCorr->SetEtaAxis(fEtaAxis);
130 results->Add(fSecCorr);
131}
563a673f 132
c8b1a7db 133//____________________________________________________________________
134Bool_t
135AliForwardMCCorrectionsTask::FinalizeVtxBin(AliBaseMCCorrectionsTask::VtxBin*
136 bin, UShort_t iVz)
137{
563a673f 138
c8b1a7db 139 AliForwardMCCorrectionsTask::VtxBin* vb =
140 static_cast<AliForwardMCCorrectionsTask::VtxBin*>(bin);
141 vb->Terminate(fList, fResults, iVz, fSecCorr);
142 return true;
563a673f 143}
144
c8b1a7db 145
563a673f 146//____________________________________________________________________
147void
f91358e1 148AliForwardMCCorrectionsTask::Print(Option_t* option) const
563a673f 149{
c8b1a7db 150 AliBaseMCCorrectionsTask::Print(option);
f91358e1 151 gROOT->IncreaseDirLevel();
f91358e1 152 fTrackDensity.Print(option);
153 gROOT->DecreaseDirLevel();
154}
7984e5f7 155
f91358e1 156//====================================================================
f91358e1 157AliForwardMCCorrectionsTask::VtxBin::VtxBin()
c8b1a7db 158 : AliBaseMCCorrectionsTask::VtxBin(),
159 fHists()
f91358e1 160{
161}
162//____________________________________________________________________
163AliForwardMCCorrectionsTask::VtxBin::VtxBin(Double_t low,
164 Double_t high,
165 const TAxis& axis)
c8b1a7db 166 : AliBaseMCCorrectionsTask::VtxBin(low, high, axis, 40),
167 fHists()
f91358e1 168{
169 fHists.Init(axis);
f91358e1 170}
563a673f 171
563a673f 172
f91358e1 173//____________________________________________________________________
c8b1a7db 174TList*
5934a3e3 175AliForwardMCCorrectionsTask::VtxBin::CreateOutputObjects(TList* l)
f91358e1 176{
c8b1a7db 177 TList* d = AliBaseMCCorrectionsTask::VtxBin::CreateOutputObjects(l);
f91358e1 178
179 d->Add(fHists.fFMD1i);
180 d->Add(fHists.fFMD2i);
181 d->Add(fHists.fFMD2o);
182 d->Add(fHists.fFMD3i);
183 d->Add(fHists.fFMD3o);
184
c8b1a7db 185 return d;
563a673f 186}
187
f91358e1 188//____________________________________________________________________
189TH2D*
190AliForwardMCCorrectionsTask::VtxBin::MakeBg(const TH2D* hits,
191 const TH2D* primary) const
192{
193 TH2D* h = static_cast<TH2D*>(hits->Clone());
194 h->SetDirectory(0);
195 TString n(h->GetName());
196 n.ReplaceAll("_cache", "");
197 h->SetName(n);
198 h->Divide(primary);
199
200 return h;
201}
202
563a673f 203//____________________________________________________________________
204void
5934a3e3 205AliForwardMCCorrectionsTask::VtxBin::Terminate(const TList* input,
f91358e1 206 TList* output,
207 UShort_t iVz,
208 AliFMDCorrSecondaryMap* map)
563a673f 209{
f91358e1 210 TList* out = new TList;
211 out->SetName(GetName());
212 out->SetOwner();
213 output->Add(out);
214
215 TList* l = static_cast<TList*>(input->FindObject(GetName()));
216 if (!l) {
217 AliError(Form("List %s not found in %s", GetName(), input->GetName()));
218 return;
219 }
220
221 TH2D* fmd1i = static_cast<TH2D*>(l->FindObject("FMD1I_cache"));
222 TH2D* fmd2i = static_cast<TH2D*>(l->FindObject("FMD2I_cache"));
223 TH2D* fmd2o = static_cast<TH2D*>(l->FindObject("FMD2O_cache"));
224 TH2D* fmd3i = static_cast<TH2D*>(l->FindObject("FMD3I_cache"));
225 TH2D* fmd3o = static_cast<TH2D*>(l->FindObject("FMD3O_cache"));
226 TH2D* primO = static_cast<TH2D*>(l->FindObject("primary"));
227 if (!fmd1i || !fmd2i || !fmd2o || !fmd3i || !fmd3o || !primO) {
228 AliError(Form("Missing histogram(s): %p,%p,%p,%p,%p,%p",
229 fmd1i, fmd2i, fmd2o, fmd3i, fmd3o, primO));
230 return;
231 }
232
233 // Half coverage in phi for inners
234 TH2D* primI = static_cast<TH2D*>(primO->Clone());
235 primI->SetDirectory(0);
236 primI->RebinY(2);
237
238 TH2D* bg1i = MakeBg(fmd1i, primI);
239 TH2D* bg2i = MakeBg(fmd2i, primI);
240 TH2D* bg2o = MakeBg(fmd2o, primO);
241 TH2D* bg3i = MakeBg(fmd3i, primI);
242 TH2D* bg3o = MakeBg(fmd3o, primO);
243 map->SetCorrection(1, 'I', iVz, bg1i);
244 map->SetCorrection(2, 'I', iVz, bg2i);
245 map->SetCorrection(2, 'O', iVz, bg2o);
246 map->SetCorrection(3, 'I', iVz, bg3i);
247 map->SetCorrection(3, 'O', iVz, bg3o);
248 out->Add(bg1i);
249 out->Add(bg2i);
250 out->Add(bg2o);
251 out->Add(bg3i);
252 out->Add(bg3o);
253
563a673f 254}
255
256//
257// EOF
258//