]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Added separate Invariant mass tab for CALO hist
authorslindal <slindal@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 14 May 2010 15:39:45 +0000 (15:39 +0000)
committerslindal <slindal@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 14 May 2010 15:39:45 +0000 (15:39 +0000)
HLT/EVE/AliHLTEveCalo.cxx
HLT/EVE/AliHLTEveCalo.h

index ccba7643553d6fd484b50f7326d8ccb27b1a9b0c..34a3f7d2f02a51944e0dbdd6dae9168350cad923 100644 (file)
@@ -46,18 +46,19 @@ AliHLTEveCalo::AliHLTEveCalo(Int_t nm, TString name) :
   fElementList(NULL),
   fNModules(nm),
   fName(name), 
-  fPadTitles(NULL)
+  fPadTitles(NULL),
+  fInvMassCanvas(NULL)
 {
   // Constructor.
 
-  SetMaxHistograms(9);
+  SetMaxHistograms(6);
 
   fPadTitles = new TString[GetMaxHistograms()];
-
   for(int i = 0; i < GetMaxHistograms(); i++) {
     fPadTitles[i] = "";
   }
-  
+
 
 }
 
@@ -77,6 +78,11 @@ AliHLTEveCalo::~AliHLTEveCalo()
     delete fElementList;
   }
   fElementList = NULL;
+
+  if(fPadTitles)
+    delete [] fPadTitles;
+  fPadTitles = NULL;
+
 }
 
 
