]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
updated cluster histo component for ITS (Gaute)
authorrichterm <richterm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 27 Jul 2009 14:13:05 +0000 (14:13 +0000)
committerrichterm <richterm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 27 Jul 2009 14:13:05 +0000 (14:13 +0000)
    - Removed 3D plot, used to much memory
    - Added documentation
    - New plot: phi vs eta

HLT/ITS/AliHLTITSClusterHistoComponent.cxx
HLT/ITS/AliHLTITSClusterHistoComponent.h

index ea7de6c1fd8ae22ba8b4434a02e4ce6b615a927e..258efc478de8a49edc7377b91bb1727398b020b8 100644 (file)
@@ -29,6 +29,7 @@ using namespace std;
 #include "AliCDBEntry.h"
 #include "AliCDBManager.h"
 #include "AliITSRecPoint.h"
+#include "TMath.h"
 #include <TFile.h>
 #include <TString.h>
 #include "TObjString.h"
@@ -43,11 +44,11 @@ ClassImp(AliHLTITSClusterHistoComponent)
 AliHLTITSClusterHistoComponent::AliHLTITSClusterHistoComponent()
 :
 fXY(NULL),                     
-  fXYZ(NULL),                   
+  fPhieta(NULL),                   
   fCharge(NULL),   
   fPlotCharge(kFALSE),   
   fPlotXY(kTRUE),
