]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSPlaneEffSSD.cxx
moving class
[u/mrichter/AliRoot.git] / ITS / AliITSPlaneEffSSD.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 module by module efficiency of the SSD,        
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 "AliITSPlaneEffSSD.h"
33 #include "AliLog.h"
34 #include "AliCDBStorage.h"
35 #include "AliCDBEntry.h"
36 #include "AliCDBManager.h"
37 //#include "AliCDBRunRange.h"
38 #include "AliITSCalibrationSSD.h"
39
40 ClassImp(AliITSPlaneEffSSD)     
41 //______________________________________________________________________
42 AliITSPlaneEffSSD::AliITSPlaneEffSSD():
43   AliITSPlaneEff(),
44   fHisResX(0),
45   fHisResZ(0),
46   fHisResXZ(0),
47   fHisClusterSize(0){
48   for (UInt_t i=0; i<kNModule; i++){
49     fFound[i]=0;
50     fTried[i]=0;
51   }
52   // default constructor
53   AliDebug(1,Form("Calling default constructor"));
54 }
55 //______________________________________________________________________
56 AliITSPlaneEffSSD::~AliITSPlaneEffSSD(){
57     // destructor
58     // Inputs:
59     //    none.
60     // Outputs:
61     //    none.
62     // Return:
63     //     none.
64     DeleteHistos();
65 }
66 //______________________________________________________________________
67 AliITSPlaneEffSSD::AliITSPlaneEffSSD(const AliITSPlaneEffSSD &s) : AliITSPlaneEff(s),
68 fHisResX(0),
69 fHisResZ(0),
70 fHisResXZ(0),
71 fHisClusterSize(0)
72 {
73     //     Copy Constructor
74     // Inputs:
75     //    AliITSPlaneEffSSD &s The original class for which
76     //                                this class is a copy of
77     // Outputs:
78     //    none.
79     // Return:
80
81  for (UInt_t i=0; i<kNModule; i++){
82     fFound[i]=s.fFound[i];
83     fTried[i]=s.fTried[i];
84  }
85  if(fHis) {
86    InitHistos();
87    for(Int_t i=0; i<kNHisto; i++) {
88       s.fHisResX[i]->Copy(*fHisResX[i]);
89       s.fHisResZ[i]->Copy(*fHisResZ[i]);
90       s.fHisResXZ[i]->Copy(*fHisResXZ[i]);
91       s.fHisClusterSize[i]->Copy(*fHisClusterSize[i]);
92    }
93  }
94 }
95 //_________________________________________________________________________
96 AliITSPlaneEffSSD& AliITSPlaneEffSSD::operator+=(const AliITSPlaneEffSSD &add){
97     //    Add-to-me operator
98     // Inputs:
99     //    const AliITSPlaneEffSSD &add  simulation class to be added
100     // Outputs:
101     //    none.
102     // Return:
103     //    none
104     for (UInt_t i=0; i<kNModule; i++){
105       fFound[i] += add.fFound[i];
106       fTried[i] += add.fTried[i];
107     }
108     if(fHis && add.fHis) {
109       for(Int_t i=0; i<kNHisto; i++) {
110         fHisResX[i]->Add(add.fHisResX[i]);
111         fHisResZ[i]->Add(add.fHisResZ[i]);
112         fHisResXZ[i]->Add(add.fHisResXZ[i]);
113         fHisClusterSize[i]->Add(add.fHisClusterSize[i]);
114       }
115     }
116     return *this;
117 }
118 //______________________________________________________________________
119 AliITSPlaneEffSSD&  AliITSPlaneEffSSD::operator=(const
120                                            AliITSPlaneEffSSD &s){
121     //    Assignment operator
122     // Inputs:
123     //    AliITSPlaneEffSSD &s The original class for which
124     //                                this class is a copy of
125     // Outputs:
126     //    none.
127     // Return:
128  
129     if(this==&s) return *this;
130     s.Copy(*this);
131     return *this;
132 }
133 //______________________________________________________________________
134 void AliITSPlaneEffSSD::Copy(TObject &obj) const {
135   // protected method. copy this to obj
136   AliITSPlaneEff::Copy(obj);
137   AliITSPlaneEffSSD& target = (AliITSPlaneEffSSD &) obj;
138   for(Int_t i=0;i<kNModule;i++) {
139       target.fFound[i] = fFound[i];
140       target.fTried[i] = fTried[i];
141   }
142   CopyHistos(target);
143   return;
144 }
145 //_______________________________________________________________________
146 void AliITSPlaneEffSSD::CopyHistos(AliITSPlaneEffSSD &target) const {
147   // protected method: copy histos from this to target
148   target.fHis  = fHis; // this is redundant only in some cases. Leave as it is.
149   if(fHis) {
150     target.fHisResX=new TH1F*[kNHisto];
151     target.fHisResZ=new TH1F*[kNHisto];
152     target.fHisResXZ=new TH2F*[kNHisto];
153     target.fHisClusterSize=new TH2I*[kNHisto];
154     for(Int_t i=0; i<kNHisto; i++) {
155       target.fHisResX[i] = new TH1F(*fHisResX[i]);
156       target.fHisResZ[i] = new TH1F(*fHisResZ[i]);
157       target.fHisResXZ[i] = new TH2F(*fHisResXZ[i]);
158       target.fHisClusterSize[i] = new TH2I(*fHisClusterSize[i]);
159     }
160   }
161 return;
162 }
163 /* Commented out by M.Masera 8/3/08
164 //______________________________________________________________________
165 AliITSPlaneEff&  AliITSPlaneEffSSD::operator=(const
166                                            AliITSPlaneEff &s){
167     //    Assignment operator
168     // Inputs:
169     //    AliITSPlaneEffSSD &s The original class for which
170     //                                this class is a copy of
171     // Outputs:
172     //    none.
173     // Return:
174
175     if(&s == this) return *this;
176     AliError("operator=: Not allowed to make a =, use default creater instead");
177     return *this;
178 }
179 */
180 //_______________________________________________________________________
181 Int_t AliITSPlaneEffSSD::GetMissingTracksForGivenEff(Double_t eff, Double_t RelErr,
182           UInt_t im) const {
183    
184   //   Estimate the number of tracks still to be collected to attain a 
185   //   given efficiency eff, with relative error RelErr
186   //   Inputs:
187   //         eff    -> Expected efficiency (e.g. those from actual estimate)
188   //         RelErr -> tollerance [0,1] 
189   //         im     -> module number [0,1697]
190   //   Outputs: none
191   //   Return: the estimated n. of tracks 
192   //
193 if (im>=kNModule) 
194  {AliError("GetMissingTracksForGivenEff: you asked for a non existing module");
195  return -1;}
196 else return GetNTracksForGivenEff(eff,RelErr)-fTried[GetKey(im)];
197 }
198 //_________________________________________________________________________
199 Double_t  AliITSPlaneEffSSD::PlaneEff(const UInt_t im) const {
200 // Compute the efficiency for a basic block, 
201 // Inputs:
202 //        im     -> module number [0,1697]
203 if (im>=kNModule) 
204  {AliError("PlaneEff(UInt_t): you asked for a non existing module"); return -1.;}
205  Int_t nf=fFound[GetKey(im)];
206  Int_t nt=fTried[GetKey(im)];
207 return AliITSPlaneEff::PlaneEff(nf,nt);
208 }
209 //_________________________________________________________________________
210 Double_t  AliITSPlaneEffSSD::ErrPlaneEff(const UInt_t im) const {
211     // Compute the statistical error on efficiency for a basic block,
212     // using binomial statistics 
213     // Inputs:
214     //        im     -> module number [0,1697]
215 if (im>=kNModule) 
216  {AliError("ErrPlaneEff(UInt_t): you asked for a non existing module"); return -1.;}
217 Int_t nf=fFound[GetKey(im)];
218 Int_t nt=fTried[GetKey(im)];
219 return AliITSPlaneEff::ErrPlaneEff(nf,nt);
220
221 //_________________________________________________________________________
222 Bool_t AliITSPlaneEffSSD::UpDatePlaneEff(const Bool_t Kfound, const UInt_t im) {
223   // Update efficiency for a basic block
224 if (im>=kNModule) 
225  {AliError("UpDatePlaneEff: you asked for a non existing module"); return kFALSE;}
226  fTried[GetKey(im)]++;
227  if(Kfound) fFound[GetKey(im)]++;
228  return kTRUE;
229 }
230 //_________________________________________________________________________
231 UInt_t AliITSPlaneEffSSD::GetKey(const UInt_t mod) const {
232   // get key given a basic block
233 if(mod>=kNModule)
234   {AliError("GetKey: you asked for a non existing block"); return 99999;}
235 return mod;
236 }
237 //__________________________________________________________________________
238 UInt_t AliITSPlaneEffSSD::GetModFromKey(const UInt_t key) const {
239   // get mod. from key
240 if(key>=kNModule)
241   {AliError("GetModFromKey: you asked for a non existing key"); return 9999;}
242 return key;
243 }
244 //__________________________________________________________________________
245 Double_t AliITSPlaneEffSSD::LivePlaneEff(UInt_t key) const {
246   // returns plane efficieny after adding the fraction of sensor which is bad
247 if(key>=kNModule)
248   {AliError("LivePlaneEff: you asked for a non existing key");
249    return -1.;}
250 Double_t leff=AliITSPlaneEff::LivePlaneEff(0); // this just for the Warning
251 leff=PlaneEff(key)+GetFracBad(key);
252 return leff>1?1:leff;
253 }
254 //____________________________________________________________________________
255 Double_t AliITSPlaneEffSSD::ErrLivePlaneEff(UInt_t key) const {
256   // returns error on live plane efficiency
257 if(key>=kNModule)
258   {AliError("ErrLivePlaneEff: you asked for a non existing key");
259    return -1.;}
260 Int_t nf=fFound[key];
261 Double_t triedInLive=GetFracLive(key)*fTried[key];
262 Int_t nt=TMath::Max(nf,TMath::Nint(triedInLive));
263 return AliITSPlaneEff::ErrPlaneEff(nf,nt); // for the time being: to be checked
264 }
265 //_____________________________________________________________________________
266 Double_t AliITSPlaneEffSSD::GetFracLive(const UInt_t key) const {
267   // returns the fraction of the sensor area which is OK (neither noisy nor dead)
268   // As for now, it computes only the fraction of good strips / total strips. 
269   // If this fraction is large, then the computation is a good approximation. 
270   // In any case, it is a lower limit of the fraction of the live area. 
271   // The next upgrades would be to add the fraction of area of superoposition 
272   // between bad N-side strips and bad P-side strips.  
273 if(key>=kNModule)
274   {AliError("GetFracLive: you asked for a non existing key");
275    return -1.;}
276 AliInfo("GetFracLive: it computes only the fraction of working strips (N+P side) / total strips");
277 UInt_t bad=0;
278 GetBadInModule(key,bad);
279 Double_t live=bad;
280 live/=(kNChip*kNSide*kNStrip);
281 return 1.-live;
282 }
283 //_____________________________________________________________________________
284 void AliITSPlaneEffSSD::GetBadInModule(const UInt_t key, UInt_t& nrBadInMod) const {
285   // returns the number of dead and noisy strips (sum of P and N sides).
286 nrBadInMod=0;
287 if(key>=kNModule)
288   {AliError("GetBadInModule: you asked for a non existing key");
289    return;}
290     // Compute the number of bad (dead+noisy) pixel in a module
291 //
292 if(!fInitCDBCalled) 
293   {AliError("GetBadInModule: CDB not inizialized: call InitCDB first");
294    return;};
295 AliCDBManager* man = AliCDBManager::Instance();
296 // retrieve map of dead Pixel 
297 AliCDBEntry *cdbSSD = man->Get("ITS/Calib/BadChannelsSSD", fRunNumber);
298 TObjArray* ssdEntry;
299 if(cdbSSD) {
300   ssdEntry = (TObjArray*)cdbSSD->GetObject();
301   if(!ssdEntry) 
302   {AliError("GetBadInChip: SSDEntry not found in CDB");
303    return;}
304 } else {
305   AliError("GetBadInChip: did not find Calib/BadChannelsSSD");
306   return;
307 }
308 //
309 UInt_t mod=GetModFromKey(key);
310 //
311 AliITSBadChannelsSSD* badchannels=(AliITSBadChannelsSSD*) ssdEntry->At(mod);
312 // count the  number of bad channels on the p side
313 nrBadInMod += (badchannels->GetBadPChannelsList()).GetSize();
314 // add the  number of bad channels on the s side
315 nrBadInMod += (badchannels->GetBadNChannelsList()).GetSize();
316 return;
317 }
318 //_____________________________________________________________________________
319 Double_t AliITSPlaneEffSSD::GetFracBad(const UInt_t key) const {
320   // returns 1-fractional live
321 if(key>=kNModule)
322   {AliError("GetFracBad: you asked for a non existing key");
323    return -1.;}
324 return 1.-GetFracLive(key);
325 }
326 //_____________________________________________________________________________
327 Bool_t AliITSPlaneEffSSD::WriteIntoCDB() const {
328 // write onto CDB
329 if(!fInitCDBCalled)
330   {AliError("WriteIntoCDB: CDB not inizialized. Call InitCDB first");
331    return kFALSE;}
332 // to be written properly: now only for debugging 
333   AliCDBMetaData *md= new AliCDBMetaData(); // metaData describing the object
334   md->SetObjectClassName("AliITSPlaneEff");
335   md->SetResponsible("Giuseppe Eugenio Bruno");
336   md->SetBeamPeriod(0);
337   md->SetAliRootVersion("head 02/01/08"); //root version
338   AliCDBId id("ITS/PlaneEff/PlaneEffSSD",0,AliCDBRunRange::Infinity()); 
339   AliITSPlaneEffSSD eff; 
340   eff=*this;
341   Bool_t r=AliCDBManager::Instance()->GetDefaultStorage()->Put(&eff,id,md);
342   delete md;
343   return r;
344 }
345 //_____________________________________________________________________________
346 Bool_t AliITSPlaneEffSSD::ReadFromCDB() {
347 // read from CDB
348 if(!fInitCDBCalled)
349   {AliError("ReadFromCDB: CDB not inizialized. Call InitCDB first");
350    return kFALSE;}
351 AliCDBEntry *cdbEntry = AliCDBManager::Instance()->Get("ITS/PlaneEff/PlaneEffSSD",fRunNumber);
352 if(!cdbEntry) return kFALSE;
353 AliITSPlaneEffSSD* eff= (AliITSPlaneEffSSD*)cdbEntry->GetObject();
354 if(this==eff) return kFALSE;
355 if(fHis) CopyHistos(*eff); // If histos already exist then copy them to eff
356 eff->Copy(*this);          // copy everything (statistics and histos) from eff to this
357 return kTRUE;
358 }
359 //_____________________________________________________________________________
360 Bool_t AliITSPlaneEffSSD::AddFromCDB(AliCDBId *cdbId) {
361 AliCDBEntry *cdbEntry=0;
362 if (!cdbId) {
363   if(!fInitCDBCalled)
364     {AliError("ReadFromCDB: CDB not inizialized. Call InitCDB first"); return kFALSE;}
365   cdbEntry = AliCDBManager::Instance()->Get("ITS/PlaneEff/PlaneEffSSD",fRunNumber);
366 } else {
367   cdbEntry = AliCDBManager::Instance()->Get(*cdbId);
368 }
369 if(!cdbEntry) return kFALSE;
370 AliITSPlaneEffSSD* eff= (AliITSPlaneEffSSD*)cdbEntry->GetObject();
371 *this+=*eff;
372 return kTRUE;
373 }
374 //_____________________________________________________________________________
375 UInt_t AliITSPlaneEffSSD::GetKeyFromDetLocCoord(Int_t ilay, Int_t idet,
376                                                 Float_t, Float_t) const {
377 // method to locate a basic block from Detector Local coordinate (to be used in tracking)
378 UInt_t key=999999;
379 if(ilay<4 || ilay>5)
380   {AliError("GetKeyFromDetLocCoord: you asked for a non existing layer");
381    return key;}
382 if(ilay==4 && (idet<0 || idet>747))
383  {AliError("GetKeyFromDetLocCoord: you asked for a non existing detector");
384    return key;}
385 if(ilay==5 && (idet<0 || idet>949))
386  {AliError("GetKeyFromDetLocCoord: you asked for a non existing detector");
387    return key;}
388
389 UInt_t mod=idet;
390 if(ilay==5) mod+=748;
391 key=GetKey(mod);
392 return key;
393 }
394 //__________________________________________________________
395 void AliITSPlaneEffSSD::InitHistos() {
396   // for the moment let's create the histograms
397   // module by  module
398   TString histnameResX="HistResX_mod_",aux;
399   TString histnameResZ="HistResZ_mod_";
400   TString histnameResXZ="HistResXZ_mod_";
401   TString histnameClusterType="HistClusterType_mod_";
402
403 //
404   fHisResX=new TH1F*[kNHisto];
405   fHisResZ=new TH1F*[kNHisto];
406   fHisResXZ=new TH2F*[kNHisto];
407   fHisClusterSize=new TH2I*[kNHisto];
408
409   for (Int_t nhist=0;nhist<kNHisto;nhist++){
410     aux=histnameResX;
411     aux+=nhist;
412     fHisResX[nhist]=new TH1F("histname","histname",500,-0.05,0.05); // +- 500 micron; 1 bin=2 micron
413     fHisResX[nhist]->SetName(aux.Data());
414     fHisResX[nhist]->SetTitle(aux.Data());
415
416     aux=histnameResZ;
417     aux+=nhist;
418     fHisResZ[nhist]=new TH1F("histname","histname",500,-0.50,0.50); // +-5000 micron; 1 bin=20 micron
419     fHisResZ[nhist]->SetName(aux.Data());
420     fHisResZ[nhist]->SetTitle(aux.Data());
421
422     aux=histnameResXZ;
423     aux+=nhist;
424     fHisResXZ[nhist]=new TH2F("histname","histname",40,-0.02,0.02,40,-0.16,0.16); // binning:
425                                                                                    // 10 micron in x;
426                                                                                    // 80 micron in z;
427     fHisResXZ[nhist]->SetName(aux.Data());
428     fHisResXZ[nhist]->SetTitle(aux.Data());
429
430     aux=histnameClusterType;
431     aux+=nhist;
432     fHisClusterSize[nhist]=new TH2I("histname","histname",6,0.5,6.5,6,0.5,6.5);
433     fHisClusterSize[nhist]->SetName(aux.Data());
434     fHisClusterSize[nhist]->SetTitle(aux.Data());
435
436   }
437 return;
438 }
439 //__________________________________________________________
440 void AliITSPlaneEffSSD::DeleteHistos() {
441   if(fHisResX) {
442     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResX[i];
443     delete [] fHisResX; fHisResX=0;
444   }
445   if(fHisResZ) {
446     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResZ[i];
447     delete [] fHisResZ; fHisResZ=0;
448   }
449   if(fHisResXZ) {
450     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResXZ[i];
451     delete [] fHisResXZ; fHisResXZ=0;
452   }
453   if(fHisClusterSize) {
454     for (Int_t i=0; i<kNHisto; i++ ) delete fHisClusterSize[i];
455     delete [] fHisClusterSize; fHisClusterSize=0;
456   }
457
458 return;
459 }
460 //__________________________________________________________
461 Bool_t AliITSPlaneEffSSD::FillHistos(UInt_t key, Bool_t found,
462                                  //    Float_t tXZ[2], Float_t cXZ[2], Int_t ctXZ[2]) {
463                                      Float_t *tr, Float_t *clu, Int_t *csize) {
464 // this method fill the histograms
465 // input: - key: unique key of the basic block
466 //        - found: Boolean to asses whether a cluster has been associated to the track or not
467 //        - tr[0],tr[1] local X and Z coordinates of the track prediction, respectively
468 //        - tr[2],tr[3] error on local X and Z coordinates of the track prediction, respectively
469 //        - clu[0],clu[1] local X and Z coordinates of the cluster associated to the track, respectively
470 //        - clu[2],clu[3] error on local X and Z coordinates of the cluster associated to the track, respectively
471 //        - csize[0][1] cluster size in X and Z, respectively
472 // output: kTRUE if filling was succesfull kFALSE otherwise
473 // side effects: updating of the histograms.
474 //
475   if (!fHis) {
476     AliWarning("FillHistos: histograms do not exist! Call SetCreateHistos(kTRUE) first");
477     return kFALSE;
478   }
479   if(key>=kNModule)
480     {AliWarning("FillHistos: you asked for a non existing key"); return kFALSE;}
481   Int_t id=GetModFromKey(key);
482   if(id>=kNHisto)
483     {AliWarning("FillHistos: you want to fill a non-existing histos"); return kFALSE;}
484   if(found) {
485     Float_t resx=tr[0]-clu[0];
486     Float_t resz=tr[1]-clu[1];
487     fHisResX[id]->Fill(resx);
488     fHisResZ[id]->Fill(resz);
489     fHisResXZ[id]->Fill(resx,resz);
490     fHisClusterSize[id]->Fill((Double_t)csize[0],(Double_t)csize[1]);
491   }
492   return kTRUE;
493 }
494 //__________________________________________________________
495 Bool_t AliITSPlaneEffSSD::WriteHistosToFile(TString filename, Option_t* option) {
496   //
497   // Saves the histograms into a tree and saves the trees into a file
498   //
499   if (!fHis) return kFALSE;
500   if (filename.Data()=="") {
501      AliWarning("WriteHistosToFile: null output filename!");
502      return kFALSE;
503   }
504 //  char branchname[30];
505   TFile *hFile=new TFile(filename.Data(),option,
506                          "The File containing the TREEs with ITS PlaneEff Histos");
507   TTree *SSDTree=new TTree("SSDTree","Tree whith Residuals and Cluster Type distributions for SSD");
508   TH1F *histZ,*histX;
509   TH2F *histXZ;
510   TH2I *histClusterType;
511
512   histZ=new TH1F();
513   histX=new TH1F();
514   histXZ=new TH2F();
515   histClusterType=new TH2I();
516
517   SSDTree->Branch("histX","TH1F",&histX,128000,0);
518   SSDTree->Branch("histZ","TH1F",&histZ,128000,0);
519   SSDTree->Branch("histXZ","TH2F",&histXZ,128000,0);
520   SSDTree->Branch("histClusterType","TH2I",&histClusterType,128000,0);
521
522   for(Int_t j=0;j<kNHisto;j++){
523     histX=fHisResX[j];
524     histZ=fHisResZ[j];
525     histXZ=fHisResXZ[j];
526     histClusterType=fHisClusterSize[j];
527
528     SSDTree->Fill();
529   }
530   hFile->Write();
531   hFile->Close();
532 return kTRUE;
533 }
534 //__________________________________________________________
535 Bool_t AliITSPlaneEffSSD::ReadHistosFromFile(TString filename) {
536   //
537   // Read histograms from an already existing file
538   //
539   if (!fHis) return kFALSE;
540   if (filename.Data()=="") {
541      AliWarning("ReadHistosFromFile: incorrect output filename!");
542      return kFALSE;
543   }
544   //char branchname[30];
545
546   TH1F *h  = 0;
547   TH2F *h2 = 0;
548   TH2I *h2i= 0;
549
550   TFile *file=TFile::Open(filename.Data(),"READONLY");
551
552   if (!file || file->IsZombie()) {
553     AliWarning(Form("Can't open %s !",filename.Data()));
554     delete file;
555     return kFALSE;
556   }
557   TTree *tree = (TTree*) file->Get("SSDTree");
558
559   TBranch *histX = (TBranch*) tree->GetBranch("histX");
560   TBranch *histZ = (TBranch*) tree->GetBranch("histZ");
561   TBranch *histXZ = (TBranch*) tree->GetBranch("histXZ");
562   TBranch *histClusterType = (TBranch*) tree->GetBranch("histClusterType");
563
564   gROOT->cd();
565
566   Int_t nevent = (Int_t)histX->GetEntries();
567   if(nevent!=kNHisto)
568     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
569   histX->SetAddress(&h);
570   for(Int_t j=0;j<kNHisto;j++){
571     delete h; h=0;
572     histX->GetEntry(j);
573     fHisResX[j]->Add(h);
574   }
575
576   nevent = (Int_t)histZ->GetEntries();
577   if(nevent!=kNHisto)
578     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
579   histZ->SetAddress(&h);
580   for(Int_t j=0;j<kNHisto;j++){
581     delete h; h=0;
582     histZ->GetEntry(j);
583     fHisResZ[j]->Add(h);
584   }
585
586   nevent = (Int_t)histXZ->GetEntries();
587   if(nevent!=kNHisto)
588     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
589   histXZ->SetAddress(&h2);
590   for(Int_t j=0;j<kNHisto;j++){
591     delete h2; h2=0;
592     histXZ->GetEntry(j);
593     fHisResXZ[j]->Add(h2);
594   }
595
596   nevent = (Int_t)histClusterType->GetEntries();
597   if(nevent!=kNHisto)
598     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
599   histClusterType->SetAddress(&h2i);
600   for(Int_t j=0;j<kNHisto;j++){
601     delete h2i; h2i=0;
602     histClusterType->GetEntry(j);
603     fHisClusterSize[j]->Add(h2i);
604   }
605
606   delete h;   h=0;
607   delete h2;  h2=0;
608   delete h2i; h2i=0;
609
610   if (file) {
611     file->Close();
612   }
613 return kTRUE;
614 }
615