]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
o) moving includes to cxx
authorjgrosseo <jgrosseo@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 29 May 2006 12:13:36 +0000 (12:13 +0000)
committerjgrosseo <jgrosseo@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 29 May 2006 12:13:36 +0000 (12:13 +0000)
o) possibility for dNdEta in different vertices in dNdEtaAnalysis
o) Adding Merge function for Outputlist of PROOF

PWG0/dNdEta/AlidNdEtaAnalysisESDSelector.cxx
PWG0/dNdEta/AlidNdEtaAnalysisMCSelector.cxx
PWG0/dNdEta/AlidNdEtaAnalysisSelector.cxx
PWG0/dNdEta/AlidNdEtaAnalysisSelector.h
PWG0/dNdEta/dNdEtaAnalysis.cxx
PWG0/dNdEta/dNdEtaAnalysis.h

index 0a32b3fc2f4d3a7bb9377e60d33608dee3cb06b0..03f140a118e1e8bfdccc2a92492f286b84711d73 100644 (file)
@@ -6,9 +6,9 @@
 #include <TSystem.h>
 #include <TCanvas.h>
 #include <TVector3.h>
+#include <TChain.h>
 
 #include <AliLog.h>
-#include <AliGenEventHeader.h>
 
 #include "esdTrackCuts/AliESDtrackCuts.h"
 #include "dNdEtaCorrection.h"
index f74e1a85f64042472b2300315e176f88dc9b2cee..b486000878b749ec03cd9aaff26f6439856b2176 100644 (file)
@@ -9,9 +9,11 @@
 #include <TParticlePDG.h>
 #include <TVector3.h>
 #include <TH3F.h>
+#include <TTree.h>
 
 #include <AliLog.h>
 #include <AliGenEventHeader.h>
+#include <AliHeader.h>
 
 #include "dNdEtaAnalysis.h"
 
index 26e2a6ea9cf34d68bcf9adfeade836f519a64505..d2fe1b64bbdb4915eb0ce57c999e8bf1ea5f6af5 100644 (file)
@@ -7,9 +7,12 @@
 #include <TCanvas.h>
 #include <TVector3.h>
 #include <TH2F.h>
+#include <TChain.h>
+#include <TFile.h>
 
 #include <AliLog.h>
 #include <AliGenEventHeader.h>
+#include <AliHeader.h>
 
 #include "dNdEtaAnalysis.h"
 
@@ -17,8 +20,7 @@ ClassImp(AlidNdEtaAnalysisSelector)
 
 AlidNdEtaAnalysisSelector::AlidNdEtaAnalysisSelector() :
   AliSelector(),
-  fdNdEtaAnalysis(0),
-  fdNdEtaAnalysisFinal(0)
+  fdNdEtaAnalysis(0)
 {
   //
   // Constructor. Initialization of pointers
@@ -43,7 +45,7 @@ void AlidNdEtaAnalysisSelector::SlaveBegin(TTree * tree)
 
   AliSelector::SlaveBegin(tree);
 
-  fdNdEtaAnalysis = new dNdEtaAnalysis("dndeta");
+  fdNdEtaAnalysis = new dNdEtaAnalysis("dndeta", "dndeta");
 }
 
 void AlidNdEtaAnalysisSelector::SlaveTerminate()
@@ -61,12 +63,7 @@ void AlidNdEtaAnalysisSelector::SlaveTerminate()
     return;
   }
 
-  fOutput->Add(fdNdEtaAnalysis->GetEtaVsVtxHistogram());
-  fOutput->Add(fdNdEtaAnalysis->GetEtaVsVtxUncorrectedHistogram());
-  fOutput->Add(fdNdEtaAnalysis->GetVtxHistogram());
-
-  fdNdEtaAnalysis->GetVtxHistogram()->Print();
-  fOutput->Print();
+  fOutput->Add(fdNdEtaAnalysis);
 }
 
 void AlidNdEtaAnalysisSelector::Terminate()
@@ -77,30 +74,22 @@ void AlidNdEtaAnalysisSelector::Terminate()
 
   AliSelector::Terminate();
 
