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