]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Setter for K
authormorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 25 Jan 2010 07:53:21 +0000 (07:53 +0000)
committermorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 25 Jan 2010 07:53:21 +0000 (07:53 +0000)
Histos for responsibilities.

JETAN/AliAnalysisTaskKMeans.cxx
JETAN/AliAnalysisTaskKMeans.h

index 24162f9b17b0b54766de46947abe46850f43c617..9ccf7d6752b0e901434fc32dba61df324f047a1e 100644 (file)
@@ -61,7 +61,8 @@
 ClassImp(AliAnalysisTaskKMeans)
 
 AliAnalysisTaskKMeans::AliAnalysisTaskKMeans() 
-    : AliAnalysisTaskSE() 
+    : AliAnalysisTaskSE()
+    ,fK(0)
     ,fHists(0)
     ,fH1CEta(0)
     ,fH1CPhi(0)
@@ -82,6 +83,8 @@ AliAnalysisTaskKMeans::AliAnalysisTaskKMeans()
     ,fH2DPhiEtaLR(0)
     ,fH2DPhiEtaC(0)
     ,fH2DPhiEtaCR(0)
+    ,fH1Resp(0)
+    ,fH1RespR(0)
     ,fCuts(0)
 {
   //
@@ -92,6 +95,7 @@ AliAnalysisTaskKMeans::AliAnalysisTaskKMeans()
 //________________________________________________________________________
 AliAnalysisTaskKMeans::AliAnalysisTaskKMeans(const char *name) 
     : AliAnalysisTaskSE(name) 
+      ,fK(0)
       ,fHists(0)
       ,fH1CEta(0)
       ,fH1CPhi(0)
@@ -112,6 +116,8 @@ AliAnalysisTaskKMeans::AliAnalysisTaskKMeans(const char *name)
       ,fH2DPhiEtaLR(0)
       ,fH2DPhiEtaC(0)
       ,fH2DPhiEtaCR(0)
+      ,fH1Resp(0)
+      ,fH1RespR(0)
       ,fCuts(0)
 {
   //
@@ -149,7 +155,8 @@ void AliAnalysisTaskKMeans::UserCreateOutputObjects()
     fH2DPhiEtaC   = new TH2F("fH2DPhiEtaC","eta phi distribution", 31, 0., TMath::Pi(), 20, 0., 2.);
     fH2DPhiEtaCR  = new TH2F("fH2DPhiEtaCR","eta phi distribution", 31, 0., TMath::Pi(), 20, 0., 2.);
     fH1DRR        = new TH1F("fH1DRR",    "dR distribution", 50, 0., 5.);
-    
+    fH1Resp       = new TH1F("fH1Resp",   "Responsibility", 50, 0., 1.);
+    fH1RespR      = new TH1F("fH1RespR",  "Responsibility", 50, 0., 1.);
     fHists->SetOwner();
 
     fHists->Add(fH1CEta);
@@ -171,9 +178,10 @@ void AliAnalysisTaskKMeans::UserCreateOutputObjects()
     fHists->Add(fH2DPhiEtaC);
     fHists->Add(fH2DPhiEtaCR);
     fHists->Add(fH1DRR);
-    
+    fHists->Add(fH1RespR);
+    fHists->Add(fH1Resp);    
     //
-    AliKMeansClustering::SetBeta(20.);
+    AliKMeansClustering::SetBeta(4.);
     
 }
 
@@ -238,10 +246,11 @@ void AliAnalysisTaskKMeans::UserExec(Option_t *)
       return;
   }
   //
-  AliKMeansClustering::SoftKMeans(2, ic, phi, eta, mPhi, mEta, rk);
+  AliKMeansClustering::SoftKMeans(fK, ic, phi, eta, mPhi, mEta, rk);
   //
   // Sort
-  TMath::Sort(2, rk, ind);
+  TMath::Sort(fK, rk, ind);
+  fH1Resp->Fill(rk[ind[0]]/(rk[0]+rk[1]));
   //
   // Analyse
   //
@@ -298,10 +307,11 @@ void AliAnalysisTaskKMeans::UserExec(Option_t *)
 
   // Randomized phi
   //
-  AliKMeansClustering::SoftKMeans(2, ic, phiR, etaR, mPhi, mEta, rk);
+  AliKMeansClustering::SoftKMeans(fK, ic, phiR, etaR, mPhi, mEta, rk);
   //
   // Sort
-  TMath::Sort(2, rk, ind);
+  TMath::Sort(fK, rk, ind);
+  fH1RespR->Fill(rk[ind[0]]/(rk[0]+rk[1]));
   //
   // Analyse
   //
index 8e38474ecabe2389c00f7bde1a171cd5f25322ad..34380c338fa0355707d1bd85eceedba377a9431d 100644 (file)
@@ -37,9 +37,10 @@ class AliAnalysisTaskKMeans : public AliAnalysisTaskSE {
   virtual void     SetCuts(AliESDtrackCuts* cuts) {fCuts = cuts;}
   virtual Double_t DeltaPhi(Double_t phi1, Double_t phi2);
   virtual Double_t DeltaR(Double_t phi1, Double_t eta1, Double_t phi2, Double_t eta2);
-  
+  virtual void     SetK(Int_t k) {fK = k;} 
  private:
   // Others
+  Int_t            fK;             // K                        
   TList*           fHists;         // Histograms
   TH1F*            fH1CEta;        // Eta distribution of clusters
   TH1F*            fH1CPhi;        // Phi distribution of clusters  
@@ -60,7 +61,8 @@ class AliAnalysisTaskKMeans : public AliAnalysisTaskSE {
   TH2F*            fH2DPhiEtaLR;   // eta-phi of leading particle
   TH2F*            fH2DPhiEtaC;    // eta-phi of Clusters
   TH2F*            fH2DPhiEtaCR;   // eta-phi of Clusters
-  
+  TH1F*            fH1Resp;        // responsibility
+  TH1F*            fH1RespR;       // responsibility
   AliESDtrackCuts* fCuts;             // List of cuts
   ClassDef(AliAnalysisTaskKMeans, 1); // A k-means clustering analysis
 };