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