]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
- adding more histgrams for the invariant mass component (more energy intervals)
authorodjuvsla <odjuvsla@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sun, 25 Apr 2010 13:54:50 +0000 (13:54 +0000)
committerodjuvsla <odjuvsla@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sun, 25 Apr 2010 13:54:50 +0000 (13:54 +0000)
- adding a 2D phi/eta plot of deposited energy

HLT/global/physics/AliHLTCaloHistoClusterEnergy.cxx
HLT/global/physics/AliHLTCaloHistoClusterEnergy.h
HLT/global/physics/AliHLTCaloHistoComponent.cxx
HLT/global/physics/AliHLTCaloHistoInvMass.cxx
HLT/global/physics/AliHLTCaloHistoInvMass.h

index 73e4510493622167012d44615a70707c92b16fe3..93b15788f80952da907045bf48c7a4ab2789f564 100644 (file)
 #include "TH1F.h"
 #include "TH2F.h"
 #include "TRefArray.h"
+#include "TVector3.h"
 
 ClassImp(AliHLTCaloHistoClusterEnergy);
 
 AliHLTCaloHistoClusterEnergy::AliHLTCaloHistoClusterEnergy(TString det) :
   AliHLTCaloHistoProducer(),
   fHistClusterEnergy(NULL),
-  fHistClusterEnergyVsNCells(NULL)
+  fHistClusterEnergyVsNCells(NULL),
+  fHistClusterEnergyDepositEtaPhi(NULL)
 {
   // See header file for documentation
   fHistClusterEnergy = new TH1F(Form("%s fHistClusterEnergy", det.Data()), Form("%s Distribution of total energy in clusters", det.Data()), 5000, 0, 100);
@@ -50,6 +52,24 @@ AliHLTCaloHistoClusterEnergy::AliHLTCaloHistoClusterEnergy(TString det) :
   fHistClusterEnergyVsNCells->GetXaxis()->SetTitle("Energy in cluster (GeV)");
   fHistClusterEnergyVsNCells->GetYaxis()->SetTitle("Number of Cells in cluster");
   fHistArray->AddLast(fHistClusterEnergyVsNCells);
+  
+  Float_t phiMin = 0.;
+  Float_t phiMax = 360.;
+  Float_t etaMin = -1.;
+  Float_t etaMax = 1.;
+  
+  if(det == "PHOS")
+  {
+     phiMin = 255.0;
+     phiMax = 325.0;
+     etaMin = -0.13;
+     etaMax = 0.13;
+  }
+  
+  fHistClusterEnergyDepositEtaPhi = new TH2F(Form("%s fHistClusterEnergyDepositedEtaPhi", det.Data()), Form("%s Amount of energy deposited in Phi vs Eta", det.Data()), 200, phiMin, phiMax, 50, etaMin , etaMax);
+  fHistClusterEnergyVsNCells->GetXaxis()->SetTitle("#phi");
+  fHistClusterEnergyVsNCells->GetYaxis()->SetTitle("#eta");
+  fHistArray->AddLast(fHistClusterEnergyDepositEtaPhi);
 
 }
 
@@ -89,5 +109,12 @@ template <class T>
 Int_t AliHLTCaloHistoClusterEnergy::FillClusterEnergyHistos(T* cluster) {
   fHistClusterEnergy->Fill(cluster->E());
   fHistClusterEnergyVsNCells->Fill(cluster->E(), cluster->GetNCells());
+  
+  Float_t pos[3];
+  cluster->GetPosition(pos);
+  TVector3 vec(pos);
+  
+  fHistClusterEnergyDepositEtaPhi->Fill(vec.Phi(), vec.Eta(), cluster->E());
+  
   return 0;
 }
index 298ba5e97e9a68228959b638ff272e495f4f6592..bf0c3d34e70ce33f982baa06ac430faa899b69af 100644 (file)
@@ -88,6 +88,9 @@ public:
   /** 2D histogram of cluster energy vs the number of cells in the cluster */
   TH2F * fHistClusterEnergyVsNCells;         //!transient
 
+/** 2D histogram of cluster energy deposit in eta vs phi */
+  TH2F * fHistClusterEnergyDepositEtaPhi;         //!transient
+  
   ClassDef(AliHLTCaloHistoClusterEnergy, 0);
 
 };
