]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Fixed bug in invariant mass histo,
authorslindal <slindal@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sat, 27 Mar 2010 14:53:34 +0000 (14:53 +0000)
committerslindal <slindal@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sat, 27 Mar 2010 14:53:34 +0000 (14:53 +0000)
increased room on calo histo canvas

HLT/EVE/AliHLTEveCalo.cxx
HLT/global/physics/AliHLTCaloHistoClusterEnergy.cxx
HLT/global/physics/AliHLTCaloHistoInvMass.cxx
HLT/global/physics/AliHLTCaloHistoInvMass.h

index fc8bfed6848af5bf8e712dbca9cb7ab22b9843c1..21917c3b4694e46bd0f76e2f51fab18416882fec 100644 (file)
@@ -89,7 +89,7 @@ void AliHLTEveCalo::ProcessHistogram(AliHLTHOMERBlockDesc * block ) {
   
   if(!fCanvas) {
     fCanvas = CreateCanvas(Form("%s QA", fName.Data()), Form("%s QA", fName.Data()));
-    fCanvas->Divide(3, 2);
+    fCanvas->Divide(3, 3);
   }
 
   AddHistogramsToCanvas(block, fCanvas, fHistoCount);
index 233bb472f112cfcba228322ca414b45fe58c6518..7eec22ae0abfa9195cfa7807836664a75ac64376 100644 (file)
@@ -67,18 +67,20 @@ AliHLTCaloHistoClusterEnergy::~AliHLTCaloHistoClusterEnergy()
 
 Int_t AliHLTCaloHistoClusterEnergy::FillHistograms(Int_t nc, TRefArray * clusterArray) {
   //See header file for documentation
+
   for(int ic = 0; ic < nc; ic++) {
     AliESDCaloCluster * cluster = static_cast<AliESDCaloCluster*>(clusterArray->At(ic));
-    return FillClusterEnergyHistos(cluster);
+    FillClusterEnergyHistos(cluster);
   }
   return 0;
 }
+
 Int_t AliHLTCaloHistoClusterEnergy::FillHistograms(Int_t nc, vector<AliHLTCaloClusterDataStruct*> &cVec) {
   //See header file for documentation
   // HLTInfo("histo");
   for(int ic = 0; ic < nc; ic++) {
     AliHLTCaloClusterDataStruct * cluster = cVec.at(ic);
-    return FillClusterEnergyHistos(cluster);
+    FillClusterEnergyHistos(cluster);
   }
   return 0;
 }
@@ -86,7 +88,6 @@ Int_t AliHLTCaloHistoClusterEnergy::FillHistograms(Int_t nc, vector<AliHLTCaloCl
 template <class T>
 Int_t AliHLTCaloHistoClusterEnergy::FillClusterEnergyHistos(T* cluster) {
   fHistClusterEnergy->Fill(cluster->E());
-  fHistClusterEnergyVsNCells->Fill(cluster->GetNCells(), cluster->E());
-
+  fHistClusterEnergyVsNCells->Fill(cluster->E(), cluster->GetNCells());
   return 0;
 }
index 37435d012d2a73368de5907cd6a19969bc8dda82..87d73af0f24abfc2a6328813243e03265b324b05 100644 (file)
 #include "TLorentzVector.h"
 
 AliHLTCaloHistoInvMass::AliHLTCaloHistoInvMass(TString det) :
-  fHistTwoClusterInvMass(NULL)
+  fHistTwoClusterInvMass(NULL),
+  fHistTwoClusterInvMass2(NULL)
 {
   // See header file for documentation
-  fHistTwoClusterInvMass = new TH1F(Form("%s fHistTwoClusterInvMass", det.Data()), Form("%s Invariant mass of two clusters PHOS", det.Data()), 200, 0, 1);
+  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("Number of counts");
+  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);
+  fHistTwoClusterInvMass2->GetXaxis()->SetTitle("m_{#gamma#gamma} GeV");
+  fHistTwoClusterInvMass2->GetYaxis()->SetTitle("Counts");
+  fHistTwoClusterInvMass2->SetMarkerStyle(21);
+  fHistArray->AddLast(fHistTwoClusterInvMass2);
+
+
 }
 
 AliHLTCaloHistoInvMass::~AliHLTCaloHistoInvMass()
@@ -54,6 +65,13 @@ AliHLTCaloHistoInvMass::~AliHLTCaloHistoInvMass()
   if(fHistTwoClusterInvMass)
     delete fHistTwoClusterInvMass;
   fHistTwoClusterInvMass = NULL;
+
+
+  if(fHistTwoClusterInvMass2)
+    delete fHistTwoClusterInvMass2;
+  fHistTwoClusterInvMass2 = NULL;
+
+
 }
 
 
@@ -69,22 +87,11 @@ Int_t AliHLTCaloHistoInvMass::FillHistograms(Int_t nc, vector<AliHLTCaloClusterD
     cEnergy[ic] = cluster->E();
   }
 
