]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibDB.cxx
Make OCDB Entry
[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);
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==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   AliCDBEntry * entry = 0;
856   if (run>= fRunList.fN){
857     fRunList.Set(run*2+1);
858     fGRPArray.Expand(run*2+1);
859     fGRPMaps.Expand(run*2+1);
860     fGoofieArray.Expand(run*2+1);
861     fVoltageArray.Expand(run*2+1); 
862     fTemperatureArray.Expand(run*2+1);
863     fVdriftArray.Expand(run*2+1);
864     fDriftCorrectionArray.Expand(run*2+1);
865     fTimeGainSplinesArray.Expand(run*2+1);
866     //
867     //
868     fALTROConfigData->Expand(run*2+1);    // ALTRO configuration data
869     fPulserData->Expand(run*2+1);         // Calibration Pulser data
870     fCEData->Expand(run*2+1);             // CE data
871     if (!fTimeGainSplines) fTimeGainSplines = new TObjArray(run*2+1);
872     fTimeGainSplines->Expand(run*2+1); // Array of AliSplineFits: at 0 MIP position in
873   }
874   if (fRunList[run]>0 &&force==kFALSE) return;
875
876   fRunList[run]=1;  // sign as used
877
878   //
879   entry = AliCDBManager::Instance()->Get("GRP/GRP/Data",run);
880   if (entry)  {
881     AliGRPObject * grpRun = dynamic_cast<AliGRPObject*>(entry->GetObject());
882     if (!grpRun){
883       TMap* map = dynamic_cast<TMap*>(entry->GetObject());
884       if (map){
885         //grpRun = new AliGRPObject; 
886         //grpRun->ReadValuesFromMap(map);
887         grpRun =  MakeGRPObjectFromMap(map);
888
889         fGRPMaps.AddAt(map,run);
890       }
891     }
892     fGRPArray.AddAt(grpRun,run);
893   }
894   entry = AliCDBManager::Instance()->Get("TPC/Calib/Goofie",run);
895   if (entry){
896     fGoofieArray.AddAt(entry->GetObject(),run);
897   }
898   //
899   entry = AliCDBManager::Instance()->Get("TPC/Calib/HighVoltage",run);
900   if (entry)  {
901     fVoltageArray.AddAt(entry->GetObject(),run);
902   }
903   //
904   entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeGain",run);
905   if (entry)  {
906     fTimeGainSplinesArray.AddAt(entry->GetObject(),run);
907   }
908   //
909   entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeDrift",run);
910   if (entry)  {
911     fDriftCorrectionArray.AddAt(entry->GetObject(),run);
912   }
913   //
914   entry = AliCDBManager::Instance()->Get("TPC/Calib/Temperature",run);
915   if (entry)  {
916     fTemperatureArray.AddAt(entry->GetObject(),run);
917   }
918   //apply fDButil filters
919
920   fDButil->UpdateFromCalibDB();
921   if (fTemperature) fDButil->FilterTemperature(fTemperature);
922
923   AliDCSSensor * press = GetPressureSensor(run,0);
924   AliTPCSensorTempArray * temp = GetTemperatureSensor(run);
925   Bool_t accept=kTRUE;
926   if (temp) {
927     accept = fDButil->FilterTemperature(temp)>0.1;
928   }
929   if (press) {
930     const Double_t kMinP=950.;
931     const Double_t kMaxP=1050.;
932     const Double_t kMaxdP=10.;
933     const Double_t kSigmaCut=4.;
934     fDButil->FilterSensor(press,kMinP,kMaxP,kMaxdP,kSigmaCut);
935     if (press->GetFit()==0) accept=kFALSE;
936   }
937   if (press && temp &&accept){
938     AliTPCCalibVdrift * vdrift = new AliTPCCalibVdrift(temp, press,0);
939     fVdriftArray.AddAt(vdrift,run);
940   }
941   fDButil->FilterCE(120., 3., 4.,0);
942   fDButil->FilterTracks(run, 10.,0);
943 }
944
945
946 Float_t AliTPCcalibDB::GetGain(Int_t sector, Int_t row, Int_t pad){
947   //
948   //
949   AliTPCCalPad *calPad = Instance()->fDedxGainFactor;;
950   if (!calPad) return 0;
951   return calPad->GetCalROC(sector)->GetValue(row,pad);
952 }
953
954 AliSplineFit* AliTPCcalibDB::GetVdriftSplineFit(const char* name, Int_t run){
955   //
956   //
957   //
958   TObjArray *arr=GetTimeVdriftSplineRun(run);
959   if (!arr) return 0;
960   return dynamic_cast<AliSplineFit*>(arr->FindObject(name));
961 }
962
963 AliSplineFit* AliTPCcalibDB::CreateVdriftSplineFit(const char* graphName, Int_t run){
964   //
965   // create spline fit from the drift time graph in TimeDrift
966   //
967   TObjArray *arr=GetTimeVdriftSplineRun(run);
968   if (!arr) return 0;
969   TGraph *graph=dynamic_cast<TGraph*>(arr->FindObject(graphName));
970   if (!graph) return 0;
971   AliSplineFit *fit = new AliSplineFit();
972   fit->SetGraph(graph);
973   fit->SetMinPoints(graph->GetN()+1);
974   fit->InitKnots(graph,2,0,0.001);
975   fit->SplineFit(0);
976   return fit;
977 }
978
979 AliGRPObject *AliTPCcalibDB::GetGRP(Int_t run){
980   //
981   // Get GRP object for given run 
982   //
983   AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>((Instance()->fGRPArray).At(run));
984   if (!grpRun) {
985     Instance()->UpdateRunInformations(run);
986     grpRun = dynamic_cast<AliGRPObject *>(Instance()->fGRPArray.At(run));
987     if (!grpRun) return 0; 
988   }
989   return grpRun;
990 }
991
992 TMap *  AliTPCcalibDB::GetGRPMap(Int_t run){
993   //
994   //
995   //
996   TMap * grpRun = dynamic_cast<TMap *>((Instance()->fGRPMaps).At(run));
997   if (!grpRun) {
998     Instance()->UpdateRunInformations(run);
999     grpRun = dynamic_cast<TMap *>(Instance()->fGRPMaps.At(run));
1000     if (!grpRun) return 0; 
1001   }
1002   return grpRun;
1003 }
1004
1005
1006 AliDCSSensor * AliTPCcalibDB::GetPressureSensor(Int_t run, Int_t type){
1007   //
1008   // Get Pressure sensor
1009   // run  = run number
1010   // type = 0 - Cavern pressure
1011   //        1 - Suface pressure
1012   // First try to get if trom map - if existing  (Old format of data storing)
1013   //
1014
1015
1016   TMap *map = GetGRPMap(run);  
1017   if (map){
1018     AliDCSSensor * sensor = 0;
1019     TObject *osensor=0;
1020     if (type==0) osensor = ((*map)("fCavernPressure"));
1021     if (type==1) osensor = ((*map)("fP2Pressure"));
1022     sensor =dynamic_cast<AliDCSSensor *>(osensor); 
1023     if (sensor) return sensor;
1024   }
1025   //
1026   // If not map try to get it from the GRPObject
1027   //
1028   AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run)); 
1029   if (!grpRun) {
1030     UpdateRunInformations(run);
1031     grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run));
1032     if (!grpRun) return 0; 
1033   }
1034   AliDCSSensor * sensor = grpRun->GetCavernAtmosPressure();
1035   if (type==1) sensor = grpRun->GetSurfaceAtmosPressure();
1036   return sensor; 
1037 }
1038
1039 AliTPCSensorTempArray * AliTPCcalibDB::GetTemperatureSensor(Int_t run){
1040   //
1041   // Get temperature sensor array
1042   //
1043   AliTPCSensorTempArray * tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
1044   if (!tempArray) {
1045     UpdateRunInformations(run);
1046     tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
1047   }
1048   return tempArray;
1049 }
1050
1051
1052 TObjArray * AliTPCcalibDB::GetTimeGainSplinesRun(Int_t run){
1053   //
1054   // Get temperature sensor array
1055   //
1056   TObjArray * gainSplines = (TObjArray *)fTimeGainSplinesArray.At(run);
1057   if (!gainSplines) {
1058     UpdateRunInformations(run);
1059     gainSplines = (TObjArray *)fTimeGainSplinesArray.At(run);
1060   }
1061   return gainSplines;
1062 }
1063
1064 TObjArray * AliTPCcalibDB::GetTimeVdriftSplineRun(Int_t run){
1065   //
1066   // Get drift spline array
1067   //
1068   TObjArray * driftSplines = (TObjArray *)fDriftCorrectionArray.At(run);
1069   if (!driftSplines) {
1070     UpdateRunInformations(run);
1071     driftSplines = (TObjArray *)fDriftCorrectionArray.At(run);
1072   }
1073   return driftSplines;
1074 }
1075
1076 AliDCSSensorArray * AliTPCcalibDB::GetVoltageSensors(Int_t run){
1077   //
1078   // Get temperature sensor array
1079   //
1080   AliDCSSensorArray * voltageArray = (AliDCSSensorArray *)fVoltageArray.At(run);
1081   if (!voltageArray) {
1082     UpdateRunInformations(run);
1083     voltageArray = (AliDCSSensorArray *)fVoltageArray.At(run);
1084   }
1085   return voltageArray;
1086 }
1087
1088 AliDCSSensorArray * AliTPCcalibDB::GetGoofieSensors(Int_t run){
1089   //
1090   // Get temperature sensor array
1091   //
1092   AliDCSSensorArray * goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
1093   if (!goofieArray) {
1094     UpdateRunInformations(run);
1095     goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
1096   }
1097   return goofieArray;
1098 }
1099
1100
1101
1102 AliTPCCalibVdrift *     AliTPCcalibDB::GetVdrift(Int_t run){
1103   //
1104   // Get the interface to the the vdrift 
1105   //
1106   AliTPCCalibVdrift  * vdrift = (AliTPCCalibVdrift*)fVdriftArray.At(run);
1107   if (!vdrift) {
1108     UpdateRunInformations(run);
1109     vdrift= (AliTPCCalibVdrift*)fVdriftArray.At(run);
1110   }
1111   return vdrift;
1112 }
1113
1114 Float_t AliTPCcalibDB::GetCEdriftTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1115 {
1116   //
1117   // GetCE drift time information for 'sector'
1118   // sector 72 is the mean drift time of the A-Side
1119   // sector 73 is the mean drift time of the C-Side
1120   // it timestamp==-1 return mean value
1121   //
1122   AliTPCcalibDB::Instance()->SetRun(run);
1123   TGraph *gr=AliTPCcalibDB::Instance()->GetCErocTgraph(sector);
1124   if (!gr||sector<0||sector>73) {
1125     if (entries) *entries=0;
1126     return 0.;
1127   }
1128   Float_t val=0.;
1129   if (timeStamp==-1.){
1130     val=gr->GetMean(2);
1131   }else{
1132     for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1133       Double_t x,y;
1134       gr->GetPoint(ipoint,x,y);
1135       if (x<timeStamp) continue;
1136       val=y;
1137       break;
1138     }
1139   }
1140   return val;
1141 }
1142   
1143 Float_t AliTPCcalibDB::GetCEchargeTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1144 {
1145   //
1146   // GetCE mean charge for 'sector'
1147   // it timestamp==-1 return mean value
1148   //
1149   AliTPCcalibDB::Instance()->SetRun(run);
1150   TGraph *gr=AliTPCcalibDB::Instance()->GetCErocQgraph(sector);
1151   if (!gr||sector<0||sector>71) {
1152     if (entries) *entries=0;
1153     return 0.;
1154   }
1155   Float_t val=0.;
1156   if (timeStamp==-1.){
1157     val=gr->GetMean(2);
1158   }else{
1159     for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1160       Double_t x,y;
1161       gr->GetPoint(ipoint,x,y);
1162       if (x<timeStamp) continue;
1163       val=y;
1164       break;
1165     }
1166   }
1167   return val;
1168 }
1169
1170 Float_t AliTPCcalibDB::GetDCSSensorValue(AliDCSSensorArray *arr, Int_t timeStamp, const char * sensorName, Int_t sigDigits)
1171 {
1172   //
1173   // Get Value for a DCS sensor 'sensorName', run 'run' at time 'timeStamp'
1174   //
1175   Float_t val=0;
1176   const TString sensorNameString(sensorName);
1177   AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1178   if (!sensor) return val;
1179   //use the dcs graph if possible
1180   TGraph *gr=sensor->GetGraph();
1181   if (gr){
1182     for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1183       Double_t x,y;
1184       gr->GetPoint(ipoint,x,y);
1185       Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1186       if (time<timeStamp) continue;
1187       val=y;
1188       break;
1189     }
1190     //if val is still 0, test if if the requested time if within 5min of the first/last
1191     //data point. If this is the case return the firs/last entry
1192     //the timestamps might not be syncronised for all calibration types, sometimes a 'pre'
1193     //and 'pos' period is requested. Especially to the HV this is not the case!
1194     //first point
1195     if (val==0 ){
1196       Double_t x,y;
1197       gr->GetPoint(0,x,y);
1198       Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1199       if ((time-timeStamp)<5*60) val=y;
1200     }
1201     //last point
1202     if (val==0 ){
1203       Double_t x,y;
1204       gr->GetPoint(gr->GetN()-1,x,y);
1205       Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1206       if ((timeStamp-time)<5*60) val=y;
1207     }
1208   } else {
1209     val=sensor->GetValue(timeStamp);
1210   }
1211   if (sigDigits>=0){
1212     val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1213   }
1214   return val;
1215 }
1216
1217 Float_t AliTPCcalibDB::GetDCSSensorMeanValue(AliDCSSensorArray *arr, const char * sensorName, Int_t sigDigits)
1218 {
1219   //
1220   // Get mean Value for a DCS sensor 'sensorName' during run 'run'
1221   //
1222   Float_t val=0;
1223   const TString sensorNameString(sensorName);
1224   AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1225   if (!sensor) return val;
1226
1227   //use dcs graph if it exists
1228   TGraph *gr=sensor->GetGraph();
1229   if (gr){
1230     val=gr->GetMean(2);
1231   } else {
1232     //if we don't have the dcs graph, try to get some meaningful information
1233     if (!sensor->GetFit()) return val;
1234     Int_t nKnots=sensor->GetFit()->GetKnots();
1235     Double_t tMid=(sensor->GetEndTime()-sensor->GetStartTime())/2.;
1236     for (Int_t iKnot=0;iKnot<nKnots;++iKnot){
1237       if (sensor->GetFit()->GetX()[iKnot]>tMid/3600.) break;
1238       val=(Float_t)sensor->GetFit()->GetY0()[iKnot];
1239     }
1240   }
1241   if (sigDigits>=0){
1242     val/=10;
1243     val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1244     val*=10;
1245   }
1246   return val;
1247 }
1248
1249 Float_t AliTPCcalibDB::GetChamberHighVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits) {
1250   //
1251   // return the chamber HV for given run and time: 0-35 IROC, 36-72 OROC
1252   // if timeStamp==-1 return mean value
1253   //
1254   Float_t val=0;
1255   TString sensorName="";
1256   TTimeStamp stamp(timeStamp);
1257   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1258   if (!voltageArray || (sector<0) || (sector>71)) return val;
1259   Char_t sideName='A';
1260   if ((sector/18)%2==1) sideName='C';
1261   if (sector<36){
1262     //IROC
1263     sensorName=Form("TPC_ANODE_I_%c%02d_VMEAS",sideName,sector%18);
1264   }else{
1265     //OROC
1266     sensorName=Form("TPC_ANODE_O_%c%02d_0_VMEAS",sideName,sector%18);
1267   }
1268   if (timeStamp==-1){
1269     val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1270   } else {
1271     val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1272   }
1273   return val;
1274 }
1275 Float_t AliTPCcalibDB::GetSkirtVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1276 {
1277   //
1278   // Get the skirt voltage for 'run' at 'timeStamp' and 'sector': 0-35 IROC, 36-72 OROC
1279   // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1280   // if timeStamp==-1 return the mean value for the run
1281   //
1282   Float_t val=0;
1283   TString sensorName="";
1284   TTimeStamp stamp(timeStamp);
1285   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1286   if (!voltageArray || (sector<0) || (sector>71)) return val;
1287   Char_t sideName='A';
1288   if ((sector/18)%2==1) sideName='C';
1289   sensorName=Form("TPC_SKIRT_%c_VMEAS",sideName);
1290   if (timeStamp==-1){
1291     val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1292   } else {
1293     val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1294   }
1295   return val;
1296 }
1297
1298 Float_t AliTPCcalibDB::GetCoverVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1299 {
1300   //
1301   // Get the cover voltage for run 'run' at time 'timeStamp'
1302   // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1303   // if timeStamp==-1 return the mean value for the run
1304   //
1305   Float_t val=0;
1306   TString sensorName="";
1307   TTimeStamp stamp(timeStamp);
1308   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1309   if (!voltageArray || (sector<0) || (sector>71)) return val;
1310   Char_t sideName='A';
1311   if ((sector/18)%2==1) sideName='C';
1312   if (sector<36){
1313     //IROC
1314     sensorName=Form("TPC_COVER_I_%c_VMEAS",sideName);
1315   }else{
1316     //OROC
1317     sensorName=Form("TPC_COVER_O_%c_VMEAS",sideName);
1318   }
1319   if (timeStamp==-1){
1320     val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1321   } else {
1322     val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1323   }
1324   return val;
1325 }
1326
1327 Float_t AliTPCcalibDB::GetGGoffsetVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1328 {
1329   //
1330   // Get the GG offset voltage for run 'run' at time 'timeStamp'
1331   // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1332   // if timeStamp==-1 return the mean value for the run
1333   //
1334   Float_t val=0;
1335   TString sensorName="";
1336   TTimeStamp stamp(timeStamp);
1337   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1338   if (!voltageArray || (sector<0) || (sector>71)) return val;
1339   Char_t sideName='A';
1340   if ((sector/18)%2==1) sideName='C';
1341   if (sector<36){
1342     //IROC
1343     sensorName=Form("TPC_GATE_I_%c_OFF_VMEAS",sideName);
1344   }else{
1345     //OROC
1346     sensorName=Form("TPC_GATE_O_%c_OFF_VMEAS",sideName);
1347   }
1348   if (timeStamp==-1){
1349     val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1350   } else {
1351     val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1352   }
1353   return val;
1354 }
1355
1356 Float_t AliTPCcalibDB::GetGGnegVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1357 {
1358   //
1359   // Get the GG offset voltage for run 'run' at time 'timeStamp'
1360   // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1361   // if timeStamp==-1 return the mean value for the run
1362   //
1363   Float_t val=0;
1364   TString sensorName="";
1365   TTimeStamp stamp(timeStamp);
1366   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1367   if (!voltageArray || (sector<0) || (sector>71)) return val;
1368   Char_t sideName='A';
1369   if ((sector/18)%2==1) sideName='C';
1370   if (sector<36){
1371     //IROC
1372     sensorName=Form("TPC_GATE_I_%c_NEG_VMEAS",sideName);
1373   }else{
1374     //OROC
1375     sensorName=Form("TPC_GATE_O_%c_NEG_VMEAS",sideName);
1376   }
1377   if (timeStamp==-1){
1378     val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1379   } else {
1380     val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1381   }
1382   return val;
1383 }
1384
1385 Float_t AliTPCcalibDB::GetGGposVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1386 {
1387   //
1388   // Get the GG offset voltage for run 'run' at time 'timeStamp'
1389   // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1390   // if timeStamp==-1 return the mean value for the run
1391   //
1392   Float_t val=0;
1393   TString sensorName="";
1394   TTimeStamp stamp(timeStamp);
1395   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1396   if (!voltageArray || (sector<0) || (sector>71)) return val;
1397   Char_t sideName='A';
1398   if ((sector/18)%2==1) sideName='C';
1399   if (sector<36){
1400     //IROC
1401     sensorName=Form("TPC_GATE_I_%c_POS_VMEAS",sideName);
1402   }else{
1403     //OROC
1404     sensorName=Form("TPC_GATE_O_%c_POS_VMEAS",sideName);
1405   }
1406   if (timeStamp==-1){
1407     val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1408   } else {
1409     val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1410   }
1411   return val;
1412 }
1413
1414 Float_t AliTPCcalibDB::GetPressure(Int_t timeStamp, Int_t run, Int_t type){
1415   //
1416   // GetPressure for given time stamp and runt
1417   //
1418   TTimeStamp stamp(timeStamp);
1419   AliDCSSensor * sensor = Instance()->GetPressureSensor(run,type);
1420   if (!sensor) return 0;
1421   return sensor->GetValue(stamp);
1422 }
1423
1424 Float_t AliTPCcalibDB::GetL3Current(Int_t run, Int_t statType){
1425   //
1426   // return L3 current
1427   // stat type is: AliGRPObject::Stats: kMean = 0, kTruncMean = 1, kMedian = 2, kSDMean = 3, kSDMedian = 4
1428   //
1429   Float_t current=-1;
1430   AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1431   if (grp) current=grp->GetL3Current((AliGRPObject::Stats)statType);
1432   return current;
1433 }
1434
1435 Float_t AliTPCcalibDB::GetBz(Int_t run){
1436   //
1437   // calculate BZ in T from L3 current
1438   //
1439   Float_t bz=-1;
1440   Float_t current=AliTPCcalibDB::GetL3Current(run);
1441   if (current>-1) bz=5*current/30000.*.1;
1442   return bz;
1443 }
1444
1445 Char_t  AliTPCcalibDB::GetL3Polarity(Int_t run) {
1446   //
1447   // get l3 polarity from GRP
1448   //
1449   Char_t pol=-100;
1450   AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1451   if (grp) pol=grp->GetL3Polarity();
1452   return pol;
1453 }
1454
1455 TString AliTPCcalibDB::GetRunType(Int_t run){
1456   //
1457   // return run type from grp
1458   //
1459
1460 //   TString type("UNKNOWN");
1461   AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1462   if (grp) return grp->GetRunType();
1463   return "UNKNOWN";
1464 }
1465
1466 Float_t AliTPCcalibDB::GetValueGoofie(Int_t timeStamp, Int_t run, Int_t type){
1467   //
1468   // GetPressure for given time stamp and runt
1469   //
1470   TTimeStamp stamp(timeStamp);
1471   AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(run);
1472   if (!goofieArray) return 0;
1473   AliDCSSensor *sensor = goofieArray->GetSensor(type);
1474   return sensor->GetValue(stamp);
1475 }
1476
1477
1478
1479
1480
1481
1482 Bool_t  AliTPCcalibDB::GetTemperatureFit(Int_t timeStamp, Int_t run, Int_t side,TVectorD& fit){
1483   //
1484   //
1485   //
1486   TTimeStamp tstamp(timeStamp);
1487   AliTPCSensorTempArray* tempArray  = Instance()->GetTemperatureSensor(run);
1488   if (! tempArray) return kFALSE;
1489   AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
1490   TLinearFitter * fitter = tempMap->GetLinearFitter(3,side,tstamp);
1491   if (fitter){
1492     fitter->Eval(); 
1493     fitter->GetParameters(fit);
1494   }
1495   delete fitter;
1496   delete tempMap;
1497   if (!fitter) return kFALSE;
1498   return kTRUE;
1499 }
1500
1501 Float_t AliTPCcalibDB::GetTemperature(Int_t timeStamp, Int_t run, Int_t side){
1502   //
1503   //
1504   //
1505   TVectorD vec(5);
1506   if (side==0) {
1507     GetTemperatureFit(timeStamp,run,0,vec);
1508     return vec[0];
1509   }
1510   if (side==1){
1511     GetTemperatureFit(timeStamp,run,0,vec);
1512     return vec[0];
1513   }
1514   return 0;
1515 }
1516
1517
1518 Double_t AliTPCcalibDB::GetPTRelative(UInt_t timeSec, Int_t run, Int_t side){
1519   //
1520   // Get relative P/T 
1521   // time - absolute time
1522   // run  - run number
1523   // side - 0 - A side   1-C side
1524   AliTPCCalibVdrift * vdrift =  Instance()->GetVdrift(run);
1525   if (!vdrift) return 0;
1526   return vdrift->GetPTRelative(timeSec,side);
1527 }
1528
1529 AliGRPObject * AliTPCcalibDB::MakeGRPObjectFromMap(TMap *map){
1530   //
1531   // Function to covert old GRP run information from TMap to GRPObject
1532   //
1533   //  TMap * map = AliTPCcalibDB::GetGRPMap(52406);
1534   if (!map) return 0;
1535   AliDCSSensor * sensor = 0;
1536   TObject *osensor=0;
1537   osensor = ((*map)("fP2Pressure"));
1538   sensor  =dynamic_cast<AliDCSSensor *>(osensor); 
1539   //
1540   if (!sensor) return 0;
1541   //
1542   AliDCSSensor * sensor2 = new AliDCSSensor(*sensor);
1543   osensor = ((*map)("fCavernPressure"));
1544   TGraph * gr = new TGraph(2);
1545   gr->GetX()[0]= -100000.;
1546   gr->GetX()[1]= 1000000.;
1547   gr->GetY()[0]= atof(osensor->GetName());
1548   gr->GetY()[1]= atof(osensor->GetName());
1549   sensor2->SetGraph(gr);
1550   sensor2->SetFit(0);
1551   
1552
1553   AliGRPObject *grpRun = new AliGRPObject; 
1554   grpRun->ReadValuesFromMap(map);
1555   grpRun->SetCavernAtmosPressure(sensor2);
1556   grpRun->SetSurfaceAtmosPressure(sensor);
1557   return grpRun;
1558 }
1559
1560 Bool_t AliTPCcalibDB::CreateGUITree(Int_t run, const char* filename)
1561 {
1562   //
1563   // Create a gui tree for run number 'run'
1564   //
1565
1566   if (!AliCDBManager::Instance()->GetDefaultStorage()){
1567     AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
1568                     MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
1569     return kFALSE;
1570   }
1571   //db instance
1572   AliTPCcalibDB *db=AliTPCcalibDB::Instance();
1573   // retrieve cal pad objects
1574   db->SetRun(run);
1575   db->CreateGUITree(filename);
1576   return kTRUE;
1577 }
1578
1579 Bool_t AliTPCcalibDB::CreateGUITree(const char* filename){
1580   //
1581   //
1582   //
1583   if (!AliCDBManager::Instance()->GetDefaultStorage()){
1584     AliError("Default Storage not set. Cannot create calibration Tree!");
1585     return kFALSE;
1586   }
1587   
1588   AliTPCPreprocessorOnline prep;
1589   //noise and pedestals
1590   if (GetPedestals()) prep.AddComponent(new AliTPCCalPad(*(GetPedestals())));
1591   if (GetPadNoise() ) prep.AddComponent(new AliTPCCalPad(*(GetPadNoise())));
1592   //pulser data
1593   if (GetPulserTmean()) prep.AddComponent(new AliTPCCalPad(*(GetPulserTmean())));
1594   if (GetPulserTrms() ) prep.AddComponent(new AliTPCCalPad(*(GetPulserTrms())));
1595   if (GetPulserQmean()) prep.AddComponent(new AliTPCCalPad(*(GetPulserQmean())));
1596   //CE data
1597   if (GetCETmean()) prep.AddComponent(new AliTPCCalPad(*(GetCETmean())));
1598   if (GetCETrms() ) prep.AddComponent(new AliTPCCalPad(*(GetCETrms())));
1599   if (GetCEQmean()) prep.AddComponent(new AliTPCCalPad(*(GetCEQmean())));
1600   //Altro data
1601   if (GetALTROAcqStart() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROAcqStart() )));
1602   if (GetALTROZsThr()    ) prep.AddComponent(new AliTPCCalPad(*(GetALTROZsThr()    )));
1603   if (GetALTROFPED()     ) prep.AddComponent(new AliTPCCalPad(*(GetALTROFPED()     )));
1604   if (GetALTROAcqStop()  ) prep.AddComponent(new AliTPCCalPad(*(GetALTROAcqStop()  )));
1605   if (GetALTROMasked()   ) prep.AddComponent(new AliTPCCalPad(*(GetALTROMasked()   )));
1606   //QA
1607   AliTPCdataQA *dataQA=GetDataQA();
1608   if (dataQA) {
1609     if (dataQA->GetNLocalMaxima())
1610       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNLocalMaxima())));
1611     if (dataQA->GetMaxCharge())
1612       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetMaxCharge())));
1613     if (dataQA->GetMeanCharge())
1614       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetMeanCharge())));
1615     if (dataQA->GetNoThreshold())
1616       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNoThreshold())));
1617     if (dataQA->GetNTimeBins())
1618       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNTimeBins())));
1619     if (dataQA->GetNPads())
1620       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNPads())));
1621     if (dataQA->GetTimePosition())
1622       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetTimePosition())));
1623   }
1624   
1625   //
1626   TString file(filename);
1627   if (file.IsNull()) file=Form("guiTreeRun_%d.root",fRun);
1628   prep.DumpToFile(file.Data());
1629   return kTRUE;
1630 }
1631
1632 Bool_t AliTPCcalibDB::CreateRefFile(Int_t run, const char* filename)
1633 {
1634   //
1635   // Create a gui tree for run number 'run'
1636   //
1637   
1638   if (!AliCDBManager::Instance()->GetDefaultStorage()){
1639     AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
1640                     MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
1641     return kFALSE;
1642   }
1643   TString file(filename);
1644   if (file.IsNull()) file=Form("RefCalPads_%d.root",run);
1645   TDirectory *currDir=gDirectory;
1646   //db instance
1647   AliTPCcalibDB *db=AliTPCcalibDB::Instance();
1648   // retrieve cal pad objects
1649   db->SetRun(run);
1650   //open file
1651   TFile f(file.Data(),"recreate");
1652   //noise and pedestals
1653   db->GetPedestals()->Write("Pedestals");
1654   db->GetPadNoise()->Write("PadNoise");
1655   //pulser data
1656   db->GetPulserTmean()->Write("PulserTmean");
1657   db->GetPulserTrms()->Write("PulserTrms");
1658   db->GetPulserQmean()->Write("PulserQmean");
1659   //CE data
1660   db->GetCETmean()->Write("CETmean");
1661   db->GetCETrms()->Write("CETrms");
1662   db->GetCEQmean()->Write("CEQmean");
1663   //Altro data
1664   db->GetALTROAcqStart() ->Write("ALTROAcqStart");
1665   db->GetALTROZsThr()    ->Write("ALTROZsThr");
1666   db->GetALTROFPED()     ->Write("ALTROFPED");
1667   db->GetALTROAcqStop()  ->Write("ALTROAcqStop");
1668   db->GetALTROMasked()   ->Write("ALTROMasked");
1669   //
1670   f.Close();
1671   currDir->cd();
1672   return kTRUE;
1673 }
1674
1675
1676
1677 Double_t AliTPCcalibDB::GetVDriftCorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
1678   //
1679   // Get time dependent drift velocity correction
1680   // multiplication factor        vd = vdnom *(1+vdriftcorr)
1681   // Arguments:
1682   // mode determines the algorith how to combine the Laser Track, LaserCE and physics tracks
1683   // timestamp - timestamp
1684   // run       - run number
1685   // side      - the drift velocity per side (possible for laser and CE)
1686   //
1687   // Notice - Extrapolation outside of calibration range  - using constant function
1688   //
1689   Double_t result;
1690   // mode 1  automatic mode - according to the distance to the valid calibration
1691   //                        -  
1692   Double_t deltaP=0,  driftP=0,      wP  = 0.;
1693   Double_t deltaITS=0,driftITS=0,    wITS= 0.;
1694   Double_t deltaLT=0, driftLT=0,     wLT = 0.;
1695   Double_t deltaCE=0, driftCE=0,     wCE = 0.;
1696   driftP  = fDButil->GetVDriftTPC(deltaP,run,timeStamp); 
1697   driftITS= fDButil->GetVDriftTPCITS(deltaITS,run,timeStamp);
1698   driftCE = fDButil->GetVDriftTPCCE(deltaCE, run,timeStamp,36000,2);
1699   driftLT = fDButil->GetVDriftTPCLaserTracks(deltaLT,run,timeStamp,36000,2);
1700   deltaITS = TMath::Abs(deltaITS);
1701   deltaP   = TMath::Abs(deltaP);
1702   deltaLT  = TMath::Abs(deltaLT);
1703   deltaCE  = TMath::Abs(deltaCE);
1704   if (mode==1) {
1705     const Double_t kEpsilon=0.00000000001;
1706     const Double_t kdeltaT=360.; // 10 minutes
1707     wITS  = 64.*kdeltaT/(deltaITS +kdeltaT);
1708     wLT   = 16.*kdeltaT/(deltaLT  +kdeltaT);
1709     wP    = 0. *kdeltaT/(deltaP   +kdeltaT);
1710     wCE   = 1. *kdeltaT/(deltaCE  +kdeltaT);
1711     //
1712     //
1713     if (TMath::Abs(driftP)<kEpsilon)  wP=0;  // invalid calibration
1714     if (TMath::Abs(driftITS)<kEpsilon)wITS=0;  // invalid calibration
1715     if (TMath::Abs(driftLT)<kEpsilon) wLT=0;  // invalid calibration
1716     if (TMath::Abs(driftCE)<kEpsilon) wCE=0;  // invalid calibration
1717     if (wP+wITS+wLT+wCE<kEpsilon) return 0;
1718     result = (driftP*wP+driftITS*wITS+driftLT*wLT+driftCE*wCE)/(wP+wITS+wLT+wCE);
1719   }
1720
1721   return result;
1722 }
1723
1724 Double_t AliTPCcalibDB::GetTime0CorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
1725   //
1726   // Get time dependent time 0 (trigger delay in cm) correction
1727   // additive correction        time0 = time0+ GetTime0CorrectionTime
1728   // Value etracted combining the vdrift correction using laser tracks and CE and the physics track matchin
1729   // Arguments:
1730   // mode determines the algorith how to combine the Laser Track and physics tracks
1731   // timestamp - timestamp
1732   // run       - run number
1733   // side      - the drift velocity per side (possible for laser and CE)
1734   //
1735   // Notice - Extrapolation outside of calibration range  - using constant function
1736   //
1737   Double_t result=0;
1738   if (mode==2) {
1739     // TPC-TPC mode
1740     result=fDButil->GetTriggerOffsetTPC(run,timeStamp);    
1741     result  *=fParam->GetZLength();
1742   }
1743   if (mode==1){
1744     // TPC-ITS mode
1745     Double_t dist=0;
1746     result= -fDButil->GetTime0TPCITS(dist, run, timeStamp)*fParam->GetDriftV()/1000000.;
1747   }
1748   return result;
1749
1750 }
1751
1752
1753
1754
1755 Double_t AliTPCcalibDB::GetVDriftCorrectionGy(Int_t timeStamp, Int_t run, Int_t side, Int_t /*mode*/){
1756   //
1757   // Get global y correction drift velocity correction factor
1758   // additive factor        vd = vdnom*(1+GetVDriftCorrectionGy *gy)
1759   // Value etracted combining the vdrift correction using laser tracks and CE
1760   // Arguments:
1761   // mode determines the algorith how to combine the Laser Track, LaserCE
1762   // timestamp - timestamp
1763   // run       - run number
1764   // side      - the drift velocity gy correction per side (CE and Laser tracks)
1765   //
1766   // Notice - Extrapolation outside of calibration range  - using constant function
1767   // 
1768   if (run<=0 && fTransform) run = fTransform->GetCurrentRunNumber();
1769   UpdateRunInformations(run,kFALSE);
1770   TObjArray *array =AliTPCcalibDB::Instance()->GetTimeVdriftSplineRun(run);
1771   if (!array) return 0;
1772   TGraphErrors *laserA= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
1773   TGraphErrors *laserC= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
1774   
1775   Double_t result=0;
1776   if (laserA && laserC){
1777    result= (laserA->Eval(timeStamp)+laserC->Eval(timeStamp))*0.5;
1778   }
1779   if (laserA && side==0){
1780     result = (laserA->Eval(timeStamp));
1781   }
1782   if (laserC &&side==1){
1783     result = (laserC->Eval(timeStamp));
1784   }
1785   return -result/250.; //normalized before
1786 }
1787
1788
1789