3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project *
5 //* ALICE Experiment at CERN, All rights reserved. *
7 //* Primary Authors: Kalliopi Kanaki <Kalliopi.Kanaki@ift.uib.no> *
8 //* for The ALICE HLT Project. *
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 //**************************************************************************
19 /** @file AliHLTTPCNoiseMapComponent.cxx
20 @author Kalliopi Kanaki
22 @brief The TPC Noise Map component
29 #include "AliHLTTPCNoiseMapComponent.h"
30 #include "AliHLTTPCDigitReaderDecoder.h"
31 //#include "AliHLTTPCDigitReaderPacked.h"
32 #include "AliHLTTPCTransform.h"
33 #include "AliHLTTPCDefinitions.h"
35 #include "AliCDBEntry.h"
36 #include "AliCDBManager.h"
37 #include "AliHLTTPCNoiseMap.h"
39 #include "AliTPCCalPad.h"
40 #include "AliTPCROC.h"
41 #include "AliTPCCalROC.h"
47 #include "TObjArray.h"
48 #include "TObjString.h"
52 ClassImp(AliHLTTPCNoiseMapComponent) //ROOT macro for the implementation of ROOT specific class methods
54 AliHLTTPCNoiseMapComponent::AliHLTTPCNoiseMapComponent()
64 fCurrentPartition(-99),
74 // see header file for class documentation
76 // refer to README to build package
78 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
81 AliHLTTPCNoiseMapComponent::~AliHLTTPCNoiseMapComponent() {
82 // see header file for class documentation
86 // Public functions to implement AliHLTComponent's interface.
87 // These functions are required for the registration process
89 const char* AliHLTTPCNoiseMapComponent::GetComponentID() {
90 // see header file for class documentation
95 void AliHLTTPCNoiseMapComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) {
96 // see header file for class documentation
99 list.push_back( kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTPC );
102 AliHLTComponentDataType AliHLTTPCNoiseMapComponent::GetOutputDataType() {
103 // see header file for class documentation
105 return kAliHLTDataTypeHistogram;
108 int AliHLTTPCNoiseMapComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList) {
109 // see header file for class documentation
112 tgtList.push_back(kAliHLTDataTypeHistogram);
113 return tgtList.size();
116 void AliHLTTPCNoiseMapComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) {
117 // see header file for class documentation
123 AliHLTComponent* AliHLTTPCNoiseMapComponent::Spawn() {
124 // see header file for class documentation
126 return new AliHLTTPCNoiseMapComponent();
129 int AliHLTTPCNoiseMapComponent::DoInit( int argc, const char** argv ) {
130 // see header file for class documentation
137 TString configuration="";
139 for (int j=0; j<argc && iResult>=0; j++) {
142 if (!configuration.IsNull()) configuration+=" ";
143 configuration+=argument;
146 if (!configuration.IsNull()) {
147 iResult=Configure(configuration.Data());
149 iResult=Reconfigure(NULL, NULL);
154 if (!strcmp( argv[i], "-apply-noisemap")) {
155 fApplyNoiseMap = strtoul( argv[i+1], &cpErr ,0);
158 HLTError("Cannot convert apply-noisemap specifier '%s'.", argv[i+1]);
165 if (!strcmp( argv[i], "-plot-side-a")) {
166 fPlotSideA = strtoul( argv[i+1], &cpErr ,0);
169 HLTError("Cannot convert plot-side-a specifier '%s'.", argv[i+1]);
176 if (!strcmp( argv[i], "-plot-side-c")) {
177 fPlotSideC = strtoul( argv[i+1], &cpErr ,0);
180 HLTError("Cannot convert plot-side-c specifier '%s'.", argv[i+1]);
187 if (!strcmp( argv[i], "-reset-histograms")) {
188 fResetHistograms = strtoul( argv[i+1], &cpErr ,0);
191 HLTError("Cannot convert reset-histograms specifier '%s'.", argv[i+1]);
198 Logging(kHLTLogError, "HLT::TPCNoiseMap::DoInit", "Unknown Option", "Unknown option '%s'", argv[i] );
204 AliHLTTPCNoiseMap *nm = AliHLTTPCNoiseMap::Instance();
205 AliTPCCalPad *noisePad = nm->ReadNoiseMap();
206 fHistCDBMap = noisePad->MakeHisto2D(1);
209 // if(fApplyNoiseMap){
210 // //TFile *f = TFile::Open("/scratch/noiseComp/Run3398_4000_v0_s72.root");
211 // TFile *f = TFile::Open("/home/kanaki/noiseComp/Run3398_4000_v0_s72.root");
212 // AliCDBEntry *pEntry = (AliCDBEntry*)f->Get("AliCDBEntry");
213 // noisePad = (AliTPCCalPad*)pEntry->GetObject();
214 // //fHistCDBMap = noisePad->MakeHisto2D(1); //side C
218 // HLTDebug("using AliHLTTPCDigitReaderDecoder");
219 // pDigitReader = new AliHLTTPCDigitReaderDecoder(); // double-loop
220 // pDigitReader = new AliHLTTPCDigitReaderPacked();
226 int AliHLTTPCNoiseMapComponent::DoDeinit() {
227 // see header file for class documentation
231 int AliHLTTPCNoiseMapComponent::DoEvent(const AliHLTComponentEventData& evtData, AliHLTComponentTriggerData& /*trigData*/){
232 // see header file for class documentation
234 HLTInfo("--- Entering DoEvent() in TPCNoiseMap ---");
236 if(GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR )) return 0;
239 fHistSideA = new TH2F("fHistSideA","TPC Side A",250,-250,250,250,-250,250);
240 fHistSideA->SetXTitle("global X (cm)"); fHistSideA->SetYTitle("global Y (cm)");
244 fHistSideC = new TH2F("fHistSideC","TPC Side C",250,-250,250,250,-250,250);
245 fHistSideC->SetXTitle("global X (cm)"); fHistSideC->SetYTitle("global Y (cm)");
248 const AliHLTComponentBlockData *iter = NULL;
251 Int_t thissector, thisrow;
253 fHistMaxSignal = new TH2F("fHistMaxSignal","maximum signal", 250,-250,250,250,-250,250);
254 fHistTotSignal = new TH2F("fHistTotSignal","total signal", 250,-250,250,250,-250,250);
255 fHistPadRMS = new TH2F("fHistPadRMS", "RMS", 250,-250,250,250,-250,250);
257 for(iter = GetFirstInputBlock(kAliHLTDataTypeDDLRaw|kAliHLTDataOriginTPC); iter != NULL; iter = GetNextInputBlock()){
259 HLTInfo("Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s",
260 evtData.fEventID, evtData.fEventID,
261 DataType2Text(iter->fDataType).c_str(),
262 DataType2Text(kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTPC).c_str());
264 if (iter->fDataType == AliHLTTPCDefinitions::fgkDDLPackedRawDataType && GetEventCount()<2){
265 HLTWarning("data type %s is depricated, use %s (kAliHLTDataTypeDDLRaw)!",
266 DataType2Text(AliHLTTPCDefinitions::fgkDDLPackedRawDataType).c_str(),
267 DataType2Text(kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTPC).c_str());
270 if (iter->fDataType != (kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTPC)) continue;
272 UInt_t slice = AliHLTTPCDefinitions::GetMinSliceNr(*iter);
273 UInt_t partition = AliHLTTPCDefinitions::GetMinPatchNr(*iter);
275 fSpecification = iter->fSpecification;
277 AliHLTTPCDigitReader *pDigitReader = new AliHLTTPCDigitReaderDecoder;
279 pDigitReader->InitBlock(iter->fPtr,iter->fSize,partition,slice);
280 if(!pDigitReader) break;
282 //sprintf(name,"hMaxSignal_slice%d_partition%d", slice, partition);
283 //fHistMaxSignal = new TH2F(name,name,250,-250,250,250,-250,250);
285 while(pDigitReader->Next()){
286 //while( pDigitReader->NextChannel()) { // pad loop
288 fCurrentRow = pDigitReader->GetRow();
289 fCurrentRow += pDigitReader->GetRowOffset();
291 AliHLTTPCTransform::Slice2Sector(slice,fCurrentRow,thissector,thisrow);
292 AliHLTTPCTransform::Raw2Local(xyz,thissector,thisrow,pDigitReader->GetPad(),0);
294 if(slice>17) xyz[1] = (-1.0)*xyz[1];
297 AliHLTTPCTransform::Local2Global(xyz,slice);
298 // temporarily the transformation Raw2Global will be broken down to 2 steps,
299 // as there is a correction necessary at the y coordinate of the local xyz.
301 //AliHLTTPCTransform::Raw2Global(xyz,thissector,thisrow,pDigitReader->GetPad(),0);
302 // transformation from pad-row coordinates to global ones
303 // time info is not taken into account
305 // AliTPCCalROC *calRoc = noisePad->GetCalROC(thissector);
306 // calRoc->GetValue(thisrow,pDigitReader->GetPad());
308 //while( pDigitReader->NextBunch()) {
310 const UInt_t *bunchData = pDigitReader->GetSignals();
311 Float_t maxSignal = 0.;
312 Float_t totalSignal = 0.;
314 fHistSignal = new TH1F("fHistSignal", "signal distribution per pad",1024,0,1024);
316 for(Int_t i=0;i<pDigitReader->GetBunchSize();i++){
318 if((Float_t)(bunchData[i])>maxSignal){ maxSignal = (Float_t)(bunchData[i]); }
319 totalSignal += bunchData[i];
320 fHistSignal->Fill(bunchData[i]);
321 } // end for loop over bunches
323 //cout << "total signal: " << totalSignal << endl;
324 //cout << " integral: " << fHistSignal->Integral() << endl;
326 //} // end of inner while loop
328 fHistMaxSignal->Fill(xyz[0],xyz[1],maxSignal);
329 fHistTotSignal->Fill(xyz[0],xyz[1],totalSignal);
330 fHistPadRMS->Fill(xyz[0],xyz[1],fHistSignal->GetRMS());
331 delete fHistSignal; fHistSignal = NULL;
333 if(fPlotSideA || fPlotSideC){
334 if(slice<18) fHistSideA->Fill(xyz[0],xyz[1],maxSignal);
335 else fHistSideC->Fill(xyz[0],xyz[1],maxSignal);
336 } // end if plotting sides
337 } // end of while loop over pads
340 } // end of for loop over data blocks
342 if(fResetHistograms) ResetHistograms();
349 void AliHLTTPCNoiseMapComponent::MakeHistosPublic() {
350 // see header file for class documentation
352 // TFile *outputfile = new TFile("test.root","RECREATE");
353 // fHistSignal->Write();
354 // outputfile->Save();
355 // outputfile->Close();
358 histos.Add(fHistMaxSignal);
359 histos.Add(fHistTotSignal);
360 histos.Add(fHistPadRMS);
361 histos.Add(fHistCDBMap);
362 if(fPlotSideA) histos.Add(fHistSideA);
363 if(fPlotSideC) histos.Add(fHistSideC);
365 TIter iterator(&histos);
366 while(TObject *pObj=iterator.Next()){ PushBack(pObj, kAliHLTDataTypeHistogram|kAliHLTDataOriginTPC, fSpecification); }
368 //PushBack( (TObject*) &histos, kAliHLTDataTypeHistogram, fSpecification);
370 if(fHistMaxSignal) delete fHistMaxSignal; fHistMaxSignal = NULL;
371 if(fHistTotSignal) delete fHistTotSignal; fHistTotSignal = NULL;
372 if(fHistPadRMS) delete fHistPadRMS; fHistPadRMS = NULL;
373 if(fHistSideA) delete fHistSideA; fHistSideA = NULL;
374 if(fHistSideC) delete fHistSideC; fHistSideC = NULL;
378 void AliHLTTPCNoiseMapComponent::ResetHistograms(){
379 // see header file for class documentation
381 //if(fHistPartition) fHistPartition->Reset();
382 if(fHistMaxSignal) fHistMaxSignal->Reset();
383 if(fHistTotSignal) fHistTotSignal->Reset();
384 if(fHistPadRMS) fHistPadRMS->Reset();
386 if(fHistSideA) fHistSideA->Reset();
387 if(fHistSideC) fHistSideC->Reset();
390 int AliHLTTPCNoiseMapComponent::Configure(const char* arguments) {
391 // see header file for class documentation
394 if (!arguments) return iResult;
395 HLTInfo("parsing configuration string \'%s\'", arguments);
397 TString allArgs=arguments;
401 TObjArray* pTokens=allArgs.Tokenize(" ");
403 for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
404 argument=((TObjString*)pTokens->At(i))->GetString();
405 if (argument.IsNull()) continue;
407 if (argument.CompareTo("-apply-noisemap")==0) {
408 if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
409 HLTInfo("got \'-apply-noisemap\': %s", ((TObjString*)pTokens->At(i))->GetString().Data());
412 else if (argument.CompareTo("-plot-side-c")==0) {
413 if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
414 HLTInfo("got \'-plot-side-c\': %s", ((TObjString*)pTokens->At(i))->GetString().Data());
417 else if (argument.CompareTo("-plot-side-a")==0) {
418 if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
419 HLTInfo("got \'-plot-side-a\': %s", ((TObjString*)pTokens->At(i))->GetString().Data());
422 else if(argument.CompareTo("-reset-histograms")==0){
423 if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
424 HLTInfo("got \'-reset-histograms\': %s", ((TObjString*)pTokens->At(i))->GetString().Data());
428 HLTError("unknown argument %s", argument.Data());
439 HLTError("missing parameter for argument %s", argument.Data());
445 int AliHLTTPCNoiseMapComponent::Reconfigure(const char* cdbEntry, const char* chainId) {
446 // see header file for class documentation
448 const char* path="HLT/ConfigTPC/TPCNoiseMapComponent";
449 const char* defaultNotify="";
452 defaultNotify=" (default)";
456 HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
457 AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
459 TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
461 HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
462 iResult=Configure(pString->GetString().Data());
464 HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
467 HLTError("cannot fetch object \"%s\" from CDB", path);