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