]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibDB.cxx
AliTPCcalibBase - use of the current event, track and seed innthe base class
[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 void AliTPCcalibDB::UpdateNonRec(){
450   //
451   // Update/Load the parameters which are important for QA studies
452   // and not used yet for the reconstruction
453   //
454    //RAW calibration data
455   AliCDBEntry * entry=0;
456   entry          = GetCDBEntry("TPC/Calib/Raw");
457   if (entry){
458     entry->SetOwner(kTRUE);
459     TObjArray *arr=(TObjArray*)(entry->GetObject());
460     if (arr) fCalibRaw=(AliTPCCalibRaw*)arr->At(0);
461   }
462   //QA calibration data
463   entry          = GetCDBEntry("TPC/Calib/QA");
464   if (entry){
465     entry->SetOwner(kTRUE);
466     fDataQA=dynamic_cast<AliTPCdataQA*>(entry->GetObject());
467   }
468   
469 }
470
471
472
473 void AliTPCcalibDB::CreateObjectList(const Char_t *filename, TObjArray *calibObjects)
474 {
475 //
476 // Create calibration objects and read contents from OCDB
477 //
478    if ( calibObjects == 0x0 ) return;
479    ifstream in;
480    in.open(filename);
481    if ( !in.is_open() ){
482       fprintf(stderr,"Error: cannot open list file '%s'", filename);
483       return;
484    }
485    
486    AliTPCCalPad *calPad=0x0;
487    
488    TString sFile;
489    sFile.ReadFile(in);
490    in.close();
491    
492    TObjArray *arrFileLine = sFile.Tokenize("\n");
493    
494    TIter nextLine(arrFileLine);
495    
496    TObjString *sObjLine=0x0;
497    while ( (sObjLine = (TObjString*)nextLine()) ){
498       TString sLine(sObjLine->GetString());
499       
500       TObjArray *arrNextCol = sLine.Tokenize("\t");
501       
502       TObjString *sObjType     = (TObjString*)(arrNextCol->At(0));
503       TObjString *sObjFileName = (TObjString*)(arrNextCol->At(1));
504       
505       if ( !sObjType || ! sObjFileName ) continue;
506       TString sType(sObjType->GetString());
507       TString sFileName(sObjFileName->GetString());
508       printf("%s\t%s\n",sType.Data(),sFileName.Data());
509       
510       TFile *fIn = TFile::Open(sFileName);
511       if ( !fIn ){
512          fprintf(stderr,"File not found: '%s'", sFileName.Data());
513          continue;
514       }
515       
516       if ( sType == "CE" ){
517          AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
518          
519          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());         
520          calPad->SetNameTitle("CETmean","CETmean");
521          calibObjects->Add(calPad);
522          
523          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());         
524          calPad->SetNameTitle("CEQmean","CEQmean");
525          calibObjects->Add(calPad);        
526          
527          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
528          calPad->SetNameTitle("CETrms","CETrms");
529          calibObjects->Add(calPad);         
530                   
531       } else if ( sType == "Pulser") {
532          AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
533          
534          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());         
535          calPad->SetNameTitle("PulserTmean","PulserTmean");
536          calibObjects->Add(calPad);
537          
538          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());         
539          calPad->SetNameTitle("PulserQmean","PulserQmean");
540          calibObjects->Add(calPad);        
541          
542          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
543          calPad->SetNameTitle("PulserTrms","PulserTrms");
544          calibObjects->Add(calPad);         
545       
546       } else if ( sType == "Pedestals") {
547          AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
548          
549          calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());         
550          calPad->SetNameTitle("Pedestals","Pedestals");
551          calibObjects->Add(calPad);
552          
553          calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());         
554          calPad->SetNameTitle("Noise","Noise");
555          calibObjects->Add(calPad);        
556      
557       } else {
558          fprintf(stderr,"Undefined Type: '%s'",sType.Data());
559          
560       }
561       delete fIn;
562    }
563 }
564
565
566
567 void AliTPCcalibDB::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
568   //
569   // Write a tree with all available information
570   // if mapFileName is specified, the Map information are also written to the tree
571   // pads specified in outlierPad are not used for calculating statistics
572   //  - the same function as AliTPCCalPad::MakeTree - 
573   //
574    AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
575
576    TObjArray* mapIROCs = 0;
577    TObjArray* mapOROCs = 0;
578    TVectorF *mapIROCArray = 0;
579    TVectorF *mapOROCArray = 0;
580    Int_t mapEntries = 0;
581    TString* mapNames = 0;
582    
583    if (mapFileName) {
584       TFile mapFile(mapFileName, "read");
585       
586       TList* listOfROCs = mapFile.GetListOfKeys();
587       mapEntries = listOfROCs->GetEntries()/2;
588       mapIROCs = new TObjArray(mapEntries*2);
589       mapOROCs = new TObjArray(mapEntries*2);
590       mapIROCArray = new TVectorF[mapEntries];
591       mapOROCArray = new TVectorF[mapEntries];
592       
593       mapNames = new TString[mapEntries];
594       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
595         TString nameROC(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
596          nameROC.Remove(nameROC.Length()-4, 4);
597          mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "IROC").Data()), ivalue);
598          mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "OROC").Data()), ivalue);
599          mapNames[ivalue].Append(nameROC);
600       }
601       
602       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
603          mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
604          mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
605       
606          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
607             (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
608          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
609             (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
610       }
611
612    } //  if (mapFileName)
613   
614    TTreeSRedirector cstream(fileName);
615    Int_t arrayEntries = array->GetEntries();
616    
617    TString* names = new TString[arrayEntries];
618    for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
619       names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
620
621    for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
622       //
623       // get statistic for given sector
624       //
625       TVectorF median(arrayEntries);
626       TVectorF mean(arrayEntries);
627       TVectorF rms(arrayEntries);
628       TVectorF ltm(arrayEntries);
629       TVectorF ltmrms(arrayEntries);
630       TVectorF medianWithOut(arrayEntries);
631       TVectorF meanWithOut(arrayEntries);
632       TVectorF rmsWithOut(arrayEntries);
633       TVectorF ltmWithOut(arrayEntries);
634       TVectorF ltmrmsWithOut(arrayEntries);
635       
636       TVectorF *vectorArray = new TVectorF[arrayEntries];
637       for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
638          vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
639       
640       for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
641          AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
642          AliTPCCalROC* calROC = calPad->GetCalROC(isector);
643          AliTPCCalROC* outlierROC = 0;
644          if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
645          if (calROC) {
646             median[ivalue] = calROC->GetMedian();
647             mean[ivalue] = calROC->GetMean();
648             rms[ivalue] = calROC->GetRMS();
649             Double_t ltmrmsValue = 0;
650             ltm[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction);
651             ltmrms[ivalue] = ltmrmsValue;
652             if (outlierROC) {
653                medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
654                meanWithOut[ivalue] = calROC->GetMean(outlierROC);
655                rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
656                ltmrmsValue = 0;
657                ltmWithOut[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction, outlierROC);
658                ltmrmsWithOut[ivalue] = ltmrmsValue;
659             }
660          }
661          else {
662             median[ivalue] = 0.;
663             mean[ivalue] = 0.;
664             rms[ivalue] = 0.;
665             ltm[ivalue] = 0.;
666             ltmrms[ivalue] = 0.;
667             medianWithOut[ivalue] = 0.;
668             meanWithOut[ivalue] = 0.;
669             rmsWithOut[ivalue] = 0.;
670             ltmWithOut[ivalue] = 0.;
671             ltmrmsWithOut[ivalue] = 0.;
672          }
673       }
674       
675       //
676       // fill vectors of variable per pad
677       //
678       TVectorF *posArray = new TVectorF[8];
679       for (Int_t ivalue = 0; ivalue < 8; ivalue++)
680          posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
681
682       Float_t posG[3] = {0};
683       Float_t posL[3] = {0};
684       Int_t ichannel = 0;
685       for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
686          for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
687             tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
688             tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
689             posArray[0][ichannel] = irow;
690             posArray[1][ichannel] = ipad;
691             posArray[2][ichannel] = posL[0];
692             posArray[3][ichannel] = posL[1];
693             posArray[4][ichannel] = posG[0];
694             posArray[5][ichannel] = posG[1];
695             posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
696             posArray[7][ichannel] = ichannel;
697             
698             // loop over array containing AliTPCCalPads
699             for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
700                AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
701                AliTPCCalROC* calROC = calPad->GetCalROC(isector);
702                if (calROC)
703                   (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
704                else
705                   (vectorArray[ivalue])[ichannel] = 0;
706             }
707             ichannel++;
708          }
709       }
710       
711       cstream << "calPads" <<
712          "sector=" << isector;
713       
714       for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
715          cstream << "calPads" <<
716             (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
717             (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
718             (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
719             (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
720             (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
721          if (outlierPad) {
722             cstream << "calPads" <<
723                (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
724                (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
725                (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
726                (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
727                (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
728          }
729       }
730
731       for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
732          cstream << "calPads" <<
733             (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
734       }
735
736       if (mapFileName) {
737          for  (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
738             if (isector < 36)
739                cstream << "calPads" <<
740                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
741             else
742                cstream << "calPads" <<
743                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
744          }
745       }
746
747       cstream << "calPads" <<
748          "row.=" << &posArray[0] <<
749          "pad.=" << &posArray[1] <<
750          "lx.=" << &posArray[2] <<
751          "ly.=" << &posArray[3] <<
752          "gx.=" << &posArray[4] <<
753          "gy.=" << &posArray[5] <<
754          "rpad.=" << &posArray[6] <<
755          "channel.=" << &posArray[7];
756          
757       cstream << "calPads" <<
758          "\n";
759
760       delete[] posArray;
761       delete[] vectorArray;
762    }
763    
764
765    delete[] names;
766    if (mapFileName) {
767       delete mapIROCs;
768       delete mapOROCs;
769       delete[] mapIROCArray;
770       delete[] mapOROCArray;
771       delete[] mapNames;
772    }
773 }
774
775 Int_t AliTPCcalibDB::GetRCUTriggerConfig() const
776 {
777   //
778   // return the RCU trigger configuration register
779   //
780   TMap *map=GetRCUconfig();
781   if (!map) return -1;
782   TVectorF *v=(TVectorF*)map->GetValue("TRGCONF_TRG_MODE");
783   Float_t mode=-1;
784   for (Int_t i=0; i<v->GetNrows(); ++i){
785     Float_t newmode=v->GetMatrixArray()[i];
786     if (newmode>-1){
787       if (mode>-1&&newmode!=mode) AliWarning("Found different RCU trigger configurations!!!");
788       mode=newmode;
789     }
790   }
791   return (Int_t)mode;
792 }
793
794 Bool_t AliTPCcalibDB::IsTrgL0() 
795 {
796   //
797   // return if the FEE readout was triggered on L0
798   //
799   Int_t mode=GetRCUTriggerConfig();
800   if (mode<0) return kFALSE;
801   return (mode==1);
802 }
803
804 Bool_t AliTPCcalibDB::IsTrgL1()
805 {
806   //
807   // return if the FEE readout was triggered on L1
808   //
809   Int_t mode=GetRCUTriggerConfig();
810   if (mode<0) return kFALSE;
811   return (mode==0);
812 }
813
814 void AliTPCcalibDB::RegisterExB(Int_t index, Float_t bz, Bool_t bdelete){
815   //
816   // Register static ExB correction map
817   // index - registration index - used for visualization
818   // bz    - bz field in kGaus
819
820   //  Float_t factor =  bz/(-5.);  // default b filed in Cheb with minus sign
821   Float_t factor =  bz/(5.);  // default b filed in Cheb with minus sign
822                               // was chenged in the Revision ???? (Ruben can you add here number)
823   
824   AliMagF*   bmap = new AliMagF("MapsExB","MapsExB", factor,TMath::Sign(1.f,factor),AliMagF::k5kG);
825   
826   AliTPCExBFirst *exb  = new  AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
827   AliTPCExB::SetInstance(exb);
828   
829   if (bdelete){
830     delete bmap;
831   }else{
832     AliTPCExB::RegisterField(index,bmap);
833   }
834   if (index>=fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
835   fgExBArray.AddAt(exb,index);
836 }
837
838
839 AliTPCExB*    AliTPCcalibDB::GetExB(Float_t bz, Bool_t deleteB) {
840   //
841   // bz filed in KGaus not in tesla
842   // Get ExB correction map
843   // if doesn't exist - create it
844   //
845   Int_t index = TMath::Nint(5+bz);
846   if (index>fgExBArray.GetEntries()) fgExBArray.Expand((index+1)*2+11);
847   if (!fgExBArray.At(index)) AliTPCcalibDB::RegisterExB(index,bz,deleteB);
848   return (AliTPCExB*)fgExBArray.At(index);
849 }
850
851
852 void  AliTPCcalibDB::SetExBField(Float_t bz){
853   //
854   // Set magnetic filed for ExB correction
855   //
856   fExB = GetExB(bz,kFALSE);
857 }
858
859 void  AliTPCcalibDB::SetExBField(const AliMagF*   bmap){
860   //
861   // Set magnetic field for ExB correction
862   //
863   AliTPCExBFirst *exb  = new  AliTPCExBFirst(bmap,0.88*2.6400e+04,50,50,50);
864   AliTPCExB::SetInstance(exb);
865   fExB=exb;
866 }
867
868
869
870
871
872 void AliTPCcalibDB::UpdateRunInformations( Int_t run, Bool_t force){
873   //
874   // - > Don't use it for reconstruction - Only for Calibration studies
875   //
876   if (run<=0) return;
877   AliCDBEntry * entry = 0;
878   if (run>= fRunList.fN){
879     fRunList.Set(run*2+1);
880     fGRPArray.Expand(run*2+1);
881     fGRPMaps.Expand(run*2+1);
882     fGoofieArray.Expand(run*2+1);
883     fVoltageArray.Expand(run*2+1); 
884     fTemperatureArray.Expand(run*2+1);
885     fVdriftArray.Expand(run*2+1);
886     fDriftCorrectionArray.Expand(run*2+1);
887     fTimeGainSplinesArray.Expand(run*2+1);
888     //
889     //
890     fALTROConfigData->Expand(run*2+1);    // ALTRO configuration data
891     fPulserData->Expand(run*2+1);         // Calibration Pulser data
892     fCEData->Expand(run*2+1);             // CE data
893     if (!fTimeGainSplines) fTimeGainSplines = new TObjArray(run*2+1);
894     fTimeGainSplines->Expand(run*2+1); // Array of AliSplineFits: at 0 MIP position in
895   }
896   if (fRunList[run]>0 &&force==kFALSE) return;
897
898   fRunList[run]=1;  // sign as used
899
900   //
901   entry = AliCDBManager::Instance()->Get("GRP/GRP/Data",run);
902   if (entry)  {
903     AliGRPObject * grpRun = dynamic_cast<AliGRPObject*>(entry->GetObject());
904     if (!grpRun){
905       TMap* map = dynamic_cast<TMap*>(entry->GetObject());
906       if (map){
907         //grpRun = new AliGRPObject; 
908         //grpRun->ReadValuesFromMap(map);
909         grpRun =  MakeGRPObjectFromMap(map);
910
911         fGRPMaps.AddAt(map,run);
912       }
913     }
914     fGRPArray.AddAt(grpRun,run);
915   }
916   entry = AliCDBManager::Instance()->Get("TPC/Calib/Goofie",run);
917   if (entry){
918     fGoofieArray.AddAt(entry->GetObject(),run);
919   }
920   //
921   entry = AliCDBManager::Instance()->Get("TPC/Calib/HighVoltage",run);
922   if (entry)  {
923     fVoltageArray.AddAt(entry->GetObject(),run);
924   }
925   //
926   entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeGain",run);
927   if (entry)  {
928     fTimeGainSplinesArray.AddAt(entry->GetObject(),run);
929   }
930   //
931   entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeDrift",run);
932   if (entry)  {
933     fDriftCorrectionArray.AddAt(entry->GetObject(),run);
934   }
935   //
936   entry = AliCDBManager::Instance()->Get("TPC/Calib/Temperature",run);
937   if (entry)  {
938     fTemperatureArray.AddAt(entry->GetObject(),run);
939   }
940   //apply fDButil filters
941
942   fDButil->UpdateFromCalibDB();
943   if (fTemperature) fDButil->FilterTemperature(fTemperature);
944
945   AliDCSSensor * press = GetPressureSensor(run,0);
946   AliTPCSensorTempArray * temp = GetTemperatureSensor(run);
947   Bool_t accept=kTRUE;
948   if (temp) {
949     accept = fDButil->FilterTemperature(temp)>0.1;
950   }
951   if (press) {
952     const Double_t kMinP=950.;
953     const Double_t kMaxP=1050.;
954     const Double_t kMaxdP=10.;
955     const Double_t kSigmaCut=4.;
956     fDButil->FilterSensor(press,kMinP,kMaxP,kMaxdP,kSigmaCut);
957     if (press->GetFit()==0) accept=kFALSE;
958   }
959   if (press && temp &&accept){
960     AliTPCCalibVdrift * vdrift = new AliTPCCalibVdrift(temp, press,0);
961     fVdriftArray.AddAt(vdrift,run);
962   }
963   fDButil->FilterCE(120., 3., 4.,0);
964   fDButil->FilterTracks(run, 10.,0);
965 }
966
967
968 Float_t AliTPCcalibDB::GetGain(Int_t sector, Int_t row, Int_t pad){
969   //
970   //
971   AliTPCCalPad *calPad = Instance()->fDedxGainFactor;;
972   if (!calPad) return 0;
973   return calPad->GetCalROC(sector)->GetValue(row,pad);
974 }
975
976 AliSplineFit* AliTPCcalibDB::GetVdriftSplineFit(const char* name, Int_t run){
977   //
978   //
979   //
980   TObjArray *arr=GetTimeVdriftSplineRun(run);
981   if (!arr) return 0;
982   return dynamic_cast<AliSplineFit*>(arr->FindObject(name));
983 }
984
985 AliSplineFit* AliTPCcalibDB::CreateVdriftSplineFit(const char* graphName, Int_t run){
986   //
987   // create spline fit from the drift time graph in TimeDrift
988   //
989   TObjArray *arr=GetTimeVdriftSplineRun(run);
990   if (!arr) return 0;
991   TGraph *graph=dynamic_cast<TGraph*>(arr->FindObject(graphName));
992   if (!graph) return 0;
993   AliSplineFit *fit = new AliSplineFit();
994   fit->SetGraph(graph);
995   fit->SetMinPoints(graph->GetN()+1);
996   fit->InitKnots(graph,2,0,0.001);
997   fit->SplineFit(0);
998   return fit;
999 }
1000
1001 AliGRPObject *AliTPCcalibDB::GetGRP(Int_t run){
1002   //
1003   // Get GRP object for given run 
1004   //
1005   AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>((Instance()->fGRPArray).At(run));
1006   if (!grpRun) {
1007     Instance()->UpdateRunInformations(run);
1008     grpRun = dynamic_cast<AliGRPObject *>(Instance()->fGRPArray.At(run));
1009     if (!grpRun) return 0; 
1010   }
1011   return grpRun;
1012 }
1013
1014 TMap *  AliTPCcalibDB::GetGRPMap(Int_t run){
1015   //
1016   //
1017   //
1018   TMap * grpRun = dynamic_cast<TMap *>((Instance()->fGRPMaps).At(run));
1019   if (!grpRun) {
1020     Instance()->UpdateRunInformations(run);
1021     grpRun = dynamic_cast<TMap *>(Instance()->fGRPMaps.At(run));
1022     if (!grpRun) return 0; 
1023   }
1024   return grpRun;
1025 }
1026
1027
1028 AliDCSSensor * AliTPCcalibDB::GetPressureSensor(Int_t run, Int_t type){
1029   //
1030   // Get Pressure sensor
1031   // run  = run number
1032   // type = 0 - Cavern pressure
1033   //        1 - Suface pressure
1034   // First try to get if trom map - if existing  (Old format of data storing)
1035   //
1036
1037
1038   TMap *map = GetGRPMap(run);  
1039   if (map){
1040     AliDCSSensor * sensor = 0;
1041     TObject *osensor=0;
1042     if (type==0) osensor = ((*map)("fCavernPressure"));
1043     if (type==1) osensor = ((*map)("fP2Pressure"));
1044     sensor =dynamic_cast<AliDCSSensor *>(osensor); 
1045     if (sensor) return sensor;
1046   }
1047   //
1048   // If not map try to get it from the GRPObject
1049   //
1050   AliGRPObject * grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run)); 
1051   if (!grpRun) {
1052     UpdateRunInformations(run);
1053     grpRun = dynamic_cast<AliGRPObject *>(fGRPArray.At(run));
1054     if (!grpRun) return 0; 
1055   }
1056   AliDCSSensor * sensor = grpRun->GetCavernAtmosPressure();
1057   if (type==1) sensor = grpRun->GetSurfaceAtmosPressure();
1058   return sensor; 
1059 }
1060
1061 AliTPCSensorTempArray * AliTPCcalibDB::GetTemperatureSensor(Int_t run){
1062   //
1063   // Get temperature sensor array
1064   //
1065   AliTPCSensorTempArray * tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
1066   if (!tempArray) {
1067     UpdateRunInformations(run);
1068     tempArray = (AliTPCSensorTempArray *)fTemperatureArray.At(run);
1069   }
1070   return tempArray;
1071 }
1072
1073
1074 TObjArray * AliTPCcalibDB::GetTimeGainSplinesRun(Int_t run){
1075   //
1076   // Get temperature sensor array
1077   //
1078   TObjArray * gainSplines = (TObjArray *)fTimeGainSplinesArray.At(run);
1079   if (!gainSplines) {
1080     UpdateRunInformations(run);
1081     gainSplines = (TObjArray *)fTimeGainSplinesArray.At(run);
1082   }
1083   return gainSplines;
1084 }
1085
1086 TObjArray * AliTPCcalibDB::GetTimeVdriftSplineRun(Int_t run){
1087   //
1088   // Get drift spline array
1089   //
1090   TObjArray * driftSplines = (TObjArray *)fDriftCorrectionArray.At(run);
1091   if (!driftSplines) {
1092     UpdateRunInformations(run);
1093     driftSplines = (TObjArray *)fDriftCorrectionArray.At(run);
1094   }
1095   return driftSplines;
1096 }
1097
1098 AliDCSSensorArray * AliTPCcalibDB::GetVoltageSensors(Int_t run){
1099   //
1100   // Get temperature sensor array
1101   //
1102   AliDCSSensorArray * voltageArray = (AliDCSSensorArray *)fVoltageArray.At(run);
1103   if (!voltageArray) {
1104     UpdateRunInformations(run);
1105     voltageArray = (AliDCSSensorArray *)fVoltageArray.At(run);
1106   }
1107   return voltageArray;
1108 }
1109
1110 AliDCSSensorArray * AliTPCcalibDB::GetGoofieSensors(Int_t run){
1111   //
1112   // Get temperature sensor array
1113   //
1114   AliDCSSensorArray * goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
1115   if (!goofieArray) {
1116     UpdateRunInformations(run);
1117     goofieArray = (AliDCSSensorArray *)fGoofieArray.At(run);
1118   }
1119   return goofieArray;
1120 }
1121
1122
1123
1124 AliTPCCalibVdrift *     AliTPCcalibDB::GetVdrift(Int_t run){
1125   //
1126   // Get the interface to the the vdrift 
1127   //
1128   AliTPCCalibVdrift  * vdrift = (AliTPCCalibVdrift*)fVdriftArray.At(run);
1129   if (!vdrift) {
1130     UpdateRunInformations(run);
1131     vdrift= (AliTPCCalibVdrift*)fVdriftArray.At(run);
1132   }
1133   return vdrift;
1134 }
1135
1136 Float_t AliTPCcalibDB::GetCEdriftTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1137 {
1138   //
1139   // GetCE drift time information for 'sector'
1140   // sector 72 is the mean drift time of the A-Side
1141   // sector 73 is the mean drift time of the C-Side
1142   // it timestamp==-1 return mean value
1143   //
1144   AliTPCcalibDB::Instance()->SetRun(run);
1145   TGraph *gr=AliTPCcalibDB::Instance()->GetCErocTgraph(sector);
1146   if (!gr||sector<0||sector>73) {
1147     if (entries) *entries=0;
1148     return 0.;
1149   }
1150   Float_t val=0.;
1151   if (timeStamp==-1.){
1152     val=gr->GetMean(2);
1153   }else{
1154     for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1155       Double_t x,y;
1156       gr->GetPoint(ipoint,x,y);
1157       if (x<timeStamp) continue;
1158       val=y;
1159       break;
1160     }
1161   }
1162   return val;
1163 }
1164   
1165 Float_t AliTPCcalibDB::GetCEchargeTime(Int_t run, Int_t sector, Double_t timeStamp, Int_t *entries)
1166 {
1167   //
1168   // GetCE mean charge for 'sector'
1169   // it timestamp==-1 return mean value
1170   //
1171   AliTPCcalibDB::Instance()->SetRun(run);
1172   TGraph *gr=AliTPCcalibDB::Instance()->GetCErocQgraph(sector);
1173   if (!gr||sector<0||sector>71) {
1174     if (entries) *entries=0;
1175     return 0.;
1176   }
1177   Float_t val=0.;
1178   if (timeStamp==-1.){
1179     val=gr->GetMean(2);
1180   }else{
1181     for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1182       Double_t x,y;
1183       gr->GetPoint(ipoint,x,y);
1184       if (x<timeStamp) continue;
1185       val=y;
1186       break;
1187     }
1188   }
1189   return val;
1190 }
1191
1192 Float_t AliTPCcalibDB::GetDCSSensorValue(AliDCSSensorArray *arr, Int_t timeStamp, const char * sensorName, Int_t sigDigits)
1193 {
1194   //
1195   // Get Value for a DCS sensor 'sensorName', run 'run' at time 'timeStamp'
1196   //
1197   Float_t val=0;
1198   const TString sensorNameString(sensorName);
1199   AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1200   if (!sensor) return val;
1201   //use the dcs graph if possible
1202   TGraph *gr=sensor->GetGraph();
1203   if (gr){
1204     for (Int_t ipoint=0;ipoint<gr->GetN();++ipoint){
1205       Double_t x,y;
1206       gr->GetPoint(ipoint,x,y);
1207       Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1208       if (time<timeStamp) continue;
1209       val=y;
1210       break;
1211     }
1212     //if val is still 0, test if if the requested time if within 5min of the first/last
1213     //data point. If this is the case return the firs/last entry
1214     //the timestamps might not be syncronised for all calibration types, sometimes a 'pre'
1215     //and 'pos' period is requested. Especially to the HV this is not the case!
1216     //first point
1217     if (val==0 ){
1218       Double_t x,y;
1219       gr->GetPoint(0,x,y);
1220       Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1221       if ((time-timeStamp)<5*60) val=y;
1222     }
1223     //last point
1224     if (val==0 ){
1225       Double_t x,y;
1226       gr->GetPoint(gr->GetN()-1,x,y);
1227       Int_t time=TMath::Nint(sensor->GetStartTime()+x*3600); //time in graph is hours
1228       if ((timeStamp-time)<5*60) val=y;
1229     }
1230   } else {
1231     val=sensor->GetValue(timeStamp);
1232   }
1233   if (sigDigits>=0){
1234     val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1235   }
1236   return val;
1237 }
1238
1239 Float_t AliTPCcalibDB::GetDCSSensorMeanValue(AliDCSSensorArray *arr, const char * sensorName, Int_t sigDigits)
1240 {
1241   //
1242   // Get mean Value for a DCS sensor 'sensorName' during run 'run'
1243   //
1244   Float_t val=0;
1245   const TString sensorNameString(sensorName);
1246   AliDCSSensor *sensor = arr->GetSensor(sensorNameString);
1247   if (!sensor) return val;
1248
1249   //use dcs graph if it exists
1250   TGraph *gr=sensor->GetGraph();
1251   if (gr){
1252     val=gr->GetMean(2);
1253   } else {
1254     //if we don't have the dcs graph, try to get some meaningful information
1255     if (!sensor->GetFit()) return val;
1256     Int_t nKnots=sensor->GetFit()->GetKnots();
1257     Double_t tMid=(sensor->GetEndTime()-sensor->GetStartTime())/2.;
1258     for (Int_t iKnot=0;iKnot<nKnots;++iKnot){
1259       if (sensor->GetFit()->GetX()[iKnot]>tMid/3600.) break;
1260       val=(Float_t)sensor->GetFit()->GetY0()[iKnot];
1261     }
1262   }
1263   if (sigDigits>=0){
1264     val/=10;
1265     val=(Float_t)TMath::Floor(val * TMath::Power(10., sigDigits) + .5) / TMath::Power(10., sigDigits);
1266     val*=10;
1267   }
1268   return val;
1269 }
1270
1271 Float_t AliTPCcalibDB::GetChamberHighVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits) {
1272   //
1273   // return the chamber HV for given run and time: 0-35 IROC, 36-72 OROC
1274   // if timeStamp==-1 return mean value
1275   //
1276   Float_t val=0;
1277   TString sensorName="";
1278   TTimeStamp stamp(timeStamp);
1279   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1280   if (!voltageArray || (sector<0) || (sector>71)) return val;
1281   Char_t sideName='A';
1282   if ((sector/18)%2==1) sideName='C';
1283   if (sector<36){
1284     //IROC
1285     sensorName=Form("TPC_ANODE_I_%c%02d_VMEAS",sideName,sector%18);
1286   }else{
1287     //OROC
1288     sensorName=Form("TPC_ANODE_O_%c%02d_0_VMEAS",sideName,sector%18);
1289   }
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 Float_t AliTPCcalibDB::GetSkirtVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1298 {
1299   //
1300   // Get the skirt voltage for 'run' at 'timeStamp' and 'sector': 0-35 IROC, 36-72 OROC
1301   // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1302   // if timeStamp==-1 return the mean value for the run
1303   //
1304   Float_t val=0;
1305   TString sensorName="";
1306   TTimeStamp stamp(timeStamp);
1307   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1308   if (!voltageArray || (sector<0) || (sector>71)) return val;
1309   Char_t sideName='A';
1310   if ((sector/18)%2==1) sideName='C';
1311   sensorName=Form("TPC_SKIRT_%c_VMEAS",sideName);
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::GetCoverVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1321 {
1322   //
1323   // Get the cover 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_COVER_I_%c_VMEAS",sideName);
1337   }else{
1338     //OROC
1339     sensorName=Form("TPC_COVER_O_%c_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::GetGGoffsetVoltage(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_OFF_VMEAS",sideName);
1366   }else{
1367     //OROC
1368     sensorName=Form("TPC_GATE_O_%c_OFF_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::GetGGnegVoltage(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_NEG_VMEAS",sideName);
1395   }else{
1396     //OROC
1397     sensorName=Form("TPC_GATE_O_%c_NEG_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::GetGGposVoltage(Int_t run, Int_t sector, Int_t timeStamp, Int_t sigDigits)
1408 {
1409   //
1410   // Get the GG offset voltage for run 'run' at time 'timeStamp'
1411   // type corresponds to the following: 0 - IROC A-Side; 1 - IROC C-Side; 2 - OROC A-Side; 3 - OROC C-Side
1412   // if timeStamp==-1 return the mean value for the run
1413   //
1414   Float_t val=0;
1415   TString sensorName="";
1416   TTimeStamp stamp(timeStamp);
1417   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(run);
1418   if (!voltageArray || (sector<0) || (sector>71)) return val;
1419   Char_t sideName='A';
1420   if ((sector/18)%2==1) sideName='C';
1421   if (sector<36){
1422     //IROC
1423     sensorName=Form("TPC_GATE_I_%c_POS_VMEAS",sideName);
1424   }else{
1425     //OROC
1426     sensorName=Form("TPC_GATE_O_%c_POS_VMEAS",sideName);
1427   }
1428   if (timeStamp==-1){
1429     val=AliTPCcalibDB::GetDCSSensorMeanValue(voltageArray, sensorName.Data(),sigDigits);
1430   } else {
1431     val=AliTPCcalibDB::GetDCSSensorValue(voltageArray, timeStamp, sensorName.Data(),sigDigits);
1432   }
1433   return val;
1434 }
1435
1436 Float_t AliTPCcalibDB::GetPressure(Int_t timeStamp, Int_t run, Int_t type){
1437   //
1438   // GetPressure for given time stamp and runt
1439   //
1440   TTimeStamp stamp(timeStamp);
1441   AliDCSSensor * sensor = Instance()->GetPressureSensor(run,type);
1442   if (!sensor) return 0;
1443   return sensor->GetValue(stamp);
1444 }
1445
1446 Float_t AliTPCcalibDB::GetL3Current(Int_t run, Int_t statType){
1447   //
1448   // return L3 current
1449   // stat type is: AliGRPObject::Stats: kMean = 0, kTruncMean = 1, kMedian = 2, kSDMean = 3, kSDMedian = 4
1450   //
1451   Float_t current=-1;
1452   AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1453   if (grp) current=grp->GetL3Current((AliGRPObject::Stats)statType);
1454   return current;
1455 }
1456
1457 Float_t AliTPCcalibDB::GetBz(Int_t run){
1458   //
1459   // calculate BZ in T from L3 current
1460   //
1461   Float_t bz=-1;
1462   Float_t current=AliTPCcalibDB::GetL3Current(run);
1463   if (current>-1) bz=5*current/30000.*.1;
1464   return bz;
1465 }
1466
1467 Char_t  AliTPCcalibDB::GetL3Polarity(Int_t run) {
1468   //
1469   // get l3 polarity from GRP
1470   //
1471   Char_t pol=-100;
1472   AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1473   if (grp) pol=grp->GetL3Polarity();
1474   return pol;
1475 }
1476
1477 TString AliTPCcalibDB::GetRunType(Int_t run){
1478   //
1479   // return run type from grp
1480   //
1481
1482 //   TString type("UNKNOWN");
1483   AliGRPObject *grp=AliTPCcalibDB::GetGRP(run);
1484   if (grp) return grp->GetRunType();
1485   return "UNKNOWN";
1486 }
1487
1488 Float_t AliTPCcalibDB::GetValueGoofie(Int_t timeStamp, Int_t run, Int_t type){
1489   //
1490   // GetPressure for given time stamp and runt
1491   //
1492   TTimeStamp stamp(timeStamp);
1493   AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(run);
1494   if (!goofieArray) return 0;
1495   AliDCSSensor *sensor = goofieArray->GetSensor(type);
1496   return sensor->GetValue(stamp);
1497 }
1498
1499
1500
1501
1502
1503
1504 Bool_t  AliTPCcalibDB::GetTemperatureFit(Int_t timeStamp, Int_t run, Int_t side,TVectorD& fit){
1505   //
1506   //
1507   //
1508   TTimeStamp tstamp(timeStamp);
1509   AliTPCSensorTempArray* tempArray  = Instance()->GetTemperatureSensor(run);
1510   if (! tempArray) return kFALSE;
1511   AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
1512   TLinearFitter * fitter = tempMap->GetLinearFitter(3,side,tstamp);
1513   if (fitter){
1514     fitter->Eval(); 
1515     fitter->GetParameters(fit);
1516   }
1517   delete fitter;
1518   delete tempMap;
1519   if (!fitter) return kFALSE;
1520   return kTRUE;
1521 }
1522
1523 Float_t AliTPCcalibDB::GetTemperature(Int_t timeStamp, Int_t run, Int_t side){
1524   //
1525   //
1526   //
1527   TVectorD vec(5);
1528   if (side==0) {
1529     GetTemperatureFit(timeStamp,run,0,vec);
1530     return vec[0];
1531   }
1532   if (side==1){
1533     GetTemperatureFit(timeStamp,run,0,vec);
1534     return vec[0];
1535   }
1536   return 0;
1537 }
1538
1539
1540 Double_t AliTPCcalibDB::GetPTRelative(UInt_t timeSec, Int_t run, Int_t side){
1541   //
1542   // Get relative P/T 
1543   // time - absolute time
1544   // run  - run number
1545   // side - 0 - A side   1-C side
1546   AliTPCCalibVdrift * vdrift =  Instance()->GetVdrift(run);
1547   if (!vdrift) return 0;
1548   return vdrift->GetPTRelative(timeSec,side);
1549 }
1550
1551 AliGRPObject * AliTPCcalibDB::MakeGRPObjectFromMap(TMap *map){
1552   //
1553   // Function to covert old GRP run information from TMap to GRPObject
1554   //
1555   //  TMap * map = AliTPCcalibDB::GetGRPMap(52406);
1556   if (!map) return 0;
1557   AliDCSSensor * sensor = 0;
1558   TObject *osensor=0;
1559   osensor = ((*map)("fP2Pressure"));
1560   sensor  =dynamic_cast<AliDCSSensor *>(osensor); 
1561   //
1562   if (!sensor) return 0;
1563   //
1564   AliDCSSensor * sensor2 = new AliDCSSensor(*sensor);
1565   osensor = ((*map)("fCavernPressure"));
1566   TGraph * gr = new TGraph(2);
1567   gr->GetX()[0]= -100000.;
1568   gr->GetX()[1]= 1000000.;
1569   gr->GetY()[0]= atof(osensor->GetName());
1570   gr->GetY()[1]= atof(osensor->GetName());
1571   sensor2->SetGraph(gr);
1572   sensor2->SetFit(0);
1573   
1574
1575   AliGRPObject *grpRun = new AliGRPObject; 
1576   grpRun->ReadValuesFromMap(map);
1577   grpRun->SetCavernAtmosPressure(sensor2);
1578   grpRun->SetSurfaceAtmosPressure(sensor);
1579   return grpRun;
1580 }
1581
1582 Bool_t AliTPCcalibDB::CreateGUITree(Int_t run, const char* filename)
1583 {
1584   //
1585   // Create a gui tree for run number 'run'
1586   //
1587
1588   if (!AliCDBManager::Instance()->GetDefaultStorage()){
1589     AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
1590                     MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
1591     return kFALSE;
1592   }
1593   //db instance
1594   AliTPCcalibDB *db=AliTPCcalibDB::Instance();
1595   // retrieve cal pad objects
1596   db->SetRun(run);
1597   db->CreateGUITree(filename);
1598   return kTRUE;
1599 }
1600
1601 Bool_t AliTPCcalibDB::CreateGUITree(const char* filename){
1602   //
1603   //
1604   //
1605   if (!AliCDBManager::Instance()->GetDefaultStorage()){
1606     AliError("Default Storage not set. Cannot create calibration Tree!");
1607     return kFALSE;
1608   }
1609   UpdateNonRec();  // load all infromation now
1610
1611   AliTPCPreprocessorOnline prep;
1612   //noise and pedestals
1613   if (GetPedestals()) prep.AddComponent(new AliTPCCalPad(*(GetPedestals())));
1614   if (GetPadNoise() ) prep.AddComponent(new AliTPCCalPad(*(GetPadNoise())));
1615   //pulser data
1616   if (GetPulserTmean()) prep.AddComponent(new AliTPCCalPad(*(GetPulserTmean())));
1617   if (GetPulserTrms() ) prep.AddComponent(new AliTPCCalPad(*(GetPulserTrms())));
1618   if (GetPulserQmean()) prep.AddComponent(new AliTPCCalPad(*(GetPulserQmean())));
1619   //CE data
1620   if (GetCETmean()) prep.AddComponent(new AliTPCCalPad(*(GetCETmean())));
1621   if (GetCETrms() ) prep.AddComponent(new AliTPCCalPad(*(GetCETrms())));
1622   if (GetCEQmean()) prep.AddComponent(new AliTPCCalPad(*(GetCEQmean())));
1623   //Altro data
1624   if (GetALTROAcqStart() ) prep.AddComponent(new AliTPCCalPad(*(GetALTROAcqStart() )));
1625   if (GetALTROZsThr()    ) prep.AddComponent(new AliTPCCalPad(*(GetALTROZsThr()    )));
1626   if (GetALTROFPED()     ) prep.AddComponent(new AliTPCCalPad(*(GetALTROFPED()     )));
1627   if (GetALTROAcqStop()  ) prep.AddComponent(new AliTPCCalPad(*(GetALTROAcqStop()  )));
1628   if (GetALTROMasked()   ) prep.AddComponent(new AliTPCCalPad(*(GetALTROMasked()   )));
1629   //QA
1630   AliTPCdataQA *dataQA=GetDataQA();
1631   if (dataQA) {
1632     if (dataQA->GetNLocalMaxima())
1633       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNLocalMaxima())));
1634     if (dataQA->GetMaxCharge())
1635       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetMaxCharge())));
1636     if (dataQA->GetMeanCharge())
1637       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetMeanCharge())));
1638     if (dataQA->GetNoThreshold())
1639       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNoThreshold())));
1640     if (dataQA->GetNTimeBins())
1641       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNTimeBins())));
1642     if (dataQA->GetNPads())
1643       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetNPads())));
1644     if (dataQA->GetTimePosition())
1645       prep.AddComponent(new AliTPCCalPad(*(dataQA->GetTimePosition())));
1646   }
1647   
1648   //
1649   TString file(filename);
1650   if (file.IsNull()) file=Form("guiTreeRun_%d.root",fRun);
1651   prep.DumpToFile(file.Data());
1652   return kTRUE;
1653 }
1654
1655 Bool_t AliTPCcalibDB::CreateRefFile(Int_t run, const char* filename)
1656 {
1657   //
1658   // Create a gui tree for run number 'run'
1659   //
1660   
1661   if (!AliCDBManager::Instance()->GetDefaultStorage()){
1662     AliLog::Message(AliLog::kError, "Default Storage not set. Cannot create Calibration Tree!",
1663                     MODULENAME(), "AliTPCcalibDB", FUNCTIONNAME(), __FILE__, __LINE__);
1664     return kFALSE;
1665   }
1666   TString file(filename);
1667   if (file.IsNull()) file=Form("RefCalPads_%d.root",run);
1668   TDirectory *currDir=gDirectory;
1669   //db instance
1670   AliTPCcalibDB *db=AliTPCcalibDB::Instance();
1671   // retrieve cal pad objects
1672   db->SetRun(run);
1673   //open file
1674   TFile f(file.Data(),"recreate");
1675   //noise and pedestals
1676   db->GetPedestals()->Write("Pedestals");
1677   db->GetPadNoise()->Write("PadNoise");
1678   //pulser data
1679   db->GetPulserTmean()->Write("PulserTmean");
1680   db->GetPulserTrms()->Write("PulserTrms");
1681   db->GetPulserQmean()->Write("PulserQmean");
1682   //CE data
1683   db->GetCETmean()->Write("CETmean");
1684   db->GetCETrms()->Write("CETrms");
1685   db->GetCEQmean()->Write("CEQmean");
1686   //Altro data
1687   db->GetALTROAcqStart() ->Write("ALTROAcqStart");
1688   db->GetALTROZsThr()    ->Write("ALTROZsThr");
1689   db->GetALTROFPED()     ->Write("ALTROFPED");
1690   db->GetALTROAcqStop()  ->Write("ALTROAcqStop");
1691   db->GetALTROMasked()   ->Write("ALTROMasked");
1692   //
1693   f.Close();
1694   currDir->cd();
1695   return kTRUE;
1696 }
1697
1698
1699
1700 Double_t AliTPCcalibDB::GetVDriftCorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
1701   //
1702   // Get time dependent drift velocity correction
1703   // multiplication factor        vd = vdnom *(1+vdriftcorr)
1704   // Arguments:
1705   // mode determines the algorith how to combine the Laser Track, LaserCE and physics tracks
1706   // timestamp - timestamp
1707   // run       - run number
1708   // side      - the drift velocity per side (possible for laser and CE)
1709   //
1710   // Notice - Extrapolation outside of calibration range  - using constant function
1711   //
1712   Double_t result;
1713   // mode 1  automatic mode - according to the distance to the valid calibration
1714   //                        -  
1715   Double_t deltaP=0,  driftP=0,      wP  = 0.;
1716   Double_t deltaITS=0,driftITS=0,    wITS= 0.;
1717   Double_t deltaLT=0, driftLT=0,     wLT = 0.;
1718   Double_t deltaCE=0, driftCE=0,     wCE = 0.;
1719   driftP  = fDButil->GetVDriftTPC(deltaP,run,timeStamp); 
1720   driftITS= fDButil->GetVDriftTPCITS(deltaITS,run,timeStamp);
1721   driftCE = fDButil->GetVDriftTPCCE(deltaCE, run,timeStamp,36000,2);
1722   driftLT = fDButil->GetVDriftTPCLaserTracks(deltaLT,run,timeStamp,36000,2);
1723   deltaITS = TMath::Abs(deltaITS);
1724   deltaP   = TMath::Abs(deltaP);
1725   deltaLT  = TMath::Abs(deltaLT);
1726   deltaCE  = TMath::Abs(deltaCE);
1727   if (mode==1) {
1728     const Double_t kEpsilon=0.00000000001;
1729     const Double_t kdeltaT=360.; // 10 minutes
1730     wITS  = 64.*kdeltaT/(deltaITS +kdeltaT);
1731     wLT   = 16.*kdeltaT/(deltaLT  +kdeltaT);
1732     wP    = 0. *kdeltaT/(deltaP   +kdeltaT);
1733     wCE   = 1. *kdeltaT/(deltaCE  +kdeltaT);
1734     //
1735     //
1736     if (TMath::Abs(driftP)<kEpsilon)  wP=0;  // invalid calibration
1737     if (TMath::Abs(driftITS)<kEpsilon)wITS=0;  // invalid calibration
1738     if (TMath::Abs(driftLT)<kEpsilon) wLT=0;  // invalid calibration
1739     if (TMath::Abs(driftCE)<kEpsilon) wCE=0;  // invalid calibration
1740     if (wP+wITS+wLT+wCE<kEpsilon) return 0;
1741     result = (driftP*wP+driftITS*wITS+driftLT*wLT+driftCE*wCE)/(wP+wITS+wLT+wCE);
1742   }
1743
1744   return result;
1745 }
1746
1747 Double_t AliTPCcalibDB::GetTime0CorrectionTime(Int_t timeStamp, Int_t run, Int_t /*side*/, Int_t mode){
1748   //
1749   // Get time dependent time 0 (trigger delay in cm) correction
1750   // additive correction        time0 = time0+ GetTime0CorrectionTime
1751   // Value etracted combining the vdrift correction using laser tracks and CE and the physics track matchin
1752   // Arguments:
1753   // mode determines the algorith how to combine the Laser Track and physics tracks
1754   // timestamp - timestamp
1755   // run       - run number
1756   // side      - the drift velocity per side (possible for laser and CE)
1757   //
1758   // Notice - Extrapolation outside of calibration range  - using constant function
1759   //
1760   Double_t result=0;
1761   if (mode==2) {
1762     // TPC-TPC mode
1763     result=fDButil->GetTriggerOffsetTPC(run,timeStamp);    
1764     result  *=fParam->GetZLength();
1765   }
1766   if (mode==1){
1767     // TPC-ITS mode
1768     Double_t dist=0;
1769     result= -fDButil->GetTime0TPCITS(dist, run, timeStamp)*fParam->GetDriftV()/1000000.;
1770   }
1771   return result;
1772
1773 }
1774
1775
1776
1777
1778 Double_t AliTPCcalibDB::GetVDriftCorrectionGy(Int_t timeStamp, Int_t run, Int_t side, Int_t /*mode*/){
1779   //
1780   // Get global y correction drift velocity correction factor
1781   // additive factor        vd = vdnom*(1+GetVDriftCorrectionGy *gy)
1782   // Value etracted combining the vdrift correction using laser tracks and CE
1783   // Arguments:
1784   // mode determines the algorith how to combine the Laser Track, LaserCE
1785   // timestamp - timestamp
1786   // run       - run number
1787   // side      - the drift velocity gy correction per side (CE and Laser tracks)
1788   //
1789   // Notice - Extrapolation outside of calibration range  - using constant function
1790   // 
1791   if (run<=0 && fTransform) run = fTransform->GetCurrentRunNumber();
1792   UpdateRunInformations(run,kFALSE);
1793   TObjArray *array =AliTPCcalibDB::Instance()->GetTimeVdriftSplineRun(run);
1794   if (!array) return 0;
1795   TGraphErrors *laserA= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
1796   TGraphErrors *laserC= (TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
1797   
1798   Double_t result=0;
1799   if (laserA && laserC){
1800    result= (laserA->Eval(timeStamp)+laserC->Eval(timeStamp))*0.5;
1801   }
1802   if (laserA && side==0){
1803     result = (laserA->Eval(timeStamp));
1804   }
1805   if (laserC &&side==1){
1806     result = (laserC->Eval(timeStamp));
1807   }
1808   return -result/250.; //normalized before
1809 }
1810
1811
1812