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>
88 #include <AliSplineFit.h>
89 #include <AliCTPTimeParams.h>
91 #include "AliTPCcalibDB.h"
92 #include "AliTPCdataQA.h"
93 #include "AliTPCcalibDButil.h"
94 #include "AliTPCAltroMapping.h"
95 #include "AliTPCExB.h"
97 #include "AliTPCCalROC.h"
98 #include "AliTPCCalPad.h"
99 #include "AliTPCSensorTempArray.h"
100 #include "AliGRPObject.h"
101 #include "AliTPCTransform.h"
102 #include "AliTPCmapper.h"
111 #include "TGraphErrors.h"
113 #include "TObjArray.h"
114 #include "TObjString.h"
116 #include "TDirectory.h"
118 #include "AliTPCCalPad.h"
119 #include "AliTPCCalibPulser.h"
120 #include "AliTPCCalibPedestal.h"
121 #include "AliTPCCalibCE.h"
122 #include "AliTPCExBFirst.h"
123 #include "AliTPCTempMap.h"
124 #include "AliTPCCalibVdrift.h"
125 #include "AliTPCCalibRaw.h"
126 #include "AliTPCParam.h"
127 #include "AliTPCCorrection.h"
128 #include "AliTPCPreprocessorOnline.h"
131 ClassImp(AliTPCcalibDB)
133 AliTPCcalibDB* AliTPCcalibDB::fgInstance = 0;
134 Bool_t AliTPCcalibDB::fgTerminated = kFALSE;
135 TObjArray AliTPCcalibDB::fgExBArray; // array of ExB corrections
138 //_ singleton implementation __________________________________________________
139 AliTPCcalibDB* AliTPCcalibDB::Instance()
142 // Singleton implementation
143 // Returns an instance of this class, it is created if neccessary
146 if (fgTerminated != kFALSE)
150 fgInstance = new AliTPCcalibDB();
155 void AliTPCcalibDB::Terminate()
158 // Singleton implementation
159 // Deletes the instance of this class and sets the terminated flag, instances cannot be requested anymore
160 // This function can be called several times.
163 fgTerminated = kTRUE;
172 //_____________________________________________________________________________
173 AliTPCcalibDB::AliTPCcalibDB():
182 fComposedCorrection(0),
183 fComposedCorrectionArray(0),
196 fTimeGainSplinesArray(100000),
197 fGRPArray(100000), //! array of GRPs - per run - JUST for calibration studies
198 fGRPMaps(100000), //! array of GRPs - per run - JUST for calibration studies
199 fGoofieArray(100000), //! array of GOOFIE values -per run - Just for calibration studies
200 fVoltageArray(100000),
201 fTemperatureArray(100000), //! array of temperature sensors - per run - Just for calibration studies
202 fVdriftArray(100000), //! array of v drift interfaces
203 fDriftCorrectionArray(100000), //! array of drift correction
204 fRunList(100000), //! run list - indicates try to get the run param
213 Update(); // temporary
216 AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& ):
225 fComposedCorrection(0),
226 fComposedCorrectionArray(0),
239 fTimeGainSplinesArray(100000),
240 fGRPArray(0), //! array of GRPs - per run - JUST for calibration studies
241 fGRPMaps(0), //! array of GRPs - per run - JUST for calibration studies
242 fGoofieArray(0), //! array of GOOFIE values -per run - Just for calibration studies
244 fTemperatureArray(0), //! array of temperature sensors - per run - Just for calibration studies
245 fVdriftArray(0), //! array of v drift interfaces
246 fDriftCorrectionArray(0), //! array of v drift interfaces
247 fRunList(0), //! run list - indicates try to get the run param
252 // Copy constructor invalid -- singleton implementation
254 Error("copy constructor","invalid -- singleton implementation");
257 AliTPCcalibDB& AliTPCcalibDB::operator= (const AliTPCcalibDB& )
260 // Singleton implementation - no assignment operator
262 Error("operator =", "assignment operator not implemented");
268 //_____________________________________________________________________________
269 AliTPCcalibDB::~AliTPCcalibDB()
276 AliTPCCalPad* AliTPCcalibDB::GetDistortionMap(Int_t i) const {
278 // get distortion map - due E field distortions
280 return (fDistortionMap) ? (AliTPCCalPad*)fDistortionMap->At(i):0;
283 //_____________________________________________________________________________
284 AliCDBEntry* AliTPCcalibDB::GetCDBEntry(const char* cdbPath)
287 // Retrieves an entry with path <cdbPath> from the CDB.
291 AliCDBEntry* entry = AliCDBManager::Instance()->Get(cdbPath, fRun);
294 sprintf(chinfo,"AliTPCcalibDB: Failed to get entry:\t%s ", cdbPath);
302 //_____________________________________________________________________________
303 void AliTPCcalibDB::SetRun(Long64_t run)
306 // Sets current run number. Calibration data is read from the corresponding file.
316 void AliTPCcalibDB::Update(){
318 // cache the OCDB entries for simulation, reconstruction, calibration
321 AliCDBEntry * entry=0;
322 Bool_t cdbCache = AliCDBManager::Instance()->GetCacheFlag(); // save cache status
323 AliCDBManager::Instance()->SetCacheFlag(kTRUE); // activate CDB cache
324 fDButil = new AliTPCcalibDButil;
326 entry = GetCDBEntry("TPC/Calib/PadGainFactor");
328 //if (fPadGainFactor) delete fPadGainFactor;
329 entry->SetOwner(kTRUE);
330 fPadGainFactor = (AliTPCCalPad*)entry->GetObject();
332 AliFatal("TPC - Missing calibration entry TPC/Calib/PadGainFactor")
335 entry = GetCDBEntry("TPC/Calib/TimeGain");
337 //if (fTimeGainSplines) delete fTimeGainSplines;
338 entry->SetOwner(kTRUE);
339 fTimeGainSplines = (TObjArray*)entry->GetObject();
341 AliFatal("TPC - Missing calibration entry TPC/Calib/Timegain")
344 entry = GetCDBEntry("TPC/Calib/GainFactorDedx");
346 entry->SetOwner(kTRUE);
347 fDedxGainFactor = (AliTPCCalPad*)entry->GetObject();
349 AliFatal("TPC - Missing calibration entry TPC/Calib/gainFactordEdx")
352 entry = GetCDBEntry("TPC/Calib/PadTime0");
354 //if (fPadTime0) delete fPadTime0;
355 entry->SetOwner(kTRUE);
356 fPadTime0 = (AliTPCCalPad*)entry->GetObject();
358 AliFatal("TPC - Missing calibration entry")
361 entry = GetCDBEntry("TPC/Calib/Distortion");
363 //if (fPadTime0) delete fPadTime0;
364 entry->SetOwner(kTRUE);
365 fDistortionMap =dynamic_cast<TObjArray*>(entry->GetObject());
367 //AliFatal("TPC - Missing calibration entry")
373 entry = GetCDBEntry("TPC/Calib/PadNoise");
375 //if (fPadNoise) delete fPadNoise;
376 entry->SetOwner(kTRUE);
377 fPadNoise = (AliTPCCalPad*)entry->GetObject();
379 AliFatal("TPC - Missing calibration entry")
382 entry = GetCDBEntry("TPC/Calib/Pedestals");
384 //if (fPedestals) delete fPedestals;
385 entry->SetOwner(kTRUE);
386 fPedestals = (AliTPCCalPad*)entry->GetObject();
389 entry = GetCDBEntry("TPC/Calib/Temperature");
391 //if (fTemperature) delete fTemperature;
392 entry->SetOwner(kTRUE);
393 fTemperature = (AliTPCSensorTempArray*)entry->GetObject();
396 entry = GetCDBEntry("TPC/Calib/Parameters");
398 //if (fPadNoise) delete fPadNoise;
399 entry->SetOwner(kTRUE);
400 fParam = (AliTPCParam*)(entry->GetObject()->Clone());
402 AliFatal("TPC - Missing calibration entry TPC/Calib/Parameters")
405 entry = GetCDBEntry("TPC/Calib/ClusterParam");
407 entry->SetOwner(kTRUE);
408 fClusterParam = (AliTPCClusterParam*)(entry->GetObject()->Clone());
410 AliFatal("TPC - Missing calibration entry")
413 //ALTRO configuration data
414 entry = GetCDBEntry("TPC/Calib/AltroConfig");
416 entry->SetOwner(kTRUE);
417 fALTROConfigData=(TObjArray*)(entry->GetObject());
419 AliFatal("TPC - Missing calibration entry")
422 //Calibration Pulser data
423 entry = GetCDBEntry("TPC/Calib/Pulser");
425 entry->SetOwner(kTRUE);
426 fPulserData=(TObjArray*)(entry->GetObject());
430 entry = GetCDBEntry("TPC/Calib/CE");
432 entry->SetOwner(kTRUE);
433 fCEData=(TObjArray*)(entry->GetObject());
435 //RAW calibration data
436 // entry = GetCDBEntry("TPC/Calib/Raw");
438 entry = GetCDBEntry("TPC/Calib/Mapping");
440 //if (fPadNoise) delete fPadNoise;
441 entry->SetOwner(kTRUE);
442 TObjArray * array = dynamic_cast<TObjArray*>(entry->GetObject());
443 if (array && array->GetEntriesFast()==6){
444 fMapping = new AliTPCAltroMapping*[6];
445 for (Int_t i=0; i<6; i++){
446 fMapping[i] = dynamic_cast<AliTPCAltroMapping*>(array->At(i));
451 //CTP calibration data
452 entry = GetCDBEntry("GRP/CTP/CTPtiming");
454 //entry->SetOwner(kTRUE);
455 fCTPTimeParams=dynamic_cast<AliCTPTimeParams*>(entry->GetObject());
457 AliError("TPC - Missing calibration entry")
459 //TPC space point correction data
460 entry = GetCDBEntry("TPC/Calib/Correction");
462 //entry->SetOwner(kTRUE);
463 fComposedCorrection=dynamic_cast<AliTPCCorrection*>(entry->GetObject());
464 if (fComposedCorrection) fComposedCorrection->Init();
465 fComposedCorrectionArray=dynamic_cast<TObjArray*>(entry->GetObject());
466 if (fComposedCorrectionArray){
467 for (Int_t i=0; i<fComposedCorrectionArray->GetEntries(); i++){
468 AliTPCCorrection* composedCorrection= dynamic_cast<AliTPCCorrection*>(fComposedCorrectionArray->At(i));
469 if (composedCorrection) composedCorrection->Init();
473 AliError("TPC - Missing calibration entry- TPC/Calib/Correction")
478 fTransform=new AliTPCTransform();
479 fTransform->SetCurrentRun(AliCDBManager::Instance()->GetRun());
483 AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
486 void AliTPCcalibDB::UpdateNonRec(){
488 // Update/Load the parameters which are important for QA studies
489 // and not used yet for the reconstruction
491 //RAW calibration data
492 AliCDBEntry * entry=0;
493 entry = GetCDBEntry("TPC/Calib/Raw");
495 entry->SetOwner(kTRUE);
496 TObjArray *arr=(TObjArray*)(entry->GetObject());
497 if (arr) fCalibRaw=(AliTPCCalibRaw*)arr->At(0);
499 //QA calibration data
500 entry = GetCDBEntry("TPC/Calib/QA");
502 entry->SetOwner(kTRUE);
503 fDataQA=dynamic_cast<AliTPCdataQA*>(entry->GetObject());
507 entry = AliCDBManager::Instance()->Get("TPC/Calib/HighVoltage",fRun);
509 fVoltageArray.AddAt(entry->GetObject(),fRun);
517 void AliTPCcalibDB::CreateObjectList(const Char_t *filename, TObjArray *calibObjects)
520 // Create calibration objects and read contents from OCDB
522 if ( calibObjects == 0x0 ) return;
525 if ( !in.is_open() ){
526 fprintf(stderr,"Error: cannot open list file '%s'", filename);
530 AliTPCCalPad *calPad=0x0;
536 TObjArray *arrFileLine = sFile.Tokenize("\n");
538 TIter nextLine(arrFileLine);
540 TObjString *sObjLine=0x0;
541 while ( (sObjLine = (TObjString*)nextLine()) ){
542 TString sLine(sObjLine->GetString());
544 TObjArray *arrNextCol = sLine.Tokenize("\t");
546 TObjString *sObjType = (TObjString*)(arrNextCol->At(0));
547 TObjString *sObjFileName = (TObjString*)(arrNextCol->At(1));
549 if ( !sObjType || ! sObjFileName ) continue;
550 TString sType(sObjType->GetString());
551 TString sFileName(sObjFileName->GetString());
552 printf("%s\t%s\n",sType.Data(),sFileName.Data());
554 TFile *fIn = TFile::Open(sFileName);
556 fprintf(stderr,"File not found: '%s'", sFileName.Data());
560 if ( sType == "CE" ){
561 AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
563 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());
564 calPad->SetNameTitle("CETmean","CETmean");
565 calibObjects->Add(calPad);
567 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());
568 calPad->SetNameTitle("CEQmean","CEQmean");
569 calibObjects->Add(calPad);
571 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
572 calPad->SetNameTitle("CETrms","CETrms");
573 calibObjects->Add(calPad);
575 } else if ( sType == "Pulser") {
576 AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
578 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());
579 calPad->SetNameTitle("PulserTmean","PulserTmean");
580 calibObjects->Add(calPad);
582 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());
583 calPad->SetNameTitle("PulserQmean","PulserQmean");
584 calibObjects->Add(calPad);
586 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
587 calPad->SetNameTitle("PulserTrms","PulserTrms");
588 calibObjects->Add(calPad);
590 } else if ( sType == "Pedestals") {
591 AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
593 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());
594 calPad->SetNameTitle("Pedestals","Pedestals");
595 calibObjects->Add(calPad);
597 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());
598 calPad->SetNameTitle("Noise","Noise");
599 calibObjects->Add(calPad);
602 fprintf(stderr,"Undefined Type: '%s'",sType.Data());
611 void AliTPCcalibDB::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
613 // Write a tree with all available information
614 // if mapFileName is specified, the Map information are also written to the tree
615 // pads specified in outlierPad are not used for calculating statistics
616 // - the same function as AliTPCCalPad::MakeTree -
618 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
620 TObjArray* mapIROCs = 0;
621 TObjArray* mapOROCs = 0;
622 TVectorF *mapIROCArray = 0;
623 TVectorF *mapOROCArray = 0;
624 Int_t mapEntries = 0;
625 TString* mapNames = 0;
628 TFile mapFile(mapFileName, "read");
630 TList* listOfROCs = mapFile.GetListOfKeys();
631 mapEntries = listOfROCs->GetEntries()/2;
632 mapIROCs = new TObjArray(mapEntries*2);
633 mapOROCs = new TObjArray(mapEntries*2);
634 mapIROCArray = new TVectorF[mapEntries];
635 mapOROCArray = new TVectorF[mapEntries];
637 mapNames = new TString[mapEntries];
638 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
639 TString nameROC(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
640 nameROC.Remove(nameROC.Length()-4, 4);
641 mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "IROC").Data()), ivalue);
642 mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "OROC").Data()), ivalue);
643 mapNames[ivalue].Append(nameROC);
646 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
647 mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
648 mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
650 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
651 (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
652 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
653 (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
656 } // if (mapFileName)
658 TTreeSRedirector cstream(fileName);
659 Int_t arrayEntries = array->GetEntries();
661 TString* names = new TString[arrayEntries];
662 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
663 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
665 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
667 // get statistic for given sector
669 TVectorF median(arrayEntries);
670 TVectorF mean(arrayEntries);
671 TVectorF rms(arrayEntries);
672 TVectorF ltm(arrayEntries);
673 TVectorF ltmrms(arrayEntries);
674 TVectorF medianWithOut(arrayEntries);
675 TVectorF meanWithOut(arrayEntries);
676 TVectorF rmsWithOut(arrayEntries);
677 TVectorF ltmWithOut(arrayEntries);
678 TVectorF ltmrmsWithOut(arrayEntries);
680 TVectorF *vectorArray = new TVectorF[arrayEntries];
681 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
682 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
684 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
685 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
686 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
687 AliTPCCalROC* outlierROC = 0;
688 if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
690 median[ivalue] = calROC->GetMedian();
691 mean[ivalue] = calROC->GetMean();
692 rms[ivalue] = calROC->GetRMS();
693 Double_t ltmrmsValue = 0;
694 ltm[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction);
695 ltmrms[ivalue] = ltmrmsValue;
697 medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
698 meanWithOut[ivalue] = calROC->GetMean(outlierROC);
699 rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
701 ltmWithOut[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction, outlierROC);
702 ltmrmsWithOut[ivalue] = ltmrmsValue;
711 medianWithOut[ivalue] = 0.;
712 meanWithOut[ivalue] = 0.;
713 rmsWithOut[ivalue] = 0.;
714 ltmWithOut[ivalue] = 0.;
715 ltmrmsWithOut[ivalue] = 0.;
720 // fill vectors of variable per pad
722 TVectorF *posArray = new TVectorF[8];
723 for (Int_t ivalue = 0; ivalue < 8; ivalue++)
724 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
726 Float_t posG[3] = {0};
727 Float_t posL[3] = {0};
729 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
730 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
731 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
732 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
733 posArray[0][ichannel] = irow;
734 posArray[1][ichannel] = ipad;
735 posArray[2][ichannel] = posL[0];
736 posArray[3][ichannel] = posL[1];
737 posArray[4][ichannel] = posG[0];
738 posArray[5][ichannel] = posG[1];
739 posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
740 posArray[7][ichannel] = ichannel;
742 // loop over array containing AliTPCCalPads
743 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
744 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
745 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
747 (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
749 (vectorArray[ivalue])[ichannel] = 0;
755 cstream << "calPads" <<
756 "sector=" << isector;
758 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
759 cstream << "calPads" <<
760 (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
761 (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
762 (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
763 (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
764 (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
766 cstream << "calPads" <<
767 (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
768 (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
769 (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
770 (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
771 (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
775 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
776 cstream << "calPads" <<
777 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
781 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
783 cstream << "calPads" <<
784 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
786 cstream << "calPads" <<
787 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
791 cstream << "calPads" <<
792 "row.=" << &posArray[0] <<
793 "pad.=" << &posArray[1] <<
794 "lx.=" << &posArray[2] <<
795 "ly.=" << &posArray[3] <<
796 "gx.=" << &posArray[4] <<
797 "gy.=" << &posArray[5] <<
798 "rpad.=" << &posArray[6] <<
799 "channel.=" << &posArray[7];
801 cstream << "calPads" <<
805 delete[] vectorArray;
813 delete[] mapIROCArray;
814 delete[] mapOROCArray;
819 Int_t AliTPCcalibDB::GetRCUTriggerConfig() const
822 // return the RCU trigger configuration register
824 TMap *map=GetRCUconfig();
826 TVectorF *v=(TVectorF*)map->GetValue("TRGCONF_TRG_MODE");
828 for (Int_t i=0; i<v->GetNrows(); ++i){
829 Float_t newmode=v->GetMatrixArray()[i];
831 if (mode>-1&&newmode!=mode) AliWarning("Found different RCU trigger configurations!!!");
838 Bool_t AliTPCcalibDB::IsTrgL0()
841 // return if the FEE readout was triggered on L0
843 Int_t mode=GetRCUTriggerConfig();
844 if (mode<0) return kFALSE;
848 Bool_t AliTPCcalibDB::IsTrgL1()
851 // return if the FEE readout was triggered on L1
853 Int_t mode=GetRCUTriggerConfig();
854 if (mode<0) return kFALSE;
858 void AliTPCcalibDB::RegisterExB(Int_t index, Float_t bz, Bool_t bdelete){
860 // Register static ExB correction map
861 // index - registration index - used for visualization
862 // bz - bz field in kGaus
864 // Float_t factor = bz/(-5.); // default b filed in Cheb with minus sign
865 Float_t factor = bz/(5.); // default b filed in Cheb with minus sign
866 // was chenged in the Revision ???? (Ruben can you add here number)
868 AliMagF* bmap = new AliMagF("MapsExB","MapsExB", factor,TMath::Sign(1.f,factor),AliMagF::k5kG);
870 AliTPCExBFirst *exb = new AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
871 AliTPCExB::SetInstance(exb);
876 AliTPCExB::RegisterField(index,bmap);
878 if (index>=fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
879 fgExBArray.AddAt(exb,index);
883 AliTPCExB* AliTPCcalibDB::GetExB(Float_t bz, Bool_t deleteB) {
885 // bz filed in KGaus not in tesla
886 // Get ExB correction map
887 // if doesn't exist - create it
889 Int_t index = TMath::Nint(5+bz);
890 if (index>fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
891 if (!fgExBArray.At(index)) AliTPCcalibDB::RegisterExB(index,bz,deleteB);
892 return (AliTPCExB*)fgExBArray.At(index);
896 void AliTPCcalibDB::SetExBField(Float_t bz){
898 // Set magnetic filed for ExB correction
900 fExB = GetExB(bz,kFALSE);
903 void AliTPCcalibDB::SetExBField(const AliMagF* bmap){
905 // Set magnetic field for ExB correction
907 AliTPCExBFirst *exb = new AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
908 AliTPCExB::SetInstance(exb);
916 void AliTPCcalibDB::UpdateRunInformations( Int_t run, Bool_t force){
918 // - > Don't use it for reconstruction - Only for Calibration studies
922 AliCDBEntry * entry = 0;
923 if (run>= fRunList.fN){
924 fRunList.Set(run*2+1);
925 fGRPArray.Expand(run*2+1);
926 fGRPMaps.Expand(run*2+1);
927 fGoofieArray.Expand(run*2+1);
928 fVoltageArray.Expand(run*2+1);
929 fTemperatureArray.Expand(run*2+1);
930 fVdriftArray.Expand(run*2+1);
931 fDriftCorrectionArray.Expand(run*2+1);
932 fTimeGainSplinesArray.Expand(run*2+1);
935 fALTROConfigData->Expand(run*2+1); // ALTRO configuration data
936 fPulserData->Expand(run*2+1); // Calibration Pulser data
937 fCEData->Expand(run*2+1); // CE data
938 if (!fTimeGainSplines) fTimeGainSplines = new TObjArray(run*2+1);
939 fTimeGainSplines->Expand(run*2+1); // Array of AliSplineFits: at 0 MIP position in
941 if (fRunList[run]>0 &&force==kFALSE) return;
943 fRunList[run]=1; // sign as used
946 entry = AliCDBManager::Instance()->Get("GRP/GRP/Data",run);
948 AliGRPObject * grpRun = dynamic_cast<AliGRPObject*>(entry->GetObject());
950 TMap* map = dynamic_cast<TMap*>(entry->GetObject());
952 //grpRun = new AliGRPObject;
953 //grpRun->ReadValuesFromMap(map);
954 grpRun = MakeGRPObjectFromMap(map);
956 fGRPMaps.AddAt(map,run);
959 fGRPArray.AddAt(grpRun,run);
961 entry = AliCDBManager::Instance()->Get("TPC/Calib/Goofie",run);
963 fGoofieArray.AddAt(entry->GetObject(),run);
968 entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeGain",run);
970 fTimeGainSplinesArray.AddAt(entry->GetObject(),run);
972 AliFatal("TPC - Missing calibration entry TimeGain")
975 entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeDrift",run);
977 fDriftCorrectionArray.AddAt(entry->GetObject(),run);
979 AliFatal("TPC - Missing calibration entry TimeDrift")
982 entry = AliCDBManager::Instance()->Get("TPC/Calib/Temperature",run);
984 fTemperatureArray.AddAt(entry->GetObject(),run);
986 //apply fDButil filters
988 fDButil->UpdateFromCalibDB();
989 if (fTemperature) fDButil->FilterTemperature(fTemperature);
991 AliDCSSensor * press = GetPressureSensor(run,0);
992 AliTPCSensorTempArray * temp = GetTemperatureSensor(run);
995 accept = fDButil->FilterTemperature(temp)>0.1;
998 const Double_t kMinP=950.;
999 const Double_t kMaxP=1050.;
1000 const Double_t kMaxdP=10.;
1001 const Double_t kSigmaCut=4.;
1002 fDButil->FilterSensor(press,kMinP,kMaxP,kMaxdP,kSigmaCut);
1003 if (press->GetFit()==0) accept=kFALSE;
1005 if (press && temp &&accept){
1006 AliTPCCalibVdrift * vdrift = new AliTPCCalibVdrift(temp, press,0);
1007 fVdriftArray.AddAt(vdrift,run);
1009 fDButil->FilterCE(120., 3., 4.,0);
1010 fDButil->FilterTracks(run, 10.,0);
1014 Float_t AliTPCcalibDB::GetGain(Int_t sector, Int_t row, Int_t pad){
1016 // Get Gain factor for given pad
1018 AliTPCCalPad *calPad = Instance()->fDedxGainFactor;;
1019 if (!calPad) return 0;
1020 return calPad->GetCalROC(sector)->GetValue(row,pad);
1023 AliSplineFit* AliTPCcalibDB::GetVdriftSplineFit(const char* name, Int_t run){
1025 // GetDrift velocity spline fit
1027 TObjArray *arr=GetTimeVdriftSplineRun(run);
1029 return dynamic_cast<AliSplineFit*>(arr->FindObject(name));
1032 AliSplineFit* AliTPCcalibDB::CreateVdriftSplineFit(const char* graphName, Int_t run){
1034 // create spline fit from the drift time graph in TimeDrift
1036 TObjArray *arr=GetTimeVdriftSplineRun(run);
1038 TGraph *graph=dynamic_cast<TGraph*>(arr->FindObject(graphName));
1039 if (!graph) return 0;
1040 AliSplineFit *fit = new AliSplineFit();
1041 fit->SetGraph(graph);
1042 fit->SetMinPoints(graph->GetN()+1);
1043 fit->InitKnots(graph,2,0,0.001);
1048 AliGRPObject *AliTPCcalibDB::GetGRP(Int_t run){
1050 // Get GRP object for given run
1052 if (run>= ((Instance()->fGRPArray)).GetEntriesFast()){
1053 Instance()->UpdateRunInformations(run);
1055 AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>((Instance()->fGRPArray).At(run));
1057 Instance()->UpdateRunInformations(run);
1058 grpRun = dynamic_cast<AliGRPObject *>(Instance()->fGRPArray.At(run));
1059 if (!grpRun) return 0;
1064 TMap * AliTPCcalibDB::GetGRPMap(Int_t run){
1066 // Get GRP map for given run
1068 TMap * grpRun = dynamic_cast<TMap *>((Instance()->fGRPMaps).At(run));
1070 Instance()->UpdateRunInformations(run);
1071 grpRun = dynamic_cast<TMap *>(Instance()->fGRPMaps.At(run));
1072 if (!grpRun) return 0;
1078 AliDCSSensor * AliTPCcalibDB::GetPressureSensor(Int_t run, Int_t type){
1080 // Get Pressure sensor
1082 // type = 0 - Cavern pressure
1083 // 1 - Suface pressure
1084 // First try to get if trom map - if existing (Old format of data storing)
1088 TMap *map = GetGRPMap(run);
1090 AliDCSSensor * sensor = 0;
1092 if (type==0) osensor = ((*map)("fCavernPressure"));
1093 if (type==1) osensor = ((*map)("fP2Pressure"));
1094 sensor =dynamic_cast<AliDCSSensor *>(osensor);
1095 if (sensor) return sensor;
1098 // If not map try to get it from the GRPObject
1100 AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run));
1102 UpdateRunInformations(run);
1103 grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run));
1104 if (!grpRun) return 0;
1106 AliDCSSensor * sensor = grpRun->GetCavernAtmosPressure();
1107 if (type==1) sensor = grpRun->GetSurfaceAtmosPressure();
1111 AliTPCSensorTempArray * AliTPCcalibDB::GetTemperatureSensor(Int_t run){
1113 // Get temperature sensor array
1115 AliTPCSensorTempArray * tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
1117 UpdateRunInformations(run);
1118 tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
1124 TObjArray * AliTPCcalibDB::GetTimeGainSplinesRun(Int_t run){
1126 // Get temperature sensor array
1128 TObjArray * gainSplines = (TObjArray *)fTimeGainSplinesArray.At(run);
1130 UpdateRunInformations(run);
1131 gainSplines = (TObjArray *)fTimeGainSplinesArray.At(run);
1136 TObjArray * AliTPCcalibDB::GetTimeVdriftSplineRun(Int_t run){
1138 // Get drift spline array
1140 TObjArray * driftSplines = (TObjArray *)fDriftCorrectionArray.At(run);
1141 if (!driftSplines) {
1142 UpdateRunInformations(run);
1143 driftSplines = (TObjArray *)fDriftCorrectionArray.At(run);
1145 return driftSplines;
1148 AliDCSSensorArray * AliTPCcalibDB::GetVoltageSensors(Int_t run){
1150 // Get temperature sensor array
1152 AliDCSSensorArray * voltageArray = (AliDCSSensorArray *)fVoltageArray.At(run);
1153 if (!voltageArray) {
1154 UpdateRunInformations(run);
1155 voltageArray = (AliDCSSensorArray *)fVoltageArray.At(run);
1157 return voltageArray;
1160 AliDCSSensorArray * AliTPCcalibDB::GetGoofieSensors(Int_t run){
1162 // Get temperature sensor array
1164 AliDCSSensorArray * goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
1166 UpdateRunInformations(run);
1167 goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
1174 AliTPCCalibVdrift * AliTPCcalibDB::GetVdrift(Int_t run){
1176 // Get the interface to the the vdrift
1178 AliTPCCalibVdrift * vdrift = (AliTPCCalibVdrift*)fVdriftArray.At(run);
1180 UpdateRunInformations(run);
1181 vdrift= (AliTPCCalibVdrift*)fVdriftArray.At(run);
1186 Float_t AliTPCcalibDB::GetCEdriftTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1189 // GetCE drift time information for 'sector'
1190 // sector 72 is the mean drift time of the A-Side
1191 // sector 73 is the mean drift time of the C-Side
1192 // it timestamp==-1 return mean value
1194 AliTPCcalibDB::Instance()->SetRun(run);
1195 TGraph *gr=AliTPCcalibDB::Instance()->GetCErocTgraph(sector);
1196 if (!gr||sector<0||sector>73) {
1197 if (entries) *entries=0;
1201 if (timeStamp==-1.){
1204 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1206 gr->GetPoint(ipoint,x,y);
1207 if (x<timeStamp) continue;
1215 Float_t AliTPCcalibDB::GetCEchargeTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1218 // GetCE mean charge for 'sector'
1219 // it timestamp==-1 return mean value
1221 AliTPCcalibDB::Instance()->SetRun(run);
1222 TGraph *gr=AliTPCcalibDB::Instance()->GetCErocQgraph(sector);
1223 if (!gr||sector<0||sector>71) {
1224 if (entries) *entries=0;
1228 if (timeStamp==-1.){
1231 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1233 gr->GetPoint(ipoint,x,y);
1234 if (x<timeStamp) continue;
1242 Float_t AliTPCcalibDB::GetDCSSensorValue(AliDCSSensorArray *arr, Int_t timeStamp, const char * sensorName, Int_t sigDigits)
1245 // Get Value for a DCS sensor 'sensorName', run 'run' at time 'timeStamp'
1248 const TString sensorNameString(sensorName);
1249 AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1250 if (!sensor) return val;
1251 //use the dcs graph if possible
1252 TGraph *gr=sensor->GetGraph();
1254 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1256 gr->GetPoint(ipoint,x,y);
1257 Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1258 if (time<timeStamp) continue;
1262 //if val is still 0, test if if the requested time if within 5min of the first/last
1263 //data point. If this is the case return the firs/last entry
1264 //the timestamps might not be syncronised for all calibration types, sometimes a 'pre'
1265 //and 'pos' period is requested. Especially to the HV this is not the case!
1269 gr->GetPoint(0,x,y);
1270 Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1271 if ((time-timeStamp)<5*60) val=y;
1276 gr->GetPoint(gr->GetN()-1,x,y);
1277 Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1278 if ((timeStamp-time)<5*60) val=y;
1281 val=sensor->GetValue(timeStamp);
1284 val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1289 Float_t AliTPCcalibDB::GetDCSSensorMeanValue(AliDCSSensorArray *arr, const char * sensorName, Int_t sigDigits)
1292 // Get mean Value for a DCS sensor 'sensorName' during run 'run'
1295 const TString sensorNameString(sensorName);
1296 AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1297 if (!sensor) return val;
1299 //use dcs graph if it exists
1300 TGraph *gr=sensor->GetGraph();
1304 //if we don't have the dcs graph, try to get some meaningful information
1305 if (!sensor->GetFit()) return val;
1306 Int_t nKnots=sensor->GetFit()->GetKnots();
1307 Double_t tMid=(sensor->GetEndTime()-sensor->GetStartTime())/2.;
1308 for (Int_t iKnot=0;iKnot<nKnots;++iKnot){
1309 if (sensor->GetFit()->GetX()[iKnot]>tMid/3600.) break;
1310 val=(Float_t)sensor->GetFit()->GetY0()[iKnot];
1315 val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1321 Float_t AliTPCcalibDB::GetChamberHighVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits) {
1323 // return the chamber HV for given run and time: 0-35 IROC, 36-72 OROC
1324 // if timeStamp==-1 return mean value
1327 TString sensorName="";
1328 TTimeStamp stamp(timeStamp);
1329 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1330 if (!voltageArray || (sector<0) || (sector>71)) return val;
1331 Char_t sideName='A';
1332 if ((sector/18)%2==1) sideName='C';
1335 sensorName=Form("TPC_ANODE_I_%c%02d_VMEAS",sideName,sector%18);
1338 sensorName=Form("TPC_ANODE_O_%c%02d_0_VMEAS",sideName,sector%18);
1341 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1343 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1347 Float_t AliTPCcalibDB::GetSkirtVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1350 // Get the skirt voltage for 'run' at 'timeStamp' and 'sector': 0-35 IROC, 36-72 OROC
1351 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1352 // if timeStamp==-1 return the mean value for the run
1355 TString sensorName="";
1356 TTimeStamp stamp(timeStamp);
1357 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1358 if (!voltageArray || (sector<0) || (sector>71)) return val;
1359 Char_t sideName='A';
1360 if ((sector/18)%2==1) sideName='C';
1361 sensorName=Form("TPC_SKIRT_%c_VMEAS",sideName);
1363 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1365 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1370 Float_t AliTPCcalibDB::GetCoverVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1373 // Get the cover voltage for run 'run' at time 'timeStamp'
1374 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1375 // if timeStamp==-1 return the mean value for the run
1378 TString sensorName="";
1379 TTimeStamp stamp(timeStamp);
1380 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1381 if (!voltageArray || (sector<0) || (sector>71)) return val;
1382 Char_t sideName='A';
1383 if ((sector/18)%2==1) sideName='C';
1386 sensorName=Form("TPC_COVER_I_%c_VMEAS",sideName);
1389 sensorName=Form("TPC_COVER_O_%c_VMEAS",sideName);
1392 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1394 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1399 Float_t AliTPCcalibDB::GetGGoffsetVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1402 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1403 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1404 // if timeStamp==-1 return the mean value for the run
1407 TString sensorName="";
1408 TTimeStamp stamp(timeStamp);
1409 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1410 if (!voltageArray || (sector<0) || (sector>71)) return val;
1411 Char_t sideName='A';
1412 if ((sector/18)%2==1) sideName='C';
1415 sensorName=Form("TPC_GATE_I_%c_OFF_VMEAS",sideName);
1418 sensorName=Form("TPC_GATE_O_%c_OFF_VMEAS",sideName);
1421 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1423 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1428 Float_t AliTPCcalibDB::GetGGnegVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1431 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1432 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1433 // if timeStamp==-1 return the mean value for the run
1436 TString sensorName="";
1437 TTimeStamp stamp(timeStamp);
1438 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1439 if (!voltageArray || (sector<0) || (sector>71)) return val;
1440 Char_t sideName='A';
1441 if ((sector/18)%2==1) sideName='C';
1444 sensorName=Form("TPC_GATE_I_%c_NEG_VMEAS",sideName);
1447 sensorName=Form("TPC_GATE_O_%c_NEG_VMEAS",sideName);
1450 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1452 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1457 Float_t AliTPCcalibDB::GetGGposVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1460 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1461 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1462 // if timeStamp==-1 return the mean value for the run
1465 TString sensorName="";
1466 TTimeStamp stamp(timeStamp);
1467 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1468 if (!voltageArray || (sector<0) || (sector>71)) return val;
1469 Char_t sideName='A';
1470 if ((sector/18)%2==1) sideName='C';
1473 sensorName=Form("TPC_GATE_I_%c_POS_VMEAS",sideName);
1476 sensorName=Form("TPC_GATE_O_%c_POS_VMEAS",sideName);
1479 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1481 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1486 Float_t AliTPCcalibDB::GetPressure(Int_t timeStamp, Int_t run, Int_t type){
1488 // GetPressure for given time stamp and runt
1490 TTimeStamp stamp(timeStamp);
1491 AliDCSSensor * sensor = Instance()->GetPressureSensor(run,type);
1492 if (!sensor) return 0;
1493 return sensor->GetValue(stamp);
1496 Float_t AliTPCcalibDB::GetL3Current(Int_t run, Int_t statType){
1498 // return L3 current
1499 // stat type is: AliGRPObject::Stats: kMean = 0, kTruncMean = 1, kMedian = 2, kSDMean = 3, kSDMedian = 4
1502 AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1503 if (grp) current=grp->GetL3Current((AliGRPObject::Stats)statType);
1507 Float_t AliTPCcalibDB::GetBz(Int_t run){
1509 // calculate BZ in T from L3 current
1512 Float_t current=AliTPCcalibDB::GetL3Current(run);
1513 if (current>-1) bz=5*current/30000.*.1;
1517 Char_t AliTPCcalibDB::GetL3Polarity(Int_t run) {
1519 // get l3 polarity from GRP
1522 AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1523 if (grp) pol=grp->GetL3Polarity();
1527 TString AliTPCcalibDB::GetRunType(Int_t run){
1529 // return run type from grp
1532 // TString type("UNKNOWN");
1533 AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1534 if (grp) return grp->GetRunType();
1538 Float_t AliTPCcalibDB::GetValueGoofie(Int_t timeStamp, Int_t run, Int_t type){
1540 // GetPressure for given time stamp and runt
1542 TTimeStamp stamp(timeStamp);
1543 AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(run);
1544 if (!goofieArray) return 0;
1545 AliDCSSensor *sensor = goofieArray->GetSensor(type);
1546 return sensor->GetValue(stamp);
1554 Bool_t AliTPCcalibDB::GetTemperatureFit(Int_t timeStamp, Int_t run, Int_t side,TVectorD& fit){
1556 // GetTmeparature fit at parameter for given time stamp
1558 TTimeStamp tstamp(timeStamp);
1559 AliTPCSensorTempArray* tempArray = Instance()->GetTemperatureSensor(run);
1560 if (! tempArray) return kFALSE;
1561 AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
1562 TLinearFitter * fitter = tempMap->GetLinearFitter(3,side,tstamp);
1565 fitter->GetParameters(fit);
1569 if (!fitter) return kFALSE;
1573 Float_t AliTPCcalibDB::GetTemperature(Int_t timeStamp, Int_t run, Int_t side){
1575 // Get mean temperature
1579 GetTemperatureFit(timeStamp,run,0,vec);
1583 GetTemperatureFit(timeStamp,run,0,vec);
1590 Double_t AliTPCcalibDB::GetPTRelative(UInt_t timeSec, Int_t run, Int_t side){
1593 // time - absolute time
1595 // side - 0 - A side 1-C side
1596 AliTPCCalibVdrift * vdrift = Instance()->GetVdrift(run);
1597 if (!vdrift) return 0;
1598 return vdrift->GetPTRelative(timeSec,side);
1601 AliGRPObject * AliTPCcalibDB::MakeGRPObjectFromMap(TMap *map){
1603 // Function to covert old GRP run information from TMap to GRPObject
1605 // TMap * map = AliTPCcalibDB::GetGRPMap(52406);
1607 AliDCSSensor * sensor = 0;
1609 osensor = ((*map)("fP2Pressure"));
1610 sensor =dynamic_cast<AliDCSSensor *>(osensor);
1612 if (!sensor) return 0;
1614 AliDCSSensor * sensor2 = new AliDCSSensor(*sensor);
1615 osensor = ((*map)("fCavernPressure"));
1616 TGraph * gr = new TGraph(2);
1617 gr->GetX()[0]= -100000.;
1618 gr->GetX()[1]= 1000000.;
1619 gr->GetY()[0]= atof(osensor->GetName());
1620 gr->GetY()[1]= atof(osensor->GetName());
1621 sensor2->SetGraph(gr);
1625 AliGRPObject *grpRun = new AliGRPObject;
1626 grpRun->ReadValuesFromMap(map);
1627 grpRun->SetCavernAtmosPressure(sensor2);
1628 grpRun->SetSurfaceAtmosPressure(sensor);
1632 Bool_t AliTPCcalibDB::CreateGUITree(Int_t run, const char* filename)
1635 // Create a gui tree for run number 'run'
1638 if (!AliCDBManager::Instance()->GetDefaultStorage()){
1639 AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
1640 MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
1644 AliTPCcalibDB *db=AliTPCcalibDB::Instance();
1645 // retrieve cal pad objects
1647 db->CreateGUITree(filename);
1651 Bool_t AliTPCcalibDB::CreateGUITree(const char* filename){
1655 if (!AliCDBManager::Instance()->GetDefaultStorage()){
1656 AliError("Default Storage not set. Cannot create calibration Tree!");
1659 UpdateNonRec(); // load all infromation now
1661 AliTPCPreprocessorOnline prep;
1662 //noise and pedestals
1663 if (GetPedestals()) prep.AddComponent(new AliTPCCalPad(*(GetPedestals())));
1664 if (GetPadNoise() ) prep.AddComponent(new AliTPCCalPad(*(GetPadNoise())));
1666 if (GetPulserTmean()) prep.AddComponent(new AliTPCCalPad(*(GetPulserTmean())));
1667 if (GetPulserTrms() ) prep.AddComponent(new AliTPCCalPad(*(GetPulserTrms())));
1668 if (GetPulserQmean()) prep.AddComponent(new AliTPCCalPad(*(GetPulserQmean())));
1670 if (GetCETmean()) prep.AddComponent(new AliTPCCalPad(*(GetCETmean())));
1671 if (GetCETrms() ) prep.AddComponent(new AliTPCCalPad(*(GetCETrms())));
1672 if (GetCEQmean()) prep.AddComponent(new AliTPCCalPad(*(GetCEQmean())));
1674 if (GetALTROAcqStart() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROAcqStart() )));
1675 if (GetALTROZsThr() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROZsThr() )));
1676 if (GetALTROFPED() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROFPED() )));
1677 if (GetALTROAcqStop() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROAcqStop() )));
1678 if (GetALTROMasked() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROMasked() )));
1680 AliTPCdataQA *dataQA=GetDataQA();
1682 if (dataQA->GetNLocalMaxima())
1683 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNLocalMaxima())));
1684 if (dataQA->GetMaxCharge())
1685 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetMaxCharge())));
1686 if (dataQA->GetMeanCharge())
1687 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetMeanCharge())));
1688 if (dataQA->GetNoThreshold())
1689 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNoThreshold())));
1690 if (dataQA->GetNTimeBins())
1691 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNTimeBins())));
1692 if (dataQA->GetNPads())
1693 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNPads())));
1694 if (dataQA->GetTimePosition())
1695 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetTimePosition())));
1699 TString file(filename);
1700 if (file.IsNull()) file=Form("guiTreeRun_%lld.root",fRun);
1701 prep.DumpToFile(file.Data());
1705 Bool_t AliTPCcalibDB::CreateRefFile(Int_t run, const char* filename)
1708 // Create a gui tree for run number 'run'
1711 if (!AliCDBManager::Instance()->GetDefaultStorage()){
1712 AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
1713 MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
1716 TString file(filename);
1717 if (file.IsNull()) file=Form("RefCalPads_%d.root",run);
1718 TDirectory *currDir=gDirectory;
1720 AliTPCcalibDB *db=AliTPCcalibDB::Instance();
1721 // retrieve cal pad objects
1724 TFile f(file.Data(),"recreate");
1725 //noise and pedestals
1726 db->GetPedestals()->Write("Pedestals");
1727 db->GetPadNoise()->Write("PadNoise");
1729 db->GetPulserTmean()->Write("PulserTmean");
1730 db->GetPulserTrms()->Write("PulserTrms");
1731 db->GetPulserQmean()->Write("PulserQmean");
1733 db->GetCETmean()->Write("CETmean");
1734 db->GetCETrms()->Write("CETrms");
1735 db->GetCEQmean()->Write("CEQmean");
1737 db->GetALTROAcqStart() ->Write("ALTROAcqStart");
1738 db->GetALTROZsThr() ->Write("ALTROZsThr");
1739 db->GetALTROFPED() ->Write("ALTROFPED");
1740 db->GetALTROAcqStop() ->Write("ALTROAcqStop");
1741 db->GetALTROMasked() ->Write("ALTROMasked");
1750 Double_t AliTPCcalibDB::GetVDriftCorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
1752 // Get time dependent drift velocity correction
1753 // multiplication factor vd = vdnom *(1+vdriftcorr)
1755 // mode determines the algorith how to combine the Laser Track, LaserCE and physics tracks
1756 // timestamp - timestamp
1758 // side - the drift velocity per side (possible for laser and CE)
1760 // Notice - Extrapolation outside of calibration range - using constant function
1763 // mode 1 automatic mode - according to the distance to the valid calibration
1765 Double_t deltaP=0, driftP=0, wP = 0.;
1766 Double_t deltaITS=0,driftITS=0, wITS= 0.;
1767 Double_t deltaLT=0, driftLT=0, wLT = 0.;
1768 Double_t deltaCE=0, driftCE=0, wCE = 0.;
1769 driftP = fDButil->GetVDriftTPC(deltaP,run,timeStamp);
1770 driftITS= fDButil->GetVDriftTPCITS(deltaITS,run,timeStamp);
1771 driftCE = fDButil->GetVDriftTPCCE(deltaCE, run,timeStamp,36000,2);
1772 driftLT = fDButil->GetVDriftTPCLaserTracks(deltaLT,run,timeStamp,36000,2);
1773 deltaITS = TMath::Abs(deltaITS);
1774 deltaP = TMath::Abs(deltaP);
1775 deltaLT = TMath::Abs(deltaLT);
1776 deltaCE = TMath::Abs(deltaCE);
1778 const Double_t kEpsilon=0.00000000001;
1779 const Double_t kdeltaT=360.; // 10 minutes
1780 wITS = 64.*kdeltaT/(deltaITS +kdeltaT);
1781 wLT = 16.*kdeltaT/(deltaLT +kdeltaT);
1782 wP = 0. *kdeltaT/(deltaP +kdeltaT);
1783 wCE = 1. *kdeltaT/(deltaCE +kdeltaT);
1786 if (TMath::Abs(driftP)<kEpsilon) wP=0; // invalid calibration
1787 if (TMath::Abs(driftITS)<kEpsilon)wITS=0; // invalid calibration
1788 if (TMath::Abs(driftLT)<kEpsilon) wLT=0; // invalid calibration
1789 if (TMath::Abs(driftCE)<kEpsilon) wCE=0; // invalid calibration
1790 if (wP+wITS+wLT+wCE<kEpsilon) return 0;
1791 result = (driftP*wP+driftITS*wITS+driftLT*wLT+driftCE*wCE)/(wP+wITS+wLT+wCE);
1797 Double_t AliTPCcalibDB::GetTime0CorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
1799 // Get time dependent time 0 (trigger delay in cm) correction
1800 // additive correction time0 = time0+ GetTime0CorrectionTime
1801 // Value etracted combining the vdrift correction using laser tracks and CE and the physics track matchin
1803 // mode determines the algorith how to combine the Laser Track and physics tracks
1804 // timestamp - timestamp
1806 // side - the drift velocity per side (possible for laser and CE)
1808 // Notice - Extrapolation outside of calibration range - using constant function
1813 result=fDButil->GetTriggerOffsetTPC(run,timeStamp);
1814 result *=fParam->GetZLength();
1819 result= -fDButil->GetTime0TPCITS(dist, run, timeStamp)*fParam->GetDriftV()/1000000.;
1828 Double_t AliTPCcalibDB::GetVDriftCorrectionGy(Int_t timeStamp, Int_t run, Int_t side, Int_t /*mode*/){
1830 // Get global y correction drift velocity correction factor
1831 // additive factor vd = vdnom*(1+GetVDriftCorrectionGy *gy)
1832 // Value etracted combining the vdrift correction using laser tracks and CE
1834 // mode determines the algorith how to combine the Laser Track, LaserCE
1835 // timestamp - timestamp
1837 // side - the drift velocity gy correction per side (CE and Laser tracks)
1839 // Notice - Extrapolation outside of calibration range - using constant function
1841 if (run<=0 && fTransform) run = fTransform->GetCurrentRunNumber();
1842 UpdateRunInformations(run,kFALSE);
1843 TObjArray *array =AliTPCcalibDB::Instance()->GetTimeVdriftSplineRun(run);
1844 if (!array) return 0;
1845 TGraphErrors *laserA= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
1846 TGraphErrors *laserC= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
1849 if (laserA && laserC){
1850 result= (laserA->Eval(timeStamp)+laserC->Eval(timeStamp))*0.5;
1852 if (laserA && side==0){
1853 result = (laserA->Eval(timeStamp));
1855 if (laserC &&side==1){
1856 result = (laserC->Eval(timeStamp));
1858 return -result/250.; //normalized before
1861 AliTPCCalPad* AliTPCcalibDB::MakeDeadMap(Double_t notInMap, const char* nameMappingFile) {
1863 // Read list of active DDLs from OCDB entry
1864 // Generate and return AliTPCCalPad containing 1 for all pads in active DDLs,
1865 // 0 for all pads in non-active DDLs.
1866 // For DDLs with missing status information (no DCS input point to Shuttle),
1867 // the value of the AliTPCCalPad entry is determined by the parameter
1868 // notInMap (default value 1)
1872 TFile *fileMapping = new TFile(nameMappingFile, "read");
1873 AliTPCmapper *mapping = (AliTPCmapper*) fileMapping->Get("tpcMapping");
1875 sprintf(chinfo,"Failed to get mapping object from %s. ...\n", nameMappingFile);
1880 AliTPCCalPad *deadMap = new AliTPCCalPad("deadMap","deadMap");
1882 AliError("Failed to allocate dead map AliTPCCalPad");
1886 /// get list of active DDLs from OCDB entry
1888 if (!fALTROConfigData ) {
1889 AliError("No ALTRO config OCDB entry available");
1892 TMap *activeDDL = (TMap*)fALTROConfigData->FindObject("DDLArray");
1893 TObjString *ddlArray=0;
1895 ddlArray = (TObjString*)activeDDL->GetValue("DDLArray");
1897 AliError("Empty list of active DDLs in OCDB entry");
1901 AliError("List of active DDLs not available in OCDB entry");
1904 TString arrDDL=ddlArray->GetString();
1905 Int_t offset = mapping->GetTpcDdlOffset();
1907 for (Int_t i=0; i<mapping->GetNumDdl(); i++) {
1909 Int_t patch = mapping->GetPatchFromEquipmentID(idDDL);
1910 Int_t roc=mapping->GetRocFromEquipmentID(idDDL);
1911 AliTPCCalROC *calRoc=deadMap->GetCalROC(roc);
1913 for ( Int_t branch = 0; branch < 2; branch++ ) {
1914 for ( Int_t fec = 0; fec < mapping->GetNfec(patch, branch); fec++ ) {
1915 for ( Int_t altro = 0; altro < 8; altro++ ) {
1916 for ( Int_t channel = 0; channel < 16; channel++ ) {
1917 Int_t hwadd = mapping->CodeHWAddress(branch, fec, altro, channel);
1918 Int_t row = mapping->GetPadRow(patch, hwadd); // row in a ROC (IROC or OROC)
1919 // Int_t globalrow = mapping.GetGlobalPadRow(patch, hwadd); // row in full sector (IROC plus OROC)
1920 Int_t pad = mapping->GetPad(patch, hwadd);
1921 if (!TString(arrDDL[i]).IsDigit()) {
1924 active=TString(arrDDL[i]).Atof();
1926 calRoc->SetValue(row,pad,active);
1927 } // end channel for loop
1928 } // end altro for loop
1929 } // end fec for loop
1930 } // end branch for loop
1932 } // end loop on active DDLs
1938 AliTPCCorrection * AliTPCcalibDB::GetTPCComposedCorrection(Float_t field) const{
1940 // GetComposed correction for given field setting
1942 if (!fComposedCorrectionArray) return 0;
1943 if (field>0.1) return (AliTPCCorrection *)fComposedCorrectionArray->At(1);
1944 if (field<-0.1) return (AliTPCCorrection *)fComposedCorrectionArray->At(2);
1945 return (AliTPCCorrection *)fComposedCorrectionArray->At(0);