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