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