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