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 <AliMagWrapCheb.h>
89 #include "AliTPCcalibDB.h"
90 #include "AliTPCAltroMapping.h"
91 #include "AliTPCExB.h"
93 #include "AliTPCCalROC.h"
94 #include "AliTPCCalPad.h"
95 #include "AliTPCSensorTempArray.h"
96 #include "AliTPCTransform.h"
106 #include "TObjArray.h"
107 #include "TObjString.h"
109 #include "AliTPCCalPad.h"
110 #include "AliTPCCalibPulser.h"
111 #include "AliTPCCalibPedestal.h"
112 #include "AliTPCCalibCE.h"
113 #include "AliTPCExBFirst.h"
119 ClassImp(AliTPCcalibDB)
121 AliTPCcalibDB* AliTPCcalibDB::fgInstance = 0;
122 Bool_t AliTPCcalibDB::fgTerminated = kFALSE;
123 TObjArray AliTPCcalibDB::fgExBArray; // array of ExB corrections
126 //_ singleton implementation __________________________________________________
127 AliTPCcalibDB* AliTPCcalibDB::Instance()
130 // Singleton implementation
131 // Returns an instance of this class, it is created if neccessary
134 if (fgTerminated != kFALSE)
138 fgInstance = new AliTPCcalibDB();
143 void AliTPCcalibDB::Terminate()
146 // Singleton implementation
147 // Deletes the instance of this class and sets the terminated flag, instances cannot be requested anymore
148 // This function can be called several times.
151 fgTerminated = kTRUE;
160 //_____________________________________________________________________________
161 AliTPCcalibDB::AliTPCcalibDB():
180 Update(); // temporary
183 AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& ):
200 // Copy constructor invalid -- singleton implementation
202 Error("copy constructor","invalid -- singleton implementation");
205 AliTPCcalibDB& AliTPCcalibDB::operator= (const AliTPCcalibDB& )
208 // Singleton implementation - no assignment operator
210 Error("operator =", "assignment operator not implemented");
216 //_____________________________________________________________________________
217 AliTPCcalibDB::~AliTPCcalibDB()
223 // don't delete anything, CDB cache is active!
224 //if (fPadGainFactor) delete fPadGainFactor;
225 //if (fPadTime0) delete fPadTime0;
226 //if (fPadNoise) delete fPadNoise;
230 //_____________________________________________________________________________
231 AliCDBEntry* AliTPCcalibDB::GetCDBEntry(const char* cdbPath)
234 // Retrieves an entry with path <cdbPath> from the CDB.
238 AliCDBEntry* entry = AliCDBManager::Instance()->Get(cdbPath, fRun);
241 sprintf(chinfo,"AliTPCcalibDB: Failed to get entry:\t%s ", cdbPath);
249 //_____________________________________________________________________________
250 void AliTPCcalibDB::SetRun(Long64_t run)
253 // Sets current run number. Calibration data is read from the corresponding file.
263 void AliTPCcalibDB::Update(){
265 AliCDBEntry * entry=0;
267 Bool_t cdbCache = AliCDBManager::Instance()->GetCacheFlag(); // save cache status
268 AliCDBManager::Instance()->SetCacheFlag(kTRUE); // activate CDB cache
271 entry = GetCDBEntry("TPC/Calib/PadGainFactor");
273 //if (fPadGainFactor) delete fPadGainFactor;
274 entry->SetOwner(kTRUE);
275 fPadGainFactor = (AliTPCCalPad*)entry->GetObject();
278 entry = GetCDBEntry("TPC/Calib/GainFactorDedx");
280 entry->SetOwner(kTRUE);
281 fDedxGainFactor = (AliTPCCalPad*)entry->GetObject();
284 entry = GetCDBEntry("TPC/Calib/PadTime0");
286 //if (fPadTime0) delete fPadTime0;
287 entry->SetOwner(kTRUE);
288 fPadTime0 = (AliTPCCalPad*)entry->GetObject();
292 entry = GetCDBEntry("TPC/Calib/PadNoise");
294 //if (fPadNoise) delete fPadNoise;
295 entry->SetOwner(kTRUE);
296 fPadNoise = (AliTPCCalPad*)entry->GetObject();
299 entry = GetCDBEntry("TPC/Calib/Pedestals");
301 //if (fPedestals) delete fPedestals;
302 entry->SetOwner(kTRUE);
303 fPedestals = (AliTPCCalPad*)entry->GetObject();
306 entry = GetCDBEntry("TPC/Calib/Temperature");
308 //if (fTemperature) delete fTemperature;
309 entry->SetOwner(kTRUE);
310 fTemperature = (AliTPCSensorTempArray*)entry->GetObject();
313 entry = GetCDBEntry("TPC/Calib/Parameters");
315 //if (fPadNoise) delete fPadNoise;
316 entry->SetOwner(kTRUE);
317 fParam = (AliTPCParam*)(entry->GetObject()->Clone());
320 entry = GetCDBEntry("TPC/Calib/ClusterParam");
322 //if (fPadNoise) delete fPadNoise;
323 entry->SetOwner(kTRUE);
324 fClusterParam = (AliTPCClusterParam*)(entry->GetObject()->Clone());
327 entry = GetCDBEntry("TPC/Calib/Mapping");
329 //if (fPadNoise) delete fPadNoise;
330 entry->SetOwner(kTRUE);
331 TObjArray * array = dynamic_cast<TObjArray*>(entry->GetObject());
332 if (array && array->GetEntriesFast()==6){
333 fMapping = new AliTPCAltroMapping*[6];
334 for (Int_t i=0; i<6; i++){
335 fMapping[i] = dynamic_cast<AliTPCAltroMapping*>(array->At(i));
342 //entry = GetCDBEntry("TPC/Calib/ExB");
344 // entry->SetOwner(kTRUE);
345 // fExB=dynamic_cast<AliTPCExB*>(entry->GetObject()->Clone());
348 // ExB - calculate during initialization
350 fExB = GetExB(-5,kTRUE);
353 fTransform=new AliTPCTransform();
357 AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
363 void AliTPCcalibDB::CreateObjectList(const Char_t *filename, TObjArray *calibObjects)
366 // Create calibration objects and read contents from OCDB
368 if ( calibObjects == 0x0 ) return;
371 if ( !in.is_open() ){
372 fprintf(stderr,"Error: cannot open list file '%s'", filename);
376 AliTPCCalPad *calPad=0x0;
382 TObjArray *arrFileLine = sFile.Tokenize("\n");
384 TIter nextLine(arrFileLine);
386 TObjString *sObjLine=0x0;
387 while ( (sObjLine = (TObjString*)nextLine()) ){
388 TString sLine(sObjLine->GetString());
390 TObjArray *arrNextCol = sLine.Tokenize("\t");
392 TObjString *sObjType = (TObjString*)(arrNextCol->At(0));
393 TObjString *sObjFileName = (TObjString*)(arrNextCol->At(1));
395 if ( !sObjType || ! sObjFileName ) continue;
396 TString sType(sObjType->GetString());
397 TString sFileName(sObjFileName->GetString());
398 printf("%s\t%s\n",sType.Data(),sFileName.Data());
400 TFile *fIn = TFile::Open(sFileName);
402 fprintf(stderr,"File not found: '%s'", sFileName.Data());
406 if ( sType == "CE" ){
407 AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
409 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());
410 calPad->SetNameTitle("CETmean","CETmean");
411 calibObjects->Add(calPad);
413 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());
414 calPad->SetNameTitle("CEQmean","CEQmean");
415 calibObjects->Add(calPad);
417 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
418 calPad->SetNameTitle("CETrms","CETrms");
419 calibObjects->Add(calPad);
421 } else if ( sType == "Pulser") {
422 AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
424 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());
425 calPad->SetNameTitle("PulserTmean","PulserTmean");
426 calibObjects->Add(calPad);
428 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());
429 calPad->SetNameTitle("PulserQmean","PulserQmean");
430 calibObjects->Add(calPad);
432 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
433 calPad->SetNameTitle("PulserTrms","PulserTrms");
434 calibObjects->Add(calPad);
436 } else if ( sType == "Pedestals") {
437 AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
439 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());
440 calPad->SetNameTitle("Pedestals","Pedestals");
441 calibObjects->Add(calPad);
443 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());
444 calPad->SetNameTitle("Noise","Noise");
445 calibObjects->Add(calPad);
448 fprintf(stderr,"Undefined Type: '%s'",sType.Data());
457 void AliTPCcalibDB::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
459 // Write a tree with all available information
460 // if mapFileName is specified, the Map information are also written to the tree
461 // pads specified in outlierPad are not used for calculating statistics
462 // - the same function as AliTPCCalPad::MakeTree -
464 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
466 TObjArray* mapIROCs = 0;
467 TObjArray* mapOROCs = 0;
468 TVectorF *mapIROCArray = 0;
469 TVectorF *mapOROCArray = 0;
470 Int_t mapEntries = 0;
471 TString* mapNames = 0;
474 TFile mapFile(mapFileName, "read");
476 TList* listOfROCs = mapFile.GetListOfKeys();
477 mapEntries = listOfROCs->GetEntries()/2;
478 mapIROCs = new TObjArray(mapEntries*2);
479 mapOROCs = new TObjArray(mapEntries*2);
480 mapIROCArray = new TVectorF[mapEntries];
481 mapOROCArray = new TVectorF[mapEntries];
483 mapNames = new TString[mapEntries];
484 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
485 TString nameROC(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
486 nameROC.Remove(nameROC.Length()-4, 4);
487 mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "IROC").Data()), ivalue);
488 mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "OROC").Data()), ivalue);
489 mapNames[ivalue].Append(nameROC);
492 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
493 mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
494 mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
496 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
497 (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
498 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
499 (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
502 } // if (mapFileName)
504 TTreeSRedirector cstream(fileName);
505 Int_t arrayEntries = array->GetEntries();
507 TString* names = new TString[arrayEntries];
508 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
509 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
511 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
513 // get statistic for given sector
515 TVectorF median(arrayEntries);
516 TVectorF mean(arrayEntries);
517 TVectorF rms(arrayEntries);
518 TVectorF ltm(arrayEntries);
519 TVectorF ltmrms(arrayEntries);
520 TVectorF medianWithOut(arrayEntries);
521 TVectorF meanWithOut(arrayEntries);
522 TVectorF rmsWithOut(arrayEntries);
523 TVectorF ltmWithOut(arrayEntries);
524 TVectorF ltmrmsWithOut(arrayEntries);
526 TVectorF *vectorArray = new TVectorF[arrayEntries];
527 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
528 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
530 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
531 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
532 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
533 AliTPCCalROC* outlierROC = 0;
534 if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
536 median[ivalue] = calROC->GetMedian();
537 mean[ivalue] = calROC->GetMean();
538 rms[ivalue] = calROC->GetRMS();
539 Double_t ltmrmsValue = 0;
540 ltm[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction);
541 ltmrms[ivalue] = ltmrmsValue;
543 medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
544 meanWithOut[ivalue] = calROC->GetMean(outlierROC);
545 rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
547 ltmWithOut[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction, outlierROC);
548 ltmrmsWithOut[ivalue] = ltmrmsValue;
557 medianWithOut[ivalue] = 0.;
558 meanWithOut[ivalue] = 0.;
559 rmsWithOut[ivalue] = 0.;
560 ltmWithOut[ivalue] = 0.;
561 ltmrmsWithOut[ivalue] = 0.;
566 // fill vectors of variable per pad
568 TVectorF *posArray = new TVectorF[8];
569 for (Int_t ivalue = 0; ivalue < 8; ivalue++)
570 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
572 Float_t posG[3] = {0};
573 Float_t posL[3] = {0};
575 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
576 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
577 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
578 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
579 posArray[0][ichannel] = irow;
580 posArray[1][ichannel] = ipad;
581 posArray[2][ichannel] = posL[0];
582 posArray[3][ichannel] = posL[1];
583 posArray[4][ichannel] = posG[0];
584 posArray[5][ichannel] = posG[1];
585 posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
586 posArray[7][ichannel] = ichannel;
588 // loop over array containing AliTPCCalPads
589 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
590 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
591 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
593 (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
595 (vectorArray[ivalue])[ichannel] = 0;
601 cstream << "calPads" <<
602 "sector=" << isector;
604 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
605 cstream << "calPads" <<
606 (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
607 (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
608 (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
609 (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
610 (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
612 cstream << "calPads" <<
613 (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
614 (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
615 (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
616 (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
617 (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
621 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
622 cstream << "calPads" <<
623 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
627 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
629 cstream << "calPads" <<
630 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
632 cstream << "calPads" <<
633 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
637 cstream << "calPads" <<
638 "row.=" << &posArray[0] <<
639 "pad.=" << &posArray[1] <<
640 "lx.=" << &posArray[2] <<
641 "ly.=" << &posArray[3] <<
642 "gx.=" << &posArray[4] <<
643 "gy.=" << &posArray[5] <<
644 "rpad.=" << &posArray[6] <<
645 "channel.=" << &posArray[7];
647 cstream << "calPads" <<
651 delete[] vectorArray;
659 delete[] mapIROCArray;
660 delete[] mapOROCArray;
667 void AliTPCcalibDB::RegisterExB(Int_t index, Float_t bz, Bool_t bdelete){
669 // Register static ExB correction map
670 // index - registration index - used for visualization
671 // bz - bz field in kGaus
673 Float_t factor = bz/(-5.); // default b filed in Cheb with minus sign
675 AliMagF* bmap = new AliMagWrapCheb("Maps","Maps", 2, factor, 10., AliMagWrapCheb::k5kG,kTRUE,"$(ALICE_ROOT)/data/maps/mfchebKGI_sym.root");
677 AliTPCExBFirst *exb = new AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
678 AliTPCExB::SetInstance(exb);
683 AliTPCExB::RegisterField(index,bmap);
685 if (index>=fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
686 fgExBArray.AddAt(exb,index);
690 AliTPCExB* AliTPCcalibDB::GetExB(Float_t bz, Bool_t deleteB) {
692 // bz filed in KGaus not in tesla
693 // Get ExB correction map
694 // if doesn't exist - create it
696 Int_t index = TMath::Nint(5+bz);
697 if (index>fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
698 if (!fgExBArray.At(index)) AliTPCcalibDB::RegisterExB(index,bz,deleteB);
699 return (AliTPCExB*)fgExBArray.At(index);
703 void AliTPCcalibDB::SetExBField(Float_t bz){
705 // Set magnetic filed for ExB correction
707 printf("Set magnetic field for ExB correction = %f\n",bz);
708 fExB = GetExB(bz,kFALSE);