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),
195 fTimeGainSplinesArray(100000),
196 fGRPArray(100000), //! array of GRPs - per run - JUST for calibration studies
197 fGRPMaps(100000), //! array of GRPs - per run - JUST for calibration studies
198 fGoofieArray(100000), //! array of GOOFIE values -per run - Just for calibration studies
199 fVoltageArray(100000),
200 fTemperatureArray(100000), //! array of temperature sensors - per run - Just for calibration studies
201 fVdriftArray(100000), //! array of v drift interfaces
202 fDriftCorrectionArray(100000), //! array of drift correction
203 fRunList(100000), //! run list - indicates try to get the run param
212 Update(); // temporary
215 AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& ):
224 fComposedCorrection(0),
237 fTimeGainSplinesArray(100000),
238 fGRPArray(0), //! array of GRPs - per run - JUST for calibration studies
239 fGRPMaps(0), //! array of GRPs - per run - JUST for calibration studies
240 fGoofieArray(0), //! array of GOOFIE values -per run - Just for calibration studies
242 fTemperatureArray(0), //! array of temperature sensors - per run - Just for calibration studies
243 fVdriftArray(0), //! array of v drift interfaces
244 fDriftCorrectionArray(0), //! array of v drift interfaces
245 fRunList(0), //! run list - indicates try to get the run param
250 // Copy constructor invalid -- singleton implementation
252 Error("copy constructor","invalid -- singleton implementation");
255 AliTPCcalibDB& AliTPCcalibDB::operator= (const AliTPCcalibDB& )
258 // Singleton implementation - no assignment operator
260 Error("operator =", "assignment operator not implemented");
266 //_____________________________________________________________________________
267 AliTPCcalibDB::~AliTPCcalibDB()
274 AliTPCCalPad* AliTPCcalibDB::GetDistortionMap(Int_t i) const {
276 // get distortion map - due E field distortions
278 return (fDistortionMap) ? (AliTPCCalPad*)fDistortionMap->At(i):0;
281 //_____________________________________________________________________________
282 AliCDBEntry* AliTPCcalibDB::GetCDBEntry(const char* cdbPath)
285 // Retrieves an entry with path <cdbPath> from the CDB.
289 AliCDBEntry* entry = AliCDBManager::Instance()->Get(cdbPath, fRun);
292 sprintf(chinfo,"AliTPCcalibDB: Failed to get entry:\t%s ", cdbPath);
300 //_____________________________________________________________________________
301 void AliTPCcalibDB::SetRun(Long64_t run)
304 // Sets current run number. Calibration data is read from the corresponding file.
314 void AliTPCcalibDB::Update(){
316 // cache the OCDB entries for simulation, reconstruction, calibration
319 AliCDBEntry * entry=0;
320 Bool_t cdbCache = AliCDBManager::Instance()->GetCacheFlag(); // save cache status
321 AliCDBManager::Instance()->SetCacheFlag(kTRUE); // activate CDB cache
322 fDButil = new AliTPCcalibDButil;
324 entry = GetCDBEntry("TPC/Calib/PadGainFactor");
326 //if (fPadGainFactor) delete fPadGainFactor;
327 entry->SetOwner(kTRUE);
328 fPadGainFactor = (AliTPCCalPad*)entry->GetObject();
330 AliFatal("TPC - Missing calibration entry TPC/Calib/PadGainFactor")
333 entry = GetCDBEntry("TPC/Calib/TimeGain");
335 //if (fTimeGainSplines) delete fTimeGainSplines;
336 entry->SetOwner(kTRUE);
337 fTimeGainSplines = (TObjArray*)entry->GetObject();
339 AliFatal("TPC - Missing calibration entry TPC/Calib/Timegain")
342 entry = GetCDBEntry("TPC/Calib/GainFactorDedx");
344 entry->SetOwner(kTRUE);
345 fDedxGainFactor = (AliTPCCalPad*)entry->GetObject();
347 AliFatal("TPC - Missing calibration entry TPC/Calib/gainFactordEdx")
350 entry = GetCDBEntry("TPC/Calib/PadTime0");
352 //if (fPadTime0) delete fPadTime0;
353 entry->SetOwner(kTRUE);
354 fPadTime0 = (AliTPCCalPad*)entry->GetObject();
356 AliFatal("TPC - Missing calibration entry")
359 entry = GetCDBEntry("TPC/Calib/Distortion");
361 //if (fPadTime0) delete fPadTime0;
362 entry->SetOwner(kTRUE);
363 fDistortionMap =dynamic_cast<TObjArray*>(entry->GetObject());
365 //AliFatal("TPC - Missing calibration entry")
371 entry = GetCDBEntry("TPC/Calib/PadNoise");
373 //if (fPadNoise) delete fPadNoise;
374 entry->SetOwner(kTRUE);
375 fPadNoise = (AliTPCCalPad*)entry->GetObject();
377 AliFatal("TPC - Missing calibration entry")
380 entry = GetCDBEntry("TPC/Calib/Pedestals");
382 //if (fPedestals) delete fPedestals;
383 entry->SetOwner(kTRUE);
384 fPedestals = (AliTPCCalPad*)entry->GetObject();
387 entry = GetCDBEntry("TPC/Calib/Temperature");
389 //if (fTemperature) delete fTemperature;
390 entry->SetOwner(kTRUE);
391 fTemperature = (AliTPCSensorTempArray*)entry->GetObject();
394 entry = GetCDBEntry("TPC/Calib/Parameters");
396 //if (fPadNoise) delete fPadNoise;
397 entry->SetOwner(kTRUE);
398 fParam = (AliTPCParam*)(entry->GetObject()->Clone());
400 AliFatal("TPC - Missing calibration entry TPC/Calib/Parameters")
403 entry = GetCDBEntry("TPC/Calib/ClusterParam");
405 entry->SetOwner(kTRUE);
406 fClusterParam = (AliTPCClusterParam*)(entry->GetObject()->Clone());
408 AliFatal("TPC - Missing calibration entry")
411 //ALTRO configuration data
412 entry = GetCDBEntry("TPC/Calib/AltroConfig");
414 entry->SetOwner(kTRUE);
415 fALTROConfigData=(TObjArray*)(entry->GetObject());
417 AliFatal("TPC - Missing calibration entry")
420 //Calibration Pulser data
421 entry = GetCDBEntry("TPC/Calib/Pulser");
423 entry->SetOwner(kTRUE);
424 fPulserData=(TObjArray*)(entry->GetObject());
428 entry = GetCDBEntry("TPC/Calib/CE");
430 entry->SetOwner(kTRUE);
431 fCEData=(TObjArray*)(entry->GetObject());
433 //RAW calibration data
434 // entry = GetCDBEntry("TPC/Calib/Raw");
436 entry = GetCDBEntry("TPC/Calib/Mapping");
438 //if (fPadNoise) delete fPadNoise;
439 entry->SetOwner(kTRUE);
440 TObjArray * array = dynamic_cast<TObjArray*>(entry->GetObject());
441 if (array && array->GetEntriesFast()==6){
442 fMapping = new AliTPCAltroMapping*[6];
443 for (Int_t i=0; i<6; i++){
444 fMapping[i] = dynamic_cast<AliTPCAltroMapping*>(array->At(i));
449 //CTP calibration data
450 entry = GetCDBEntry("GRP/CTP/CTPtiming");
452 //entry->SetOwner(kTRUE);
453 fCTPTimeParams=dynamic_cast<AliCTPTimeParams*>(entry->GetObject());
455 AliError("TPC - Missing calibration entry")
457 //TPC space point correction data
458 entry = GetCDBEntry("TPC/Calib/Correction");
460 //entry->SetOwner(kTRUE);
461 fComposedCorrection=dynamic_cast<AliTPCCorrection*>(entry->GetObject());
462 fComposedCorrection->Init();
464 AliError("TPC - Missing calibration entry- TPC/Calib/Correction")
469 fTransform=new AliTPCTransform();
470 fTransform->SetCurrentRun(AliCDBManager::Instance()->GetRun());
474 AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
477 void AliTPCcalibDB::UpdateNonRec(){
479 // Update/Load the parameters which are important for QA studies
480 // and not used yet for the reconstruction
482 //RAW calibration data
483 AliCDBEntry * entry=0;
484 entry = GetCDBEntry("TPC/Calib/Raw");
486 entry->SetOwner(kTRUE);
487 TObjArray *arr=(TObjArray*)(entry->GetObject());
488 if (arr) fCalibRaw=(AliTPCCalibRaw*)arr->At(0);
490 //QA calibration data
491 entry = GetCDBEntry("TPC/Calib/QA");
493 entry->SetOwner(kTRUE);
494 fDataQA=dynamic_cast<AliTPCdataQA*>(entry->GetObject());
498 entry = AliCDBManager::Instance()->Get("TPC/Calib/HighVoltage",fRun);
500 fVoltageArray.AddAt(entry->GetObject(),fRun);
508 void AliTPCcalibDB::CreateObjectList(const Char_t *filename, TObjArray *calibObjects)
511 // Create calibration objects and read contents from OCDB
513 if ( calibObjects == 0x0 ) return;
516 if ( !in.is_open() ){
517 fprintf(stderr,"Error: cannot open list file '%s'", filename);
521 AliTPCCalPad *calPad=0x0;
527 TObjArray *arrFileLine = sFile.Tokenize("\n");
529 TIter nextLine(arrFileLine);
531 TObjString *sObjLine=0x0;
532 while ( (sObjLine = (TObjString*)nextLine()) ){
533 TString sLine(sObjLine->GetString());
535 TObjArray *arrNextCol = sLine.Tokenize("\t");
537 TObjString *sObjType = (TObjString*)(arrNextCol->At(0));
538 TObjString *sObjFileName = (TObjString*)(arrNextCol->At(1));
540 if ( !sObjType || ! sObjFileName ) continue;
541 TString sType(sObjType->GetString());
542 TString sFileName(sObjFileName->GetString());
543 printf("%s\t%s\n",sType.Data(),sFileName.Data());
545 TFile *fIn = TFile::Open(sFileName);
547 fprintf(stderr,"File not found: '%s'", sFileName.Data());
551 if ( sType == "CE" ){
552 AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
554 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());
555 calPad->SetNameTitle("CETmean","CETmean");
556 calibObjects->Add(calPad);
558 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());
559 calPad->SetNameTitle("CEQmean","CEQmean");
560 calibObjects->Add(calPad);
562 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
563 calPad->SetNameTitle("CETrms","CETrms");
564 calibObjects->Add(calPad);
566 } else if ( sType == "Pulser") {
567 AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
569 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());
570 calPad->SetNameTitle("PulserTmean","PulserTmean");
571 calibObjects->Add(calPad);
573 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());
574 calPad->SetNameTitle("PulserQmean","PulserQmean");
575 calibObjects->Add(calPad);
577 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
578 calPad->SetNameTitle("PulserTrms","PulserTrms");
579 calibObjects->Add(calPad);
581 } else if ( sType == "Pedestals") {
582 AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
584 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());
585 calPad->SetNameTitle("Pedestals","Pedestals");
586 calibObjects->Add(calPad);
588 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());
589 calPad->SetNameTitle("Noise","Noise");
590 calibObjects->Add(calPad);
593 fprintf(stderr,"Undefined Type: '%s'",sType.Data());
602 void AliTPCcalibDB::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
604 // Write a tree with all available information
605 // if mapFileName is specified, the Map information are also written to the tree
606 // pads specified in outlierPad are not used for calculating statistics
607 // - the same function as AliTPCCalPad::MakeTree -
609 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
611 TObjArray* mapIROCs = 0;
612 TObjArray* mapOROCs = 0;
613 TVectorF *mapIROCArray = 0;
614 TVectorF *mapOROCArray = 0;
615 Int_t mapEntries = 0;
616 TString* mapNames = 0;
619 TFile mapFile(mapFileName, "read");
621 TList* listOfROCs = mapFile.GetListOfKeys();
622 mapEntries = listOfROCs->GetEntries()/2;
623 mapIROCs = new TObjArray(mapEntries*2);
624 mapOROCs = new TObjArray(mapEntries*2);
625 mapIROCArray = new TVectorF[mapEntries];
626 mapOROCArray = new TVectorF[mapEntries];
628 mapNames = new TString[mapEntries];
629 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
630 TString nameROC(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
631 nameROC.Remove(nameROC.Length()-4, 4);
632 mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "IROC").Data()), ivalue);
633 mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "OROC").Data()), ivalue);
634 mapNames[ivalue].Append(nameROC);
637 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
638 mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
639 mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
641 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
642 (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
643 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
644 (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
647 } // if (mapFileName)
649 TTreeSRedirector cstream(fileName);
650 Int_t arrayEntries = array->GetEntries();
652 TString* names = new TString[arrayEntries];
653 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
654 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
656 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
658 // get statistic for given sector
660 TVectorF median(arrayEntries);
661 TVectorF mean(arrayEntries);
662 TVectorF rms(arrayEntries);
663 TVectorF ltm(arrayEntries);
664 TVectorF ltmrms(arrayEntries);
665 TVectorF medianWithOut(arrayEntries);
666 TVectorF meanWithOut(arrayEntries);
667 TVectorF rmsWithOut(arrayEntries);
668 TVectorF ltmWithOut(arrayEntries);
669 TVectorF ltmrmsWithOut(arrayEntries);
671 TVectorF *vectorArray = new TVectorF[arrayEntries];
672 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
673 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
675 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
676 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
677 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
678 AliTPCCalROC* outlierROC = 0;
679 if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
681 median[ivalue] = calROC->GetMedian();
682 mean[ivalue] = calROC->GetMean();
683 rms[ivalue] = calROC->GetRMS();
684 Double_t ltmrmsValue = 0;
685 ltm[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction);
686 ltmrms[ivalue] = ltmrmsValue;
688 medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
689 meanWithOut[ivalue] = calROC->GetMean(outlierROC);
690 rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
692 ltmWithOut[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction, outlierROC);
693 ltmrmsWithOut[ivalue] = ltmrmsValue;
702 medianWithOut[ivalue] = 0.;
703 meanWithOut[ivalue] = 0.;
704 rmsWithOut[ivalue] = 0.;
705 ltmWithOut[ivalue] = 0.;
706 ltmrmsWithOut[ivalue] = 0.;
711 // fill vectors of variable per pad
713 TVectorF *posArray = new TVectorF[8];
714 for (Int_t ivalue = 0; ivalue < 8; ivalue++)
715 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
717 Float_t posG[3] = {0};
718 Float_t posL[3] = {0};
720 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
721 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
722 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
723 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
724 posArray[0][ichannel] = irow;
725 posArray[1][ichannel] = ipad;
726 posArray[2][ichannel] = posL[0];
727 posArray[3][ichannel] = posL[1];
728 posArray[4][ichannel] = posG[0];
729 posArray[5][ichannel] = posG[1];
730 posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
731 posArray[7][ichannel] = ichannel;
733 // loop over array containing AliTPCCalPads
734 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
735 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
736 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
738 (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
740 (vectorArray[ivalue])[ichannel] = 0;
746 cstream << "calPads" <<
747 "sector=" << isector;
749 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
750 cstream << "calPads" <<
751 (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
752 (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
753 (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
754 (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
755 (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
757 cstream << "calPads" <<
758 (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
759 (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
760 (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
761 (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
762 (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
766 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
767 cstream << "calPads" <<
768 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
772 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
774 cstream << "calPads" <<
775 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
777 cstream << "calPads" <<
778 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
782 cstream << "calPads" <<
783 "row.=" << &posArray[0] <<
784 "pad.=" << &posArray[1] <<
785 "lx.=" << &posArray[2] <<
786 "ly.=" << &posArray[3] <<
787 "gx.=" << &posArray[4] <<
788 "gy.=" << &posArray[5] <<
789 "rpad.=" << &posArray[6] <<
790 "channel.=" << &posArray[7];
792 cstream << "calPads" <<
796 delete[] vectorArray;
804 delete[] mapIROCArray;
805 delete[] mapOROCArray;
810 Int_t AliTPCcalibDB::GetRCUTriggerConfig() const
813 // return the RCU trigger configuration register
815 TMap *map=GetRCUconfig();
817 TVectorF *v=(TVectorF*)map->GetValue("TRGCONF_TRG_MODE");
819 for (Int_t i=0; i<v->GetNrows(); ++i){
820 Float_t newmode=v->GetMatrixArray()[i];
822 if (mode>-1&&newmode!=mode) AliWarning("Found different RCU trigger configurations!!!");
829 Bool_t AliTPCcalibDB::IsTrgL0()
832 // return if the FEE readout was triggered on L0
834 Int_t mode=GetRCUTriggerConfig();
835 if (mode<0) return kFALSE;
839 Bool_t AliTPCcalibDB::IsTrgL1()
842 // return if the FEE readout was triggered on L1
844 Int_t mode=GetRCUTriggerConfig();
845 if (mode<0) return kFALSE;
849 void AliTPCcalibDB::RegisterExB(Int_t index, Float_t bz, Bool_t bdelete){
851 // Register static ExB correction map
852 // index - registration index - used for visualization
853 // bz - bz field in kGaus
855 // Float_t factor = bz/(-5.); // default b filed in Cheb with minus sign
856 Float_t factor = bz/(5.); // default b filed in Cheb with minus sign
857 // was chenged in the Revision ???? (Ruben can you add here number)
859 AliMagF* bmap = new AliMagF("MapsExB","MapsExB", factor,TMath::Sign(1.f,factor),AliMagF::k5kG);
861 AliTPCExBFirst *exb = new AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
862 AliTPCExB::SetInstance(exb);
867 AliTPCExB::RegisterField(index,bmap);
869 if (index>=fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
870 fgExBArray.AddAt(exb,index);
874 AliTPCExB* AliTPCcalibDB::GetExB(Float_t bz, Bool_t deleteB) {
876 // bz filed in KGaus not in tesla
877 // Get ExB correction map
878 // if doesn't exist - create it
880 Int_t index = TMath::Nint(5+bz);
881 if (index>fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
882 if (!fgExBArray.At(index)) AliTPCcalibDB::RegisterExB(index,bz,deleteB);
883 return (AliTPCExB*)fgExBArray.At(index);
887 void AliTPCcalibDB::SetExBField(Float_t bz){
889 // Set magnetic filed for ExB correction
891 fExB = GetExB(bz,kFALSE);
894 void AliTPCcalibDB::SetExBField(const AliMagF* bmap){
896 // Set magnetic field for ExB correction
898 AliTPCExBFirst *exb = new AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
899 AliTPCExB::SetInstance(exb);
907 void AliTPCcalibDB::UpdateRunInformations( Int_t run, Bool_t force){
909 // - > Don't use it for reconstruction - Only for Calibration studies
913 AliCDBEntry * entry = 0;
914 if (run>= fRunList.fN){
915 fRunList.Set(run*2+1);
916 fGRPArray.Expand(run*2+1);
917 fGRPMaps.Expand(run*2+1);
918 fGoofieArray.Expand(run*2+1);
919 fVoltageArray.Expand(run*2+1);
920 fTemperatureArray.Expand(run*2+1);
921 fVdriftArray.Expand(run*2+1);
922 fDriftCorrectionArray.Expand(run*2+1);
923 fTimeGainSplinesArray.Expand(run*2+1);
926 fALTROConfigData->Expand(run*2+1); // ALTRO configuration data
927 fPulserData->Expand(run*2+1); // Calibration Pulser data
928 fCEData->Expand(run*2+1); // CE data
929 if (!fTimeGainSplines) fTimeGainSplines = new TObjArray(run*2+1);
930 fTimeGainSplines->Expand(run*2+1); // Array of AliSplineFits: at 0 MIP position in
932 if (fRunList[run]>0 &&force==kFALSE) return;
934 fRunList[run]=1; // sign as used
937 entry = AliCDBManager::Instance()->Get("GRP/GRP/Data",run);
939 AliGRPObject * grpRun = dynamic_cast<AliGRPObject*>(entry->GetObject());
941 TMap* map = dynamic_cast<TMap*>(entry->GetObject());
943 //grpRun = new AliGRPObject;
944 //grpRun->ReadValuesFromMap(map);
945 grpRun = MakeGRPObjectFromMap(map);
947 fGRPMaps.AddAt(map,run);
950 fGRPArray.AddAt(grpRun,run);
952 entry = AliCDBManager::Instance()->Get("TPC/Calib/Goofie",run);
954 fGoofieArray.AddAt(entry->GetObject(),run);
959 entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeGain",run);
961 fTimeGainSplinesArray.AddAt(entry->GetObject(),run);
963 AliFatal("TPC - Missing calibration entry TimeGain")
966 entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeDrift",run);
968 fDriftCorrectionArray.AddAt(entry->GetObject(),run);
970 AliFatal("TPC - Missing calibration entry TimeDrift")
973 entry = AliCDBManager::Instance()->Get("TPC/Calib/Temperature",run);
975 fTemperatureArray.AddAt(entry->GetObject(),run);
977 //apply fDButil filters
979 fDButil->UpdateFromCalibDB();
980 if (fTemperature) fDButil->FilterTemperature(fTemperature);
982 AliDCSSensor * press = GetPressureSensor(run,0);
983 AliTPCSensorTempArray * temp = GetTemperatureSensor(run);
986 accept = fDButil->FilterTemperature(temp)>0.1;
989 const Double_t kMinP=950.;
990 const Double_t kMaxP=1050.;
991 const Double_t kMaxdP=10.;
992 const Double_t kSigmaCut=4.;
993 fDButil->FilterSensor(press,kMinP,kMaxP,kMaxdP,kSigmaCut);
994 if (press->GetFit()==0) accept=kFALSE;
996 if (press && temp &&accept){
997 AliTPCCalibVdrift * vdrift = new AliTPCCalibVdrift(temp, press,0);
998 fVdriftArray.AddAt(vdrift,run);
1000 fDButil->FilterCE(120., 3., 4.,0);
1001 fDButil->FilterTracks(run, 10.,0);
1005 Float_t AliTPCcalibDB::GetGain(Int_t sector, Int_t row, Int_t pad){
1007 // Get Gain factor for given pad
1009 AliTPCCalPad *calPad = Instance()->fDedxGainFactor;;
1010 if (!calPad) return 0;
1011 return calPad->GetCalROC(sector)->GetValue(row,pad);
1014 AliSplineFit* AliTPCcalibDB::GetVdriftSplineFit(const char* name, Int_t run){
1016 // GetDrift velocity spline fit
1018 TObjArray *arr=GetTimeVdriftSplineRun(run);
1020 return dynamic_cast<AliSplineFit*>(arr->FindObject(name));
1023 AliSplineFit* AliTPCcalibDB::CreateVdriftSplineFit(const char* graphName, Int_t run){
1025 // create spline fit from the drift time graph in TimeDrift
1027 TObjArray *arr=GetTimeVdriftSplineRun(run);
1029 TGraph *graph=dynamic_cast<TGraph*>(arr->FindObject(graphName));
1030 if (!graph) return 0;
1031 AliSplineFit *fit = new AliSplineFit();
1032 fit->SetGraph(graph);
1033 fit->SetMinPoints(graph->GetN()+1);
1034 fit->InitKnots(graph,2,0,0.001);
1039 AliGRPObject *AliTPCcalibDB::GetGRP(Int_t run){
1041 // Get GRP object for given run
1043 if (run>= ((Instance()->fGRPArray)).GetEntriesFast()){
1044 Instance()->UpdateRunInformations(run);
1046 AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>((Instance()->fGRPArray).At(run));
1048 Instance()->UpdateRunInformations(run);
1049 grpRun = dynamic_cast<AliGRPObject *>(Instance()->fGRPArray.At(run));
1050 if (!grpRun) return 0;
1055 TMap * AliTPCcalibDB::GetGRPMap(Int_t run){
1057 // Get GRP map for given run
1059 TMap * grpRun = dynamic_cast<TMap *>((Instance()->fGRPMaps).At(run));
1061 Instance()->UpdateRunInformations(run);
1062 grpRun = dynamic_cast<TMap *>(Instance()->fGRPMaps.At(run));
1063 if (!grpRun) return 0;
1069 AliDCSSensor * AliTPCcalibDB::GetPressureSensor(Int_t run, Int_t type){
1071 // Get Pressure sensor
1073 // type = 0 - Cavern pressure
1074 // 1 - Suface pressure
1075 // First try to get if trom map - if existing (Old format of data storing)
1079 TMap *map = GetGRPMap(run);
1081 AliDCSSensor * sensor = 0;
1083 if (type==0) osensor = ((*map)("fCavernPressure"));
1084 if (type==1) osensor = ((*map)("fP2Pressure"));
1085 sensor =dynamic_cast<AliDCSSensor *>(osensor);
1086 if (sensor) return sensor;
1089 // If not map try to get it from the GRPObject
1091 AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run));
1093 UpdateRunInformations(run);
1094 grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run));
1095 if (!grpRun) return 0;
1097 AliDCSSensor * sensor = grpRun->GetCavernAtmosPressure();
1098 if (type==1) sensor = grpRun->GetSurfaceAtmosPressure();
1102 AliTPCSensorTempArray * AliTPCcalibDB::GetTemperatureSensor(Int_t run){
1104 // Get temperature sensor array
1106 AliTPCSensorTempArray * tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
1108 UpdateRunInformations(run);
1109 tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
1115 TObjArray * AliTPCcalibDB::GetTimeGainSplinesRun(Int_t run){
1117 // Get temperature sensor array
1119 TObjArray * gainSplines = (TObjArray *)fTimeGainSplinesArray.At(run);
1121 UpdateRunInformations(run);
1122 gainSplines = (TObjArray *)fTimeGainSplinesArray.At(run);
1127 TObjArray * AliTPCcalibDB::GetTimeVdriftSplineRun(Int_t run){
1129 // Get drift spline array
1131 TObjArray * driftSplines = (TObjArray *)fDriftCorrectionArray.At(run);
1132 if (!driftSplines) {
1133 UpdateRunInformations(run);
1134 driftSplines = (TObjArray *)fDriftCorrectionArray.At(run);
1136 return driftSplines;
1139 AliDCSSensorArray * AliTPCcalibDB::GetVoltageSensors(Int_t run){
1141 // Get temperature sensor array
1143 AliDCSSensorArray * voltageArray = (AliDCSSensorArray *)fVoltageArray.At(run);
1144 if (!voltageArray) {
1145 UpdateRunInformations(run);
1146 voltageArray = (AliDCSSensorArray *)fVoltageArray.At(run);
1148 return voltageArray;
1151 AliDCSSensorArray * AliTPCcalibDB::GetGoofieSensors(Int_t run){
1153 // Get temperature sensor array
1155 AliDCSSensorArray * goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
1157 UpdateRunInformations(run);
1158 goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
1165 AliTPCCalibVdrift * AliTPCcalibDB::GetVdrift(Int_t run){
1167 // Get the interface to the the vdrift
1169 AliTPCCalibVdrift * vdrift = (AliTPCCalibVdrift*)fVdriftArray.At(run);
1171 UpdateRunInformations(run);
1172 vdrift= (AliTPCCalibVdrift*)fVdriftArray.At(run);
1177 Float_t AliTPCcalibDB::GetCEdriftTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1180 // GetCE drift time information for 'sector'
1181 // sector 72 is the mean drift time of the A-Side
1182 // sector 73 is the mean drift time of the C-Side
1183 // it timestamp==-1 return mean value
1185 AliTPCcalibDB::Instance()->SetRun(run);
1186 TGraph *gr=AliTPCcalibDB::Instance()->GetCErocTgraph(sector);
1187 if (!gr||sector<0||sector>73) {
1188 if (entries) *entries=0;
1192 if (timeStamp==-1.){
1195 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1197 gr->GetPoint(ipoint,x,y);
1198 if (x<timeStamp) continue;
1206 Float_t AliTPCcalibDB::GetCEchargeTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1209 // GetCE mean charge for 'sector'
1210 // it timestamp==-1 return mean value
1212 AliTPCcalibDB::Instance()->SetRun(run);
1213 TGraph *gr=AliTPCcalibDB::Instance()->GetCErocQgraph(sector);
1214 if (!gr||sector<0||sector>71) {
1215 if (entries) *entries=0;
1219 if (timeStamp==-1.){
1222 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1224 gr->GetPoint(ipoint,x,y);
1225 if (x<timeStamp) continue;
1233 Float_t AliTPCcalibDB::GetDCSSensorValue(AliDCSSensorArray *arr, Int_t timeStamp, const char * sensorName, Int_t sigDigits)
1236 // Get Value for a DCS sensor 'sensorName', run 'run' at time 'timeStamp'
1239 const TString sensorNameString(sensorName);
1240 AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1241 if (!sensor) return val;
1242 //use the dcs graph if possible
1243 TGraph *gr=sensor->GetGraph();
1245 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1247 gr->GetPoint(ipoint,x,y);
1248 Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1249 if (time<timeStamp) continue;
1253 //if val is still 0, test if if the requested time if within 5min of the first/last
1254 //data point. If this is the case return the firs/last entry
1255 //the timestamps might not be syncronised for all calibration types, sometimes a 'pre'
1256 //and 'pos' period is requested. Especially to the HV this is not the case!
1260 gr->GetPoint(0,x,y);
1261 Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1262 if ((time-timeStamp)<5*60) val=y;
1267 gr->GetPoint(gr->GetN()-1,x,y);
1268 Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1269 if ((timeStamp-time)<5*60) val=y;
1272 val=sensor->GetValue(timeStamp);
1275 val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1280 Float_t AliTPCcalibDB::GetDCSSensorMeanValue(AliDCSSensorArray *arr, const char * sensorName, Int_t sigDigits)
1283 // Get mean Value for a DCS sensor 'sensorName' during run 'run'
1286 const TString sensorNameString(sensorName);
1287 AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1288 if (!sensor) return val;
1290 //use dcs graph if it exists
1291 TGraph *gr=sensor->GetGraph();
1295 //if we don't have the dcs graph, try to get some meaningful information
1296 if (!sensor->GetFit()) return val;
1297 Int_t nKnots=sensor->GetFit()->GetKnots();
1298 Double_t tMid=(sensor->GetEndTime()-sensor->GetStartTime())/2.;
1299 for (Int_t iKnot=0;iKnot<nKnots;++iKnot){
1300 if (sensor->GetFit()->GetX()[iKnot]>tMid/3600.) break;
1301 val=(Float_t)sensor->GetFit()->GetY0()[iKnot];
1306 val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1312 Float_t AliTPCcalibDB::GetChamberHighVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits) {
1314 // return the chamber HV for given run and time: 0-35 IROC, 36-72 OROC
1315 // if timeStamp==-1 return mean value
1318 TString sensorName="";
1319 TTimeStamp stamp(timeStamp);
1320 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1321 if (!voltageArray || (sector<0) || (sector>71)) return val;
1322 Char_t sideName='A';
1323 if ((sector/18)%2==1) sideName='C';
1326 sensorName=Form("TPC_ANODE_I_%c%02d_VMEAS",sideName,sector%18);
1329 sensorName=Form("TPC_ANODE_O_%c%02d_0_VMEAS",sideName,sector%18);
1332 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1334 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1338 Float_t AliTPCcalibDB::GetSkirtVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1341 // Get the skirt voltage for 'run' at 'timeStamp' and 'sector': 0-35 IROC, 36-72 OROC
1342 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1343 // if timeStamp==-1 return the mean value for the run
1346 TString sensorName="";
1347 TTimeStamp stamp(timeStamp);
1348 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1349 if (!voltageArray || (sector<0) || (sector>71)) return val;
1350 Char_t sideName='A';
1351 if ((sector/18)%2==1) sideName='C';
1352 sensorName=Form("TPC_SKIRT_%c_VMEAS",sideName);
1354 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1356 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1361 Float_t AliTPCcalibDB::GetCoverVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1364 // Get the cover voltage for run 'run' at time 'timeStamp'
1365 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1366 // if timeStamp==-1 return the mean value for the run
1369 TString sensorName="";
1370 TTimeStamp stamp(timeStamp);
1371 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1372 if (!voltageArray || (sector<0) || (sector>71)) return val;
1373 Char_t sideName='A';
1374 if ((sector/18)%2==1) sideName='C';
1377 sensorName=Form("TPC_COVER_I_%c_VMEAS",sideName);
1380 sensorName=Form("TPC_COVER_O_%c_VMEAS",sideName);
1383 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1385 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1390 Float_t AliTPCcalibDB::GetGGoffsetVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1393 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1394 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1395 // if timeStamp==-1 return the mean value for the run
1398 TString sensorName="";
1399 TTimeStamp stamp(timeStamp);
1400 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1401 if (!voltageArray || (sector<0) || (sector>71)) return val;
1402 Char_t sideName='A';
1403 if ((sector/18)%2==1) sideName='C';
1406 sensorName=Form("TPC_GATE_I_%c_OFF_VMEAS",sideName);
1409 sensorName=Form("TPC_GATE_O_%c_OFF_VMEAS",sideName);
1412 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1414 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1419 Float_t AliTPCcalibDB::GetGGnegVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1422 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1423 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1424 // if timeStamp==-1 return the mean value for the run
1427 TString sensorName="";
1428 TTimeStamp stamp(timeStamp);
1429 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1430 if (!voltageArray || (sector<0) || (sector>71)) return val;
1431 Char_t sideName='A';
1432 if ((sector/18)%2==1) sideName='C';
1435 sensorName=Form("TPC_GATE_I_%c_NEG_VMEAS",sideName);
1438 sensorName=Form("TPC_GATE_O_%c_NEG_VMEAS",sideName);
1441 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1443 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1448 Float_t AliTPCcalibDB::GetGGposVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1451 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1452 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1453 // if timeStamp==-1 return the mean value for the run
1456 TString sensorName="";
1457 TTimeStamp stamp(timeStamp);
1458 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1459 if (!voltageArray || (sector<0) || (sector>71)) return val;
1460 Char_t sideName='A';
1461 if ((sector/18)%2==1) sideName='C';
1464 sensorName=Form("TPC_GATE_I_%c_POS_VMEAS",sideName);
1467 sensorName=Form("TPC_GATE_O_%c_POS_VMEAS",sideName);
1470 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1472 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1477 Float_t AliTPCcalibDB::GetPressure(Int_t timeStamp, Int_t run, Int_t type){
1479 // GetPressure for given time stamp and runt
1481 TTimeStamp stamp(timeStamp);
1482 AliDCSSensor * sensor = Instance()->GetPressureSensor(run,type);
1483 if (!sensor) return 0;
1484 return sensor->GetValue(stamp);
1487 Float_t AliTPCcalibDB::GetL3Current(Int_t run, Int_t statType){
1489 // return L3 current
1490 // stat type is: AliGRPObject::Stats: kMean = 0, kTruncMean = 1, kMedian = 2, kSDMean = 3, kSDMedian = 4
1493 AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1494 if (grp) current=grp->GetL3Current((AliGRPObject::Stats)statType);
1498 Float_t AliTPCcalibDB::GetBz(Int_t run){
1500 // calculate BZ in T from L3 current
1503 Float_t current=AliTPCcalibDB::GetL3Current(run);
1504 if (current>-1) bz=5*current/30000.*.1;
1508 Char_t AliTPCcalibDB::GetL3Polarity(Int_t run) {
1510 // get l3 polarity from GRP
1513 AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1514 if (grp) pol=grp->GetL3Polarity();
1518 TString AliTPCcalibDB::GetRunType(Int_t run){
1520 // return run type from grp
1523 // TString type("UNKNOWN");
1524 AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1525 if (grp) return grp->GetRunType();
1529 Float_t AliTPCcalibDB::GetValueGoofie(Int_t timeStamp, Int_t run, Int_t type){
1531 // GetPressure for given time stamp and runt
1533 TTimeStamp stamp(timeStamp);
1534 AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(run);
1535 if (!goofieArray) return 0;
1536 AliDCSSensor *sensor = goofieArray->GetSensor(type);
1537 return sensor->GetValue(stamp);
1545 Bool_t AliTPCcalibDB::GetTemperatureFit(Int_t timeStamp, Int_t run, Int_t side,TVectorD& fit){
1547 // GetTmeparature fit at parameter for given time stamp
1549 TTimeStamp tstamp(timeStamp);
1550 AliTPCSensorTempArray* tempArray = Instance()->GetTemperatureSensor(run);
1551 if (! tempArray) return kFALSE;
1552 AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
1553 TLinearFitter * fitter = tempMap->GetLinearFitter(3,side,tstamp);
1556 fitter->GetParameters(fit);
1560 if (!fitter) return kFALSE;
1564 Float_t AliTPCcalibDB::GetTemperature(Int_t timeStamp, Int_t run, Int_t side){
1566 // Get mean temperature
1570 GetTemperatureFit(timeStamp,run,0,vec);
1574 GetTemperatureFit(timeStamp,run,0,vec);
1581 Double_t AliTPCcalibDB::GetPTRelative(UInt_t timeSec, Int_t run, Int_t side){
1584 // time - absolute time
1586 // side - 0 - A side 1-C side
1587 AliTPCCalibVdrift * vdrift = Instance()->GetVdrift(run);
1588 if (!vdrift) return 0;
1589 return vdrift->GetPTRelative(timeSec,side);
1592 AliGRPObject * AliTPCcalibDB::MakeGRPObjectFromMap(TMap *map){
1594 // Function to covert old GRP run information from TMap to GRPObject
1596 // TMap * map = AliTPCcalibDB::GetGRPMap(52406);
1598 AliDCSSensor * sensor = 0;
1600 osensor = ((*map)("fP2Pressure"));
1601 sensor =dynamic_cast<AliDCSSensor *>(osensor);
1603 if (!sensor) return 0;
1605 AliDCSSensor * sensor2 = new AliDCSSensor(*sensor);
1606 osensor = ((*map)("fCavernPressure"));
1607 TGraph * gr = new TGraph(2);
1608 gr->GetX()[0]= -100000.;
1609 gr->GetX()[1]= 1000000.;
1610 gr->GetY()[0]= atof(osensor->GetName());
1611 gr->GetY()[1]= atof(osensor->GetName());
1612 sensor2->SetGraph(gr);
1616 AliGRPObject *grpRun = new AliGRPObject;
1617 grpRun->ReadValuesFromMap(map);
1618 grpRun->SetCavernAtmosPressure(sensor2);
1619 grpRun->SetSurfaceAtmosPressure(sensor);
1623 Bool_t AliTPCcalibDB::CreateGUITree(Int_t run, const char* filename)
1626 // Create a gui tree for run number 'run'
1629 if (!AliCDBManager::Instance()->GetDefaultStorage()){
1630 AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
1631 MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
1635 AliTPCcalibDB *db=AliTPCcalibDB::Instance();
1636 // retrieve cal pad objects
1638 db->CreateGUITree(filename);
1642 Bool_t AliTPCcalibDB::CreateGUITree(const char* filename){
1646 if (!AliCDBManager::Instance()->GetDefaultStorage()){
1647 AliError("Default Storage not set. Cannot create calibration Tree!");
1650 UpdateNonRec(); // load all infromation now
1652 AliTPCPreprocessorOnline prep;
1653 //noise and pedestals
1654 if (GetPedestals()) prep.AddComponent(new AliTPCCalPad(*(GetPedestals())));
1655 if (GetPadNoise() ) prep.AddComponent(new AliTPCCalPad(*(GetPadNoise())));
1657 if (GetPulserTmean()) prep.AddComponent(new AliTPCCalPad(*(GetPulserTmean())));
1658 if (GetPulserTrms() ) prep.AddComponent(new AliTPCCalPad(*(GetPulserTrms())));
1659 if (GetPulserQmean()) prep.AddComponent(new AliTPCCalPad(*(GetPulserQmean())));
1661 if (GetCETmean()) prep.AddComponent(new AliTPCCalPad(*(GetCETmean())));
1662 if (GetCETrms() ) prep.AddComponent(new AliTPCCalPad(*(GetCETrms())));
1663 if (GetCEQmean()) prep.AddComponent(new AliTPCCalPad(*(GetCEQmean())));
1665 if (GetALTROAcqStart() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROAcqStart() )));
1666 if (GetALTROZsThr() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROZsThr() )));
1667 if (GetALTROFPED() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROFPED() )));
1668 if (GetALTROAcqStop() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROAcqStop() )));
1669 if (GetALTROMasked() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROMasked() )));
1671 AliTPCdataQA *dataQA=GetDataQA();
1673 if (dataQA->GetNLocalMaxima())
1674 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNLocalMaxima())));
1675 if (dataQA->GetMaxCharge())
1676 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetMaxCharge())));
1677 if (dataQA->GetMeanCharge())
1678 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetMeanCharge())));
1679 if (dataQA->GetNoThreshold())
1680 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNoThreshold())));
1681 if (dataQA->GetNTimeBins())
1682 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNTimeBins())));
1683 if (dataQA->GetNPads())
1684 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNPads())));
1685 if (dataQA->GetTimePosition())
1686 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetTimePosition())));
1690 TString file(filename);
1691 if (file.IsNull()) file=Form("guiTreeRun_%d.root",fRun);
1692 prep.DumpToFile(file.Data());
1696 Bool_t AliTPCcalibDB::CreateRefFile(Int_t run, const char* filename)
1699 // Create a gui tree for run number 'run'
1702 if (!AliCDBManager::Instance()->GetDefaultStorage()){
1703 AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
1704 MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
1707 TString file(filename);
1708 if (file.IsNull()) file=Form("RefCalPads_%d.root",run);
1709 TDirectory *currDir=gDirectory;
1711 AliTPCcalibDB *db=AliTPCcalibDB::Instance();
1712 // retrieve cal pad objects
1715 TFile f(file.Data(),"recreate");
1716 //noise and pedestals
1717 db->GetPedestals()->Write("Pedestals");
1718 db->GetPadNoise()->Write("PadNoise");
1720 db->GetPulserTmean()->Write("PulserTmean");
1721 db->GetPulserTrms()->Write("PulserTrms");
1722 db->GetPulserQmean()->Write("PulserQmean");
1724 db->GetCETmean()->Write("CETmean");
1725 db->GetCETrms()->Write("CETrms");
1726 db->GetCEQmean()->Write("CEQmean");
1728 db->GetALTROAcqStart() ->Write("ALTROAcqStart");
1729 db->GetALTROZsThr() ->Write("ALTROZsThr");
1730 db->GetALTROFPED() ->Write("ALTROFPED");
1731 db->GetALTROAcqStop() ->Write("ALTROAcqStop");
1732 db->GetALTROMasked() ->Write("ALTROMasked");
1741 Double_t AliTPCcalibDB::GetVDriftCorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
1743 // Get time dependent drift velocity correction
1744 // multiplication factor vd = vdnom *(1+vdriftcorr)
1746 // mode determines the algorith how to combine the Laser Track, LaserCE and physics tracks
1747 // timestamp - timestamp
1749 // side - the drift velocity per side (possible for laser and CE)
1751 // Notice - Extrapolation outside of calibration range - using constant function
1754 // mode 1 automatic mode - according to the distance to the valid calibration
1756 Double_t deltaP=0, driftP=0, wP = 0.;
1757 Double_t deltaITS=0,driftITS=0, wITS= 0.;
1758 Double_t deltaLT=0, driftLT=0, wLT = 0.;
1759 Double_t deltaCE=0, driftCE=0, wCE = 0.;
1760 driftP = fDButil->GetVDriftTPC(deltaP,run,timeStamp);
1761 driftITS= fDButil->GetVDriftTPCITS(deltaITS,run,timeStamp);
1762 driftCE = fDButil->GetVDriftTPCCE(deltaCE, run,timeStamp,36000,2);
1763 driftLT = fDButil->GetVDriftTPCLaserTracks(deltaLT,run,timeStamp,36000,2);
1764 deltaITS = TMath::Abs(deltaITS);
1765 deltaP = TMath::Abs(deltaP);
1766 deltaLT = TMath::Abs(deltaLT);
1767 deltaCE = TMath::Abs(deltaCE);
1769 const Double_t kEpsilon=0.00000000001;
1770 const Double_t kdeltaT=360.; // 10 minutes
1771 wITS = 64.*kdeltaT/(deltaITS +kdeltaT);
1772 wLT = 16.*kdeltaT/(deltaLT +kdeltaT);
1773 wP = 0. *kdeltaT/(deltaP +kdeltaT);
1774 wCE = 1. *kdeltaT/(deltaCE +kdeltaT);
1777 if (TMath::Abs(driftP)<kEpsilon) wP=0; // invalid calibration
1778 if (TMath::Abs(driftITS)<kEpsilon)wITS=0; // invalid calibration
1779 if (TMath::Abs(driftLT)<kEpsilon) wLT=0; // invalid calibration
1780 if (TMath::Abs(driftCE)<kEpsilon) wCE=0; // invalid calibration
1781 if (wP+wITS+wLT+wCE<kEpsilon) return 0;
1782 result = (driftP*wP+driftITS*wITS+driftLT*wLT+driftCE*wCE)/(wP+wITS+wLT+wCE);
1788 Double_t AliTPCcalibDB::GetTime0CorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
1790 // Get time dependent time 0 (trigger delay in cm) correction
1791 // additive correction time0 = time0+ GetTime0CorrectionTime
1792 // Value etracted combining the vdrift correction using laser tracks and CE and the physics track matchin
1794 // mode determines the algorith how to combine the Laser Track and physics tracks
1795 // timestamp - timestamp
1797 // side - the drift velocity per side (possible for laser and CE)
1799 // Notice - Extrapolation outside of calibration range - using constant function
1804 result=fDButil->GetTriggerOffsetTPC(run,timeStamp);
1805 result *=fParam->GetZLength();
1810 result= -fDButil->GetTime0TPCITS(dist, run, timeStamp)*fParam->GetDriftV()/1000000.;
1819 Double_t AliTPCcalibDB::GetVDriftCorrectionGy(Int_t timeStamp, Int_t run, Int_t side, Int_t /*mode*/){
1821 // Get global y correction drift velocity correction factor
1822 // additive factor vd = vdnom*(1+GetVDriftCorrectionGy *gy)
1823 // Value etracted combining the vdrift correction using laser tracks and CE
1825 // mode determines the algorith how to combine the Laser Track, LaserCE
1826 // timestamp - timestamp
1828 // side - the drift velocity gy correction per side (CE and Laser tracks)
1830 // Notice - Extrapolation outside of calibration range - using constant function
1832 if (run<=0 && fTransform) run = fTransform->GetCurrentRunNumber();
1833 UpdateRunInformations(run,kFALSE);
1834 TObjArray *array =AliTPCcalibDB::Instance()->GetTimeVdriftSplineRun(run);
1835 if (!array) return 0;
1836 TGraphErrors *laserA= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
1837 TGraphErrors *laserC= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
1840 if (laserA && laserC){
1841 result= (laserA->Eval(timeStamp)+laserC->Eval(timeStamp))*0.5;
1843 if (laserA && side==0){
1844 result = (laserA->Eval(timeStamp));
1846 if (laserC &&side==1){
1847 result = (laserC->Eval(timeStamp));
1849 return -result/250.; //normalized before
1852 AliTPCCalPad* AliTPCcalibDB::MakeDeadMap(Double_t notInMap, const char* nameMappingFile) {
1854 // Read list of active DDLs from OCDB entry
1855 // Generate and return AliTPCCalPad containing 1 for all pads in active DDLs,
1856 // 0 for all pads in non-active DDLs.
1857 // For DDLs with missing status information (no DCS input point to Shuttle),
1858 // the value of the AliTPCCalPad entry is determined by the parameter
1859 // notInMap (default value 1)
1863 TFile *fileMapping = new TFile(nameMappingFile, "read");
1864 AliTPCmapper *mapping = (AliTPCmapper*) fileMapping->Get("tpcMapping");
1866 sprintf(chinfo,"Failed to get mapping object from %s. ...\n", nameMappingFile);
1871 AliTPCCalPad *deadMap = new AliTPCCalPad("deadMap","deadMap");
1873 AliError("Failed to allocate dead map AliTPCCalPad");
1877 /// get list of active DDLs from OCDB entry
1879 if (!fALTROConfigData ) {
1880 AliError("No ALTRO config OCDB entry available");
1883 TMap *activeDDL = (TMap*)fALTROConfigData->FindObject("DDLArray");
1884 TObjString *ddlArray=0;
1886 ddlArray = (TObjString*)activeDDL->GetValue("DDLArray");
1888 AliError("Empty list of active DDLs in OCDB entry");
1892 AliError("List of active DDLs not available in OCDB entry");
1895 TString arrDDL=ddlArray->GetString();
1896 Int_t offset = mapping->GetTpcDdlOffset();
1898 for (Int_t i=0; i<mapping->GetNumDdl(); i++) {
1900 Int_t patch = mapping->GetPatchFromEquipmentID(idDDL);
1901 Int_t roc=mapping->GetRocFromEquipmentID(idDDL);
1902 AliTPCCalROC *calRoc=deadMap->GetCalROC(roc);
1904 for ( Int_t branch = 0; branch < 2; branch++ ) {
1905 for ( Int_t fec = 0; fec < mapping->GetNfec(patch, branch); fec++ ) {
1906 for ( Int_t altro = 0; altro < 8; altro++ ) {
1907 for ( Int_t channel = 0; channel < 16; channel++ ) {
1908 Int_t hwadd = mapping->CodeHWAddress(branch, fec, altro, channel);
1909 Int_t row = mapping->GetPadRow(patch, hwadd); // row in a ROC (IROC or OROC)
1910 // Int_t globalrow = mapping.GetGlobalPadRow(patch, hwadd); // row in full sector (IROC plus OROC)
1911 Int_t pad = mapping->GetPad(patch, hwadd);
1912 if (!TString(arrDDL[i]).IsDigit()) {
1915 active=TString(arrDDL[i]).Atof();
1917 calRoc->SetValue(row,pad,active);
1918 } // end channel for loop
1919 } // end altro for loop
1920 } // end fec for loop
1921 } // end branch for loop
1923 } // end loop on active DDLs