-  for(Int_t ipho = 0; ipho<(nc-1); ipho++) { 
-    for(Int_t jpho = ipho+1; jpho<nc; jpho++) { 
-      
-      // Calculate the theta angle between two photons
-      Double_t theta = (2* asin(0.5*TMath::Sqrt((cPos[ipho][0]-cPos[jpho][0])*(cPos[ipho][0]-cPos[jpho][0]) +(cPos[ipho][1]-cPos[jpho][1])*(cPos[ipho][1]-cPos[jpho][1]))/460));
-      
-      // Calculate the mass m of the pion candidate
-      Double_t m =(TMath::Sqrt(2 * cEnergy[ipho]* cEnergy[jpho]*(1-TMath::Cos(theta))));
-      
-      fHistTwoClusterInvMass->Fill(m);
-    }
-  }
-
-  return 0;
+  return FillInvariantMassHistograms(nc, cPos, cEnergy);
 }
 
+
+
 Int_t AliHLTCaloHistoInvMass::FillHistograms(Int_t nc, TRefArray * clusterArray) {
   //See header file for documentation
   
@@ -97,29 +104,89 @@ Int_t AliHLTCaloHistoInvMass::FillHistograms(Int_t nc, TRefArray * clusterArray)
     cEnergy[ic] = cluster->E();
   }
 
+  return FillInvariantMassHistograms(nc, cPos, cEnergy);
+
+}
+
+
+Int_t AliHLTCaloHistoInvMass::FillInvariantMassHistograms(Int_t nc, Float_t cPos[][3], Float_t  cEnergy[]){ 
+
+  Int_t iResult = 0;
+
   for(Int_t ic = 0; ic<(nc-1); ic++) { 
+     
+    //BALLE hardcoded variable
+    if(cEnergy[ic] < 0.8)
+      continue;
 
     //Get momentum vector
     TVector3 iVec(cPos[ic]);
     iVec = iVec.Unit();
     iVec = cEnergy[ic] * iVec;
-
+    
     
     for(Int_t jc = ic+1; jc<nc; jc++) { 
+     
+    //BALLE hardcoded variable
+      if(cEnergy[jc] < 0.8)
+       continue;
 
+
+      
       //Get second momentum vector
       TVector3 jVec(cPos[jc]);
       jVec = jVec.Unit();
       jVec = cEnergy[jc] * jVec;
-
+      
       //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);
+      
     }
   }
 
+  return iResult;
 
-  return 0;
 }
+
+// template <class T>
+// Int_t AliHLTCaloHistoInvMass::FillHistograms(Int_t nc, vector<T*> clusterVec ) {
+//   Float_t cPos[nc][3];
+//   Float_t cEnergy[nc];
+  
+//   for(int ic = 0; ic < nc; ic++) {
+//     T * cluster = cVec.at(ic);
+//     cluster->GetPosition(cPos[ic]);
+//     cEnergy[ic] = cluster->E();
+//   }
+  
+
+//   for(Int_t ic = 0; ic<(nc-1); ic++) { 
+
+//     //Get momentum vector
+//     TVector3 iVec(cPos[ic]);
+//     iVec = iVec.Unit();
+//     iVec = cEnergy[ic] * iVec;
+
+    
+//     for(Int_t jc = ic+1; jc<nc; jc++) { 
+
+//       //Get second momentum vector
+//       TVector3 jVec(cPos[jc]);
+//       jVec = jVec.Unit();
+//       jVec = cEnergy[jc] * jVec;
+
+//       //Calculate inv mass
+//       Double_t m = TMath::Sqrt( 2 *(cEnergy[ic]* cEnergy[jc] - iVec.Dot(jVec) ) );
+
+//       //Fill histogram
+//       fHistTwoClusterInvMass->Fill(m);
+//     }
+//   }
+// }
index 214b9543ec2880d1a89d4f96a79eacd448ff80c2..600eb3c03716c2c19903d5c6cd431a076fc0c141 100644 (file)
@@ -67,6 +67,9 @@ public:
   //** Loops of the calo clusters and fills histos
   Int_t FillHistograms(Int_t nc, TRefArray * clusterArray);
   Int_t FillHistograms(Int_t nc, vector<AliHLTCaloClusterDataStruct*> &cVec);
+  //  template <class T> Int_t FillHistograms(Int_t nc, vector<T*> clusterVec);
+
+
   //   template <class T>
   //   Int_t FillHistograms(Int_t nc, vector<T*> &cVec);
   
@@ -81,8 +84,12 @@ private:
   /** Assignment operator */
   AliHLTCaloHistoInvMass & operator= (const AliHLTCaloHistoInvMass &);
   
+  /** Calculate 2 cluster inv mass and fill histograms */
+  Int_t FillInvariantMassHistograms(Int_t nc, Float_t cPos[][3], Float_t cEnergy[]);
+
   /** Histogram of the 2 cluster invariant mass */
   TH1F *fHistTwoClusterInvMass;                 //!transient
+  TH1F *fHistTwoClusterInvMass2;                 //!transient
 
   ClassDef(AliHLTCaloHistoInvMass, 0);