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