index 202b72fa112076625454164f2732ea4a130cb119..d7c35671b9baca8efcdd7109797c474d0a1b4cc7 100644 (file)
@@ -315,9 +315,27 @@ Int_t AliHLTCaloHistoComponent::ProcessBlocks(const AliHLTComponentBlockData * p
   
   clustersVector.resize((int) (clusterHeader->fNClusters)); 
   Int_t nClusters = 0;
-  
+  Double_t ampFrac;
+  UShort_t cellId;
+  Bool_t cutCluster = false;
   while( (clusterStruct = fClusterReader->NextCluster()) != 0) {
-    clustersVector[nClusters++] = clusterStruct;  
+     cutCluster = false;
+     if(clusterStruct->fEnergy > 0.5)
+     {
+        for(UInt_t i = 0; i < clusterStruct->fNCells; i++)
+        {
+           fClusterReader->GetCell(clusterStruct, cellId, ampFrac, i);
+           if(ampFrac > 0.9)
+           {
+              cutCluster = true;
+              break;
+           }
+        }
+     if(!cutCluster)
+     {
+        clustersVector[nClusters++] = clusterStruct;  
+     }
+     }
   }
   
   nClusters = clusterHeader->fNClusters;
index 87d73af0f24abfc2a6328813243e03265b324b05..e60209adbbc073634939944edb7156711eaf6873 100644 (file)
 #include "TLorentzVector.h"
 
 AliHLTCaloHistoInvMass::AliHLTCaloHistoInvMass(TString det) :
-  fHistTwoClusterInvMass(NULL),
-  fHistTwoClusterInvMass2(NULL)
+  fHistTwoClusterInvMass0(NULL),
+  fHistTwoClusterInvMass1(NULL),
+  fHistTwoClusterInvMass2(NULL),
+  fHistTwoClusterInvMass3(NULL),
+  fHistTwoClusterInvMass4(NULL)
 {
   // See header file for documentation
-  fHistTwoClusterInvMass = new TH1F(Form("%s fHistTwoClusterInvMass", det.Data()), Form("Invariant mass of two clusters in %s, E > 0.8 GeV", det.Data()), 200, 0, 1);
-  fHistTwoClusterInvMass->GetXaxis()->SetTitle("m_{#gamma#gamma} GeV");
-  fHistTwoClusterInvMass->GetYaxis()->SetTitle("Counts");
-  fHistTwoClusterInvMass->SetMarkerStyle(21);
-  fHistArray->AddLast(fHistTwoClusterInvMass);
-
-
-  fHistTwoClusterInvMass2 = new TH1F(Form("%s fHistTwoClusterInvMass2", det.Data()), Form("Invariant mass of two clusters in %s E > 2.4 GeV", det.Data()), 200, 0, 1);
+  fHistTwoClusterInvMass0 = new TH1F(Form("%s fHistTwoClusterInvMass0", det.Data()), Form("Invariant mass of two clusters in %s, 0.8 GeV < E < 1.2", det.Data()), 200, 0, 1);
+  fHistTwoClusterInvMass0->GetXaxis()->SetTitle("m_{#gamma#gamma} GeV");
+  fHistTwoClusterInvMass0->GetYaxis()->SetTitle("Counts");
+  fHistTwoClusterInvMass0->SetMarkerStyle(21);
+  fHistArray->AddLast(fHistTwoClusterInvMass0);
+
+  fHistTwoClusterInvMass1 = new TH1F(Form("%s fHistTwoClusterInvMass1", det.Data()), Form("Invariant mass of two clusters in %s, 1.2 GeV < E < 1.6", det.Data()), 200, 0, 1);
+  fHistTwoClusterInvMass1->GetXaxis()->SetTitle("m_{#gamma#gamma} GeV");
+  fHistTwoClusterInvMass1->GetYaxis()->SetTitle("Counts");
+  fHistTwoClusterInvMass1->SetMarkerStyle(21);
+  fHistArray->AddLast(fHistTwoClusterInvMass1);
+
+  fHistTwoClusterInvMass2 = new TH1F(Form("%s fHistTwoClusterInvMass2", det.Data()), Form("Invariant mass of two clusters in %s, 1.6 GeV < E < 2.0", det.Data()), 200, 0, 1);
   fHistTwoClusterInvMass2->GetXaxis()->SetTitle("m_{#gamma#gamma} GeV");
   fHistTwoClusterInvMass2->GetYaxis()->SetTitle("Counts");
   fHistTwoClusterInvMass2->SetMarkerStyle(21);
   fHistArray->AddLast(fHistTwoClusterInvMass2);
 
+  fHistTwoClusterInvMass3 = new TH1F(Form("%s fHistTwoClusterInvMass3", det.Data()), Form("Invariant mass of two clusters in %s, 2.0 GeV < E < 4.0", det.Data()), 200, 0, 1);
+  fHistTwoClusterInvMass3->GetXaxis()->SetTitle("m_{#gamma#gamma} GeV");
+  fHistTwoClusterInvMass3->GetYaxis()->SetTitle("Counts");
+  fHistTwoClusterInvMass3->SetMarkerStyle(21);
+  fHistArray->AddLast(fHistTwoClusterInvMass3);
+
+  fHistTwoClusterInvMass4 = new TH1F(Form("%s fHistTwoClusterInvMass4", det.Data()), Form("Invariant mass of two clusters in %s E > 4.0 GeV", det.Data()), 200, 0, 1);
+  fHistTwoClusterInvMass4->GetXaxis()->SetTitle("m_{#gamma#gamma} GeV");
+  fHistTwoClusterInvMass4->GetYaxis()->SetTitle("Counts");
+  fHistTwoClusterInvMass4->SetMarkerStyle(21);
+  fHistArray->AddLast(fHistTwoClusterInvMass4);
+
 
 }
 
 AliHLTCaloHistoInvMass::~AliHLTCaloHistoInvMass()
 {
-  if(fHistTwoClusterInvMass)
-    delete fHistTwoClusterInvMass;
-  fHistTwoClusterInvMass = NULL;
 
+  if(fHistTwoClusterInvMass0)
+    delete fHistTwoClusterInvMass0;
+  fHistTwoClusterInvMass0 = NULL;
+
+  if(fHistTwoClusterInvMass1)
+    delete fHistTwoClusterInvMass1;
+  fHistTwoClusterInvMass1 = NULL;
 
   if(fHistTwoClusterInvMass2)
     delete fHistTwoClusterInvMass2;
   fHistTwoClusterInvMass2 = NULL;
 
+  if(fHistTwoClusterInvMass3)
+    delete fHistTwoClusterInvMass3;
+  fHistTwoClusterInvMass3 = NULL;
+
+  if(fHistTwoClusterInvMass4)
+    delete fHistTwoClusterInvMass4;
+  fHistTwoClusterInvMass4 = NULL;
+
 
 }
 
@@ -116,7 +147,7 @@ Int_t AliHLTCaloHistoInvMass::FillInvariantMassHistograms(Int_t nc, Float_t cPos
   for(Int_t ic = 0; ic<(nc-1); ic++) { 
      
     //BALLE hardcoded variable
-    if(cEnergy[ic] < 0.8)
+    if(cEnergy[ic] < 0.4)
       continue;
 
     //Get momentum vector
@@ -128,9 +159,10 @@ Int_t AliHLTCaloHistoInvMass::FillInvariantMassHistograms(Int_t nc, Float_t cPos
     for(Int_t jc = ic+1; jc<nc; jc++) { 
      
     //BALLE hardcoded variable
-      if(cEnergy[jc] < 0.8)
+      if(cEnergy[jc] < 0.4)
        continue;
 
+      
 
       
       //Get second momentum vector
@@ -141,13 +173,38 @@ Int_t AliHLTCaloHistoInvMass::FillInvariantMassHistograms(Int_t nc, Float_t cPos
       //Calculate inv mass
       Double_t m = TMath::Sqrt( 2 *(cEnergy[ic]* cEnergy[jc] - iVec.Dot(jVec) ) );
       
-      //Fill histogram
-      fHistTwoClusterInvMass->Fill(m);
-      
-      //BALLE hardcoded variable
-      if( (cEnergy[ic] > 2.4) && (cEnergy[jc] > 2.4))
-        fHistTwoClusterInvMass2->Fill(m);
+      //Fill histograms
       
+      Float_t sum = cEnergy[ic]+cEnergy[ic];
+      if(sum > 1.2)
+      {
+        if(sum > 1.6)
+        {
+           if(sum > 2.0)
+           {
+              if(sum > 4.0)
+              {
+                 fHistTwoClusterInvMass4->Fill(m);
+              }
+              else
+              {
+                 fHistTwoClusterInvMass3->Fill(m);
+              }
+           }
+           else
+           {
+              fHistTwoClusterInvMass2->Fill(m);
+           }
+        }
+        else
+        {
+           fHistTwoClusterInvMass1->Fill(m);
+        }
+      }
+      else
+      {
+        fHistTwoClusterInvMass0->Fill(m);
+      }
     }
   }
 
index 600eb3c03716c2c19903d5c6cd431a076fc0c141..2f80372c66b82a9022d47613387e20f40ec959ec 100644 (file)
@@ -88,8 +88,11 @@ private:
   Int_t FillInvariantMassHistograms(Int_t nc, Float_t cPos[][3], Float_t cEnergy[]);
 
   /** Histogram of the 2 cluster invariant mass */
-  TH1F *fHistTwoClusterInvMass;                 //!transient
+  TH1F *fHistTwoClusterInvMass0;                 //!transient
+  TH1F *fHistTwoClusterInvMass1;                 //!transient
   TH1F *fHistTwoClusterInvMass2;                 //!transient
+  TH1F *fHistTwoClusterInvMass3;                 //!transient
+  TH1F *fHistTwoClusterInvMass4;                 //!transient
 
   ClassDef(AliHLTCaloHistoInvMass, 0);