]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - PWG0/AliCorrectionMatrix.cxx
added correction for events with vertex but 0 tracks
[u/mrichter/AliRoot.git] / PWG0 / AliCorrectionMatrix.cxx
index 5137cbc91a000cad4098b7864565ae9e236163d0..66b704865523838c2f5192dd415e90f3964784f2 100644 (file)
@@ -157,11 +157,19 @@ Long64_t AliCorrectionMatrix::Merge(TCollection* list)
 void AliCorrectionMatrix::Divide()
 {
   //
-  // divide the histograms to get the correction
-  // 
+  // divides generated by measured to get the correction
+  //
+
+//   if (!fhCorr) {
+//     fhCorr = (TH1*)fhGene->Clone("correction");
+//     fhCorr->SetTitle(Form("%s correction",GetTitle()));
+//     fhCorr->Reset();
+//   }
 
-  if (!fhMeas || !fhGene)
+  if (!fhMeas || !fhGene || !fhCorr) {
+    AliDebug(AliLog::kError, "measured or generated histograms not available");
     return;
+  }
 
   fhCorr->Divide(fhGene, fhMeas, 1, 1, "B");
 
@@ -173,40 +181,86 @@ void AliCorrectionMatrix::Divide()
           ++emptyBins;
 
   if (emptyBins > 0)
-    printf("INFO: In %s we have %d empty bins (of %d) in the correction map\n", GetName(), emptyBins, fhCorr->GetNbinsX() * fhCorr->GetNbinsY() * fhCorr->GetNbinsZ());
+    printf("INFO: In %s we have %d empty bins (of %d) in the correction map\n", GetTitle(), emptyBins, fhCorr->GetNbinsX() * fhCorr->GetNbinsY() * fhCorr->GetNbinsZ());
+}
+
+//____________________________________________________________________
+void AliCorrectionMatrix::Multiply()
+{
+  //
+  // multiplies measured with correction to get the generated
+  //
+
+  if (!fhMeas || !fhGene || !fhCorr)
+    return;
+
+  fhGene->Multiply(fhMeas, fhCorr, 1, 1, "B");
+}
+
+//____________________________________________________________________
+void AliCorrectionMatrix::Add(AliCorrectionMatrix* aMatrixToAdd, Float_t c) {
+  //
+  // adds the measured and generated of aMatrixToAdd to measured and generated of this
+  // 
+  // NB: the correction will naturally stay the same
+  
+  fhMeas->Add(aMatrixToAdd->GetMeasuredHistogram(), c);
+  fhGene->Add(aMatrixToAdd->GetGeneratedHistogram(), c);
 }
 
+
 //____________________________________________________________________
