3 /**************************************************************************
4 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
6 * Authors: Matthias Richter <Matthias.Richter@ift.uib.no> *
7 * Timm Steinbeck <timm@kip.uni-heidelberg.de> *
8 * for The ALICE Off-line 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 AliHLTTRDCalibrationComponent.cxx
20 @author Timm Steinbeck, Matthias Richter
22 @brief A TRDCalibration processing component for the HLT. */
33 #include "TProfile2D.h"
35 #include "AliHLTReadoutList.h"
37 #include "AliHLTTRDCalibrationComponent.h"
38 #include "AliHLTTRDDefinitions.h"
39 #include "AliHLTTRDUtils.h"
41 #include "AliCDBManager.h"
42 #include "AliCDBStorage.h"
43 #include "AliRawReaderMemory.h"
45 #include "AliTRDCalPad.h"
46 #include "AliTRDCalDet.h"
48 #include "AliTRDCalibraFillHisto.h"
49 #include "AliTRDtrackV1.h"
51 #include "AliTRDCalibraFit.h"
52 #include "AliTRDCalibraMode.h"
53 #include "AliTRDCalibraVector.h"
54 #include "AliTRDCalibraVdriftLinearFit.h"
60 ClassImp(AliHLTTRDCalibrationComponent);
62 AliHLTTRDCalibrationComponent::AliHLTTRDCalibrationComponent()
63 : AliHLTCalibrationProcessor(),
64 fTRDCalibraFillHisto(NULL),
70 fRecievedTimeBins(kFALSE),
74 // Default constructor
77 AliHLTTRDCalibrationComponent::~AliHLTTRDCalibrationComponent()
82 const char* AliHLTTRDCalibrationComponent::GetComponentID()
84 // Return the component ID const char *
85 return "TRDCalibration"; // The ID of this component
88 void AliHLTTRDCalibrationComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
90 // Get the list of input data
91 list.clear(); // We do not have any requirements for our input data type(s).
92 list.push_back(AliHLTTRDDefinitions::fgkTracksDataType);
95 AliHLTComponentDataType AliHLTTRDCalibrationComponent::GetOutputDataType()
97 // Get the output data type
98 return kAliHLTMultipleDataType;
99 // return AliHLTTRDDefinitions::fgkCalibrationDataType;
103 int AliHLTTRDCalibrationComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
105 // Get the output data type
107 tgtList.push_back(AliHLTTRDDefinitions::fgkCalibrationDataType);
108 tgtList.push_back(AliHLTTRDDefinitions::fgkEORCalibrationDataType);
109 return tgtList.size();
112 void AliHLTTRDCalibrationComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
114 // Get the output data size
115 constBase = fOutputSize;
119 AliHLTComponent* AliHLTTRDCalibrationComponent::Spawn()
121 // Spawn function, return new instance of this class
122 return new AliHLTTRDCalibrationComponent;
125 Int_t AliHLTTRDCalibrationComponent::ScanArgument( int argc, const char** argv )
127 // perform initialization. We check whether our relative output size is specified in the arguments.
131 fTrgStrings = new TObjArray();
135 HLTDebug("argv[%d] == %s", i, argv[i] );
136 if ( !strcmp( argv[i], "output_size" ) )
140 HLTError("Missing output_size parameter");
143 HLTDebug("argv[%d+1] == %s", i, argv[i+1] );
144 fOutputSize = strtoul( argv[i+1], &cpErr, 0 );
147 HLTError("Cannot convert output_size parameter '%s'", argv[i+1] );
150 HLTInfo("Output size set to %lu %%", fOutputSize );
154 if ( !strcmp( argv[i], "-everyNevent" ) )
158 HLTError("Missing everyNevent parameter");
161 HLTDebug("argv[%d+1] == %s", i, argv[i+1] );
162 HLTInfo("Option -everyNevent depreceated");
166 if ( !strcmp( argv[i], "-TrgStr" ) )
170 HLTError("Missing parameter for mbTriggerString");
173 HLTDebug("argv[%d+1] == %s", i, argv[i+1] );
174 fTrgStrings->Add(new TObjString(argv[i+1]));
179 if ( !strcmp( argv[i], "-acceptTrgStr" ) )
185 if ( !strcmp( argv[i], "-rejectTrgStr" ) )
193 HLTError("Unknown option '%s'", argv[i] );
200 Int_t AliHLTTRDCalibrationComponent::InitCalibration()
204 fTrgStrings = new TObjArray();
206 if(!AliCDBManager::Instance()->IsDefaultStorageSet()){
207 HLTError("DefaultStorage is not set in CDBManager");
210 if(AliCDBManager::Instance()->GetRun()<0){
211 HLTError("Run Number is not set in CDBManager");
214 HLTInfo("CDB default storage: %s; RunNo: %i", (AliCDBManager::Instance()->GetDefaultStorage()->GetBaseFolder()).Data(), AliCDBManager::Instance()->GetRun());
216 if(fTrgStrings->GetEntriesFast()>0 && !fAccRejTrg){
217 HLTError("Trigger string(s) given, but acceptTrgStr or rejectTrgStr not selected");
221 fTRDCalibraFillHisto = AliTRDCalibraFillHisto::Instance();
222 fTRDCalibraFillHisto->SetIsHLT(kTRUE);
223 fTRDCalibraFillHisto->SetHisto2d(); // choose to use histograms
224 fTRDCalibraFillHisto->SetCH2dOn(); // choose to calibrate the gain
225 fTRDCalibraFillHisto->SetPH2dOn(); // choose to calibrate the drift velocity
226 fTRDCalibraFillHisto->SetPRF2dOn(); // choose to look at the PRF
227 fTRDCalibraFillHisto->SetIsHLT(); // per detector
228 //fTRDCalibraFillHisto->SetDebugLevel(1);// debug
229 fTRDCalibraFillHisto->SetFillWithZero(kTRUE);
230 fTRDCalibraFillHisto->SetLinearFitterOn(kTRUE);
231 fTRDCalibraFillHisto->SetNumberBinCharge(100);
233 fTracksArray = new TClonesArray("AliTRDtrackV1");
234 fOutArray = new TObjArray(4);
235 fAfterRunArray=new TObjArray(5);
236 fDisplayArray=new TObjArray(4);
238 HLTDebug("run SetupCTPData");
244 Int_t AliHLTTRDCalibrationComponent::DeinitCalibration()
247 // Deinitialization of the component
249 HLTDebug("DeinitCalibration");
250 delete fTracksArray; fTracksArray=0;
251 //fTRDCalibraFillHisto->Destroy();
252 //fOutArray->Delete();
253 delete fOutArray; fOutArray=0;
254 fAfterRunArray->Delete();
255 delete fAfterRunArray; fAfterRunArray=0;
256 fDisplayArray->Delete();
257 delete fDisplayArray; fDisplayArray=0;
258 fTrgStrings->Delete();
259 delete fTrgStrings; fTrgStrings=0;
263 Int_t AliHLTTRDCalibrationComponent::ProcessCalibration(const AliHLTComponent_EventData& evtData,
264 const AliHLTComponent_BlockData* blocks,
265 AliHLTComponent_TriggerData& trigData,
266 AliHLTUInt8_t* /*outputPtr*/,
267 AliHLTUInt32_t& /*size*/,
268 vector<AliHLTComponent_BlockData>& /*outputBlocks*/)
270 HLTDebug("NofBlocks %lu", evtData.fBlockCnt );
274 // Loop over all input blocks in the event
275 vector<AliHLTComponent_DataType> expectedDataTypes;
276 GetInputDataTypes(expectedDataTypes);
277 for ( unsigned long iBlock = 0; iBlock < evtData.fBlockCnt; iBlock++ )
279 const AliHLTComponentBlockData &block = blocks[iBlock];
280 AliHLTComponentDataType inputDataType = block.fDataType;
281 Bool_t correctDataType = kFALSE;
283 for(UInt_t i = 0; i < expectedDataTypes.size(); i++)
284 if( expectedDataTypes.at(i) == inputDataType)
285 correctDataType = kTRUE;
286 if (!correctDataType) {
287 HLTDebug( "Block # %i/%i; Event 0x%08LX (%Lu) Wrong received datatype: %s - Skipping",
288 iBlock, evtData.fBlockCnt,
289 evtData.fEventID, evtData.fEventID,
290 DataType2Text(inputDataType).c_str());
294 HLTDebug("We get the right data type: Block # %i/%i; Event 0x%08LX (%Lu) Received datatype: %s; Block Size: %i",
295 iBlock, evtData.fBlockCnt-1,
296 evtData.fEventID, evtData.fEventID,
297 DataType2Text(inputDataType).c_str(),
302 AliHLTTRDUtils::ReadTracks(fTracksArray, block.fPtr, block.fSize, &nTimeBins);
305 if(!fRecievedTimeBins){
306 HLTDebug("Reading number of time bins from input block. Value is: %d", nTimeBins);
307 fTRDCalibraFillHisto->Init2Dhistos(nTimeBins); // initialise the histos
308 fTRDCalibraFillHisto->SetNumberClusters(0); // At least 1 clusters
309 fTRDCalibraFillHisto->SetNumberClustersf(nTimeBins); // Not more than %d clusters
310 fRecievedTimeBins=kTRUE;
314 Bool_t TriggerPassed=kFALSE;
318 TriggerPassed=kFALSE;
319 for(int i = 0; i < fTrgStrings->GetEntriesFast(); i++){
320 const TObjString *const obString=(TObjString*)fTrgStrings->At(i);
321 const TString tString=obString->GetString();
322 //printf("Trigger Output: %i\n",EvaluateCTPTriggerClass(tString.Data(),trigData));
323 if(EvaluateCTPTriggerClass(tString.Data(),trigData)){TriggerPassed=kTRUE; break;}
328 for(int i = 0; i < fTrgStrings->GetEntriesFast(); i++){
329 const TObjString *const obString=(TObjString*)fTrgStrings->At(i);
330 const TString tString=obString->GetString();
331 if(EvaluateCTPTriggerClass(tString.Data(),trigData)){TriggerPassed=kFALSE; break;}
337 Int_t nbEntries = fTracksArray->GetEntries();
338 HLTDebug(" %i TRDtracks in tracksArray", nbEntries);
339 AliTRDtrackV1* trdTrack = 0x0;
340 for (Int_t i = 0; i < nbEntries; i++){
341 HLTDebug("%i/%i: ", i+1, nbEntries);
342 trdTrack = (AliTRDtrackV1*)fTracksArray->At(i);
344 fTRDCalibraFillHisto->SetCH2dOn(TriggerPassed);
345 fTRDCalibraFillHisto->UpdateHistogramsV1(trdTrack);
346 fTRDCalibraFillHisto->SetCH2dOn(kTRUE);
350 if(!fOutArray->At(0))FormOutput(0);
351 if(!fDisplayArray->At(0))FormOutput(1);
352 PushBack(fDisplayArray, AliHLTTRDDefinitions::fgkCalibrationDataType);
354 fTracksArray->Delete();
363 * Form output array of histrograms
365 //============================================================================
366 void AliHLTTRDCalibrationComponent::FormOutput(Int_t param)
369 TH2I *hCH2d = fTRDCalibraFillHisto->GetCH2d();
370 if(!param)fOutArray->Add(hCH2d);
371 else fDisplayArray->Add(hCH2d);
373 // drift velocity histo
374 TProfile2D *hPH2d = fTRDCalibraFillHisto->GetPH2d();
375 if(!param)fOutArray->Add(hPH2d);
376 else fDisplayArray->Add(hPH2d);
379 TProfile2D *hPRF2d = fTRDCalibraFillHisto->GetPRF2d();
380 if(!param)fOutArray->Add(hPRF2d);
381 else fDisplayArray->Add(hPRF2d);
385 AliTRDCalibraVdriftLinearFit *hVdriftLinearFitOne=(AliTRDCalibraVdriftLinearFit *)fTRDCalibraFillHisto->GetVdriftLinearFit();
386 fOutArray->Add(hVdriftLinearFitOne);
389 TH2S *hVdriftLinearFitOne = (TH2S *)(((AliTRDCalibraVdriftLinearFit *)fTRDCalibraFillHisto->GetVdriftLinearFit())->GetLinearFitterHisto(10,kTRUE));
390 fDisplayArray->Add(hVdriftLinearFitOne);
393 HLTDebug("GetCH2d = 0x%x; NEntries = %i; size = %i", hCH2d, hCH2d->GetEntries(), sizeof(*hCH2d));
395 HLTDebug("GetPH2d = 0x%x; NEntries = %i; size = %i", hPH2d, hPH2d->GetEntries(), sizeof(*hPH2d));
397 HLTDebug("GetPRF2d = 0x%x; NEntries = %i; size = %i", hPRF2d, hPRF2d->GetEntries(), sizeof(*hPRF2d));
399 //HLTDebug("GetVdriftLinearFit = 0x%x; size = %i", hVdriftLinearFitOne, sizeof(hVdriftLinearFitOne));
401 HLTDebug("output Array: pointer = 0x%x; NEntries = %i; size = %i", fOutArray, fOutArray->GetEntries(), sizeof(fOutArray));
405 Int_t AliHLTTRDCalibrationComponent::ShipDataToFXS(const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/)
407 //fTRDCalibraFillHisto->DestroyDebugStreamer();
409 AliHLTReadoutList rdList(AliHLTReadoutList::kTRD);
413 fOutArray->Remove(fOutArray->FindObject("AliTRDCalibraVdriftLinearFit"));
414 //fOutArray->Remove(fOutArray->FindObject("PRF2d"));
415 //fOutArray->Remove(fOutArray->FindObject("PH2d"));
416 //fOutArray->Remove(fOutArray->FindObject("CH2d"));
418 //if(!(fOutArray->FindObject("CH2d"))) {
419 // TH2I * ch2d = new TH2I("CH2d","Nz0Nrphi0",100,0.0,300.0,540,0,540);
420 // fOutArray->Add(ch2d);
423 //if(!(fOutArray->FindObject("PH2d"))) {
424 // TProfile2D * ph2d = new TProfile2D("PH2d","Nz0Nrphi0",30,-0.05,2.95,540,0,540);
425 // fOutArray->Add(ph2d);
428 //if(!(fOutArray->FindObject("PRF2d"))) {
429 // TProfile2D * prf2d = new TProfile2D("PRF2d","Nz0Nrphi0Ngp3",60,-9.0,9.0,540,0,540);
430 // fOutArray->Add(prf2d);
434 HLTDebug("Size of the fOutArray is %d\n",fOutArray->GetEntriesFast());
437 TString fileName="$ALIHLT_TOPDIR/build-debug/output/CalibHistoDump_run";
439 HLTInfo("Dumping Histogram file to %s",fileName.Data());
440 TFile* file = TFile::Open(fileName, "RECREATE");
441 //fAfterRunArray->Write();
444 HLTInfo("Histogram file dumped");
447 PushToFXS((TObject*)fOutArray, "TRD", "GAINDRIFTPRF", rdList.Buffer() );
448 //PushToFXS((TObject*)fOutArray->FindObject("CH2d"), "TRD", "GAINDRIFTPRF", rdList.Buffer() );
453 Int_t AliHLTTRDCalibrationComponent::EORCalibration()
455 //Also Fill histograms for the online display
456 TH2I *hCH2d=(TH2I*)fOutArray->FindObject("CH2d");
457 TProfile2D *hPH2d=(TProfile2D*)fOutArray->FindObject("PH2d");
458 TProfile2D *hPRF2d= (TProfile2D*)fOutArray->FindObject("PRF2d");
459 AliTRDCalibraVdriftLinearFit* hVdriftLinearFit = (AliTRDCalibraVdriftLinearFit*)fOutArray->FindObject("AliTRDCalibraVdriftLinearFit");
462 if(!hCH2d || !hPH2d || !hPRF2d || !hVdriftLinearFit) return 0;
465 AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
468 calibra->SetMinEntries(100);
469 calibra->AnalyseCH(hCH2d);
470 Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(0))
471 + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(0));
472 Int_t nbfit = calibra->GetNumberFit();
473 Int_t nbE = calibra->GetNumberEnt();
474 TH1F *coefgain = 0x0;
477 // (nbfit >= 0.2*nbE)) {
478 // create the cal objects
479 //calibra->PutMeanValueOtherVectorFit(1,kTRUE);
480 TObjArray object = calibra->GetVectorFit();
481 AliTRDCalDet *objgaindet = calibra->CreateDetObjectGain(&object,kFALSE);
482 coefgain = objgaindet->MakeHisto1DAsFunctionOfDet();
484 calibra->ResetVectorFit();
486 // vdrift second method
487 calibra->SetMinEntries(100); // If there is less than 100
488 hVdriftLinearFit->FillPEArray();
489 calibra->AnalyseLinearFitters(hVdriftLinearFit);
491 nbfit = calibra->GetNumberFit();
492 nbE = calibra->GetNumberEnt();
493 TH1F *coefdriftsecond = 0x0;
496 // (nbfit >= 0.1*nbE)) {
497 // create the cal objects
498 //calibra->PutMeanValueOtherVectorFit(1,kTRUE);
499 object = calibra->GetVectorFit();
500 AliTRDCalDet *objdriftvelocitydetsecond = calibra->CreateDetObjectVdrift(&object,kTRUE);
501 objdriftvelocitydetsecond->SetTitle("secondmethodvdrift");
502 coefdriftsecond = objdriftvelocitydetsecond->MakeHisto1DAsFunctionOfDet();
504 calibra->ResetVectorFit();
506 // vdrift first method
507 calibra->SetMinEntries(100*20); // If there is less than 20000
508 calibra->AnalysePH(hPH2d);
509 nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(1))
510 + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(1));
511 nbfit = calibra->GetNumberFit();
512 nbE = calibra->GetNumberEnt();
513 TH1F *coefdrift = 0x0;
517 // (nbfit >= 0.2*nbE)) {
518 // create the cal objects
519 //calibra->PutMeanValueOtherVectorFit(1,kTRUE);
520 //calibra->PutMeanValueOtherVectorFit2(1,kTRUE);
521 object = calibra->GetVectorFit();
522 AliTRDCalDet *objdriftvelocitydet = calibra->CreateDetObjectVdrift(&object,kTRUE);
523 coefdrift = objdriftvelocitydet->MakeHisto1DAsFunctionOfDet();
524 object = calibra->GetVectorFit2();
525 AliTRDCalDet *objtime0det = calibra->CreateDetObjectT0(&object,kTRUE);
526 coeft0 = objtime0det->MakeHisto1DAsFunctionOfDet();
528 calibra->ResetVectorFit();
532 calibra->SetMinEntries(200);
533 calibra->AnalysePRFMarianFit(hPRF2d);
534 nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(2))
535 + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(2));
536 nbfit = calibra->GetNumberFit();
537 nbE = calibra->GetNumberEnt();
541 // (nbfit >= 0.95*nbE)) {
542 // create cal pad objects
543 object = calibra->GetVectorFit();
544 TObject *objPRFpad = calibra->CreatePadObjectPRF(&object);
545 coefprf = ((AliTRDCalPad *) objPRFpad)->MakeHisto1D();
547 calibra->ResetVectorFit();
550 coefgain->SetName("coefgain");
551 coefprf->SetName("coefprf");
552 coefdrift->SetName("coefdrift");
553 coefdriftsecond->SetName("coefdriftsecond");
554 coeft0->SetName("coeft0");
555 if(coefgain) fAfterRunArray->Add(coefgain);
556 if(coefprf) fAfterRunArray->Add(coefprf);
557 if(coefdrift) fAfterRunArray->Add(coefdrift);
558 if(coefdriftsecond) fAfterRunArray->Add(coefdriftsecond);
559 if(coeft0) fAfterRunArray->Add(coeft0);
562 if(coefgain||coefprf||coefdrift||coeft0||coefdriftsecond) {
563 PushBack(fAfterRunArray, AliHLTTRDDefinitions::fgkEORCalibrationDataType);
567 TString fileName="$ALIHLT_TOPDIR/build-debug/output/CalibHistoDump_run";
569 HLTInfo("Dumping Histogram file to %s",fileName.Data());
570 TFile* file = TFile::Open(fileName, "RECREATE");
571 fAfterRunArray->Write();
574 HLTInfo("Histogram file dumped");