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