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