]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSPlaneEffSPD.cxx
Restore backward compatibility for SDD map objects - Savannah bug 49895 (F. Prino)
[u/mrichter/AliRoot.git] / ITS / AliITSPlaneEffSPD.cxx
1 /**************************************************************************
2  * Copyright(c) 2007-2009, 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 //  Plane Efficiency class for ITS                      
17 //  It is used for chip by chip efficiency of the SPD,        
18 //  evaluated by tracks
19 //  (Inherits from AliITSPlaneEff)
20 //  Author: G.E. Bruno 
21 //          giuseppe.bruno@ba.infn.it
22 //
23 ///////////////////////////////////////////////////////////////////////////
24
25 /* $Id$ */
26
27 #include <TMath.h>
28 #include <TH1F.h>
29 #include <TFile.h>
30 #include <TTree.h>
31 #include <TROOT.h>
32 #include "AliITSPlaneEffSPD.h"
33 #include "AliLog.h"
34 #include "AliCDBStorage.h"
35 #include "AliCDBEntry.h"
36 #include "AliCDBManager.h"
37 //#include "AliCDBRunRange.h"
38 #include "AliITSsegmentationSPD.h"
39 #include "AliITSCalibrationSPD.h"
40
41 ClassImp(AliITSPlaneEffSPD)     
42 //______________________________________________________________________
43 AliITSPlaneEffSPD::AliITSPlaneEffSPD():
44   AliITSPlaneEff(),
45   fHisResX(0),
46   fHisResZ(0),
47   fHisResXZ(0),
48   fHisClusterSize(0),
49   fHisResXclu(0),
50   fHisResZclu(0),
51   fHisResXchip(0),
52   fHisResZchip(0),
53   fHisTrackErrX(0),
54   fHisTrackErrZ(0),
55   fHisClusErrX(0),
56   fHisClusErrZ(0){
57   for (UInt_t i=0; i<kNModule*kNChip; i++){
58     fFound[i]=0;
59     fTried[i]=0;
60   }
61   // default constructor
62   AliDebug(1,Form("Calling default constructor"));
63 }
64 //______________________________________________________________________
65 AliITSPlaneEffSPD::~AliITSPlaneEffSPD(){
66     // destructor
67     // Inputs:
68     //    none.
69     // Outputs:
70     //    none.
71     // Return:
72     //     none.
73     DeleteHistos();
74 }
75 //______________________________________________________________________
76 AliITSPlaneEffSPD::AliITSPlaneEffSPD(const AliITSPlaneEffSPD &s) : AliITSPlaneEff(s), 
77 //fHis(s.fHis),
78 fHisResX(0),
79 fHisResZ(0),
80 fHisResXZ(0),
81 fHisClusterSize(0),
82 fHisResXclu(0),
83 fHisResZclu(0),
84 fHisResXchip(0),
85 fHisResZchip(0),
86 fHisTrackErrX(0),
87 fHisTrackErrZ(0),
88 fHisClusErrX(0),
89 fHisClusErrZ(0)
90 {
91     //     Copy Constructor
92     // Inputs:
93     //    AliITSPlaneEffSPD &s The original class for which
94     //                                this class is a copy of
95     // Outputs:
96     //    none.
97     // Return:
98
99  for (UInt_t i=0; i<kNModule*kNChip; i++){
100     fFound[i]=s.fFound[i];
101     fTried[i]=s.fTried[i];
102  }
103  if(fHis) { 
104    InitHistos();
105    for(Int_t i=0; i<kNHisto; i++) {
106       s.fHisResX[i]->Copy(*fHisResX[i]);
107       s.fHisResZ[i]->Copy(*fHisResZ[i]);
108       s.fHisResXZ[i]->Copy(*fHisResXZ[i]);
109       s.fHisClusterSize[i]->Copy(*fHisClusterSize[i]);
110       for(Int_t clu=0; clu<kNclu; clu++) {  // clu=0 --> cluster size 1
111         s.fHisResXclu[i][clu]->Copy(*fHisResXclu[i][clu]);
112         s.fHisResZclu[i][clu]->Copy(*fHisResZclu[i][clu]);
113       }
114       for(Int_t chip=0; chip<kNChip; chip++) { 
115         s.fHisResXchip[i][chip]->Copy(*fHisResXchip[i][chip]);
116         s.fHisResZchip[i][chip]->Copy(*fHisResZchip[i][chip]);
117       }
118       s.fHisTrackErrX[i]->Copy(*fHisTrackErrX[i]);
119       s.fHisTrackErrZ[i]->Copy(*fHisTrackErrZ[i]);
120       s.fHisClusErrX[i]->Copy(*fHisClusErrX[i]);
121       s.fHisClusErrZ[i]->Copy(*fHisClusErrZ[i]);
122    }
123  }
124 }
125 //_________________________________________________________________________
126 AliITSPlaneEffSPD& AliITSPlaneEffSPD::operator+=(const AliITSPlaneEffSPD &add){
127     //    Add-to-me operator
128     // Inputs:
129     //    const AliITSPlaneEffSPD &add  simulation class to be added
130     // Outputs:
131     //    none.
132     // Return:
133     //    none
134     for (UInt_t i=0; i<kNModule*kNChip; i++){
135       fFound[i] += add.fFound[i];
136       fTried[i] += add.fTried[i];
137     }
138     if(fHis && add.fHis) {
139       for(Int_t i=0; i<kNHisto; i++) {
140         fHisResX[i]->Add(add.fHisResX[i]); 
141         fHisResZ[i]->Add(add.fHisResZ[i]); 
142         fHisResXZ[i]->Add(add.fHisResXZ[i]); 
143         fHisClusterSize[i]->Add(add.fHisClusterSize[i]); 
144         for(Int_t clu=0; clu<kNclu; clu++) {  // clu=0 --> cluster size 1
145           fHisResXclu[i][clu]->Add(add.fHisResXclu[i][clu]); 
146           fHisResZclu[i][clu]->Add(add.fHisResZclu[i][clu]); 
147         }
148         for(Int_t chip=0; chip<kNChip; chip++) {  
149           fHisResXchip[i][chip]->Add(add.fHisResXchip[i][chip]); 
150           fHisResZchip[i][chip]->Add(add.fHisResZchip[i][chip]); 
151         }
152         fHisTrackErrX[i]->Add(add.fHisTrackErrX[i]);
153         fHisTrackErrZ[i]->Add(add.fHisTrackErrZ[i]);
154         fHisClusErrX[i]->Add(add.fHisClusErrX[i]);
155         fHisClusErrZ[i]->Add(add.fHisClusErrZ[i]);
156       }
157     }
158     return *this;
159 }
160 //______________________________________________________________________
161 AliITSPlaneEffSPD&  AliITSPlaneEffSPD::operator=(const
162                                            AliITSPlaneEffSPD &s){
163     //    Assignment operator
164     // Inputs:
165     //    AliITSPlaneEffSPD &s The original class for which
166     //                                this class is a copy of
167     // Outputs:
168     //    none.
169     // Return:
170  
171     if(this==&s) return *this;
172     s.Copy(*this);
173     return *this;
174 }
175 //______________________________________________________________________
176 void AliITSPlaneEffSPD::Copy(TObject &obj) const {
177   // protected method. copy this to obj
178   AliITSPlaneEff::Copy(obj);
179   AliITSPlaneEffSPD& target = (AliITSPlaneEffSPD &) obj;
180   for(Int_t i=0;i<kNModule*kNChip;i++) {
181       target.fFound[i] = fFound[i];
182       target.fTried[i] = fTried[i];
183   }
184   CopyHistos(target);
185   return;
186 }
187 //_______________________________________________________________________
188 void AliITSPlaneEffSPD::CopyHistos(AliITSPlaneEffSPD &target) const {
189   // protected method: copy histos from this to target
190   target.fHis  = fHis; // this is redundant only in some cases. Leave as it is.
191   if(fHis) {
192     target.fHisResX=new TH1F*[kNHisto];
193     target.fHisResZ=new TH1F*[kNHisto];
194     target.fHisResXZ=new TH2F*[kNHisto];
195     target.fHisClusterSize=new TH2I*[kNHisto];
196     target.fHisResXclu=new TH1F**[kNHisto];
197     target.fHisResZclu=new TH1F**[kNHisto];
198     target.fHisResXchip=new TH1F**[kNHisto];
199     target.fHisResZchip=new TH1F**[kNHisto];
200     target.fHisTrackErrX=new TH1F*[kNHisto];
201     target.fHisTrackErrZ=new TH1F*[kNHisto];
202     target.fHisClusErrX=new TH1F*[kNHisto];
203     target.fHisClusErrZ=new TH1F*[kNHisto];
204     for(Int_t i=0; i<kNHisto; i++) {
205       target.fHisResX[i] = new TH1F(*fHisResX[i]);
206       target.fHisResZ[i] = new TH1F(*fHisResZ[i]);
207       target.fHisResXZ[i] = new TH2F(*fHisResXZ[i]);
208       target.fHisClusterSize[i] = new TH2I(*fHisClusterSize[i]);
209       target.fHisResXclu[i]=new TH1F*[kNclu];
210       target.fHisResZclu[i]=new TH1F*[kNclu];
211       for(Int_t clu=0; clu<kNclu; clu++) {  // clu=0 --> cluster size 1
212         target.fHisResXclu[i][clu] = new TH1F(*fHisResXclu[i][clu]);
213         target.fHisResZclu[i][clu] = new TH1F(*fHisResZclu[i][clu]);
214       }
215       target.fHisResXchip[i]=new TH1F*[kNChip];
216       target.fHisResZchip[i]=new TH1F*[kNChip];
217       for(Int_t chip=0; chip<kNChip; chip++) {  
218         target.fHisResXchip[i][chip] = new TH1F(*fHisResXchip[i][chip]);
219         target.fHisResZchip[i][chip] = new TH1F(*fHisResZchip[i][chip]);
220       }
221       target.fHisTrackErrX[i] = new TH1F(*fHisTrackErrX[i]);
222       target.fHisTrackErrZ[i] = new TH1F(*fHisTrackErrZ[i]);
223       target.fHisClusErrX[i] = new TH1F(*fHisClusErrX[i]);
224       target.fHisClusErrZ[i] = new TH1F(*fHisClusErrZ[i]);
225     }
226   }
227 return;
228 }
229 /*    commented out by M.Masera 8/3/08
230 //______________________________________________________________________
231 AliITSPlaneEff&  AliITSPlaneEffSPD::operator=(const
232                                            AliITSPlaneEff &s){
233     //    Assignment operator
234     // Inputs:
235     //    AliITSPlaneEffSPD &s The original class for which
236     //                                this class is a copy of
237     // Outputs:
238     //    none.
239     // Return:
240
241     if(&s == this) return *this;
242     AliError("operator=: Not allowed to make a =, use default creater instead");
243     return *this;
244 }
245 */
246 //_______________________________________________________________________
247 Int_t AliITSPlaneEffSPD::GetMissingTracksForGivenEff(Double_t eff, Double_t RelErr,
248           UInt_t im, UInt_t ic) const {
249    
250   //   Estimate the number of tracks still to be collected to attain a 
251   //   given efficiency eff, with relative error RelErr
252   //   Inputs:
253   //         eff    -> Expected efficiency (e.g. those from actual estimate)
254   //         RelErr -> tollerance [0,1] 
255   //         im     -> module number [0,249]
256   //         ic     -> chip number [0,4]
257   //   Outputs: none
258   //   Return: the estimated n. of tracks 
259   //
260 if (im>=kNModule || ic>=kNChip) 
261  {AliError("GetMissingTracksForGivenEff: you asked for a non existing chip");
262  return -1;}
263 else return GetNTracksForGivenEff(eff,RelErr)-fTried[GetKey(im,ic)];
264 }
265 //_________________________________________________________________________
266 Double_t  AliITSPlaneEffSPD::PlaneEff(const UInt_t im,const UInt_t ic) const {
267 // Compute the efficiency for a basic block, 
268 // Inputs:
269 //        im     -> module number [0,249]
270 //        ic     -> chip number [0,4] 
271 if (im>=kNModule || ic>=kNChip) 
272  {AliError("PlaneEff(Uint_t,Uint_t): you asked for a non existing chip"); return -1.;}
273  Int_t nf=fFound[GetKey(im,ic)];
274  Int_t nt=fTried[GetKey(im,ic)];
275 return AliITSPlaneEff::PlaneEff(nf,nt);
276 }
277 //_________________________________________________________________________
278 Double_t  AliITSPlaneEffSPD::ErrPlaneEff(const UInt_t im,const UInt_t ic) const {
279     // Compute the statistical error on efficiency for a basic block,
280     // using binomial statistics 
281     // Inputs:
282     //        im     -> module number [0,249]
283     //        ic     -> chip number [0,4] 
284 if (im>=kNModule || ic>=kNChip) 
285  {AliError("ErrPlaneEff(Uint_t,Uint_t): you asked for a non existing chip"); return -1.;}
286 Int_t nf=fFound[GetKey(im,ic)];
287 Int_t nt=fTried[GetKey(im,ic)];
288 return AliITSPlaneEff::ErrPlaneEff(nf,nt);
289
290 //_________________________________________________________________________
291 Bool_t AliITSPlaneEffSPD::UpDatePlaneEff(const Bool_t Kfound,
292                                          const UInt_t im, const UInt_t ic) {
293   // Update efficiency for a basic block
294 if (im>=kNModule || ic>=kNChip) 
295  {AliError("UpDatePlaneEff: you asked for a non existing chip"); return kFALSE;}
296  fTried[GetKey(im,ic)]++;
297  if(Kfound) fFound[GetKey(im,ic)]++;
298  return kTRUE;
299 }
300 //_________________________________________________________________________
301 UInt_t AliITSPlaneEffSPD::GetChipFromCol(const UInt_t col) const {
302   // get chip given the column
303 if(col>=kNCol*kNChip) 
304  {AliDebug(1,Form("GetChipFromCol: you asked for a non existing column %d",col)); return 10;}
305 return col/kNCol;
306 }
307 //__________________________________________________________________________
308 UInt_t AliITSPlaneEffSPD::GetKey(const UInt_t mod, const UInt_t chip) const {
309   // get key given a basic block
310 if(mod>=kNModule || chip>=kNChip)
311   {AliWarning("GetKey: you asked for a non existing block"); return 99999;}
312 return mod*kNChip+chip;
313 }
314 //__________________________________________________________________________
315 UInt_t AliITSPlaneEffSPD::GetModFromKey(const UInt_t key) const {
316   // get mod. from key
317 if(key>=kNModule*kNChip)
318   {AliError("GetModFromKey: you asked for a non existing key"); return 9999;}
319 return key/kNChip;
320 }
321 //__________________________________________________________________________
322 UInt_t AliITSPlaneEffSPD::GetChipFromKey(const UInt_t key) const {
323   // retrieves chip from key
324 if(key>=kNModule*kNChip)
325   {AliError("GetChipFromKey: you asked for a non existing key"); return 999;}
326 return (key%(kNModule*kNChip))%kNChip;
327 }
328 //__________________________________________________________________________
329 void AliITSPlaneEffSPD::GetModAndChipFromKey(const UInt_t key,UInt_t& mod,UInt_t& chip) const {
330   // get module and chip from a key
331 if(key>=kNModule*kNChip)
332   {AliError("GetModAndChipFromKey: you asked for a non existing key"); 
333   mod=9999;
334   chip=999;
335   return;}
336 mod=key/kNChip;
337 chip=(key%(kNModule*kNChip))%kNChip;
338 return;
339 }
340 //____________________________________________________________________________
341 Double_t AliITSPlaneEffSPD::LivePlaneEff(UInt_t key) const {
342   // returns plane efficieny after adding the fraction of sensor which is bad
343 if(key>=kNModule*kNChip)
344   {AliError("LivePlaneEff: you asked for a non existing key");
345    return -1.;}
346 Double_t leff=AliITSPlaneEff::LivePlaneEff(0); // this just for the Warning
347 leff=PlaneEff(key)+GetFracBad(key);
348 return leff>1?1:leff;
349 }
350 //____________________________________________________________________________
351 Double_t AliITSPlaneEffSPD::ErrLivePlaneEff(UInt_t key) const {
352   // returns error on live plane efficiency
353 if(key>=kNModule*kNChip)
354   {AliError("ErrLivePlaneEff: you asked for a non existing key");
355    return -1.;}
356 Int_t nf=fFound[key];
357 Double_t triedInLive=GetFracLive(key)*fTried[key];
358 Int_t nt=TMath::Max(nf,TMath::Nint(triedInLive));
359 return AliITSPlaneEff::ErrPlaneEff(nf,nt); // for the time being: to be checked
360 }
361 //_____________________________________________________________________________
362 Double_t AliITSPlaneEffSPD::GetFracLive(const UInt_t key) const {
363   // returns the fraction of the sensor which is OK
364 if(key>=kNModule*kNChip)
365   {AliError("GetFracLive: you asked for a non existing key");
366    return -1.;}
367     // Compute the fraction of bad (dead+noisy) detector 
368 UInt_t dead=0,noisy=0;
369 GetDeadAndNoisyInChip(key,dead,noisy);
370 Double_t live=dead+noisy;
371 live/=(kNRow*kNCol);
372 return 1.-live;
373 }
374 //_____________________________________________________________________________
375 void AliITSPlaneEffSPD::GetDeadAndNoisyInChip(const UInt_t key,
376       UInt_t& nrDeadInChip, UInt_t& nrNoisyInChip) const {
377   // returns the number of dead and noisy pixels
378 nrDeadInChip=0;
379 nrNoisyInChip=0;
380 if(key>=kNModule*kNChip)
381   {AliError("GetDeadAndNoisyInChip: you asked for a non existing key");
382    return;}
383     // Compute the number of bad (dead+noisy) pixel in a chip
384 //
385 if(!fInitCDBCalled) 
386   {AliError("GetDeadAndNoisyInChip: CDB not inizialized: call InitCDB first");
387    return;};
388 AliCDBManager* man = AliCDBManager::Instance();
389 // retrieve map of dead Pixel 
390 AliCDBEntry *cdbSPDDead = man->Get("ITS/Calib/SPDDead", fRunNumber);
391 TObjArray* spdDead;
392 if(cdbSPDDead) {
393   spdDead = (TObjArray*)cdbSPDDead->GetObject();
394   if(!spdDead) 
395   {AliError("GetDeadAndNoisyInChip: SPDDead not found in CDB");
396    return;}
397 } else {
398   AliError("GetDeadAndNoisyInChip: did not find Calib/SPDDead.");
399   return;
400 }
401 // retrieve map of noisy Pixel 
402 AliCDBEntry *cdbSPDNoisy = man->Get("ITS/Calib/SPDNoisy", fRunNumber);
403 TObjArray* spdNoisy;
404 if(cdbSPDNoisy) {
405   spdNoisy = (TObjArray*)cdbSPDNoisy->GetObject();
406   if(!spdNoisy) 
407   {AliError("GetDeadAndNoisyInChip: SPDNoisy not found in CDB");
408    return;}
409 } else {
410   AliError("GetDeadAndNoisyInChip: did not find Calib/SPDNoisy.");
411   return;
412 }
413 //
414 UInt_t mod=GetModFromKey(key);
415 UInt_t chip=GetChipFromKey(key);
416 // count number of dead
417 AliITSCalibrationSPD* calibSPD=(AliITSCalibrationSPD*) spdDead->At(mod);
418 UInt_t nrDead = calibSPD->GetNrBad();
419 for (UInt_t index=0; index<nrDead; index++) {
420   if(GetChipFromCol(calibSPD->GetBadColAt(index))==chip) nrDeadInChip++;
421 }
422 calibSPD=(AliITSCalibrationSPD*) spdNoisy->At(mod);
423 UInt_t nrNoisy = calibSPD->GetNrBad();
424 for (UInt_t index=0; index<nrNoisy; index++) {
425   if(GetChipFromCol(calibSPD->GetBadColAt(index))==chip) nrNoisyInChip++;
426 }
427 return;
428 }
429 //_____________________________________________________________________________
430 Double_t AliITSPlaneEffSPD::GetFracBad(const UInt_t key) const {
431   // returns 1-fractional live
432 if(key>=kNModule*kNChip)
433   {AliError("GetFracBad: you asked for a non existing key");
434    return -1.;}
435 return 1.-GetFracLive(key);
436 }
437 //_____________________________________________________________________________
438 Bool_t AliITSPlaneEffSPD::WriteIntoCDB() const {
439 // write onto CDB
440 if(!fInitCDBCalled)
441   {AliError("WriteIntoCDB: CDB not inizialized. Call InitCDB first");
442    return kFALSE;}
443 // to be written properly: now only for debugging 
444   AliCDBMetaData *md= new AliCDBMetaData(); // metaData describing the object
445   //md->SetObjectClassName("AliITSPlaneEff");
446   md->SetResponsible("Giuseppe Eugenio Bruno");
447   md->SetBeamPeriod(0);
448   md->SetAliRootVersion("head 19/11/07"); //root version
449   AliCDBId id("ITS/PlaneEff/PlaneEffSPD",0,AliCDBRunRange::Infinity()); 
450   AliITSPlaneEffSPD eff; 
451   eff=*this;
452   Bool_t r=AliCDBManager::Instance()->GetDefaultStorage()->Put(&eff,id,md);
453   delete md;
454   return r;
455 }
456 //_____________________________________________________________________________
457 Bool_t AliITSPlaneEffSPD::ReadFromCDB() {
458 // read from CDB
459 if(!fInitCDBCalled)
460   {AliError("ReadFromCDB: CDB not inizialized. Call InitCDB first");
461    return kFALSE;}
462 AliCDBEntry *cdbEntry = AliCDBManager::Instance()->Get("ITS/PlaneEff/PlaneEffSPD",fRunNumber);
463 if(!cdbEntry) return kFALSE;
464 AliITSPlaneEffSPD* eff= (AliITSPlaneEffSPD*)cdbEntry->GetObject();
465 if(this==eff) return kFALSE;
466 if(fHis) CopyHistos(*eff); // If histos already exist then copy them to eff
467 eff->Copy(*this);          // copy everything (statistics and histos) from eff to this
468 return kTRUE;
469 }
470 //_____________________________________________________________________________
471 Bool_t AliITSPlaneEffSPD::AddFromCDB(AliCDBId *cdbId) {
472 AliCDBEntry *cdbEntry=0;
473 if (!cdbId) {
474   if(!fInitCDBCalled)  
475     {AliError("ReadFromCDB: CDB not inizialized. Call InitCDB first"); return kFALSE;}
476   cdbEntry = AliCDBManager::Instance()->Get("ITS/PlaneEff/PlaneEffSPD",fRunNumber);
477 } else {
478   cdbEntry = AliCDBManager::Instance()->Get(*cdbId);
479 }
480 if(!cdbEntry) return kFALSE;
481 AliITSPlaneEffSPD* eff= (AliITSPlaneEffSPD*)cdbEntry->GetObject();
482 *this+=*eff;
483 return kTRUE;
484 }
485 //_____________________________________________________________________________
486 UInt_t AliITSPlaneEffSPD::GetKeyFromDetLocCoord(Int_t ilay, Int_t idet, 
487                                                 Float_t, Float_t locz) const {
488 // method to locate a basic block from Detector Local coordinate (to be used in tracking)
489 UInt_t key=999999;
490 if(ilay<0 || ilay>1) 
491   {AliError("GetKeyFromDetLocCoord: you asked for a non existing layer");
492    return key;}
493 if(ilay==0 && (idet<0 || idet>79))
494  {AliError("GetKeyFromDetLocCoord: you asked for a non existing detector");
495    return key;}
496 if(ilay==1 && (idet<0 || idet>159))
497  {AliError("GetKeyFromDetLocCoord: you asked for a non existing detector");
498    return key;}
499
500 UInt_t mod=idet;
501 if(ilay==1) mod+=80;
502 key=GetKey(mod,GetChipFromCol(GetColFromLocZ(locz)));
503 return key;
504 }
505 //_____________________________________________________________________________
506 UInt_t AliITSPlaneEffSPD::GetColFromLocZ(Float_t zloc) const {
507 // method to retrieve column number from the local z coordinate
508   UInt_t col=0;
509   AliITSsegmentationSPD spd;
510   Int_t ix,iz;
511   if(spd.LocalToDet(0,zloc,ix,iz)) col+=iz;
512   else {
513     AliDebug(1,Form("cannot compute column number from local z=%f",zloc));
514     col=99999;}
515   return col;
516 /*
517 const Float_t kconv = 1.0E-04; // converts microns to cm.
518 Float_t bz[160];
519 for(Int_t i=000;i<160;i++) bz[i] = 425.0; // most are 425 microns except below
520 bz[ 31] = bz[ 32] = 625.0; // first chip boundry
521 bz[ 63] = bz[ 64] = 625.0; // first chip boundry
522 bz[ 95] = bz[ 96] = 625.0; // first chip boundry
523 bz[127] = bz[128] = 625.0; // first chip boundry
524 //
525 Int_t j=-1;
526 Float_t dz=0;
527 for(Int_t i=000;i<160;i++) dz+=bz[i];
528 dz = -0.5*kconv*dz;
529 if(zloc<dz || zloc>-1*dz) { // outside z range
530   AliDebug(1,Form("GetColFromLocZ: cannot compute column number from local z=%f",zloc));
531   return 99999;}
532 for(j=0;j<160;j++){
533   dz += kconv*bz[j];
534   if(zloc<dz) break;
535 } // end for j
536 col+=j;
537 //
538 return col;
539 */
540 }
541 //________________________________________________________
542 Bool_t AliITSPlaneEffSPD::GetBlockBoundaries(const UInt_t key, Float_t& xmn,Float_t& xmx,
543                                              Float_t& zmn,Float_t& zmx) const {
544 //
545 //  This method return the geometrical boundaries of the active volume of a given 
546 //  basic block, in the detector reference system.
547 //  Input: unique key to locate a basic block.
548 //  
549 //  Output: Ymin, Ymax, Zmin, Zmax of a basic block (chip for SPD)
550 //  Return: kTRUE if computation was succesfully, kFALSE otherwise
551 //
552 if(key>=kNModule*kNChip)
553   {AliWarning("GetBlockBoundaries: you asked for a non existing key"); return kFALSE;}
554 UInt_t chip=GetChipFromKey(key);
555 zmn=GetLocZFromCol(chip*kNCol);
556 zmx=GetLocZFromCol((chip+1)*kNCol);
557 xmn=GetLocXFromRow(0);
558 xmx=GetLocXFromRow(kNRow);
559 //
560 Float_t tmp=zmn;
561 if(zmx<zmn) {zmn=zmx; zmx=tmp;}
562 tmp=xmn;
563 if(xmx<xmn) {xmn=xmx; xmx=tmp;}
564 return kTRUE;
565 }
566 //________________________________________________________
567 Float_t AliITSPlaneEffSPD::GetLocXFromRow(const UInt_t row) const {
568 // 
569 //  This method return the local (i.e. detector reference system) lower x coordinate 
570 //  of the row. To get the central value of a given row, you can do 
571 //  1/2*[LocXFromRow(row)+LocXFromRow(row+1)].
572 //
573 //  Input: row number in the range [0,kNRow] 
574 //  Output: lower local X coordinate of this row.
575 //
576 if(row>kNRow)  // not >= ! allow also computation of upper limit of the last row. 
577   {AliError("LocYFromRow: you asked for a non existing row"); return 9999999.;}
578 // Use only AliITSsegmentationSPD
579 AliITSsegmentationSPD spd;
580 Double_t dummy,x;
581 if(row==kNRow) spd.CellBoundries((Int_t)row-1,0,dummy,x,dummy,dummy);
582 else spd.CellBoundries((Int_t)row,0,x,dummy,dummy,dummy);
583 return (Float_t)x;
584
585 }
586 //________________________________________________________
587 Float_t AliITSPlaneEffSPD::GetLocZFromCol(const UInt_t col) const {
588 //
589 //  This method return the local (i.e. detector reference system) lower Z coordinate
590 //  of the column. To get the central value of a given column, you can do
591 //  1/2*[LocZFromCol(col)+LocZFromCol(col+1)].
592 //
593 //  Input: col number in the range [0,kNChip*kNCol]
594 //  Output: lower local Y coordinate of this row.
595 //
596 if(col>kNChip*kNCol) // not >= ! allow also computation of upper limit of the last column
597   {AliError("LocZFromCol: you asked for a non existing column"); return 9999999.;}
598 // Use only AliITSsegmentationSPD
599 AliITSsegmentationSPD spd;
600 Double_t dummy,y;
601 if(col==kNChip*kNCol) spd.CellBoundries(0,(Int_t)col-1,dummy,dummy,dummy,y);
602 else spd.CellBoundries(0,(Int_t)col,dummy,dummy,y,dummy);
603 return (Float_t)y;
604
605 }
606 //__________________________________________________________
607 void AliITSPlaneEffSPD::InitHistos() {
608   // for the moment let's create the histograms 
609   // module by  module
610   TString histnameResX="HistResX_mod_",aux;
611   TString histnameResZ="HistResZ_mod_";
612   TString histnameResXZ="HistResXZ_mod_";
613   TString histnameClusterType="HistClusterType_mod_";
614   TString histnameResXclu="HistResX_mod_";
615   TString histnameResZclu="HistResZ_mod_";
616   TString histnameResXchip="HistResX_mod_";
617   TString histnameResZchip="HistResZ_mod_";
618   TString histnameTrackErrX="HistTrackErrX_mod_";
619   TString histnameTrackErrZ="HistTrackErrZ_mod_";
620   TString histnameClusErrX="HistClusErrX_mod_";
621   TString histnameClusErrZ="HistClusErrZ_mod_";
622 //
623
624   TH1::AddDirectory(kFALSE);
625
626   fHisResX=new TH1F*[kNHisto];
627   fHisResZ=new TH1F*[kNHisto];
628   fHisResXZ=new TH2F*[kNHisto];
629   fHisClusterSize=new TH2I*[kNHisto];
630   fHisResXclu=new TH1F**[kNHisto];
631   fHisResZclu=new TH1F**[kNHisto];
632   fHisResXchip=new TH1F**[kNHisto];
633   fHisResZchip=new TH1F**[kNHisto];
634   fHisTrackErrX=new TH1F*[kNHisto];
635   fHisTrackErrZ=new TH1F*[kNHisto];
636   fHisClusErrX=new TH1F*[kNHisto];
637   fHisClusErrZ=new TH1F*[kNHisto];
638
639   for (Int_t nhist=0;nhist<kNHisto;nhist++){
640     aux=histnameResX;
641     aux+=nhist;
642     fHisResX[nhist]=new TH1F("histname","histname",800,-0.16,0.16); // +- 1600 micron; 1 bin=4 micron
643     fHisResX[nhist]->SetName(aux.Data());
644     fHisResX[nhist]->SetTitle(aux.Data());
645
646     aux=histnameResZ;
647     aux+=nhist;
648     fHisResZ[nhist]=new TH1F("histname","histname",800,-0.32,0.32); // +-3200 micron; 1 bin=8 micron
649     fHisResZ[nhist]->SetName(aux.Data());
650     fHisResZ[nhist]->SetTitle(aux.Data());
651
652     aux=histnameResXZ;
653     aux+=nhist;
654     fHisResXZ[nhist]=new TH2F("histname","histname",40,-0.08,0.08,40,-0.16,0.16); // binning:
655     fHisResXZ[nhist]->SetName(aux.Data());                                         // 40 micron in x;
656     fHisResXZ[nhist]->SetTitle(aux.Data());                                        // 80 micron in z;
657
658     aux=histnameClusterType;
659     aux+=nhist;
660     fHisClusterSize[nhist]=new TH2I("histname","histname",10,0.5,10.5,10,0.5,10.5);
661     fHisClusterSize[nhist]->SetName(aux.Data());
662     fHisClusterSize[nhist]->SetTitle(aux.Data());
663
664     fHisResXclu[nhist]=new TH1F*[kNclu];
665     fHisResZclu[nhist]=new TH1F*[kNclu];
666     for(Int_t clu=0; clu<kNclu; clu++) {  // clu=0 --> cluster size 1
667       aux=histnameResXclu;
668       aux+=nhist;
669       aux+="_clu_";
670       aux+=clu+1; // clu=0 --> cluster size 1
671       fHisResXclu[nhist][clu]=new TH1F("histname","histname",800,-0.16,0.16); // +- 1600 micron; 1 bin=4 micron
672       fHisResXclu[nhist][clu]->SetName(aux.Data());
673       fHisResXclu[nhist][clu]->SetTitle(aux.Data());
674
675       aux=histnameResZclu;
676       aux+=nhist;
677       aux+="_clu_";
678       aux+=clu+1; // clu=0 --> cluster size 1
679       fHisResZclu[nhist][clu]=new TH1F("histname","histname",800,-0.32,0.32); // +-3200 micron; 1 bin=8 micron
680       fHisResZclu[nhist][clu]->SetName(aux.Data());
681       fHisResZclu[nhist][clu]->SetTitle(aux.Data());
682     }
683
684     fHisResXchip[nhist]=new TH1F*[kNChip];
685     fHisResZchip[nhist]=new TH1F*[kNChip];
686     for(Int_t chip=0; chip<kNChip; chip++) { 
687       aux=histnameResXchip;
688       aux+=nhist;
689       aux+="_chip_";
690       aux+=chip; 
691       fHisResXchip[nhist][chip]=new TH1F("histname","histname",200,-0.08,0.08); // +- 800 micron; 1 bin=8 micron
692       fHisResXchip[nhist][chip]->SetName(aux.Data());
693       fHisResXchip[nhist][chip]->SetTitle(aux.Data());
694
695       aux=histnameResZchip;
696       aux+=nhist;
697       aux+="_chip_";
698       aux+=chip;
699       fHisResZchip[nhist][chip]=new TH1F("histname","histname",200,-0.32,0.32); // +-3200 micron; 1 bin=32 micron
700       fHisResZchip[nhist][chip]->SetName(aux.Data());
701       fHisResZchip[nhist][chip]->SetTitle(aux.Data());
702     }
703
704     aux=histnameTrackErrX;
705     aux+=nhist;
706     fHisTrackErrX[nhist]=new TH1F("histname","histname",200,0.,0.16); // 0-1600 micron; 1 bin=8 micron
707     fHisTrackErrX[nhist]->SetName(aux.Data());
708     fHisTrackErrX[nhist]->SetTitle(aux.Data());
709
710     aux=histnameTrackErrZ;
711     aux+=nhist;
712     fHisTrackErrZ[nhist]=new TH1F("histname","histname",200,0.,0.32); // 0-3200 micron; 1 bin=16 micron
713     fHisTrackErrZ[nhist]->SetName(aux.Data());
714     fHisTrackErrZ[nhist]->SetTitle(aux.Data());
715
716     aux=histnameClusErrX;
717     aux+=nhist;
718     fHisClusErrX[nhist]=new TH1F("histname","histname",200,0.,0.04); //  0-400 micron; 1 bin=2 micron
719     fHisClusErrX[nhist]->SetName(aux.Data());
720     fHisClusErrX[nhist]->SetTitle(aux.Data());
721
722     aux=histnameClusErrZ;
723     aux+=nhist;
724     fHisClusErrZ[nhist]=new TH1F("histname","histname",200,0.,0.16); //  0-1600 micron; 1 bin=8 micron
725     fHisClusErrZ[nhist]->SetName(aux.Data());
726     fHisClusErrZ[nhist]->SetTitle(aux.Data());
727
728   }
729
730   TH1::AddDirectory(kTRUE);
731
732 return;
733 }
734 //__________________________________________________________
735 void AliITSPlaneEffSPD::DeleteHistos() {
736   if(fHisResX) {
737     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResX[i];
738     delete [] fHisResX; fHisResX=0;
739   }
740   if(fHisResZ) {
741     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResZ[i];
742     delete [] fHisResZ; fHisResZ=0;
743   }
744   if(fHisResXZ) {
745     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResXZ[i];
746     delete [] fHisResXZ; fHisResXZ=0;
747   }
748   if(fHisClusterSize) {
749     for (Int_t i=0; i<kNHisto; i++ ) delete fHisClusterSize[i];
750     delete [] fHisClusterSize; fHisClusterSize=0;
751   }
752   if(fHisResXclu) {
753     for (Int_t i=0; i<kNHisto; i++ ) {
754       for (Int_t clu=0; clu<kNclu; clu++) if (fHisResXclu[i][clu]) delete fHisResXclu[i][clu];
755       delete [] fHisResXclu[i];
756     }
757     delete [] fHisResXclu;
758     fHisResXclu = 0;
759   }
760   if(fHisResZclu) {
761     for (Int_t i=0; i<kNHisto; i++ ) {
762       for (Int_t clu=0; clu<kNclu; clu++) if (fHisResZclu[i][clu]) delete fHisResZclu[i][clu];
763       delete [] fHisResZclu[i];
764     }
765     delete [] fHisResZclu;
766     fHisResZclu = 0;
767   }
768   if(fHisResXchip) {
769     for (Int_t i=0; i<kNHisto; i++ ) {
770       for (Int_t chip=0; chip<kNChip; chip++) if (fHisResXchip[i][chip]) delete fHisResXchip[i][chip];
771       delete [] fHisResXchip[i];
772     }
773     delete [] fHisResXchip;
774     fHisResXchip = 0;
775   }
776   if(fHisResZchip) {
777     for (Int_t i=0; i<kNHisto; i++ ) {
778       for (Int_t chip=0; chip<kNChip; chip++) if (fHisResZchip[i][chip]) delete fHisResZchip[i][chip];
779       delete [] fHisResZchip[i];
780     }
781     delete [] fHisResZchip;
782     fHisResZchip = 0;
783   }
784   if(fHisTrackErrX) {
785     for (Int_t i=0; i<kNHisto; i++ ) delete fHisTrackErrX[i];
786     delete [] fHisTrackErrX; fHisTrackErrX=0;
787   }
788   if(fHisTrackErrZ) {
789     for (Int_t i=0; i<kNHisto; i++ ) delete fHisTrackErrZ[i];
790     delete [] fHisTrackErrZ; fHisTrackErrZ=0;
791   }
792   if(fHisClusErrX) {
793     for (Int_t i=0; i<kNHisto; i++ ) delete fHisClusErrX[i];
794     delete [] fHisClusErrX; fHisClusErrX=0;
795   }
796   if(fHisClusErrZ) {
797     for (Int_t i=0; i<kNHisto; i++ ) delete fHisClusErrZ[i];
798     delete [] fHisClusErrZ; fHisClusErrZ=0;
799   }
800
801 return;
802 }
803 //__________________________________________________________
804 Bool_t AliITSPlaneEffSPD::FillHistos(UInt_t key, Bool_t found, 
805                                     // Float_t tXZ[2], Float_t cXZ[2], Int_t ctXZ[2]) {
806                                      Float_t *tr, Float_t *clu, Int_t *csize) {
807 // this method fill the histograms
808 // input: - key: unique key of the basic block 
809 //        - found: Boolean to asses whether a cluster has been associated to the track or not 
810 //        - tr[0],tr[1] local X and Z coordinates of the track prediction, respectively
811 //        - tr[2],tr[3] error on local X and Z coordinates of the track prediction, respectively
812 //        - clu[0],clu[1] local X and Z coordinates of the cluster associated to the track, respectively
813 //        - clu[2],clu[3] error on local X and Z coordinates of the cluster associated to the track, respectively
814 //        - csize[0][1] cluster size in X and Z, respectively
815 // output: kTRUE if filling was succesfull kFALSE otherwise
816 // side effects: updating of the histograms. 
817 //
818   if (!fHis) {
819     AliWarning("FillHistos: histograms do not exist! Call SetCreateHistos(kTRUE) first");
820     return kFALSE;
821   }
822   if(key>=kNModule*kNChip)
823     {AliWarning("FillHistos: you asked for a non existing key"); return kFALSE;}
824   Int_t id=GetModFromKey(key);
825   Int_t chip=GetChipFromKey(key);
826   if(id>=kNHisto) 
827     {AliWarning("FillHistos: you want to fill a non-existing histos"); return kFALSE;}
828   if(found) {
829     Float_t resx=tr[0]-clu[0];
830     Float_t resz=tr[1]-clu[1];
831     fHisResX[id]->Fill(resx);
832     fHisResZ[id]->Fill(resz);
833     fHisResXZ[id]->Fill(resx,resz);
834     fHisClusterSize[id]->Fill((Double_t)csize[0],(Double_t)csize[1]);
835     if(csize[0]>0 &&  csize[0]<=kNclu) fHisResXclu[id][csize[0]-1]->Fill(resx);
836     if(csize[1]>0 &&  csize[1]<=kNclu) fHisResZclu[id][csize[1]-1]->Fill(resz);
837     fHisResXchip[id][chip]->Fill(resx);
838     fHisResZchip[id][chip]->Fill(resz);
839   }
840   fHisTrackErrX[id]->Fill(tr[2]);
841   fHisTrackErrZ[id]->Fill(tr[3]);
842   fHisClusErrX[id]->Fill(clu[2]);
843   fHisClusErrZ[id]->Fill(clu[3]);
844   return kTRUE;
845 }
846 //__________________________________________________________
847 Bool_t AliITSPlaneEffSPD::WriteHistosToFile(TString filename, Option_t* option) {
848   //
849   // Saves the histograms into a tree and saves the trees into a file
850   //
851   if (!fHis) return kFALSE;
852   if (filename.IsNull() || filename.IsWhitespace()) {
853      AliWarning("WriteHistosToFile: null output filename!");
854      return kFALSE;
855   }
856   char branchname[30];
857   TFile *hFile=new TFile(filename.Data(),option,
858                          "The File containing the TREEs with ITS PlaneEff Histos");
859   TTree *SPDTree=new TTree("SPDTree","Tree whith Residuals and Cluster Type distributions for SPD");
860   TH1F *histZ,*histX;
861   TH2F *histXZ;
862   TH2I *histClusterType;
863   TH1F *histXclu[kNclu];
864   TH1F *histZclu[kNclu];
865   TH1F *histXchip[kNChip];
866   TH1F *histZchip[kNChip];
867   TH1F *histTrErrZ,*histTrErrX;
868   TH1F *histClErrZ,*histClErrX;
869
870   histZ=new TH1F();
871   histX=new TH1F();
872   histXZ=new TH2F();
873   histClusterType=new TH2I();
874   for(Int_t clu=0;clu<kNclu;clu++) {
875     histXclu[clu]=new TH1F();
876     histZclu[clu]=new TH1F();
877   }
878   for(Int_t chip=0;chip<kNChip;chip++) {
879     histXchip[chip]=new TH1F();
880     histZchip[chip]=new TH1F();
881   }
882   histTrErrX=new TH1F();
883   histTrErrZ=new TH1F();
884   histClErrX=new TH1F();
885   histClErrZ=new TH1F();
886
887   SPDTree->Branch("histX","TH1F",&histX,128000,0);
888   SPDTree->Branch("histZ","TH1F",&histZ,128000,0);
889   SPDTree->Branch("histXZ","TH2F",&histXZ,128000,0);
890   SPDTree->Branch("histClusterType","TH2I",&histClusterType,128000,0);
891   for(Int_t clu=0;clu<kNclu;clu++) {
892     sprintf(branchname,"histXclu_%d",clu+1);
893     SPDTree->Branch(branchname,"TH1F",&histXclu[clu],128000,0);
894     sprintf(branchname,"histZclu_%d",clu+1);
895     SPDTree->Branch(branchname,"TH1F",&histZclu[clu],128000,0);
896   }
897   for(Int_t chip=0;chip<kNChip;chip++) {
898     sprintf(branchname,"histXchip_%d",chip);
899     SPDTree->Branch(branchname,"TH1F",&histXchip[chip],128000,0);
900     sprintf(branchname,"histZchip_%d",chip);
901     SPDTree->Branch(branchname,"TH1F",&histZchip[chip],128000,0);
902   }
903   SPDTree->Branch("histTrErrX","TH1F",&histTrErrX,128000,0);
904   SPDTree->Branch("histTrErrZ","TH1F",&histTrErrZ,128000,0);
905   SPDTree->Branch("histClErrX","TH1F",&histClErrX,128000,0);
906   SPDTree->Branch("histClErrZ","TH1F",&histClErrZ,128000,0);
907
908   for(Int_t j=0;j<kNHisto;j++){
909     histX=fHisResX[j];
910     histZ=fHisResZ[j];
911     histXZ=fHisResXZ[j];
912     histClusterType=fHisClusterSize[j];
913     for(Int_t clu=0;clu<kNclu;clu++) {
914       histXclu[clu]=fHisResXclu[j][clu];
915       histZclu[clu]=fHisResZclu[j][clu];
916     }
917     for(Int_t chip=0;chip<kNChip;chip++) {
918       histXchip[chip]=fHisResXchip[j][chip];
919       histZchip[chip]=fHisResZchip[j][chip];
920     }
921     histTrErrX=fHisTrackErrX[j];
922     histTrErrZ=fHisTrackErrZ[j];
923     histClErrX=fHisClusErrX[j];
924     histClErrZ=fHisClusErrZ[j];
925     SPDTree->Fill();
926   }
927   hFile->Write();
928   hFile->Close();
929 return kTRUE;
930 }
931 //__________________________________________________________
932 Bool_t AliITSPlaneEffSPD::ReadHistosFromFile(TString filename) {
933   //
934   // Read histograms from an already existing file 
935   //
936   if (!fHis) return kFALSE;
937   if (filename.IsNull() || filename.IsWhitespace()) {
938      AliWarning("ReadHistosFromFile: incorrect output filename!");
939      return kFALSE;
940   }
941   char branchname[30];
942
943   TH1F *h  = 0;
944   TH2F *h2 = 0;
945   TH2I *h2i= 0;
946
947   TFile *file=TFile::Open(filename.Data(),"READONLY");
948
949   if (!file || file->IsZombie()) {
950     AliWarning(Form("Can't open %s !",filename.Data()));
951     delete file;
952     return kFALSE;
953   }
954   TTree *tree = (TTree*) file->Get("SPDTree");
955
956   TBranch *histX = (TBranch*) tree->GetBranch("histX");
957   TBranch *histZ = (TBranch*) tree->GetBranch("histZ");
958   TBranch *histXZ = (TBranch*) tree->GetBranch("histXZ");
959   TBranch *histClusterType = (TBranch*) tree->GetBranch("histClusterType");
960    
961   TBranch *histXclu[kNclu], *histZclu[kNclu];
962   for(Int_t clu=0; clu<kNclu; clu++) {
963     sprintf(branchname,"histXclu_%d",clu+1);
964     histXclu[clu]= (TBranch*) tree->GetBranch(branchname);
965     sprintf(branchname,"histZclu_%d",clu+1);
966     histZclu[clu]= (TBranch*) tree->GetBranch(branchname);
967   }
968
969   TBranch *histXchip[kNChip], *histZchip[kNChip];
970   for(Int_t chip=0; chip<kNChip; chip++) {
971     sprintf(branchname,"histXchip_%d",chip);
972     histXchip[chip]= (TBranch*) tree->GetBranch(branchname);
973     sprintf(branchname,"histZchip_%d",chip);
974     histZchip[chip]= (TBranch*) tree->GetBranch(branchname);
975   }
976
977   TBranch *histTrErrX = (TBranch*) tree->GetBranch("histTrErrX");
978   TBranch *histTrErrZ = (TBranch*) tree->GetBranch("histTrErrZ");
979   TBranch *histClErrX = (TBranch*) tree->GetBranch("histClErrX");
980   TBranch *histClErrZ = (TBranch*) tree->GetBranch("histClErrZ");
981
982   gROOT->cd();
983
984   Int_t nevent = (Int_t)histX->GetEntries();
985   if(nevent!=kNHisto) 
986     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
987   histX->SetAddress(&h);
988   for(Int_t j=0;j<kNHisto;j++){
989     delete h; h=0;
990     histX->GetEntry(j);
991     fHisResX[j]->Add(h);
992   }
993
994   nevent = (Int_t)histZ->GetEntries();
995   if(nevent!=kNHisto) 
996     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
997   histZ->SetAddress(&h);
998   for(Int_t j=0;j<kNHisto;j++){
999     delete h; h=0;
1000     histZ->GetEntry(j);
1001     fHisResZ[j]->Add(h);
1002   }
1003
1004   nevent = (Int_t)histXZ->GetEntries();
1005   if(nevent!=kNHisto) 
1006     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1007   histXZ->SetAddress(&h2);
1008   for(Int_t j=0;j<kNHisto;j++){
1009     delete h2; h2=0;
1010     histXZ->GetEntry(j);
1011     fHisResXZ[j]->Add(h2);
1012   }
1013
1014   nevent = (Int_t)histClusterType->GetEntries();
1015   if(nevent!=kNHisto) 
1016     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1017   histClusterType->SetAddress(&h2i);
1018   for(Int_t j=0;j<kNHisto;j++){
1019     delete h2i; h2i=0;
1020     histClusterType->GetEntry(j);
1021     fHisClusterSize[j]->Add(h2i);
1022   }
1023
1024   for(Int_t clu=0; clu<kNclu; clu++) {
1025
1026     nevent = (Int_t)histXclu[clu]->GetEntries();
1027     if(nevent!=kNHisto)
1028       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1029     histXclu[clu]->SetAddress(&h);
1030     for(Int_t j=0;j<kNHisto;j++){
1031       delete h; h=0;
1032       histXclu[clu]->GetEntry(j);
1033       fHisResXclu[j][clu]->Add(h);
1034     }
1035
1036    nevent = (Int_t)histZclu[clu]->GetEntries();
1037     if(nevent!=kNHisto)
1038       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1039     histZclu[clu]->SetAddress(&h);
1040     for(Int_t j=0;j<kNHisto;j++){
1041       delete h; h=0;
1042       histZclu[clu]->GetEntry(j);
1043       fHisResZclu[j][clu]->Add(h);
1044     }
1045   }
1046
1047
1048     for(Int_t chip=0; chip<kNChip; chip++) {
1049
1050     nevent = (Int_t)histXchip[chip]->GetEntries();
1051     if(nevent!=kNHisto)
1052       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1053     histXchip[chip]->SetAddress(&h);
1054     for(Int_t j=0;j<kNHisto;j++){
1055       delete h; h=0;
1056       histXchip[chip]->GetEntry(j);
1057       fHisResXchip[j][chip]->Add(h);
1058     }
1059
1060     nevent = (Int_t)histZchip[chip]->GetEntries();
1061     if(nevent!=kNHisto)
1062       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1063     histZchip[chip]->SetAddress(&h);
1064     for(Int_t j=0;j<kNHisto;j++){
1065       delete h; h=0;
1066       histZchip[chip]->GetEntry(j);
1067       fHisResZchip[j][chip]->Add(h);
1068     }
1069   }
1070
1071   nevent = (Int_t)histTrErrX->GetEntries(); 
1072   if(nevent!=kNHisto)
1073     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1074   histTrErrX->SetAddress(&h);
1075   for(Int_t j=0;j<kNHisto;j++){
1076     delete h; h=0;
1077     histTrErrX->GetEntry(j);
1078     fHisTrackErrX[j]->Add(h);
1079   }
1080
1081   nevent = (Int_t)histTrErrZ->GetEntries();
1082   if(nevent!=kNHisto)
1083     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1084   histTrErrZ->SetAddress(&h);
1085   for(Int_t j=0;j<kNHisto;j++){
1086     delete h; h=0;
1087     histTrErrZ->GetEntry(j);
1088     fHisTrackErrZ[j]->Add(h);
1089   }
1090
1091   nevent = (Int_t)histClErrX->GetEntries();
1092   if(nevent!=kNHisto)
1093     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1094   histClErrX->SetAddress(&h);
1095   for(Int_t j=0;j<kNHisto;j++){
1096     delete h; h=0;
1097     histClErrX->GetEntry(j);
1098     fHisClusErrX[j]->Add(h);
1099   }
1100
1101   nevent = (Int_t)histClErrZ->GetEntries();
1102   if(nevent!=kNHisto)
1103     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1104   histClErrZ->SetAddress(&h);
1105   for(Int_t j=0;j<kNHisto;j++){
1106     delete h; h=0;
1107     histClErrZ->GetEntry(j);
1108     fHisClusErrZ[j]->Add(h);
1109   }
1110
1111   delete h;   h=0;
1112   delete h2;  h2=0;
1113   delete h2i; h2i=0;
1114
1115   if (file) {
1116     file->Close();
1117   }
1118 return kTRUE;
1119 }