]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCNoiseMapComponent.cxx
added new helper components to libAliHLTUtil (EsdCollector and AliHLTOUTPublisher...
[u/mrichter/AliRoot.git] / HLT / TPCLib / AliHLTTPCNoiseMapComponent.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   AliHLTTPCNoiseMapComponent.cxx
20     @author Kalliopi Kanaki
21     @date   
22     @brief  The TPC Noise Map component
23 */
24
25 #if __GNUC__>= 3
26 using namespace std;
27 #endif
28
29 #include "AliHLTTPCNoiseMapComponent.h"
30 #include "AliHLTTPCDigitReaderDecoder.h"
31 #include "AliHLTTPCDigitReaderPacked.h"
32 #include "AliHLTTPCTransform.h"
33 #include "AliHLTTPCDefinitions.h"
34
35 #include "AliCDBEntry.h"
36 #include "AliCDBManager.h"
37
38 #include "AliTPCCalPad.h"
39 #include "AliTPCROC.h"
40 #include "AliTPCCalROC.h"
41
42 #include <cstdlib>
43 #include <cerrno>
44 #include "TString.h"
45 #include "TFile.h"
46 #include "TObjArray.h"
47 #include "TObjString.h"
48 #include <sys/time.h>
49 #include "TH2.h"
50
51
52 AliHLTTPCNoiseMapComponent gAliHLTTPCNoiseMapComponent;
53
54 ClassImp(AliHLTTPCNoiseMapComponent) //ROOT macro for the implementation of ROOT specific class methods
55
56 AliHLTTPCNoiseMapComponent::AliHLTTPCNoiseMapComponent()
57     :    
58     fSpecification(0),
59     noisePad(NULL),
60     //pDigitReader(0),
61     fPlotSideA(0),
62     fPlotSideC(0),    
63     fApplyNoiseMap(0),
64     fResetHistograms(0),
65     fIsPacked(0),
66     fIsUnpacked(0),
67     fCurrentSlice(-99),
68     fCurrentPartition(-99),
69     fCurrentRow(-99),
70     fHistPartition(NULL),
71     fHistSideA(NULL),  
72     fHistSideC(NULL),
73     fHistCDBMap(NULL)    
74 {
75   // see header file for class documentation
76   // or
77   // refer to README to build package
78   // or
79   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
80 }
81
82 AliHLTTPCNoiseMapComponent::~AliHLTTPCNoiseMapComponent() { 
83 // see header file for class documentation
84
85 }
86
87 // Public functions to implement AliHLTComponent's interface.
88 // These functions are required for the registration process
89
90 const char* AliHLTTPCNoiseMapComponent::GetComponentID() { 
91 // see header file for class documentation
92
93   return "TPCNoiseMap";
94 }
95
96 void AliHLTTPCNoiseMapComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) { 
97 // see header file for class documentation
98
99   list.clear(); 
100   list.push_back( kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTPC );
101 }
102
103 AliHLTComponentDataType AliHLTTPCNoiseMapComponent::GetOutputDataType() { 
104 // see header file for class documentation
105
106   return kAliHLTDataTypeHistogram;
107 }
108
109 int AliHLTTPCNoiseMapComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList) { 
110 // see header file for class documentation
111
112   tgtList.clear();
113   tgtList.push_back(kAliHLTDataTypeHistogram);
114   return tgtList.size();
115 }
116
117 void AliHLTTPCNoiseMapComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) { 
118 // see header file for class documentation
119
120   constBase=800000;
121   inputMultiplier=0.0;
122 }
123
124 AliHLTComponent* AliHLTTPCNoiseMapComponent::Spawn() { 
125 // see header file for class documentation
126
127   return new AliHLTTPCNoiseMapComponent();
128 }
129         
130 int AliHLTTPCNoiseMapComponent::DoInit( int argc, const char** argv ) { 
131 // see header file for class documentation
132  
133   Int_t i = 0;
134   Char_t* cpErr;
135   
136   int iResult=0;
137   
138   TString configuration="";
139   TString argument="";
140   for (int j=0; j<argc && iResult>=0; j++) {
141     
142     argument=argv[j];
143     if (!configuration.IsNull()) configuration+=" ";
144     configuration+=argument;    
145   }
146    
147   if (!configuration.IsNull()) {
148     iResult=Configure(configuration.Data());
149   } else {
150     iResult=Reconfigure(NULL, NULL);
151   }
152
153  
154   while ( i < argc ) {      
155     if (!strcmp( argv[i], "-apply-noisemap")) {
156         fApplyNoiseMap = strtoul( argv[i+1], &cpErr ,0);
157             
158     if ( *cpErr ) {
159         HLTError("Cannot convert apply-noisemap specifier '%s'.", argv[i+1]);
160         return EINVAL;
161     }
162       i+=2;
163       continue;
164     }
165     
166     if (!strcmp( argv[i], "-plot-side-a")) {
167         fPlotSideA = strtoul( argv[i+1], &cpErr ,0);
168             
169     if ( *cpErr ) {
170         HLTError("Cannot convert plot-side-a specifier '%s'.", argv[i+1]);
171         return EINVAL;
172     }
173       i+=2;
174       continue;
175     }
176     
177     if (!strcmp( argv[i], "-plot-side-c")) {
178         fPlotSideC = strtoul( argv[i+1], &cpErr ,0);
179     
180     if ( *cpErr ) {
181         HLTError("Cannot convert plot-side-c specifier '%s'.", argv[i+1]);
182         return EINVAL;
183     }
184       i+=2;
185       continue;
186     }
187
188     if (!strcmp( argv[i], "-reset-histograms")) {
189         fResetHistograms = strtoul( argv[i+1], &cpErr ,0);
190     
191     if ( *cpErr ) {
192         HLTError("Cannot convert reset-histograms specifier '%s'.", argv[i+1]);
193         return EINVAL;
194     }
195       i+=2;
196       continue;
197     }
198                    
199     Logging(kHLTLogError, "HLT::TPCNoiseMap::DoInit", "Unknown Option", "Unknown option '%s'", argv[i] );
200     return EINVAL;
201
202   } // end while
203   
204
205 //   if(fApplyNoiseMap){
206 //     //TFile *f = TFile::Open("/scratch/noiseComp/Run3398_4000_v0_s72.root");
207 //     TFile *f = TFile::Open("/home/kanaki/noiseComp/Run3398_4000_v0_s72.root");
208 //     AliCDBEntry *pEntry = (AliCDBEntry*)f->Get("AliCDBEntry"); 
209 //     noisePad = (AliTPCCalPad*)pEntry->GetObject();
210 //     //fHistCDBMap = noisePad->MakeHisto2D(1); //side C
211 //   }
212    
213  
214 //   HLTDebug("using AliHLTTPCDigitReaderDecoder");
215 //   pDigitReader = new AliHLTTPCDigitReaderDecoder(); // double-loop
216 //   pDigitReader = new AliHLTTPCDigitReaderPacked();
217   
218   return 0;
219
220 } // end DoInit()
221
222 int AliHLTTPCNoiseMapComponent::DoDeinit() { 
223 // see header file for class documentation 
224       
225     return 0;
226 }
227
228 int AliHLTTPCNoiseMapComponent::DoEvent(const AliHLTComponentEventData& evtData, AliHLTComponentTriggerData& /*trigData*/){
229 // see header file for class documentation
230  
231   HLTInfo("--- Entering DoEvent() in TPCNoiseMap ---");
232  
233   if(GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR )) return 0;
234   
235   if(fPlotSideA){
236      fHistSideA = new TH2F("fHistSideA","TPC Side A",250,-250,250,250,-250,250);                
237      fHistSideA->SetXTitle("global X (cm)"); fHistSideA->SetYTitle("global Y (cm)");
238   }   
239   
240   if(fPlotSideC){    
241      fHistSideC = new TH2F("fHistSideC","TPC Side C",250,-250,250,250,-250,250);
242      fHistSideC->SetXTitle("global X (cm)"); fHistSideC->SetYTitle("global Y (cm)");
243   }
244   
245   const AliHLTComponentBlockData *iter = NULL;
246
247   Float_t xyz[3]; 
248   Int_t thissector, thisrow;
249  
250   fHistPartition = new TH2F("fHistPartition","fHistPartition",250,-250,250,250,-250,250);
251  
252   for(iter = GetFirstInputBlock(kAliHLTDataTypeDDLRaw|kAliHLTDataOriginTPC); iter != NULL; iter = GetNextInputBlock()){
253       
254      HLTInfo("Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s", 
255               evtData.fEventID, evtData.fEventID,
256               DataType2Text(iter->fDataType).c_str(), 
257               DataType2Text(kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTPC).c_str());
258
259      if (iter->fDataType == AliHLTTPCDefinitions::fgkDDLPackedRawDataType && GetEventCount()<2){
260          HLTWarning("data type %s is depricated, use %s (kAliHLTDataTypeDDLRaw)!", 
261          DataType2Text(AliHLTTPCDefinitions::fgkDDLPackedRawDataType).c_str(),
262          DataType2Text(kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTPC).c_str());
263      }      
264      
265      if (iter->fDataType != (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTPC)) continue;
266       
267      UInt_t slice     = AliHLTTPCDefinitions::GetMinSliceNr(*iter); 
268      UInt_t partition = AliHLTTPCDefinitions::GetMinPatchNr(*iter);
269      
270      fSpecification = iter->fSpecification;
271      
272      AliHLTTPCDigitReader *pDigitReader = new AliHLTTPCDigitReaderDecoder;
273
274      pDigitReader->InitBlock(iter->fPtr,iter->fSize,partition,slice);
275      if(!pDigitReader) break;
276        
277      //sprintf(name,"hMaxSignal_slice%d_partition%d", slice, partition);
278      //fHistPartition = new TH2F(name,name,250,-250,250,250,-250,250);
279             
280      while( pDigitReader->Next() ){ 
281      //while( pDigitReader->NextChannel()) { // pad loop 
282       
283       fCurrentRow  = pDigitReader->GetRow();  
284       fCurrentRow += pDigitReader->GetRowOffset();
285
286       AliHLTTPCTransform::Slice2Sector(slice,fCurrentRow,thissector,thisrow);
287       //AliHLTTPCTransform::Raw2Global(xyz,thissector,thisrow,pDigitReader->GetPad(),0);
288       
289 //       AliTPCCalROC *calRoc = noisePad->GetCalROC(thissector);
290 //       calRoc->GetValue(thisrow,pDigitReader->GetPad());
291       
292       Float_t maxSignal = 0.;
293       //while( pDigitReader->NextBunch()) {
294     
295       const UInt_t *bunchData = pDigitReader->GetSignals();
296       for(Int_t i=0;i<pDigitReader->GetBunchSize();i++) {
297           
298           if((Float_t)(bunchData[i])>maxSignal){
299               maxSignal = (Float_t)(bunchData[i]);
300           }
301           
302 //        if((Float_t)(bunchData[i])>maxSignal){          
303 //              if(fApplyNoiseMap) { //still in local coordinates
304 // 
305 //               if(calRoc->GetValue(thisrow,pDigitReader->GetPad())>0.) maxSignal = 0.;
306 //               else maxSignal = bunchData[i];
307 //               
308 //              } else maxSignal = bunchData[i];
309 //        } // end if
310       } // end for loop
311       //} // end of inner while loop
312             
313       AliHLTTPCTransform::Raw2Global(xyz,thissector,thisrow,pDigitReader->GetPad(),0); 
314       // transformation from pad-row coordinates to global ones
315       // time info is not taken into account
316            
317       fHistPartition->Fill(xyz[0],xyz[1],maxSignal);
318       
319       if(fPlotSideA || fPlotSideC){
320          if(slice<18) fHistSideA->Fill(xyz[0],xyz[1],maxSignal);
321          else         fHistSideC->Fill(xyz[0],xyz[1],maxSignal);                             
322       } // end if plotting sides            
323      } // end of while loop  
324      delete pDigitReader;
325   } // end of for loop over data blocks
326  
327   if(fResetHistograms) ResetHistograms();
328   MakeHistosPublic();
329   
330   return 0;
331 } // end DoEvent()
332
333 void AliHLTTPCNoiseMapComponent::MakeHistosPublic() {
334 // see header file for class documentation
335   
336 //   TFile *outputfile = new TFile("test.root","RECREATE");
337 //   fHistSideC->Write();
338 //   fHistPartition->Write();
339 //   fHistCDBMap->Write();
340 //   outputfile->Save();
341 //   outputfile->Close();
342   
343   TObjArray histos;
344   histos.Add(fHistPartition);
345   if(fPlotSideA) histos.Add(fHistSideA);
346   if(fPlotSideC) histos.Add(fHistSideC);
347   //histos.Add(fHistCDBMap);
348   
349   TIter iterator(&histos);
350   while(TObject *pObj=iterator.Next()){
351   
352         PushBack(pObj, kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC, fSpecification);
353   }
354     
355   //PushBack( (TObject*) &histos, kAliHLTDataTypeHistogram, fSpecification);    
356  
357   delete fHistPartition;
358   if(fHistSideA) delete fHistSideA; fHistSideA=NULL;
359   if(fHistSideC) delete fHistSideC; fHistSideC=NULL;
360 }
361
362 void AliHLTTPCNoiseMapComponent::ResetHistograms(){
363 // see header file for class documentation
364
365   fHistPartition->Reset();
366   if(fHistSideA) fHistSideA->Reset();
367   if(fHistSideC) fHistSideC->Reset();
368 }
369
370 int AliHLTTPCNoiseMapComponent::Configure(const char* arguments) { 
371 // see header file for class documentation
372   
373   int iResult=0;
374   if (!arguments) return iResult;
375   HLTInfo("parsing configuration string \'%s\'", arguments);
376
377   TString allArgs=arguments;
378   TString argument;
379   int bMissingParam=0;
380
381   TObjArray* pTokens=allArgs.Tokenize(" ");
382   if (pTokens) {
383     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
384       argument=((TObjString*)pTokens->At(i))->GetString();
385       if (argument.IsNull()) continue;
386      
387       if (argument.CompareTo("-apply-noisemap")==0) {
388         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
389         HLTInfo("got \'-apply-noisemap\': %s", ((TObjString*)pTokens->At(i))->GetString().Data());
390         
391       } 
392       else if (argument.CompareTo("-plot-side-c")==0) {
393         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
394         HLTInfo("got \'-plot-side-c\': %s", ((TObjString*)pTokens->At(i))->GetString().Data());
395         
396       } 
397       else if (argument.CompareTo("-plot-side-a")==0) {
398         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
399         HLTInfo("got \'-plot-side-a\': %s", ((TObjString*)pTokens->At(i))->GetString().Data());
400         
401       } 
402       else {
403         HLTError("unknown argument %s", argument.Data());
404         iResult=-EINVAL;
405         break;
406       }
407     } // end for
408   
409     delete pTokens;
410   
411   } // end if pTokens
412   
413   if (bMissingParam) {
414     HLTError("missing parameter for argument %s", argument.Data());
415     iResult=-EINVAL;
416   }
417   return iResult;  
418 }
419
420 int AliHLTTPCNoiseMapComponent::Reconfigure(const char* cdbEntry, const char* chainId) { 
421 // see header file for class documentation
422   int iResult=0;
423   const char* path="HLT/ConfigTPC/TPCNoiseMapComponent";
424   const char* defaultNotify="";
425   if (cdbEntry) {
426       path=cdbEntry;
427       defaultNotify=" (default)";
428   }
429   
430   if (path) {
431     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
432     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
433     if (pEntry) {
434       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
435       if (pString) {
436         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
437         iResult=Configure(pString->GetString().Data());
438       } else {
439         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
440       }
441     } else {
442       HLTError("cannot fetch object \"%s\" from CDB", path);
443     }
444   }
445   
446   const char* pathNoiseMap="TPC/Config/NoiseMap";
447
448   if (pathNoiseMap) {
449     HLTInfo("reconfigure noise map from entry %s, chain id %s", path,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
450     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(pathNoiseMap/*,GetRunNo()*/);
451     if (pEntry) {
452       TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
453       if (pString) {
454         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
455         iResult=Configure(pString->GetString().Data());
456       } else {
457         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
458       }
459     } else {
460       HLTError("cannot fetch object \"%s\" from CDB", path);
461     }
462   }
463
464   return iResult;
465 }