-Bool_t AliCorrectionMatrix::LoadHistograms(const Char_t* fileName, const Char_t* dir)
+Bool_t AliCorrectionMatrix::LoadHistograms(const Char_t* dir)
 {
   //
   // loads the histograms from a file
+  // if dir is empty a directory with the name of this object is taken (like in SaveHistogram)
   //
 
-  TFile* fin = TFile::Open(fileName);
+  if (!dir)
+    dir = GetName();
 
-  if(!fin) {
-    //Info("LoadHistograms",Form(" %s file does not exist",fileName));
+  if (!gDirectory->cd(dir))
     return kFALSE;
+
+  if (fhGene)
+  {
+    delete fhGene;
+    fhGene=0;
   }
-  
-  if(fhGene)  {delete fhGene;  fhGene=0;}
-  if(fhCorr)  {delete fhCorr;  fhCorr=0;}
-  if(fhMeas)  {delete fhMeas;  fhMeas=0;}
-  
-  fhMeas  = dynamic_cast<TH1*> (fin->Get(Form("%s/meas_%s", dir,GetName())));
-  if(!fhMeas)  Info("LoadHistograms","No meas. (%s) hist available",Form("%s/meas_%s", dir,GetName()));
 
-  fhGene  = dynamic_cast<TH1*> (fin->Get(Form("%s/gene_%s",dir, GetName())));
-  if(!fhGene)  Info("LoadHistograms","No gene. (%s) hist available",Form("%s/gene_%s",dir, GetName()));
+  if (fhCorr)
+  {
+    delete fhCorr;
+    fhCorr=0;
+  }
 
-  fhCorr  = dynamic_cast<TH1*> (fin->Get(Form("%s/corr_%s",dir, GetName())));
-  if(!fhCorr) {
-    Info("LoadHistograms","No corr.(%s) hist available",Form("%s/corr_%s",dir, GetName()));
-    return kFALSE;
+  if (fhMeas)
+  {
+    delete fhMeas;
+    fhMeas=0;
   }
-      
-  return kTRUE;
+
+  fhMeas  = dynamic_cast<TH1*> (gDirectory->Get("measured"));
+  if (!fhMeas)
+    Info("LoadHistograms", "No measured hist available");
+
+  fhGene  = dynamic_cast<TH1*> (gDirectory->Get("generated"));
+  if (!fhGene)
+    Info("LoadHistograms", "No generated hist available");
+
+  fhCorr  = dynamic_cast<TH1*> (gDirectory->Get("correction"));
+
+  Bool_t success = kTRUE;
+  if (!fhCorr)
+  {
+    Info("LoadHistograms", "No correction hist available");
+    success = kFALSE;
+  }
+
+  gDirectory->cd("..");
+
+  return success;
 }
 
 //____________________________________________________________________
@@ -216,6 +270,9 @@ void AliCorrectionMatrix::SaveHistograms()
   // saves the histograms
   //
 
+  gDirectory->mkdir(GetName());
+  gDirectory->cd(GetName());
+
   if (fhMeas)
     fhMeas ->Write();
 
@@ -224,34 +281,40 @@ void AliCorrectionMatrix::SaveHistograms()
 
   if (fhCorr)
     fhCorr->Write();
+
+  gDirectory->cd("..");
 }
 
 //____________________________________________________________________
-void AliCorrectionMatrix::DrawHistograms()
+void AliCorrectionMatrix::DrawHistograms(const Char_t* canvasName)
 {
   //
-  // draws all the four histograms on one TCanvas
+  // draws all histograms on one TCanvas
+  // if canvasName is 0 the name of this object is taken
   //
 
-  TCanvas* canvas = new TCanvas(Form("correction_%s",fName.Data()), 
-                               Form("correction_%s",fName.Data()), 800, 800);
-  canvas->Divide(2, 2);
+  if (!canvasName)
+    canvasName = Form("%s_canvas", GetName());
+
+  TCanvas* canvas = new TCanvas(canvasName, GetTitle(), 1200, 400);
+  canvas->Divide(3, 1);
 
   canvas->cd(1);
   if (fhMeas)
     fhMeas->Draw("COLZ");
-  
+
   canvas->cd(2);
   if (fhGene)
+  {
+    // work around ROOT bug #22011
+    if (fhGene->GetEntries() == 0)
+      fhGene->SetEntries(1);
     fhGene->Draw("COLZ");
+  }
 
   canvas->cd(3);
   if (fhCorr)
     fhCorr->Draw("COLZ");
-
-  canvas->cd(4);
-
-  // add: draw here the stat. errors of the correction histogram
 }
 
 //____________________________________________________________________
@@ -272,3 +335,46 @@ void AliCorrectionMatrix::ReduceInformation()
     fhGene = 0;
   }
 }
+
+//____________________________________________________________________
+void AliCorrectionMatrix::Reset(Option_t* option)
+{
+  // resets the histograms
+
+  if (fhGene)
+    fhGene->Reset(option);
+
+  if (fhMeas)
+    fhMeas->Reset(option);
+
+  if (fhCorr)
+    fhCorr->Reset(option);
+}
+
+//____________________________________________________________________
+void AliCorrectionMatrix::SetCorrectionToUnity()
+{
+  // sets the correction matrix to unity
+
+  if (!fhCorr)
+    return;
+
+  for (Int_t x=1; x<=fhCorr->GetNbinsX(); ++x)
+    for (Int_t y=1; y<=fhCorr->GetNbinsY(); ++y)
+      for (Int_t z=1; z<=fhCorr->GetNbinsZ(); ++z)
+      {
+        fhCorr->SetBinContent(x, y, z, 1);
+        fhCorr->SetBinError(x, y, z, 0);
+      }
+}
+
+//____________________________________________________________________
+void AliCorrectionMatrix::Scale(Double_t factor)
+{
+  // scales the generated and measured histogram with the given factor
+
+  Printf("Scaling histograms with %f", factor);
+
+  fhMeas->Scale(factor);
+  fhGene->Scale(factor);
+}