-  TH2F* etaVsVtxHistogram = dynamic_cast<TH2F*> (fOutput->FindObject("eta_vs_vtx"));
-  TH2F* etaVsVtxUncorrectedHistogram = dynamic_cast<TH2F*> (fOutput->FindObject("eta_vs_vtx_uncorrected"));
-  TH1D* vtxHistogram = dynamic_cast<TH1D*> (fOutput->FindObject("vtx"));
+  fdNdEtaAnalysis = dynamic_cast<dNdEtaAnalysis*> (fOutput->FindObject("dndeta"));
 
-  if (!etaVsVtxHistogram || !vtxHistogram || !etaVsVtxUncorrectedHistogram)
+  if (!fdNdEtaAnalysis)
   {
-     AliDebug(AliLog::kError, Form("ERROR: Histograms not available %p %p %p", (void*) etaVsVtxHistogram, (void*) etaVsVtxUncorrectedHistogram, (void*) vtxHistogram));
+    AliDebug(AliLog::kError, Form("ERROR: Histograms not available %p", (void*) fdNdEtaAnalysis));
     return;
   }
 
-  fdNdEtaAnalysisFinal = new dNdEtaAnalysis("dNdEtaResult");
-
-  fdNdEtaAnalysisFinal->SetEtaVsVtxHistogram(etaVsVtxHistogram);
-  fdNdEtaAnalysisFinal->SetEtaVsVtxUncorrectedHistogram(etaVsVtxUncorrectedHistogram);
-  fdNdEtaAnalysisFinal->SetVtxHistogram(vtxHistogram);
-
-  fdNdEtaAnalysisFinal->Finish();
+  fdNdEtaAnalysis->Finish();
 
   TFile* fout = new TFile("out.root","RECREATE");
   WriteObjects();
   fout->Write();
   fout->Close();
 
-  fdNdEtaAnalysisFinal->DrawHistograms();
+  fdNdEtaAnalysis->DrawHistograms();
 }
 
 void AlidNdEtaAnalysisSelector::WriteObjects()
@@ -109,5 +98,5 @@ void AlidNdEtaAnalysisSelector::WriteObjects()
   // this is an extra function to be overloaded...
   //
 
-  fdNdEtaAnalysisFinal->SaveHistograms();
+  fdNdEtaAnalysis->SaveHistograms();
 }
