1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
17 ///////////////////////////////////////////////////////////////////////////////
19 // Class providing the calibration parameters by accessing the CDB //
21 // Request an instance with AliTPCcalibDB::Instance() //
22 // If a new event is processed set the event number with SetRun //
23 // Then request the calibration data ////
28 // Simulation - not yet
29 // Reconstruction - AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
31 // 1.) pad by pad calibration - AliTPCCalPad
34 // Simulation: AliTPCDigitizer::ExecFast - Multiply by gain
35 // Reconstruction : AliTPCclustererMI::Digits2Clusters - Divide by gain
38 // Simulation: AliTPCDigitizer::ExecFast
39 // Reconstruction: AliTPCclustererMI::FindClusters(AliTPCCalROC * noiseROC)
40 // Noise depending cut on clusters charge (n sigma)
42 // Simulation: Not used yet - To be impleneted - Rounding to the nearest integer
43 // Reconstruction: Used in AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
44 // if data taken without zero suppression
45 // Currently switch in fRecoParam->GetCalcPedestal();
48 // Simulation: applied in the AliTPC::MakeSector - adding offset
49 // Reconstruction: AliTPCTransform::Transform() - remove offset
50 // AliTPCTransform::Transform() - to be called
51 // in AliTPCtracker::Transform()
54 // 2.) Space points transformation:
56 // a.) General coordinate tranformation - AliTPCtransform (see $ALICE_ROOT/TPC/AliTPCtransform.cxx)
57 // Created on fly - use the other calibration components
58 // Unisochronity - (substract time0 - pad by pad)
59 // Drift velocity - Currently common drift velocity - functionality of AliTPCParam
61 // Simulation - Not used directly (the effects are applied one by one (see AliTPC::MakeSector)
63 // AliTPCclustererMI::AddCluster
64 // AliTPCtrackerMI::Transform
65 // b.) ExB effect calibration -
66 // classes (base class AliTPCExB, implementation- AliTPCExBExact.h AliTPCExBFirst.h)
67 // a.a) Simulation: applied in the AliTPC::MakeSector -
68 // calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
69 // a.b) Reconstruction -
71 // in AliTPCtransform::Correct() - called calib->GetExB()->Correct(dxyz0,dxyz1)
73 // 3.) cluster error, shape and Q parameterization
77 ///////////////////////////////////////////////////////////////////////////////
83 #include <AliCDBManager.h>
84 #include <AliCDBEntry.h>
87 #include "AliTPCcalibDB.h"
88 #include "AliTPCAltroMapping.h"
89 #include "AliTPCExB.h"
91 #include "AliTPCCalROC.h"
92 #include "AliTPCCalPad.h"
93 #include "AliTPCSensorTempArray.h"
94 #include "AliTPCTransform.h"
104 #include "TObjArray.h"
105 #include "TObjString.h"
107 #include "AliTPCCalPad.h"
108 #include "AliTPCCalibPulser.h"
109 #include "AliTPCCalibPedestal.h"
110 #include "AliTPCCalibCE.h"
116 ClassImp(AliTPCcalibDB)
118 AliTPCcalibDB* AliTPCcalibDB::fgInstance = 0;
119 Bool_t AliTPCcalibDB::fgTerminated = kFALSE;
122 //_ singleton implementation __________________________________________________
123 AliTPCcalibDB* AliTPCcalibDB::Instance()
126 // Singleton implementation
127 // Returns an instance of this class, it is created if neccessary
130 if (fgTerminated != kFALSE)
134 fgInstance = new AliTPCcalibDB();
139 void AliTPCcalibDB::Terminate()
142 // Singleton implementation
143 // Deletes the instance of this class and sets the terminated flag, instances cannot be requested anymore
144 // This function can be called several times.
147 fgTerminated = kTRUE;
156 //_____________________________________________________________________________
157 AliTPCcalibDB::AliTPCcalibDB():
175 Update(); // temporary
178 AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& ):
194 // Copy constructor invalid -- singleton implementation
196 Error("copy constructor","invalid -- singleton implementation");
199 AliTPCcalibDB& AliTPCcalibDB::operator= (const AliTPCcalibDB& )
202 // Singleton implementation - no assignment operator
204 Error("operator =", "assignment operator not implemented");
210 //_____________________________________________________________________________
211 AliTPCcalibDB::~AliTPCcalibDB()
217 // don't delete anything, CDB cache is active!
218 //if (fPadGainFactor) delete fPadGainFactor;
219 //if (fPadTime0) delete fPadTime0;
220 //if (fPadNoise) delete fPadNoise;
224 //_____________________________________________________________________________
225 AliCDBEntry* AliTPCcalibDB::GetCDBEntry(const char* cdbPath)
228 // Retrieves an entry with path <cdbPath> from the CDB.
232 AliCDBEntry* entry = AliCDBManager::Instance()->Get(cdbPath, fRun);
235 sprintf(chinfo,"AliTPCcalibDB: Failed to get entry:\t%s ", cdbPath);
243 //_____________________________________________________________________________
244 void AliTPCcalibDB::SetRun(Long64_t run)
247 // Sets current run number. Calibration data is read from the corresponding file.
257 void AliTPCcalibDB::Update(){
259 AliCDBEntry * entry=0;
261 Bool_t cdbCache = AliCDBManager::Instance()->GetCacheFlag(); // save cache status
262 AliCDBManager::Instance()->SetCacheFlag(kTRUE); // activate CDB cache
265 entry = GetCDBEntry("TPC/Calib/PadGainFactor");
267 //if (fPadGainFactor) delete fPadGainFactor;
268 entry->SetOwner(kTRUE);
269 fPadGainFactor = (AliTPCCalPad*)entry->GetObject();
272 entry = GetCDBEntry("TPC/Calib/PadTime0");
274 //if (fPadTime0) delete fPadTime0;
275 entry->SetOwner(kTRUE);
276 fPadTime0 = (AliTPCCalPad*)entry->GetObject();
280 entry = GetCDBEntry("TPC/Calib/PadNoise");
282 //if (fPadNoise) delete fPadNoise;
283 entry->SetOwner(kTRUE);
284 fPadNoise = (AliTPCCalPad*)entry->GetObject();
287 entry = GetCDBEntry("TPC/Calib/Pedestals");
289 //if (fPedestals) delete fPedestals;
290 entry->SetOwner(kTRUE);
291 fPedestals = (AliTPCCalPad*)entry->GetObject();
294 entry = GetCDBEntry("TPC/Calib/Temperature");
296 //if (fTemperature) delete fTemperature;
297 entry->SetOwner(kTRUE);
298 fTemperature = (AliTPCSensorTempArray*)entry->GetObject();
301 entry = GetCDBEntry("TPC/Calib/Parameters");
303 //if (fPadNoise) delete fPadNoise;
304 entry->SetOwner(kTRUE);
305 fParam = (AliTPCParam*)(entry->GetObject()->Clone());
308 entry = GetCDBEntry("TPC/Calib/ClusterParam");
310 //if (fPadNoise) delete fPadNoise;
311 entry->SetOwner(kTRUE);
312 fClusterParam = (AliTPCClusterParam*)(entry->GetObject()->Clone());
315 entry = GetCDBEntry("TPC/Calib/Mapping");
317 //if (fPadNoise) delete fPadNoise;
318 entry->SetOwner(kTRUE);
319 TObjArray * array = dynamic_cast<TObjArray*>(entry->GetObject());
320 if (array && array->GetEntriesFast()==6){
321 fMapping = new AliTPCAltroMapping*[6];
322 for (Int_t i=0; i<6; i++){
323 fMapping[i] = dynamic_cast<AliTPCAltroMapping*>(array->At(i));
330 entry = GetCDBEntry("TPC/Calib/ExB");
332 entry->SetOwner(kTRUE);
333 fExB=dynamic_cast<AliTPCExB*>(entry->GetObject()->Clone());
337 fTransform=new AliTPCTransform();
341 AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
347 void AliTPCcalibDB::CreateObjectList(const Char_t *filename, TObjArray *calibObjects)
350 // Create calibration objects and read contents from OCDB
352 if ( calibObjects == 0x0 ) return;
355 if ( !in.is_open() ){
356 fprintf(stderr,"Error: cannot open list file '%s'", filename);
360 AliTPCCalPad *calPad=0x0;
366 TObjArray *arrFileLine = sFile.Tokenize("\n");
368 TIter nextLine(arrFileLine);
370 TObjString *sObjLine=0x0;
371 while ( (sObjLine = (TObjString*)nextLine()) ){
372 TString sLine(sObjLine->GetString());
374 TObjArray *arrNextCol = sLine.Tokenize("\t");
376 TObjString *sObjType = (TObjString*)(arrNextCol->At(0));
377 TObjString *sObjFileName = (TObjString*)(arrNextCol->At(1));
379 if ( !sObjType || ! sObjFileName ) continue;
380 TString sType(sObjType->GetString());
381 TString sFileName(sObjFileName->GetString());
382 printf("%s\t%s\n",sType.Data(),sFileName.Data());
384 TFile *fIn = TFile::Open(sFileName);
386 fprintf(stderr,"File not found: '%s'", sFileName.Data());
390 if ( sType == "CE" ){
391 AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
393 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());
394 calPad->SetNameTitle("CETmean","CETmean");
395 calibObjects->Add(calPad);
397 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());
398 calPad->SetNameTitle("CEQmean","CEQmean");
399 calibObjects->Add(calPad);
401 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
402 calPad->SetNameTitle("CETrms","CETrms");
403 calibObjects->Add(calPad);
405 } else if ( sType == "Pulser") {
406 AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
408 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());
409 calPad->SetNameTitle("PulserTmean","PulserTmean");
410 calibObjects->Add(calPad);
412 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());
413 calPad->SetNameTitle("PulserQmean","PulserQmean");
414 calibObjects->Add(calPad);
416 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
417 calPad->SetNameTitle("PulserTrms","PulserTrms");
418 calibObjects->Add(calPad);
420 } else if ( sType == "Pedestals") {
421 AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
423 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());
424 calPad->SetNameTitle("Pedestals","Pedestals");
425 calibObjects->Add(calPad);
427 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());
428 calPad->SetNameTitle("Noise","Noise");
429 calibObjects->Add(calPad);
432 fprintf(stderr,"Undefined Type: '%s'",sType.Data());
441 void AliTPCcalibDB::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
443 // Write a tree with all available information
444 // if mapFileName is specified, the Map information are also written to the tree
445 // pads specified in outlierPad are not used for calculating statistics
446 // - the same function as AliTPCCalPad::MakeTree -
448 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
450 TObjArray* mapIROCs = 0;
451 TObjArray* mapOROCs = 0;
452 TVectorF *mapIROCArray = 0;
453 TVectorF *mapOROCArray = 0;
454 Int_t mapEntries = 0;
455 TString* mapNames = 0;
458 TFile mapFile(mapFileName, "read");
460 TList* listOfROCs = mapFile.GetListOfKeys();
461 mapEntries = listOfROCs->GetEntries()/2;
462 mapIROCs = new TObjArray(mapEntries*2);
463 mapOROCs = new TObjArray(mapEntries*2);
464 mapIROCArray = new TVectorF[mapEntries];
465 mapOROCArray = new TVectorF[mapEntries];
467 mapNames = new TString[mapEntries];
468 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
469 TString nameROC(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
470 nameROC.Remove(nameROC.Length()-4, 4);
471 mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "IROC").Data()), ivalue);
472 mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "OROC").Data()), ivalue);
473 mapNames[ivalue].Append(nameROC);
476 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
477 mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
478 mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
480 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
481 (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
482 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
483 (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
486 } // if (mapFileName)
488 TTreeSRedirector cstream(fileName);
489 Int_t arrayEntries = array->GetEntries();
491 TString* names = new TString[arrayEntries];
492 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
493 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
495 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
497 // get statistic for given sector
499 TVectorF median(arrayEntries);
500 TVectorF mean(arrayEntries);
501 TVectorF rms(arrayEntries);
502 TVectorF ltm(arrayEntries);
503 TVectorF ltmrms(arrayEntries);
504 TVectorF medianWithOut(arrayEntries);
505 TVectorF meanWithOut(arrayEntries);
506 TVectorF rmsWithOut(arrayEntries);
507 TVectorF ltmWithOut(arrayEntries);
508 TVectorF ltmrmsWithOut(arrayEntries);
510 TVectorF *vectorArray = new TVectorF[arrayEntries];
511 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
512 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
514 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
515 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
516 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
517 AliTPCCalROC* outlierROC = 0;
518 if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
520 median[ivalue] = calROC->GetMedian();
521 mean[ivalue] = calROC->GetMean();
522 rms[ivalue] = calROC->GetRMS();
523 Double_t ltmrmsValue = 0;
524 ltm[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction);
525 ltmrms[ivalue] = ltmrmsValue;
527 medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
528 meanWithOut[ivalue] = calROC->GetMean(outlierROC);
529 rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
531 ltmWithOut[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction, outlierROC);
532 ltmrmsWithOut[ivalue] = ltmrmsValue;
541 medianWithOut[ivalue] = 0.;
542 meanWithOut[ivalue] = 0.;
543 rmsWithOut[ivalue] = 0.;
544 ltmWithOut[ivalue] = 0.;
545 ltmrmsWithOut[ivalue] = 0.;
550 // fill vectors of variable per pad
552 TVectorF *posArray = new TVectorF[8];
553 for (Int_t ivalue = 0; ivalue < 8; ivalue++)
554 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
556 Float_t posG[3] = {0};
557 Float_t posL[3] = {0};
559 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
560 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
561 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
562 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
563 posArray[0][ichannel] = irow;
564 posArray[1][ichannel] = ipad;
565 posArray[2][ichannel] = posL[0];
566 posArray[3][ichannel] = posL[1];
567 posArray[4][ichannel] = posG[0];
568 posArray[5][ichannel] = posG[1];
569 posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
570 posArray[7][ichannel] = ichannel;
572 // loop over array containing AliTPCCalPads
573 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
574 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
575 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
577 (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
579 (vectorArray[ivalue])[ichannel] = 0;
585 cstream << "calPads" <<
586 "sector=" << isector;
588 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
589 cstream << "calPads" <<
590 (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
591 (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
592 (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
593 (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
594 (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
596 cstream << "calPads" <<
597 (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
598 (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
599 (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
600 (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
601 (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
605 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
606 cstream << "calPads" <<
607 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
611 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
613 cstream << "calPads" <<
614 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
616 cstream << "calPads" <<
617 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
621 cstream << "calPads" <<
622 "row.=" << &posArray[0] <<
623 "pad.=" << &posArray[1] <<
624 "lx.=" << &posArray[2] <<
625 "ly.=" << &posArray[3] <<
626 "gx.=" << &posArray[4] <<
627 "gy.=" << &posArray[5] <<
628 "rpad.=" << &posArray[6] <<
629 "channel.=" << &posArray[7];
631 cstream << "calPads" <<
635 delete[] vectorArray;
643 delete[] mapIROCArray;
644 delete[] mapOROCArray;