]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
o) Added DrawHistograms function
authorjgrosseo <jgrosseo@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 18 May 2006 10:08:32 +0000 (10:08 +0000)
committerjgrosseo <jgrosseo@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 18 May 2006 10:08:32 +0000 (10:08 +0000)
o) Added histogram that stores uncorrected data
o) Added Getters, Setters for the histograms

PWG0/dNdEta/dNdEtaAnalysis.cxx
PWG0/dNdEta/dNdEtaAnalysis.h
PWG0/dNdEta/testAnalysis.C

index 06501f1e34abe180e6215f62045790d50028f41d..076baa2e5a15645d31f6042f2fad3c9adbe0d802 100644 (file)
@@ -1,5 +1,11 @@
 #include "dNdEtaAnalysis.h"
 
+#include <TFile.h>
+#include <TH2F.h>
+#include <TH1D.h>
+#include <TMath.h>
+#include <TCanvas.h>
+
 //____________________________________________________________________
 ClassImp(dNdEtaAnalysis)
 
@@ -8,9 +14,10 @@ dNdEtaAnalysis::dNdEtaAnalysis(Char_t* name) {
 
   fName = TString(name);
 
-  hEtaVsVtx  = new TH2F("eta_vs_vtx","",80,-20,20,120,-6,6);   
-  hVtx       = (TH1F*)hEtaVsVtx->ProjectionX("vtx");
-  hdNdEta    = (TH1F*)hEtaVsVtx->ProjectionY("dNdEta");
+  hEtaVsVtx  = new TH2F("eta_vs_vtx","",80,-20,20,120,-6,6);
+  hEtaVsVtxUncorrected = new TH2F("eta_vs_vtx_uncorrected","",80,-20,20,120,-6,6);
+  hVtx       = hEtaVsVtx->ProjectionX("vtx");
+  hdNdEta    = hEtaVsVtx->ProjectionY("dNdEta");
  
   hEtaVsVtx->SetXTitle("vtx z [cm]");
   hEtaVsVtx->SetYTitle("#eta");
@@ -20,13 +27,13 @@ dNdEtaAnalysis::dNdEtaAnalysis(Char_t* name) {
 
   hEtaVsVtx->Sumw2();
   hVtx->Sumw2();
-
 }
 
 //____________________________________________________________________
 void
 dNdEtaAnalysis::FillTrack(Float_t vtx, Float_t eta, Float_t weight) {
   hEtaVsVtx->Fill(vtx, eta, weight);
+  hEtaVsVtxUncorrected->Fill(vtx,eta);
 }
 
 //____________________________________________________________________
@@ -93,8 +100,8 @@ dNdEtaAnalysis::Finish() {
       dndeta = dndeta/hdNdEta->GetBinWidth(i_eta);
       error  = error/hdNdEta->GetBinWidth(i_eta);
 
-      hdNdEta->SetBinContent(i_eta, dndeta);    
-      hdNdEta->SetBinError(i_eta, error);    
+      hdNdEta->SetBinContent(i_eta, dndeta);
+      hdNdEta->SetBinError(i_eta, error);
     }
 
   }
@@ -116,3 +123,25 @@ dNdEtaAnalysis::SaveHistograms() {
   gDirectory->cd("../");
 }
 
+//____________________________________________________________________
+void dNdEtaAnalysis::DrawHistograms()
+{
+  TCanvas* canvas = new TCanvas("dNdEtaAnalysis", "dNdEtaAnalysis", 800, 800);
+  canvas->Divide(2, 2);
+
+  canvas->cd(1);
+  if (hEtaVsVtx)
+    hEtaVsVtx->Draw("COLZ");
+
+  canvas->cd(2);
+  if (hEtaVsVtxUncorrected)
+    hEtaVsVtxUncorrected->Draw("COLZ");
+
+  canvas->cd(3);
+  if (hVtx)
+    hVtx->Draw();
+
+  canvas->cd(4);
+  if (hdNdEta)
+    hdNdEta->Draw();
+}
index 818ba196aba693176a22d281d21ca16b379782f5..067abf8bf3d530aafc56dcb2ae16f81591279b0e 100644 (file)
 // - figure out correct way to treat the errors
 // - add functionality to make dn/deta for different mult classes?
 
