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 <AliSplineFit.h>
89 #include "AliTPCcalibDB.h"
90 #include "AliTPCcalibDButil.h"
91 #include "AliTPCAltroMapping.h"
92 #include "AliTPCExB.h"
94 #include "AliTPCCalROC.h"
95 #include "AliTPCCalPad.h"
96 #include "AliTPCSensorTempArray.h"
97 #include "AliGRPObject.h"
98 #include "AliTPCTransform.h"
107 #include "TGraphErrors.h"
109 #include "TObjArray.h"
110 #include "TObjString.h"
112 #include "TDirectory.h"
113 #include "AliTPCCalPad.h"
114 #include "AliTPCCalibPulser.h"
115 #include "AliTPCCalibPedestal.h"
116 #include "AliTPCCalibCE.h"
117 #include "AliTPCExBFirst.h"
118 #include "AliTPCTempMap.h"
119 #include "AliTPCCalibVdrift.h"
120 #include "AliTPCCalibRaw.h"
121 #include "AliTPCParam.h"
123 #include "AliTPCPreprocessorOnline.h"
126 ClassImp(AliTPCcalibDB)
128 AliTPCcalibDB* AliTPCcalibDB::fgInstance = 0;
129 Bool_t AliTPCcalibDB::fgTerminated = kFALSE;
130 TObjArray AliTPCcalibDB::fgExBArray; // array of ExB corrections
133 //_ singleton implementation __________________________________________________
134 AliTPCcalibDB* AliTPCcalibDB::Instance()
137 // Singleton implementation
138 // Returns an instance of this class, it is created if neccessary
141 if (fgTerminated != kFALSE)
145 fgInstance = new AliTPCcalibDB();
150 void AliTPCcalibDB::Terminate()
153 // Singleton implementation
154 // Deletes the instance of this class and sets the terminated flag, instances cannot be requested anymore
155 // This function can be called several times.
158 fgTerminated = kTRUE;
167 //_____________________________________________________________________________
168 AliTPCcalibDB::AliTPCcalibDB():
187 fTimeGainSplinesArray(100000),
188 fGRPArray(100000), //! array of GRPs - per run - JUST for calibration studies
189 fGRPMaps(100000), //! array of GRPs - per run - JUST for calibration studies
190 fGoofieArray(100000), //! array of GOOFIE values -per run - Just for calibration studies
191 fVoltageArray(100000),
192 fTemperatureArray(100000), //! array of temperature sensors - per run - Just for calibration studies
193 fVdriftArray(100000), //! array of v drift interfaces
194 fDriftCorrectionArray(100000), //! array of drift correction
195 fRunList(100000) //! run list - indicates try to get the run param
202 Update(); // temporary
205 AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& ):
224 fTimeGainSplinesArray(100000),
225 fGRPArray(0), //! array of GRPs - per run - JUST for calibration studies
226 fGRPMaps(0), //! array of GRPs - per run - JUST for calibration studies
227 fGoofieArray(0), //! array of GOOFIE values -per run - Just for calibration studies
229 fTemperatureArray(0), //! array of temperature sensors - per run - Just for calibration studies
230 fVdriftArray(0), //! array of v drift interfaces
231 fDriftCorrectionArray(0), //! array of v drift interfaces
232 fRunList(0) //! run list - indicates try to get the run param
235 // Copy constructor invalid -- singleton implementation
237 Error("copy constructor","invalid -- singleton implementation");
240 AliTPCcalibDB& AliTPCcalibDB::operator= (const AliTPCcalibDB& )
243 // Singleton implementation - no assignment operator
245 Error("operator =", "assignment operator not implemented");
251 //_____________________________________________________________________________
252 AliTPCcalibDB::~AliTPCcalibDB()
258 // don't delete anything, CDB cache is active!
259 //if (fPadGainFactor) delete fPadGainFactor;
260 //if (fPadTime0) delete fPadTime0;
261 //if (fPadNoise) delete fPadNoise;
265 //_____________________________________________________________________________
266 AliCDBEntry* AliTPCcalibDB::GetCDBEntry(const char* cdbPath)
269 // Retrieves an entry with path <cdbPath> from the CDB.
273 AliCDBEntry* entry = AliCDBManager::Instance()->Get(cdbPath, fRun);
276 sprintf(chinfo,"AliTPCcalibDB: Failed to get entry:\t%s ", cdbPath);
284 //_____________________________________________________________________________
285 void AliTPCcalibDB::SetRun(Long64_t run)
288 // Sets current run number. Calibration data is read from the corresponding file.
298 void AliTPCcalibDB::Update(){
300 AliCDBEntry * entry=0;
301 Bool_t cdbCache = AliCDBManager::Instance()->GetCacheFlag(); // save cache status
302 AliCDBManager::Instance()->SetCacheFlag(kTRUE); // activate CDB cache
305 entry = GetCDBEntry("TPC/Calib/PadGainFactor");
307 //if (fPadGainFactor) delete fPadGainFactor;
308 entry->SetOwner(kTRUE);
309 fPadGainFactor = (AliTPCCalPad*)entry->GetObject();
312 entry = GetCDBEntry("TPC/Calib/TimeGain");
314 //if (fTimeGainSplines) delete fTimeGainSplines;
315 entry->SetOwner(kTRUE);
316 fTimeGainSplines = (TObjArray*)entry->GetObject();
319 entry = GetCDBEntry("TPC/Calib/GainFactorDedx");
321 entry->SetOwner(kTRUE);
322 fDedxGainFactor = (AliTPCCalPad*)entry->GetObject();
325 entry = GetCDBEntry("TPC/Calib/PadTime0");
327 //if (fPadTime0) delete fPadTime0;
328 entry->SetOwner(kTRUE);
329 fPadTime0 = (AliTPCCalPad*)entry->GetObject();
333 entry = GetCDBEntry("TPC/Calib/PadNoise");
335 //if (fPadNoise) delete fPadNoise;
336 entry->SetOwner(kTRUE);
337 fPadNoise = (AliTPCCalPad*)entry->GetObject();
340 entry = GetCDBEntry("TPC/Calib/Pedestals");
342 //if (fPedestals) delete fPedestals;
343 entry->SetOwner(kTRUE);
344 fPedestals = (AliTPCCalPad*)entry->GetObject();
347 entry = GetCDBEntry("TPC/Calib/Temperature");
349 //if (fTemperature) delete fTemperature;
350 entry->SetOwner(kTRUE);
351 fTemperature = (AliTPCSensorTempArray*)entry->GetObject();
354 entry = GetCDBEntry("TPC/Calib/Parameters");
356 //if (fPadNoise) delete fPadNoise;
357 entry->SetOwner(kTRUE);
358 fParam = (AliTPCParam*)(entry->GetObject()->Clone());
361 entry = GetCDBEntry("TPC/Calib/ClusterParam");
363 entry->SetOwner(kTRUE);
364 fClusterParam = (AliTPCClusterParam*)(entry->GetObject()->Clone());
367 //ALTRO configuration data
368 entry = GetCDBEntry("TPC/Calib/AltroConfig");
370 entry->SetOwner(kTRUE);
371 fALTROConfigData=(TObjArray*)(entry->GetObject());
374 //Calibration Pulser data
375 entry = GetCDBEntry("TPC/Calib/Pulser");
377 entry->SetOwner(kTRUE);
378 fPulserData=(TObjArray*)(entry->GetObject());
382 entry = GetCDBEntry("TPC/Calib/CE");
384 entry->SetOwner(kTRUE);
385 fCEData=(TObjArray*)(entry->GetObject());
387 //RAW calibration data
388 entry = GetCDBEntry("TPC/Calib/Raw");
390 entry->SetOwner(kTRUE);
391 TObjArray *arr=(TObjArray*)(entry->GetObject());
392 if (arr) fCalibRaw=(AliTPCCalibRaw*)arr->At(0);
395 entry = GetCDBEntry("TPC/Calib/Mapping");
397 //if (fPadNoise) delete fPadNoise;
398 entry->SetOwner(kTRUE);
399 TObjArray * array = dynamic_cast<TObjArray*>(entry->GetObject());
400 if (array && array->GetEntriesFast()==6){
401 fMapping = new AliTPCAltroMapping*[6];
402 for (Int_t i=0; i<6; i++){
403 fMapping[i] = dynamic_cast<AliTPCAltroMapping*>(array->At(i));
410 //entry = GetCDBEntry("TPC/Calib/ExB");
412 // entry->SetOwner(kTRUE);
413 // fExB=dynamic_cast<AliTPCExB*>(entry->GetObject()->Clone());
416 // ExB - calculate during initialization - in simulation /reconstruction
417 // - not invoked here anymore
418 //fExB = GetExB(-5,kTRUE);
421 fTransform=new AliTPCTransform();
422 fTransform->SetCurrentRun(AliCDBManager::Instance()->GetRun());
426 AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
431 void AliTPCcalibDB::CreateObjectList(const Char_t *filename, TObjArray *calibObjects)
434 // Create calibration objects and read contents from OCDB
436 if ( calibObjects == 0x0 ) return;
439 if ( !in.is_open() ){
440 fprintf(stderr,"Error: cannot open list file '%s'", filename);
444 AliTPCCalPad *calPad=0x0;
450 TObjArray *arrFileLine = sFile.Tokenize("\n");
452 TIter nextLine(arrFileLine);
454 TObjString *sObjLine=0x0;
455 while ( (sObjLine = (TObjString*)nextLine()) ){
456 TString sLine(sObjLine->GetString());
458 TObjArray *arrNextCol = sLine.Tokenize("\t");
460 TObjString *sObjType = (TObjString*)(arrNextCol->At(0));
461 TObjString *sObjFileName = (TObjString*)(arrNextCol->At(1));
463 if ( !sObjType || ! sObjFileName ) continue;
464 TString sType(sObjType->GetString());
465 TString sFileName(sObjFileName->GetString());
466 printf("%s\t%s\n",sType.Data(),sFileName.Data());
468 TFile *fIn = TFile::Open(sFileName);
470 fprintf(stderr,"File not found: '%s'", sFileName.Data());
474 if ( sType == "CE" ){
475 AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
477 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());
478 calPad->SetNameTitle("CETmean","CETmean");
479 calibObjects->Add(calPad);
481 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());
482 calPad->SetNameTitle("CEQmean","CEQmean");
483 calibObjects->Add(calPad);
485 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
486 calPad->SetNameTitle("CETrms","CETrms");
487 calibObjects->Add(calPad);
489 } else if ( sType == "Pulser") {
490 AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
492 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());
493 calPad->SetNameTitle("PulserTmean","PulserTmean");
494 calibObjects->Add(calPad);
496 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());
497 calPad->SetNameTitle("PulserQmean","PulserQmean");
498 calibObjects->Add(calPad);
500 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
501 calPad->SetNameTitle("PulserTrms","PulserTrms");
502 calibObjects->Add(calPad);
504 } else if ( sType == "Pedestals") {
505 AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
507 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());
508 calPad->SetNameTitle("Pedestals","Pedestals");
509 calibObjects->Add(calPad);
511 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());
512 calPad->SetNameTitle("Noise","Noise");
513 calibObjects->Add(calPad);
516 fprintf(stderr,"Undefined Type: '%s'",sType.Data());
525 void AliTPCcalibDB::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
527 // Write a tree with all available information
528 // if mapFileName is specified, the Map information are also written to the tree
529 // pads specified in outlierPad are not used for calculating statistics
530 // - the same function as AliTPCCalPad::MakeTree -
532 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
534 TObjArray* mapIROCs = 0;
535 TObjArray* mapOROCs = 0;
536 TVectorF *mapIROCArray = 0;
537 TVectorF *mapOROCArray = 0;
538 Int_t mapEntries = 0;
539 TString* mapNames = 0;
542 TFile mapFile(mapFileName, "read");
544 TList* listOfROCs = mapFile.GetListOfKeys();
545 mapEntries = listOfROCs->GetEntries()/2;
546 mapIROCs = new TObjArray(mapEntries*2);
547 mapOROCs = new TObjArray(mapEntries*2);
548 mapIROCArray = new TVectorF[mapEntries];
549 mapOROCArray = new TVectorF[mapEntries];
551 mapNames = new TString[mapEntries];
552 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
553 TString nameROC(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
554 nameROC.Remove(nameROC.Length()-4, 4);
555 mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "IROC").Data()), ivalue);
556 mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "OROC").Data()), ivalue);
557 mapNames[ivalue].Append(nameROC);
560 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
561 mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
562 mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
564 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
565 (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
566 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
567 (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
570 } // if (mapFileName)
572 TTreeSRedirector cstream(fileName);
573 Int_t arrayEntries = array->GetEntries();
575 TString* names = new TString[arrayEntries];
576 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
577 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
579 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
581 // get statistic for given sector
583 TVectorF median(arrayEntries);
584 TVectorF mean(arrayEntries);
585 TVectorF rms(arrayEntries);
586 TVectorF ltm(arrayEntries);
587 TVectorF ltmrms(arrayEntries);
588 TVectorF medianWithOut(arrayEntries);
589 TVectorF meanWithOut(arrayEntries);
590 TVectorF rmsWithOut(arrayEntries);
591 TVectorF ltmWithOut(arrayEntries);
592 TVectorF ltmrmsWithOut(arrayEntries);
594 TVectorF *vectorArray = new TVectorF[arrayEntries];
595 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
596 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
598 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
599 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
600 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
601 AliTPCCalROC* outlierROC = 0;
602 if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
604 median[ivalue] = calROC->GetMedian();
605 mean[ivalue] = calROC->GetMean();
606 rms[ivalue] = calROC->GetRMS();
607 Double_t ltmrmsValue = 0;
608 ltm[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction);
609 ltmrms[ivalue] = ltmrmsValue;
611 medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
612 meanWithOut[ivalue] = calROC->GetMean(outlierROC);
613 rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
615 ltmWithOut[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction, outlierROC);
616 ltmrmsWithOut[ivalue] = ltmrmsValue;
625 medianWithOut[ivalue] = 0.;
626 meanWithOut[ivalue] = 0.;
627 rmsWithOut[ivalue] = 0.;
628 ltmWithOut[ivalue] = 0.;
629 ltmrmsWithOut[ivalue] = 0.;
634 // fill vectors of variable per pad
636 TVectorF *posArray = new TVectorF[8];
637 for (Int_t ivalue = 0; ivalue < 8; ivalue++)
638 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
640 Float_t posG[3] = {0};
641 Float_t posL[3] = {0};
643 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
644 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
645 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
646 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
647 posArray[0][ichannel] = irow;
648 posArray[1][ichannel] = ipad;
649 posArray[2][ichannel] = posL[0];
650 posArray[3][ichannel] = posL[1];
651 posArray[4][ichannel] = posG[0];
652 posArray[5][ichannel] = posG[1];
653 posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
654 posArray[7][ichannel] = ichannel;
656 // loop over array containing AliTPCCalPads
657 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
658 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
659 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
661 (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
663 (vectorArray[ivalue])[ichannel] = 0;
669 cstream << "calPads" <<
670 "sector=" << isector;
672 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
673 cstream << "calPads" <<
674 (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
675 (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
676 (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
677 (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
678 (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
680 cstream << "calPads" <<
681 (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
682 (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
683 (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
684 (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
685 (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
689 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
690 cstream << "calPads" <<
691 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
695 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
697 cstream << "calPads" <<
698 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
700 cstream << "calPads" <<
701 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
705 cstream << "calPads" <<
706 "row.=" << &posArray[0] <<
707 "pad.=" << &posArray[1] <<
708 "lx.=" << &posArray[2] <<
709 "ly.=" << &posArray[3] <<
710 "gx.=" << &posArray[4] <<
711 "gy.=" << &posArray[5] <<
712 "rpad.=" << &posArray[6] <<
713 "channel.=" << &posArray[7];
715 cstream << "calPads" <<
719 delete[] vectorArray;
727 delete[] mapIROCArray;
728 delete[] mapOROCArray;
735 void AliTPCcalibDB::RegisterExB(Int_t index, Float_t bz, Bool_t bdelete){
737 // Register static ExB correction map
738 // index - registration index - used for visualization
739 // bz - bz field in kGaus
741 // Float_t factor = bz/(-5.); // default b filed in Cheb with minus sign
742 Float_t factor = bz/(5.); // default b filed in Cheb with minus sign
743 // was chenged in the Revision ???? (Ruben can you add here number)
745 AliMagF* bmap = new AliMagF("MapsExB","MapsExB", factor,TMath::Sign(1.f,factor),AliMagF::k5kG);
747 AliTPCExBFirst *exb = new AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
748 AliTPCExB::SetInstance(exb);
753 AliTPCExB::RegisterField(index,bmap);
755 if (index>=fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
756 fgExBArray.AddAt(exb,index);
760 AliTPCExB* AliTPCcalibDB::GetExB(Float_t bz, Bool_t deleteB) {
762 // bz filed in KGaus not in tesla
763 // Get ExB correction map
764 // if doesn't exist - create it
766 Int_t index = TMath::Nint(5+bz);
767 if (index>fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
768 if (!fgExBArray.At(index)) AliTPCcalibDB::RegisterExB(index,bz,deleteB);
769 return (AliTPCExB*)fgExBArray.At(index);
773 void AliTPCcalibDB::SetExBField(Float_t bz){
775 // Set magnetic filed for ExB correction
777 fExB = GetExB(bz,kFALSE);
780 void AliTPCcalibDB::SetExBField(const AliMagF* bmap){
782 // Set magnetic field for ExB correction
784 AliTPCExBFirst *exb = new AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
785 AliTPCExB::SetInstance(exb);
793 void AliTPCcalibDB::UpdateRunInformations( Int_t run, Bool_t force){
795 // - > Don't use it for reconstruction - Only for Calibration studies
797 AliCDBEntry * entry = 0;
798 if (run>= fRunList.GetSize()){
799 fRunList.Set(run*2+1);
800 fGRPArray.Expand(run*2+1);
801 fGRPMaps.Expand(run*2+1);
802 fGoofieArray.Expand(run*2+1);
803 fVoltageArray.Expand(run*2+1);
804 fTemperatureArray.Expand(run*2+1);
805 fVdriftArray.Expand(run*2+1);
806 fDriftCorrectionArray.Expand(run*2+1);
807 fTimeGainSplinesArray.Expand(run*2+1);
809 if (fRunList[run]>0 &&force==kFALSE) return;
811 entry = AliCDBManager::Instance()->Get("GRP/GRP/Data",run);
813 AliGRPObject * grpRun = dynamic_cast<AliGRPObject*>(entry->GetObject());
815 TMap* map = dynamic_cast<TMap*>(entry->GetObject());
817 //grpRun = new AliGRPObject;
818 //grpRun->ReadValuesFromMap(map);
819 grpRun = MakeGRPObjectFromMap(map);
821 fGRPMaps.AddAt(map,run);
824 fGRPArray.AddAt(grpRun,run);
826 entry = AliCDBManager::Instance()->Get("TPC/Calib/Goofie",run);
828 fGoofieArray.AddAt(entry->GetObject(),run);
831 entry = AliCDBManager::Instance()->Get("TPC/Calib/HighVoltage",run);
833 fVoltageArray.AddAt(entry->GetObject(),run);
836 entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeGain",run);
838 fTimeGainSplinesArray.AddAt(entry->GetObject(),run);
841 entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeDrift",run);
843 fDriftCorrectionArray.AddAt(entry->GetObject(),run);
846 entry = AliCDBManager::Instance()->Get("TPC/Calib/Temperature",run);
848 fTemperatureArray.AddAt(entry->GetObject(),run);
850 fRunList[run]=1; // sign as used
852 AliDCSSensor * press = GetPressureSensor(run,0);
853 AliTPCSensorTempArray * temp = GetTemperatureSensor(run);
855 AliTPCCalibVdrift * vdrift = new AliTPCCalibVdrift(temp, press,0);
856 fVdriftArray.AddAt(vdrift,run);
861 Float_t AliTPCcalibDB::GetGain(Int_t sector, Int_t row, Int_t pad){
864 AliTPCCalPad *calPad = Instance()->fDedxGainFactor;;
865 if (!calPad) return 0;
866 return calPad->GetCalROC(sector)->GetValue(row,pad);
869 AliSplineFit* AliTPCcalibDB::GetVdriftSplineFit(const char* name, Int_t run){
873 TObjArray *arr=GetTimeVdriftSplineRun(run);
875 return dynamic_cast<AliSplineFit*>(arr->FindObject(name));
878 AliSplineFit* AliTPCcalibDB::CreateVdriftSplineFit(const char* graphName, Int_t run){
880 // create spline fit from the drift time graph in TimeDrift
882 TObjArray *arr=GetTimeVdriftSplineRun(run);
884 TGraph *graph=dynamic_cast<TGraph*>(arr->FindObject(graphName));
885 if (!graph) return 0;
886 AliSplineFit *fit = new AliSplineFit();
887 fit->SetGraph(graph);
888 fit->SetMinPoints(graph->GetN()+1);
889 fit->InitKnots(graph,2,0,0.001);
894 AliGRPObject *AliTPCcalibDB::GetGRP(Int_t run){
896 // Get GRP object for given run
898 AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>((Instance()->fGRPArray).At(run));
900 Instance()->UpdateRunInformations(run);
901 grpRun = dynamic_cast<AliGRPObject *>(Instance()->fGRPArray.At(run));
902 if (!grpRun) return 0;
907 TMap * AliTPCcalibDB::GetGRPMap(Int_t run){
911 TMap * grpRun = dynamic_cast<TMap *>((Instance()->fGRPMaps).At(run));
913 Instance()->UpdateRunInformations(run);
914 grpRun = dynamic_cast<TMap *>(Instance()->fGRPMaps.At(run));
915 if (!grpRun) return 0;
921 AliDCSSensor * AliTPCcalibDB::GetPressureSensor(Int_t run, Int_t type){
923 // Get Pressure sensor
925 // type = 0 - Cavern pressure
926 // 1 - Suface pressure
927 // First try to get if trom map - if existing (Old format of data storing)
931 TMap *map = GetGRPMap(run);
933 AliDCSSensor * sensor = 0;
935 if (type==0) osensor = ((*map)("fCavernPressure"));
936 if (type==1) osensor = ((*map)("fP2Pressure"));
937 sensor =dynamic_cast<AliDCSSensor *>(osensor);
938 if (sensor) return sensor;
941 // If not map try to get it from the GRPObject
943 AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run));
945 UpdateRunInformations(run);
946 grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run));
947 if (!grpRun) return 0;
949 AliDCSSensor * sensor = grpRun->GetCavernAtmosPressure();
950 if (type==1) sensor = grpRun->GetSurfaceAtmosPressure();
954 AliTPCSensorTempArray * AliTPCcalibDB::GetTemperatureSensor(Int_t run){
956 // Get temperature sensor array
958 AliTPCSensorTempArray * tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
960 UpdateRunInformations(run);
961 tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
967 TObjArray * AliTPCcalibDB::GetTimeGainSplinesRun(Int_t run){
969 // Get temperature sensor array
971 TObjArray * gainSplines = (TObjArray *)fTimeGainSplinesArray.At(run);
973 UpdateRunInformations(run);
974 gainSplines = (TObjArray *)fTimeGainSplinesArray.At(run);
979 TObjArray * AliTPCcalibDB::GetTimeVdriftSplineRun(Int_t run){
981 // Get drift spline array
983 TObjArray * driftSplines = (TObjArray *)fDriftCorrectionArray.At(run);
985 UpdateRunInformations(run);
986 driftSplines = (TObjArray *)fDriftCorrectionArray.At(run);
991 AliDCSSensorArray * AliTPCcalibDB::GetVoltageSensors(Int_t run){
993 // Get temperature sensor array
995 AliDCSSensorArray * voltageArray = (AliDCSSensorArray *)fVoltageArray.At(run);
997 UpdateRunInformations(run);
998 voltageArray = (AliDCSSensorArray *)fVoltageArray.At(run);
1000 return voltageArray;
1003 AliDCSSensorArray * AliTPCcalibDB::GetGoofieSensors(Int_t run){
1005 // Get temperature sensor array
1007 AliDCSSensorArray * goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
1009 UpdateRunInformations(run);
1010 goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
1017 AliTPCCalibVdrift * AliTPCcalibDB::GetVdrift(Int_t run){
1019 // Get the interface to the the vdrift
1021 AliTPCCalibVdrift * vdrift = (AliTPCCalibVdrift*)fVdriftArray.At(run);
1023 UpdateRunInformations(run);
1024 vdrift= (AliTPCCalibVdrift*)fVdriftArray.At(run);
1029 Float_t AliTPCcalibDB::GetCEdriftTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1032 // GetCE drift time information for 'sector'
1033 // sector 72 is the mean drift time of the A-Side
1034 // sector 73 is the mean drift time of the C-Side
1035 // it timestamp==-1 return mean value
1037 AliTPCcalibDB::Instance()->SetRun(run);
1038 TGraph *gr=AliTPCcalibDB::Instance()->GetCErocTgraph(sector);
1039 if (!gr||sector<0||sector>73) {
1040 if (entries) *entries=0;
1044 if (timeStamp==-1.){
1047 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1049 gr->GetPoint(ipoint,x,y);
1050 if (x<timeStamp) continue;
1058 Float_t AliTPCcalibDB::GetCEchargeTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1061 // GetCE mean charge for 'sector'
1062 // it timestamp==-1 return mean value
1064 AliTPCcalibDB::Instance()->SetRun(run);
1065 TGraph *gr=AliTPCcalibDB::Instance()->GetCErocQgraph(sector);
1066 if (!gr||sector<0||sector>71) {
1067 if (entries) *entries=0;
1071 if (timeStamp==-1.){
1074 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1076 gr->GetPoint(ipoint,x,y);
1077 if (x<timeStamp) continue;
1085 Float_t AliTPCcalibDB::GetDCSSensorValue(AliDCSSensorArray *arr, Int_t timeStamp, const char * sensorName, Int_t sigDigits)
1088 // Get Value for a DCS sensor 'sensorName', run 'run' at time 'timeStamp'
1091 const TString sensorNameString(sensorName);
1092 AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1093 if (!sensor) return val;
1094 //use the dcs graph if possible
1095 TGraph *gr=sensor->GetGraph();
1097 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1099 gr->GetPoint(ipoint,x,y);
1100 Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1101 if (time<timeStamp) continue;
1105 //if val is still 0, test if if the requested time if within 5min of the first/last
1106 //data point. If this is the case return the firs/last entry
1107 //the timestamps might not be syncronised for all calibration types, sometimes a 'pre'
1108 //and 'pos' period is requested. Especially to the HV this is not the case!
1112 gr->GetPoint(0,x,y);
1113 Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1114 if ((time-timeStamp)<5*60) val=y;
1119 gr->GetPoint(gr->GetN()-1,x,y);
1120 Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1121 if ((timeStamp-time)<5*60) val=y;
1124 val=sensor->GetValue(timeStamp);
1127 val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1132 Float_t AliTPCcalibDB::GetDCSSensorMeanValue(AliDCSSensorArray *arr, const char * sensorName, Int_t sigDigits)
1135 // Get mean Value for a DCS sensor 'sensorName' during run 'run'
1138 const TString sensorNameString(sensorName);
1139 AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1140 if (!sensor) return val;
1142 //use dcs graph if it exists
1143 TGraph *gr=sensor->GetGraph();
1147 //if we don't have the dcs graph, try to get some meaningful information
1148 if (!sensor->GetFit()) return val;
1149 Int_t nKnots=sensor->GetFit()->GetKnots();
1150 Double_t tMid=(sensor->GetEndTime()-sensor->GetStartTime())/2.;
1151 for (Int_t iKnot=0;iKnot<nKnots;++iKnot){
1152 if (sensor->GetFit()->GetX()[iKnot]>tMid/3600.) break;
1153 val=(Float_t)sensor->GetFit()->GetY0()[iKnot];
1158 val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1164 Float_t AliTPCcalibDB::GetChamberHighVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits) {
1166 // return the chamber HV for given run and time: 0-35 IROC, 36-72 OROC
1167 // if timeStamp==-1 return mean value
1170 TString sensorName="";
1171 TTimeStamp stamp(timeStamp);
1172 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1173 if (!voltageArray || (sector<0) || (sector>71)) return val;
1174 Char_t sideName='A';
1175 if ((sector/18)%2==1) sideName='C';
1178 sensorName=Form("TPC_ANODE_I_%c%02d_VMEAS",sideName,sector%18);
1181 sensorName=Form("TPC_ANODE_O_%c%02d_0_VMEAS",sideName,sector%18);
1184 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1186 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1190 Float_t AliTPCcalibDB::GetSkirtVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1193 // Get the skirt voltage for 'run' at 'timeStamp' and 'sector': 0-35 IROC, 36-72 OROC
1194 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1195 // if timeStamp==-1 return the mean value for the run
1198 TString sensorName="";
1199 TTimeStamp stamp(timeStamp);
1200 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1201 if (!voltageArray || (sector<0) || (sector>71)) return val;
1202 Char_t sideName='A';
1203 if ((sector/18)%2==1) sideName='C';
1204 sensorName=Form("TPC_SKIRT_%c_VMEAS",sideName);
1206 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1208 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1213 Float_t AliTPCcalibDB::GetCoverVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1216 // Get the cover voltage for run 'run' at time 'timeStamp'
1217 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1218 // if timeStamp==-1 return the mean value for the run
1221 TString sensorName="";
1222 TTimeStamp stamp(timeStamp);
1223 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1224 if (!voltageArray || (sector<0) || (sector>71)) return val;
1225 Char_t sideName='A';
1226 if ((sector/18)%2==1) sideName='C';
1229 sensorName=Form("TPC_COVER_I_%c_VMEAS",sideName);
1232 sensorName=Form("TPC_COVER_O_%c_VMEAS",sideName);
1235 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1237 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1242 Float_t AliTPCcalibDB::GetGGoffsetVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1245 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1246 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1247 // if timeStamp==-1 return the mean value for the run
1250 TString sensorName="";
1251 TTimeStamp stamp(timeStamp);
1252 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1253 if (!voltageArray || (sector<0) || (sector>71)) return val;
1254 Char_t sideName='A';
1255 if ((sector/18)%2==1) sideName='C';
1258 sensorName=Form("TPC_GATE_I_%c_OFF_VMEAS",sideName);
1261 sensorName=Form("TPC_GATE_O_%c_OFF_VMEAS",sideName);
1264 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1266 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1271 Float_t AliTPCcalibDB::GetGGnegVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1274 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1275 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1276 // if timeStamp==-1 return the mean value for the run
1279 TString sensorName="";
1280 TTimeStamp stamp(timeStamp);
1281 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1282 if (!voltageArray || (sector<0) || (sector>71)) return val;
1283 Char_t sideName='A';
1284 if ((sector/18)%2==1) sideName='C';
1287 sensorName=Form("TPC_GATE_I_%c_NEG_VMEAS",sideName);
1290 sensorName=Form("TPC_GATE_O_%c_NEG_VMEAS",sideName);
1293 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1295 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1300 Float_t AliTPCcalibDB::GetGGposVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1303 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1304 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1305 // if timeStamp==-1 return the mean value for the run
1308 TString sensorName="";
1309 TTimeStamp stamp(timeStamp);
1310 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1311 if (!voltageArray || (sector<0) || (sector>71)) return val;
1312 Char_t sideName='A';
1313 if ((sector/18)%2==1) sideName='C';
1316 sensorName=Form("TPC_GATE_I_%c_POS_VMEAS",sideName);
1319 sensorName=Form("TPC_GATE_O_%c_POS_VMEAS",sideName);
1322 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1324 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1329 Float_t AliTPCcalibDB::GetPressure(Int_t timeStamp, Int_t run, Int_t type){
1331 // GetPressure for given time stamp and runt
1333 TTimeStamp stamp(timeStamp);
1334 AliDCSSensor * sensor = Instance()->GetPressureSensor(run,type);
1335 if (!sensor) return 0;
1336 return sensor->GetValue(stamp);
1339 Float_t AliTPCcalibDB::GetL3Current(Int_t run, Int_t statType){
1341 // return L3 current
1342 // stat type is: AliGRPObject::Stats: kMean = 0, kTruncMean = 1, kMedian = 2, kSDMean = 3, kSDMedian = 4
1345 AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1346 if (grp) current=grp->GetL3Current((AliGRPObject::Stats)statType);
1350 Float_t AliTPCcalibDB::GetBz(Int_t run){
1352 // calculate BZ in T from L3 current
1355 Float_t current=AliTPCcalibDB::GetL3Current(run);
1356 if (current>-1) bz=5*current/30000.*.1;
1360 Char_t AliTPCcalibDB::GetL3Polarity(Int_t run) {
1362 // get l3 polarity from GRP
1365 AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1366 if (grp) pol=grp->GetL3Polarity();
1370 TString AliTPCcalibDB::GetRunType(Int_t run){
1372 // return run type from grp
1375 // TString type("UNKNOWN");
1376 AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1377 if (grp) return grp->GetRunType();
1381 Float_t AliTPCcalibDB::GetValueGoofie(Int_t timeStamp, Int_t run, Int_t type){
1383 // GetPressure for given time stamp and runt
1385 TTimeStamp stamp(timeStamp);
1386 AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(run);
1387 if (!goofieArray) return 0;
1388 AliDCSSensor *sensor = goofieArray->GetSensor(type);
1389 return sensor->GetValue(stamp);
1397 Bool_t AliTPCcalibDB::GetTemperatureFit(Int_t timeStamp, Int_t run, Int_t side,TVectorD& fit){
1401 TTimeStamp tstamp(timeStamp);
1402 AliTPCSensorTempArray* tempArray = Instance()->GetTemperatureSensor(run);
1403 if (! tempArray) return kFALSE;
1404 AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
1405 TLinearFitter * fitter = tempMap->GetLinearFitter(3,side,tstamp);
1408 fitter->GetParameters(fit);
1412 if (!fitter) return kFALSE;
1416 Float_t AliTPCcalibDB::GetTemperature(Int_t timeStamp, Int_t run, Int_t side){
1422 GetTemperatureFit(timeStamp,run,0,vec);
1426 GetTemperatureFit(timeStamp,run,0,vec);
1433 Double_t AliTPCcalibDB::GetPTRelative(UInt_t timeSec, Int_t run, Int_t side){
1436 // time - absolute time
1438 // side - 0 - A side 1-C side
1439 AliTPCCalibVdrift * vdrift = Instance()->GetVdrift(run);
1440 if (!vdrift) return 0;
1441 return vdrift->GetPTRelative(timeSec,side);
1444 AliGRPObject * AliTPCcalibDB::MakeGRPObjectFromMap(TMap *map){
1446 // Function to covert old GRP run information from TMap to GRPObject
1448 // TMap * map = AliTPCcalibDB::GetGRPMap(52406);
1450 AliDCSSensor * sensor = 0;
1452 osensor = ((*map)("fP2Pressure"));
1453 sensor =dynamic_cast<AliDCSSensor *>(osensor);
1455 if (!sensor) return 0;
1457 AliDCSSensor * sensor2 = new AliDCSSensor(*sensor);
1458 osensor = ((*map)("fCavernPressure"));
1459 TGraph * gr = new TGraph(2);
1460 gr->GetX()[0]= -100000.;
1461 gr->GetX()[1]= 1000000.;
1462 gr->GetY()[0]= atof(osensor->GetName());
1463 gr->GetY()[1]= atof(osensor->GetName());
1464 sensor2->SetGraph(gr);
1468 AliGRPObject *grpRun = new AliGRPObject;
1469 grpRun->ReadValuesFromMap(map);
1470 grpRun->SetCavernAtmosPressure(sensor2);
1471 grpRun->SetSurfaceAtmosPressure(sensor);
1475 Bool_t AliTPCcalibDB::CreateGUITree(Int_t run, const char* filename)
1478 // Create a gui tree for run number 'run'
1481 if (!AliCDBManager::Instance()->GetDefaultStorage()){
1482 AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
1483 MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
1487 AliTPCcalibDB *db=AliTPCcalibDB::Instance();
1488 // retrieve cal pad objects
1490 AliTPCPreprocessorOnline prep;
1491 //noise and pedestals
1492 prep.AddComponent(db->GetPedestals());
1493 prep.AddComponent(db->GetPadNoise());
1495 prep.AddComponent(db->GetPulserTmean());
1496 prep.AddComponent(db->GetPulserTrms());
1497 prep.AddComponent(db->GetPulserQmean());
1499 prep.AddComponent(db->GetCETmean());
1500 prep.AddComponent(db->GetCETrms());
1501 prep.AddComponent(db->GetCEQmean());
1503 prep.AddComponent(db->GetALTROAcqStart() );
1504 prep.AddComponent(db->GetALTROZsThr() );
1505 prep.AddComponent(db->GetALTROFPED() );
1506 prep.AddComponent(db->GetALTROAcqStop() );
1507 prep.AddComponent(db->GetALTROMasked() );
1509 TString file(filename);
1510 if (file.IsNull()) file=Form("guiTreeRun_%d.root",run);
1511 prep.DumpToFile(file.Data());
1515 Bool_t AliTPCcalibDB::CreateRefFile(Int_t run, const char* filename)
1518 // Create a gui tree for run number 'run'
1521 if (!AliCDBManager::Instance()->GetDefaultStorage()){
1522 AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
1523 MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
1526 TString file(filename);
1527 if (file.IsNull()) file=Form("RefCalPads_%d.root",run);
1528 TDirectory *currDir=gDirectory;
1530 AliTPCcalibDB *db=AliTPCcalibDB::Instance();
1531 // retrieve cal pad objects
1534 TFile f(file.Data(),"recreate");
1535 //noise and pedestals
1536 db->GetPedestals()->Write("Pedestals");
1537 db->GetPadNoise()->Write("PadNoise");
1539 db->GetPulserTmean()->Write("PulserTmean");
1540 db->GetPulserTrms()->Write("PulserTrms");
1541 db->GetPulserQmean()->Write("PulserQmean");
1543 db->GetCETmean()->Write("CETmean");
1544 db->GetCETrms()->Write("CETrms");
1545 db->GetCEQmean()->Write("CEQmean");
1547 db->GetALTROAcqStart() ->Write("ALTROAcqStart");
1548 db->GetALTROZsThr() ->Write("ALTROZsThr");
1549 db->GetALTROFPED() ->Write("ALTROFPED");
1550 db->GetALTROAcqStop() ->Write("ALTROAcqStop");
1551 db->GetALTROMasked() ->Write("ALTROMasked");
1560 Double_t AliTPCcalibDB::GetVDriftCorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
1562 // Get time dependent drift velocity correction
1563 // multiplication factor vd = vdnom *(1+vdriftcorr)
1565 // mode determines the algorith how to combine the Laser Track, LaserCE and physics tracks
1566 // timestamp - timestamp
1568 // side - the drift velocity per side (possible for laser and CE)
1570 // Notice - Extrapolation outside of calibration range - using constant function
1573 // mode TPC crossing and laser
1575 result=AliTPCcalibDButil::GetVDriftTPC(run,timeStamp);
1582 Double_t AliTPCcalibDB::GetTime0CorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
1584 // Get time dependent time 0 (trigger delay in cm) correction
1585 // additive correction time0 = time0+ GetTime0CorrectionTime
1586 // Value etracted combining the vdrift correction using laser tracks and CE and the physics track matchin
1588 // mode determines the algorith how to combine the Laser Track and physics tracks
1589 // timestamp - timestamp
1591 // side - the drift velocity per side (possible for laser and CE)
1593 // Notice - Extrapolation outside of calibration range - using constant function
1596 if (mode==1) result=AliTPCcalibDButil::GetTriggerOffsetTPC(run,timeStamp);
1597 result *=fParam->GetZLength();
1606 Double_t AliTPCcalibDB::GetVDriftCorrectionGy(Int_t timeStamp, Int_t run, Int_t side, Int_t /*mode*/){
1608 // Get global y correction drift velocity correction factor
1609 // additive factor vd = vdnom*(1+GetVDriftCorrectionGy *gy)
1610 // Value etracted combining the vdrift correction using laser tracks and CE
1612 // mode determines the algorith how to combine the Laser Track, LaserCE
1613 // timestamp - timestamp
1615 // side - the drift velocity gy correction per side (CE and Laser tracks)
1617 // Notice - Extrapolation outside of calibration range - using constant function
1619 if (run<=0 && fTransform) run = fTransform->GetCurrentRunNumber();
1620 UpdateRunInformations(run,kFALSE);
1621 TObjArray *array =AliTPCcalibDB::Instance()->GetTimeVdriftSplineRun(run);
1622 if (!array) return 0;
1623 TGraphErrors *laserA= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
1624 TGraphErrors *laserC= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
1627 if (laserA && laserC){
1628 result= (laserA->Eval(timeStamp)+laserC->Eval(timeStamp))*0.5;
1630 if (laserA && side==0){
1631 result = (laserA->Eval(timeStamp));
1633 if (laserC &&side==1){
1634 result = (laserC->Eval(timeStamp));
1636 return -result/250.; //normalized before