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