]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibDB.cxx
AliTPCcalibDB.cxx.difF Add possibility to get RCU configuration, RCU Trigger...
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibDB.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 // Class providing the calibration parameters by accessing the CDB           //
20 //                                                                           //
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                                         ////
24 //
25 //
26 // Calibration data:
27 // 0.)  Altro mapping
28 //          Simulation      - not yet 
29 //          Reconstruction  - AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
30 //
31 // 1.)  pad by pad calibration -  AliTPCCalPad
32 //      
33 //      a.) fPadGainFactor
34 //          Simulation: AliTPCDigitizer::ExecFast - Multiply by gain
35 //          Reconstruction : AliTPCclustererMI::Digits2Clusters - Divide by gain  
36 //
37 //      b.) fPadNoise -
38 //          Simulation:        AliTPCDigitizer::ExecFast
39 //          Reconstruction:    AliTPCclustererMI::FindClusters(AliTPCCalROC * noiseROC)
40 //                             Noise depending cut on clusters charge (n sigma)
41 //      c.) fPedestal:
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();
46 //      
47 //      d.) fPadTime0
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()      
52 //
53 // 
54 // 2.)  Space points transformation:
55 //
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
60 //                 ExB effect    
61 //          Simulation     - Not used directly (the effects are applied one by one (see AliTPC::MakeSector)
62 //          Reconstruction - 
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 -  
70 //                  
71 //                  in AliTPCtransform::Correct() - called calib->GetExB()->Correct(dxyz0,dxyz1)
72 //
73 //  3.)   cluster error, shape and Q parameterization
74 //
75 //
76 //
77 ///////////////////////////////////////////////////////////////////////////////
78
79 #include <iostream>
80 #include <fstream>
81
82
83 #include <AliCDBManager.h>
84 #include <AliCDBEntry.h>
85 #include <AliCDBId.h>
86 #include <AliLog.h>
87 #include <AliMagF.h>
88 #include <AliSplineFit.h>
89 #include <AliCTPTimeParams.h>
90
91 #include "AliTPCcalibDB.h"
92 #include "AliTPCdataQA.h"
93 #include "AliTPCcalibDButil.h"
94 #include "AliTPCAltroMapping.h"
95 #include "AliTPCExB.h"
96
97 #include "AliTPCCalROC.h"
98 #include "AliTPCCalPad.h"
99 #include "AliTPCSensorTempArray.h"
100 #include "AliGRPObject.h"
101 #include "AliTPCTransform.h"
102
103 class AliCDBStorage;
104 class AliTPCCalDet;
105 //
106 //
107
108 #include "TFile.h"
109 #include "TKey.h"
110 #include "TGraphErrors.h"
111
112 #include "TObjArray.h"
113 #include "TObjString.h"
114 #include "TString.h"
115 #include "TDirectory.h"
116 #include "AliTPCCalPad.h"
117 #include "AliTPCCalibPulser.h"
118 #include "AliTPCCalibPedestal.h"
119 #include "AliTPCCalibCE.h"
120 #include "AliTPCExBFirst.h"
121 #include "AliTPCTempMap.h"
122 #include "AliTPCCalibVdrift.h"
123 #include "AliTPCCalibRaw.h"
124 #include "AliTPCParam.h"
125
126 #include "AliTPCPreprocessorOnline.h"
127
128
129 ClassImp(AliTPCcalibDB)
130
131 AliTPCcalibDB* AliTPCcalibDB::fgInstance = 0;
132 Bool_t AliTPCcalibDB::fgTerminated = kFALSE;
133 TObjArray    AliTPCcalibDB::fgExBArray;    // array of ExB corrections
134
135
136 //_ singleton implementation __________________________________________________
137 AliTPCcalibDB* AliTPCcalibDB::Instance()
138 {
139   //
140   // Singleton implementation
141   // Returns an instance of this class, it is created if neccessary
142   //
143   
144   if (fgTerminated != kFALSE)
145     return 0;
146
147   if (fgInstance == 0)
148     fgInstance = new AliTPCcalibDB();
149   
150   return fgInstance;
151 }
152
153 void AliTPCcalibDB::Terminate()
154 {
155   //
156   // Singleton implementation
157   // Deletes the instance of this class and sets the terminated flag, instances cannot be requested anymore
158   // This function can be called several times.
159   //
160   
161   fgTerminated = kTRUE;
162   
163   if (fgInstance != 0)
164   {
165     delete fgInstance;
166     fgInstance = 0;
167   }
168 }
169
170 //_____________________________________________________________________________
171 AliTPCcalibDB::AliTPCcalibDB():
172   TObject(),
173   fRun(-1),
174   fTransform(0),
175   fExB(0),
176   fPadGainFactor(0),
177   fDedxGainFactor(0),
178   fPadTime0(0),
179   fPadNoise(0),
180   fPedestals(0),
181   fCalibRaw(0),
182   fDataQA(0),
183   fALTROConfigData(0),
184   fPulserData(0),
185   fCEData(0),
186   fTemperature(0),
187   fMapping(0),
188   fParam(0),
189   fClusterParam(0),
190   fTimeGainSplines(0),
191   fTimeGainSplinesArray(100000),
192   fGRPArray(100000),            //! array of GRPs  -  per run  - JUST for calibration studies
193   fGRPMaps(100000),            //! array of GRPs  -  per run  - JUST for calibration studies
194   fGoofieArray(100000),         //! array of GOOFIE values -per run - Just for calibration studies
195   fVoltageArray(100000),
196   fTemperatureArray(100000),    //! array of temperature sensors - per run - Just for calibration studies
197   fVdriftArray(100000),                 //! array of v drift interfaces
198   fDriftCorrectionArray(100000),  //! array of drift correction
199   fRunList(100000),              //! run list - indicates try to get the run param 
200   fDButil(0),
201   fCTPTimeParams(0)
202 {
203   //
204   // constructor
205   //  
206   //
207   Update();    // temporary
208 }
209
210 AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& ):
211   TObject(),
212   fRun(-1),
213   fTransform(0),
214   fExB(0),
215   fPadGainFactor(0),
216   fDedxGainFactor(0),
217   fPadTime0(0),
218   fPadNoise(0),
219   fPedestals(0),
220   fCalibRaw(0),
221   fDataQA(0),
222   fALTROConfigData(0),
223   fPulserData(0),
224   fCEData(0),
225   fTemperature(0),
226   fMapping(0),
227   fParam(0),
228   fClusterParam(0),
229   fTimeGainSplines(0),
230   fTimeGainSplinesArray(100000),
231   fGRPArray(0),          //! array of GRPs  -  per run  - JUST for calibration studies
232   fGRPMaps(0),          //! array of GRPs  -  per run  - JUST for calibration studies
233   fGoofieArray(0),        //! array of GOOFIE values -per run - Just for calibration studies
234   fVoltageArray(0),
235   fTemperatureArray(0),   //! array of temperature sensors - per run - Just for calibration studies
236   fVdriftArray(0),         //! array of v drift interfaces
237   fDriftCorrectionArray(0),         //! array of v drift interfaces
238   fRunList(0),              //! run list - indicates try to get the run param 
239   fDButil(0),
240   fCTPTimeParams(0)
241 {
242   //
243   // Copy constructor invalid -- singleton implementation
244   //
245    Error("copy constructor","invalid -- singleton implementation");
246 }
247
248 AliTPCcalibDB& AliTPCcalibDB::operator= (const AliTPCcalibDB& )
249 {
250 //
251 // Singleton implementation - no assignment operator
252 //
253   Error("operator =", "assignment operator not implemented");
254   return *this;
255 }
256
257
258
259 //_____________________________________________________________________________
260 AliTPCcalibDB::~AliTPCcalibDB() 
261 {
262   //
263   // destructor
264   //
265   
266   // don't delete anything, CDB cache is active!
267   //if (fPadGainFactor) delete fPadGainFactor;
268   //if (fPadTime0) delete fPadTime0;
269   //if (fPadNoise) delete fPadNoise;
270 }
271
272
273 //_____________________________________________________________________________
274 AliCDBEntry* AliTPCcalibDB::GetCDBEntry(const char* cdbPath)
275 {
276   // 
277   // Retrieves an entry with path <cdbPath> from the CDB.
278   //
279   char chinfo[1000];
280     
281   AliCDBEntry* entry = AliCDBManager::Instance()->Get(cdbPath, fRun); 
282   if (!entry) 
283   { 
284     sprintf(chinfo,"AliTPCcalibDB: Failed to get entry:\t%s ", cdbPath);
285     AliError(chinfo); 
286     return 0; 
287   }
288   return entry;
289 }
290
291
292 //_____________________________________________________________________________
293 void AliTPCcalibDB::SetRun(Long64_t run)
294 {
295   //
296   // Sets current run number. Calibration data is read from the corresponding file. 
297   //  
298   if (fRun == run)
299     return;  
300         fRun = run;
301   Update();
302 }
303   
304
305
306 void AliTPCcalibDB::Update(){
307         //
308   AliCDBEntry * entry=0;
309   Bool_t cdbCache = AliCDBManager::Instance()->GetCacheFlag(); // save cache status
310   AliCDBManager::Instance()->SetCacheFlag(kTRUE); // activate CDB cache
311   fDButil = new AliTPCcalibDButil;   
312   //
313   entry          = GetCDBEntry("TPC/Calib/PadGainFactor");
314   if (entry){
315     //if (fPadGainFactor) delete fPadGainFactor;
316     entry->SetOwner(kTRUE);
317     fPadGainFactor = (AliTPCCalPad*)entry->GetObject();
318   }
319   //
320   entry          = GetCDBEntry("TPC/Calib/TimeGain");
321   if (entry){
322     //if (fTimeGainSplines) delete fTimeGainSplines;
323     entry->SetOwner(kTRUE);
324     fTimeGainSplines = (TObjArray*)entry->GetObject();
325   }
326   //
327   entry          = GetCDBEntry("TPC/Calib/GainFactorDedx");
328   if (entry){
329     entry->SetOwner(kTRUE);
330     fDedxGainFactor = (AliTPCCalPad*)entry->GetObject();
331   }
332   //
333   entry          = GetCDBEntry("TPC/Calib/PadTime0");
334   if (entry){
335     //if (fPadTime0) delete fPadTime0;
336     entry->SetOwner(kTRUE);
337     fPadTime0 = (AliTPCCalPad*)entry->GetObject();
338   }
339   //
340   //
341   entry          = GetCDBEntry("TPC/Calib/PadNoise");
342   if (entry){
343     //if (fPadNoise) delete fPadNoise;
344     entry->SetOwner(kTRUE);
345     fPadNoise = (AliTPCCalPad*)entry->GetObject();
346   }
347
348   entry          = GetCDBEntry("TPC/Calib/Pedestals");
349   if (entry){
350     //if (fPedestals) delete fPedestals;
351     entry->SetOwner(kTRUE);
352     fPedestals = (AliTPCCalPad*)entry->GetObject();
353   }
354
355   entry          = GetCDBEntry("TPC/Calib/Temperature");
356   if (entry){
357     //if (fTemperature) delete fTemperature;
358     entry->SetOwner(kTRUE);
359     fTemperature = (AliTPCSensorTempArray*)entry->GetObject();
360   }
361
362   entry          = GetCDBEntry("TPC/Calib/Parameters");
363   if (entry){
364     //if (fPadNoise) delete fPadNoise;
365     entry->SetOwner(kTRUE);
366     fParam = (AliTPCParam*)(entry->GetObject()->Clone());
367   }
368
369   entry          = GetCDBEntry("TPC/Calib/ClusterParam");
370   if (entry){
371     entry->SetOwner(kTRUE);
372     fClusterParam = (AliTPCClusterParam*)(entry->GetObject()->Clone());
373   }
374
375   //ALTRO configuration data
376   entry          = GetCDBEntry("TPC/Calib/AltroConfig");
377   if (entry){
378     entry->SetOwner(kTRUE);
379     fALTROConfigData=(TObjArray*)(entry->GetObject());
380   }
381   
382   //Calibration Pulser data
383   entry          = GetCDBEntry("TPC/Calib/Pulser");
384   if (entry){
385     entry->SetOwner(kTRUE);
386     fPulserData=(TObjArray*)(entry->GetObject());
387   }
388   
389   //CE data
390   entry          = GetCDBEntry("TPC/Calib/CE");
391   if (entry){
392     entry->SetOwner(kTRUE);
393     fCEData=(TObjArray*)(entry->GetObject());
394   }
395   //RAW calibration data
396   entry          = GetCDBEntry("TPC/Calib/Raw");
397   if (entry){
398     entry->SetOwner(kTRUE);
399     TObjArray *arr=(TObjArray*)(entry->GetObject());
400     if (arr) fCalibRaw=(AliTPCCalibRaw*)arr->At(0);
401   }
402   //QA calibration data
403   entry          = GetCDBEntry("TPC/Calib/QA");
404   if (entry){
405     entry->SetOwner(kTRUE);
406     fDataQA=dynamic_cast<AliTPCdataQA*>(entry->GetObject());
407   }
408   
409   entry          = GetCDBEntry("TPC/Calib/Mapping");
410   if (entry){
411     //if (fPadNoise) delete fPadNoise;
412     entry->SetOwner(kTRUE);
413     TObjArray * array = dynamic_cast<TObjArray*>(entry->GetObject());
414     if (array && array->GetEntriesFast()==6){
415       fMapping = new AliTPCAltroMapping*[6];
416       for (Int_t i=0; i<6; i++){
417         fMapping[i] =  dynamic_cast<AliTPCAltroMapping*>(array->At(i));
418       }
419     }
420   }
421
422   //QA calibration data
423   entry          = GetCDBEntry("GRP/CTP/CTPtiming");
424   if (entry){
425     entry->SetOwner(kTRUE);
426     fCTPTimeParams=dynamic_cast<AliCTPTimeParams*>(entry->GetObject());
427   }
428   
429
430   //entry          = GetCDBEntry("TPC/Calib/ExB");
431   //if (entry) {
432   //  entry->SetOwner(kTRUE);
433   //  fExB=dynamic_cast<AliTPCExB*>(entry->GetObject()->Clone());
434   //}
435   //
436   // ExB  - calculate during initialization - in simulation /reconstruction
437   //      - not invoked here anymore
438   //fExB =  GetExB(-5,kTRUE);
439      //
440   if (!fTransform) {
441     fTransform=new AliTPCTransform(); 
442     fTransform->SetCurrentRun(AliCDBManager::Instance()->GetRun());
443   }
444
445   //
446   AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
447 }
448
449
450
451 void AliTPCcalibDB::CreateObjectList(const Char_t *filename, TObjArray *calibObjects)
452 {
453 //
454 // Create calibration objects and read contents from OCDB
455 //
456    if ( calibObjects == 0x0 ) return;
457    ifstream in;
458    in.open(filename);
459    if ( !in.is_open() ){
460       fprintf(stderr,"Error: cannot open list file '%s'", filename);
461       return;
462    }
463    
464    AliTPCCalPad *calPad=0x0;
465    
466    TString sFile;
467    sFile.ReadFile(in);
468    in.close();
469    
470    TObjArray *arrFileLine = sFile.Tokenize("\n");
471    
472    TIter nextLine(arrFileLine);
473    
474    TObjString *sObjLine=0x0;
475    while ( (sObjLine = (TObjString*)nextLine()) ){
476       TString sLine(sObjLine->GetString());
477       
478       TObjArray *arrNextCol = sLine.Tokenize("\t");
479       
480       TObjString *sObjType     = (TObjString*)(arrNextCol->At(0));
481       TObjString *sObjFileName = (TObjString*)(arrNextCol->At(1));
482       
483       if ( !sObjType || ! sObjFileName ) continue;
484       TString sType(sObjType->GetString());
485       TString sFileName(sObjFileName->GetString());
486       printf("%s\t%s\n",sType.Data(),sFileName.Data());
487       
488       TFile *fIn = TFile::Open(sFileName);
489       if ( !fIn ){
490          fprintf(stderr,"File not found: '%s'", sFileName.Data());
491          continue;
492       }
493       
494       if ( sType == "CE" ){
495          AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
496          
497          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());         
498          calPad->SetNameTitle("CETmean","CETmean");
499          calibObjects->Add(calPad);
500          
501          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());         
502          calPad->SetNameTitle("CEQmean","CEQmean");
503          calibObjects->Add(calPad);        
504          
505          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
506          calPad->SetNameTitle("CETrms","CETrms");
507          calibObjects->Add(calPad);         
508                   
509       } else if ( sType == "Pulser") {
510          AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
511          
512          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());         
513          calPad->SetNameTitle("PulserTmean","PulserTmean");
514          calibObjects->Add(calPad);
515          
516          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());         
517          calPad->SetNameTitle("PulserQmean","PulserQmean");
518          calibObjects->Add(calPad);        
519          
520          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
521          calPad->SetNameTitle("PulserTrms","PulserTrms");
522          calibObjects->Add(calPad);         
523       
524       } else if ( sType == "Pedestals") {
525          AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
526          
527          calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());         
528          calPad->SetNameTitle("Pedestals","Pedestals");
529          calibObjects->Add(calPad);
530          
531          calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());         
532          calPad->SetNameTitle("Noise","Noise");
533          calibObjects->Add(calPad);        
534      
535       } else {
536          fprintf(stderr,"Undefined Type: '%s'",sType.Data());
537          
538       }
539       delete fIn;
540    }
541 }
542
543
544
545 void AliTPCcalibDB::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
546   //
547   // Write a tree with all available information
548   // if mapFileName is specified, the Map information are also written to the tree
549   // pads specified in outlierPad are not used for calculating statistics
550   //  - the same function as AliTPCCalPad::MakeTree - 
551   //
552    AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
553
554    TObjArray* mapIROCs = 0;
555    TObjArray* mapOROCs = 0;
556    TVectorF *mapIROCArray = 0;
557    TVectorF *mapOROCArray = 0;
558    Int_t mapEntries = 0;
559    TString* mapNames = 0;
560    
561    if (mapFileName) {
562       TFile mapFile(mapFileName, "read");
563       
564       TList* listOfROCs = mapFile.GetListOfKeys();
565       mapEntries = listOfROCs->GetEntries()/2;
566       mapIROCs = new TObjArray(mapEntries*2);
567       mapOROCs = new TObjArray(mapEntries*2);
568       mapIROCArray = new TVectorF[mapEntries];
569       mapOROCArray = new TVectorF[mapEntries];
570       
571       mapNames = new TString[mapEntries];
572       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
573         TString nameROC(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
574          nameROC.Remove(nameROC.Length()-4, 4);
575          mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "IROC").Data()), ivalue);
576          mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "OROC").Data()), ivalue);
577          mapNames[ivalue].Append(nameROC);
578       }
579       
580       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
581          mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
582          mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
583       
584          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
585             (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
586          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
587             (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
588       }
589
590    } //  if (mapFileName)
591   
592    TTreeSRedirector cstream(fileName);
593    Int_t arrayEntries = array->GetEntries();
594    
595    TString* names = new TString[arrayEntries];
596    for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
597       names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
598
599    for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
600       //
601       // get statistic for given sector
602       //
603       TVectorF median(arrayEntries);
604       TVectorF mean(arrayEntries);
605       TVectorF rms(arrayEntries);
606       TVectorF ltm(arrayEntries);
607       TVectorF ltmrms(arrayEntries);
608       TVectorF medianWithOut(arrayEntries);
609       TVectorF meanWithOut(arrayEntries);
610       TVectorF rmsWithOut(arrayEntries);
611       TVectorF ltmWithOut(arrayEntries);
612       TVectorF ltmrmsWithOut(arrayEntries);
613       
614       TVectorF *vectorArray = new TVectorF[arrayEntries];
615       for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
616          vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
617       
618       for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
619          AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
620          AliTPCCalROC* calROC = calPad->GetCalROC(isector);
621          AliTPCCalROC* outlierROC = 0;
622          if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
623          if (calROC) {
624             median[ivalue] = calROC->GetMedian();
625             mean[ivalue] = calROC->GetMean();
626             rms[ivalue] = calROC->GetRMS();
627             Double_t ltmrmsValue = 0;
628             ltm[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction);
629             ltmrms[ivalue] = ltmrmsValue;
630             if (outlierROC) {
631                medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
632                meanWithOut[ivalue] = calROC->GetMean(outlierROC);
633                rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
634                ltmrmsValue = 0;
635                ltmWithOut[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction, outlierROC);
636                ltmrmsWithOut[ivalue] = ltmrmsValue;
637             }
638          }
639          else {
640             median[ivalue] = 0.;
641             mean[ivalue] = 0.;
642             rms[ivalue] = 0.;
643             ltm[ivalue] = 0.;
644             ltmrms[ivalue] = 0.;
645             medianWithOut[ivalue] = 0.;
646             meanWithOut[ivalue] = 0.;
647             rmsWithOut[ivalue] = 0.;
648             ltmWithOut[ivalue] = 0.;
649             ltmrmsWithOut[ivalue] = 0.;
650          }
651       }
652       
653       //
654       // fill vectors of variable per pad
655       //
656       TVectorF *posArray = new TVectorF[8];
657       for (Int_t ivalue = 0; ivalue < 8; ivalue++)
658          posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
659
660       Float_t posG[3] = {0};
661       Float_t posL[3] = {0};
662       Int_t ichannel = 0;
663       for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
664          for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
665             tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
666             tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
667             posArray[0][ichannel] = irow;
668             posArray[1][ichannel] = ipad;
669             posArray[2][ichannel] = posL[0];
670             posArray[3][ichannel] = posL[1];
671             posArray[4][ichannel] = posG[0];
672             posArray[5][ichannel] = posG[1];
673             posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
674             posArray[7][ichannel] = ichannel;
675             
676             // loop over array containing AliTPCCalPads
677             for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
678                AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
679                AliTPCCalROC* calROC = calPad->GetCalROC(isector);
680                if (calROC)
681                   (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
682                else
683                   (vectorArray[ivalue])[ichannel] = 0;
684             }
685             ichannel++;
686          }
687       }
688       
689       cstream << "calPads" <<
690          "sector=" << isector;
691       
692       for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
693          cstream << "calPads" <<
694             (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
695             (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
696             (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
697             (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
698             (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
699          if (outlierPad) {
700             cstream << "calPads" <<
701                (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
702                (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
703                (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
704                (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
705                (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
706          }
707       }
708
709       for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
710          cstream << "calPads" <<
711             (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
712       }
713
714       if (mapFileName) {
715          for  (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
716             if (isector < 36)
717                cstream << "calPads" <<
718                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
719             else
720                cstream << "calPads" <<
721                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
722          }
723       }
724
725       cstream << "calPads" <<
726          "row.=" << &posArray[0] <<
727          "pad.=" << &posArray[1] <<
728          "lx.=" << &posArray[2] <<
729          "ly.=" << &posArray[3] <<
730          "gx.=" << &posArray[4] <<
731          "gy.=" << &posArray[5] <<
732          "rpad.=" << &posArray[6] <<
733          "channel.=" << &posArray[7];
734          
735       cstream << "calPads" <<
736          "\n";
737
738       delete[] posArray;
739       delete[] vectorArray;
740    }
741    
742
743    delete[] names;
744    if (mapFileName) {
745       delete mapIROCs;
746       delete mapOROCs;
747       delete[] mapIROCArray;
748       delete[] mapOROCArray;
749       delete[] mapNames;
750    }
751 }
752
753 Int_t AliTPCcalibDB::GetRCUTriggerConfig() const
754 {
755   //
756   // return the RCU trigger configuration register
757   //
758   TMap *map=GetRCUconfig();
759   if (!map) return -1;
760   TVectorF *v=(TVectorF*)map->GetValue("TRGCONF_TRG_MODE");
761   Float_t mode=-1;
762   for (Int_t i=0; i<v->GetNrows(); ++i){
763     Float_t newmode=v->GetMatrixArray()[i];
764     if (newmode>-1){
765       if (mode>-1&&newmode!=mode) AliWarning("Found different RCU trigger configurations!!!");
766       mode=newmode;
767     }
768   }
769   return (Int_t)mode;
770 }
771
772 Bool_t AliTPCcalibDB::IsTrgL0() 
773 {
774   //
775   // return if the FEE readout was triggered on L0
776   //
777   Int_t mode=GetRCUTriggerConfig();
778   if (mode<0) return kFALSE;
779   return ((mode&(1<<13))>>13)==1;
780 }
781
782 Bool_t AliTPCcalibDB::IsTrgL1()
783 {
784   //
785   // return if the FEE readout was triggered on L1
786   //
787   Int_t mode=GetRCUTriggerConfig();
788   if (mode<0) return kFALSE;
789   return ((mode&(1<<13))>>13)==0;
790 }
791
792 void AliTPCcalibDB::RegisterExB(Int_t index, Float_t bz, Bool_t bdelete){
793   //
794   // Register static ExB correction map
795   // index - registration index - used for visualization
796   // bz    - bz field in kGaus
797
798   //  Float_t factor =  bz/(-5.);  // default b filed in Cheb with minus sign
799   Float_t factor =  bz/(5.);  // default b filed in Cheb with minus sign
800                               // was chenged in the Revision ???? (Ruben can you add here number)
801   
802   AliMagF*   bmap = new AliMagF("MapsExB","MapsExB", factor,TMath::Sign(1.f,factor),AliMagF::k5kG);
803   
804   AliTPCExBFirst *exb  = new  AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
805   AliTPCExB::SetInstance(exb);
806   
807   if (bdelete){
808     delete bmap;
809   }else{
810     AliTPCExB::RegisterField(index,bmap);
811   }
812   if (index>=fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
813   fgExBArray.AddAt(exb,index);
814 }
815
816
817 AliTPCExB*    AliTPCcalibDB::GetExB(Float_t bz, Bool_t deleteB) {
818   //
819   // bz filed in KGaus not in tesla
820   // Get ExB correction map
821   // if doesn't exist - create it
822   //
823   Int_t index = TMath::Nint(5+bz);
824   if (index>fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
825   if (!fgExBArray.At(index)) AliTPCcalibDB::RegisterExB(index,bz,deleteB);
826   return (AliTPCExB*)fgExBArray.At(index);
827 }
828
829
830 void  AliTPCcalibDB::SetExBField(Float_t bz){
831   //
832   // Set magnetic filed for ExB correction
833   //
834   fExB = GetExB(bz,kFALSE);
835 }
836
837 void  AliTPCcalibDB::SetExBField(const AliMagF*   bmap){
838   //
839   // Set magnetic field for ExB correction
840   //
841   AliTPCExBFirst *exb  = new  AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
842   AliTPCExB::SetInstance(exb);
843   fExB=exb;
844 }
845
846
847
848
849
850 void AliTPCcalibDB::UpdateRunInformations( Int_t run, Bool_t force){
851   //
852   // - > Don't use it for reconstruction - Only for Calibration studies
853   //
854   if (run<0) return;
855   if (fRunList[run]>0 &&force==kFALSE) return;
856   AliCDBEntry * entry = 0;
857   if (run>= fRunList.GetSize()){
858     fRunList.Set(run*2+1);
859     fGRPArray.Expand(run*2+1);
860     fGRPMaps.Expand(run*2+1);
861     fGoofieArray.Expand(run*2+1);
862     fVoltageArray.Expand(run*2+1); 
863     fTemperatureArray.Expand(run*2+1);
864     fVdriftArray.Expand(run*2+1);
865     fDriftCorrectionArray.Expand(run*2+1);
866     fTimeGainSplinesArray.Expand(run*2+1);
867   }
868
869   fRunList[run]=1;  // sign as used
870
871   //
872   entry = AliCDBManager::Instance()->Get("GRP/GRP/Data",run);
873   if (entry)  {
874     AliGRPObject * grpRun = dynamic_cast<AliGRPObject*>(entry->GetObject());
875     if (!grpRun){
876       TMap* map = dynamic_cast<TMap*>(entry->GetObject());
877       if (map){
878         //grpRun = new AliGRPObject; 
879         //grpRun->ReadValuesFromMap(map);
880         grpRun =  MakeGRPObjectFromMap(map);
881
882         fGRPMaps.AddAt(map,run);
883       }
884     }
885     fGRPArray.AddAt(grpRun,run);
886   }
887   entry = AliCDBManager::Instance()->Get("TPC/Calib/Goofie",run);
888   if (entry){
889     fGoofieArray.AddAt(entry->GetObject(),run);
890   }
891   //
892   entry = AliCDBManager::Instance()->Get("TPC/Calib/HighVoltage",run);
893   if (entry)  {
894     fVoltageArray.AddAt(entry->GetObject(),run);
895   }
896   //
897   entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeGain",run);
898   if (entry)  {
899     fTimeGainSplinesArray.AddAt(entry->GetObject(),run);
900   }
901   //
902   entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeDrift",run);
903   if (entry)  {
904     fDriftCorrectionArray.AddAt(entry->GetObject(),run);
905   }
906   //
907   entry = AliCDBManager::Instance()->Get("TPC/Calib/Temperature",run);
908   if (entry)  {
909     fTemperatureArray.AddAt(entry->GetObject(),run);
910   }
911   //apply fDButil filters
912
913   fDButil->UpdateFromCalibDB();
914   if (fTemperature) fDButil->FilterTemperature(fTemperature);
915
916   AliDCSSensor * press = GetPressureSensor(run,0);
917   AliTPCSensorTempArray * temp = GetTemperatureSensor(run);
918   Bool_t accept=kTRUE;
919   if (temp) {
920     accept = fDButil->FilterTemperature(temp)>0.1;
921   }
922   if (press) {
923     const Double_t kMinP=950.;
924     const Double_t kMaxP=1050.;
925     const Double_t kMaxdP=10.;
926     const Double_t kSigmaCut=4.;
927     fDButil->FilterSensor(press,kMinP,kMaxP,kMaxdP,kSigmaCut);
928     if (press->GetFit()==0) accept=kFALSE;
929   }
930   if (press && temp &&accept){
931     AliTPCCalibVdrift * vdrift = new AliTPCCalibVdrift(temp, press,0);
932     fVdriftArray.AddAt(vdrift,run);
933   }
934   fDButil->FilterCE(120., 3., 4.,0);
935   fDButil->FilterTracks(run, 10.,0);
936 }
937
938
939 Float_t AliTPCcalibDB::GetGain(Int_t sector, Int_t row, Int_t pad){
940   //
941   //
942   AliTPCCalPad *calPad = Instance()->fDedxGainFactor;;
943   if (!calPad) return 0;
944   return calPad->GetCalROC(sector)->GetValue(row,pad);
945 }
946
947 AliSplineFit* AliTPCcalibDB::GetVdriftSplineFit(const char* name, Int_t run){
948   //
949   //
950   //
951   TObjArray *arr=GetTimeVdriftSplineRun(run);
952   if (!arr) return 0;
953   return dynamic_cast<AliSplineFit*>(arr->FindObject(name));
954 }
955
956 AliSplineFit* AliTPCcalibDB::CreateVdriftSplineFit(const char* graphName, Int_t run){
957   //
958   // create spline fit from the drift time graph in TimeDrift
959   //
960   TObjArray *arr=GetTimeVdriftSplineRun(run);
961   if (!arr) return 0;
962   TGraph *graph=dynamic_cast<TGraph*>(arr->FindObject(graphName));
963   if (!graph) return 0;
964   AliSplineFit *fit = new AliSplineFit();
965   fit->SetGraph(graph);
966   fit->SetMinPoints(graph->GetN()+1);
967   fit->InitKnots(graph,2,0,0.001);
968   fit->SplineFit(0);
969   return fit;
970 }
971
972 AliGRPObject *AliTPCcalibDB::GetGRP(Int_t run){
973   //
974   // Get GRP object for given run 
975   //
976   AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>((Instance()->fGRPArray).At(run));
977   if (!grpRun) {
978     Instance()->UpdateRunInformations(run);
979     grpRun = dynamic_cast<AliGRPObject *>(Instance()->fGRPArray.At(run));
980     if (!grpRun) return 0; 
981   }
982   return grpRun;
983 }
984
985 TMap *  AliTPCcalibDB::GetGRPMap(Int_t run){
986   //
987   //
988   //
989   TMap * grpRun = dynamic_cast<TMap *>((Instance()->fGRPMaps).At(run));
990   if (!grpRun) {
991     Instance()->UpdateRunInformations(run);
992     grpRun = dynamic_cast<TMap *>(Instance()->fGRPMaps.At(run));
993     if (!grpRun) return 0; 
994   }
995   return grpRun;
996 }
997
998
999 AliDCSSensor * AliTPCcalibDB::GetPressureSensor(Int_t run, Int_t type){
1000   //
1001   // Get Pressure sensor
1002   // run  = run number
1003   // type = 0 - Cavern pressure
1004   //        1 - Suface pressure
1005   // First try to get if trom map - if existing  (Old format of data storing)
1006   //
1007
1008
1009   TMap *map = GetGRPMap(run);  
1010   if (map){
1011     AliDCSSensor * sensor = 0;
1012     TObject *osensor=0;
1013     if (type==0) osensor = ((*map)("fCavernPressure"));
1014     if (type==1) osensor = ((*map)("fP2Pressure"));
1015     sensor =dynamic_cast<AliDCSSensor *>(osensor); 
1016     if (sensor) return sensor;
1017   }
1018   //
1019   // If not map try to get it from the GRPObject
1020   //
1021   AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run)); 
1022   if (!grpRun) {
1023     UpdateRunInformations(run);
1024     grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run));
1025     if (!grpRun) return 0; 
1026   }
1027   AliDCSSensor * sensor = grpRun->GetCavernAtmosPressure();
1028   if (type==1) sensor = grpRun->GetSurfaceAtmosPressure();
1029   return sensor; 
1030 }
1031
1032 AliTPCSensorTempArray * AliTPCcalibDB::GetTemperatureSensor(Int_t run){
1033   //
1034   // Get temperature sensor array
1035   //
1036   AliTPCSensorTempArray * tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
1037   if (!tempArray) {
1038     UpdateRunInformations(run);
1039     tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
1040   }
1041   return tempArray;
1042 }
1043
1044
1045 TObjArray * AliTPCcalibDB::GetTimeGainSplinesRun(Int_t run){
1046   //
1047   // Get temperature sensor array
1048   //
1049   TObjArray * gainSplines = (TObjArray *)fTimeGainSplinesArray.At(run);
1050   if (!gainSplines) {
1051     UpdateRunInformations(run);
1052     gainSplines = (TObjArray *)fTimeGainSplinesArray.At(run);
1053   }
1054   return gainSplines;
1055 }
1056
1057 TObjArray * AliTPCcalibDB::GetTimeVdriftSplineRun(Int_t run){
1058   //
1059   // Get drift spline array
1060   //
1061   TObjArray * driftSplines = (TObjArray *)fDriftCorrectionArray.At(run);
1062   if (!driftSplines) {
1063     UpdateRunInformations(run);
1064     driftSplines = (TObjArray *)fDriftCorrectionArray.At(run);
1065   }
1066   return driftSplines;
1067 }
1068
1069 AliDCSSensorArray * AliTPCcalibDB::GetVoltageSensors(Int_t run){
1070   //
1071   // Get temperature sensor array
1072   //
1073   AliDCSSensorArray * voltageArray = (AliDCSSensorArray *)fVoltageArray.At(run);
1074   if (!voltageArray) {
1075     UpdateRunInformations(run);
1076     voltageArray = (AliDCSSensorArray *)fVoltageArray.At(run);
1077   }
1078   return voltageArray;
1079 }
1080
1081 AliDCSSensorArray * AliTPCcalibDB::GetGoofieSensors(Int_t run){
1082   //
1083   // Get temperature sensor array
1084   //
1085   AliDCSSensorArray * goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
1086   if (!goofieArray) {
1087     UpdateRunInformations(run);
1088     goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
1089   }
1090   return goofieArray;
1091 }
1092
1093
1094
1095 AliTPCCalibVdrift *     AliTPCcalibDB::GetVdrift(Int_t run){
1096   //
1097   // Get the interface to the the vdrift 
1098   //
1099   AliTPCCalibVdrift  * vdrift = (AliTPCCalibVdrift*)fVdriftArray.At(run);
1100   if (!vdrift) {
1101     UpdateRunInformations(run);
1102     vdrift= (AliTPCCalibVdrift*)fVdriftArray.At(run);
1103   }
1104   return vdrift;
1105 }
1106
1107 Float_t AliTPCcalibDB::GetCEdriftTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1108 {
1109   //
1110   // GetCE drift time information for 'sector'
1111   // sector 72 is the mean drift time of the A-Side
1112   // sector 73 is the mean drift time of the C-Side
1113   // it timestamp==-1 return mean value
1114   //
1115   AliTPCcalibDB::Instance()->SetRun(run);
1116   TGraph *gr=AliTPCcalibDB::Instance()->GetCErocTgraph(sector);
1117   if (!gr||sector<0||sector>73) {
1118     if (entries) *entries=0;
1119     return 0.;
1120   }
1121   Float_t val=0.;
1122   if (timeStamp==-1.){
1123     val=gr->GetMean(2);
1124   }else{
1125     for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1126       Double_t x,y;
1127       gr->GetPoint(ipoint,x,y);
1128       if (x<timeStamp) continue;
1129       val=y;
1130       break;
1131     }
1132   }
1133   return val;
1134 }
1135   
1136 Float_t AliTPCcalibDB::GetCEchargeTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1137 {
1138   //
1139   // GetCE mean charge for 'sector'
1140   // it timestamp==-1 return mean value
1141   //
1142   AliTPCcalibDB::Instance()->SetRun(run);
1143   TGraph *gr=AliTPCcalibDB::Instance()->GetCErocQgraph(sector);
1144   if (!gr||sector<0||sector>71) {
1145     if (entries) *entries=0;
1146     return 0.;
1147   }
1148   Float_t val=0.;
1149   if (timeStamp==-1.){
1150     val=gr->GetMean(2);
1151   }else{
1152     for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1153       Double_t x,y;
1154       gr->GetPoint(ipoint,x,y);
1155       if (x<timeStamp) continue;
1156       val=y;
1157       break;
1158     }
1159   }
1160   return val;
1161 }
1162
1163 Float_t AliTPCcalibDB::GetDCSSensorValue(AliDCSSensorArray *arr, Int_t timeStamp, const char * sensorName, Int_t sigDigits)
1164 {
1165   //
1166   // Get Value for a DCS sensor 'sensorName', run 'run' at time 'timeStamp'
1167   //
1168   Float_t val=0;
1169   const TString sensorNameString(sensorName);
1170   AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1171   if (!sensor) return val;
1172   //use the dcs graph if possible
1173   TGraph *gr=sensor->GetGraph();
1174   if (gr){
1175     for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1176       Double_t x,y;
1177       gr->GetPoint(ipoint,x,y);
1178       Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1179       if (time<timeStamp) continue;
1180       val=y;
1181       break;
1182     }
1183     //if val is still 0, test if if the requested time if within 5min of the first/last
1184     //data point. If this is the case return the firs/last entry
1185     //the timestamps might not be syncronised for all calibration types, sometimes a 'pre'
1186     //and 'pos' period is requested. Especially to the HV this is not the case!
1187     //first point
1188     if (val==0 ){
1189       Double_t x,y;
1190       gr->GetPoint(0,x,y);
1191       Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1192       if ((time-timeStamp)<5*60) val=y;
1193     }
1194     //last point
1195     if (val==0 ){
1196       Double_t x,y;
1197       gr->GetPoint(gr->GetN()-1,x,y);
1198       Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1199       if ((timeStamp-time)<5*60) val=y;
1200     }
1201   } else {
1202     val=sensor->GetValue(timeStamp);
1203   }
1204   if (sigDigits>=0){
1205     val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1206   }
1207   return val;
1208 }
1209
1210 Float_t AliTPCcalibDB::GetDCSSensorMeanValue(AliDCSSensorArray *arr, const char * sensorName, Int_t sigDigits)
1211 {
1212   //
1213   // Get mean Value for a DCS sensor 'sensorName' during run 'run'
1214   //
1215   Float_t val=0;
1216   const TString sensorNameString(sensorName);
1217   AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1218   if (!sensor) return val;
1219
1220   //use dcs graph if it exists
1221   TGraph *gr=sensor->GetGraph();
1222   if (gr){
1223     val=gr->GetMean(2);
1224   } else {
1225     //if we don't have the dcs graph, try to get some meaningful information
1226     if (!sensor->GetFit()) return val;
1227     Int_t nKnots=sensor->GetFit()->GetKnots();
1228     Double_t tMid=(sensor->GetEndTime()-sensor->GetStartTime())/2.;
1229     for (Int_t iKnot=0;iKnot<nKnots;++iKnot){
1230       if (sensor->GetFit()->GetX()[iKnot]>tMid/3600.) break;
1231       val=(Float_t)sensor->GetFit()->GetY0()[iKnot];
1232     }
1233   }
1234   if (sigDigits>=0){
1235     val/=10;
1236     val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1237     val*=10;
1238   }
1239   return val;
1240 }
1241
1242 Float_t AliTPCcalibDB::GetChamberHighVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits) {
1243   //
1244   // return the chamber HV for given run and time: 0-35 IROC, 36-72 OROC
1245   // if timeStamp==-1 return mean value
1246   //
1247   Float_t val=0;
1248   TString sensorName="";
1249   TTimeStamp stamp(timeStamp);
1250   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1251   if (!voltageArray || (sector<0) || (sector>71)) return val;
1252   Char_t sideName='A';
1253   if ((sector/18)%2==1) sideName='C';
1254   if (sector<36){
1255     //IROC
1256     sensorName=Form("TPC_ANODE_I_%c%02d_VMEAS",sideName,sector%18);
1257   }else{
1258     //OROC
1259     sensorName=Form("TPC_ANODE_O_%c%02d_0_VMEAS",sideName,sector%18);
1260   }
1261   if (timeStamp==-1){
1262     val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1263   } else {
1264     val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1265   }
1266   return val;
1267 }
1268 Float_t AliTPCcalibDB::GetSkirtVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1269 {
1270   //
1271   // Get the skirt voltage for 'run' at 'timeStamp' and 'sector': 0-35 IROC, 36-72 OROC
1272   // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1273   // if timeStamp==-1 return the mean value for the run
1274   //
1275   Float_t val=0;
1276   TString sensorName="";
1277   TTimeStamp stamp(timeStamp);
1278   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1279   if (!voltageArray || (sector<0) || (sector>71)) return val;
1280   Char_t sideName='A';
1281   if ((sector/18)%2==1) sideName='C';
1282   sensorName=Form("TPC_SKIRT_%c_VMEAS",sideName);
1283   if (timeStamp==-1){
1284     val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1285   } else {
1286     val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1287   }
1288   return val;
1289 }
1290
1291 Float_t AliTPCcalibDB::GetCoverVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1292 {
1293   //
1294   // Get the cover voltage for run 'run' at time 'timeStamp'
1295   // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1296   // if timeStamp==-1 return the mean value for the run
1297   //
1298   Float_t val=0;
1299   TString sensorName="";
1300   TTimeStamp stamp(timeStamp);
1301   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1302   if (!voltageArray || (sector<0) || (sector>71)) return val;
1303   Char_t sideName='A';
1304   if ((sector/18)%2==1) sideName='C';
1305   if (sector<36){
1306     //IROC
1307     sensorName=Form("TPC_COVER_I_%c_VMEAS",sideName);
1308   }else{
1309     //OROC
1310     sensorName=Form("TPC_COVER_O_%c_VMEAS",sideName);
1311   }
1312   if (timeStamp==-1){
1313     val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1314   } else {
1315     val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1316   }
1317   return val;
1318 }
1319
1320 Float_t AliTPCcalibDB::GetGGoffsetVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1321 {
1322   //
1323   // Get the GG offset voltage for run 'run' at time 'timeStamp'
1324   // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1325   // if timeStamp==-1 return the mean value for the run
1326   //
1327   Float_t val=0;
1328   TString sensorName="";
1329   TTimeStamp stamp(timeStamp);
1330   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1331   if (!voltageArray || (sector<0) || (sector>71)) return val;
1332   Char_t sideName='A';
1333   if ((sector/18)%2==1) sideName='C';
1334   if (sector<36){
1335     //IROC
1336     sensorName=Form("TPC_GATE_I_%c_OFF_VMEAS",sideName);
1337   }else{
1338     //OROC
1339     sensorName=Form("TPC_GATE_O_%c_OFF_VMEAS",sideName);
1340   }
1341   if (timeStamp==-1){
1342     val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1343   } else {
1344     val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1345   }
1346   return val;
1347 }
1348
1349 Float_t AliTPCcalibDB::GetGGnegVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1350 {
1351   //
1352   // Get the GG offset voltage for run 'run' at time 'timeStamp'
1353   // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1354   // if timeStamp==-1 return the mean value for the run
1355   //
1356   Float_t val=0;
1357   TString sensorName="";
1358   TTimeStamp stamp(timeStamp);
1359   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1360   if (!voltageArray || (sector<0) || (sector>71)) return val;
1361   Char_t sideName='A';
1362   if ((sector/18)%2==1) sideName='C';
1363   if (sector<36){
1364     //IROC
1365     sensorName=Form("TPC_GATE_I_%c_NEG_VMEAS",sideName);
1366   }else{
1367     //OROC
1368     sensorName=Form("TPC_GATE_O_%c_NEG_VMEAS",sideName);
1369   }
1370   if (timeStamp==-1){
1371     val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1372   } else {
1373     val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1374   }
1375   return val;
1376 }
1377
1378 Float_t AliTPCcalibDB::GetGGposVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1379 {
1380   //
1381   // Get the GG offset voltage for run 'run' at time 'timeStamp'
1382   // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1383   // if timeStamp==-1 return the mean value for the run
1384   //
1385   Float_t val=0;
1386   TString sensorName="";
1387   TTimeStamp stamp(timeStamp);
1388   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1389   if (!voltageArray || (sector<0) || (sector>71)) return val;
1390   Char_t sideName='A';
1391   if ((sector/18)%2==1) sideName='C';
1392   if (sector<36){
1393     //IROC
1394     sensorName=Form("TPC_GATE_I_%c_POS_VMEAS",sideName);
1395   }else{
1396     //OROC
1397     sensorName=Form("TPC_GATE_O_%c_POS_VMEAS",sideName);
1398   }
1399   if (timeStamp==-1){
1400     val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1401   } else {
1402     val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1403   }
1404   return val;
1405 }
1406
1407 Float_t AliTPCcalibDB::GetPressure(Int_t timeStamp, Int_t run, Int_t type){
1408   //
1409   // GetPressure for given time stamp and runt
1410   //
1411   TTimeStamp stamp(timeStamp);
1412   AliDCSSensor * sensor = Instance()->GetPressureSensor(run,type);
1413   if (!sensor) return 0;
1414   return sensor->GetValue(stamp);
1415 }
1416
1417 Float_t AliTPCcalibDB::GetL3Current(Int_t run, Int_t statType){
1418   //
1419   // return L3 current
1420   // stat type is: AliGRPObject::Stats: kMean = 0, kTruncMean = 1, kMedian = 2, kSDMean = 3, kSDMedian = 4
1421   //
1422   Float_t current=-1;
1423   AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1424   if (grp) current=grp->GetL3Current((AliGRPObject::Stats)statType);
1425   return current;
1426 }
1427
1428 Float_t AliTPCcalibDB::GetBz(Int_t run){
1429   //
1430   // calculate BZ in T from L3 current
1431   //
1432   Float_t bz=-1;
1433   Float_t current=AliTPCcalibDB::GetL3Current(run);
1434   if (current>-1) bz=5*current/30000.*.1;
1435   return bz;
1436 }
1437
1438 Char_t  AliTPCcalibDB::GetL3Polarity(Int_t run) {
1439   //
1440   // get l3 polarity from GRP
1441   //
1442   Char_t pol=-100;
1443   AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1444   if (grp) pol=grp->GetL3Polarity();
1445   return pol;
1446 }
1447
1448 TString AliTPCcalibDB::GetRunType(Int_t run){
1449   //
1450   // return run type from grp
1451   //
1452
1453 //   TString type("UNKNOWN");
1454   AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1455   if (grp) return grp->GetRunType();
1456   return "UNKNOWN";
1457 }
1458
1459 Float_t AliTPCcalibDB::GetValueGoofie(Int_t timeStamp, Int_t run, Int_t type){
1460   //
1461   // GetPressure for given time stamp and runt
1462   //
1463   TTimeStamp stamp(timeStamp);
1464   AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(run);
1465   if (!goofieArray) return 0;
1466   AliDCSSensor *sensor = goofieArray->GetSensor(type);
1467   return sensor->GetValue(stamp);
1468 }
1469
1470
1471
1472
1473
1474
1475 Bool_t  AliTPCcalibDB::GetTemperatureFit(Int_t timeStamp, Int_t run, Int_t side,TVectorD& fit){
1476   //
1477   //
1478   //
1479   TTimeStamp tstamp(timeStamp);
1480   AliTPCSensorTempArray* tempArray  = Instance()->GetTemperatureSensor(run);
1481   if (! tempArray) return kFALSE;
1482   AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
1483   TLinearFitter * fitter = tempMap->GetLinearFitter(3,side,tstamp);
1484   if (fitter){
1485     fitter->Eval(); 
1486     fitter->GetParameters(fit);
1487   }
1488   delete fitter;
1489   delete tempMap;
1490   if (!fitter) return kFALSE;
1491   return kTRUE;
1492 }
1493
1494 Float_t AliTPCcalibDB::GetTemperature(Int_t timeStamp, Int_t run, Int_t side){
1495   //
1496   //
1497   //
1498   TVectorD vec(5);
1499   if (side==0) {
1500     GetTemperatureFit(timeStamp,run,0,vec);
1501     return vec[0];
1502   }
1503   if (side==1){
1504     GetTemperatureFit(timeStamp,run,0,vec);
1505     return vec[0];
1506   }
1507   return 0;
1508 }
1509
1510
1511 Double_t AliTPCcalibDB::GetPTRelative(UInt_t timeSec, Int_t run, Int_t side){
1512   //
1513   // Get relative P/T 
1514   // time - absolute time
1515   // run  - run number
1516   // side - 0 - A side   1-C side
1517   AliTPCCalibVdrift * vdrift =  Instance()->GetVdrift(run);
1518   if (!vdrift) return 0;
1519   return vdrift->GetPTRelative(timeSec,side);
1520 }
1521
1522 AliGRPObject * AliTPCcalibDB::MakeGRPObjectFromMap(TMap *map){
1523   //
1524   // Function to covert old GRP run information from TMap to GRPObject
1525   //
1526   //  TMap * map = AliTPCcalibDB::GetGRPMap(52406);
1527   if (!map) return 0;
1528   AliDCSSensor * sensor = 0;
1529   TObject *osensor=0;
1530   osensor = ((*map)("fP2Pressure"));
1531   sensor  =dynamic_cast<AliDCSSensor *>(osensor); 
1532   //
1533   if (!sensor) return 0;
1534   //
1535   AliDCSSensor * sensor2 = new AliDCSSensor(*sensor);
1536   osensor = ((*map)("fCavernPressure"));
1537   TGraph * gr = new TGraph(2);
1538   gr->GetX()[0]= -100000.;
1539   gr->GetX()[1]= 1000000.;
1540   gr->GetY()[0]= atof(osensor->GetName());
1541   gr->GetY()[1]= atof(osensor->GetName());
1542   sensor2->SetGraph(gr);
1543   sensor2->SetFit(0);
1544   
1545
1546   AliGRPObject *grpRun = new AliGRPObject; 
1547   grpRun->ReadValuesFromMap(map);
1548   grpRun->SetCavernAtmosPressure(sensor2);
1549   grpRun->SetSurfaceAtmosPressure(sensor);
1550   return grpRun;
1551 }
1552
1553 Bool_t AliTPCcalibDB::CreateGUITree(Int_t run, const char* filename)
1554 {
1555   //
1556   // Create a gui tree for run number 'run'
1557   //
1558
1559   if (!AliCDBManager::Instance()->GetDefaultStorage()){
1560     AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
1561                     MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
1562     return kFALSE;
1563   }
1564   //db instance
1565   AliTPCcalibDB *db=AliTPCcalibDB::Instance();
1566   // retrieve cal pad objects
1567   db->SetRun(run);
1568   db->CreateGUITree(filename);
1569 }
1570
1571 Bool_t AliTPCcalibDB::CreateGUITree(const char* filename){
1572   //
1573   //
1574   //
1575   if (!AliCDBManager::Instance()->GetDefaultStorage()){
1576     AliError("Default Storage not set. Cannot create calibration Tree!");
1577     return kFALSE;
1578   }
1579   
1580   AliTPCPreprocessorOnline prep;
1581   //noise and pedestals
1582   if (GetPedestals()) prep.AddComponent(new AliTPCCalPad(*(GetPedestals())));
1583   if (GetPadNoise() ) prep.AddComponent(new AliTPCCalPad(*(GetPadNoise())));
1584   //pulser data
1585   if (GetPulserTmean()) prep.AddComponent(new AliTPCCalPad(*(GetPulserTmean())));
1586   if (GetPulserTrms() ) prep.AddComponent(new AliTPCCalPad(*(GetPulserTrms())));
1587   if (GetPulserQmean()) prep.AddComponent(new AliTPCCalPad(*(GetPulserQmean())));
1588   //CE data
1589   if (GetCETmean()) prep.AddComponent(new AliTPCCalPad(*(GetCETmean())));
1590   if (GetCETrms() ) prep.AddComponent(new AliTPCCalPad(*(GetCETrms())));
1591   if (GetCEQmean()) prep.AddComponent(new AliTPCCalPad(*(GetCEQmean())));
1592   //Altro data
1593   if (GetALTROAcqStart() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROAcqStart() )));
1594   if (GetALTROZsThr()    ) prep.AddComponent(new AliTPCCalPad(*(GetALTROZsThr()    )));
1595   if (GetALTROFPED()     ) prep.AddComponent(new AliTPCCalPad(*(GetALTROFPED()     )));
1596   if (GetALTROAcqStop()  ) prep.AddComponent(new AliTPCCalPad(*(GetALTROAcqStop()  )));
1597   if (GetALTROMasked()   ) prep.AddComponent(new AliTPCCalPad(*(GetALTROMasked()   )));
1598   //QA
1599   AliTPCdataQA *dataQA=GetDataQA();
1600   if (dataQA) {
1601     if (dataQA->GetNLocalMaxima())
1602       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNLocalMaxima())));
1603     if (dataQA->GetMaxCharge())
1604       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetMaxCharge())));
1605     if (dataQA->GetMeanCharge())
1606       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetMeanCharge())));
1607     if (dataQA->GetNoThreshold())
1608       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNoThreshold())));
1609     if (dataQA->GetNTimeBins())
1610       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNTimeBins())));
1611     if (dataQA->GetNPads())
1612       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNPads())));
1613     if (dataQA->GetTimePosition())
1614       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetTimePosition())));
1615   }
1616   
1617   //
1618   TString file(filename);
1619   if (file.IsNull()) file=Form("guiTreeRun_%d.root",fRun);
1620   prep.DumpToFile(file.Data());
1621   return kTRUE;
1622 }
1623
1624 Bool_t AliTPCcalibDB::CreateRefFile(Int_t run, const char* filename)
1625 {
1626   //
1627   // Create a gui tree for run number 'run'
1628   //
1629   
1630   if (!AliCDBManager::Instance()->GetDefaultStorage()){
1631     AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
1632                     MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
1633     return kFALSE;
1634   }
1635   TString file(filename);
1636   if (file.IsNull()) file=Form("RefCalPads_%d.root",run);
1637   TDirectory *currDir=gDirectory;
1638   //db instance
1639   AliTPCcalibDB *db=AliTPCcalibDB::Instance();
1640   // retrieve cal pad objects
1641   db->SetRun(run);
1642   //open file
1643   TFile f(file.Data(),"recreate");
1644   //noise and pedestals
1645   db->GetPedestals()->Write("Pedestals");
1646   db->GetPadNoise()->Write("PadNoise");
1647   //pulser data
1648   db->GetPulserTmean()->Write("PulserTmean");
1649   db->GetPulserTrms()->Write("PulserTrms");
1650   db->GetPulserQmean()->Write("PulserQmean");
1651   //CE data
1652   db->GetCETmean()->Write("CETmean");
1653   db->GetCETrms()->Write("CETrms");
1654   db->GetCEQmean()->Write("CEQmean");
1655   //Altro data
1656   db->GetALTROAcqStart() ->Write("ALTROAcqStart");
1657   db->GetALTROZsThr()    ->Write("ALTROZsThr");
1658   db->GetALTROFPED()     ->Write("ALTROFPED");
1659   db->GetALTROAcqStop()  ->Write("ALTROAcqStop");
1660   db->GetALTROMasked()   ->Write("ALTROMasked");
1661   //
1662   f.Close();
1663   currDir->cd();
1664   return kTRUE;
1665 }
1666
1667
1668
1669 Double_t AliTPCcalibDB::GetVDriftCorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
1670   //
1671   // Get time dependent drift velocity correction
1672   // multiplication factor        vd = vdnom *(1+vdriftcorr)
1673   // Arguments:
1674   // mode determines the algorith how to combine the Laser Track, LaserCE and physics tracks
1675   // timestamp - timestamp
1676   // run       - run number
1677   // side      - the drift velocity per side (possible for laser and CE)
1678   //
1679   // Notice - Extrapolation outside of calibration range  - using constant function
1680   //
1681   Double_t result;
1682   // mode 1  automatic mode - according to the distance to the valid calibration
1683   //                        -  
1684   Double_t deltaP=0,  driftP=0,  wP  = 0.;
1685   Double_t deltaLT=0, driftLT=0, wLT = 0.;
1686   Double_t deltaCE=0, driftCE=0, wCE = 0.;
1687   driftP  = fDButil->GetVDriftTPC(deltaP,run,timeStamp); 
1688   driftCE = fDButil->GetVDriftTPCCE(deltaCE, run,timeStamp,36000,2);
1689   driftLT = fDButil->GetVDriftTPCLaserTracks(deltaLT,run,timeStamp,36000,2);
1690   deltaP   = TMath::Abs(deltaP);
1691   deltaLT  = TMath::Abs(deltaLT);
1692   deltaCE  = TMath::Abs(deltaCE);
1693   if (mode==1) {
1694     const Double_t kEpsilon=0.0000000001;
1695     Double_t meanDist= (deltaP+deltaLT+deltaCE)*0.3;
1696     if (meanDist<1.) return driftLT;
1697     wP  = meanDist/(deltaP +0.005*meanDist);
1698     wLT = meanDist/(deltaLT+0.005*meanDist);
1699     wCE = meanDist/(deltaCE+0.001*meanDist);
1700     if (TMath::Abs(driftCE)<kEpsilon) wCE=0;  // invalid calibration
1701     result = (driftP*wP+driftLT*wLT+driftCE*wCE)/(wP+wLT+wCE);
1702   }
1703
1704   return result;
1705 }
1706
1707 Double_t AliTPCcalibDB::GetTime0CorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
1708   //
1709   // Get time dependent time 0 (trigger delay in cm) correction
1710   // additive correction        time0 = time0+ GetTime0CorrectionTime
1711   // Value etracted combining the vdrift correction using laser tracks and CE and the physics track matchin
1712   // Arguments:
1713   // mode determines the algorith how to combine the Laser Track and physics tracks
1714   // timestamp - timestamp
1715   // run       - run number
1716   // side      - the drift velocity per side (possible for laser and CE)
1717   //
1718   // Notice - Extrapolation outside of calibration range  - using constant function
1719   //
1720   Double_t result=0;
1721   if (mode==1) result=fDButil->GetTriggerOffsetTPC(run,timeStamp);    
1722   result  *=fParam->GetZLength();
1723
1724   return result;
1725
1726 }
1727
1728
1729
1730
1731 Double_t AliTPCcalibDB::GetVDriftCorrectionGy(Int_t timeStamp, Int_t run, Int_t side, Int_t /*mode*/){
1732   //
1733   // Get global y correction drift velocity correction factor
1734   // additive factor        vd = vdnom*(1+GetVDriftCorrectionGy *gy)
1735   // Value etracted combining the vdrift correction using laser tracks and CE
1736   // Arguments:
1737   // mode determines the algorith how to combine the Laser Track, LaserCE
1738   // timestamp - timestamp
1739   // run       - run number
1740   // side      - the drift velocity gy correction per side (CE and Laser tracks)
1741   //
1742   // Notice - Extrapolation outside of calibration range  - using constant function
1743   // 
1744   if (run<=0 && fTransform) run = fTransform->GetCurrentRunNumber();
1745   UpdateRunInformations(run,kFALSE);
1746   TObjArray *array =AliTPCcalibDB::Instance()->GetTimeVdriftSplineRun(run);
1747   if (!array) return 0;
1748   TGraphErrors *laserA= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
1749   TGraphErrors *laserC= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
1750   
1751   Double_t result=0;
1752   if (laserA && laserC){
1753    result= (laserA->Eval(timeStamp)+laserC->Eval(timeStamp))*0.5;
1754   }
1755   if (laserA && side==0){
1756     result = (laserA->Eval(timeStamp));
1757   }
1758   if (laserC &&side==1){
1759     result = (laserC->Eval(timeStamp));
1760   }
1761   return -result/250.; //normalized before
1762 }
1763
1764
1765