1 /**************************************************************************
2 * Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 //---------------------------------------------------------------------------------------
18 // Analysis Task that uses the Soft K-Means Algorithm to find clusters in
19 // the eta-phi space of Minimum Bias. No pt information is used for the clustering.
22 // Author: Andreas Morsch (CERN)
23 // andreas.morsch@cern.ch
24 //---------------------------------------------------------------------------------------
36 #include "TParticle.h"
37 #include "TParticlePDG.h"
41 #include "TObjArray.h"
43 #include "AliAnalysisTask.h"
44 #include "AliAnalysisManager.h"
46 #include "AliVEvent.h"
47 #include "AliESDEvent.h"
48 #include "AliExternalTrackParam.h"
50 #include "AliESDVertex.h"
51 #include "AliESDInputHandler.h"
52 #include "AliESDtrackCuts.h"
53 #include "AliMultiplicity.h"
55 #include "AliMCParticle.h"
56 #include "AliMCEvent.h"
57 #include "AliAnalysisTaskKMeans.h"
58 #include "AliTrackReference.h"
60 #include "AliHeader.h"
61 #include "AliKMeansClustering.h"
65 ClassImp(AliAnalysisTaskKMeans)
67 AliAnalysisTaskKMeans::AliAnalysisTaskKMeans()
112 //________________________________________________________________________
113 AliAnalysisTaskKMeans::AliAnalysisTaskKMeans(const char *name)
114 : AliAnalysisTaskSE(name)
156 DefineOutput(1, TList::Class());
159 AliAnalysisTaskKMeans::AliAnalysisTaskKMeans(const AliAnalysisTaskKMeans &res)
160 : AliAnalysisTaskSE(res)
198 // Dummy copy constructor
201 AliAnalysisTaskKMeans& AliAnalysisTaskKMeans::operator=(const AliAnalysisTaskKMeans& /*trclass*/)
203 // Dummy assignment operator
209 //________________________________________________________________________
210 void AliAnalysisTaskKMeans::UserCreateOutputObjects()
214 fHists = new TList();
215 fH1CEta = new TH1F("fH1CEta", "eta distribution of clusters", 90, -1.5, 1.5);
216 fH1CEtaR = new TH1F("fH1CEtaR", "eta distribution of clusters", 90, -1.5, 1.5);
217 fH1CPhi = new TH1F("fH1CPhi", "phi distribution of clusters", 157, 0.0, 2. * TMath::Pi());
218 fH2N1N2 = new TH2F("fH2N1N2", "multiplicity distribution", 50, 0., 50., 50, 0., 50.);
220 fH1Pt = new TH1F("fH1Pt", "pt distribution",50, 0., 10.);
221 fH1PtC = new TH1F("fH1PtC", "pt distribution",50, 0., 10.);
222 fH1PtC1 = new TH1F("fH1PtC1", "pt distribution",50, 0., 10.);
223 fH1PtC2 = new TH1F("fH1PtC2", "pt distribution",50, 0., 10.);
224 fH1PtAS = new TH1F("fH1PtAS", "pt distribution",50, 0., 10.);
225 fH1PtR = new TH1F("fH1PtR", "pt distribution",50, 0., 10.);
227 fH1SPt = new TH1F("fH1SPt", "sum pt distribution",50, 0., 10.);
228 fH1SPtC = new TH1F("fH1SPtC", "sum pt distribution",50, 0., 10.);
230 fH1DR = new TH1F("fH1DR", "dR distribution", 50, 0., 5.);
231 fH1DPhi = new TH1F("fH1DPhi", "dPhi distribution", 31, 0., TMath::Pi());
232 fH2DPhiEta = new TH2F("fH2DPhiEta","eta phi distribution", 31, 0., TMath::Pi(), 20, 0., 2.);
233 fH2DPhiEtaR = new TH2F("fH2DPhiEtaR","eta phi distribution", 31, 0., TMath::Pi(), 20, 0., 2.);
234 fH2DPhiEtaL = new TH2F("fH2DPhiEtaL","eta phi distribution", 31, 0., TMath::Pi(), 20, 0., 2.);
235 fH2DPhiEtaLR = new TH2F("fH2DPhiEtaLR","eta phi distribution", 31, 0., TMath::Pi(), 20, 0., 2.);
236 fH2DPhiEtaC = new TH2F("fH2DPhiEtaC","eta phi distribution", 31, 0., TMath::Pi(), 20, 0., 2.);
237 fH2DPhiEtaCR = new TH2F("fH2DPhiEtaCR","eta phi distribution", 31, 0., TMath::Pi(), 20, 0., 2.);
238 fH1DRR = new TH1F("fH1DRR", "dR distribution", 50, 0., 5.);
239 fH1Resp = new TH1F("fH1Resp", "Responsibility", 50, 0., 1.);
240 fH1RespR = new TH1F("fH1RespR", "Responsibility", 50, 0., 1.);
241 fH2Sigma = new TH2F("fH2Sigma", "Sigma", 31, 0., TMath::Pi(), 20, 0., 2.);
242 fH2SigmaR = new TH2F("fH2SigmaR", "Sigma", 31, 0., TMath::Pi(), 20, 0., 2.);
245 fHDensity = new TH1F("fHDensity", "density distribution", 100, 0., 20.);
246 fHCSize = new TH1F("fHCSize", "cluster size", 20, -0.5, 19.5);
248 fHNCluster = new TH1F("fHNCluster", "Number of clusters", 11, -0.5, 10.5);
249 fHPtDensity = new TH2F("fHPtDensity", "Pt vs density", 100, 0., 20., 50, 0., 10.);
250 fHDPhi = new TH1F("fHDPhi", "phi correlation", 100, 0., 0.5);
251 fH2EtaPhi = new TH2F("fH2EtaPhi", "Eta Phi", 200, -1., 1., 628, 0., 2. * TMath::Pi());
254 fHists->Add(fH1CEta);
255 fHists->Add(fH1CEtaR);
256 fHists->Add(fH1CPhi);
257 fHists->Add(fH2N1N2);
261 fHists->Add(fH1PtC1);
262 fHists->Add(fH1PtC2);
263 fHists->Add(fH1PtAS);
265 fHists->Add(fH1SPtC);
267 fHists->Add(fH1DPhi);
268 fHists->Add(fH2DPhiEta);
269 fHists->Add(fH2DPhiEtaR);
270 fHists->Add(fH2DPhiEtaL);
271 fHists->Add(fH2DPhiEtaLR);
272 fHists->Add(fH2DPhiEtaC);
273 fHists->Add(fH2DPhiEtaCR);
275 fHists->Add(fH1RespR);
276 fHists->Add(fH1Resp);
277 fHists->Add(fH2Sigma);
278 fHists->Add(fH2SigmaR);
279 fHists->Add(fHCSize);
280 fHists->Add(fHDensity);
281 fHists->Add(fHNCluster);
282 fHists->Add(fHPtDensity);
284 fHists->Add(fH2EtaPhi);
286 for (Int_t i = 0; i < 10; i++) {
287 fA[i] = new AliKMeansResult(i+1);
288 fB[i] = new AliKMeansResult(i+1);
292 //________________________________________________________________________
293 void AliAnalysisTaskKMeans::UserExec(Option_t *)
296 // Called for each event
304 Printf("ERROR: fESD not available");
312 // Fill eta-phi positions
318 for (Int_t iTracks = 0; iTracks < fInputEvent->GetNumberOfTracks(); iTracks++) {
319 AliVParticle* track = fInputEvent->GetTrack(iTracks);
320 if ((fCuts->AcceptTrack((AliESDtrack*)track)))
322 const AliExternalTrackParam * tpcT = ((AliESDtrack*) track)->GetTPCInnerParam();
324 if (TMath::Abs(tpcT->Eta()) > 0.9) continue;
326 phi [ic] = tpcT->Phi();
327 eta [ic] = tpcT->Eta();
328 pt [ic] = tpcT->Pt();
331 fH2REtaPhi->GetRandom2(etaR[ic], phiR[ic]);
333 phiR[ic] = 2. * TMath::Pi() * gRandom->Rndm();
334 etaR[ic] = 1.8 * gRandom->Rndm() - 0.9;
338 if (pt[ic] > ptmax) {
343 fH2EtaPhi->Fill(eta[ic], phi[ic]);
348 for (Int_t i = 0; i < ic; i++) {
349 for (Int_t j = i+1; j < ic; j++) {
350 Double_t dphi = TMath::Abs(phi[i] - phi[j]);
364 AliKMeansResult* res = 0;
365 AliKMeansResult best(10);
368 for (Int_t k = 0; k < 20; k++) {
371 for (Int_t i = 0; i < fK; i++) {
373 AliKMeansClustering::SoftKMeans2(i+1, ic, phi, eta, res->GetMx(),res->GetMy(), res->GetSigma2(), res->GetRk());
374 res->Sort(ic, phi, eta);
375 Int_t j = (res->GetInd())[0];
376 rk0[i] = (res->GetTarget())[j];
384 best.CopyResults(fA[imax]);
387 Double_t* mPhi = best.GetMx();
388 Double_t* mEta = best.GetMy();
389 Double_t* sigma2 = best.GetSigma2();
390 Int_t nk = best.GetK();
391 Int_t im = (best.GetInd())[0];
392 Double_t etaC = mEta[im];
394 fHDensity->Fill(rmaxG / TMath::Pi());
395 fHCSize->Fill(2.28 * rmaxG * sigma2[im]);
396 fHNCluster->Fill(Float_t(nk));
398 Double_t dphic, detac;
400 if (rmaxG > 0. && TMath::Abs(etaC) < 0.4) {
403 // Cluster Multiplicity
404 Int_t mult[2] = {0, 0};
407 dphic = DeltaPhi(mPhi[0], mPhi[1]);
408 detac = TMath::Abs(mEta[0] - mEta[1]);
409 fH2DPhiEtaC->Fill(dphic, detac);
412 // Random cluster position
413 Int_t ir = Int_t(Float_t(ic) * gRandom->Rndm());
415 Double_t crPhi = phi[ir];
416 Double_t crEta = eta[ir];
421 for (Int_t i = 0; i < 1; i++) {
422 fH1CEta->Fill(mEta[im]);
423 fH1CPhi->Fill(mPhi[im]);
424 for (Int_t j = 0; j < ic; j++) {
425 Double_t r = DeltaR(mPhi[im], mEta[im], phi[j], eta[j]);
426 Double_t dphi = DeltaPhi(mPhi[im], phi[j]);
427 Double_t deta = mEta[im] - eta[j];
428 Double_t rr = DeltaR(crPhi, crEta, phi[j], eta[j]);
432 fH2DPhiEta->Fill(dphi, TMath::Abs(deta));
433 if (j == icmax) fH2DPhiEtaL->Fill(dphi, TMath::Abs(deta));
436 fH1PtC2->Fill(pt[j]);
439 fH1PtC1->Fill(pt[j]);
440 fHPtDensity->Fill(rmaxG/TMath::Pi(), pt[j]);
451 if (r > 0.7 && dphi < (TMath::Pi() - 0.3)) {
455 if (r > 0.7 && r < 1.) {
459 if (dphi > (TMath::Pi() - 0.3)) {
460 fH1PtAS->Fill(pt[j], 1.);
465 fH2N1N2->Fill(Float_t(mult[0]), Float_t(mult[1]));
466 fH1SPt ->Fill(sumPt);
467 fH1SPtC->Fill(sumPtC);
473 for (Int_t k = 0; k < 20; k++) {
476 for (Int_t i = 0; i < fK; i++) {
478 AliKMeansClustering::SoftKMeans2(i+1, ic, phiR, etaR, res->GetMx(),res->GetMy(), res->GetSigma2(), res->GetRk());
479 res->Sort(ic, phiR, etaR);
480 Int_t j = (res->GetInd())[0];
481 rk0[i] = (res->GetTarget())[j];
489 best.CopyResults(fB[imax]);
495 sigma2 = best.GetSigma2();
497 im = (best.GetInd())[0];
501 // Cluster Multiplicity
502 if (rmaxG > 0. && TMath::Abs(etaC) < 0.4) {
503 for (Int_t i = 0; i < 1; i++) {
504 im = (best.GetInd())[i];
505 fH1CEtaR->Fill(mEta[im]);
508 for (Int_t i = 0; i < 1; i++) {
509 im = (best.GetInd())[i];
510 for (Int_t j = 0; j < ic; j++) {
511 Double_t dphi = DeltaPhi(mPhi[im], phiR[j]);
512 Double_t deta = mEta[im] - etaR[j];
513 Double_t r = DeltaR(mPhi[im], mEta[im], phiR[j], etaR[j]);
515 fH2DPhiEtaR->Fill(dphi, TMath::Abs(deta));
516 if (j == icmax) fH2DPhiEtaLR->Fill(dphi, TMath::Abs(deta));
520 dphic = DeltaPhi(mPhi[0], mPhi[1]);
521 detac = TMath::Abs(mEta[0] - mEta[1]);
522 fH2DPhiEtaCR->Fill(dphic, detac);
528 Double_t AliAnalysisTaskKMeans::DeltaPhi(Double_t phi1, Double_t phi2)
530 Double_t dphi = TMath::Abs(phi1 - phi2);
531 if (dphi > TMath::Pi()) dphi = 2. * TMath::Pi() - dphi;
535 Double_t AliAnalysisTaskKMeans::DeltaR(Double_t phi1, Double_t eta1, Double_t phi2, Double_t eta2)
537 Double_t dphi = DeltaPhi(phi1, phi2);
538 Double_t deta = eta1 - eta2;
539 return (TMath::Sqrt(dphi * dphi + deta * deta));
543 //________________________________________________________________________
544 void AliAnalysisTaskKMeans::Terminate(Option_t *)