-#ifndef ROOT_TObject
-#include "TObject.h"
-#endif
-#ifndef ROOT_TFile
-#include "TFile.h"
-#endif
-#ifndef ROOT_TH2
-#include "TH2.h"
-#endif
-#ifndef ROOT_TMath
-#include "TMath.h"
-#endif
+#include <TObject.h>
+#include <TString.h>
 
-class dNdEtaAnalysis : public TObject 
-{
-protected:    
-
-  TString  fName; 
-
-  TH2F* hEtaVsVtx;
-  TH1F* hVtx;
-  TH1F* hdNdEta;
+class TH2F;
+class TH1D;
 
+class dNdEtaAnalysis : public TObject
+{
 public:
   dNdEtaAnalysis(Char_t* name="dndeta_correction");
 
   void FillTrack(Float_t vtx, Float_t eta, Float_t weight);
   void FillEvent(Float_t vtx);
-  
+
   void Finish();
-  
+
+  void DrawHistograms();
   void SaveHistograms();
-  
+
+  TH2F* GetEtaVsVtxHistogram() { return hEtaVsVtx; }
+  TH2F* GetEtaVsVtxUncorrectedHistogram() { return hEtaVsVtxUncorrected; }
+  TH1D* GetVtxHistogram() { return hVtx; }
+  TH1D* GetdNdEtaHistogram() { return hdNdEta; }
+
+  void SetEtaVsVtxHistogram(TH2F* aHist) { hEtaVsVtx = aHist; }
+  void SetEtaVsVtxUncorrectedHistogram(TH2F* aHist) { hEtaVsVtxUncorrected = aHist; }
+  void SetVtxHistogram(TH1D* aHist) { hVtx = aHist; }
+  void SetdNdEtaHistogram(TH1D* aHist) { hdNdEta = aHist; }
+
+protected:
+  TString  fName;
+
+  TH2F* hEtaVsVtx;
+  TH2F* hEtaVsVtxUncorrected;
+  TH1D* hVtx;
+  TH1D* hdNdEta;
+
   ClassDef(dNdEtaAnalysis,0)
 };
 
index c509646a9ce4869cf3e4f6798db621897fd0f7ad..fcb266f9fa7d11025c610a241820090799a2badf 100644 (file)
@@ -64,9 +64,10 @@ testAnalysis(Char_t* dataDir, Int_t nRuns=20) {
     if (!presentDir->IsDirectory())
       continue;
     // check that the files are there
-    sprintf(str,"%s/%s",dataDir, presentDir->GetName());
-    if ((!gSystem->Which(str,"galice.root")) ||
-        (!gSystem->Which(str,"AliESDs.root"))) 
+    TString currentDataDir;
+    currentDataDir.Form("%s/%s",dataDir, presentDir->GetName());
+    if ((!gSystem->Which(currentDataDir.Data(),"galice.root")) ||
+        (!gSystem->Which(currentDataDir.Data(),"AliESDs.root")))
       continue;
     
     if (nRunCounter++ >= nRuns)
@@ -82,7 +83,7 @@ testAnalysis(Char_t* dataDir, Int_t nRuns=20) {
       gAlice=0;
     }
 
-    sprintf(str,"%s/run%d/galice.root",dataDir,r);
+    sprintf(str,"%s/galice.root",currentDataDir.Data());
     AliRunLoader* runLoader = AliRunLoader::Open(str);
 
     runLoader->LoadgAlice();
@@ -93,7 +94,7 @@ testAnalysis(Char_t* dataDir, Int_t nRuns=20) {
     // #########################################################
     // open esd file and get the tree
 
-    sprintf(str,"%s/run%d/AliESDs.root",dataDir,r);
+    sprintf(str,"%s/AliESDs.root",currentDataDir.Data());
     // close it first to avoid memory leak
     if (esdFile)
       if (esdFile->IsOpen())