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