]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCHistogramHandlerComponent.cxx
added new helper components to libAliHLTUtil (EsdCollector and AliHLTOUTPublisher...
[u/mrichter/AliRoot.git] / HLT / TPCLib / AliHLTTPCHistogramHandlerComponent.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        *
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Kalliopi Kanaki <Kalliopi.Kanaki@ift.uib.no>          *
8 //*                  for The ALICE HLT Project.                            *
9 //*                                                                        *
10 //* Permission to use, copy, modify and distribute this software and its   *
11 //* documentation strictly for non-commercial purposes is hereby granted   *
12 //* without fee, provided that the above copyright notice appears in all   *
13 //* copies and that both the copyright notice and this permission notice   *
14 //* appear in the supporting documentation. The authors make no claims     *
15 //* about the suitability of this software for any purpose. It is          *
16 //* provided "as is" without express or implied warranty.                  *
17 //**************************************************************************
18
19 /** @file   AliHLTTPCHistogramHandlerComponent.cxx
20     @author Kalliopi Kanaki
21     @date   
22     @brief  The Histogram Handler component
23 */
24
25 #if __GNUC__>= 3
26 using namespace std;
27 #endif
28 #include "AliHLTTPCHistogramHandlerComponent.h"
29 #include "AliHLTTPCDefinitions.h"
30 #include "AliCDBEntry.h"
31 #include "AliCDBManager.h"
32 #include "AliHLTTPCTransform.h"
33
34 #include <cstdlib>
35 #include <cerrno>
36 #include "TString.h"
37 #include "TFile.h"
38 #include "TObjArray.h"
39 #include "TObjString.h"
40 #include <sys/time.h>
41 #include "TH1.h"
42 #include "TH2.h"
43 #include "TLine.h"
44 #include "TMath.h"
45
46
47 AliHLTTPCHistogramHandlerComponent gAliHLTTPCHistogramHandlerComponent;
48
49 ClassImp(AliHLTTPCHistogramHandlerComponent) //ROOT macro for the implementation of ROOT specific class methods
50
51 AliHLTTPCHistogramHandlerComponent::AliHLTTPCHistogramHandlerComponent()
52     :    
53     fSpecification(0),
54     fNoiseHistograms(0),
55     fKryptonHistograms(0),
56     fSlice(-99),
57     
58     fHistTH1Tmp(NULL),
59     fTotalClusterChargeIROCAll(NULL),
60     fTotalClusterChargeOROCAll(NULL),
61     fQMaxPartitionAll(NULL),
62     fPlotQmaxROCAll(NULL),
63     fNumberOfClusters(NULL),
64     
65     fHistTH2Tmp(NULL),
66     fHistTPCSideA(NULL),  
67     fHistTPCSideC(NULL)    
68 {
69   // see header file for class documentation
70   // or
71   // refer to README to build package
72   // or
73   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
74 }
75
76 AliHLTTPCHistogramHandlerComponent::~AliHLTTPCHistogramHandlerComponent() { 
77 // see header file for class documentation
78
79 }
80
81 // Public functions to implement AliHLTComponent's interface.
82 // These functions are required for the registration process
83
84 const char* AliHLTTPCHistogramHandlerComponent::GetComponentID() { 
85 // see header file for class documentation
86
87   return "TPCHistogramHandler";
88 }
89
90 void AliHLTTPCHistogramHandlerComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) { 
91 // see header file for class documentation
92
93   list.clear(); 
94   list.push_back( kAliHLTDataTypeHistogram );
95 }
96
97 AliHLTComponentDataType AliHLTTPCHistogramHandlerComponent::GetOutputDataType() { 
98 // see header file for class documentation
99
100   return kAliHLTDataTypeHistogram;
101 }
102
103 int AliHLTTPCHistogramHandlerComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList) { 
104 // see header file for class documentation
105
106   tgtList.clear();
107   tgtList.push_back(kAliHLTDataTypeHistogram);
108   return tgtList.size();
109 }
110
111 void AliHLTTPCHistogramHandlerComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) { 
112 // see header file for class documentation
113
114   constBase=0;
115   inputMultiplier=2.0;
116 }
117
118 AliHLTComponent* AliHLTTPCHistogramHandlerComponent::Spawn() { 
119 // see header file for class documentation
120
121   return new AliHLTTPCHistogramHandlerComponent();
122 }
123         
124 int AliHLTTPCHistogramHandlerComponent::DoInit( int argc, const char** argv ) { 
125 // see header file for class documentation
126
127   Int_t i = 0;
128   Char_t* cpErr;
129   
130   int iResult=0;
131   
132   TString configuration="";
133   TString argument="";
134   for (int j=0; j<argc && iResult>=0; j++) {
135     
136     argument=argv[j];
137     if (!configuration.IsNull()) configuration+=" ";
138     configuration+=argument;    
139   }
140    
141   if (!configuration.IsNull()) {
142     iResult=Configure(configuration.Data());
143   } else {
144     iResult=Reconfigure(NULL, NULL);
145   }
146
147  
148   while ( i < argc ) {      
149     if (!strcmp( argv[i], "-sum-noise-histograms")) {
150         fNoiseHistograms = strtoul( argv[i+1], &cpErr ,0);
151             
152     if ( *cpErr ) {
153         HLTError("Cannot convert sum-noise-histograms specifier '%s'.", argv[i+1]);
154         return EINVAL;
155     }
156       i+=2;
157       continue;
158     }
159     
160     if (!strcmp( argv[i], "-sum-krypton-histograms")) {
161         fKryptonHistograms = strtoul( argv[i+1], &cpErr ,0);
162             
163     if ( *cpErr ) {
164         HLTError("Cannot convert sum-krypton-histograms specifier '%s'.", argv[i+1]);
165         return EINVAL;
166     }
167       i+=2;
168       continue;
169     }    
170                   
171     Logging(kHLTLogError, "HLT::TPCHistogramHandler::DoInit", "Unknown Option", "Unknown option '%s'", argv[i] );
172     return EINVAL;
173
174   } // end while
175   
176   return 0;
177 } // end DoInit()
178
179 int AliHLTTPCHistogramHandlerComponent::DoDeinit() { 
180 // see header file for class documentation 
181    
182    return 0;
183 }
184
185 int AliHLTTPCHistogramHandlerComponent::DoEvent(const AliHLTComponentEventData&/* evtData*/, AliHLTComponentTriggerData& /*trigData*/){
186 // see header file for class documentation
187
188   HLTInfo("--- Entering DoEvent() in TPCHistogramHandler ---");
189   
190   if(GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR )) return 0;
191  
192   fTotalClusterChargeIROCAll = new TH1F("fTotalClusterChargeIROCAll","Total Charge of clusters in all IROC",4000,0,4000);  
193   fTotalClusterChargeOROCAll = new TH1F("fTotalClusterChargeOROCAll","Total Charge of clusters in all OROC",4000,0,4000);
194   fQMaxPartitionAll          = new TH1F("fQMaxPartitionAll",         "QMax for All Partitions",             216,0,216);
195   fPlotQmaxROCAll            = new TH1F("fQMaxROCAll",               "QMax for All ROC",                    72,0,72);
196   fNumberOfClusters          = new TH1F("fNumberOfClusters",         "Total Number of Clusters",            1,0,1);
197     
198   fHistTH2Tmp = new TH2F("fHistTH2Tmp","fHistTH2Tmp",250,-250,250,250,-250,250);    
199   fHistTPCSideA = new TH2F("fHistTPCSideA","TPC side A (max signal)",250,-250,250,250,-250,250);
200   fHistTPCSideA->SetXTitle("global X (cm)"); fHistTPCSideA->SetYTitle("global Y (cm)");
201   fHistTPCSideC = new TH2F("fHistTPCSideC","TPC side C (max signal)",250,-250,250,250,-250,250);
202   fHistTPCSideC->SetXTitle("global X (cm)"); fHistTPCSideC->SetYTitle("global Y (cm)");
203   
204   const TObject *iter = NULL;  
205         
206   for(iter = GetFirstInputObject(kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC); iter != NULL; iter = GetNextInputObject()){
207    
208 //      HLTInfo("Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s", 
209 //               evtData.fEventID, evtData.fEventID,
210 //               DataType2Text(GetDataType(iter)).c_str(), 
211 //               DataType2Text(kAliHLTDataTypeHistogram | kAliHLTDataOriginTPC).c_str());
212    
213 //      if (GetDataType(iter) == (kAliHLTDataTypeHistogram | kAliHLTDataOriginTPC) && GetEventCount()<2){
214 //          HLTWarning("data type %s is depricated, use %s (kAliHLTDataTypeHistogram)!", 
215 //       DataType2Text(kAliHLTDataTypeHistogram).c_str(),
216 //          DataType2Text(kAliHLTDataTypeHistogram | kAliHLTDataOriginTPC).c_str());
217 //      }      
218      
219      if (GetDataType(iter) != (kAliHLTDataTypeHistogram | kAliHLTDataOriginTPC)) continue;
220          
221      // Summing the output histograms of the AliHLTTPCNoiseMapComponent (from partition to TPC sides)
222      if(fNoiseHistograms){  
223        
224         fHistTH2Tmp = (TH2F*)iter;
225         UInt_t minSlice     = AliHLTTPCDefinitions::GetMinSliceNr(GetSpecification(iter)); 
226         UInt_t maxSlice     = AliHLTTPCDefinitions::GetMaxSliceNr(GetSpecification(iter)); 
227         UInt_t minPartition = AliHLTTPCDefinitions::GetMinPatchNr(GetSpecification(iter)); 
228         UInt_t maxPartition = AliHLTTPCDefinitions::GetMaxPatchNr(GetSpecification(iter)); 
229        
230         if((minSlice!=maxSlice) || (minPartition!=maxPartition)){
231             HLTWarning("TPCHistogramHandler::The Noise Map component is not running on partition level!");
232         }
233
234         // minSlice=maxSlice, when the Noise Map component runs on partition level (as it should)
235         if(minSlice<18) fHistTPCSideA->Add(fHistTPCSideA,fHistTH2Tmp,1,1);
236         else            fHistTPCSideC->Add(fHistTPCSideC,fHistTH2Tmp,1,1);
237      } // endif fNoiseHistograms==kTRUE   
238      
239      
240      // Summing the output of AliHLTTPCClusterHistoComponent
241      if(fKryptonHistograms){
242        Int_t thisrow=-1,thissector=-1,row=-1;
243        
244        AliHLTUInt8_t slice = AliHLTTPCDefinitions::GetMinSliceNr(GetSpecification(iter));
245        AliHLTUInt8_t patch = AliHLTTPCDefinitions::GetMinPatchNr(GetSpecification(iter));
246        row = AliHLTTPCTransform::GetFirstRow(patch); 
247        AliHLTTPCTransform::Slice2Sector(slice,row,thissector,thisrow);
248        
249        fHistTH1Tmp = (TH1F*)iter;       
250        //cout << fHistTH1Tmp->GetName() << "\t" << fHistTH1Tmp->GetEntries() << endl;
251        
252         TString name = fHistTH1Tmp->GetName();
253                 
254         if(name=="fTotalClusterChargeIROCAll"){
255            fTotalClusterChargeIROCAll->Add(fTotalClusterChargeIROCAll,fHistTH1Tmp,1,1);
256         } 
257         else if(name=="fTotalClusterChargeOROCAll"){
258            fTotalClusterChargeOROCAll->Add(fTotalClusterChargeOROCAll,fHistTH1Tmp,1,1);
259         } 
260         else if(name=="fQMaxPartitionAll"){
261           for(Int_t t=0;t<216;t++){
262             if(fHistTH1Tmp->GetBinContent(t)>fQMaxPartitionAll->GetBinContent(t)){
263               fQMaxPartitionAll->SetBinContent(t,fHistTH1Tmp->GetBinContent(t));
264             }
265           } 
266         } 
267         else if(name=="fQMaxROCAll"){
268           for(Int_t t=0;t<72;t++){
269             if(fHistTH1Tmp->GetBinContent(t)>fPlotQmaxROCAll->GetBinContent(t)){
270               fPlotQmaxROCAll->SetBinContent(t,fHistTH1Tmp->GetBinContent(t));
271             }
272           }
273         } 
274         else if(name=="fNumberOfClusters"){ 
275           fNumberOfClusters->Add(fNumberOfClusters,fHistTH1Tmp,1,1);
276         } 
277         else{
278           HLTWarning("No histogram names match. %s",name.Data());
279           continue;
280         }     
281      } //endif fKryptonHistograms==kTRUE                         
282   } // end for loop over histogram blocks
283   
284   MakeHistosPublic();
285   return 0;
286 } // end DoEvent()
287
288 void AliHLTTPCHistogramHandlerComponent::MakeHistosPublic() {
289 // see header file for class documentation
290  
291   if(fNoiseHistograms){ 
292     PushBack((TObject*)fHistTPCSideA,kAliHLTDataTypeHistogram,AliHLTTPCDefinitions::EncodeDataSpecification( 0,17,0,5));
293     PushBack((TObject*)fHistTPCSideC,kAliHLTDataTypeHistogram,AliHLTTPCDefinitions::EncodeDataSpecification(18,35,0,5));
294
295     if(fHistTH2Tmp){
296       delete fHistTH2Tmp;
297     }
298     if(fHistTPCSideA){
299       delete fHistTPCSideA;
300     }
301     if(fHistTPCSideC){
302       delete fHistTPCSideC;
303     }
304     
305     fHistTH2Tmp=NULL;
306     fHistTPCSideA=NULL;
307     fHistTPCSideC=NULL;
308   }  
309   
310   if(fKryptonHistograms){
311      PushBack((TObject*)fTotalClusterChargeIROCAll,kAliHLTDataTypeHistogram,AliHLTTPCDefinitions::EncodeDataSpecification(0,17,0,1));
312      PushBack((TObject*)fTotalClusterChargeOROCAll,kAliHLTDataTypeHistogram,AliHLTTPCDefinitions::EncodeDataSpecification(0,17,2,5));
313      PushBack((TObject*)fQMaxPartitionAll,         kAliHLTDataTypeHistogram,AliHLTTPCDefinitions::EncodeDataSpecification(0,35,0,5));
314      PushBack((TObject*)fPlotQmaxROCAll,           kAliHLTDataTypeHistogram,AliHLTTPCDefinitions::EncodeDataSpecification(0,35,0,5));
315      PushBack((TObject*)fNumberOfClusters,         kAliHLTDataTypeHistogram,AliHLTTPCDefinitions::EncodeDataSpecification(0,35,0,5));
316           
317      if(fTotalClusterChargeIROCAll){
318        delete fTotalClusterChargeIROCAll;
319        fTotalClusterChargeIROCAll=NULL;
320      }
321      if(fTotalClusterChargeOROCAll){
322        delete fTotalClusterChargeOROCAll;
323        fTotalClusterChargeOROCAll=NULL;
324      }
325      if(fQMaxPartitionAll){
326        delete fQMaxPartitionAll;
327        fQMaxPartitionAll=NULL;
328      }
329      if(fPlotQmaxROCAll){
330        delete fPlotQmaxROCAll;
331        fPlotQmaxROCAll=NULL;
332      }
333      if(fNumberOfClusters){
334        delete fNumberOfClusters;
335        fNumberOfClusters=NULL;
336      }
337      if(fHistTH1Tmp){
338        delete fHistTH1Tmp;
339        fHistTH1Tmp=NULL;
340      }
341   }
342  
343 //  TObjArray histos;
344    
345 //   if(fPlotSideA) histos.Add(fHistSideA);
346 //   if(fPlotSideC) histos.Add(fHistSideC);
347 //   if(fApplyNoiseMap) histos.Add(fHistCDBMap);
348 //   
349 //   TIter iterator(&histos);
350 //   while(TObject *pObj=iterator.Next()){
351 //   
352 //         PushBack(pObj, kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC, fSpecification);
353 //   }
354 //   
355 //   
356 //   //PushBack( (TObject*) &histos, kAliHLTDataTypeHistogram, fSpecification);    
357 }
358
359 int AliHLTTPCHistogramHandlerComponent::Configure(const char* arguments) { 
360 // see header file for class documentation
361   
362   int iResult=0;
363   if (!arguments) return iResult;
364   HLTInfo("parsing configuration string \'%s\'", arguments);
365
366   TString allArgs=arguments;
367   TString argument;
368   int bMissingParam=0;
369
370   TObjArray* pTokens=allArgs.Tokenize(" ");
371   if (pTokens) {
372     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
373       argument=((TObjString*)pTokens->At(i))->GetString();
374       if (argument.IsNull()) continue;
375      
376       if (argument.CompareTo("-sum-noise-histograms")==0) {
377         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
378         HLTInfo("got \'-sum-noise-histograms\': %s", ((TObjString*)pTokens->At(i))->GetString().Data());
379         
380       } else if (argument.CompareTo("-sum-krypton-histograms")==0) {
381         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
382         HLTInfo("got \'-sum-krypton-histograms\': %s", ((TObjString*)pTokens->At(i))->GetString().Data());
383         
384       } 
385       else {
386         HLTError("unknown argument %s", argument.Data());
387         iResult=-EINVAL;
388         break;
389       }
390     } // end for
391   
392     delete pTokens;
393   
394   } // end if pTokens
395   
396   if (bMissingParam) {
397     HLTError("missing parameter for argument %s", argument.Data());
398     iResult=-EINVAL;
399   }
400   return iResult;
401 }
402
403 int AliHLTTPCHistogramHandlerComponent::Reconfigure(const char* cdbEntry, const char* chainId) { 
404 // see header file for class documentation
405
406   int iResult=0;
407   const char* path="HLT/ConfigTPC/TPCHistogramHandlerComponent";
408   const char* defaultNotify="";
409   if (cdbEntry) {
410       path=cdbEntry;
411       defaultNotify=" (default)";
412   }
413   
414   if (path) {
415     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
416     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
417     if (pEntry) {
418       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
419       if (pString) {
420         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
421         iResult=Configure(pString->GetString().Data());
422       } else {
423         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
424       }
425     } else {
426       HLTError("cannot fetch object \"%s\" from CDB", path);
427     }
428   }
429   return iResult;
430 }