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