]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibDB.cxx
Minor issues solved. Additional printouts for debugging purposes (H. Tydesjo)
[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 <AliLog.h>
86
87 #include "AliTPCcalibDB.h"
88 #include "AliTPCAltroMapping.h"
89 #include "AliTPCExB.h"
90
91 #include "AliTPCCalROC.h"
92 #include "AliTPCCalPad.h"
93 #include "AliTPCSensorTempArray.h"
94 #include "AliTPCTransform.h"
95
96 class AliCDBStorage;
97 class AliTPCCalDet;
98 //
99 //
100
101 #include "TFile.h"
102 #include "TKey.h"
103
104 #include "TObjArray.h"
105 #include "TObjString.h"
106 #include "TString.h"
107 #include "AliTPCCalPad.h"
108 #include "AliTPCCalibPulser.h"
109 #include "AliTPCCalibPedestal.h"
110 #include "AliTPCCalibCE.h"
111
112
113
114
115
116 ClassImp(AliTPCcalibDB)
117
118 AliTPCcalibDB* AliTPCcalibDB::fgInstance = 0;
119 Bool_t AliTPCcalibDB::fgTerminated = kFALSE;
120
121
122 //_ singleton implementation __________________________________________________
123 AliTPCcalibDB* AliTPCcalibDB::Instance()
124 {
125   //
126   // Singleton implementation
127   // Returns an instance of this class, it is created if neccessary
128   //
129   
130   if (fgTerminated != kFALSE)
131     return 0;
132
133   if (fgInstance == 0)
134     fgInstance = new AliTPCcalibDB();
135   
136   return fgInstance;
137 }
138
139 void AliTPCcalibDB::Terminate()
140 {
141   //
142   // Singleton implementation
143   // Deletes the instance of this class and sets the terminated flag, instances cannot be requested anymore
144   // This function can be called several times.
145   //
146   
147   fgTerminated = kTRUE;
148   
149   if (fgInstance != 0)
150   {
151     delete fgInstance;
152     fgInstance = 0;
153   }
154 }
155
156 //_____________________________________________________________________________
157 AliTPCcalibDB::AliTPCcalibDB():
158   fRun(-1),
159   fTransform(0),
160   fExB(0),
161   fPadGainFactor(0),
162   fPadTime0(0),
163   fPadNoise(0),
164   fPedestals(0),
165   fTemperature(0),
166   fMapping(0),
167   fRecoParamArray(0),
168   fParam(0),
169   fClusterParam(0)
170 {
171   //
172   // constructor
173   //  
174   //
175   Update();    // temporary
176 }
177
178 //_____________________________________________________________________________
179 AliTPCcalibDB::~AliTPCcalibDB() 
180 {
181   //
182   // destructor
183   //
184   
185   // don't delete anything, CDB cache is active!
186   //if (fPadGainFactor) delete fPadGainFactor;
187   //if (fPadTime0) delete fPadTime0;
188   //if (fPadNoise) delete fPadNoise;
189 }
190
191
192 //_____________________________________________________________________________
193 AliCDBEntry* AliTPCcalibDB::GetCDBEntry(const char* cdbPath)
194 {
195   // 
196   // Retrieves an entry with path <cdbPath> from the CDB.
197   //
198   char chinfo[1000];
199     
200   AliCDBEntry* entry = AliCDBManager::Instance()->Get(cdbPath, fRun); 
201   if (!entry) 
202   { 
203     sprintf(chinfo,"AliTPCcalibDB: Failed to get entry:\t%s ", cdbPath);
204     AliError(chinfo); 
205     return 0; 
206   }
207   return entry;
208 }
209
210
211 //_____________________________________________________________________________
212 void AliTPCcalibDB::SetRun(Long64_t run)
213 {
214   //
215   // Sets current run number. Calibration data is read from the corresponding file. 
216   //  
217   if (fRun == run)
218     return;  
219   fRun = run;
220   Update();
221 }
222   
223
224
225 void AliTPCcalibDB::Update(){
226   //
227   AliCDBEntry * entry=0;
228   
229   Bool_t cdbCache = AliCDBManager::Instance()->GetCacheFlag(); // save cache status
230   AliCDBManager::Instance()->SetCacheFlag(kTRUE); // activate CDB cache
231   
232   //
233   entry          = GetCDBEntry("TPC/Calib/PadGainFactor");
234   if (entry){
235     //if (fPadGainFactor) delete fPadGainFactor;
236     entry->SetOwner(kTRUE);
237     fPadGainFactor = (AliTPCCalPad*)entry->GetObject();
238   }
239   //
240   entry          = GetCDBEntry("TPC/Calib/PadTime0");
241   if (entry){
242     //if (fPadTime0) delete fPadTime0;
243     entry->SetOwner(kTRUE);
244     fPadTime0 = (AliTPCCalPad*)entry->GetObject();
245   }
246   //
247   //
248   entry          = GetCDBEntry("TPC/Calib/PadNoise");
249   if (entry){
250     //if (fPadNoise) delete fPadNoise;
251     entry->SetOwner(kTRUE);
252     fPadNoise = (AliTPCCalPad*)entry->GetObject();
253   }
254
255   entry          = GetCDBEntry("TPC/Calib/Pedestals");
256   if (entry){
257     //if (fPedestals) delete fPedestals;
258     entry->SetOwner(kTRUE);
259     fPedestals = (AliTPCCalPad*)entry->GetObject();
260   }
261
262   entry          = GetCDBEntry("TPC/Calib/Temperature");
263   if (entry){
264     //if (fTemperature) delete fTemperature;
265     entry->SetOwner(kTRUE);
266     fTemperature = (AliTPCSensorTempArray*)entry->GetObject();
267   }
268
269
270   entry          = GetCDBEntry("TPC/Calib/RecoParam");
271   if (entry){
272     entry->SetOwner(kTRUE);
273     fRecoParamArray = (TObjArray*)(entry->GetObject());
274   }
275
276
277   entry          = GetCDBEntry("TPC/Calib/Parameters");
278   if (entry){
279     //if (fPadNoise) delete fPadNoise;
280     entry->SetOwner(kTRUE);
281     fParam = (AliTPCParam*)(entry->GetObject()->Clone());
282   }
283
284   entry          = GetCDBEntry("TPC/Calib/ClusterParam");
285   if (entry){
286     //if (fPadNoise) delete fPadNoise;
287     entry->SetOwner(kTRUE);
288     fClusterParam = (AliTPCClusterParam*)(entry->GetObject()->Clone());
289   }
290
291   entry          = GetCDBEntry("TPC/Calib/Mapping");
292   if (entry){
293     //if (fPadNoise) delete fPadNoise;
294     entry->SetOwner(kTRUE);
295     TObjArray * array = dynamic_cast<TObjArray*>(entry->GetObject());
296     if (array && array->GetEntriesFast()==6){
297       fMapping = new AliTPCAltroMapping*[6];
298       for (Int_t i=0; i<6; i++){
299         fMapping[i] =  dynamic_cast<AliTPCAltroMapping*>(array->At(i));
300       }
301     }
302   }
303
304
305
306   entry          = GetCDBEntry("TPC/Calib/ExB");
307   if (entry) {
308     entry->SetOwner(kTRUE);
309     fExB=dynamic_cast<AliTPCExB*>(entry->GetObject()->Clone());
310   }
311
312   if (!fTransform) {
313     fTransform=new AliTPCTransform(); 
314   }
315
316   //
317   AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
318   
319 }
320 AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& org)
321 {
322   //
323   // Copy constructor invalid -- singleton implementation
324   //
325    Error("copy constructor","invalid -- singleton implementation");
326 }
327
328 AliTPCcalibDB& AliTPCcalibDB::operator= (const AliTPCcalibDB& rhs)
329 {
330 //
331 // Singleton implementation - no assignment operator
332 //
333   Error("operator =", "assignment operator not implemented");
334   return *this;
335 }
336
337
338
339 void AliTPCcalibDB::CreateObjectList(const Char_t *filename, TObjArray *calibObjects)
340 {
341 //
342 // Create calibration objects and read contents from OCDB
343 //
344    if ( calibObjects == 0x0 ) return;
345    ifstream in;
346    in.open(filename);
347    if ( !in.is_open() ){
348       fprintf(stderr,"Error: cannot open list file '%s'", filename);
349       return;
350    }
351    
352    AliTPCCalPad *calPad=0x0;
353    
354    TString sFile;
355    sFile.ReadFile(in);
356    in.close();
357    
358    TObjArray *arrFileLine = sFile.Tokenize("\n");
359    
360    TIter nextLine(arrFileLine);
361    
362    TObjString *sObjLine=0x0;
363    while ( sObjLine = (TObjString*)nextLine() ){
364       TString sLine(sObjLine->GetString());
365       
366       TObjArray *arrNextCol = sLine.Tokenize("\t");
367       
368       TObjString *sObjType     = (TObjString*)(arrNextCol->At(0));
369       TObjString *sObjFileName = (TObjString*)(arrNextCol->At(1));
370       
371       if ( !sObjType || ! sObjFileName ) continue;
372       TString sType(sObjType->GetString());
373       TString sFileName(sObjFileName->GetString());
374       printf("%s\t%s\n",sType.Data(),sFileName.Data());
375       
376       TFile *fIn = TFile::Open(sFileName);
377       if ( !fIn ){
378          fprintf(stderr,"File not found: '%s'", sFileName.Data());
379          continue;
380       }
381       
382       if ( sType == "CE" ){
383          AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
384          
385          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());         
386          calPad->SetNameTitle("CETmean","CETmean");
387          calibObjects->Add(calPad);
388          
389          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());         
390          calPad->SetNameTitle("CEQmean","CEQmean");
391          calibObjects->Add(calPad);        
392          
393          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
394          calPad->SetNameTitle("CETrms","CETrms");
395          calibObjects->Add(calPad);         
396                   
397       } else if ( sType == "Pulser") {
398          AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
399          
400          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());         
401          calPad->SetNameTitle("PulserTmean","PulserTmean");
402          calibObjects->Add(calPad);
403          
404          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());         
405          calPad->SetNameTitle("PulserQmean","PulserQmean");
406          calibObjects->Add(calPad);        
407          
408          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
409          calPad->SetNameTitle("PulserTrms","PulserTrms");
410          calibObjects->Add(calPad);         
411       
412       } else if ( sType == "Pedestals") {
413          AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
414          
415          calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());         
416          calPad->SetNameTitle("Pedestals","Pedestals");
417          calibObjects->Add(calPad);
418          
419          calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());         
420          calPad->SetNameTitle("Noise","Noise");
421          calibObjects->Add(calPad);        
422      
423       } else {
424          fprintf(stderr,"Undefined Type: '%s'",sType.Data());
425          
426       }
427       delete fIn;
428    }
429 }
430
431
432
433 void AliTPCcalibDB::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
434   //
435   // Write a tree with all available information
436   // if mapFileName is specified, the Map information are also written to the tree
437   // pads specified in outlierPad are not used for calculating statistics
438   //  - the same function as AliTPCCalPad::MakeTree - 
439   //
440    AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
441
442    TObjArray* mapIROCs = 0;
443    TObjArray* mapOROCs = 0;
444    TVectorF *mapIROCArray = 0;
445    TVectorF *mapOROCArray = 0;
446    Int_t mapEntries = 0;
447    TString* mapNames = 0;
448    
449    if (mapFileName) {
450       TFile mapFile(mapFileName, "read");
451       
452       TList* listOfROCs = mapFile.GetListOfKeys();
453       mapEntries = listOfROCs->GetEntries()/2;
454       mapIROCs = new TObjArray(mapEntries*2);
455       mapOROCs = new TObjArray(mapEntries*2);
456       mapIROCArray = new TVectorF[mapEntries];
457       mapOROCArray = new TVectorF[mapEntries];
458       
459       mapNames = new TString[mapEntries];
460       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
461         TString nameROC(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
462          nameROC.Remove(nameROC.Length()-4, 4);
463          mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "IROC").Data()), ivalue);
464          mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "OROC").Data()), ivalue);
465          mapNames[ivalue].Append(nameROC);
466       }
467       
468       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
469          mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
470          mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
471       
472          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
473             (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
474          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
475             (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
476       }
477
478    } //  if (mapFileName)
479   
480    TTreeSRedirector cstream(fileName);
481    Int_t arrayEntries = array->GetEntries();
482    
483    TString* names = new TString[arrayEntries];
484    for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
485       names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
486
487    for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
488       //
489       // get statistic for given sector
490       //
491       TVectorF median(arrayEntries);
492       TVectorF mean(arrayEntries);
493       TVectorF rms(arrayEntries);
494       TVectorF ltm(arrayEntries);
495       TVectorF ltmrms(arrayEntries);
496       TVectorF medianWithOut(arrayEntries);
497       TVectorF meanWithOut(arrayEntries);
498       TVectorF rmsWithOut(arrayEntries);
499       TVectorF ltmWithOut(arrayEntries);
500       TVectorF ltmrmsWithOut(arrayEntries);
501       
502       TVectorF *vectorArray = new TVectorF[arrayEntries];
503       for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
504          vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
505       
506       for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
507          AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
508          AliTPCCalROC* calROC = calPad->GetCalROC(isector);
509          AliTPCCalROC* outlierROC = 0;
510          if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
511          if (calROC) {
512             median[ivalue] = calROC->GetMedian();
513             mean[ivalue] = calROC->GetMean();
514             rms[ivalue] = calROC->GetRMS();
515             Double_t ltmrmsValue = 0;
516             ltm[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction);
517             ltmrms[ivalue] = ltmrmsValue;
518             if (outlierROC) {
519                medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
520                meanWithOut[ivalue] = calROC->GetMean(outlierROC);
521                rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
522                ltmrmsValue = 0;
523                ltmWithOut[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction, outlierROC);
524                ltmrmsWithOut[ivalue] = ltmrmsValue;
525             }
526          }
527          else {
528             median[ivalue] = 0.;
529             mean[ivalue] = 0.;
530             rms[ivalue] = 0.;
531             ltm[ivalue] = 0.;
532             ltmrms[ivalue] = 0.;
533             medianWithOut[ivalue] = 0.;
534             meanWithOut[ivalue] = 0.;
535             rmsWithOut[ivalue] = 0.;
536             ltmWithOut[ivalue] = 0.;
537             ltmrmsWithOut[ivalue] = 0.;
538          }
539       }
540       
541       //
542       // fill vectors of variable per pad
543       //
544       TVectorF *posArray = new TVectorF[8];
545       for (Int_t ivalue = 0; ivalue < 8; ivalue++)
546          posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
547
548       Float_t posG[3] = {0};
549       Float_t posL[3] = {0};
550       Int_t ichannel = 0;
551       for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
552          for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
553             tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
554             tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
555             posArray[0][ichannel] = irow;
556             posArray[1][ichannel] = ipad;
557             posArray[2][ichannel] = posL[0];
558             posArray[3][ichannel] = posL[1];
559             posArray[4][ichannel] = posG[0];
560             posArray[5][ichannel] = posG[1];
561             posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
562             posArray[7][ichannel] = ichannel;
563             
564             // loop over array containing AliTPCCalPads
565             for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
566                AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
567                AliTPCCalROC* calROC = calPad->GetCalROC(isector);
568                if (calROC)
569                   (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
570                else
571                   (vectorArray[ivalue])[ichannel] = 0;
572             }
573             ichannel++;
574          }
575       }
576       
577       cstream << "calPads" <<
578          "sector=" << isector;
579       
580       for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
581          cstream << "calPads" <<
582             (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
583             (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
584             (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
585             (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
586             (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
587          if (outlierPad) {
588             cstream << "calPads" <<
589                (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
590                (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
591                (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
592                (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
593                (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
594          }
595       }
596
597       for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
598          cstream << "calPads" <<
599             (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
600       }
601
602       if (mapFileName) {
603          for  (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
604             if (isector < 36)
605                cstream << "calPads" <<
606                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
607             else
608                cstream << "calPads" <<
609                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
610          }
611       }
612
613       cstream << "calPads" <<
614          "row.=" << &posArray[0] <<
615          "pad.=" << &posArray[1] <<
616          "lx.=" << &posArray[2] <<
617          "ly.=" << &posArray[3] <<
618          "gx.=" << &posArray[4] <<
619          "gy.=" << &posArray[5] <<
620          "rpad.=" << &posArray[6] <<
621          "channel.=" << &posArray[7];
622          
623       cstream << "calPads" <<
624          "\n";
625
626       delete[] posArray;
627       delete[] vectorArray;
628    }
629    
630
631    delete[] names;
632    if (mapFileName) {
633       delete mapIROCs;
634       delete mapOROCs;
635       delete[] mapIROCArray;
636       delete[] mapOROCArray;
637       delete[] mapNames;
638    }
639 }
640
641
642 AliTPCRecoParam *  AliTPCcalibDB::GetRecoParam(Int_t */*eventtype*/){
643   //
644   // 
645   //
646   if (!fRecoParamArray){
647     return 0;  // back compatible sollution
648   };
649   
650   AliTPCRecoParam * param = (AliTPCRecoParam*)fRecoParamArray->At(0);
651   return param;
652
653 }