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