-  fPlotXYZ(kFALSE) 
+  fPlotPhieta(kFALSE) 
 {
   // see header file for class documentation
   // or
@@ -106,11 +107,11 @@ int AliHLTITSClusterHistoComponent::DoInit( int argc, const char** argv )
 {
   fPlotCharge=kFALSE;   
   fPlotXY=kTRUE;
-  fPlotXYZ=kFALSE; 
+  fPlotPhieta=kFALSE; 
      
-  if(fPlotCharge){fCharge = new TH1F("fCharge","Total Charge of clusters",4000,0,4000);}
+  if(fPlotCharge){fCharge = new TH1F("fCharge","Total Charge of clusters",2000,0,2000);}
   if(fPlotXY){fXY = new TH2F("fXY","Global XY of ITS clusters",1600,-80,80,1600,-80,80);}
-  if(fPlotXYZ){fXYZ = new TH3F("fXYZ","Global XYZ of ITS clusters",1600,-80,80,1600,-80,80,2000,-100,100);}
+  if(fPlotPhieta){fPhieta = new TH2F("fPhieta","Global Phieta of ITS clusters",30,-1.5,1.5,60,0,2*TMath::Pi());}
   
   int iResult=0;
   TString configuration="";
@@ -133,7 +134,7 @@ int AliHLTITSClusterHistoComponent::DoDeinit()
   // see header file for class documentation
   if(fCharge!=NULL) delete fCharge;
   if(fXY!=NULL) delete fXY;     
-  if(fXYZ!=NULL) delete fXYZ;
+  if(fPhieta!=NULL) delete fPhieta;
   return 0;
 }
 
@@ -184,8 +185,14 @@ int AliHLTITSClusterHistoComponent::DoEvent(const AliHLTComponentEventData& /*ev
       if(fPlotXY){
        fXY->Fill(xyz[0],xyz[1]);
       }
-      if(fPlotXYZ){
-       fXYZ->Fill(xyz[0],xyz[1],xyz[2]);
+      if(fPlotPhieta){
+       Float_t rad=TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]); 
+       Float_t theta=TMath::ATan2(rad,xyz[2]);
+       Float_t eta=-1*TMath::Log(TMath::Tan(theta/2.0));
+       Float_t phi=TMath::ATan2(xyz[1],xyz[0]);
+       if(phi<0.0){phi=2 * TMath::Pi() - TMath::Abs(phi);} 
+       //fPhieta->Fill(theta,phi);
+       fPhieta->Fill(eta,phi);
       }
       if(fPlotCharge){
        fCharge->Fill(recpoint.GetQ());
@@ -201,9 +208,9 @@ int AliHLTITSClusterHistoComponent::DoEvent(const AliHLTComponentEventData& /*ev
     AliHLTUInt32_t fSpecification = 0x0;
     PushBack( (TObject*) fXY,kAliHLTDataTypeHistogram,fSpecification);
   }
-  if(fPlotXYZ){
+  if(fPlotPhieta){
     AliHLTUInt32_t fSpecification = 0x0;
-    PushBack( (TObject*) fXYZ,kAliHLTDataTypeHistogram,fSpecification);
+    PushBack( (TObject*) fPhieta,kAliHLTDataTypeHistogram,fSpecification);
   }
   
   HLTInfo("ITSClusterHisto found %d Total Spacepoints", TotalSpacePoint);
@@ -231,7 +238,7 @@ int AliHLTITSClusterHistoComponent::Configure(const char* arguments)
       if (argument.CompareTo("-plot-all")==0) {
        HLTInfo("Ploting all historgams");
        fPlotXY = kTRUE;
-       fPlotXYZ = kTRUE;
+       fPlotPhieta = kTRUE;
        fPlotCharge = kTRUE;
        continue;
       }
@@ -242,9 +249,9 @@ int AliHLTITSClusterHistoComponent::Configure(const char* arguments)
        continue;
       }
 
-      else if (argument.CompareTo("-plot-xyz")==0) {
-       HLTInfo("Ploting Global XYZ");
-       //fPlotXYZ = kTRUE;
+      else if (argument.CompareTo("-plot-phieta")==0) {
+       HLTInfo("Ploting Global Phieta");
+       fPlotPhieta = kTRUE;
        continue;
       }
       else if (argument.CompareTo("-plot-charge")==0) {
@@ -262,9 +269,9 @@ int AliHLTITSClusterHistoComponent::Configure(const char* arguments)
     delete pTokens;
   }
   
-  if(!fCharge && fPlotCharge){fCharge = new TH1F("fCharge","Total Charge of clusters",4000,0,4000);}
+  if(!fCharge && fPlotCharge){fCharge = new TH1F("fCharge","Total Charge of clusters",2000,0,2000);}
   if(!fXY && fPlotXY){fXY = new TH2F("fXY","Global XY of ITS clusters",1600,-80,80,1600,-80,80);}
-  if(!fXYZ && fPlotXYZ){fXYZ = new TH3F("fXYZ","Global XYZ of ITS clusters",1600,-80,80,1600,-80,80,2000,-100,100);}
+  if(!fPhieta && fPlotPhieta){fPhieta = new TH2F("fPhieta","Global Phieta of ITS clusters",30,-1.5,1.5,60,0,2*TMath::Pi());}
   
   return iResult;
 }
index 4037378242855da64479458821f91cfba4cba8a6..415b6c401d63336fa65dadfc3e7c826eb6079b2c 100644 (file)
@@ -6,7 +6,7 @@
 //* ALICE Experiment at CERN, All rights reserved.                         *
 //* See cxx source for full Copyright notice                               *
 
-/** @file   AliHLTITSQHistoComponent.h
+/** @file   AliHLTITSClusterHistoComponent.h
     @author Gaute Ovrebekk
     @brief  Component for ploting clusters
 */
 #include "AliHLTProcessor.h"
 #include "TH1F.h"
 #include "TH2F.h"
-#include "TH3F.h"
 #include "AliHLTITSSpacePointData.h"
 #include "TClonesArray.h"
 #include "AliITSRecPoint.h"
 
-class AliHLTTPCConfMapper;
-
 /**
  * @class AliHLTITSQHistoComponent
  * Component for ploting charge in clusters
  * 
- * Component ID: \b ITSQHisto <br>
- * Library: \b libAliHLTITS.
+ * Component ID: \b ITSClusterHisto <br>
+ * Library: \b libAliHLTITS.so
  *
  * Mandatory arguments: <br>
  * 
  * 
  * Optional arguments: <br>
- * 
- *
+ * \li -plot-all <br>  
+ *      will plot all histograms in the component.
+ * \li -plot-xy <br>  
+ *      will plot a historgam of the x and y projection of the clusters.
+ * \li -plot-charge     <i> multiplicity  </i> <br> 
+ *      will plot the charge of the clusters.
+ * \li -plot-phieta <br>
+ *      will plot the phi vs eta distrubution of the clusters.
  * @ingroup alihlt_tpc_components
  */
 class AliHLTITSClusterHistoComponent : public AliHLTProcessor
@@ -88,14 +91,14 @@ private:
   int Configure(const char* arguments);
   
   TH2F * fXY;                              //! transient
-  TH3F * fXYZ;                             //! transient
+  TH2F * fPhieta;                          //! transient
   TH1F * fCharge;                          //! transient
     
   Bool_t fPlotCharge;                      //! transient
   Bool_t fPlotXY;                          //! transient
-  Bool_t fPlotXYZ;                         //! transient
+  Bool_t fPlotPhieta;                      //! transient
    
-  ClassDef(AliHLTITSClusterHistoComponent, 0);
+  ClassDef(AliHLTITSClusterHistoComponent, 1);
 
 };
 #endif