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 - AliTPCclusterer::Digits2Clusters(AliRawReader* rawReader)
31 // 1.) pad by pad calibration - AliTPCCalPad
34 // Simulation: AliTPCDigitizer::ExecFast - Multiply by gain
35 // Reconstruction : AliTPCclusterer::Digits2Clusters - Divide by gain
38 // Simulation: AliTPCDigitizer::ExecFast
39 // Reconstruction: AliTPCclusterer::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 AliTPCclusterer::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 // AliTPCclusterer::AddCluster
64 // AliTPCtracker::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 "AliTPCComposedCorrection.h"
129 #include "AliTPCPreprocessorOnline.h"
130 #include "AliTimeStamp.h"
131 #include "AliTriggerRunScalers.h"
132 #include "AliTriggerScalers.h"
133 #include "AliTriggerScalersRecord.h"
135 ClassImp(AliTPCcalibDB)
137 AliTPCcalibDB* AliTPCcalibDB::fgInstance = 0;
138 Bool_t AliTPCcalibDB::fgTerminated = kFALSE;
139 TObjArray AliTPCcalibDB::fgExBArray; // array of ExB corrections
142 //_ singleton implementation __________________________________________________
143 AliTPCcalibDB* AliTPCcalibDB::Instance()
146 // Singleton implementation
147 // Returns an instance of this class, it is created if necessary
150 if (fgTerminated != kFALSE)
154 fgInstance = new AliTPCcalibDB();
159 void AliTPCcalibDB::Terminate()
162 // Singleton implementation
163 // Deletes the instance of this class and sets the terminated flag, instances cannot be requested anymore
164 // This function can be called several times.
167 fgTerminated = kTRUE;
176 //_____________________________________________________________________________
177 AliTPCcalibDB::AliTPCcalibDB():
183 fActiveChannelMap(0),
187 fComposedCorrection(0),
188 fComposedCorrectionArray(0),
204 fTimeGainSplinesArray(1),
205 fGRPArray(1), //! array of GRPs - per run - JUST for calibration studies
206 fGRPMaps(1), //! array of GRPs - per run - JUST for calibration studies
207 fGoofieArray(1), //! array of GOOFIE values -per run - Just for calibration studies
209 fTemperatureArray(1), //! array of temperature sensors - per run - Just for calibration studies
210 fVdriftArray(1), //! array of v drift interfaces
211 fDriftCorrectionArray(1), //! array of drift correction
212 fRunList(1), //! run list - indicates try to get the run param
213 fBHasAlignmentOCDB(kFALSE), // Flag - has the alignment on the composed correction ?
223 for (Int_t i=0;i<72;++i){
224 fChamberHVStatus[i]=kTRUE;
225 fChamberHVmedian[i]=-1;
226 fCurrentNominalVoltage[i]=0.;
227 fChamberHVgoodFraction[i]=0.;
229 Update(); // temporary
230 fTimeGainSplinesArray.SetOwner(); //own the keys
231 fGRPArray.SetOwner(); //own the keys
232 fGRPMaps.SetOwner(); //own the keys
233 fGoofieArray.SetOwner(); //own the keys
234 fVoltageArray.SetOwner(); //own the keys
235 fTemperatureArray.SetOwner(); //own the keys
236 fVdriftArray.SetOwner(); //own the keys
237 fDriftCorrectionArray.SetOwner(); //own the keys
240 AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& ):
246 fActiveChannelMap(0),
250 fComposedCorrection(0),
251 fComposedCorrectionArray(0),
267 fTimeGainSplinesArray(1),
268 fGRPArray(0), //! array of GRPs - per run - JUST for calibration studies
269 fGRPMaps(0), //! array of GRPs - per run - JUST for calibration studies
270 fGoofieArray(0), //! array of GOOFIE values -per run - Just for calibration studies
272 fTemperatureArray(0), //! array of temperature sensors - per run - Just for calibration studies
273 fVdriftArray(0), //! array of v drift interfaces
274 fDriftCorrectionArray(0), //! array of v drift corrections
275 fRunList(0), //! run list - indicates try to get the run param
276 fBHasAlignmentOCDB(kFALSE), // Flag - has the alignment on the composed correction ?
282 // Copy constructor invalid -- singleton implementation
284 Error("copy constructor","invalid -- singleton implementation");
285 for (Int_t i=0;i<72;++i){
286 fChamberHVStatus[i]=kTRUE;
287 fChamberHVmedian[i]=-1;
288 fCurrentNominalVoltage[i]=0.;
289 fChamberHVgoodFraction[i]=0.;
291 fTimeGainSplinesArray.SetOwner(); //own the keys
292 fGRPArray.SetOwner(); //own the keys
293 fGRPMaps.SetOwner(); //own the keys
294 fGoofieArray.SetOwner(); //own the keys
295 fVoltageArray.SetOwner(); //own the keys
296 fTemperatureArray.SetOwner(); //own the keys
297 fVdriftArray.SetOwner(); //own the keys
298 fDriftCorrectionArray.SetOwner(); //own the keys
301 AliTPCcalibDB& AliTPCcalibDB::operator= (const AliTPCcalibDB& )
304 // Singleton implementation - no assignment operator
306 Error("operator =", "assignment operator not implemented");
312 //_____________________________________________________________________________
313 AliTPCcalibDB::~AliTPCcalibDB()
319 delete fActiveChannelMap;
322 AliTPCCalPad* AliTPCcalibDB::GetDistortionMap(Int_t i) const {
324 // get distortion map - due E field distortions
326 return (fDistortionMap) ? (AliTPCCalPad*)fDistortionMap->At(i):0;
329 //_____________________________________________________________________________
330 AliCDBEntry* AliTPCcalibDB::GetCDBEntry(const char* cdbPath)
333 // Retrieves an entry with path <cdbPath> from the CDB.
337 AliCDBEntry* entry = AliCDBManager::Instance()->Get(cdbPath, fRun);
340 snprintf(chinfo,1000,"AliTPCcalibDB: Failed to get entry:\t%s ", cdbPath);
348 //_____________________________________________________________________________
349 void AliTPCcalibDB::SetRun(Long64_t run)
352 // Sets current run number. Calibration data is read from the corresponding file.
362 void AliTPCcalibDB::Update(){
364 // cache the OCDB entries for simulation, reconstruction, calibration
367 AliCDBEntry * entry=0;
368 Bool_t cdbCache = AliCDBManager::Instance()->GetCacheFlag(); // save cache status
369 AliCDBManager::Instance()->SetCacheFlag(kTRUE); // activate CDB cache
370 fDButil = new AliTPCcalibDButil;
372 fRun = AliCDBManager::Instance()->GetRun();
374 entry = GetCDBEntry("TPC/Calib/PadGainFactor");
376 //if (fPadGainFactor) delete fPadGainFactor;
377 entry->SetOwner(kTRUE);
378 fPadGainFactor = (AliTPCCalPad*)entry->GetObject();
380 AliFatal("TPC - Missing calibration entry TPC/Calib/PadGainFactor");
383 entry = GetCDBEntry("TPC/Calib/TimeGain");
385 //if (fTimeGainSplines) delete fTimeGainSplines;
386 entry->SetOwner(kTRUE);
387 fTimeGainSplines = (TObjArray*)entry->GetObject();
389 AliFatal("TPC - Missing calibration entry TPC/Calib/Timegain");
392 entry = GetCDBEntry("TPC/Calib/GainFactorDedx");
394 entry->SetOwner(kTRUE);
395 fDedxGainFactor = (AliTPCCalPad*)entry->GetObject();
397 AliFatal("TPC - Missing calibration entry TPC/Calib/gainFactordEdx");
400 entry = GetCDBEntry("TPC/Calib/PadTime0");
402 //if (fPadTime0) delete fPadTime0;
403 entry->SetOwner(kTRUE);
404 fPadTime0 = (AliTPCCalPad*)entry->GetObject();
406 AliFatal("TPC - Missing calibration entry");
409 entry = GetCDBEntry("TPC/Calib/Distortion");
411 //if (fPadTime0) delete fPadTime0;
412 entry->SetOwner(kTRUE);
413 fDistortionMap =dynamic_cast<TObjArray*>(entry->GetObject());
415 //AliFatal("TPC - Missing calibration entry")
421 entry = GetCDBEntry("TPC/Calib/PadNoise");
423 //if (fPadNoise) delete fPadNoise;
424 entry->SetOwner(kTRUE);
425 fPadNoise = (AliTPCCalPad*)entry->GetObject();
427 AliFatal("TPC - Missing calibration entry");
430 entry = GetCDBEntry("TPC/Calib/Pedestals");
432 //if (fPedestals) delete fPedestals;
433 entry->SetOwner(kTRUE);
434 fPedestals = (AliTPCCalPad*)entry->GetObject();
437 entry = GetCDBEntry("TPC/Calib/Temperature");
439 //if (fTemperature) delete fTemperature;
440 entry->SetOwner(kTRUE);
441 fTemperature = (AliTPCSensorTempArray*)entry->GetObject();
444 entry = GetCDBEntry("TPC/Calib/Parameters");
446 //if (fPadNoise) delete fPadNoise;
447 entry->SetOwner(kTRUE);
448 fParam = (AliTPCParam*)(entry->GetObject());
450 AliFatal("TPC - Missing calibration entry TPC/Calib/Parameters");
453 entry = GetCDBEntry("TPC/Calib/ClusterParam");
455 entry->SetOwner(kTRUE);
456 fClusterParam = (AliTPCClusterParam*)(entry->GetObject());
458 AliFatal("TPC - Missing calibration entry");
461 //ALTRO configuration data
462 entry = GetCDBEntry("TPC/Calib/AltroConfig");
464 entry->SetOwner(kTRUE);
465 fALTROConfigData=(TObjArray*)(entry->GetObject());
467 AliFatal("TPC - Missing calibration entry");
470 //Calibration Pulser data
471 entry = GetCDBEntry("TPC/Calib/Pulser");
473 entry->SetOwner(kTRUE);
474 fPulserData=(TObjArray*)(entry->GetObject());
477 //Calibration ION tail data
478 // entry = GetCDBEntry("TPC/Calib/IonTail");
480 // entry->SetOwner(kTRUE);
481 // fIonTailArray=(TObjArray*)(entry->GetObject());
486 entry = GetCDBEntry("TPC/Calib/CE");
488 entry->SetOwner(kTRUE);
489 fCEData=(TObjArray*)(entry->GetObject());
491 //RAW calibration data
492 // entry = GetCDBEntry("TPC/Calib/Raw");
494 entry = GetCDBEntry("TPC/Calib/Mapping");
496 //if (fPadNoise) delete fPadNoise;
497 entry->SetOwner(kTRUE);
498 TObjArray * array = dynamic_cast<TObjArray*>(entry->GetObject());
499 if (array && array->GetEntriesFast()==6){
500 fMapping = new AliTPCAltroMapping*[6];
501 for (Int_t i=0; i<6; i++){
502 fMapping[i] = dynamic_cast<AliTPCAltroMapping*>(array->At(i));
507 //CTP calibration data
508 entry = GetCDBEntry("GRP/CTP/CTPtiming");
510 //entry->SetOwner(kTRUE);
511 fCTPTimeParams=dynamic_cast<AliCTPTimeParams*>(entry->GetObject());
513 AliError("TPC - Missing calibration entry");
515 //TPC space point correction data
516 entry = GetCDBEntry("TPC/Calib/Correction");
518 //entry->SetOwner(kTRUE);
519 fComposedCorrection=dynamic_cast<AliTPCCorrection*>(entry->GetObject());
520 if (fComposedCorrection) fComposedCorrection->Init();
521 fComposedCorrectionArray=dynamic_cast<TObjArray*>(entry->GetObject());
522 if (fComposedCorrectionArray){
523 for (Int_t i=0; i<fComposedCorrectionArray->GetEntries(); i++){
524 AliTPCComposedCorrection* composedCorrection= dynamic_cast<AliTPCComposedCorrection*>(fComposedCorrectionArray->At(i));
525 if (composedCorrection) {
526 composedCorrection->Init();
527 if (composedCorrection->GetCorrections()){
528 if (composedCorrection->GetCorrections()->FindObject("FitAlignTPC")){
529 fBHasAlignmentOCDB=kTRUE;
536 AliError("TPC - Missing calibration entry- TPC/Calib/Correction");
538 //RCU trigger config mode
539 fMode=GetRCUTriggerConfig();
542 fTransform=new AliTPCTransform();
543 fTransform->SetCurrentRun(AliCDBManager::Instance()->GetRun());
547 // needs to be called before InitDeadMap
548 UpdateChamberHighVoltageData();
550 // Create Dead Channel Map
554 AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
557 void AliTPCcalibDB::UpdateNonRec(){
559 // Update/Load the parameters which are important for QA studies
560 // and not used yet for the reconstruction
562 //RAW calibration data
563 AliCDBEntry * entry=0;
564 entry = GetCDBEntry("TPC/Calib/Raw");
566 entry->SetOwner(kTRUE);
567 TObjArray *arr=(TObjArray*)(entry->GetObject());
568 if (arr) fCalibRaw=(AliTPCCalibRaw*)arr->At(0);
570 //QA calibration data
571 entry = GetCDBEntry("TPC/Calib/QA");
573 entry->SetOwner(kTRUE);
574 fDataQA=dynamic_cast<AliTPCdataQA*>(entry->GetObject());
577 if (fRun>=0 && !fVoltageArray.GetValue(Form("%i",fRun))){
578 entry = AliCDBManager::Instance()->Get("TPC/Calib/HighVoltage",fRun);
580 fVoltageArray.Add(new TObjString(Form("%i",fRun)),entry->GetObject());
588 void AliTPCcalibDB::CreateObjectList(const Char_t *filename, TObjArray *calibObjects)
591 // Create calibration objects and read contents from OCDB
593 if ( calibObjects == 0x0 ) return;
596 if ( !in.is_open() ){
597 fprintf(stderr,"Error: cannot open list file '%s'", filename);
601 AliTPCCalPad *calPad=0x0;
607 TObjArray *arrFileLine = sFile.Tokenize("\n");
609 TIter nextLine(arrFileLine);
611 TObjString *sObjLine=0x0;
612 while ( (sObjLine = (TObjString*)nextLine()) ){
613 TString sLine(sObjLine->GetString());
615 TObjArray *arrNextCol = sLine.Tokenize("\t");
617 TObjString *sObjType = (TObjString*)(arrNextCol->At(0));
618 TObjString *sObjFileName = (TObjString*)(arrNextCol->At(1));
621 if ( !sObjType || ! sObjFileName ) continue;
622 TString sType(sObjType->GetString());
623 TString sFileName(sObjFileName->GetString());
624 // printf("%s\t%s\n",sType.Data(),sFileName.Data());
626 TFile *fIn = TFile::Open(sFileName);
628 fprintf(stderr,"File not found: '%s'", sFileName.Data());
632 if ( sType == "CE" ){
633 AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
635 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());
636 calPad->SetNameTitle("CETmean","CETmean");
637 calibObjects->Add(calPad);
639 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());
640 calPad->SetNameTitle("CEQmean","CEQmean");
641 calibObjects->Add(calPad);
643 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
644 calPad->SetNameTitle("CETrms","CETrms");
645 calibObjects->Add(calPad);
647 } else if ( sType == "Pulser") {
648 AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
650 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());
651 calPad->SetNameTitle("PulserTmean","PulserTmean");
652 calibObjects->Add(calPad);
654 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());
655 calPad->SetNameTitle("PulserQmean","PulserQmean");
656 calibObjects->Add(calPad);
658 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
659 calPad->SetNameTitle("PulserTrms","PulserTrms");
660 calibObjects->Add(calPad);
662 } else if ( sType == "Pedestals") {
663 AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
665 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());
666 calPad->SetNameTitle("Pedestals","Pedestals");
667 calibObjects->Add(calPad);
669 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());
670 calPad->SetNameTitle("Noise","Noise");
671 calibObjects->Add(calPad);
674 fprintf(stderr,"Undefined Type: '%s'",sType.Data());
682 Int_t AliTPCcalibDB::InitDeadMap() {
683 // Initialize DeadChannel Map
684 // Source of information:
685 // - HV (see UpdateChamberHighVoltageData())
686 // - Altro disabled channels. Noisy channels.
689 // check necessary information
690 const Int_t run=GetRun();
692 AliError("run not set in CDB manager. Cannot create active channel map");
695 AliDCSSensorArray* voltageArray = GetVoltageSensors(run);
696 AliTPCCalPad* altroMap = GetALTROMasked();
697 TMap* mapddl = GetDDLMap();
699 if (!voltageArray && !altroMap && !mapddl) {
700 AliError("All necessary information to create the activate channel are map missing.");
704 //=============================================================
707 Bool_t ddlMap[216]={0};
708 for (Int_t iddl=0; iddl<216; ++iddl) ddlMap[iddl]=1;
710 TObjString *s = (TObjString*)mapddl->GetValue("DDLArray");
712 for (Int_t iddl=0; iddl<216; ++iddl) ddlMap[iddl]=TString(s->GetString()(iddl))!="0";
715 AliError("DDL map missing. ActiveChannelMap can only be created with parts of the information.");
717 // Setup DDL map done
718 // ============================================================
720 //=============================================================
721 // Setup active chnnel map
724 if (!fActiveChannelMap) fActiveChannelMap=new AliTPCCalPad("ActiveChannelMap","ActiveChannelMap");
726 AliTPCmapper map(gSystem->ExpandPathName("$ALICE_ROOT/TPC/mapping/"));
728 if (!altroMap) AliError("ALTRO dead channel map missing. ActiveChannelMap can only be created with parts of the information.");
730 for (Int_t iROC=0;iROC<AliTPCCalPad::kNsec;++iROC){
731 AliTPCCalROC *roc=fActiveChannelMap->GetCalROC(iROC);
733 AliError(Form("No ROC %d in active channel map",iROC));
737 // check for bad voltage
738 // see UpdateChamberHighVoltageData()
739 if (!fChamberHVStatus[iROC]){
744 AliTPCCalROC *masked=0x0;
745 if (altroMap) masked=altroMap->GetCalROC(iROC);
747 for (UInt_t irow=0; irow<roc->GetNrows(); ++irow){
748 for (UInt_t ipad=0; ipad<roc->GetNPads(irow); ++ipad){
749 //per default the channel is on
750 roc->SetValue(irow,ipad,1);
751 // apply altro dead channel mask (inverse logik, it is not active, but inactive channles)
752 if (masked && masked->GetValue(irow, ipad)) roc->SetValue(irow, ipad ,0);
753 // mask channels if a DDL is inactive
754 Int_t ddlId=map.GetEquipmentID(iROC, irow, ipad)-768;
755 if (ddlId>=0 && !ddlMap[ddlId]) roc->SetValue(irow, ipad ,0);
763 void AliTPCcalibDB::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
765 // Write a tree with all available information
766 // if mapFileName is specified, the Map information are also written to the tree
767 // pads specified in outlierPad are not used for calculating statistics
768 // - the same function as AliTPCCalPad::MakeTree -
770 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
772 TObjArray* mapIROCs = 0;
773 TObjArray* mapOROCs = 0;
774 TVectorF *mapIROCArray = 0;
775 TVectorF *mapOROCArray = 0;
776 Int_t mapEntries = 0;
777 TString* mapNames = 0;
780 TFile mapFile(mapFileName, "read");
782 TList* listOfROCs = mapFile.GetListOfKeys();
783 mapEntries = listOfROCs->GetEntries()/2;
784 mapIROCs = new TObjArray(mapEntries*2);
785 mapOROCs = new TObjArray(mapEntries*2);
786 mapIROCArray = new TVectorF[mapEntries];
787 mapOROCArray = new TVectorF[mapEntries];
789 mapNames = new TString[mapEntries];
790 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
791 TString nameROC(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
792 nameROC.Remove(nameROC.Length()-4, 4);
793 mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "IROC").Data()), ivalue);
794 mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "OROC").Data()), ivalue);
795 mapNames[ivalue].Append(nameROC);
798 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
799 mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
800 mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
802 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
803 (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
804 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
805 (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
808 } // if (mapFileName)
810 TTreeSRedirector cstream(fileName);
811 Int_t arrayEntries = array->GetEntries();
813 TString* names = new TString[arrayEntries];
814 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
815 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
817 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
819 // get statistic for given sector
821 TVectorF median(arrayEntries);
822 TVectorF mean(arrayEntries);
823 TVectorF rms(arrayEntries);
824 TVectorF ltm(arrayEntries);
825 TVectorF ltmrms(arrayEntries);
826 TVectorF medianWithOut(arrayEntries);
827 TVectorF meanWithOut(arrayEntries);
828 TVectorF rmsWithOut(arrayEntries);
829 TVectorF ltmWithOut(arrayEntries);
830 TVectorF ltmrmsWithOut(arrayEntries);
832 TVectorF *vectorArray = new TVectorF[arrayEntries];
833 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
834 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
836 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
837 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
838 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
839 AliTPCCalROC* outlierROC = 0;
840 if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
842 median[ivalue] = calROC->GetMedian();
843 mean[ivalue] = calROC->GetMean();
844 rms[ivalue] = calROC->GetRMS();
845 Double_t ltmrmsValue = 0;
846 ltm[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction);
847 ltmrms[ivalue] = ltmrmsValue;
849 medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
850 meanWithOut[ivalue] = calROC->GetMean(outlierROC);
851 rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
853 ltmWithOut[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction, outlierROC);
854 ltmrmsWithOut[ivalue] = ltmrmsValue;
863 medianWithOut[ivalue] = 0.;
864 meanWithOut[ivalue] = 0.;
865 rmsWithOut[ivalue] = 0.;
866 ltmWithOut[ivalue] = 0.;
867 ltmrmsWithOut[ivalue] = 0.;
872 // fill vectors of variable per pad
874 TVectorF *posArray = new TVectorF[8];
875 for (Int_t ivalue = 0; ivalue < 8; ivalue++)
876 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
878 Float_t posG[3] = {0};
879 Float_t posL[3] = {0};
881 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
882 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
883 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
884 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
885 posArray[0][ichannel] = irow;
886 posArray[1][ichannel] = ipad;
887 posArray[2][ichannel] = posL[0];
888 posArray[3][ichannel] = posL[1];
889 posArray[4][ichannel] = posG[0];
890 posArray[5][ichannel] = posG[1];
891 posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
892 posArray[7][ichannel] = ichannel;
894 // loop over array containing AliTPCCalPads
895 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
896 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
897 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
899 (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
901 (vectorArray[ivalue])[ichannel] = 0;
907 cstream << "calPads" <<
908 "sector=" << isector;
910 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
911 cstream << "calPads" <<
912 (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
913 (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
914 (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
915 (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
916 (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
918 cstream << "calPads" <<
919 (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
920 (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
921 (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
922 (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
923 (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
927 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
928 cstream << "calPads" <<
929 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
933 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
935 cstream << "calPads" <<
936 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
938 cstream << "calPads" <<
939 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
943 cstream << "calPads" <<
944 "row.=" << &posArray[0] <<
945 "pad.=" << &posArray[1] <<
946 "lx.=" << &posArray[2] <<
947 "ly.=" << &posArray[3] <<
948 "gx.=" << &posArray[4] <<
949 "gy.=" << &posArray[5] <<
950 "rpad.=" << &posArray[6] <<
951 "channel.=" << &posArray[7];
953 cstream << "calPads" <<
957 delete[] vectorArray;
965 delete[] mapIROCArray;
966 delete[] mapOROCArray;
971 Int_t AliTPCcalibDB::GetRCUTriggerConfig() const
974 // return the RCU trigger configuration register
976 TMap *map=GetRCUconfig();
978 TVectorF *v=(TVectorF*)map->GetValue("TRGCONF_TRG_MODE");
980 for (Int_t i=0; i<v->GetNrows(); ++i){
981 Float_t newmode=v->GetMatrixArray()[i];
983 if (mode>-1&&newmode!=mode) AliWarning("Found different RCU trigger configurations!!!");
990 Bool_t AliTPCcalibDB::IsTrgL0()
993 // return if the FEE readout was triggered on L0
995 if (fMode<0) return kFALSE;
999 Bool_t AliTPCcalibDB::IsTrgL1()
1002 // return if the FEE readout was triggered on L1
1004 if (fMode<0) return kFALSE;
1008 void AliTPCcalibDB::RegisterExB(Int_t index, Float_t bz, Bool_t bdelete){
1010 // Register static ExB correction map
1011 // index - registration index - used for visualization
1012 // bz - bz field in kGaus
1014 // Float_t factor = bz/(-5.); // default b filed in Cheb with minus sign
1015 Float_t factor = bz/(5.); // default b filed in Cheb with minus sign
1016 // was chenged in the Revision ???? (Ruben can you add here number)
1018 AliMagF* bmap = new AliMagF("MapsExB","MapsExB", factor,TMath::Sign(1.f,factor),AliMagF::k5kG);
1020 AliTPCExBFirst *exb = new AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
1021 AliTPCExB::SetInstance(exb);
1026 AliTPCExB::RegisterField(index,bmap);
1028 if (index>=fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
1029 fgExBArray.AddAt(exb,index);
1033 AliTPCExB* AliTPCcalibDB::GetExB(Float_t bz, Bool_t deleteB) {
1035 // bz filed in KGaus not in tesla
1036 // Get ExB correction map
1037 // if doesn't exist - create it
1039 Int_t index = TMath::Nint(5+bz);
1040 if (index>fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
1041 if (!fgExBArray.At(index)) AliTPCcalibDB::RegisterExB(index,bz,deleteB);
1042 return (AliTPCExB*)fgExBArray.At(index);
1046 void AliTPCcalibDB::SetExBField(Float_t bz){
1048 // Set magnetic filed for ExB correction
1050 fExB = GetExB(bz,kFALSE);
1053 void AliTPCcalibDB::SetExBField(const AliMagF* bmap){
1055 // Set magnetic field for ExB correction
1057 AliTPCExBFirst *exb = new AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
1058 AliTPCExB::SetInstance(exb);
1064 void AliTPCcalibDB::UpdateRunInformations( Int_t run, Bool_t force){
1066 // - > Don't use it for reconstruction - Only for Calibration studies
1069 TObjString runstr(Form("%i",run));
1071 AliCDBEntry * entry = 0;
1072 if (run>= fRunList.fN){
1073 fRunList.Set(run*2+1);
1076 fALTROConfigData->Expand(run*2+1); // ALTRO configuration data
1077 fPulserData->Expand(run*2+1); // Calibration Pulser data
1078 fCEData->Expand(run*2+1); // CE data
1079 if (!fTimeGainSplines) fTimeGainSplines = new TObjArray(run*2+1);
1080 fTimeGainSplines->Expand(run*2+1); // Array of AliSplineFits: at 0 MIP position in
1082 if (fRunList[run]>0 &&force==kFALSE) return;
1084 fRunList[run]=1; // sign as used
1087 entry = AliCDBManager::Instance()->Get("GRP/GRP/Data",run);
1089 AliGRPObject * grpRun = dynamic_cast<AliGRPObject*>(entry->GetObject());
1091 TMap* map = dynamic_cast<TMap*>(entry->GetObject());
1093 //grpRun = new AliGRPObject;
1094 //grpRun->ReadValuesFromMap(map);
1095 grpRun = MakeGRPObjectFromMap(map);
1097 fGRPMaps.Add(new TObjString(runstr),map);
1100 fGRPArray.Add(new TObjString(runstr),grpRun);
1102 entry = AliCDBManager::Instance()->Get("TPC/Calib/Goofie",run);
1104 fGoofieArray.Add(new TObjString(runstr),entry->GetObject());
1109 entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeGain",run);
1111 fTimeGainSplinesArray.Add(new TObjString(runstr),entry->GetObject());
1113 AliFatal("TPC - Missing calibration entry TimeGain");
1116 entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeDrift",run);
1118 TObjArray * timeArray = (TObjArray*)entry->GetObject();
1119 fDriftCorrectionArray.Add(new TObjString(runstr),entry->GetObject());
1120 AliTPCCorrection * correctionTime = (AliTPCCorrection *)timeArray->FindObject("FitCorrectionTime");
1121 if (correctionTime && fComposedCorrectionArray){
1122 correctionTime->Init();
1123 if (fComposedCorrectionArray->GetEntriesFast()<4) fComposedCorrectionArray->Expand(40);
1124 fComposedCorrectionArray->AddAt(correctionTime,4); //add time dependent correction to the list of available corrections
1127 AliFatal("TPC - Missing calibration entry TimeDrift");
1130 entry = AliCDBManager::Instance()->Get("TPC/Calib/Temperature",run);
1132 fTemperatureArray.Add(new TObjString(runstr),entry->GetObject());
1136 entry = AliCDBManager::Instance()->Get("TPC/Calib/HighVoltage",run);
1137 if (!fVoltageArray.GetValue(runstr.GetName()) && entry) {
1138 fVoltageArray.Add(new TObjString(runstr),entry->GetObject());
1141 //apply fDButil filters
1143 fDButil->UpdateFromCalibDB();
1144 if (fTemperature) fDButil->FilterTemperature(fTemperature);
1146 AliDCSSensor * press = GetPressureSensor(run,0);
1147 AliTPCSensorTempArray * temp = GetTemperatureSensor(run);
1148 Bool_t accept=kTRUE;
1150 accept = fDButil->FilterTemperature(temp)>0.1;
1153 const Double_t kMinP=900.;
1154 const Double_t kMaxP=1050.;
1155 const Double_t kMaxdP=10.;
1156 const Double_t kSigmaCut=4.;
1157 fDButil->FilterSensor(press,kMinP,kMaxP,kMaxdP,kSigmaCut);
1158 if (press->GetFit()==0) accept=kFALSE;
1161 if (press && temp &&accept){
1162 AliTPCCalibVdrift * vdrift = new AliTPCCalibVdrift(temp, press,0);
1163 fVdriftArray.Add(new TObjString(runstr),vdrift);
1166 fDButil->FilterCE(120., 3., 4.,0);
1167 fDButil->FilterTracks(run, 10.,0);
1172 Float_t AliTPCcalibDB::GetGain(Int_t sector, Int_t row, Int_t pad){
1174 // Get Gain factor for given pad
1176 AliTPCCalPad *calPad = Instance()->fDedxGainFactor;;
1177 if (!calPad) return 0;
1178 return calPad->GetCalROC(sector)->GetValue(row,pad);
1181 AliSplineFit* AliTPCcalibDB::GetVdriftSplineFit(const char* name, Int_t run){
1183 // GetDrift velocity spline fit
1185 TObjArray *arr=GetTimeVdriftSplineRun(run);
1187 return dynamic_cast<AliSplineFit*>(arr->FindObject(name));
1190 AliSplineFit* AliTPCcalibDB::CreateVdriftSplineFit(const char* graphName, Int_t run){
1192 // create spline fit from the drift time graph in TimeDrift
1194 TObjArray *arr=GetTimeVdriftSplineRun(run);
1196 TGraph *graph=dynamic_cast<TGraph*>(arr->FindObject(graphName));
1197 if (!graph) return 0;
1198 AliSplineFit *fit = new AliSplineFit();
1199 fit->SetGraph(graph);
1200 fit->SetMinPoints(graph->GetN()+1);
1201 fit->InitKnots(graph,2,0,0.001);
1206 AliGRPObject *AliTPCcalibDB::GetGRP(Int_t run){
1208 // Get GRP object for given run
1210 AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>((Instance()->fGRPArray).GetValue(Form("%i",run)));
1212 Instance()->UpdateRunInformations(run);
1213 grpRun = dynamic_cast<AliGRPObject *>(Instance()->fGRPArray.GetValue(Form("%i",run)));
1214 if (!grpRun) return 0;
1219 TMap * AliTPCcalibDB::GetGRPMap(Int_t run){
1221 // Get GRP map for given run
1223 TMap * grpRun = dynamic_cast<TMap *>((Instance()->fGRPMaps).GetValue(Form("%i",run)));
1225 Instance()->UpdateRunInformations(run);
1226 grpRun = dynamic_cast<TMap *>(Instance()->fGRPMaps.GetValue(Form("%i",run)));
1227 if (!grpRun) return 0;
1233 AliDCSSensor * AliTPCcalibDB::GetPressureSensor(Int_t run, Int_t type){
1235 // Get Pressure sensor
1237 // type = 0 - Cavern pressure
1238 // 1 - Suface pressure
1239 // First try to get if trom map - if existing (Old format of data storing)
1243 TMap *map = GetGRPMap(run);
1245 AliDCSSensor * sensor = 0;
1247 if (type==0) osensor = ((*map)("fCavernPressure"));
1248 if (type==1) osensor = ((*map)("fP2Pressure"));
1249 sensor =dynamic_cast<AliDCSSensor *>(osensor);
1250 if (sensor) return sensor;
1253 // If not map try to get it from the GRPObject
1255 AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.GetValue(Form("%i",run)));
1257 UpdateRunInformations(run);
1258 grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.GetValue(Form("%i",run)));
1259 if (!grpRun) return 0;
1261 AliDCSSensor * sensor = grpRun->GetCavernAtmosPressure();
1262 if (type==1) sensor = grpRun->GetSurfaceAtmosPressure();
1266 AliTPCSensorTempArray * AliTPCcalibDB::GetTemperatureSensor(Int_t run){
1268 // Get temperature sensor array
1270 AliTPCSensorTempArray * tempArray = (AliTPCSensorTempArray *)fTemperatureArray.GetValue(Form("%i",run));
1272 UpdateRunInformations(run);
1273 tempArray = (AliTPCSensorTempArray *)fTemperatureArray.GetValue(Form("%i",run));
1279 TObjArray * AliTPCcalibDB::GetTimeGainSplinesRun(Int_t run){
1281 // Get temperature sensor array
1283 TObjArray * gainSplines = (TObjArray *)fTimeGainSplinesArray.GetValue(Form("%i",run));
1285 UpdateRunInformations(run);
1286 gainSplines = (TObjArray *)fTimeGainSplinesArray.GetValue(Form("%i",run));
1291 TObjArray * AliTPCcalibDB::GetTimeVdriftSplineRun(Int_t run){
1293 // Get drift spline array
1295 TObjArray * driftSplines = (TObjArray *)fDriftCorrectionArray.GetValue(Form("%i",run));
1296 if (!driftSplines) {
1297 UpdateRunInformations(run);
1298 driftSplines = (TObjArray *)fDriftCorrectionArray.GetValue(Form("%i",run));
1300 return driftSplines;
1303 AliDCSSensorArray * AliTPCcalibDB::GetVoltageSensors(Int_t run){
1305 // Get temperature sensor array
1307 AliDCSSensorArray * voltageArray = (AliDCSSensorArray *)fVoltageArray.GetValue(Form("%i",run));
1308 if (!voltageArray) {
1309 UpdateRunInformations(run);
1310 voltageArray = (AliDCSSensorArray *)fVoltageArray.GetValue(Form("%i",run));
1312 return voltageArray;
1315 AliDCSSensorArray * AliTPCcalibDB::GetGoofieSensors(Int_t run){
1317 // Get temperature sensor array
1319 AliDCSSensorArray * goofieArray = (AliDCSSensorArray *)fGoofieArray.GetValue(Form("%i",run));
1321 UpdateRunInformations(run);
1322 goofieArray = (AliDCSSensorArray *)fGoofieArray.GetValue(Form("%i",run));
1329 AliTPCCalibVdrift * AliTPCcalibDB::GetVdrift(Int_t run){
1331 // Get the interface to the the vdrift
1333 AliTPCCalibVdrift * vdrift = (AliTPCCalibVdrift*)fVdriftArray.GetValue(Form("%i",run));
1335 UpdateRunInformations(run);
1336 vdrift= (AliTPCCalibVdrift*)fVdriftArray.GetValue(Form("%i",run));
1341 Float_t AliTPCcalibDB::GetCEdriftTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1344 // GetCE drift time information for 'sector'
1345 // sector 72 is the mean drift time of the A-Side
1346 // sector 73 is the mean drift time of the C-Side
1347 // it timestamp==-1 return mean value
1349 AliTPCcalibDB::Instance()->SetRun(run);
1350 TGraph *gr=AliTPCcalibDB::Instance()->GetCErocTgraph(sector);
1351 if (!gr||sector<0||sector>73) {
1352 if (entries) *entries=0;
1356 if (timeStamp==-1.){
1359 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1361 gr->GetPoint(ipoint,x,y);
1362 if (x<timeStamp) continue;
1370 Float_t AliTPCcalibDB::GetCEchargeTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1373 // GetCE mean charge for 'sector'
1374 // it timestamp==-1 return mean value
1376 AliTPCcalibDB::Instance()->SetRun(run);
1377 TGraph *gr=AliTPCcalibDB::Instance()->GetCErocQgraph(sector);
1378 if (!gr||sector<0||sector>71) {
1379 if (entries) *entries=0;
1383 if (timeStamp==-1.){
1386 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1388 gr->GetPoint(ipoint,x,y);
1389 if (x<timeStamp) continue;
1397 Float_t AliTPCcalibDB::GetDCSSensorValue(AliDCSSensorArray *arr, Int_t timeStamp, const char * sensorName, Int_t sigDigits)
1400 // Get Value for a DCS sensor 'sensorName', run 'run' at time 'timeStamp'
1403 const TString sensorNameString(sensorName);
1404 AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1405 if (!sensor) return val;
1406 //use the dcs graph if possible
1407 TGraph *gr=sensor->GetGraph();
1409 for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1411 gr->GetPoint(ipoint,x,y);
1412 Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1413 if (time<timeStamp) continue;
1417 //if val is still 0, test if if the requested time if within 5min of the first/last
1418 //data point. If this is the case return the firs/last entry
1419 //the timestamps might not be syncronised for all calibration types, sometimes a 'pre'
1420 //and 'pos' period is requested. Especially to the HV this is not the case!
1424 gr->GetPoint(0,x,y);
1425 const Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1426 const Int_t dtime=time-timeStamp;
1427 if ( (dtime>0) && (dtime<5*60) ) val=y;
1432 gr->GetPoint(gr->GetN()-1,x,y);
1433 const Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1434 const Int_t dtime=timeStamp-time;
1435 if ( (dtime>0) && (dtime<5*60) ) val=y;
1438 val=sensor->GetValue(timeStamp);
1441 val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1446 Float_t AliTPCcalibDB::GetDCSSensorMeanValue(AliDCSSensorArray *arr, const char * sensorName, Int_t sigDigits)
1449 // Get mean Value for a DCS sensor 'sensorName' during run 'run'
1452 const TString sensorNameString(sensorName);
1453 AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1454 if (!sensor) return val;
1456 //use dcs graph if it exists
1457 TGraph *gr=sensor->GetGraph();
1461 //if we don't have the dcs graph, try to get some meaningful information
1462 if (!sensor->GetFit()) return val;
1463 Int_t nKnots=sensor->GetFit()->GetKnots();
1464 Double_t tMid=(sensor->GetEndTime()-sensor->GetStartTime())/2.;
1465 for (Int_t iKnot=0;iKnot<nKnots;++iKnot){
1466 if (sensor->GetFit()->GetX()[iKnot]>tMid/3600.) break;
1467 val=(Float_t)sensor->GetFit()->GetY0()[iKnot];
1472 val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1478 Bool_t AliTPCcalibDB::IsDataTakingActive(time_t timeStamp)
1480 if (!fGrRunState) return kFALSE;
1481 Double_t time=Double_t(timeStamp);
1482 Int_t currentPoint=0;
1483 Bool_t currentVal=fGrRunState->GetY()[currentPoint]>0.5;
1484 Bool_t retVal=currentVal;
1485 Double_t currentTime=fGrRunState->GetX()[currentPoint];
1487 while (time>currentTime){
1489 if (currentPoint==fGrRunState->GetN()) break;
1490 currentVal=fGrRunState->GetY()[currentPoint]>0.5;
1491 currentTime=fGrRunState->GetX()[currentPoint];
1498 void AliTPCcalibDB::UpdateChamberHighVoltageData()
1501 // set chamber high voltage data
1502 // 1. Robust median (sampling the hv graphs over time)
1503 // 2. Current nominal voltages (nominal voltage corrected for common HV offset)
1504 // 3. Fraction of good HV values over time (deviation from robust median)
1505 // 4. HV status, based on the above
1508 // start and end time of the run
1509 const Int_t run=GetRun();
1512 // if no valid run information - return
1513 AliGRPObject* grp = GetGRP(run);
1516 const Int_t startTimeGRP = grp->GetTimeStart();
1517 const Int_t stopTimeGRP = grp->GetTimeEnd();
1520 // check active state by analysing the scalers
1522 // initialise graph with active running
1523 AliCDBEntry *entry = GetCDBEntry("GRP/CTP/Scalers");
1525 // entry->SetOwner(kTRUE);
1526 AliTriggerRunScalers *sca = (AliTriggerRunScalers*)entry->GetObject();
1527 Int_t nchannels = sca->GetNumClasses(); // number of scaler channels (i.e. trigger classes)
1528 Int_t npoints = sca->GetScalersRecords()->GetEntries(); // number of samples
1531 fGrRunState=new TGraph;
1532 fGrRunState->SetPoint(fGrRunState->GetN(),Double_t(startTimeGRP)-.001,0);
1533 fGrRunState->SetPoint(fGrRunState->GetN(),Double_t(startTimeGRP),1);
1534 ULong64_t lastSum=0;
1535 Double_t timeLast=0.;
1536 Bool_t active=kTRUE;
1537 for (int i=0; i<npoints; i++) {
1538 AliTriggerScalersRecord *rec = (AliTriggerScalersRecord *) sca->GetScalersRecord(i);
1539 Double_t time = ((AliTimeStamp*) rec->GetTimeStamp())->GetSeconds();
1541 for (int j=0; j<nchannels; j++) sum += ((AliTriggerScalers*) rec->GetTriggerScalers()->At(j))->GetL2CA();
1542 if (TMath::Abs(time-timeLast)<.001 && sum==lastSum ) continue;
1543 if (active && sum==lastSum){
1544 fGrRunState->SetPoint(fGrRunState->GetN(),timeLast-.01,1);
1545 fGrRunState->SetPoint(fGrRunState->GetN(),timeLast,0);
1547 } else if (!active && sum>lastSum ){
1548 fGrRunState->SetPoint(fGrRunState->GetN(),timeLast-.01,0);
1549 fGrRunState->SetPoint(fGrRunState->GetN(),timeLast,1);
1555 fGrRunState->SetPoint(fGrRunState->GetN(),Double_t(stopTimeGRP),active);
1556 fGrRunState->SetPoint(fGrRunState->GetN(),Double_t(stopTimeGRP)+.001,0);
1561 for (Int_t iROC=0;iROC<72;++iROC) {
1562 fChamberHVmedian[iROC] = -1;
1563 fChamberHVgoodFraction[iROC] = 0.;
1564 fCurrentNominalVoltage[iROC] = -999.;
1565 fChamberHVStatus[iROC] = kFALSE;
1568 AliDCSSensorArray* voltageArray = GetVoltageSensors(run);
1569 if (!voltageArray) {
1570 AliError("Voltage Array missing. Cannot calculate HV information!");
1574 // max HV diffs before a chamber is masked
1575 const Float_t maxVdiff = fParam->GetMaxVoltageDeviation();
1576 const Float_t maxDipVoltage = fParam->GetMaxDipVoltage();
1577 const Float_t maxFracHVbad = fParam->GetMaxFractionHVbad();
1579 const Int_t samplingPeriod=1;
1581 // array with sampled voltages
1582 const Int_t maxSamples=(stopTimeGRP-startTimeGRP)/samplingPeriod + 10*samplingPeriod;
1583 Float_t *vSampled = new Float_t[maxSamples];
1585 // deviation of the median from the nominal voltage
1586 Double_t chamberMedianDeviation[72]={0.};
1588 for (Int_t iROC=0; iROC<72; ++iROC){
1589 chamberMedianDeviation[iROC]=0.;
1590 TString sensorName="";
1591 Char_t sideName='A';
1592 if ((iROC/18)%2==1) sideName='C';
1593 if (iROC<36) sensorName=Form("TPC_ANODE_I_%c%02d_VMEAS",sideName,iROC%18);
1594 else sensorName=Form("TPC_ANODE_O_%c%02d_0_VMEAS",sideName,iROC%18);
1596 AliDCSSensor *sensor = voltageArray->GetSensor(sensorName);
1598 fHVsensors[iROC]=sensor;
1599 if (!sensor) continue;
1601 Int_t nPointsSampled=0;
1603 TGraph *gr=sensor->GetGraph();
1604 if ( gr && gr->GetN()>1 ){
1605 //1. sample voltage over time
1606 // get a robust median
1607 // buffer sampled voltages
1609 // current sampling time
1610 Int_t time=startTimeGRP;
1612 // input graph sampling point
1613 const Int_t nGraph=gr->GetN();
1616 //initialise graph information
1617 Int_t timeGraph=TMath::Nint(gr->GetX()[pointGraph+1]*3600+sensor->GetStartTime());
1618 Double_t sampledHV=gr->GetY()[pointGraph++];
1620 while (time<stopTimeGRP){
1621 while (timeGraph<=time && pointGraph+1<nGraph){
1622 timeGraph=TMath::Nint(gr->GetX()[pointGraph+1]*3600+sensor->GetStartTime());
1623 sampledHV=gr->GetY()[pointGraph++];
1625 time+=samplingPeriod;
1626 if (!IsDataTakingActive(time-samplingPeriod)) continue;
1627 vSampled[nPointsSampled++]=sampledHV;
1630 if (nPointsSampled<1) continue;
1632 fChamberHVmedian[iROC]=TMath::Median(nPointsSampled,vSampled);
1633 chamberMedianDeviation[iROC]=fChamberHVmedian[iROC]-fParam->GetNominalVoltage(iROC);
1635 //2. calculate good HV fraction
1637 for (Int_t ipoint=0; ipoint<nPointsSampled; ++ipoint) {
1638 if (TMath::Abs(vSampled[ipoint]-fChamberHVmedian[iROC])<maxDipVoltage) ++ngood;
1641 fChamberHVgoodFraction[iROC]=Float_t(ngood)/Float_t(nPointsSampled);
1643 AliError(Form("No Graph or too few points found for HV sensor of ROC %d",iROC));
1650 // get median deviation from all chambers (detect e.g. -50V)
1651 const Double_t medianIROC=TMath::Median( 36, chamberMedianDeviation );
1652 const Double_t medianOROC=TMath::Median( 36, chamberMedianDeviation+36 );
1654 // Define current default voltages
1655 for (Int_t iROC=0;iROC<72/*AliTPCCalPad::kNsec*/;++iROC){
1656 const Float_t averageDeviation=(iROC<36)?medianIROC:medianOROC;
1657 fCurrentNominalVoltage[iROC]=fParam->GetNominalVoltage(iROC)+averageDeviation;
1663 for (Int_t iROC=0;iROC<72/*AliTPCCalPad::kNsec*/;++iROC){
1664 fChamberHVStatus[iROC]=kTRUE;
1666 //a. Deviation of median from current nominal voltage
1667 // allow larger than nominal voltages
1668 if (fCurrentNominalVoltage[iROC]-fChamberHVmedian[iROC] > maxVdiff) fChamberHVStatus[iROC]=kFALSE;
1670 //b. Fraction of bad hv values
1671 if ( 1-fChamberHVgoodFraction[iROC] > maxFracHVbad ) fChamberHVStatus[iROC]=kFALSE;
1675 Float_t AliTPCcalibDB::GetChamberHighVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits, Bool_t current) {
1677 // return the chamber HV for given run and time: 0-35 IROC, 36-72 OROC
1678 // if timeStamp==-1 return mean value
1681 TString sensorName="";
1682 TTimeStamp stamp(timeStamp);
1683 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1684 if (!voltageArray || (sector<0) || (sector>71)) return val;
1685 Char_t sideName='A';
1686 if ((sector/18)%2==1) sideName='C';
1689 sensorName=Form("TPC_ANODE_I_%c%02d_VMEAS",sideName,sector%18);
1692 sensorName=Form("TPC_ANODE_O_%c%02d_0_VMEAS",sideName,sector%18);
1697 sensorName=Form("TPC_ANODE_I_%c%02d_IMEAS",sideName,sector%18);
1700 sensorName=Form("TPC_ANODE_O_%c%02d_0_IMEAS",sideName,sector%18);
1705 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1707 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1711 Float_t AliTPCcalibDB::GetSkirtVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1714 // Get the skirt voltage for 'run' at 'timeStamp' and 'sector': 0-35 IROC, 36-72 OROC
1715 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1716 // if timeStamp==-1 return the mean value for the run
1719 TString sensorName="";
1720 TTimeStamp stamp(timeStamp);
1721 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1722 if (!voltageArray || (sector<0) || (sector>71)) return val;
1723 Char_t sideName='A';
1724 if ((sector/18)%2==1) sideName='C';
1725 sensorName=Form("TPC_SKIRT_%c_VMEAS",sideName);
1727 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1729 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1734 Float_t AliTPCcalibDB::GetCoverVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1737 // Get the cover voltage for run 'run' at time 'timeStamp'
1738 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1739 // if timeStamp==-1 return the mean value for the run
1742 TString sensorName="";
1743 TTimeStamp stamp(timeStamp);
1744 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1745 if (!voltageArray || (sector<0) || (sector>71)) return val;
1746 Char_t sideName='A';
1747 if ((sector/18)%2==1) sideName='C';
1750 sensorName=Form("TPC_COVER_I_%c_VMEAS",sideName);
1753 sensorName=Form("TPC_COVER_O_%c_VMEAS",sideName);
1756 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1758 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1763 Float_t AliTPCcalibDB::GetGGoffsetVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1766 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1767 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1768 // if timeStamp==-1 return the mean value for the run
1771 TString sensorName="";
1772 TTimeStamp stamp(timeStamp);
1773 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1774 if (!voltageArray || (sector<0) || (sector>71)) return val;
1775 Char_t sideName='A';
1776 if ((sector/18)%2==1) sideName='C';
1779 sensorName=Form("TPC_GATE_I_%c_OFF_VMEAS",sideName);
1782 sensorName=Form("TPC_GATE_O_%c_OFF_VMEAS",sideName);
1785 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1787 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1792 Float_t AliTPCcalibDB::GetGGnegVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1795 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1796 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1797 // if timeStamp==-1 return the mean value for the run
1800 TString sensorName="";
1801 TTimeStamp stamp(timeStamp);
1802 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1803 if (!voltageArray || (sector<0) || (sector>71)) return val;
1804 Char_t sideName='A';
1805 if ((sector/18)%2==1) sideName='C';
1808 sensorName=Form("TPC_GATE_I_%c_NEG_VMEAS",sideName);
1811 sensorName=Form("TPC_GATE_O_%c_NEG_VMEAS",sideName);
1814 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1816 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1821 Float_t AliTPCcalibDB::GetGGposVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1824 // Get the GG offset voltage for run 'run' at time 'timeStamp'
1825 // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1826 // if timeStamp==-1 return the mean value for the run
1829 TString sensorName="";
1830 TTimeStamp stamp(timeStamp);
1831 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1832 if (!voltageArray || (sector<0) || (sector>71)) return val;
1833 Char_t sideName='A';
1834 if ((sector/18)%2==1) sideName='C';
1837 sensorName=Form("TPC_GATE_I_%c_POS_VMEAS",sideName);
1840 sensorName=Form("TPC_GATE_O_%c_POS_VMEAS",sideName);
1843 val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1845 val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1850 Float_t AliTPCcalibDB::GetPressure(Int_t timeStamp, Int_t run, Int_t type){
1852 // GetPressure for given time stamp and runt
1854 TTimeStamp stamp(timeStamp);
1855 AliDCSSensor * sensor = Instance()->GetPressureSensor(run,type);
1856 if (!sensor) return 0;
1857 return sensor->GetValue(stamp);
1860 Float_t AliTPCcalibDB::GetL3Current(Int_t run, Int_t statType){
1862 // return L3 current
1863 // stat type is: AliGRPObject::Stats: kMean = 0, kTruncMean = 1, kMedian = 2, kSDMean = 3, kSDMedian = 4
1866 AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1867 if (grp) current=grp->GetL3Current((AliGRPObject::Stats)statType);
1871 Float_t AliTPCcalibDB::GetBz(Int_t run){
1873 // calculate BZ in T from L3 current
1876 Float_t current=AliTPCcalibDB::GetL3Current(run);
1877 if (current>-1) bz=5*current/30000.*.1;
1881 Char_t AliTPCcalibDB::GetL3Polarity(Int_t run) {
1883 // get l3 polarity from GRP
1886 AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1887 if (grp) pol=grp->GetL3Polarity();
1891 TString AliTPCcalibDB::GetRunType(Int_t run){
1893 // return run type from grp
1896 // TString type("UNKNOWN");
1897 AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1898 if (grp) return grp->GetRunType();
1902 Float_t AliTPCcalibDB::GetValueGoofie(Int_t timeStamp, Int_t run, Int_t type){
1904 // GetPressure for given time stamp and runt
1906 TTimeStamp stamp(timeStamp);
1907 AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(run);
1908 if (!goofieArray) return 0;
1909 AliDCSSensor *sensor = goofieArray->GetSensor(type);
1910 return sensor->GetValue(stamp);
1918 Bool_t AliTPCcalibDB::GetTemperatureFit(Int_t timeStamp, Int_t run, Int_t side,TVectorD& fit){
1920 // GetTmeparature fit at parameter for given time stamp
1922 TTimeStamp tstamp(timeStamp);
1923 AliTPCSensorTempArray* tempArray = Instance()->GetTemperatureSensor(run);
1924 if (! tempArray) return kFALSE;
1925 AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
1926 TLinearFitter * fitter = tempMap->GetLinearFitter(3,side,tstamp);
1929 fitter->GetParameters(fit);
1933 if (!fitter) return kFALSE;
1937 Float_t AliTPCcalibDB::GetTemperature(Int_t timeStamp, Int_t run, Int_t side){
1939 // Get mean temperature
1943 GetTemperatureFit(timeStamp,run,0,vec);
1947 GetTemperatureFit(timeStamp,run,0,vec);
1954 Double_t AliTPCcalibDB::GetPTRelative(UInt_t timeSec, Int_t run, Int_t side){
1957 // time - absolute time
1959 // side - 0 - A side 1-C side
1960 AliTPCCalibVdrift * vdrift = Instance()->GetVdrift(run);
1961 if (!vdrift) return 0;
1962 return vdrift->GetPTRelative(timeSec,side);
1965 AliGRPObject * AliTPCcalibDB::MakeGRPObjectFromMap(TMap *map){
1967 // Function to covert old GRP run information from TMap to GRPObject
1969 // TMap * map = AliTPCcalibDB::GetGRPMap(52406);
1971 AliDCSSensor * sensor = 0;
1973 osensor = ((*map)("fP2Pressure"));
1974 sensor =dynamic_cast<AliDCSSensor *>(osensor);
1976 if (!sensor) return 0;
1978 AliDCSSensor * sensor2 = new AliDCSSensor(*sensor);
1979 osensor = ((*map)("fCavernPressure"));
1980 TGraph * gr = new TGraph(2);
1981 gr->GetX()[0]= -100000.;
1982 gr->GetX()[1]= 1000000.;
1983 gr->GetY()[0]= atof(osensor->GetName());
1984 gr->GetY()[1]= atof(osensor->GetName());
1985 sensor2->SetGraph(gr);
1989 AliGRPObject *grpRun = new AliGRPObject;
1990 grpRun->ReadValuesFromMap(map);
1991 grpRun->SetCavernAtmosPressure(sensor2);
1992 grpRun->SetCavernAtmosPressure(sensor2);
1993 grpRun->SetSurfaceAtmosPressure(sensor);
1997 Bool_t AliTPCcalibDB::CreateGUITree(Int_t run, const char* filename)
2000 // Create a gui tree for run number 'run'
2003 if (!AliCDBManager::Instance()->GetDefaultStorage()){
2004 AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
2005 MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
2009 AliTPCcalibDB *db=AliTPCcalibDB::Instance();
2010 // retrieve cal pad objects
2012 db->CreateGUITree(filename);
2016 Bool_t AliTPCcalibDB::CreateGUITree(const char* filename){
2020 if (!AliCDBManager::Instance()->GetDefaultStorage()){
2021 AliError("Default Storage not set. Cannot create calibration Tree!");
2024 UpdateNonRec(); // load all infromation now
2026 AliTPCPreprocessorOnline prep;
2027 //noise and pedestals
2028 if (GetPedestals()) prep.AddComponent(new AliTPCCalPad(*(GetPedestals())));
2029 if (GetPadNoise() ) prep.AddComponent(new AliTPCCalPad(*(GetPadNoise())));
2031 if (GetPulserTmean()) prep.AddComponent(new AliTPCCalPad(*(GetPulserTmean())));
2032 if (GetPulserTrms() ) prep.AddComponent(new AliTPCCalPad(*(GetPulserTrms())));
2033 if (GetPulserQmean()) prep.AddComponent(new AliTPCCalPad(*(GetPulserQmean())));
2035 if (GetCETmean()) prep.AddComponent(new AliTPCCalPad(*(GetCETmean())));
2036 if (GetCETrms() ) prep.AddComponent(new AliTPCCalPad(*(GetCETrms())));
2037 if (GetCEQmean()) prep.AddComponent(new AliTPCCalPad(*(GetCEQmean())));
2039 if (GetALTROAcqStart() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROAcqStart() )));
2040 if (GetALTROZsThr() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROZsThr() )));
2041 if (GetALTROFPED() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROFPED() )));
2042 if (GetALTROAcqStop() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROAcqStop() )));
2043 if (GetALTROMasked() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROMasked() )));
2045 AliTPCdataQA *dataQA=GetDataQA();
2047 if (dataQA->GetNLocalMaxima())
2048 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNLocalMaxima())));
2049 if (dataQA->GetMaxCharge())
2050 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetMaxCharge())));
2051 if (dataQA->GetMeanCharge())
2052 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetMeanCharge())));
2053 if (dataQA->GetNoThreshold())
2054 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNoThreshold())));
2055 if (dataQA->GetNTimeBins())
2056 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNTimeBins())));
2057 if (dataQA->GetNPads())
2058 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNPads())));
2059 if (dataQA->GetTimePosition())
2060 prep.AddComponent(new AliTPCCalPad(*(dataQA->GetTimePosition())));
2064 TString file(filename);
2065 if (file.IsNull()) file=Form("guiTreeRun_%i.root",fRun);
2066 prep.DumpToFile(file.Data());
2070 Bool_t AliTPCcalibDB::CreateRefFile(Int_t run, const char* filename)
2073 // Create a gui tree for run number 'run'
2076 if (!AliCDBManager::Instance()->GetDefaultStorage()){
2077 AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
2078 MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
2081 TString file(filename);
2082 if (file.IsNull()) file=Form("RefCalPads_%d.root",run);
2083 TDirectory *currDir=gDirectory;
2085 AliTPCcalibDB *db=AliTPCcalibDB::Instance();
2086 // retrieve cal pad objects
2089 TFile f(file.Data(),"recreate");
2090 //noise and pedestals
2091 db->GetPedestals()->Write("Pedestals");
2092 db->GetPadNoise()->Write("PadNoise");
2094 db->GetPulserTmean()->Write("PulserTmean");
2095 db->GetPulserTrms()->Write("PulserTrms");
2096 db->GetPulserQmean()->Write("PulserQmean");
2098 db->GetCETmean()->Write("CETmean");
2099 db->GetCETrms()->Write("CETrms");
2100 db->GetCEQmean()->Write("CEQmean");
2102 db->GetALTROAcqStart() ->Write("ALTROAcqStart");
2103 db->GetALTROZsThr() ->Write("ALTROZsThr");
2104 db->GetALTROFPED() ->Write("ALTROFPED");
2105 db->GetALTROAcqStop() ->Write("ALTROAcqStop");
2106 db->GetALTROMasked() ->Write("ALTROMasked");
2115 Double_t AliTPCcalibDB::GetVDriftCorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
2117 // Get time dependent drift velocity correction
2118 // multiplication factor vd = vdnom *(1+vdriftcorr)
2120 // mode determines the algorith how to combine the Laser Track, LaserCE and physics tracks
2121 // timestamp - timestamp
2123 // side - the drift velocity per side (possible for laser and CE)
2125 // Notice - Extrapolation outside of calibration range - using constant function
2128 // mode 1 automatic mode - according to the distance to the valid calibration
2130 Double_t deltaP=0, driftP=0, wP = 0.;
2131 Double_t deltaITS=0,driftITS=0, wITS= 0.;
2132 Double_t deltaLT=0, driftLT=0, wLT = 0.;
2133 Double_t deltaCE=0, driftCE=0, wCE = 0.;
2134 driftP = fDButil->GetVDriftTPC(deltaP,run,timeStamp);
2135 driftITS= fDButil->GetVDriftTPCITS(deltaITS,run,timeStamp);
2136 driftCE = fDButil->GetVDriftTPCCE(deltaCE, run,timeStamp,36000,2);
2137 driftLT = fDButil->GetVDriftTPCLaserTracks(deltaLT,run,timeStamp,36000,2);
2138 deltaITS = TMath::Abs(deltaITS);
2139 deltaP = TMath::Abs(deltaP);
2140 deltaLT = TMath::Abs(deltaLT);
2141 deltaCE = TMath::Abs(deltaCE);
2143 const Double_t kEpsilon=0.00000000001;
2144 const Double_t kdeltaT=360.; // 10 minutes
2145 if(TMath::Abs(deltaITS) < 12*kdeltaT) {
2148 wITS = 64.*kdeltaT/(deltaITS +kdeltaT);
2149 wLT = 16.*kdeltaT/(deltaLT +kdeltaT);
2150 wP = 0. *kdeltaT/(deltaP +kdeltaT);
2151 wCE = 1. *kdeltaT/(deltaCE +kdeltaT);
2154 if (TMath::Abs(driftP)<kEpsilon) wP=0; // invalid calibration
2155 if (TMath::Abs(driftITS)<kEpsilon)wITS=0; // invalid calibration
2156 if (TMath::Abs(driftLT)<kEpsilon) wLT=0; // invalid calibration
2157 if (TMath::Abs(driftCE)<kEpsilon) wCE=0; // invalid calibration
2158 if (wP+wITS+wLT+wCE<kEpsilon) return 0;
2159 result = (driftP*wP+driftITS*wITS+driftLT*wLT+driftCE*wCE)/(wP+wITS+wLT+wCE);
2168 Double_t AliTPCcalibDB::GetTime0CorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
2170 // Get time dependent time 0 (trigger delay in cm) correction
2171 // additive correction time0 = time0+ GetTime0CorrectionTime
2172 // Value etracted combining the vdrift correction using laser tracks and CE and the physics track matchin
2174 // mode determines the algorith how to combine the Laser Track and physics tracks
2175 // timestamp - timestamp
2177 // side - the drift velocity per side (possible for laser and CE)
2179 // Notice - Extrapolation outside of calibration range - using constant function
2184 result=fDButil->GetTriggerOffsetTPC(run,timeStamp);
2185 result *=fParam->GetZLength();
2190 result= -fDButil->GetTime0TPCITS(dist, run, timeStamp)*fParam->GetDriftV()/1000000.;
2199 Double_t AliTPCcalibDB::GetVDriftCorrectionGy(Int_t timeStamp, Int_t run, Int_t side, Int_t /*mode*/){
2201 // Get global y correction drift velocity correction factor
2202 // additive factor vd = vdnom*(1+GetVDriftCorrectionGy *gy)
2203 // Value etracted combining the vdrift correction using laser tracks and CE or TPC-ITS
2205 // mode determines the algorith how to combine the Laser Track, LaserCE or TPC-ITS
2206 // timestamp - timestamp
2208 // side - the drift velocity gy correction per side (CE and Laser tracks)
2210 // Notice - Extrapolation outside of calibration range - using constant function
2212 if (run<=0 && fTransform) run = fTransform->GetCurrentRunNumber();
2213 UpdateRunInformations(run,kFALSE);
2214 TObjArray *array =AliTPCcalibDB::Instance()->GetTimeVdriftSplineRun(run);
2215 if (!array) return 0;
2218 // use TPC-ITS if present
2219 TGraphErrors *gr= (TGraphErrors*)array->FindObject("ALIGN_ITSB_TPC_VDGY");
2220 if (!gr) gr = (TGraphErrors*)array->FindObject("ALIGN_TOFB_TPC_VDGY");
2222 result = AliTPCcalibDButil::EvalGraphConst(gr,timeStamp);
2224 // transform from [(cm/mus)/ m] to [1/cm]
2225 result /= (fParam->GetDriftV()/1000000.);
2228 //printf("result %e \n", result);
2232 // use laser if ITS-TPC not present
2233 TGraphErrors *laserA= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
2234 TGraphErrors *laserC= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
2236 if (laserA && laserC){
2237 result= (laserA->Eval(timeStamp)+laserC->Eval(timeStamp))*0.5;
2239 if (laserA && side==0){
2240 result = (laserA->Eval(timeStamp));
2242 if (laserC &&side==1){
2243 result = (laserC->Eval(timeStamp));
2245 //printf("laser result %e \n", -result/250.);
2247 return -result/250.; //normalized before
2250 AliTPCCalPad* AliTPCcalibDB::MakeDeadMap(Double_t notInMap, const char* nameMappingFile) {
2252 // Read list of active DDLs from OCDB entry
2253 // Generate and return AliTPCCalPad containing 1 for all pads in active DDLs,
2254 // 0 for all pads in non-active DDLs.
2255 // For DDLs with missing status information (no DCS input point to Shuttle),
2256 // the value of the AliTPCCalPad entry is determined by the parameter
2257 // notInMap (default value 1)
2261 TFile *fileMapping = new TFile(nameMappingFile, "read");
2262 AliTPCmapper *mapping = (AliTPCmapper*) fileMapping->Get("tpcMapping");
2264 snprintf(chinfo,1000,"Failed to get mapping object from %s. ...\n", nameMappingFile);
2269 AliTPCCalPad *deadMap = new AliTPCCalPad("deadMap","deadMap");
2271 AliError("Failed to allocate dead map AliTPCCalPad");
2275 /// get list of active DDLs from OCDB entry
2277 if (!fALTROConfigData ) {
2278 AliError("No ALTRO config OCDB entry available");
2281 TMap *activeDDL = (TMap*)fALTROConfigData->FindObject("DDLArray");
2282 TObjString *ddlArray=0;
2284 ddlArray = (TObjString*)activeDDL->GetValue("DDLArray");
2286 AliError("Empty list of active DDLs in OCDB entry");
2290 AliError("List of active DDLs not available in OCDB entry");
2293 TString arrDDL=ddlArray->GetString();
2294 Int_t offset = mapping->GetTpcDdlOffset();
2296 for (Int_t i=0; i<mapping->GetNumDdl(); i++) {
2298 if (idDDL<0) continue;
2299 Int_t patch = mapping->GetPatchFromEquipmentID(idDDL);
2300 if (patch<0) continue;
2301 Int_t roc=mapping->GetRocFromEquipmentID(idDDL);
2302 if (roc<0) continue;
2303 AliTPCCalROC *calRoc=deadMap->GetCalROC(roc);
2305 for ( Int_t branch = 0; branch < 2; branch++ ) {
2306 for ( Int_t fec = 0; fec < mapping->GetNfec(patch, branch); fec++ ) {
2307 for ( Int_t altro = 0; altro < 8; altro++ ) {
2308 for ( Int_t channel = 0; channel < 16; channel++ ) {
2309 Int_t hwadd = mapping->CodeHWAddress(branch, fec, altro, channel);
2310 Int_t row = mapping->GetPadRow(patch, hwadd); // row in a ROC (IROC or OROC)
2311 // Int_t globalrow = mapping.GetGlobalPadRow(patch, hwadd); // row in full sector (IROC plus OROC)
2312 Int_t pad = mapping->GetPad(patch, hwadd);
2313 if (!TString(arrDDL[i]).IsDigit()) {
2316 active=TString(arrDDL[i]).Atof();
2318 calRoc->SetValue(row,pad,active);
2319 } // end channel for loop
2320 } // end altro for loop
2321 } // end fec for loop
2322 } // end branch for loop
2324 } // end loop on active DDLs
2330 AliTPCCorrection * AliTPCcalibDB::GetTPCComposedCorrection(Float_t field) const{
2332 // GetComposed correction for given field setting
2333 // If not specific correction for field used return correction for all field
2334 // - Complication needed to gaurantee OCDB back compatibility
2335 // - Not neeeded for the new space point correction
2336 if (!fComposedCorrectionArray) return 0;
2337 if (field>0.1 && fComposedCorrectionArray->At(1)) {
2338 return (AliTPCCorrection *)fComposedCorrectionArray->At(1);
2340 if (field<-0.1 &&fComposedCorrectionArray->At(2)) {
2341 return (AliTPCCorrection *)fComposedCorrectionArray->At(2);
2343 return (AliTPCCorrection *)fComposedCorrectionArray->At(0);
2348 AliTPCCorrection * AliTPCcalibDB::GetTPCComposedCorrectionDelta() const{
2350 // GetComposedCorrection delta
2351 // Delta is time dependent - taken form the CalibTime OCDB entry
2353 if (!fComposedCorrectionArray) return 0;
2354 if (fRun<0) return 0;
2355 if (fDriftCorrectionArray.GetValue(Form("%i",fRun))==0) return 0;
2356 if (fComposedCorrectionArray->GetEntriesFast()<=4) {
2357 fComposedCorrectionArray->Expand(5);
2358 TObjArray * timeArray =(TObjArray*)(fDriftCorrectionArray.GetValue(Form("%i",fRun)));
2359 AliTPCCorrection * correctionTime = (AliTPCCorrection *)timeArray->FindObject("FitCorrectionTime");
2360 if (correctionTime){
2361 correctionTime->Init();
2362 fComposedCorrectionArray->AddAt(correctionTime,4); //add time dependent c
2365 return (AliTPCCorrection *)fComposedCorrectionArray->At(4); //