]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCalPad.cxx
The AliDCSSensor classes were recently upgraded to include start and end time entries...
[u/mrichter/AliRoot.git] / TPC / AliTPCCalPad.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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 //  TPC calibration class for parameters which saved per pad                 //
21 //                                                                           //
22 ///////////////////////////////////////////////////////////////////////////////
23
24 #include "AliTPCCalPad.h"
25 #include "AliTPCCalROC.h"
26 #include <TObjArray.h>
27 #include <TAxis.h>
28 #include <TGraph.h>
29 #include <TGraph2D.h>
30 #include <TH2F.h>
31 #include "TTreeStream.h"
32 #include "TFile.h"
33 #include "TKey.h"
34
35 ClassImp(AliTPCCalPad)
36
37 //_____________________________________________________________________________
38 AliTPCCalPad::AliTPCCalPad():TNamed()
39 {
40   //
41   // AliTPCCalPad default constructor
42   //
43
44   for (Int_t isec = 0; isec < kNsec; isec++) {
45     fROC[isec] = 0;
46   }
47
48 }
49
50 //_____________________________________________________________________________
51 AliTPCCalPad::AliTPCCalPad(const Text_t *name, const Text_t *title)
52                 :TNamed(name,title)
53 {
54   //
55   // AliTPCCalPad constructor
56   //
57   for (Int_t isec = 0; isec < kNsec; isec++) {
58     fROC[isec] = new AliTPCCalROC(isec);
59   }
60 }
61
62
63 //_____________________________________________________________________________
64 AliTPCCalPad::AliTPCCalPad(const AliTPCCalPad &c):TNamed(c)
65 {
66   //
67   // AliTPCCalPad copy constructor
68   //
69
70   for (Int_t isec = 0; isec < kNsec; isec++) {
71          fROC[isec] = 0;
72      if (c.fROC[isec])
73        fROC[isec] = new AliTPCCalROC(*(c.fROC[isec]));
74   }
75 }
76
77 //_____________________________________________________________________________
78 AliTPCCalPad::AliTPCCalPad(TObjArray * array):TNamed()
79 {
80   //
81   // AliTPCCalPad default constructor
82   //
83
84   for (Int_t isec = 0; isec < kNsec; isec++) {
85     fROC[isec] = (AliTPCCalROC *)array->At(isec);
86   }
87
88 }
89
90
91 ///_____________________________________________________________________________
92 AliTPCCalPad::~AliTPCCalPad()
93 {
94   //
95   // AliTPCCalPad destructor
96   //
97
98   for (Int_t isec = 0; isec < kNsec; isec++) {
99     if (fROC[isec]) {
100       delete fROC[isec];
101       fROC[isec] = 0;
102     }
103   }
104
105 }
106
107 //_____________________________________________________________________________
108 AliTPCCalPad &AliTPCCalPad::operator=(const AliTPCCalPad &c)
109 {
110   //
111   // Assignment operator
112   //
113
114   if (this != &c) ((AliTPCCalPad &) c).Copy(*this);
115   return *this;
116
117 }
118
119 //_____________________________________________________________________________
120 void AliTPCCalPad::Copy(TObject &c) const
121 {
122   //
123   // Copy function
124   //
125
126   for (Int_t isec = 0; isec < kNsec; isec++) {
127     if (fROC[isec]) {
128       fROC[isec]->Copy(*((AliTPCCalPad &) c).fROC[isec]);
129     }
130   }
131   TObject::Copy(c);
132 }
133
134 //_____________________________________________________________________________
135 void AliTPCCalPad::Add(Float_t c1)
136 {
137     //
138     // add constant for all channels of all ROCs
139     //
140
141     for (Int_t isec = 0; isec < kNsec; isec++) {
142         if (fROC[isec]){
143             fROC[isec]->Add(c1);
144         }
145     }
146 }
147
148 //_____________________________________________________________________________
149 void AliTPCCalPad::Multiply(Float_t c1)
150 {
151     //
152     // multiply constant for all channels of all ROCs
153     //
154     for (Int_t isec = 0; isec < kNsec; isec++) {
155         if (fROC[isec]){
156             fROC[isec]->Multiply(c1);
157         }
158     }
159 }
160
161 //_____________________________________________________________________________
162 void AliTPCCalPad::Add(const AliTPCCalPad * pad, Double_t c1)
163 {
164     //
165     // add calpad channel by channel multiplied by c1 - all ROCs
166     //
167     for (Int_t isec = 0; isec < kNsec; isec++) {
168         if (fROC[isec]){
169             fROC[isec]->Add(pad->GetCalROC(isec),c1);
170         }
171     }
172 }
173
174 //_____________________________________________________________________________
175 void AliTPCCalPad::Multiply(const AliTPCCalPad * pad)
176 {
177     //
178     // multiply calpad channel by channel - all ROCs
179     //
180     for (Int_t isec = 0; isec < kNsec; isec++) {
181         if (fROC[isec]){
182             fROC[isec]->Multiply(pad->GetCalROC(isec));
183         }
184     }
185 }
186
187 //_____________________________________________________________________________
188 void AliTPCCalPad::Divide(const AliTPCCalPad * pad)
189 {
190     //
191     // divide calpad channel by channel - all ROCs
192     //
193     for (Int_t isec = 0; isec < kNsec; isec++) {
194         if (fROC[isec]){
195             fROC[isec]->Divide(pad->GetCalROC(isec));
196         }
197     }
198 }
199
200 //_____________________________________________________________________________
201 TGraph  *  AliTPCCalPad::MakeGraph(Int_t type, Float_t ratio){
202   //
203   //   type=1 - mean
204   //        2 - median
205   //        3 - LTM
206   Int_t npoints = 0;
207   for (Int_t i=0;i<72;i++) if (fROC[i]) npoints++;
208   TGraph * graph = new TGraph(npoints);
209   npoints=0;   
210   for (Int_t isec=0;isec<72;isec++){
211     if (!fROC[isec]) continue;
212     if (type==0)  graph->SetPoint(npoints,isec,fROC[isec]->GetMean());      
213     if (type==1)  graph->SetPoint(npoints,isec,fROC[isec]->GetMedian());
214     if (type==2)  graph->SetPoint(npoints,isec,fROC[isec]->GetLTM(0,ratio));    
215     npoints++;
216   }
217
218   graph->GetXaxis()->SetTitle("Sector"); 
219   if (type==0) {
220     graph->GetYaxis()->SetTitle("Mean");   
221     graph->SetMarkerStyle(22);    
222   }
223   if (type==1) {
224     graph->GetYaxis()->SetTitle("Median");   
225     graph->SetMarkerStyle(22);    
226   }
227   if (type==2) {
228       graph->GetYaxis()->SetTitle(Form("Mean%f",ratio));      
229       graph->SetMarkerStyle(24);
230   }
231
232   return graph;
233 }
234
235 //_____________________________________________________________________________
236 Double_t AliTPCCalPad::GetMeanRMS(Double_t &rms)
237 {
238     //
239     // Calculate mean an RMS of all rocs
240     //
241     Double_t sum = 0, sum2 = 0, n=0, val=0;
242     for (Int_t isec = 0; isec < kNsec; isec++) {
243         AliTPCCalROC *calRoc = fROC[isec];
244         if ( calRoc ){
245             for (UInt_t irow=0; irow<calRoc->GetNrows(); irow++){
246                 for (UInt_t ipad=0; ipad<calRoc->GetNPads(irow); ipad++){
247                     val = calRoc->GetValue(irow,ipad);
248                     sum+=val;
249                     sum2+=val*val;
250                     n++;
251                 }
252             }
253
254         }
255     }
256     Double_t n1 = 1./n;
257     Double_t mean = sum*n1;
258     rms  = TMath::Sqrt(TMath::Abs(sum2*n1-mean*mean));
259     return mean;
260 }
261
262
263 //_____________________________________________________________________________
264 Double_t AliTPCCalPad::GetMean(AliTPCCalPad* outlierPad)
265 {
266     //
267     // return mean of the mean of all ROCs
268     //
269     Double_t arr[kNsec];
270     Int_t n=0;
271     for (Int_t isec = 0; isec < kNsec; isec++) {
272        AliTPCCalROC *calRoc = fROC[isec];
273        if ( calRoc ){
274           AliTPCCalROC* outlierROC = 0;
275           if (outlierPad) outlierROC = outlierPad->GetCalROC(isec);
276                arr[n] = calRoc->GetMean(outlierROC);
277           n++;
278        }
279     }
280     return TMath::Mean(n,arr);
281 }
282
283 //_____________________________________________________________________________
284 Double_t AliTPCCalPad::GetRMS(AliTPCCalPad* outlierPad)
285 {
286     //
287     // return mean of the RMS of all ROCs
288     //
289     Double_t arr[kNsec];
290     Int_t n=0;
291     for (Int_t isec = 0; isec < kNsec; isec++) {
292        AliTPCCalROC *calRoc = fROC[isec];
293        if ( calRoc ){
294           AliTPCCalROC* outlierROC = 0;
295           if (outlierPad) outlierROC = outlierPad->GetCalROC(isec);
296           arr[n] = calRoc->GetRMS(outlierROC);
297           n++;
298        }
299     }
300     return TMath::Mean(n,arr);
301 }
302
303 //_____________________________________________________________________________
304 Double_t AliTPCCalPad::GetMedian(AliTPCCalPad* outlierPad)
305 {
306     //
307     // return mean of the median of all ROCs
308     //
309     Double_t arr[kNsec];
310     Int_t n=0;
311     for (Int_t isec = 0; isec < kNsec; isec++) {
312        AliTPCCalROC *calRoc = fROC[isec];
313        if ( calRoc ){
314           AliTPCCalROC* outlierROC = 0;
315           if (outlierPad) outlierROC = outlierPad->GetCalROC(isec);
316           arr[n] = calRoc->GetMedian(outlierROC);
317           n++;
318        }
319     }
320     return TMath::Mean(n,arr);
321 }
322
323 //_____________________________________________________________________________
324 Double_t AliTPCCalPad::GetLTM(Double_t *sigma, Double_t fraction, AliTPCCalPad* outlierPad)
325 {
326     //
327     // return mean of the LTM and sigma of all ROCs
328     //
329     Double_t arrm[kNsec];
330     Double_t arrs[kNsec];
331     Double_t *sTemp=0x0;
332     Int_t n=0;
333
334     for (Int_t isec = 0; isec < kNsec; isec++) {
335         AliTPCCalROC *calRoc = fROC[isec];
336         if ( calRoc ){
337             if ( sigma ) sTemp=arrs+n;
338        AliTPCCalROC* outlierROC = 0;
339        if (outlierPad) outlierROC = outlierPad->GetCalROC(isec);
340             arrm[n] = calRoc->GetLTM(sTemp,fraction, outlierROC);
341             n++;
342         }
343     }
344     if ( sigma ) *sigma = TMath::Mean(n,arrs);
345     return TMath::Mean(n,arrm);
346 }
347
348 //_____________________________________________________________________________
349 TH1F * AliTPCCalPad::MakeHisto1D(Float_t min, Float_t max,Int_t type){
350   //
351   // make 1D histo
352   // type -1 = user defined range
353   //       0 = nsigma cut nsigma=min
354   if (type>=0){
355     if (type==0){
356       // nsigma range
357       Float_t mean  = GetMean();
358       Float_t sigma = GetRMS();
359       Float_t nsigma = TMath::Abs(min);
360       min = mean-nsigma*sigma;
361       max = mean+nsigma*sigma;
362     }
363     if (type==1){
364       // fixed range
365       Float_t mean   = GetMedian();
366       Float_t  delta = min;
367       min = mean-delta;
368       max = mean+delta;
369     }
370     if (type==2){
371       //
372       // LTM mean +- nsigma
373       //
374       Double_t sigma;
375       Float_t mean  = GetLTM(&sigma,max);
376       sigma*=min;
377       min = mean-sigma;
378       max = mean+sigma;
379     }
380   }
381   char  name[1000];
382   sprintf(name,"%s Pad 1D",GetTitle());
383   TH1F * his = new TH1F(name,name,100, min,max);
384     for (Int_t isec = 0; isec < kNsec; isec++) {
385         if (fROC[isec]){
386             for (UInt_t irow=0; irow<fROC[isec]->GetNrows(); irow++){
387                 UInt_t npads = (Int_t)fROC[isec]->GetNPads(irow);
388                 for (UInt_t ipad=0; ipad<npads; ipad++){
389                     his->Fill(fROC[isec]->GetValue(irow,ipad));
390                 }
391             }
392         }
393     }
394   return his;
395 }
396
397 //_____________________________________________________________________________
398 TH2F *AliTPCCalPad::MakeHisto2D(Int_t side){
399   //
400   // Make 2D graph
401   // side  -  specify the side A = 0 C = 1
402   // type  -  used types of determination of boundaries in z
403   Float_t kEpsilon = 0.000000000001;
404   TH2F * his = new TH2F(GetName(), GetName(), 250,-250,250,250,-250,250);
405   AliTPCROC * roc  = AliTPCROC::Instance(); 
406   for (Int_t isec=0; isec<72; isec++){
407     if (side==0 && isec%36>=18) continue;
408     if (side>0 && isec%36<18) continue;
409     if (fROC[isec]){
410       AliTPCCalROC * calRoc = fROC[isec];
411       for (UInt_t irow=0; irow<calRoc->GetNrows(); irow++)
412         for (UInt_t ipad=0; ipad<calRoc->GetNPads(irow); ipad++)
413           if (TMath::Abs(calRoc->GetValue(irow,ipad))>kEpsilon){
414             Float_t xyz[3];
415             roc->GetPositionGlobal(isec,irow,ipad,xyz);
416             Int_t binx = 1+TMath::Nint((xyz[0]+250.)*0.5);
417             Int_t biny = 1+TMath::Nint((xyz[1]+250.)*0.5);
418             Float_t value = calRoc->GetValue(irow,ipad);            
419             his->SetBinContent(binx,biny,value);
420           }
421     }
422   }
423   his->SetXTitle("x (cm)");
424   his->SetYTitle("y (cm)");
425   return his;
426 }
427
428
429
430
431 void AliTPCCalPad::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
432   //
433   // Write tree with all available information
434   //
435    AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
436
437    TObjArray* mapIROCs = 0;
438    TObjArray* mapOROCs = 0;
439    TVectorF *mapIROCArray = 0;
440    TVectorF *mapOROCArray = 0;
441    Int_t mapEntries = 0;
442    TString* mapNames = 0;
443    
444    if (mapFileName) {
445       TFile mapFile(mapFileName, "read");
446       
447       TList* listOfROCs = mapFile.GetListOfKeys();
448       mapEntries = listOfROCs->GetEntries()/2;
449       mapIROCs = new TObjArray(mapEntries*2);
450       mapOROCs = new TObjArray(mapEntries*2);
451       mapIROCArray = new TVectorF[mapEntries];
452       mapOROCArray = new TVectorF[mapEntries];
453       
454       mapNames = new TString[mapEntries];
455       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
456          TString ROCname(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
457          ROCname.Remove(ROCname.Length()-4, 4);
458          mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((ROCname + "IROC").Data()), ivalue);
459          mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((ROCname + "OROC").Data()), ivalue);
460          mapNames[ivalue].Append(ROCname);
461       }
462       
463       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
464          mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
465          mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
466       
467          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
468             (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
469          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
470             (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
471       }
472
473    } //  if (mapFileName)
474   
475    TTreeSRedirector cstream(fileName);
476    Int_t arrayEntries = array->GetEntries();
477    
478    TString* names = new TString[arrayEntries];
479    for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
480       names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
481
482    for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
483       //
484       // get statistic for given sector
485       //
486       TVectorF median(arrayEntries);
487       TVectorF mean(arrayEntries);
488       TVectorF rms(arrayEntries);
489       TVectorF ltm(arrayEntries);
490       TVectorF ltmrms(arrayEntries);
491       TVectorF medianWithOut(arrayEntries);
492       TVectorF meanWithOut(arrayEntries);
493       TVectorF rmsWithOut(arrayEntries);
494       TVectorF ltmWithOut(arrayEntries);
495       TVectorF ltmrmsWithOut(arrayEntries);
496       
497       TVectorF *vectorArray = new TVectorF[arrayEntries];
498       for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
499          vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
500       
501       for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
502          AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
503          AliTPCCalROC* calROC = calPad->GetCalROC(isector);
504          AliTPCCalROC* outlierROC = 0;
505          if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
506          if (calROC) {
507             median[ivalue] = calROC->GetMedian();
508             mean[ivalue] = calROC->GetMean();
509             rms[ivalue] = calROC->GetRMS();
510             Double_t ltmrmsValue = 0;
511             ltm[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction);
512             ltmrms[ivalue] = ltmrmsValue;
513             if (outlierROC) {
514                medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
515                meanWithOut[ivalue] = calROC->GetMean(outlierROC);
516                rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
517                ltmrmsValue = 0;
518                ltmWithOut[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction, outlierROC);
519                ltmrmsWithOut[ivalue] = ltmrmsValue;
520             }
521          }
522          else {
523             median[ivalue] = 0.;
524             mean[ivalue] = 0.;
525             rms[ivalue] = 0.;
526             ltm[ivalue] = 0.;
527             ltmrms[ivalue] = 0.;
528             medianWithOut[ivalue] = 0.;
529             meanWithOut[ivalue] = 0.;
530             rmsWithOut[ivalue] = 0.;
531             ltmWithOut[ivalue] = 0.;
532             ltmrmsWithOut[ivalue] = 0.;
533          }
534       }
535       
536       //
537       // fill vectors of variable per pad
538       //
539       TVectorF *posArray = new TVectorF[8];
540       for (Int_t ivalue = 0; ivalue < 8; ivalue++)
541          posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
542
543       Float_t posG[3] = {0};
544       Float_t posL[3] = {0};
545       Int_t ichannel = 0;
546       for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
547          for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
548             tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
549             tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
550             posArray[0][ichannel] = irow;
551             posArray[1][ichannel] = ipad;
552             posArray[2][ichannel] = posL[0];
553             posArray[3][ichannel] = posL[1];
554             posArray[4][ichannel] = posG[0];
555             posArray[5][ichannel] = posG[1];
556             posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
557             posArray[7][ichannel] = ichannel;
558             
559             // loop over array containing AliTPCCalPads
560             for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
561                AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
562                AliTPCCalROC* calROC = calPad->GetCalROC(isector);
563                if (calROC)
564                   (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
565                else
566                   (vectorArray[ivalue])[ichannel] = 0;
567             }
568             ichannel++;
569          }
570       }
571       
572       cstream << "calPads" <<
573          "sector=" << isector;
574       
575       for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
576          cstream << "calPads" <<
577             (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
578             (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
579             (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
580             (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
581             (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
582          if (outlierPad) {
583             cstream << "calPads" <<
584                (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
585                (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
586                (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
587                (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
588                (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
589          }
590       }
591
592       for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
593          cstream << "calPads" <<
594             (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
595       }
596
597       if (mapFileName) {
598          for  (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
599             if (isector < 36)
600                cstream << "calPads" <<
601                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
602             else
603                cstream << "calPads" <<
604                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
605          }
606       }
607
608       cstream << "calPads" <<
609          "row.=" << &posArray[0] <<
610          "pad.=" << &posArray[1] <<
611          "lx.=" << &posArray[2] <<
612          "ly.=" << &posArray[3] <<
613          "gx.=" << &posArray[4] <<
614          "gy.=" << &posArray[5] <<
615          "rpad.=" << &posArray[6] <<
616          "channel.=" << &posArray[7];
617          
618       cstream << "calPads" <<
619          "\n";
620
621       delete[] posArray;
622       delete[] vectorArray;
623    }
624    
625    delete[] names;
626    if (mapFileName) {
627       delete mapIROCs;
628       delete mapOROCs;
629       delete[] mapIROCArray;
630       delete[] mapOROCArray;
631       delete[] mapNames;
632    }
633 }
634
635