]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/global/physics/AliHLTCaloHistoMatchedTracks.cxx
adding PrimaryVertexFinder and V0Finder component to registration, to be consolidated...
[u/mrichter/AliRoot.git] / HLT / global / physics / AliHLTCaloHistoMatchedTracks.cxx
1 //-*- Mode: C++ -*-
2 /**************************************************************************
3  * This file is property of and copyright by the ALICE HLT Project        * 
4  * All rights reserved.                                                   *
5  *                                                                        *
6  * Primary Authors: Albin Gaignette                                       *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          * 
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 /** 
18  * @file   AliHLTPHOSMatchedclustershistoProducer
19  * @author Albin Gaignette
20  * @date 
21  * @brief  Base Class for the Calo Matched track histograms
22  */
23
24 // see header file for class documentation
25 // or
26 // refer to README to build package
27 // or
28 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
29
30 #include "AliHLTCaloHistoMatchedTracks.h"
31 #include "AliESDCaloCluster.h"
32 #include "AliHLTCaloClusterDataStruct.h"
33 #include "TObjArray.h"
34 #include "TH1F.h"
35 #include "TH2F.h"
36 #include "TObjArray.h"
37 #include "TRefArray.h"
38 #include "TString.h"
39
40 AliHLTCaloHistoMatchedTracks::AliHLTCaloHistoMatchedTracks(TString det) :
41   fHistDxy(NULL),
42   fHistDz(NULL),
43   fHistDxyDz(NULL),
44   fHistMatchedEnergy(NULL),
45   fHistUnMatchedEnergy(NULL)
46 {
47
48   fHistMatchedEnergy = new TH1F( Form("%s_fHistMatchedEnergy", det.Data()), Form("%s Energy distribution of clusters with matching tracks", det.Data()), 5000, 0, 100);
49   fHistMatchedEnergy->GetXaxis()->SetTitle("Cluster Energy (GeV)");
50   fHistMatchedEnergy->GetYaxis()->SetTitle("Number of clusters");
51   fHistMatchedEnergy->SetMarkerStyle(21);
52   fHistArray->AddLast(fHistMatchedEnergy);
53
54   fHistUnMatchedEnergy = new TH1F( Form("%s_fHistUnMatchedEnergy", det.Data()), Form("%s Energy distribution of clusters with no matching track", det.Data()), 5000, 0, 100);
55   fHistUnMatchedEnergy->GetXaxis()->SetTitle("Cluster Energy (GeV)");
56   fHistUnMatchedEnergy->GetYaxis()->SetTitle("Number of clusters");
57   fHistUnMatchedEnergy->SetMarkerStyle(21);
58   fHistArray->AddLast(fHistUnMatchedEnergy);
59
60
61   fHistDxyDz = new TH2F( Form("%s_fHistdXYdZ", det.Data()), Form("%s dXY - dZ distribution of track - cluster residuals", det.Data()), 50, -50, 50, 50, -50, 50);
62   fHistDxyDz->GetXaxis()->SetTitle("sqrt(dx^2 + dy^2)  (cm)");
63   fHistDxyDz->GetYaxis()->SetTitle("dz (cm)");
64   fHistArray->AddLast(fHistDxyDz);
65
66   fHistDxy = new TH1F( Form("%s_fHistdXY", det.Data()), Form("%s #sqrt(dx^2 + dy^2)", det.Data()), 100, -50, 50);
67   fHistDxy->GetXaxis()->SetTitle("sqrt(dx^2 + dy^2)  (cm)");
68   fHistArray->AddLast(fHistDxy);
69                        
70   fHistDz = new TH1F( Form("%s_fHistdZ", det.Data()), Form("%s dZ", det.Data()),100, -50, 50);
71   fHistDz->GetXaxis()->SetTitle("dZ (cm)");
72   fHistArray->AddLast(fHistDz);
73
74 }
75
76
77 AliHLTCaloHistoMatchedTracks::~AliHLTCaloHistoMatchedTracks()
78 {
79
80   if(fHistMatchedEnergy) 
81     delete fHistMatchedEnergy;
82   fHistMatchedEnergy = NULL;
83
84   if(fHistUnMatchedEnergy) 
85     delete fHistUnMatchedEnergy;
86   fHistUnMatchedEnergy = NULL;
87
88   if (fHistDxyDz) 
89     delete fHistDxyDz;
90   fHistDxyDz = NULL;
91
92   if (fHistDxy) 
93     delete fHistDxy;
94   fHistDxy = NULL;
95
96   if (fHistDz) 
97     delete fHistDz;
98   fHistDz = NULL;
99
100 }
101   
102
103 Int_t AliHLTCaloHistoMatchedTracks::FillHistograms(Int_t nc, TRefArray * clusterArray) {
104   //See header file for documentation
105   for(int ic = 0; ic < nc; ic++) {
106     AliESDCaloCluster * cluster = static_cast<AliESDCaloCluster*>(clusterArray->At(ic));
107     FillMatchedTracks(cluster);
108   }
109   return 0;
110 }
111
112 Int_t AliHLTCaloHistoMatchedTracks::FillHistograms(Int_t nc, vector<AliHLTCaloClusterDataStruct*> &cVec) {
113   for(int ic = 0; ic < nc; ic++) {
114     AliHLTCaloClusterDataStruct * cluster = cVec.at(ic);
115     FillMatchedTracks(cluster);
116   }
117   return 0;
118 }
119
120 template <class T>
121 Int_t AliHLTCaloHistoMatchedTracks::FillMatchedTracks(T* cluster){
122   // HLTInfo("Filling track-matching histograms");
123
124   if(cluster->GetNTracksMatched() > 0) {
125     fHistMatchedEnergy->Fill(cluster->E());
126     fHistDz->Fill(cluster->GetTrackDz());
127     fHistDxy->Fill(cluster->GetTrackDx());
128     fHistDxyDz->Fill(cluster->GetTrackDx(), cluster->GetTrackDz());
129   } else {
130     fHistUnMatchedEnergy->Fill(cluster->E());
131   }
132   
133   return 0;
134 }
135