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