index 552ba3a54e87e170f3a26abfdbba3b0a33067a9d..3bceec66f00e439b8949ea42be8b51b453e77ac4 100644 (file)
@@ -20,7 +20,6 @@ class AlidNdEtaAnalysisSelector : public AliSelector {
     virtual void WriteObjects();
 
     dNdEtaAnalysis* fdNdEtaAnalysis;      // contains the intermediate histograms (on each slave)
-    dNdEtaAnalysis* fdNdEtaAnalysisFinal; // contains the final histograms
 
  private:
     ClassDef(AlidNdEtaAnalysisSelector, 0);
index ae704487101dbacd3f5aa7a1b3d108dbec312f3f..31397c9dcc768ee67a8ad025007af4663d9bfb32 100644 (file)
@@ -7,24 +7,28 @@
 #include <TH1D.h>
 #include <TMath.h>
 #include <TCanvas.h>
+#include <TCollection.h>
+#include <TIterator.h>
+#include <TList.h>
 
 //____________________________________________________________________
 ClassImp(dNdEtaAnalysis)
 
 //____________________________________________________________________
-dNdEtaAnalysis::dNdEtaAnalysis(Char_t* name) {
-
-  fName = TString(name);
-
-  hEtaVsVtx  = new TH2F("eta_vs_vtx","",80,-20,20,120,-6,6);
+dNdEtaAnalysis::dNdEtaAnalysis(Char_t* name, Char_t* title) :
+TNamed(name, title)
+{
+  hEtaVsVtx  = new TH2F(Form("%s_eta_vs_vtx", name),"",80,-20,20,120,-6,6);
   hEtaVsVtx->SetXTitle("vtx z [cm]");
   hEtaVsVtx->SetYTitle("#eta");
 
-  hEtaVsVtxUncorrected = dynamic_cast<TH2F*> (hEtaVsVtx->Clone("eta_vs_vtx_uncorrected"));
-  hVtx       = hEtaVsVtx->ProjectionX("vtx");
-  hdNdEta    = hEtaVsVtx->ProjectionY("dNdEta");
-
-  hdNdEta->SetYTitle("dN/d#eta");
+  hEtaVsVtxUncorrected = dynamic_cast<TH2F*> (hEtaVsVtx->Clone(Form("%s_eta_vs_vtx_uncorrected", name)));
+  hVtx       = hEtaVsVtx->ProjectionX(Form("%s_vtx", name));
+  for (Int_t i=0; i<kVertexBinning; ++i)
+  {
+    hdNdEta[i]    = hEtaVsVtx->ProjectionY(Form("%s_dNdEta_%d", name, i));
+    hdNdEta[i]->SetYTitle("dN/d#eta");
+  }
 
   hEtaVsVtx->Sumw2();
   hVtx->Sumw2();
@@ -76,50 +80,57 @@ dNdEtaAnalysis::Finish() {
     Float_t sum_xw = 0;
     Float_t sum_w  = 0;
     
-    for (Int_t i_vtx=0; i_vtx<=hVtx->GetNbinsX(); i_vtx++) {
-      if (hVtx->GetBinContent(i_vtx)==0)             continue;
-      if (hEtaVsVtx->GetBinContent(i_vtx, i_eta)==0) continue;
-
-      Float_t w = 1/TMath::Power(hEtaVsVtx->GetBinError(i_vtx, i_eta),2);
-      sum_xw = sum_xw + hEtaVsVtx->GetBinContent(i_vtx, i_eta)*w;
-      sum_w  = sum_w + w;
-
-      sum = sum + hEtaVsVtx->GetBinContent(i_vtx, i_eta);
-      sumError2 = sumError2 + TMath::Power(hEtaVsVtx->GetBinError(i_vtx, i_eta),2);      
-      nMeasurements++;
-    }
-    Float_t dndeta = 0;
-    Float_t error  = 0;
-
-    if (nMeasurements!=0) {
-      dndeta = sum/Float_t(nMeasurements);
-      error  = TMath::Sqrt(sumError2)/Float_t(nMeasurements);
-      
-      dndeta = sum_xw/sum_w;
-      error  = 1/TMath::Sqrt(sum_w);
-
-      dndeta = dndeta/hdNdEta->GetBinWidth(i_eta);
-      error  = error/hdNdEta->GetBinWidth(i_eta);
-
-      hdNdEta->SetBinContent(i_eta, dndeta);
-      hdNdEta->SetBinError(i_eta, error);
+    // do we have several histograms for different vertex positions?
+    Int_t vertexBinWidth = hVtx->GetNbinsX() / kVertexBinning;
+    for (Int_t vertexPos=0; vertexPos<kVertexBinning; ++vertexPos)
+    {
+      Int_t vertexBinBegin = vertexBinWidth * vertexPos;
+      Int_t vertexBinEnd = vertexBinBegin + vertexBinWidth;
+
+      for (Int_t i_vtx=vertexBinBegin; i_vtx<=vertexBinEnd; i_vtx++) {
+        if (hVtx->GetBinContent(i_vtx)==0)             continue;
+        if (hEtaVsVtx->GetBinContent(i_vtx, i_eta)==0) continue;
+
+        Float_t w = 1/TMath::Power(hEtaVsVtx->GetBinError(i_vtx, i_eta),2);
+        sum_xw = sum_xw + hEtaVsVtx->GetBinContent(i_vtx, i_eta)*w;
+        sum_w  = sum_w + w;
+
+        sum = sum + hEtaVsVtx->GetBinContent(i_vtx, i_eta);
+        sumError2 = sumError2 + TMath::Power(hEtaVsVtx->GetBinError(i_vtx, i_eta),2);      
+        nMeasurements++;
+      }
+      Float_t dndeta = 0;
+      Float_t error  = 0;
+
+      if (nMeasurements!=0) {
+        dndeta = sum/Float_t(nMeasurements);
+        error  = TMath::Sqrt(sumError2)/Float_t(nMeasurements);
+
+        dndeta = sum_xw/sum_w;
+        error  = 1/TMath::Sqrt(sum_w);
+
+        dndeta = dndeta/hdNdEta[vertexPos]->GetBinWidth(i_eta);
+        error  = error/hdNdEta[vertexPos]->GetBinWidth(i_eta);
+
+        hdNdEta[vertexPos]->SetBinContent(i_eta, dndeta);
+        hdNdEta[vertexPos]->SetBinError(i_eta, error);
+      }
     }
 
   }
 }
 
-
-
 //____________________________________________________________________
 void
 dNdEtaAnalysis::SaveHistograms() {
 
-  gDirectory->mkdir(fName.Data());
-  gDirectory->cd(fName.Data());
+  gDirectory->mkdir(GetName());
+  gDirectory->cd(GetName());
   
   hEtaVsVtx  ->Write();
   hVtx       ->Write();
-  hdNdEta    ->Write();
+  for (Int_t i=0; i<kVertexBinning; ++i)
+    hdNdEta[i]    ->Write();
 
   gDirectory->cd("../");
 }
@@ -143,6 +154,56 @@ void dNdEtaAnalysis::DrawHistograms()
     hVtx->Draw();
 
   canvas->cd(4);
-  if (hdNdEta)
-    hdNdEta->Draw();
+  if (hdNdEta[0])
+    hdNdEta[0]->Draw();
+}
+
+Long64_t dNdEtaAnalysis::Merge(TCollection* list)
+{
+  // Merges a list of dNdEtaAnalysis objects with this one.
+  // This is needed for PROOF.
+  // Returns the number of merged objects (including this)
+
+  if (!list)
+    return 0;
+
+  if (list->IsEmpty())
+    return 1;
+
+  TIterator* iter = list->MakeIterator();
+  TObject* obj;
+
+  // sub collections
+  const Int_t nCollections = kVertexBinning + 3;
+  TList* collections[nCollections];
+  for (Int_t i=0; i<nCollections; ++i)
+    collections[i] = new TList;
+
+  Int_t count = 0;
+  while ((obj = iter->Next()))
+  {
+    dNdEtaAnalysis* entry = dynamic_cast<dNdEtaAnalysis*> (obj);
+    if (entry == 0)
+      continue;
+
+    collections[0]->Add(entry->hEtaVsVtx);
+    collections[1]->Add(entry->hEtaVsVtxUncorrected);
+    collections[2]->Add(entry->hVtx);
+
+    for (Int_t i=0; i<kVertexBinning; ++i)
+      collections[3+i]->Add(entry->hdNdEta[i]);
+
+    ++count;
+  }
+
+  hEtaVsVtx->Merge(collections[0]);
+  hEtaVsVtxUncorrected->Merge(collections[1]);
+  hVtx->Merge(collections[2]);
+  for (Int_t i=0; i<kVertexBinning; ++i)
+    hdNdEta[i]->Merge(collections[3+i]);
+
+  for (Int_t i=0; i<nCollections; ++i)
+    delete collections[i];
+
+  return count+1;
 }
index 86783f4549157715d127ee0bd5597a62c6b8a57e..bc93ad6bd3a5e61d9fd3ef91ecf8a76f243007a6 100644 (file)
 // - add functionality to make dn/deta for different mult classes?
 // - implement destructor
 
-#include <TObject.h>
-#include <TString.h>
+#include <TNamed.h>
 
 class TH2F;
 class TH1D;
+class TCollection;
 
-class dNdEtaAnalysis : public TObject
+class dNdEtaAnalysis : public TNamed
 {
 public:
-  dNdEtaAnalysis(Char_t* name="dndeta_correction");
+  enum { kVertexBinning = 1 };
+
+  dNdEtaAnalysis(Char_t* name, Char_t* title);
 
   void FillTrack(Float_t vtx, Float_t eta, Float_t weight);
   void FillEvent(Float_t vtx);
@@ -37,23 +39,18 @@ public:
   void DrawHistograms();
   void SaveHistograms();
 
+  virtual Long64_t Merge(TCollection* list);
+
   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; }
+  TH1D* GetdNdEtaHistogram(Int_t i = 0) { return hdNdEta[i]; }
 
 protected:
-  TString  fName;
-
   TH2F* hEtaVsVtx;
   TH2F* hEtaVsVtxUncorrected;
   TH1D* hVtx;
-  TH1D* hdNdEta;
+  TH1D* hdNdEta[kVertexBinning];
 
   ClassDef(dNdEtaAnalysis,0)
 };