]>
Commit | Line | Data |
---|---|---|
ea3fd2d5 | 1 | // $Id$ |
2 | ||
ea3fd2d5 | 3 | #include "AliAnalysisTaskEMCALPi0PbPb.h" |
fa443410 | 4 | #include <TAxis.h> |
717fe7de | 5 | #include <TChain.h> |
6 | #include <TClonesArray.h> | |
f5d4ab70 | 7 | #include <TFile.h> |
717fe7de | 8 | #include <TH1F.h> |
9 | #include <TH2F.h> | |
10 | #include <TList.h> | |
11 | #include <TLorentzVector.h> | |
f5d4ab70 | 12 | #include <TNtuple.h> |
717fe7de | 13 | #include "AliAODEvent.h" |
14 | #include "AliAODVertex.h" | |
15 | #include "AliAnalysisManager.h" | |
16 | #include "AliAnalysisTaskEMCALClusterizeFast.h" | |
ea3fd2d5 | 17 | #include "AliCentrality.h" |
ea3fd2d5 | 18 | #include "AliEMCALGeoUtils.h" |
717fe7de | 19 | #include "AliESDEvent.h" |
20 | #include "AliESDVertex.h" | |
43cfaa06 | 21 | #include "AliLog.h" |
ea3fd2d5 | 22 | |
23 | ClassImp(AliAnalysisTaskEMCALPi0PbPb) | |
717fe7de | 24 | |
ea3fd2d5 | 25 | //________________________________________________________________________ |
26 | AliAnalysisTaskEMCALPi0PbPb::AliAnalysisTaskEMCALPi0PbPb(const char *name) | |
27 | : AliAnalysisTaskSE(name), | |
717fe7de | 28 | fCentVar("V0M"), |
29 | fCentFrom(0), | |
30 | fCentTo(100), | |
76332037 | 31 | fVtxZMin(-10), |
32 | fVtxZMax(+10), | |
43cfaa06 | 33 | fUseQualFlag(1), |
717fe7de | 34 | fClusName(), |
f5d4ab70 | 35 | fDoNtuple(0), |
a49742b5 | 36 | fDoAfterburner(0), |
f224d35b | 37 | fAsymMax(1), |
a49742b5 | 38 | fNminCells(1), |
f224d35b | 39 | fMinErat(0), |
40 | fMinEcc(0), | |
d9f26424 | 41 | fNEvs(0), |
d595acbb | 42 | fGeom(0), |
717fe7de | 43 | fOutput(0), |
44 | fEsdEv(0), | |
45 | fAodEv(0), | |
43cfaa06 | 46 | fRecPoints(0), |
717fe7de | 47 | fEsdClusters(0), |
48 | fEsdCells(0), | |
49 | fAodClusters(0), | |
286b47a5 | 50 | fAodCells(0), |
fa443410 | 51 | fPtRanges(0), |
f5d4ab70 | 52 | fNtuple(0), |
fa443410 | 53 | fHCuts(0x0), |
54 | fHVertexZ(0x0), | |
76332037 | 55 | fHVertexZ2(0x0), |
fa443410 | 56 | fHCent(0x0), |
57 | fHCentQual(0x0), | |
d595acbb | 58 | fHColuRow(0x0), |
59 | fHColuRowE(0x0), | |
60 | fHCellMult(0x0), | |
61 | fHCellE(0x0), | |
62 | fHCellH(0x0), | |
6eb6260e | 63 | fHCellM(0x0), |
a49742b5 | 64 | fHCellM2(0x0), |
6eb6260e | 65 | fHClustEccentricity(0), |
fa443410 | 66 | fHClustEtaPhi(0x0), |
67 | fHClustEnergyPt(0x0), | |
68 | fHClustEnergySigma(0x0), | |
d595acbb | 69 | fHClustSigmaSigma(0x0), |
6eb6260e | 70 | fHClustNCellEnergyRatio(0x0), |
fa443410 | 71 | fHPionEtaPhi(0x0), |
72 | fHPionMggPt(0x0), | |
6eb6260e | 73 | fHPionMggAsym(0x0), |
74 | fHPionMggDgg(0x0) | |
ea3fd2d5 | 75 | { |
717fe7de | 76 | // Constructor. |
ea3fd2d5 | 77 | |
d595acbb | 78 | if (!name) |
79 | return; | |
80 | SetName(name); | |
ea3fd2d5 | 81 | DefineInput(0, TChain::Class()); |
ea3fd2d5 | 82 | DefineOutput(1, TList::Class()); |
fa443410 | 83 | fBranchNames="ESD:AliESDRun.,AliESDHeader.,PrimaryVertex.,SPDVertex.,TPCVertex.,EMCALCells. " |
29b496d5 | 84 | "AOD:header,vertices,emcalCells"; |
ea3fd2d5 | 85 | } |
ea3fd2d5 | 86 | |
717fe7de | 87 | //________________________________________________________________________ |
88 | AliAnalysisTaskEMCALPi0PbPb::~AliAnalysisTaskEMCALPi0PbPb() | |
89 | { | |
90 | // Destructor. | |
ea3fd2d5 | 91 | |
717fe7de | 92 | delete fOutput; fOutput = 0; |
fa443410 | 93 | delete fPtRanges; fPtRanges = 0; |
d595acbb | 94 | delete fGeom; fGeom = 0; |
95 | delete [] fHColuRow; | |
96 | delete [] fHColuRowE; | |
97 | delete [] fHCellMult; | |
ea3fd2d5 | 98 | } |
717fe7de | 99 | |
ea3fd2d5 | 100 | //________________________________________________________________________ |
101 | void AliAnalysisTaskEMCALPi0PbPb::UserCreateOutputObjects() | |
102 | { | |
717fe7de | 103 | // Create user objects here. |
ea3fd2d5 | 104 | |
717fe7de | 105 | fOutput = new TList(); |
106 | fOutput->SetOwner(); | |
ea3fd2d5 | 107 | |
f5d4ab70 | 108 | if (fDoNtuple) { |
109 | TFile *f = OpenFile(1); | |
110 | if (f) | |
6eb6260e | 111 | fNtuple = new TNtuple(Form("nt%.0fto%.0f",fCentFrom,fCentTo),"nt", |
a49742b5 | 112 | "run:evt:cent:pt:e:emax:n:db:disp:mn:ms:chi:cpv:ecc:sig:eta:phi"); |
f5d4ab70 | 113 | } |
114 | ||
d595acbb | 115 | // histograms |
a49742b5 | 116 | TH1::SetDefaultSumw2(kTRUE); |
117 | TH2::SetDefaultSumw2(kTRUE); | |
fa443410 | 118 | fHCuts = new TH1F("hEventCuts","",4,0.5,4.5); |
119 | fHCuts->GetXaxis()->SetBinLabel(1,"All (PS)"); | |
120 | fHCuts->GetXaxis()->SetBinLabel(2,Form("%s: %.0f-%.0f",fCentVar.Data(),fCentFrom,fCentTo)); | |
121 | fHCuts->GetXaxis()->SetBinLabel(3,"QFlag"); | |
122 | fHCuts->GetXaxis()->SetBinLabel(4,Form("zvtx: %.0f-%.0f",fVtxZMin,fVtxZMax)); | |
123 | fOutput->Add(fHCuts); | |
124 | fHVertexZ = new TH1F("hVertexZBeforeCut","",100,-25,25); | |
125 | fHVertexZ->SetXTitle("z [cm]"); | |
126 | fOutput->Add(fHVertexZ); | |
76332037 | 127 | fHVertexZ2 = new TH1F("hVertexZAfterCut","",100,-25,25); |
128 | fHVertexZ2->SetXTitle("z [cm]"); | |
129 | fOutput->Add(fHVertexZ2); | |
fa443410 | 130 | fHCent = new TH1F("hCentBeforeCut","",101,-1,100); |
131 | fHCent->SetXTitle(fCentVar.Data()); | |
76332037 | 132 | fOutput->Add(fHCent); |
fa443410 | 133 | fHCentQual = new TH1F("hCentAfterCut","",101,-1,100); |
134 | fHCentQual->SetXTitle(fCentVar.Data()); | |
135 | fOutput->Add(fHCentQual); | |
90d5b88b | 136 | |
d595acbb | 137 | // histograms for cells |
138 | fHColuRow = new TH2F*[4]; | |
139 | fHColuRowE = new TH2F*[4]; | |
140 | fHCellMult = new TH1F*[4]; | |
141 | for (Int_t i = 0; i<4; ++i) { | |
142 | fHColuRow[i] = new TH2F(Form("hColRow_Mod%d",i),"",49,0,49,25,0,25); | |
143 | fHColuRow[i]->SetTitle(Form("Module %d: Occupancy", i)); | |
144 | fHColuRow[i]->SetXTitle("col (i#eta)"); | |
145 | fHColuRow[i]->SetYTitle("row (i#phi)"); | |
146 | fHColuRowE[i] = new TH2F(Form("hColRowE_Mod%d", i),"",49,0,49,25,0,25); | |
90d5b88b | 147 | fHColuRowE[i]->SetTitle(Form("Module %d: Cell energy",i)); |
d595acbb | 148 | fHColuRowE[i]->SetXTitle("col (i#eta)"); |
149 | fHColuRowE[i]->SetYTitle("row (i#phi)"); | |
150 | fHCellMult[i] = new TH1F(Form("hCellMult_Mod%d",i),"",1000,0,1000); | |
151 | fHCellMult[i]->SetTitle(Form("Module %d: Cell multiplicity",i)); | |
90d5b88b | 152 | fHCellMult[i]->SetXTitle("# of cells"); |
d595acbb | 153 | fOutput->Add(fHColuRow[i]); |
154 | fOutput->Add(fHColuRowE[i]); | |
155 | fOutput->Add(fHCellMult[i]); | |
156 | } | |
a49742b5 | 157 | fHCellE = new TH1F("hCellE","",150,0.,15.); |
d595acbb | 158 | fHCellE->SetXTitle("E_{cell} [GeV]"); |
159 | fOutput->Add(fHCellE); | |
a49742b5 | 160 | fHCellH = new TH1F ("fHCellHighestE","",150,0.,15.); |
6eb6260e | 161 | fHCellH->SetXTitle("E^{max}_{cell} [GeV]"); |
d595acbb | 162 | fOutput->Add(fHCellH); |
a49742b5 | 163 | fHCellM = new TH1F ("fHCellMeanEperHitCell","",250,0.,2.5); |
6eb6260e | 164 | fHCellM->SetXTitle("#LT E_{cell}#GT [GeV]"); |
165 | fOutput->Add(fHCellM); | |
a49742b5 | 166 | fHCellM2 = new TH1F ("fHCellMeanEperAllCells","",250,0.,1); |
f4ec884e | 167 | fHCellM2->SetXTitle("1/N_{cells} #Sigma E_{cell} [GeV]"); |
a49742b5 | 168 | fOutput->Add(fHCellM2); |
90d5b88b | 169 | |
d595acbb | 170 | // histograms for clusters |
6eb6260e | 171 | fHClustEccentricity = new TH1F("hClustEccentricity","",100,-0.1,1.1); |
172 | fHClustEccentricity->SetXTitle("#epsilon_{C}"); | |
173 | fOutput->Add(fHClustEccentricity); | |
a49742b5 | 174 | fHClustEtaPhi = new TH2F("hClustEtaPhi","",500,-0.8,0.8,500,1.2,2.2); |
fa443410 | 175 | fHClustEtaPhi->SetXTitle("#eta"); |
176 | fHClustEtaPhi->SetYTitle("#varphi"); | |
177 | fOutput->Add(fHClustEtaPhi); | |
178 | fHClustEnergyPt = new TH2F("hClustEnergyPt","",250,0,50,250,0,50); | |
179 | fHClustEnergyPt->SetXTitle("E [GeV]"); | |
6eb6260e | 180 | fHClustEnergyPt->SetYTitle("p_{T} [GeV/c]"); |
fa443410 | 181 | fOutput->Add(fHClustEnergyPt); |
a49742b5 | 182 | fHClustEnergySigma = new TH2F("hClustEnergySigma","",100,0,100,500,0,50); |
183 | fHClustEnergySigma->SetXTitle("E_{C} * #sigma_{max} [GeV*cm]"); | |
6eb6260e | 184 | fHClustEnergySigma->SetYTitle("E_{C} [GeV]"); |
d595acbb | 185 | fOutput->Add(fHClustEnergySigma); |
a49742b5 | 186 | fHClustSigmaSigma = new TH2F("hClustSigmaSigma","",500,0,50,500,0,50); |
6eb6260e | 187 | fHClustSigmaSigma->SetXTitle("#lambda_{0} [cm]"); |
188 | fHClustSigmaSigma->SetYTitle("#sigma_{max} [cm]"); | |
d595acbb | 189 | fOutput->Add(fHClustSigmaSigma); |
6eb6260e | 190 | fHClustNCellEnergyRatio = new TH2F("hClustNCellEnergyRatio","",27,-0.5,26.5,101,-0.05,1.05); |
191 | fHClustNCellEnergyRatio->SetXTitle("N_{cells}"); | |
192 | fHClustNCellEnergyRatio->SetYTitle("E^{max}_{cell}/E_{clus}"); | |
193 | fOutput->Add(fHClustNCellEnergyRatio); | |
90d5b88b | 194 | |
d595acbb | 195 | // histograms for pion candidates |
fa443410 | 196 | fHPionEtaPhi = new TH2F("hPionEtaPhi","",100,-0.8,0.8,100,1.2,2.2); |
a49742b5 | 197 | fHPionEtaPhi->SetXTitle("#eta_{#gamma#gamma}"); |
198 | fHPionEtaPhi->SetYTitle("#varphi_{#gamma#gamma}"); | |
fa443410 | 199 | fOutput->Add(fHPionEtaPhi); |
a49742b5 | 200 | fHPionMggPt = new TH2F("hPionMggPt","",1000,0,2,100,0,20.0); |
fa443410 | 201 | fHPionMggPt->SetXTitle("M_{#gamma#gamma} [GeV/c^{2}]"); |
202 | fHPionMggPt->SetYTitle("p_{T}^{#gamma#gamma} [GeV/c]"); | |
203 | fOutput->Add(fHPionMggPt); | |
a49742b5 | 204 | fHPionMggAsym = new TH2F("hPionMggAsym","",1000,0,2,100,0,1); |
fa443410 | 205 | fHPionMggAsym->SetXTitle("M_{#gamma#gamma} [GeV/c^{2}]"); |
206 | fHPionMggAsym->SetYTitle("Z_{#gamma#gamma} [GeV]"); | |
207 | fOutput->Add(fHPionMggAsym); | |
a49742b5 | 208 | fHPionMggDgg = new TH2F("hPionMggDgg","",1000,0,2,100,0,10); |
6eb6260e | 209 | fHPionMggDgg->SetXTitle("M_{#gamma#gamma} [GeV/c^{2}]"); |
210 | fHPionMggDgg->SetYTitle("opening angle [grad]"); | |
211 | fOutput->Add(fHPionMggDgg); | |
78204d3d | 212 | const Int_t nbins = 20; |
213 | Double_t xbins[nbins] = {0.5,1,1.5,2,2.5,3,3.5,4,4.5,5,6,7,8,9,10,12.5,15,20,25,50}; | |
fa443410 | 214 | fPtRanges = new TAxis(nbins-1,xbins); |
215 | for (Int_t i = 0; i<=nbins; ++i) { | |
a49742b5 | 216 | fHPionInvMasses[i] = new TH1F(Form("hPionInvMass%d",i),"",1000,0,2); |
fa443410 | 217 | fHPionInvMasses[i]->SetXTitle("M_{#gamma#gamma} [GeV/c^{2}]"); |
218 | if (i==0) | |
219 | fHPionInvMasses[i]->SetTitle(Form("0 < p_{T}^{#gamma#gamma} <%.1f",xbins[0])); | |
220 | else if (i==nbins) | |
221 | fHPionInvMasses[i]->SetTitle(Form("p_{T}^{#gamma#gamma} > 50")); | |
222 | else | |
223 | fHPionInvMasses[i]->SetTitle(Form("%.1f < p_{T}^{#gamma#gamma} <%.1f",xbins[i-1],xbins[i])); | |
fa443410 | 224 | fOutput->Add(fHPionInvMasses[i]); |
225 | } | |
d595acbb | 226 | |
717fe7de | 227 | PostData(1, fOutput); |
ea3fd2d5 | 228 | } |
229 | ||
230 | //________________________________________________________________________ | |
231 | void AliAnalysisTaskEMCALPi0PbPb::UserExec(Option_t *) | |
232 | { | |
717fe7de | 233 | // Called for each event. |
234 | ||
235 | if (!InputEvent()) | |
236 | return; | |
237 | ||
43cfaa06 | 238 | AliAnalysisManager *am = AliAnalysisManager::GetAnalysisManager(); |
239 | fEsdEv = dynamic_cast<AliESDEvent*>(InputEvent()); | |
240 | if (fEsdEv) { | |
241 | am->LoadBranch("AliESDRun."); | |
242 | am->LoadBranch("AliESDHeader."); | |
243 | } else { | |
244 | fAodEv = dynamic_cast<AliAODEvent*>(InputEvent()); | |
245 | am->LoadBranch("header"); | |
246 | } | |
247 | ||
90d5b88b | 248 | if (!fGeom) { // set misalignment matrices (stored in first event) |
d595acbb | 249 | fGeom = new AliEMCALGeoUtils("EMCAL_FIRSTYEARV1","EMCAL"); |
250 | if (fEsdEv) { | |
251 | for (Int_t i=0; i<4; ++i) | |
252 | fGeom->SetMisalMatrix(fEsdEv->GetESDRun()->GetEMCALMatrix(i),i); | |
253 | } else { | |
254 | for (Int_t i=0; i<4; ++i) | |
255 | fGeom->SetMisalMatrix(fAodEv->GetHeader()->GetEMCALMatrix(i),i); | |
256 | } | |
257 | fGeom->GetEMCGeometry(); | |
258 | } | |
259 | ||
286b47a5 | 260 | Int_t cut = 1; |
fa443410 | 261 | fHCuts->Fill(cut++); |
286b47a5 | 262 | |
717fe7de | 263 | const AliCentrality *centP = InputEvent()->GetCentrality(); |
264 | Double_t cent = centP->GetCentralityPercentileUnchecked(fCentVar); | |
fa443410 | 265 | fHCent->Fill(cent); |
717fe7de | 266 | if (cent<fCentFrom||cent>fCentTo) |
286b47a5 | 267 | return; |
43cfaa06 | 268 | |
fa443410 | 269 | fHCuts->Fill(cut++); |
286b47a5 | 270 | |
43cfaa06 | 271 | if (fUseQualFlag) { |
272 | if (centP->GetQuality()>0) | |
273 | return; | |
274 | } | |
286b47a5 | 275 | |
fa443410 | 276 | fHCentQual->Fill(cent); |
277 | fHCuts->Fill(cut++); | |
717fe7de | 278 | |
d9f26424 | 279 | // count number of accepted events |
280 | ++fNEvs; | |
281 | ||
43cfaa06 | 282 | if (fEsdEv) { |
fa443410 | 283 | am->LoadBranch("PrimaryVertex."); |
284 | am->LoadBranch("SPDVertex."); | |
285 | am->LoadBranch("TPCVertex."); | |
43cfaa06 | 286 | } else { |
287 | fAodEv = dynamic_cast<AliAODEvent*>(InputEvent()); | |
288 | am->LoadBranch("vertices"); | |
289 | } | |
290 | ||
291 | const AliVVertex *vertex = InputEvent()->GetPrimaryVertex(); | |
717fe7de | 292 | if (!vertex) |
293 | return; | |
294 | ||
fa443410 | 295 | fHVertexZ->Fill(vertex->GetZ()); |
286b47a5 | 296 | |
717fe7de | 297 | if(vertex->GetZ()<fVtxZMin||vertex->GetZ()>fVtxZMax) |
298 | return; | |
286b47a5 | 299 | |
fa443410 | 300 | fHCuts->Fill(cut++); |
76332037 | 301 | fHVertexZ2->Fill(vertex->GetZ()); |
717fe7de | 302 | |
43cfaa06 | 303 | fRecPoints = 0; // will be set if fClusName is given and AliAnalysisTaskEMCALClusterizeFast is used |
304 | fEsdClusters = 0; // will be set if ESD input used and if fRecPoints are not set of if clusters are attached | |
305 | fEsdCells = 0; // will be set if ESD input used | |
306 | fAodClusters = 0; // will be set if AOD input used and if fRecPoints are not set of if clusters are attached | |
307 | // or if fClusName is given and AliAnalysisTaskEMCALClusterizeFast in AOD output mode | |
308 | fAodCells = 0; // will be set if AOD input used | |
717fe7de | 309 | |
43cfaa06 | 310 | // deal with special output from AliAnalysisTaskEMCALClusterizeFast first |
311 | Bool_t clusattached = 0; | |
312 | Bool_t recalibrated = 0; | |
313 | if (1 && !fClusName.IsNull()) { | |
314 | AliAnalysisTaskEMCALClusterizeFast *cltask = 0; | |
315 | TObjArray *ts = am->GetTasks(); | |
316 | cltask = dynamic_cast<AliAnalysisTaskEMCALClusterizeFast*>(ts->FindObject(fClusName)); | |
317 | if (cltask && cltask->GetClusters()) { | |
d595acbb | 318 | fRecPoints = const_cast<TObjArray*>(cltask->GetClusters()); |
43cfaa06 | 319 | clusattached = cltask->GetAttachClusters(); |
320 | if (cltask->GetCalibData()!=0) | |
321 | recalibrated = kTRUE; | |
322 | } | |
323 | } | |
324 | if (1 && AODEvent() && !fClusName.IsNull()) { | |
717fe7de | 325 | TList *l = AODEvent()->GetList(); |
326 | TClonesArray *clus = 0; | |
327 | if (l) { | |
328 | clus = dynamic_cast<TClonesArray*>(l->FindObject(fClusName)); | |
329 | fAodClusters = clus; | |
330 | } | |
ea3fd2d5 | 331 | } |
717fe7de | 332 | |
333 | if (fEsdEv) { // ESD input mode | |
43cfaa06 | 334 | if (1 && (!fRecPoints||clusattached)) { |
335 | if (!clusattached) | |
336 | am->LoadBranch("CaloClusters"); | |
717fe7de | 337 | TList *l = fEsdEv->GetList(); |
338 | TClonesArray *clus = 0; | |
339 | if (l) { | |
340 | clus = dynamic_cast<TClonesArray*>(l->FindObject("CaloClusters")); | |
341 | fEsdClusters = clus; | |
ea3fd2d5 | 342 | } |
343 | } | |
43cfaa06 | 344 | if (1) { |
345 | if (!recalibrated) | |
346 | am->LoadBranch("EMCALCells."); | |
717fe7de | 347 | fEsdCells = fEsdEv->GetEMCALCells(); |
348 | } | |
349 | } else if (fAodEv) { // AOD input mode | |
43cfaa06 | 350 | if (1 && (!fAodClusters || clusattached)) { |
351 | if (!clusattached) | |
352 | am->LoadBranch("caloClusters"); | |
717fe7de | 353 | TList *l = fAodEv->GetList(); |
354 | TClonesArray *clus = 0; | |
355 | if (l) { | |
356 | clus = dynamic_cast<TClonesArray*>(l->FindObject("caloClusters")); | |
357 | fAodClusters = clus; | |
ea3fd2d5 | 358 | } |
359 | } | |
43cfaa06 | 360 | if (1) { |
361 | if (!recalibrated) | |
362 | am->LoadBranch("emcalCells"); | |
717fe7de | 363 | fAodCells = fAodEv->GetEMCALCells(); |
364 | } | |
365 | } else { | |
366 | AliFatal("Impossible to not have either pointer to ESD or AOD event"); | |
367 | } | |
ea3fd2d5 | 368 | |
43cfaa06 | 369 | if (1) { |
370 | AliDebug(2,Form("fRecPoints set: %p", fRecPoints)); | |
371 | AliDebug(2,Form("fEsdClusters set: %p", fEsdClusters)); | |
372 | AliDebug(2,Form("fEsdCells set: %p", fEsdCells)); | |
373 | AliDebug(2,Form("fAodClusters set: %p", fAodClusters)); | |
374 | AliDebug(2,Form("fAodCells set: %p", fAodCells)); | |
375 | } | |
376 | ||
a49742b5 | 377 | if (fDoAfterburner) |
378 | ClusterAfterburner(); | |
6eb6260e | 379 | |
286b47a5 | 380 | FillCellHists(); |
381 | FillClusHists(); | |
382 | FillPionHists(); | |
ea3fd2d5 | 383 | |
717fe7de | 384 | PostData(1, fOutput); |
ea3fd2d5 | 385 | } |
386 | ||
387 | //________________________________________________________________________ | |
388 | void AliAnalysisTaskEMCALPi0PbPb::Terminate(Option_t *) | |
389 | { | |
717fe7de | 390 | // Terminate called at the end of analysis. |
f5d4ab70 | 391 | |
392 | if (fNtuple) { | |
393 | TFile *f = OpenFile(1); | |
394 | if (f) | |
395 | fNtuple->Write(); | |
396 | } | |
d9f26424 | 397 | |
f224d35b | 398 | AliInfo(Form("\n%s: Accepted %lld events", GetName(), fNEvs)); |
ea3fd2d5 | 399 | } |
ea3fd2d5 | 400 | |
286b47a5 | 401 | //________________________________________________________________________ |
402 | void AliAnalysisTaskEMCALPi0PbPb::FillCellHists() | |
403 | { | |
404 | // Fill histograms related to cell properties. | |
d595acbb | 405 | |
90d5b88b | 406 | AliVCaloCells *cells = fEsdCells; |
407 | if (!cells) | |
408 | cells = fAodCells; | |
409 | ||
410 | if (cells) { | |
411 | Int_t cellModCount[4] = {0,0,0,0}; | |
412 | Double_t cellMaxE = 0; | |
6eb6260e | 413 | Double_t cellMeanE = 0; |
90d5b88b | 414 | Int_t ncells = cells->GetNumberOfCells(); |
6eb6260e | 415 | if (ncells==0) |
416 | return; | |
90d5b88b | 417 | |
d595acbb | 418 | for (Int_t i = 0; i<ncells; ++i ) { |
27422847 | 419 | Int_t absID = TMath::Abs(cells->GetCellNumber(i)); |
90d5b88b | 420 | Double_t cellE = cells->GetAmplitude(i); |
d595acbb | 421 | fHCellE->Fill(cellE); |
90d5b88b | 422 | if (cellE>cellMaxE) |
423 | cellMaxE = cellE; | |
6eb6260e | 424 | cellMeanE += cellE; |
d595acbb | 425 | |
90d5b88b | 426 | Int_t iSM=-1, iTower=-1, nIphi=-1, nIeta=-1; |
427 | Bool_t ret = fGeom->GetCellIndex(absID, iSM, iTower, nIphi, nIeta); | |
428 | if (!ret) { | |
429 | AliError(Form("Could not get cell index for %d", absID)); | |
430 | continue; | |
431 | } | |
432 | ++cellModCount[iSM]; | |
433 | Int_t iPhi=-1, iEta=-1; | |
434 | fGeom->GetCellPhiEtaIndexInSModule(iSM, iTower, nIphi, nIeta, iPhi, iEta); | |
435 | fHColuRow[iSM]->Fill(iEta,iPhi,1); | |
436 | fHColuRowE[iSM]->Fill(iEta,iPhi,cellE); | |
437 | } | |
438 | fHCellH->Fill(cellMaxE); | |
6eb6260e | 439 | cellMeanE /= ncells; |
440 | fHCellM->Fill(cellMeanE); | |
a49742b5 | 441 | fHCellM2->Fill(cellMeanE*ncells/24/48/4); //hard-coded but there is a way to figure out from geometry |
90d5b88b | 442 | for (Int_t i=0; i<4; ++i) |
443 | fHCellMult[i]->Fill(cellModCount[i]); | |
444 | } | |
286b47a5 | 445 | } |
446 | ||
447 | //________________________________________________________________________ | |
448 | void AliAnalysisTaskEMCALPi0PbPb::FillClusHists() | |
449 | { | |
90d5b88b | 450 | // Fill histograms related to cluster properties. |
6eb6260e | 451 | Double_t clusterEcc = 0; |
452 | ||
fa443410 | 453 | Double_t vertex[3] = {0,0,0}; |
454 | InputEvent()->GetPrimaryVertex()->GetXYZ(vertex); | |
455 | ||
90d5b88b | 456 | TObjArray *clusters = fEsdClusters; |
457 | if (!clusters) | |
458 | clusters = fAodClusters; | |
459 | ||
460 | if (clusters) { | |
461 | TLorentzVector clusterVec; | |
462 | Int_t nclus = clusters->GetEntries(); | |
fa443410 | 463 | for(Int_t i = 0; i<nclus; ++i) { |
90d5b88b | 464 | AliVCluster *clus = static_cast<AliVCluster*>(clusters->At(i)); |
fa443410 | 465 | if (!clus) |
466 | continue; | |
467 | if (!clus->IsEMCAL()) | |
468 | continue; | |
469 | clus->GetMomentum(clusterVec,vertex); | |
6eb6260e | 470 | Double_t maxAxis = 0, minAxis = 0; |
471 | GetSigma(clus,maxAxis,minAxis); | |
472 | if (maxAxis > 0) | |
473 | clusterEcc = TMath::Sqrt(1.0 - minAxis*minAxis/(maxAxis*maxAxis)); | |
f224d35b | 474 | clus->SetChi2(clusterEcc); // store ecc in chi2 |
6eb6260e | 475 | fHClustEccentricity->Fill(clusterEcc); |
fa443410 | 476 | fHClustEtaPhi->Fill(clusterVec.Eta(),clusterVec.Phi()); |
477 | fHClustEnergyPt->Fill(clusterVec.E(),clusterVec.Pt()); | |
6eb6260e | 478 | fHClustEnergySigma->Fill(clus->E()*maxAxis,clus->E()); |
479 | fHClustSigmaSigma->Fill(max(clus->GetM02(),clus->GetM20()),clus->E()*maxAxis); | |
480 | fHClustNCellEnergyRatio->Fill(clus->GetNCells(),GetMaxCellEnergy(clus)/clus->E()); | |
a49742b5 | 481 | if (fNtuple) { |
f4ec884e | 482 | Float_t vals[17]; |
a49742b5 | 483 | vals[0] = InputEvent()->GetRunNumber(); |
881a35d9 | 484 | vals[1] = (((UInt_t)InputEvent()->GetOrbitNumber() << 12) | (UInt_t)InputEvent()->GetBunchCrossNumber()); |
d9f26424 | 485 | if (vals[1]<1) |
486 | vals[1] = fNEvs; | |
a49742b5 | 487 | vals[2] = InputEvent()->GetCentrality()->GetCentralityPercentileUnchecked(fCentVar); |
488 | vals[3] = clusterVec.Pt(); | |
489 | vals[4] = clusterVec.E(); | |
490 | vals[5] = GetMaxCellEnergy(clus); | |
491 | vals[6] = clus->GetNCells(); | |
492 | vals[7] = clus->GetDistanceToBadChannel(); | |
493 | vals[8] = clus->GetDispersion(); | |
494 | vals[9] = clus->GetM20(); | |
495 | vals[10] = clus->GetM02(); | |
496 | vals[11] = clus->Chi2(); | |
497 | vals[12] = clus->GetEmcCpvDistance(); | |
498 | vals[13] = clusterEcc; | |
f4ec884e | 499 | vals[14] = maxAxis; |
500 | vals[15] = clusterVec.Eta(); | |
501 | vals[16] = clusterVec.Phi(); | |
a49742b5 | 502 | fNtuple->Fill(vals); |
503 | } | |
fa443410 | 504 | } |
505 | } | |
286b47a5 | 506 | } |
507 | ||
508 | //________________________________________________________________________ | |
509 | void AliAnalysisTaskEMCALPi0PbPb::FillPionHists() | |
510 | { | |
511 | // Fill histograms related to pions. | |
286b47a5 | 512 | |
fa443410 | 513 | Double_t vertex[3] = {0,0,0}; |
514 | InputEvent()->GetPrimaryVertex()->GetXYZ(vertex); | |
515 | ||
90d5b88b | 516 | TObjArray *clusters = fEsdClusters; |
517 | if (!clusters) | |
518 | clusters = fAodClusters; | |
fa443410 | 519 | |
90d5b88b | 520 | if (clusters) { |
521 | TLorentzVector clusterVec1; | |
522 | TLorentzVector clusterVec2; | |
523 | TLorentzVector pionVec; | |
524 | ||
525 | Int_t nclus = clusters->GetEntries(); | |
fa443410 | 526 | for (Int_t i = 0; i<nclus; ++i) { |
90d5b88b | 527 | AliVCluster *clus1 = static_cast<AliVCluster*>(clusters->At(i)); |
fa443410 | 528 | if (!clus1) |
529 | continue; | |
530 | if (!clus1->IsEMCAL()) | |
531 | continue; | |
6eb6260e | 532 | if (clus1->E()<0.010) |
533 | continue; | |
a49742b5 | 534 | if (clus1->GetNCells()<fNminCells) |
535 | continue; | |
f224d35b | 536 | if (GetMaxCellEnergy(clus1)/clus1->E()<fMinErat) |
537 | continue; | |
538 | if (clus1->Chi2()<fMinEcc) | |
539 | continue; | |
fa443410 | 540 | clus1->GetMomentum(clusterVec1,vertex); |
541 | for (Int_t j = i+1; j<nclus; ++j) { | |
90d5b88b | 542 | AliVCluster *clus2 = static_cast<AliVCluster*>(clusters->At(j)); |
fa443410 | 543 | if (!clus2) |
544 | continue; | |
545 | if (!clus2->IsEMCAL()) | |
546 | continue; | |
6eb6260e | 547 | if (clus2->E()<0.010) |
548 | continue; | |
a49742b5 | 549 | if (clus2->GetNCells()<fNminCells) |
550 | continue; | |
f224d35b | 551 | if (GetMaxCellEnergy(clus2)/clus2->E()<fMinErat) |
552 | continue; | |
553 | if (clus2->Chi2()<fMinEcc) | |
554 | continue; | |
fa443410 | 555 | clus2->GetMomentum(clusterVec2,vertex); |
556 | pionVec = clusterVec1 + clusterVec2; | |
557 | Double_t pionZgg = TMath::Abs(clusterVec1.E()-clusterVec2.E())/pionVec.E(); | |
6eb6260e | 558 | Double_t pionOpeningAngle = clusterVec1.Angle(clusterVec2.Vect()); |
d595acbb | 559 | if (pionZgg < fAsymMax) { |
560 | fHPionEtaPhi->Fill(pionVec.Eta(),pionVec.Phi()); | |
561 | fHPionMggPt->Fill(pionVec.M(),pionVec.Pt()); | |
562 | fHPionMggAsym->Fill(pionVec.M(),pionZgg); | |
6eb6260e | 563 | fHPionMggDgg->Fill(pionVec.M(),pionOpeningAngle); |
d595acbb | 564 | Int_t bin = fPtRanges->FindBin(pionVec.Pt()); |
565 | fHPionInvMasses[bin]->Fill(pionVec.M()); | |
566 | } | |
fa443410 | 567 | } |
568 | } | |
90d5b88b | 569 | } |
fa443410 | 570 | } |
d595acbb | 571 | |
6eb6260e | 572 | //________________________________________________________________________ |
573 | void AliAnalysisTaskEMCALPi0PbPb::ClusterAfterburner() | |
574 | { | |
575 | // | |
576 | ||
577 | AliVCaloCells *cells = fEsdCells; | |
578 | if (!cells) | |
579 | cells = fAodCells; | |
580 | ||
581 | if (!cells) | |
582 | return; | |
583 | ||
584 | Int_t ncells = cells->GetNumberOfCells(); | |
585 | if (ncells<=0) | |
586 | return; | |
587 | ||
588 | Double_t cellMeanE = 0, cellSigE = 0; | |
589 | for (Int_t i = 0; i<ncells; ++i ) { | |
590 | Double_t cellE = cells->GetAmplitude(i); | |
591 | cellMeanE += cellE; | |
592 | cellSigE += cellE*cellE; | |
593 | } | |
594 | cellMeanE /= ncells; | |
595 | cellSigE /= ncells; | |
596 | cellSigE -= cellMeanE*cellMeanE; | |
597 | if (cellSigE<0) | |
598 | cellSigE = 0; | |
599 | cellSigE = TMath::Sqrt(cellSigE / ncells); | |
600 | ||
601 | Double_t subE = cellMeanE - 7*cellSigE; | |
602 | if (subE<0) | |
603 | return; | |
604 | ||
605 | for (Int_t i = 0; i<ncells; ++i) { | |
606 | Short_t id=-1; | |
607 | Double_t amp=0,time=0; | |
608 | if (!cells->GetCell(i, id, amp, time)) | |
609 | continue; | |
610 | amp -= cellMeanE; | |
611 | if (amp<0.001) | |
612 | amp = 0; | |
613 | cells->SetCell(i, id, amp, time); | |
614 | } | |
615 | ||
616 | TObjArray *clusters = fEsdClusters; | |
617 | if (!clusters) | |
618 | clusters = fAodClusters; | |
619 | ||
620 | if (!clusters) | |
621 | return; | |
622 | ||
623 | Int_t nclus = clusters->GetEntries(); | |
624 | for (Int_t i = 0; i<nclus; ++i) { | |
625 | AliVCluster *clus = static_cast<AliVCluster*>(clusters->At(i)); | |
626 | Int_t nc = clus->GetNCells(); | |
6eb6260e | 627 | Double_t clusE = 0; |
628 | Short_t nc2 = 0; | |
629 | UShort_t ids[100] = {0}; | |
630 | Double_t fra[100] = {0}; | |
631 | for (Int_t j = 0; j<nc; ++j) { | |
27422847 | 632 | Short_t id = TMath::Abs(clus->GetCellAbsId(j)); |
6eb6260e | 633 | Double_t cen = cells->GetCellAmplitude(id); |
634 | clusE += cen; | |
635 | if (cen>0) { | |
636 | ids[nc] = id; | |
637 | ++nc; | |
638 | } | |
639 | } | |
640 | if (clusE<=0) { | |
641 | clusters->RemoveAt(i); | |
642 | continue; | |
643 | } | |
644 | for (Int_t j = 0; j<nc2; ++j) { | |
645 | Short_t id = ids[j]; | |
646 | Double_t cen = cells->GetCellAmplitude(id); | |
647 | fra[j] = cen/clusE; | |
648 | } | |
649 | clus->SetE(clusE); | |
650 | AliAODCaloCluster *aodclus = dynamic_cast<AliAODCaloCluster*>(clus); | |
651 | if (aodclus) { | |
652 | aodclus->Clear(""); | |
653 | aodclus->SetNCells(nc2); | |
654 | aodclus->SetCellsAmplitudeFraction(fra); | |
655 | aodclus->SetCellsAbsId(ids); | |
656 | } | |
6eb6260e | 657 | } |
6eb6260e | 658 | clusters->Compress(); |
659 | } | |
660 | ||
d595acbb | 661 | //________________________________________________________________________ |
90d5b88b | 662 | Double_t AliAnalysisTaskEMCALPi0PbPb::GetMaxCellEnergy(AliVCluster *cluster) |
d595acbb | 663 | { |
90d5b88b | 664 | // Get maximum energy of attached cell. |
665 | ||
666 | Double_t maxe = 0; | |
90d5b88b | 667 | Int_t ncells = cluster->GetNCells(); |
668 | if (fEsdCells) { | |
669 | for (Int_t i=0; i<ncells; i++) { | |
27422847 | 670 | Double_t e = fEsdCells->GetCellAmplitude(TMath::Abs(cluster->GetCellAbsId(i))); |
90d5b88b | 671 | if (e>maxe) |
672 | maxe = e; | |
673 | } | |
674 | } else { | |
675 | for (Int_t i=0; i<ncells; i++) { | |
27422847 | 676 | Double_t e = fAodCells->GetCellAmplitude(TMath::Abs(cluster->GetCellAbsId(i))); |
90d5b88b | 677 | if (e>maxe) |
678 | maxe = e; | |
679 | } | |
680 | } | |
6eb6260e | 681 | return maxe; |
d595acbb | 682 | } |
683 | ||
684 | //________________________________________________________________________ | |
6eb6260e | 685 | void AliAnalysisTaskEMCALPi0PbPb::GetSigma(AliVCluster *cluster, Double_t& sigmaMax, Double_t &sigmaMin) |
d595acbb | 686 | { |
90d5b88b | 687 | // Calculate the (E) weighted variance along the longer (eigen) axis. |
688 | ||
6eb6260e | 689 | sigmaMax = 0; // cluster variance along its longer axis |
690 | sigmaMin = 0; // cluster variance along its shorter axis | |
90d5b88b | 691 | Double_t Ec = cluster->E(); // cluster energy |
692 | Double_t Xc = 0 ; // cluster first moment along X | |
693 | Double_t Yc = 0 ; // cluster first moment along Y | |
694 | Double_t Sxx = 0 ; // cluster second central moment along X (variance_X^2) | |
695 | Double_t Sxy = 0 ; // cluster second central moment along Y (variance_Y^2) | |
696 | Double_t Syy = 0 ; // cluster covariance^2 | |
697 | ||
698 | AliVCaloCells *cells = fEsdCells; | |
699 | if (!cells) | |
700 | cells = fAodCells; | |
701 | ||
6eb6260e | 702 | if (!cells) |
703 | return; | |
704 | ||
705 | Int_t ncells = cluster->GetNCells(); | |
706 | if (ncells==1) | |
707 | return; | |
708 | ||
709 | TVector3 pos; | |
710 | for(Int_t j=0; j<ncells; ++j) { | |
27422847 | 711 | Int_t id = TMath::Abs(cluster->GetCellAbsId(j)); |
712 | fGeom->GetGlobal(id,pos); | |
6eb6260e | 713 | Double_t cellen = cells->GetCellAmplitude(id); |
714 | Xc += cellen*pos.X(); | |
715 | Yc += cellen*pos.Y(); | |
716 | Sxx += cellen*pos.X()*pos.X(); | |
717 | Syy += cellen*pos.Y()*pos.Y(); | |
718 | Sxy += cellen*pos.X()*pos.Y(); | |
719 | } | |
720 | Xc /= Ec; | |
721 | Yc /= Ec; | |
722 | Sxx /= Ec; | |
723 | Syy /= Ec; | |
724 | Sxy /= Ec; | |
725 | Sxx -= Xc*Xc; | |
726 | Syy -= Yc*Yc; | |
727 | Sxy -= Xc*Yc; | |
728 | sigmaMax = (Sxx + Syy + TMath::Sqrt((Sxx-Syy)*(Sxx-Syy)+4.0*Sxy*Sxy))/2.0; | |
729 | sigmaMax = TMath::Sqrt(sigmaMax); | |
730 | sigmaMin = TMath::Abs(Sxx + Syy - TMath::Sqrt((Sxx-Syy)*(Sxx-Syy)+4.0*Sxy*Sxy))/2.0; | |
731 | sigmaMin = TMath::Sqrt(sigmaMin); | |
d595acbb | 732 | } |
90d5b88b | 733 |