@@ -111,9 +117,15 @@ void AliHLTEveCalo::ProcessHistogram(AliHLTHOMERBlockDesc * block ) {
   
   if(!fCanvas) {
     fCanvas = CreateCanvas(Form("%s QA", fName.Data()), Form("%s QA", fName.Data()));
-    fCanvas->Divide(3, 3);
+    fCanvas->Divide(3, 2);
   }
 
+  if(!fInvMassCanvas) {
+    fInvMassCanvas = CreateCanvas(Form("%s IM", fName.Data()), Form("%s IM", fName.Data()));
+    fInvMassCanvas->Divide(3, 2);
+  }
+
+
   AddHistogramsToCanvas(block, fCanvas, fHistoCount);
 
 
@@ -164,6 +176,8 @@ void AliHLTEveCalo::ProcessClusters(AliHLTHOMERBlockDesc* block) {
 void AliHLTEveCalo::UpdateElements() {
   //See header file for documentation
   if(fCanvas) fCanvas->Update();
+  if(fInvMassCanvas) fInvMassCanvas->Update();
+
 
   if(fBoxSetDigits) {
     for(int im = 0; im < fNModules; im++) {
@@ -220,59 +234,72 @@ Int_t AliHLTEveCalo::GetPadNumber(TString name) {
 
 }
 
+void AliHLTEveCalo::DrawInvMassHistogram(TH1F * histo) {
+  
+  fInvMassCanvas->cd(++fHistoCount);
+
+  histo->SetAxisRange(histo->GetXaxis()->GetBinLowEdge(histo->FindFirstBinAbove(0, 1) - 3), histo->GetXaxis()->GetBinUpEdge(histo->FindLastBinAbove(0, 1) + 3), "X");
+  histo->Fit("gaus", "", "", histo->GetXaxis()->GetBinLowEdge(histo->FindFirstBinAbove(0, 1)), histo->GetXaxis()->GetBinUpEdge(histo->FindLastBinAbove(0, 1)));
+
+  histo->Draw();
+
+}
+
 void AliHLTEveCalo::AddHistogramsToCanvas(AliHLTHOMERBlockDesc * block, TCanvas * canvas, Int_t &/*cdCount*/) {
   //See header file for documentation
 
   if ( ! block->GetClassName().CompareTo("TObjArray")) {
     TIter next((TObjArray*)(block->GetTObject()));
     TObject *object;
-   
-    while (( object = (TObject*) next())) {
-      
-      Int_t iPad = GetPadNumber((static_cast<TH1*>(object))->GetName());
-      canvas->cd(iPad);
 
-  
 
-      //Check if histo is 2D histo
-      TH2F* histo2 = dynamic_cast<TH2F*>(object);
-      if(histo2){
+   
+    while (( object = (TObject*) next())) {
 
-       Int_t lb = histo2->FindLastBinAbove(0,1);
-       if(lb > -1) {
-         histo2->SetAxisRange(0, histo2->GetXaxis()->GetBinUpEdge(histo2->FindLastBinAbove(0, 1) + 3), "X");
-         histo2->SetAxisRange(0, histo2->GetYaxis()->GetBinUpEdge(histo2->FindLastBinAbove(0, 2) + 3), "Y");
-       }
+      TString name = static_cast<TH1*>(object)->GetName();
+      if(name.Contains("InvMass")){
+       DrawInvMassHistogram(static_cast<TH1F*>(object));
        
-       histo2->Draw("COLZ");
-      }
-      
+      } else {
 
-      //Must be 1D histo
-      else {
-       TH1F* histo = dynamic_cast<TH1F*>(object);
-       if (histo) {
-
-         TString name = histo->GetName();
+       
+       Int_t iPad = GetPadNumber(name);
+       canvas->cd(iPad);
+       
+       
+       //Check if histo is 2D histo
+       TH2F* histo2 = dynamic_cast<TH2F*>(object);
+       if(histo2){
          
-         if(name.Contains("Energy")) {
-           histo->SetAxisRange(0, histo->GetXaxis()->GetBinUpEdge(histo->FindLastBinAbove(0, 1) + 3), "X");
+         Int_t lb = histo2->FindLastBinAbove(0,1);
+         if(lb > -1) {
+           histo2->SetAxisRange(0, histo2->GetXaxis()->GetBinUpEdge(histo2->FindLastBinAbove(0, 1) + 3), "X");
+           histo2->SetAxisRange(0, histo2->GetYaxis()->GetBinUpEdge(histo2->FindLastBinAbove(0, 2) + 3), "Y");
          }
-
-         else if(name.Contains("InvMass")) {
-           histo->SetAxisRange(histo->GetXaxis()->GetBinLowEdge(histo->FindFirstBinAbove(0, 1) - 3), histo->GetXaxis()->GetBinUpEdge(histo->FindLastBinAbove(0, 1) + 3), "X");
-           histo->Fit("gaus", "", "", histo->GetXaxis()->GetBinLowEdge(histo->FindFirstBinAbove(0, 1)), histo->GetXaxis()->GetBinUpEdge(histo->FindLastBinAbove(0, 1)));
-           //histo->Fit("gaus", "", "", histo->GetXaxis()->GetBinLowEdge(histo->FindFirstBinAbove(0, 1)), 0.3);
-           //histo->Fit("gaus", )
+         
+         histo2->Draw("COLZ");
+       }
+       
+       
+       //Must be 1D histo
+       else {
+         TH1F* histo = dynamic_cast<TH1F*>(object);
+         if (histo) {
+           
+           TString name = histo->GetName();
+           
+           if(name.Contains("Energy")) {
+             histo->SetAxisRange(0, histo->GetXaxis()->GetBinUpEdge(histo->FindLastBinAbove(0, 1) + 3), "X");
+           }
+           
+           
+           histo->Draw();
+         } else {
+           cout <<"AliHLTEveCaloBase::AddHistogramsTocCanvas: Histogram neither TH1F nor TH2F"<<endl;
          }
-
-         histo->Draw();
-       } else {
-         cout <<"AliHLTEveCaloBase::AddHistogramsTocCanvas: Histogram neither TH1F nor TH2F"<<endl;
        }
       }
     }
-  
 
   } else if ( ! block->GetClassName().CompareTo("TH1F")) {
 
index 904f340fd421660aa42ed102ca22ae483abbd967..2cb5b853d4c944b72f2e6fffff164ffeb60d9c33 100644 (file)
@@ -16,6 +16,8 @@
 class TEveElementList;
 class TEveBoxSet;
 class AliHLTHOMERBlockDesc;
+class TH1F;
+
 
 class AliHLTEveCalo : public AliHLTEveBase {
 
@@ -69,8 +71,6 @@ protected :
 
   const Int_t fNModules;          //Number of modules in calorimeter
 
-
-
 private:
   
   /** default constructor prohibited */
@@ -80,10 +80,16 @@ private:
   /** assignment operator prohibited */
   AliHLTEveCalo& operator = (const AliHLTEveCalo &);
 
+  void DrawInvMassHistogram(TH1F * histo);
+
   TString fName;  //PHOS or EMCAL
  
   TString * fPadTitles;
 
+
+  TCanvas * fInvMassCanvas;
+
+
   ClassDef(AliHLTEveCalo, 0);
 };