]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDClusterHistoComponent.cxx
- remove memory leak from TRD histo component (Theodor)
[u/mrichter/AliRoot.git] / HLT / TRD / AliHLTTRDClusterHistoComponent.cxx
1 //**************************************************************************
2 //* This file is property of and copyright by the ALICE HLT Project        * 
3 //* ALICE Experiment at CERN, All rights reserved.                         *
4 //*                                                                        *
5 //* Primary Authors: Sylwester Radomski radomski@physi.uni-heidelberg.de    *
6 //*                  for The ALICE HLT Project.                            *
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 /** @file   AliHLTTRDClusterHistoComponent.cxx
18     @author Sylwester Radomski
19     @brief  Component for ploting charge in clusters
20 */
21
22 #if __GNUC__>= 3
23 using namespace std;
24 #endif
25
26 #include "AliHLTTRDClusterHistoComponent.h"
27 #include "AliHLTTRDDefinitions.h"
28 #include "AliHLTTRDCluster.h"
29 #include "AliTRDcluster.h"
30 #include "AliCDBEntry.h"
31 #include "AliCDBManager.h"
32 #include <TFile.h>
33 #include <TString.h>
34 #include "TObjString.h"
35 #include "TClonesArray.h"
36 #include "AliHLTTRDUtils.h"
37
38 //#include "AliHLTTRD.h"
39 //#include <stdlib.h>
40 //#include <cerrno>
41
42 /** ROOT macro for the implementation of ROOT specific class methods */
43 ClassImp(AliHLTTRDClusterHistoComponent)
44
45 AliHLTTRDClusterHistoComponent::AliHLTTRDClusterHistoComponent()
46 : fNClsDet(0),
47   fClsAmp(0),
48   fClsAmpDrift(0),
49   fClsTB(0),
50   fClsAmpDist(0)
51 {
52   // see header file for class documentation
53   // or
54   // refer to README to build package
55   // or
56   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
57
58 }
59
60 AliHLTTRDClusterHistoComponent::~AliHLTTRDClusterHistoComponent()
61 {
62   // see header file for class documentation
63 }
64
65 // Public functions to implement AliHLTComponent's interface.
66 // These functions are required for the registration process
67
68 const char* AliHLTTRDClusterHistoComponent::GetComponentID()
69 {
70   // see header file for class documentation
71   
72   return "TRDClusterHisto";
73 }
74
75 void AliHLTTRDClusterHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
76 {
77   // see header file for class documentation
78   list.clear();
79   list.push_back( AliHLTTRDDefinitions::fgkClusterDataType );
80 }
81
82 AliHLTComponentDataType AliHLTTRDClusterHistoComponent::GetOutputDataType()
83 {
84   // see header file for class documentation
85   return kAliHLTDataTypeHistogram;
86
87 }
88
89 void AliHLTTRDClusterHistoComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
90 {
91   // see header file for class documentation
92   // XXX TODO: Find more realistic values.
93   constBase = 80000;
94   inputMultiplier = 1;
95 }
96
97 AliHLTComponent* AliHLTTRDClusterHistoComponent::Spawn()
98 {
99   // see header file for class documentation
100   return new AliHLTTRDClusterHistoComponent;
101 }
102
103 int AliHLTTRDClusterHistoComponent::DoInit( int argc, const char** argv )
104 {
105   // Initialize histograms
106
107   fNClsDet = new TH1D("trdClsDet", ";detector", 540, -0.5, 539.5);
108   fClsAmp  = new TH1D("trdClsAmp", ";amplitude", 200, -0.5, 199.5);
109   fClsAmpDrift = new TH1D("trdClsAmpDrift", ";amplitude", 200, -0.5, 199.5) ;
110   fClsTB = new TH1D("trdClsTB", ";time bin", 35, -0.5, 34.5);
111   fClsAmpDist = new TH1D("trdClsAmpDist", "mean amplitude", 200, 0, 100);
112
113   for(int i=0; i<540; i++)
114     fClsAmpDriftDet[i] = new TH1D(Form("trdClsDriftDet_%d",i), "", 200, -0.5, 199.5);
115
116   /*
117   // configure
118
119   int iResult=0;
120   TString configuration="";
121   TString argument="";
122   for (int i=0; i<argc && iResult>=0; i++) {
123     argument=argv[i];
124     if (!configuration.IsNull()) configuration+=" ";
125     configuration+=argument;
126   }
127   
128   if (!configuration.IsNull()) {
129     iResult=Configure(configuration.Data());
130   }  
131
132   return iResult; 
133     */
134     
135     return 0;
136 }
137   
138 int AliHLTTRDClusterHistoComponent::DoDeinit()
139 {
140   // see header file for class documentation
141
142   // delete histograms
143   if (fNClsDet) delete fNClsDet;
144   if (fClsAmp) delete fClsAmp;
145   if (fClsAmpDrift) delete fClsAmpDrift;
146   if (fClsTB) delete fClsTB;
147   if (fClsAmpDist) delete fClsAmpDist;
148
149   for(int i=0; i<540; i++)
150     if (fClsAmpDriftDet[i]) delete fClsAmpDriftDet[i];
151
152   return 0;
153 }
154
155 int AliHLTTRDClusterHistoComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/, 
156                                             AliHLTComponentTriggerData& /*trigData*/)
157 {
158
159   if ( GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR ) )
160     return 0;
161   
162   const AliHLTComponentBlockData* iter = NULL;
163   
164   for ( iter = GetFirstInputBlock(AliHLTTRDDefinitions::fgkClusterDataType); 
165         iter != NULL; iter = GetNextInputBlock() ) {
166
167     HLTDebug("We get the right data type: Block Ptr: 0x%x; Block Size: %i",
168              iter->fPtr, iter->fSize);
169
170     TClonesArray* clusterArray = new TClonesArray("AliTRDcluster");
171     AliHLTTRDUtils::ReadClusters(clusterArray, iter->fPtr, iter->fSize);
172     HLTDebug("TClonesArray of clusters: nbEntries = %i", clusterArray->GetEntriesFast());
173
174     AliTRDcluster *cls;
175         
176     // loop over clusters 
177     for(int i=0;i<clusterArray->GetEntriesFast();i++) {
178
179       cls=(AliTRDcluster*)clusterArray->At(i);
180       
181       fNClsDet->Fill(cls->GetDetector());
182       fClsAmp->Fill(cls->GetQ());
183       
184       int tb = cls->GetPadTime();
185       fClsTB->Fill(tb);
186       if (tb > 5 && tb <25)
187         fClsAmpDrift->Fill(cls->GetQ()); 
188       
189       fClsAmpDriftDet[cls->GetDetector()]->Fill(cls->GetQ());
190     }
191     
192     clusterArray->Delete();
193     delete clusterArray;
194     
195   }
196    
197   fClsAmpDist->Reset();
198   for(int det=0; det<540; det++)
199     if (fClsAmpDriftDet[det]->GetSum() > 0) 
200       fClsAmpDist->Fill(fClsAmpDriftDet[det]->GetMean());
201
202   PushBack((TObject*)fNClsDet, kAliHLTDataTypeHistogram, 0);   
203   PushBack((TObject*)fClsAmp, kAliHLTDataTypeHistogram, 0);  
204   PushBack((TObject*)fClsAmpDrift, kAliHLTDataTypeHistogram, 0);   
205   PushBack((TObject*)fClsTB, kAliHLTDataTypeHistogram, 0);  
206   PushBack((TObject*)fClsAmpDist, kAliHLTDataTypeHistogram, 0);  
207
208   //delete til dodeinit
209   // if(fPlotChargeOROCAll){
210   //  AliHLTUInt32_t fSpecification = AliHLTTRDDefinitions::EncodeDataSpecification(0,35,2,5);
211   //  PushBack( (TObject*) fTotalClusterChargeOROCAll,kAliHLTDataTypeHistogram,fSpecification);
212   //}
213   
214   return 0;
215 }
216
217 int AliHLTTRDClusterHistoComponent::Configure(const char* arguments)
218 {
219   return 0;
220 }
221
222 int AliHLTTRDClusterHistoComponent::Reconfigure(const char* cdbEntry, const char* chainId)
223 {
224   return 0;
225 }