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