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