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 ////
29 // Simulation - not yet
30 // Reconstruction - AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
32 // 1.) pad by pad calibration - AliTPCCalPad
35 // Simulation: AliTPCDigitizer::ExecFast - Multiply by gain
36 // Reconstruction : AliTPCclustererMI::Digits2Clusters - Divide by gain
39 // Simulation: AliTPCDigitizer::ExecFast
40 // Reconstruction: AliTPCclustererMI::FindClusters(AliTPCCalROC * noiseROC)
41 // Noise depending cut on clusters charge (n sigma)
43 // Simulation: Not used yet - To be impleneted - Rounding to the nearest integer
44 // Reconstruction: Used in AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
45 // if data taken without zero suppression
46 // Currently switch in fRecoParam->GetCalcPedestal();
49 // Simulation: applied in the AliTPC::MakeSector - adding offset
50 // Reconstruction: AliTPCTransform::Transform() - remove offset
51 // AliTPCTransform::Transform() - to be called
52 // in AliTPCtracker::Transform()
55 // 2.) Space points transformation:
57 // a.) General coordinate tranformation - AliTPCtransform (see $ALICE_ROOT/TPC/AliTPCtransform.cxx)
58 // Created on fly - use the other calibration components
59 // Unisochronity - (substract time0 - pad by pad)
60 // Drift velocity - Currently common drift velocity - functionality of AliTPCParam
62 // Simulation - Not used directly (the effects are applied one by one (see AliTPC::MakeSector)
64 // AliTPCclustererMI::AddCluster
65 // AliTPCtrackerMI::Transform
66 // b.) ExB effect calibration -
67 // classes (base class AliTPCExB, implementation- AliTPCExBExact.h AliTPCExBFirst.h)
68 // a.a) Simulation: applied in the AliTPC::MakeSector -
69 // calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
70 // a.b) Reconstruction -
72 // in AliTPCtransform::Correct() - called calib->GetExB()->Correct(dxyz0,dxyz1)
75 ///////////////////////////////////////////////////////////////////////////////
81 #include <AliCDBManager.h>
82 #include <AliCDBEntry.h>
85 #include "AliTPCcalibDB.h"
86 #include "AliTPCAltroMapping.h"
87 #include "AliTPCExB.h"
89 #include "AliTPCCalROC.h"
90 #include "AliTPCCalPad.h"
91 #include "AliTPCSensorTempArray.h"
92 #include "AliTPCTransform.h"
102 #include "TObjArray.h"
103 #include "TObjString.h"
105 #include "AliTPCCalPad.h"
106 #include "AliTPCCalibPulser.h"
107 #include "AliTPCCalibPedestal.h"
108 #include "AliTPCCalibCE.h"
114 ClassImp(AliTPCcalibDB)
116 AliTPCcalibDB* AliTPCcalibDB::fgInstance = 0;
117 Bool_t AliTPCcalibDB::fgTerminated = kFALSE;
120 //_ singleton implementation __________________________________________________
121 AliTPCcalibDB* AliTPCcalibDB::Instance()
124 // Singleton implementation
125 // Returns an instance of this class, it is created if neccessary
128 if (fgTerminated != kFALSE)
132 fgInstance = new AliTPCcalibDB();
137 void AliTPCcalibDB::Terminate()
140 // Singleton implementation
141 // Deletes the instance of this class and sets the terminated flag, instances cannot be requested anymore
142 // This function can be called several times.
145 fgTerminated = kTRUE;
154 //_____________________________________________________________________________
155 AliTPCcalibDB::AliTPCcalibDB():
171 Update(); // temporary
174 //_____________________________________________________________________________
175 AliTPCcalibDB::~AliTPCcalibDB()
181 // don't delete anything, CDB cache is active!
182 //if (fPadGainFactor) delete fPadGainFactor;
183 //if (fPadTime0) delete fPadTime0;
184 //if (fPadNoise) delete fPadNoise;
188 //_____________________________________________________________________________
189 AliCDBEntry* AliTPCcalibDB::GetCDBEntry(const char* cdbPath)
192 // Retrieves an entry with path <cdbPath> from the CDB.
196 AliCDBEntry* entry = AliCDBManager::Instance()->Get(cdbPath, fRun);
199 sprintf(chinfo,"AliTPCcalibDB: Failed to get entry:\t%s ", cdbPath);
207 //_____________________________________________________________________________
208 void AliTPCcalibDB::SetRun(Long64_t run)
211 // Sets current run number. Calibration data is read from the corresponding file.
221 void AliTPCcalibDB::Update(){
223 AliCDBEntry * entry=0;
225 Bool_t cdbCache = AliCDBManager::Instance()->GetCacheFlag(); // save cache status
226 AliCDBManager::Instance()->SetCacheFlag(kTRUE); // activate CDB cache
229 entry = GetCDBEntry("TPC/Calib/PadGainFactor");
231 //if (fPadGainFactor) delete fPadGainFactor;
232 entry->SetOwner(kTRUE);
233 fPadGainFactor = (AliTPCCalPad*)entry->GetObject();
236 entry = GetCDBEntry("TPC/Calib/PadTime0");
238 //if (fPadTime0) delete fPadTime0;
239 entry->SetOwner(kTRUE);
240 fPadTime0 = (AliTPCCalPad*)entry->GetObject();
244 entry = GetCDBEntry("TPC/Calib/PadNoise");
246 //if (fPadNoise) delete fPadNoise;
247 entry->SetOwner(kTRUE);
248 fPadNoise = (AliTPCCalPad*)entry->GetObject();
251 entry = GetCDBEntry("TPC/Calib/Pedestals");
253 //if (fPedestals) delete fPedestals;
254 entry->SetOwner(kTRUE);
255 fPedestals = (AliTPCCalPad*)entry->GetObject();
258 entry = GetCDBEntry("TPC/Calib/Temperature");
260 //if (fTemperature) delete fTemperature;
261 entry->SetOwner(kTRUE);
262 fTemperature = (AliTPCSensorTempArray*)entry->GetObject();
265 entry = GetCDBEntry("TPC/Calib/Parameters");
267 //if (fPadNoise) delete fPadNoise;
268 entry->SetOwner(kTRUE);
269 fParam = (AliTPCParam*)(entry->GetObject()->Clone());
272 entry = GetCDBEntry("TPC/Calib/Mapping");
274 //if (fPadNoise) delete fPadNoise;
275 entry->SetOwner(kTRUE);
276 TObjArray * array = dynamic_cast<TObjArray*>(entry->GetObject());
277 if (array && array->GetEntriesFast()==6){
278 fMapping = new AliTPCAltroMapping*[6];
279 for (Int_t i=0; i<6; i++){
280 fMapping[i] = dynamic_cast<AliTPCAltroMapping*>(array->At(i));
287 entry = GetCDBEntry("TPC/Calib/ExB");
289 entry->SetOwner(kTRUE);
290 fExB=dynamic_cast<AliTPCExB*>(entry->GetObject()->Clone());
294 fTransform=new AliTPCTransform();
298 AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
301 AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& org)
304 // Copy constructor invalid -- singleton implementation
306 Error("copy constructor","invalid -- singleton implementation");
309 AliTPCcalibDB& AliTPCcalibDB::operator= (const AliTPCcalibDB& rhs)
312 // Singleton implementation - no assignment operator
314 Error("operator =", "assignment operator not implemented");
320 void AliTPCcalibDB::CreateObjectList(const Char_t *filename, TObjArray *calibObjects)
323 // Create calibration objects and read contents from OCDB
325 if ( calibObjects == 0x0 ) return;
328 if ( !in.is_open() ){
329 fprintf(stderr,"Error: cannot open list file '%s'", filename);
333 AliTPCCalPad *calPad=0x0;
339 TObjArray *arrFileLine = sFile.Tokenize("\n");
341 TIter nextLine(arrFileLine);
343 TObjString *sObjLine=0x0;
344 while ( sObjLine = (TObjString*)nextLine() ){
345 TString sLine(sObjLine->GetString());
347 TObjArray *arrNextCol = sLine.Tokenize("\t");
349 TObjString *sObjType = (TObjString*)(arrNextCol->At(0));
350 TObjString *sObjFileName = (TObjString*)(arrNextCol->At(1));
352 if ( !sObjType || ! sObjFileName ) continue;
353 TString sType(sObjType->GetString());
354 TString sFileName(sObjFileName->GetString());
355 printf("%s\t%s\n",sType.Data(),sFileName.Data());
357 TFile *fIn = TFile::Open(sFileName);
359 fprintf(stderr,"File not found: '%s'", sFileName.Data());
363 if ( sType == "CE" ){
364 AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
366 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());
367 calPad->SetNameTitle("CETmean","CETmean");
368 calibObjects->Add(calPad);
370 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());
371 calPad->SetNameTitle("CEQmean","CEQmean");
372 calibObjects->Add(calPad);
374 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
375 calPad->SetNameTitle("CETrms","CETrms");
376 calibObjects->Add(calPad);
378 } else if ( sType == "Pulser") {
379 AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
381 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());
382 calPad->SetNameTitle("PulserTmean","PulserTmean");
383 calibObjects->Add(calPad);
385 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());
386 calPad->SetNameTitle("PulserQmean","PulserQmean");
387 calibObjects->Add(calPad);
389 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
390 calPad->SetNameTitle("PulserTrms","PulserTrms");
391 calibObjects->Add(calPad);
393 } else if ( sType == "Pedestals") {
394 AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
396 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());
397 calPad->SetNameTitle("Pedestals","Pedestals");
398 calibObjects->Add(calPad);
400 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());
401 calPad->SetNameTitle("Noise","Noise");
402 calibObjects->Add(calPad);
405 fprintf(stderr,"Undefined Type: '%s'",sType.Data());
414 void AliTPCcalibDB::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
416 // Write a tree with all available information
417 // if mapFileName is specified, the Map information are also written to the tree
418 // pads specified in outlierPad are not used for calculating statistics
419 // - the same function as AliTPCCalPad::MakeTree -
421 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
423 TObjArray* mapIROCs = 0;
424 TObjArray* mapOROCs = 0;
425 TVectorF *mapIROCArray = 0;
426 TVectorF *mapOROCArray = 0;
427 Int_t mapEntries = 0;
428 TString* mapNames = 0;
431 TFile mapFile(mapFileName, "read");
433 TList* listOfROCs = mapFile.GetListOfKeys();
434 mapEntries = listOfROCs->GetEntries()/2;
435 mapIROCs = new TObjArray(mapEntries*2);
436 mapOROCs = new TObjArray(mapEntries*2);
437 mapIROCArray = new TVectorF[mapEntries];
438 mapOROCArray = new TVectorF[mapEntries];
440 mapNames = new TString[mapEntries];
441 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
442 TString nameROC(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
443 nameROC.Remove(nameROC.Length()-4, 4);
444 mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "IROC").Data()), ivalue);
445 mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "OROC").Data()), ivalue);
446 mapNames[ivalue].Append(nameROC);
449 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
450 mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
451 mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
453 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
454 (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
455 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
456 (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
459 } // if (mapFileName)
461 TTreeSRedirector cstream(fileName);
462 Int_t arrayEntries = array->GetEntries();
464 TString* names = new TString[arrayEntries];
465 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
466 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
468 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
470 // get statistic for given sector
472 TVectorF median(arrayEntries);
473 TVectorF mean(arrayEntries);
474 TVectorF rms(arrayEntries);
475 TVectorF ltm(arrayEntries);
476 TVectorF ltmrms(arrayEntries);
477 TVectorF medianWithOut(arrayEntries);
478 TVectorF meanWithOut(arrayEntries);
479 TVectorF rmsWithOut(arrayEntries);
480 TVectorF ltmWithOut(arrayEntries);
481 TVectorF ltmrmsWithOut(arrayEntries);
483 TVectorF *vectorArray = new TVectorF[arrayEntries];
484 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
485 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
487 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
488 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
489 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
490 AliTPCCalROC* outlierROC = 0;
491 if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
493 median[ivalue] = calROC->GetMedian();
494 mean[ivalue] = calROC->GetMean();
495 rms[ivalue] = calROC->GetRMS();
496 Double_t ltmrmsValue = 0;
497 ltm[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction);
498 ltmrms[ivalue] = ltmrmsValue;
500 medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
501 meanWithOut[ivalue] = calROC->GetMean(outlierROC);
502 rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
504 ltmWithOut[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction, outlierROC);
505 ltmrmsWithOut[ivalue] = ltmrmsValue;
514 medianWithOut[ivalue] = 0.;
515 meanWithOut[ivalue] = 0.;
516 rmsWithOut[ivalue] = 0.;
517 ltmWithOut[ivalue] = 0.;
518 ltmrmsWithOut[ivalue] = 0.;
523 // fill vectors of variable per pad
525 TVectorF *posArray = new TVectorF[8];
526 for (Int_t ivalue = 0; ivalue < 8; ivalue++)
527 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
529 Float_t posG[3] = {0};
530 Float_t posL[3] = {0};
532 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
533 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
534 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
535 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
536 posArray[0][ichannel] = irow;
537 posArray[1][ichannel] = ipad;
538 posArray[2][ichannel] = posL[0];
539 posArray[3][ichannel] = posL[1];
540 posArray[4][ichannel] = posG[0];
541 posArray[5][ichannel] = posG[1];
542 posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
543 posArray[7][ichannel] = ichannel;
545 // loop over array containing AliTPCCalPads
546 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
547 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
548 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
550 (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
552 (vectorArray[ivalue])[ichannel] = 0;
558 cstream << "calPads" <<
559 "sector=" << isector;
561 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
562 cstream << "calPads" <<
563 (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
564 (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
565 (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
566 (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
567 (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
569 cstream << "calPads" <<
570 (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
571 (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
572 (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
573 (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
574 (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
578 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
579 cstream << "calPads" <<
580 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
584 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
586 cstream << "calPads" <<
587 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
589 cstream << "calPads" <<
590 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
594 cstream << "calPads" <<
595 "row.=" << &posArray[0] <<
596 "pad.=" << &posArray[1] <<
597 "lx.=" << &posArray[2] <<
598 "ly.=" << &posArray[3] <<
599 "gx.=" << &posArray[4] <<
600 "gy.=" << &posArray[5] <<
601 "rpad.=" << &posArray[6] <<
602 "channel.=" << &posArray[7];
604 cstream << "calPads" <<
608 delete[] vectorArray;
616 delete[] mapIROCArray;
617 delete[] mapOROCArray;