]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDClusterHistoComponent.cxx
unify system of the input multiplier for clusterizer and tracker and
[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 <time.h>
27
28 #include "AliHLTTRDClusterHistoComponent.h"
29 #include "AliHLTTRDDefinitions.h"
30 #include "AliTRDcluster.h"
31 #include "AliCDBEntry.h"
32 #include "AliCDBManager.h"
33 #include <TFile.h>
34 #include <TString.h>
35 #include "TObjString.h"
36 #include "TClonesArray.h"
37 #include "TTimeStamp.h"
38 #include "AliHLTTRDUtils.h"
39 #include "TH2F.h"
40
41 //#include "AliHLTTRD.h"
42 //#include <stdlib.h>
43 //#include <cerrno>
44
45 /** ROOT macro for the implementation of ROOT specific class methods */
46 ClassImp(AliHLTTRDClusterHistoComponent)
47
48 AliHLTTRDClusterHistoComponent::AliHLTTRDClusterHistoComponent()
49 : AliHLTProcessor(),
50   fOutputSize(100000),
51   fClusterArray(NULL),
52   fNClsDet(NULL),
53   fClsAmp(NULL),
54   fClsAmpDrift(NULL),
55   fClsTB(NULL),
56   fClsAmpDist(NULL),
57   fSClsDist(NULL),
58   fNScls(NULL)
59 {
60   // see header file for class documentation
61   // or
62   // refer to README to build package
63   // or
64   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
65
66 }
67
68 AliHLTTRDClusterHistoComponent::~AliHLTTRDClusterHistoComponent()
69 {
70   // see header file for class documentation
71 }
72
73 // Public functions to implement AliHLTComponent's interface.
74 // These functions are required for the registration process
75
76 const char* AliHLTTRDClusterHistoComponent::GetComponentID()
77 {
78   // see header file for class documentation
79   
80   return "TRDClusterHisto";
81 }
82
83 void AliHLTTRDClusterHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
84 {
85   // see header file for class documentation
86   list.clear();
87   list.push_back( AliHLTTRDDefinitions::fgkClusterDataType );
88 }
89
90 AliHLTComponentDataType AliHLTTRDClusterHistoComponent::GetOutputDataType()
91 {
92   // see header file for class documentation
93   return kAliHLTDataTypeHistogram  | kAliHLTDataOriginTRD;
94
95 }
96
97 void AliHLTTRDClusterHistoComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
98 {
99   // see header file for class documentation
100   constBase = fOutputSize;
101   inputMultiplier = 0;
102 }
103
104 AliHLTComponent* AliHLTTRDClusterHistoComponent::Spawn()
105 {
106   // see header file for class documentation
107   return new AliHLTTRDClusterHistoComponent;
108 }
109
110 int AliHLTTRDClusterHistoComponent::DoInit(int argc, const char** argv)
111 {
112   // Initialize histograms
113   int iResult=0;
114   
115   TString configuration="";
116   TString argument="";
117   for (int i=0; i<argc && iResult>=0; i++) {
118     argument=argv[i];
119     if (!configuration.IsNull()) configuration+=" ";
120     configuration+=argument;
121   }
122
123   if (!configuration.IsNull()) {
124     iResult=Configure(configuration.Data());
125   } 
126
127   fClusterArray = new TClonesArray("AliTRDcluster");
128
129   fNClsDet = new TH1D("trdClsDet", ";detector", 540, -0.5, 539.5);
130   fClsAmp  = new TH1D("trdClsAmp", ";amplitude", 200, -0.5, 1999.5);
131   fClsAmpDrift = new TH1D("trdClsAmpDrift", ";amplitude", 200, -0.5, 199.5) ;
132   fClsTB = new TH1D("trdClsTB", ";time bin", 35, -0.5, 34.5);
133   fClsAmpDist = new TH1D("trdClsAmpDist", "mean amplitude", 200, 0, 1000);
134   fSClsDist = new TH1D("sclsdist", "Super cluster spectrum", 200, 0, 8000);
135   fNScls = new TH1D("nscls", "No. of Kr clusters per event", 540, 0, 540);
136
137   for(int i=0; i<540; i++) {
138     fClsAmpDriftDet[i] = new TH1D(Form("trdClsDriftDet_%d",i), "", 200, -0.5, 199.5);
139   }
140   
141   return 0;
142 }
143   
144 int AliHLTTRDClusterHistoComponent::DoDeinit()
145 {
146   // see header file for class documentation
147
148   fClusterArray->Delete();
149   delete fClusterArray;
150
151   // delete histograms
152   if (fNClsDet) delete fNClsDet;
153   if (fClsAmp) delete fClsAmp;
154   if (fClsAmpDrift) delete fClsAmpDrift;
155   if (fClsTB) delete fClsTB;
156   if (fClsAmpDist) delete fClsAmpDist;
157   if (fSClsDist) delete fSClsDist;
158   if (fNScls) delete fNScls;
159
160   for(int i=0; i<540; i++)
161     if (fClsAmpDriftDet[i]) delete fClsAmpDriftDet[i];
162
163   return 0;
164 }
165
166 int AliHLTTRDClusterHistoComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/,
167                                             AliHLTComponentTriggerData& /*trigData*/)
168 {
169
170   // if (GetFirstInputBlock(kAliHLTDataTypeSOR)) return 0;
171   // else if (GetFirstInputBlock(kAliHLTDataTypeEOR))
172   //   {
173   //     TString fileName="/tmp/ClusterHistoDump_run";
174   //     fileName+=AliCDBManager::Instance()->GetRun();
175   //     fileName+=".root";
176   //     HLTInfo("Dumping Histogram file to %s",fileName.Data());
177   //     TFile* file = TFile::Open(fileName, "RECREATE");
178   //     fNClsDet->Write();
179   //     fClsAmp->Write();
180   //     fClsAmpDrift->Write();
181   //     fClsTB->Write();
182   //     fClsAmpDist->Write(); 
183   //     fSClsDist->Write();
184   //     fNScls->Write();
185   //     file->Close();
186   //     HLTInfo("Histogram file dumped");
187   //     return 0;
188   //   }
189
190   if (GetFirstInputBlock(kAliHLTDataTypeSOR) || GetFirstInputBlock(kAliHLTDataTypeEOR)) return 0;
191
192   const AliHLTComponentBlockData* iter = NULL;
193   
194   Float_t sClusterCharge[540] = { 0 };
195
196   for ( iter = GetFirstInputBlock(AliHLTTRDDefinitions::fgkClusterDataType); 
197         iter != NULL; iter = GetNextInputBlock() ) {
198
199     HLTDebug("We get the right data type: Block Ptr: 0x%x; Block Size: %i",
200              iter->fPtr, iter->fSize);
201
202     AliHLTTRDUtils::ReadClusters(fClusterArray, iter->fPtr, iter->fSize);
203     HLTDebug("TClonesArray of clusters: nbEntries = %i", fClusterArray->GetEntriesFast());
204
205     AliTRDcluster *cls;
206
207     // loop over clusters 
208     for(int i=0;i<fClusterArray->GetEntriesFast();i++) {
209
210       cls=(AliTRDcluster*)fClusterArray->At(i);
211       
212       fNClsDet->Fill(cls->GetDetector());
213       fClsAmp->Fill(cls->GetQ());
214       
215       int tb = cls->GetPadTime();
216       fClsTB->Fill(tb);
217       if (tb > 5 && tb <25)
218         fClsAmpDrift->Fill(cls->GetQ()); 
219       
220       fClsAmpDriftDet[cls->GetDetector()]->Fill(cls->GetQ());
221
222       Int_t det = cls->GetDetector();
223       sClusterCharge[det] += cls->GetQ();
224
225     }
226     
227     fClusterArray->Delete();
228     
229   }
230    
231   fClsAmpDist->Reset();
232   Int_t nSClusters = 0;
233   for(int det=0; det<540; det++) {
234     if (fClsAmpDriftDet[det]->GetSum() > 0) 
235       fClsAmpDist->Fill(fClsAmpDriftDet[det]->GetMean());
236     fSClsDist->Fill(sClusterCharge[det]);
237
238   }
239
240   fNScls->Fill(nSClusters);
241
242   PushBack((TObject*)fNClsDet, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);   
243   PushBack((TObject*)fClsAmp, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);  
244   //  PushBack((TObject*)fClsAmpDrift, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);   
245   //  PushBack((TObject*)fClsTB, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);  
246   //  PushBack((TObject*)fClsAmpDist, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);  
247   PushBack((TObject*)fNScls, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);  
248   PushBack((TObject*)fSClsDist, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);
249
250   return 0;
251 }
252
253 int AliHLTTRDClusterHistoComponent::Configure(const char* arguments){
254   int iResult=0;
255   if (!arguments) return iResult;
256   
257   TString allArgs=arguments;
258   TString argument;
259   int bMissingParam=0;
260
261   TObjArray* pTokens=allArgs.Tokenize(" ");
262   if (pTokens) {
263     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
264       argument=((TObjString*)pTokens->At(i))->GetString();
265       if (argument.IsNull()) continue;
266       
267       if (argument.CompareTo("output_size")==0) {
268         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
269         HLTInfo("Setting output size to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
270         fOutputSize=((TObjString*)pTokens->At(i))->GetString().Atoi();
271         continue;
272       } 
273       if (argument.CompareTo("-everyNevent")==0) {
274         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
275         HLTInfo("Option -everyNevent depreceated");
276         continue;
277       } 
278       else {
279         HLTError("unknown argument: %s", argument.Data());
280         iResult=-EINVAL;
281         break;
282       }
283     }
284     delete pTokens;
285   }
286   if (bMissingParam) {
287     HLTError("missing parameter for argument %s", argument.Data());
288     iResult=-EINVAL;
289   }
290   return iResult;
291 }