]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSPlaneEffSPD.cxx
Adding a reminder for coders
[u/mrichter/AliRoot.git] / ITS / AliITSPlaneEffSPD.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 of the SPD,        
18 //  evaluated by tracks
19 //  (Inherits from AliITSPlaneEff)
20 //  Author: G.E. Bruno 
21 //          giuseppe.bruno@ba.infn.it
22 //
23 ///////////////////////////////////////////////////////////////////////////
24
25 /* $Id$ */
26
27 #include <TMath.h>
28 #include <TH1F.h>
29 #include <TFile.h>
30 #include <TTree.h>
31 #include <TROOT.h>
32 #include "AliITSPlaneEffSPD.h"
33 #include "AliLog.h"
34 #include "AliCDBStorage.h"
35 #include "AliCDBEntry.h"
36 #include "AliCDBManager.h"
37 //#include "AliCDBRunRange.h"
38 #include "AliITSsegmentationSPD.h"
39 #include "AliITSCalibrationSPD.h"
40
41 ClassImp(AliITSPlaneEffSPD)     
42 //______________________________________________________________________
43 AliITSPlaneEffSPD::AliITSPlaneEffSPD():
44   AliITSPlaneEff(),
45   fHisResX(0),
46   fHisResZ(0),
47   fHisResXZ(0),
48   fHisClusterSize(0),
49   fHisResXclu(0),
50   fHisResZclu(0),
51   fHisResXchip(0),
52   fHisResZchip(0),
53   fProfResXvsPhi(0),
54   fProfResZvsDip(0),
55   fProfResXvsPhiclu(0), 
56   fProfResZvsDipclu(0),
57   fHisTrackErrX(0),
58   fHisTrackErrZ(0),
59   fHisClusErrX(0),
60   fHisClusErrZ(0),
61   fHisTrackXFOtrue(0),
62   fHisTrackZFOtrue(0),
63   fHisTrackXFOfalse(0),
64   fHisTrackZFOfalse(0),
65   fHisTrackXZFOtrue(0),
66   fHisTrackXZFOfalse(0){
67   for (UInt_t i=0; i<kNModule*kNChip*(kNClockPhase+1); i++){
68     fFound[i]=0;
69     fTried[i]=0;
70   }
71   // default constructor
72   AliDebug(1,Form("Calling default constructor"));
73 }
74 //______________________________________________________________________
75 AliITSPlaneEffSPD::~AliITSPlaneEffSPD(){
76     // destructor
77     // Inputs:
78     //    none.
79     // Outputs:
80     //    none.
81     // Return:
82     //     none.
83     DeleteHistos();
84 }
85 //______________________________________________________________________
86 AliITSPlaneEffSPD::AliITSPlaneEffSPD(const AliITSPlaneEffSPD &s) : AliITSPlaneEff(s), 
87 //fHis(s.fHis),
88 fHisResX(0),
89 fHisResZ(0),
90 fHisResXZ(0),
91 fHisClusterSize(0),
92 fHisResXclu(0),
93 fHisResZclu(0),
94 fHisResXchip(0),
95 fHisResZchip(0),
96 fProfResXvsPhi(0),
97 fProfResZvsDip(0),
98 fProfResXvsPhiclu(0),
99 fProfResZvsDipclu(0),
100 fHisTrackErrX(0),
101 fHisTrackErrZ(0),
102 fHisClusErrX(0),
103 fHisClusErrZ(0),
104 fHisTrackXFOtrue(0),
105 fHisTrackZFOtrue(0),
106 fHisTrackXFOfalse(0),
107 fHisTrackZFOfalse(0),
108 fHisTrackXZFOtrue(0),
109 fHisTrackXZFOfalse(0)
110 {
111     //     Copy Constructor
112     // Inputs:
113     //    AliITSPlaneEffSPD &s The original class for which
114     //                                this class is a copy of
115     // Outputs:
116     //    none.
117     // Return:
118
119  for (UInt_t i=0; i<kNModule*kNChip*(kNClockPhase+1); i++){
120     fFound[i]=s.fFound[i];
121     fTried[i]=s.fTried[i];
122  }
123  if(fHis) { 
124    InitHistos();
125    for(Int_t i=0; i<kNHisto; i++) {
126       s.fHisResX[i]->Copy(*fHisResX[i]);
127       s.fHisResZ[i]->Copy(*fHisResZ[i]);
128       s.fHisResXZ[i]->Copy(*fHisResXZ[i]);
129       s.fHisClusterSize[i]->Copy(*fHisClusterSize[i]);
130       for(Int_t clu=0; clu<kNclu; clu++) {  // clu=0 --> cluster size 1
131         s.fHisResXclu[i][clu]->Copy(*fHisResXclu[i][clu]);
132         s.fHisResZclu[i][clu]->Copy(*fHisResZclu[i][clu]);
133         s.fProfResXvsPhiclu[i][clu]->Copy(*fProfResXvsPhiclu[i][clu]);
134         s.fProfResZvsDipclu[i][clu]->Copy(*fProfResZvsDipclu[i][clu]);
135       }
136       for(Int_t chip=0; chip<kNChip; chip++) { 
137         s.fHisResXchip[i][chip]->Copy(*fHisResXchip[i][chip]);
138         s.fHisResZchip[i][chip]->Copy(*fHisResZchip[i][chip]);
139       }
140       s.fProfResXvsPhi[i]->Copy(*fProfResXvsPhi[i]);
141       s.fProfResZvsDip[i]->Copy(*fProfResZvsDip[i]);
142       s.fHisTrackErrX[i]->Copy(*fHisTrackErrX[i]);
143       s.fHisTrackErrZ[i]->Copy(*fHisTrackErrZ[i]);
144       s.fHisClusErrX[i]->Copy(*fHisClusErrX[i]);
145       s.fHisClusErrZ[i]->Copy(*fHisClusErrZ[i]);
146       for(Int_t phas=0; phas<kNClockPhase;phas++){
147         s.fHisTrackXFOtrue[i][phas]->Copy(*fHisTrackXFOtrue[i][phas]);
148         s.fHisTrackZFOtrue[i][phas]->Copy(*fHisTrackXFOtrue[i][phas]);
149         s.fHisTrackXFOfalse[i][phas]->Copy(*fHisTrackXFOtrue[i][phas]);
150         s.fHisTrackZFOfalse[i][phas]->Copy(*fHisTrackXFOtrue[i][phas]);
151         s.fHisTrackXZFOtrue[i][phas]->Copy(*fHisTrackXFOtrue[i][phas]);
152         s.fHisTrackXZFOfalse[i][phas]->Copy(*fHisTrackXFOtrue[i][phas]);
153       }
154    }
155  }
156 }
157 //_________________________________________________________________________
158 AliITSPlaneEffSPD& AliITSPlaneEffSPD::operator+=(const AliITSPlaneEffSPD &add){
159     //    Add-to-me operator
160     // Inputs:
161     //    const AliITSPlaneEffSPD &add  simulation class to be added
162     // Outputs:
163     //    none.
164     // Return:
165     //    none
166     for (UInt_t i=0; i<kNModule*kNChip*(kNClockPhase+1); i++){
167       fFound[i] += add.fFound[i];
168       fTried[i] += add.fTried[i];
169     }
170     if(fHis && add.fHis) {
171       for(Int_t i=0; i<kNHisto; i++) {
172         fHisResX[i]->Add(add.fHisResX[i]); 
173         fHisResZ[i]->Add(add.fHisResZ[i]); 
174         fHisResXZ[i]->Add(add.fHisResXZ[i]); 
175         fHisClusterSize[i]->Add(add.fHisClusterSize[i]); 
176         for(Int_t clu=0; clu<kNclu; clu++) {  // clu=0 --> cluster size 1
177           fHisResXclu[i][clu]->Add(add.fHisResXclu[i][clu]); 
178           fHisResZclu[i][clu]->Add(add.fHisResZclu[i][clu]); 
179           fProfResXvsPhiclu[i][clu]->Add(add.fProfResXvsPhiclu[i][clu]);
180           fProfResZvsDipclu[i][clu]->Add(add.fProfResZvsDipclu[i][clu]);
181         }
182         for(Int_t chip=0; chip<kNChip; chip++) {  
183           fHisResXchip[i][chip]->Add(add.fHisResXchip[i][chip]); 
184           fHisResZchip[i][chip]->Add(add.fHisResZchip[i][chip]); 
185         }
186         fProfResXvsPhi[i]->Add(add.fProfResXvsPhi[i]);
187         fProfResZvsDip[i]->Add(add.fProfResZvsDip[i]);
188         fHisTrackErrX[i]->Add(add.fHisTrackErrX[i]);
189         fHisTrackErrZ[i]->Add(add.fHisTrackErrZ[i]);
190         fHisClusErrX[i]->Add(add.fHisClusErrX[i]);
191         fHisClusErrZ[i]->Add(add.fHisClusErrZ[i]);
192         for(Int_t phas=0; phas<kNClockPhase;phas++){
193           fHisTrackXFOtrue[i][phas]->Add(add.fHisTrackXFOtrue[i][phas]);
194           fHisTrackZFOtrue[i][phas]->Add(add.fHisTrackXFOtrue[i][phas]);
195           fHisTrackXFOfalse[i][phas]->Add(add.fHisTrackXFOtrue[i][phas]);
196           fHisTrackZFOfalse[i][phas]->Add(add.fHisTrackXFOtrue[i][phas]);
197           fHisTrackXZFOtrue[i][phas]->Add(add.fHisTrackXFOtrue[i][phas]);
198           fHisTrackXZFOfalse[i][phas]->Add(add.fHisTrackXFOtrue[i][phas]);
199         }
200       }
201     }
202     return *this;
203 }
204 //______________________________________________________________________
205 AliITSPlaneEffSPD&  AliITSPlaneEffSPD::operator=(const
206                                            AliITSPlaneEffSPD &s){
207     //    Assignment operator
208     // Inputs:
209     //    AliITSPlaneEffSPD &s The original class for which
210     //                                this class is a copy of
211     // Outputs:
212     //    none.
213     // Return:
214  
215     if(this==&s) return *this;
216     s.Copy(*this);
217     return *this;
218 }
219 //______________________________________________________________________
220 void AliITSPlaneEffSPD::Copy(TObject &obj) const {
221   // protected method. copy this to obj
222   AliITSPlaneEff::Copy(obj);
223   AliITSPlaneEffSPD& target = (AliITSPlaneEffSPD &) obj;
224   for(Int_t i=0;i<kNModule*kNChip*(kNClockPhase+1);i++) {
225       target.fFound[i] = fFound[i];
226       target.fTried[i] = fTried[i];
227   }
228   CopyHistos(target);
229   return;
230 }
231 //_______________________________________________________________________
232 void AliITSPlaneEffSPD::CopyHistos(AliITSPlaneEffSPD &target) const {
233   // protected method: copy histos from this to target
234   target.fHis  = fHis; // this is redundant only in some cases. Leave as it is.
235   if(fHis) {
236     target.fHisResX=new TH1F*[kNHisto];
237     target.fHisResZ=new TH1F*[kNHisto];
238     target.fHisResXZ=new TH2F*[kNHisto];
239     target.fHisClusterSize=new TH2I*[kNHisto];
240     target.fHisResXclu=new TH1F**[kNHisto];
241     target.fHisResZclu=new TH1F**[kNHisto];
242     target.fHisResXchip=new TH1F**[kNHisto];
243     target.fHisResZchip=new TH1F**[kNHisto];
244     target.fProfResXvsPhi=new TProfile*[kNHisto];
245     target.fProfResZvsDip=new TProfile*[kNHisto];
246     target.fProfResXvsPhiclu=new TProfile**[kNHisto];
247     target.fProfResZvsDipclu=new TProfile**[kNHisto];
248     target.fHisTrackErrX=new TH1F*[kNHisto];
249     target.fHisTrackErrZ=new TH1F*[kNHisto];
250     target.fHisClusErrX=new TH1F*[kNHisto];
251     target.fHisClusErrZ=new TH1F*[kNHisto];
252     target.fHisTrackXFOtrue=new TH1F**[kNHisto];
253     target.fHisTrackZFOtrue=new TH1F**[kNHisto];
254     target.fHisTrackXFOfalse=new TH1F**[kNHisto];
255     target.fHisTrackZFOfalse=new TH1F**[kNHisto];
256     target.fHisTrackXZFOtrue=new TH2F**[kNHisto];
257     target.fHisTrackXZFOfalse=new TH2F**[kNHisto];
258     for(Int_t i=0; i<kNHisto; i++) {
259       target.fHisResX[i] = new TH1F(*fHisResX[i]);
260       target.fHisResZ[i] = new TH1F(*fHisResZ[i]);
261       target.fHisResXZ[i] = new TH2F(*fHisResXZ[i]);
262       target.fHisClusterSize[i] = new TH2I(*fHisClusterSize[i]);
263       target.fHisResXclu[i]=new TH1F*[kNclu];
264       target.fHisResZclu[i]=new TH1F*[kNclu];
265       target.fProfResXvsPhiclu[i]=new TProfile*[kNclu];
266       target.fProfResZvsDipclu[i]=new TProfile*[kNclu];
267       for(Int_t clu=0; clu<kNclu; clu++) {  // clu=0 --> cluster size 1
268         target.fHisResXclu[i][clu] = new TH1F(*fHisResXclu[i][clu]);
269         target.fHisResZclu[i][clu] = new TH1F(*fHisResZclu[i][clu]);
270         target.fProfResXvsPhiclu[i][clu] = new TProfile(*fProfResXvsPhiclu[i][clu]);
271         target.fProfResZvsDipclu[i][clu] = new TProfile(*fProfResZvsDipclu[i][clu]);
272       }
273       target.fHisResXchip[i]=new TH1F*[kNChip];
274       target.fHisResZchip[i]=new TH1F*[kNChip];
275       for(Int_t chip=0; chip<kNChip; chip++) {  
276         target.fHisResXchip[i][chip] = new TH1F(*fHisResXchip[i][chip]);
277         target.fHisResZchip[i][chip] = new TH1F(*fHisResZchip[i][chip]);
278       }
279       target.fProfResXvsPhi[i] = new TProfile(*fProfResXvsPhi[i]);
280       target.fProfResZvsDip[i] = new TProfile(*fProfResZvsDip[i]);
281       target.fHisTrackErrX[i] = new TH1F(*fHisTrackErrX[i]);
282       target.fHisTrackErrZ[i] = new TH1F(*fHisTrackErrZ[i]);
283       target.fHisClusErrX[i] = new TH1F(*fHisClusErrX[i]);
284       target.fHisClusErrZ[i] = new TH1F(*fHisClusErrZ[i]);
285
286       target.fHisTrackXFOtrue[i]=new TH1F*[kNClockPhase];
287       target.fHisTrackZFOtrue[i]=new TH1F*[kNClockPhase];
288       target.fHisTrackXFOfalse[i]=new TH1F*[kNClockPhase];
289       target.fHisTrackZFOfalse[i]=new TH1F*[kNClockPhase];
290       target.fHisTrackXZFOtrue[i]=new TH2F*[kNClockPhase];
291       target.fHisTrackXZFOfalse[i]=new TH2F*[kNClockPhase];
292       for(Int_t phas=0; phas<kNClockPhase;phas++){
293       target.fHisTrackXFOtrue[i][phas]=new TH1F(*fHisTrackXFOtrue[i][phas]);
294       target.fHisTrackZFOtrue[i][phas]=new TH1F(*fHisTrackZFOtrue[i][phas]);
295       target.fHisTrackXFOfalse[i][phas]=new TH1F(*fHisTrackXFOfalse[i][phas]);
296       target.fHisTrackZFOfalse[i][phas]=new TH1F(*fHisTrackZFOfalse[i][phas]);
297       target.fHisTrackXZFOtrue[i][phas]=new TH2F(*fHisTrackXZFOtrue[i][phas]);
298       target.fHisTrackXZFOfalse[i][phas]=new TH2F(*fHisTrackXZFOfalse[i][phas]);
299       }
300     }
301   }
302 return;
303 }
304
305 //_______________________________________________________________________
306 Int_t AliITSPlaneEffSPD::GetMissingTracksForGivenEff(Double_t eff, Double_t RelErr,
307           UInt_t im, UInt_t ic) const {
308    
309   //   Estimate the number of tracks still to be collected to attain a 
310   //   given efficiency eff, with relative error RelErr
311   //   Inputs:
312   //         eff    -> Expected efficiency (e.g. those from actual estimate)
313   //         RelErr -> tollerance [0,1] 
314   //         im     -> module number [0,239]
315   //         ic     -> chip number [0,4]
316   //   Outputs: none
317   //   Return: the estimated n. of tracks 
318   //
319 if (im>=kNModule || ic>=kNChip) 
320  {AliError("GetMissingTracksForGivenEff: you asked for a non existing chip");
321  return -1;}
322 else { 
323   UInt_t key=GetKey(im,ic);
324   if(key<kNModule*kNChip) return GetNTracksForGivenEff(eff,RelErr)-fTried[key];
325   else return -1;
326 }
327 }
328 //_________________________________________________________________________
329 Double_t  AliITSPlaneEffSPD::PlaneEff(const UInt_t im,const UInt_t ic, const Bool_t fo, const UInt_t bcm4) const {
330 // Compute the efficiency for a basic block, 
331 // Inputs:
332 //        im     -> module number [0,239]
333 //        ic     -> chip number [0,4] 
334 //        fo     -> boolean, true in case of Fast Or studies
335 //        bcm4   -> for Fast Or: bunch crossing % 4
336 if (im>=kNModule || ic>=kNChip) 
337  {AliError("PlaneEff(Uint_t,Uint_t): you asked for a non existing chip"); return -1.;}
338 if(fo && bcm4>=kNClockPhase)
339  {AliError("PlaneEff(Uint_t,Uint_t): you asked for Fast Or in a wrong phase"); return -1.;}
340 Int_t nf=-1;
341 Int_t nt=-1;
342 if(fo) {
343  AliWarning("PlaneEff: you asked for FO efficiency");
344  UInt_t key=GetKey(im,ic,fo,bcm4);
345  if(key<kNModule*kNChip*(kNClockPhase+1)) {
346    nf=fFound[key];
347    nt=fTried[key];
348  }
349 } else {
350  UInt_t key=GetKey(im,ic);
351  if (key<kNModule*kNChip) {
352   nf=fFound[key];
353   nt=fTried[key];
354  }
355 }
356 return AliITSPlaneEff::PlaneEff(nf,nt);
357 }
358 //_________________________________________________________________________
359 Double_t  AliITSPlaneEffSPD::ErrPlaneEff(const UInt_t im,const UInt_t ic, const Bool_t fo, const UInt_t bcm4) const {
360     // Compute the statistical error on efficiency for a basic block,
361     // using binomial statistics 
362     // Inputs:
363     //        im     -> module number [0,239]
364     //        ic     -> chip number [0,4] 
365 //        fo     -> boolean, true in case of Fast Or studies
366 //        bcm4   -> for Fast Or: bunch crossing % 4
367 if (im>=kNModule || ic>=kNChip) 
368  {AliError("ErrPlaneEff(Uint_t,Uint_t): you asked for a non existing chip"); return -1.;}
369 if(fo && bcm4>=kNClockPhase)
370  {AliError("PlaneEff(Uint_t,Uint_t): you asked for Fast Or in a wrong phase"); return -1.;}
371 Int_t nf=-1;
372 Int_t nt=-1;
373 if(fo) {
374  AliWarning("ErrPlaneEff: you asked for FO efficiency");
375  UInt_t key=GetKey(im,ic,fo,bcm4);
376  if(key<kNModule*kNChip*(kNClockPhase+1)) {
377    nf=fFound[key];
378    nt=fTried[key];
379  }
380 } else {
381  UInt_t key=GetKey(im,ic);
382  if (key<kNModule*kNChip) {
383    nf=fFound[key];
384    nt=fTried[key];
385  }
386 }
387 return AliITSPlaneEff::ErrPlaneEff(nf,nt);
388
389 //_________________________________________________________________________
390 Bool_t AliITSPlaneEffSPD::UpDatePlaneEff(const Bool_t Kfound,
391                                          const UInt_t im, const UInt_t ic, const Bool_t fo, const UInt_t bcm4) {
392   // Update efficiency for a basic block
393 if (im>=kNModule || ic>=kNChip) 
394  {AliError("UpDatePlaneEff: you asked for a non existing chip"); return kFALSE;}
395 if(fo && bcm4>=kNClockPhase)
396  {AliError("UpDatePlaneEff: you asked for Fast Or in a wrong phase"); return kFALSE;}
397 if (!fo) {
398  UInt_t key=GetKey(im,ic);
399  if(key<kNModule*kNChip) {
400    fTried[key]++;
401    if(Kfound) fFound[key]++;
402    return kTRUE;
403  }
404 }
405 else {
406  UInt_t key=GetKey(im,ic,fo,bcm4);
407  if(key<kNModule*kNChip*(kNClockPhase+1)) {
408    fTried[key]++;
409    if(Kfound) fFound[key]++;
410    return kTRUE;
411  }
412 }
413 return kFALSE;
414 }
415 //_________________________________________________________________________
416 UInt_t AliITSPlaneEffSPD::GetChipFromCol(const UInt_t col) const {
417   // get chip given the column
418 if(col>=kNCol*kNChip) 
419  {AliDebug(1,Form("GetChipFromCol: you asked for a non existing column %d",col)); return 10;}
420 return col/kNCol;
421 }
422 //__________________________________________________________________________
423 UInt_t AliITSPlaneEffSPD::GetKey(const UInt_t mod, const UInt_t chip, const Bool_t FO, const UInt_t BCm4) const {
424   // get key given a basic block
425 UInt_t key=99999;
426 if(mod>=kNModule || chip>=kNChip)
427   {AliWarning("GetKey: you asked for a non existing block"); return 99999;}
428 key = mod*kNChip+chip;
429 if(FO) { 
430   if(BCm4>= kNClockPhase) {AliWarning("GetKey: you have asked Fast OR and a non exisiting BC modulo 4"); return 99999;}
431   key += kNModule*kNChip*(BCm4+1);
432 }
433 return key;
434 }
435 //__________________________________________________________________________
436 UInt_t AliITSPlaneEffSPD::SwitchChipKeyNumbering(UInt_t key) const {
437
438 // methods to switch from offline chip key numbering 
439 // to online Raw Stream chip numbering and viceversa. 
440 // Used for Fast-Or studies.
441 // Implemented by valerio.altini@ba.infn.it
442
443 if(key>=kNModule*kNChip*(kNClockPhase+1))
444   {AliWarning("SwitchChipKeyNumbering: you asked for a non existing key"); return 99999;}
445 UInt_t mod=9999,chip=9999,phase=9999;
446 GetModAndChipFromKey(key,mod,chip);
447 if(mod<kNModuleLy1) chip = kNChip-(chip+1);
448 if(IsForFO(key))phase = GetBCm4FromKey(key);
449
450 return GetKey(mod,chip,IsForFO(key),phase);
451
452 }
453 //__________________________________________________________________________
454 UInt_t AliITSPlaneEffSPD::GetModFromKey(const UInt_t key) const {
455   // get mod. from key
456 if(key>=kNModule*kNChip*(kNClockPhase+1))
457   {AliError("GetModFromKey: you asked for a non existing key"); return 9999;}
458 return (key%(kNModule*kNChip))/kNChip;
459 }
460 //__________________________________________________________________________
461 UInt_t AliITSPlaneEffSPD::GetChipFromKey(const UInt_t key) const {
462   // retrieves chip from key
463 if(key>=kNModule*kNChip*(kNClockPhase+1))
464   {AliError("GetChipFromKey: you asked for a non existing key"); return 999;}
465 return ((key%(kNModule*kNChip))%(kNModule*kNChip))%kNChip;
466 }
467 //__________________________________________________________________________
468 UInt_t AliITSPlaneEffSPD::GetBCm4FromKey(const UInt_t key) const {
469   // retrieves the "Bunch Crossing modulo 4" (for Fast Or studies)
470 if(key>=kNModule*kNChip*(kNClockPhase+1))
471   {AliError("GetBCm4FromKey: you asked for a non existing key"); return 999;}
472 if(key<kNModule*kNChip) 
473   {AliWarning("GetBCm4FromKey: key is below 1200, why are you asking for FO related stuff"); return 999;}
474
475 return key/(kNModule*kNChip) - 1 ;
476 }
477 //__________________________________________________________________________
478 Bool_t AliITSPlaneEffSPD::IsForFO(const UInt_t key) const {
479 if(key>=kNModule*kNChip) return kTRUE;
480 else return kFALSE;
481 }
482 //__________________________________________________________________________
483 void AliITSPlaneEffSPD::GetModAndChipFromKey(const UInt_t key,UInt_t& mod,UInt_t& chip) const {
484   // get module and chip from a key
485 if(key>=kNModule*kNChip*(kNClockPhase+1))
486   {AliError("GetModAndChipFromKey: you asked for a non existing key"); 
487   mod=9999;
488   chip=999;
489   return;}
490 mod=GetModFromKey(key);
491 chip=GetChipFromKey(key);
492 return;
493 }
494 //____________________________________________________________________________
495 Double_t AliITSPlaneEffSPD::LivePlaneEff(UInt_t key) const {
496   // returns plane efficieny after adding the fraction of sensor which is bad
497 if(key>=kNModule*kNChip)
498   {AliError("LivePlaneEff: you asked for a non existing key");
499    return -1.;}
500 Double_t leff=AliITSPlaneEff::LivePlaneEff(0); // this just for the Warning
501 leff=PlaneEff(key)+GetFracBad(key);
502 return leff>1?1:leff;
503 }
504 //____________________________________________________________________________
505 Double_t AliITSPlaneEffSPD::ErrLivePlaneEff(UInt_t key) const {
506   // returns error on live plane efficiency
507 if(key>=kNModule*kNChip)
508   {AliError("ErrLivePlaneEff: you asked for a non existing key");
509    return -1.;}
510 Int_t nf=fFound[key];
511 Double_t triedInLive=GetFracLive(key)*fTried[key];
512 Int_t nt=TMath::Max(nf,TMath::Nint(triedInLive));
513 return AliITSPlaneEff::ErrPlaneEff(nf,nt); // for the time being: to be checked
514 }
515 //_____________________________________________________________________________
516 Double_t AliITSPlaneEffSPD::GetFracLive(const UInt_t key) const {
517   // returns the fraction of the sensor which is OK
518 if(key>=kNModule*kNChip)
519   {AliError("GetFracLive: you asked for a non existing key");
520    return -1.;}
521     // Compute the fraction of bad (dead+noisy) detector 
522 UInt_t dead=0,noisy=0;
523 GetDeadAndNoisyInChip(key,dead,noisy);
524 Double_t live=dead+noisy;
525 live/=(kNRow*kNCol);
526 return 1.-live;
527 }
528 //_____________________________________________________________________________
529 void AliITSPlaneEffSPD::GetDeadAndNoisyInChip(const UInt_t key,
530       UInt_t& nrDeadInChip, UInt_t& nrNoisyInChip) const {
531   // returns the number of dead and noisy pixels
532 nrDeadInChip=0;
533 nrNoisyInChip=0;
534 if(key>=kNModule*kNChip)
535   {AliError("GetDeadAndNoisyInChip: you asked for a non existing key");
536    return;}
537     // Compute the number of bad (dead+noisy) pixel in a chip
538 //
539 if(!fInitCDBCalled) 
540   {AliError("GetDeadAndNoisyInChip: CDB not inizialized: call InitCDB first");
541    return;};
542 AliCDBManager* man = AliCDBManager::Instance();
543 // retrieve map of dead Pixel 
544 AliCDBEntry *cdbSPDDead = man->Get("ITS/Calib/SPDDead", fRunNumber);
545 TObjArray* spdDead;
546 if(cdbSPDDead) {
547   spdDead = (TObjArray*)cdbSPDDead->GetObject();
548   if(!spdDead) 
549   {AliError("GetDeadAndNoisyInChip: SPDDead not found in CDB");
550    return;}
551 } else {
552   AliError("GetDeadAndNoisyInChip: did not find Calib/SPDDead.");
553   return;
554 }
555 // retrieve map of sparse dead Pixel 
556 AliCDBEntry *cdbSPDSparseDead = man->Get("ITS/Calib/SPDSparseDead", fRunNumber);
557 TObjArray* spdSparseDead;
558 if(cdbSPDSparseDead) {
559   spdSparseDead = (TObjArray*)cdbSPDSparseDead->GetObject();
560   if(!spdSparseDead) 
561   {AliError("GetDeadAndNoisyInChip: SPDSparseDead not found in CDB");
562    return;}
563 } else {
564   AliError("GetDeadAndNoisyInChip: did not find Calib/SPDSparseDead.");
565   return;
566 }
567
568 // retrieve map of noisy Pixel 
569 AliCDBEntry *cdbSPDNoisy = man->Get("ITS/Calib/SPDNoisy", fRunNumber);
570 TObjArray* spdNoisy;
571 if(cdbSPDNoisy) {
572   spdNoisy = (TObjArray*)cdbSPDNoisy->GetObject();
573   if(!spdNoisy) 
574   {AliError("GetDeadAndNoisyInChip: SPDNoisy not found in CDB");
575    return;}
576 } else {
577   AliError("GetDeadAndNoisyInChip: did not find Calib/SPDNoisy.");
578   return;
579 }
580 //
581 UInt_t mod=GetModFromKey(key);
582 UInt_t chip=GetChipFromKey(key);
583 // count number of dead
584 AliITSCalibrationSPD* calibSPD=(AliITSCalibrationSPD*) spdDead->At(mod);
585 UInt_t nrDead = calibSPD->GetNrBad();
586 for (UInt_t index=0; index<nrDead; index++) {
587   if(GetChipFromCol(calibSPD->GetBadColAt(index))==chip) nrDeadInChip++;
588 }
589 // add the number of sparse dead to the previous dead
590 calibSPD=(AliITSCalibrationSPD*) spdSparseDead->At(mod);
591 UInt_t nrSparseDead = calibSPD->GetNrBad();
592 for (UInt_t index=0; index<nrSparseDead; index++) {
593   if(GetChipFromCol(calibSPD->GetBadColAt(index))==chip) nrDeadInChip++;
594 }
595 calibSPD=(AliITSCalibrationSPD*) spdNoisy->At(mod);
596 UInt_t nrNoisy = calibSPD->GetNrBad();
597 for (UInt_t index=0; index<nrNoisy; index++) {
598   if(GetChipFromCol(calibSPD->GetBadColAt(index))==chip) nrNoisyInChip++;
599 }
600 return;
601 }
602 //_____________________________________________________________________________
603 Double_t AliITSPlaneEffSPD::GetFracBad(const UInt_t key) const {
604   // returns 1-fractional live
605 if(key>=kNModule*kNChip)
606   {AliError("GetFracBad: you asked for a non existing key");
607    return -1.;}
608 return 1.-GetFracLive(key);
609 }
610 //_____________________________________________________________________________
611 Bool_t AliITSPlaneEffSPD::WriteIntoCDB() const {
612 // write onto CDB
613 if(!fInitCDBCalled)
614   {AliError("WriteIntoCDB: CDB not inizialized. Call InitCDB first");
615    return kFALSE;}
616 // to be written properly: now only for debugging 
617   AliCDBMetaData *md= new AliCDBMetaData(); // metaData describing the object
618   //md->SetObjectClassName("AliITSPlaneEff");
619   md->SetResponsible("Giuseppe Eugenio Bruno");
620   md->SetBeamPeriod(0);
621   md->SetAliRootVersion("head 19/11/07"); //root version
622   AliCDBId id("ITS/PlaneEff/PlaneEffSPD",0,AliCDBRunRange::Infinity()); 
623   AliITSPlaneEffSPD eff; 
624   eff=*this;
625   Bool_t r=AliCDBManager::Instance()->GetDefaultStorage()->Put(&eff,id,md);
626   delete md;
627   return r;
628 }
629 //_____________________________________________________________________________
630 Bool_t AliITSPlaneEffSPD::ReadFromCDB() {
631 // read from CDB
632 if(!fInitCDBCalled)
633   {AliError("ReadFromCDB: CDB not inizialized. Call InitCDB first");
634    return kFALSE;}
635 AliCDBEntry *cdbEntry = AliCDBManager::Instance()->Get("ITS/PlaneEff/PlaneEffSPD",fRunNumber);
636 if(!cdbEntry) return kFALSE;
637 AliITSPlaneEffSPD* eff= (AliITSPlaneEffSPD*)cdbEntry->GetObject();
638 if(this==eff) return kFALSE;
639 if(fHis) CopyHistos(*eff); // If histos already exist then copy them to eff
640 eff->Copy(*this);          // copy everything (statistics and histos) from eff to this
641 return kTRUE;
642 }
643 //_____________________________________________________________________________
644 Bool_t AliITSPlaneEffSPD::AddFromCDB(AliCDBId *cdbId) {
645 AliCDBEntry *cdbEntry=0;
646 if (!cdbId) {
647   if(!fInitCDBCalled)  
648     {AliError("ReadFromCDB: CDB not inizialized. Call InitCDB first"); return kFALSE;}
649   cdbEntry = AliCDBManager::Instance()->Get("ITS/PlaneEff/PlaneEffSPD",fRunNumber);
650 } else {
651   cdbEntry = AliCDBManager::Instance()->Get(*cdbId);
652 }
653 if(!cdbEntry) return kFALSE;
654 AliITSPlaneEffSPD* eff= (AliITSPlaneEffSPD*)cdbEntry->GetObject();
655 *this+=*eff;
656 return kTRUE;
657 }
658 //_____________________________________________________________________________
659 UInt_t AliITSPlaneEffSPD::GetKeyFromDetLocCoord(Int_t ilay, Int_t idet, 
660                                                 Float_t, Float_t locz) const {
661 // method to locate a basic block from Detector Local coordinate (to be used in tracking)
662 UInt_t key=999999;
663 if(ilay<0 || ilay>1) 
664   {AliError("GetKeyFromDetLocCoord: you asked for a non existing layer");
665    return key;}
666 if(ilay==0 && (idet<0 || idet>79))
667  {AliError("GetKeyFromDetLocCoord: you asked for a non existing detector");
668    return key;}
669 if(ilay==1 && (idet<0 || idet>159))
670  {AliError("GetKeyFromDetLocCoord: you asked for a non existing detector");
671    return key;}
672
673 UInt_t mod=idet;
674 if(ilay==1) mod+=80;
675 key=GetKey(mod,GetChipFromCol(GetColFromLocZ(locz)));
676 return key;
677 }
678 //_____________________________________________________________________________
679 UInt_t AliITSPlaneEffSPD::GetColFromLocZ(Float_t zloc) const {
680 // method to retrieve column number from the local z coordinate
681   UInt_t col=0;
682   AliITSsegmentationSPD spd;
683   Int_t ix,iz;
684   if(spd.LocalToDet(0,zloc,ix,iz)) col+=iz;
685   else {
686     AliDebug(1,Form("cannot compute column number from local z=%f",zloc));
687     col=99999;}
688   return col;
689 /*
690 const Float_t kconv = 1.0E-04; // converts microns to cm.
691 Float_t bz[160];
692 for(Int_t i=000;i<160;i++) bz[i] = 425.0; // most are 425 microns except below
693 bz[ 31] = bz[ 32] = 625.0; // first chip boundry
694 bz[ 63] = bz[ 64] = 625.0; // first chip boundry
695 bz[ 95] = bz[ 96] = 625.0; // first chip boundry
696 bz[127] = bz[128] = 625.0; // first chip boundry
697 //
698 Int_t j=-1;
699 Float_t dz=0;
700 for(Int_t i=000;i<160;i++) dz+=bz[i];
701 dz = -0.5*kconv*dz;
702 if(zloc<dz || zloc>-1*dz) { // outside z range
703   AliDebug(1,Form("GetColFromLocZ: cannot compute column number from local z=%f",zloc));
704   return 99999;}
705 for(j=0;j<160;j++){
706   dz += kconv*bz[j];
707   if(zloc<dz) break;
708 } // end for j
709 col+=j;
710 //
711 return col;
712 */
713 }
714 //________________________________________________________
715 Bool_t AliITSPlaneEffSPD::GetBlockBoundaries(const UInt_t key, Float_t& xmn,Float_t& xmx,
716                                              Float_t& zmn,Float_t& zmx) const {
717 //
718 //  This method return the geometrical boundaries of the active volume of a given 
719 //  basic block, in the detector reference system.
720 //  Input: unique key to locate a basic block.
721 //  
722 //  Output: Ymin, Ymax, Zmin, Zmax of a basic block (chip for SPD)
723 //  Return: kTRUE if computation was succesfully, kFALSE otherwise
724 //
725 if(key>=kNModule*kNChip)
726   {AliWarning("GetBlockBoundaries: you asked for a non existing key"); return kFALSE;}
727 UInt_t chip=GetChipFromKey(key);
728 zmn=GetLocZFromCol(chip*kNCol);
729 zmx=GetLocZFromCol((chip+1)*kNCol);
730 xmn=GetLocXFromRow(0);
731 xmx=GetLocXFromRow(kNRow);
732 //
733 Float_t tmp=zmn;
734 if(zmx<zmn) {zmn=zmx; zmx=tmp;}
735 tmp=xmn;
736 if(xmx<xmn) {xmn=xmx; xmx=tmp;}
737 return kTRUE;
738 }
739 //________________________________________________________
740 Float_t AliITSPlaneEffSPD::GetLocXFromRow(const UInt_t row) const {
741 // 
742 //  This method return the local (i.e. detector reference system) lower x coordinate 
743 //  of the row. To get the central value of a given row, you can do 
744 //  1/2*[LocXFromRow(row)+LocXFromRow(row+1)].
745 //
746 //  Input: row number in the range [0,kNRow] 
747 //  Output: lower local X coordinate of this row.
748 //
749 if(row>kNRow)  // not >= ! allow also computation of upper limit of the last row. 
750   {AliError("LocYFromRow: you asked for a non existing row"); return 9999999.;}
751 // Use only AliITSsegmentationSPD
752 AliITSsegmentationSPD spd;
753 Double_t dummy,x;
754 if(row==kNRow) spd.CellBoundries((Int_t)row-1,0,dummy,x,dummy,dummy);
755 else spd.CellBoundries((Int_t)row,0,x,dummy,dummy,dummy);
756 return (Float_t)x;
757
758 }
759 //________________________________________________________
760 Float_t AliITSPlaneEffSPD::GetLocZFromCol(const UInt_t col) const {
761 //
762 //  This method return the local (i.e. detector reference system) lower Z coordinate
763 //  of the column. To get the central value of a given column, you can do
764 //  1/2*[LocZFromCol(col)+LocZFromCol(col+1)].
765 //
766 //  Input: col number in the range [0,kNChip*kNCol]
767 //  Output: lower local Y coordinate of this row.
768 //
769 if(col>kNChip*kNCol) // not >= ! allow also computation of upper limit of the last column
770   {AliError("LocZFromCol: you asked for a non existing column"); return 9999999.;}
771 // Use only AliITSsegmentationSPD
772 AliITSsegmentationSPD spd;
773 Double_t dummy,y;
774 if(col==kNChip*kNCol) spd.CellBoundries(0,(Int_t)col-1,dummy,dummy,dummy,y);
775 else spd.CellBoundries(0,(Int_t)col,dummy,dummy,y,dummy);
776 return (Float_t)y;
777
778 }
779 //__________________________________________________________
780 void AliITSPlaneEffSPD::InitHistos() {
781   // for the moment let's create the histograms 
782   // module by  module
783   TString histnameResX="HistResX_mod_",aux;
784   TString histnameResZ="HistResZ_mod_";
785   TString histnameResXZ="HistResXZ_mod_";
786   TString histnameClusterType="HistClusterType_mod_";
787   TString histnameResXclu="HistResX_mod_";
788   TString histnameResZclu="HistResZ_mod_";
789   TString histnameResXchip="HistResX_mod_";
790   TString histnameResZchip="HistResZ_mod_";
791   TString profnameResXvsPhi="ProfResXvsPhi_mod_";
792   TString profnameResZvsDip="ProfResZvsDip_mod_";
793   TString profnameResXvsPhiclu="ProfResXvsPhi_mod_";
794   TString profnameResZvsDipclu="ProfResZvsDip_mod_";
795   TString histnameTrackErrX="HistTrackErrX_mod_";
796   TString histnameTrackErrZ="HistTrackErrZ_mod_";
797   TString histnameClusErrX="HistClusErrX_mod_";
798   TString histnameClusErrZ="HistClusErrZ_mod_";
799   TString histnameTrackXFOtrue="HistTrackXFOok_mod_";
800   TString histnameTrackZFOtrue="HistTrackZFOok_mod_";
801   TString histnameTrackXFOfalse="HistTrackXFOko_mod_";
802   TString histnameTrackZFOfalse="HistTrackZFOko_mod_";
803   TString histnameTrackXZFOtrue="HistTrackZvsXFOok_mod_";
804   TString histnameTrackXZFOfalse="HistTrackZvsXFOko_mod_";
805 //
806
807   TH1::AddDirectory(kFALSE);
808
809   fHisResX=new TH1F*[kNHisto];
810   fHisResZ=new TH1F*[kNHisto];
811   fHisResXZ=new TH2F*[kNHisto];
812   fHisClusterSize=new TH2I*[kNHisto];
813   fHisResXclu=new TH1F**[kNHisto];
814   fHisResZclu=new TH1F**[kNHisto];
815   fHisResXchip=new TH1F**[kNHisto];
816   fHisResZchip=new TH1F**[kNHisto];
817   fProfResXvsPhi=new TProfile*[kNHisto];
818   fProfResZvsDip=new TProfile*[kNHisto];
819   fProfResXvsPhiclu=new TProfile**[kNHisto];
820   fProfResZvsDipclu=new TProfile**[kNHisto];
821   fHisTrackErrX=new TH1F*[kNHisto];
822   fHisTrackErrZ=new TH1F*[kNHisto];
823   fHisClusErrX=new TH1F*[kNHisto];
824   fHisClusErrZ=new TH1F*[kNHisto];
825   fHisTrackXFOtrue=new TH1F**[kNHisto];
826   fHisTrackZFOtrue=new TH1F**[kNHisto];
827   fHisTrackXFOfalse=new TH1F**[kNHisto];
828   fHisTrackZFOfalse=new TH1F**[kNHisto];
829   fHisTrackXZFOtrue=new TH2F**[kNHisto];
830   fHisTrackXZFOfalse=new TH2F**[kNHisto];
831
832   for (Int_t nhist=0;nhist<kNHisto;nhist++){
833     aux=histnameResX;
834     aux+=nhist;
835     fHisResX[nhist]=new TH1F("histname","histname",1600,-0.32,0.32); // +- 3200 micron; 1 bin=4 micron
836     fHisResX[nhist]->SetName(aux.Data());
837     fHisResX[nhist]->SetTitle(aux.Data());
838
839     aux=histnameResZ;
840     aux+=nhist;
841     fHisResZ[nhist]=new TH1F("histname","histname",1200,-0.48,0.48); // +-4800 micron; 1 bin=8 micron
842     fHisResZ[nhist]->SetName(aux.Data());
843     fHisResZ[nhist]->SetTitle(aux.Data());
844
845     aux=histnameResXZ;
846     aux+=nhist;
847     fHisResXZ[nhist]=new TH2F("histname","histname",80,-0.16,0.16,80,-0.32,0.32); // binning:
848     fHisResXZ[nhist]->SetName(aux.Data());                                         // 40 micron in x;
849     fHisResXZ[nhist]->SetTitle(aux.Data());                                        // 80 micron in z;
850
851     aux=histnameClusterType;
852     aux+=nhist;
853     fHisClusterSize[nhist]=new TH2I("histname","histname",10,0.5,10.5,10,0.5,10.5);
854     fHisClusterSize[nhist]->SetName(aux.Data());
855     fHisClusterSize[nhist]->SetTitle(aux.Data());
856
857     fHisResXclu[nhist]=new TH1F*[kNclu];
858     fHisResZclu[nhist]=new TH1F*[kNclu];
859     fHisTrackXFOtrue[nhist]=new TH1F*[kNClockPhase];
860     fHisTrackZFOtrue[nhist]=new TH1F*[kNClockPhase];
861     fHisTrackXFOfalse[nhist]=new TH1F*[kNClockPhase];
862     fHisTrackZFOfalse[nhist]=new TH1F*[kNClockPhase];
863     fHisTrackXZFOtrue[nhist]=new TH2F*[kNClockPhase];
864     fHisTrackXZFOfalse[nhist]=new TH2F*[kNClockPhase];
865
866     for(Int_t clu=0; clu<kNclu; clu++) {  // clu=0 --> cluster size 1
867       aux=histnameResXclu;
868       aux+=nhist;
869       aux+="_clu_";
870       aux+=clu+1; // clu=0 --> cluster size 1
871       fHisResXclu[nhist][clu]=new TH1F("histname","histname",1600,-0.32,0.32); // +- 3200 micron; 1 bin=4 micron
872       fHisResXclu[nhist][clu]->SetName(aux.Data());
873       fHisResXclu[nhist][clu]->SetTitle(aux.Data());
874
875       aux=histnameResZclu;
876       aux+=nhist;
877       aux+="_clu_";
878       aux+=clu+1; // clu=0 --> cluster size 1
879       fHisResZclu[nhist][clu]=new TH1F("histname","histname",1200,-0.48,0.48); // +-4800 micron; 1 bin=8 micron
880       fHisResZclu[nhist][clu]->SetName(aux.Data());
881       fHisResZclu[nhist][clu]->SetTitle(aux.Data());
882     }
883
884     fHisResXchip[nhist]=new TH1F*[kNChip];
885     fHisResZchip[nhist]=new TH1F*[kNChip];
886     for(Int_t chip=0; chip<kNChip; chip++) { 
887       aux=histnameResXchip;
888       aux+=nhist;
889       aux+="_chip_";
890       aux+=chip; 
891       fHisResXchip[nhist][chip]=new TH1F("histname","histname",800,-0.32,0.32); // +- 3200 micron; 1 bin=8 micron
892       fHisResXchip[nhist][chip]->SetName(aux.Data());
893       fHisResXchip[nhist][chip]->SetTitle(aux.Data());
894
895       aux=histnameResZchip;
896       aux+=nhist;
897       aux+="_chip_";
898       aux+=chip;
899       fHisResZchip[nhist][chip]=new TH1F("histname","histname",300,-0.48,0.48); // +-4800 micron; 1 bin=32 micron
900       fHisResZchip[nhist][chip]->SetName(aux.Data());
901       fHisResZchip[nhist][chip]->SetTitle(aux.Data());
902     }
903
904     aux=histnameTrackErrX;
905     aux+=nhist;
906     fHisTrackErrX[nhist]=new TH1F("histname","histname",400,0.,0.32); // 0-3200 micron; 1 bin=8 micron
907     fHisTrackErrX[nhist]->SetName(aux.Data());
908     fHisTrackErrX[nhist]->SetTitle(aux.Data());
909
910     aux=histnameTrackErrZ;
911     aux+=nhist;
912     fHisTrackErrZ[nhist]=new TH1F("histname","histname",200,0.,0.32); // 0-3200 micron; 1 bin=16 micron
913     fHisTrackErrZ[nhist]->SetName(aux.Data());
914     fHisTrackErrZ[nhist]->SetTitle(aux.Data());
915
916     aux=histnameClusErrX;
917     aux+=nhist;
918     fHisClusErrX[nhist]=new TH1F("histname","histname",400,0.,0.08); //  0-800 micron; 1 bin=2 micron
919     fHisClusErrX[nhist]->SetName(aux.Data());
920     fHisClusErrX[nhist]->SetTitle(aux.Data());
921
922     aux=histnameClusErrZ;
923     aux+=nhist;
924     fHisClusErrZ[nhist]=new TH1F("histname","histname",400,0.,0.32); //  0-3200 micron; 1 bin=8 micron
925     fHisClusErrZ[nhist]->SetName(aux.Data());
926     fHisClusErrZ[nhist]->SetTitle(aux.Data());
927
928     aux=profnameResXvsPhi;
929     aux+=nhist;
930     fProfResXvsPhi[nhist]=new TProfile("histname","histname",40,-40.,40.0); // binning: range:  -40°- 40°
931     fProfResXvsPhi[nhist]->SetName(aux.Data());                             //          bin width: 2°
932     fProfResXvsPhi[nhist]->SetTitle(aux.Data());
933
934     aux=profnameResZvsDip;
935     aux+=nhist;
936     fProfResZvsDip[nhist]=new TProfile("histname","histname",48,-72.,72.0); // binning: range:  -70°-4°
937     fProfResZvsDip[nhist]->SetName(aux.Data());                             //          bin width: 3°
938     fProfResZvsDip[nhist]->SetTitle(aux.Data());
939
940     fProfResXvsPhiclu[nhist]=new TProfile*[kNclu];
941     fProfResZvsDipclu[nhist]=new TProfile*[kNclu];
942     for(Int_t clu=0; clu<kNclu; clu++) {  // clu=0 --> cluster size 1
943       aux=profnameResXvsPhiclu;
944       aux+=nhist;
945       aux+="_clu_";
946       aux+=clu+1; // clu=0 --> cluster size 1
947       fProfResXvsPhiclu[nhist][clu]=new TProfile("histname","histname",40,-40.,40.0); // binning: range:  -40°- 40
948       fProfResXvsPhiclu[nhist][clu]->SetName(aux.Data());                             //          bin width: 2°
949       fProfResXvsPhiclu[nhist][clu]->SetTitle(aux.Data());
950
951       aux=profnameResZvsDipclu;
952       aux+=nhist;
953       aux+="_clu_";
954       aux+=clu+1; // clu=0 --> cluster size 1
955       fProfResZvsDipclu[nhist][clu]= new TProfile("histname","histname",48,-72.,72.0); // binning: range:  -70°-7°
956       fProfResZvsDipclu[nhist][clu]->SetName(aux.Data());                              //      bin width: 3°
957       fProfResZvsDipclu[nhist][clu]->SetTitle(aux.Data());
958     }
959
960     fHisTrackXFOtrue[nhist]=new TH1F*[kNClockPhase];
961     fHisTrackZFOtrue[nhist]=new TH1F*[kNClockPhase];
962     fHisTrackXFOfalse[nhist]=new TH1F*[kNClockPhase];
963     fHisTrackZFOfalse[nhist]=new TH1F*[kNClockPhase];
964     fHisTrackXZFOtrue[nhist]=new TH2F*[kNClockPhase];
965     fHisTrackXZFOfalse[nhist]=new TH2F*[kNClockPhase];
966     for(Int_t phas=0; phas<kNClockPhase;phas++){
967       aux=histnameTrackXFOtrue;
968       aux+=nhist;
969       aux+="_BCmod4_";
970       aux+=phas;
971       fHisTrackXFOtrue[nhist][phas]=new TH1F("histname","histname",128,-0.64,0.64); // +- 6.4 mm; 1 bin=0.1 mm
972       fHisTrackXFOtrue[nhist][phas]->SetName(aux.Data());
973       fHisTrackXFOtrue[nhist][phas]->SetTitle(aux.Data());
974       
975       aux=histnameTrackZFOtrue;
976       aux+=nhist;
977       aux+="_BCmod4_";
978       aux+=phas;
979       fHisTrackZFOtrue[nhist][phas]=new TH1F("histname","histname",350,-3.5,3.5); // +- 35. mm; 1 bin=0.2 mm
980       fHisTrackZFOtrue[nhist][phas]->SetName(aux.Data());
981       fHisTrackZFOtrue[nhist][phas]->SetTitle(aux.Data());
982       
983       aux=histnameTrackXFOfalse;
984       aux+=nhist;
985       aux+="_BCmod4_";
986       aux+=phas;
987       fHisTrackXFOfalse[nhist][phas]=new TH1F("histname","histname",128,-0.64,0.64); // +- 6.4 mm; 1 bin=0.1 mm
988       fHisTrackXFOfalse[nhist][phas]->SetName(aux.Data());
989       fHisTrackXFOfalse[nhist][phas]->SetTitle(aux.Data());
990
991       aux=histnameTrackZFOfalse;
992       aux+=nhist;
993       aux+="_BCmod4_";
994       aux+=phas;
995       fHisTrackZFOfalse[nhist][phas]=new TH1F("histname","histname",350,-3.5,3.5); // +- 35. mm; 1 bin=0.2 mm
996       fHisTrackZFOfalse[nhist][phas]->SetName(aux.Data());
997       fHisTrackZFOfalse[nhist][phas]->SetTitle(aux.Data());
998     
999       aux=histnameTrackXZFOtrue;
1000       aux+=nhist;
1001       aux+="_BCmod4_";
1002       aux+=phas;
1003       fHisTrackXZFOtrue[nhist][phas]=new TH2F("histname","histname",22,-3.5,3.5,32,-0.64,0.64); //  localZ +- 35. mm; 1 bin=3.2 mm
1004       fHisTrackXZFOtrue[nhist][phas]->SetName(aux.Data());                                      //  localX +- 6.4 mm; 1 bin=0.4 mm
1005       fHisTrackXZFOtrue[nhist][phas]->SetTitle(aux.Data());
1006
1007       aux=histnameTrackXZFOfalse;
1008       aux+=nhist;
1009       aux+="_BCmod4_";
1010       aux+=phas;
1011       fHisTrackXZFOfalse[nhist][phas]=new TH2F("histname","histname",22,-3.5,3.5,32,-0.64,0.64); //  localZ +- 35. mm; 1 bin=3.2 mm
1012       fHisTrackXZFOfalse[nhist][phas]->SetName(aux.Data());                                      //  localX +- 6.4 mm; 1 bin=0.4 mm
1013       fHisTrackXZFOfalse[nhist][phas]->SetTitle(aux.Data());
1014       } 
1015   } // end loop on module
1016
1017   TH1::AddDirectory(kTRUE);
1018
1019 return;
1020 }
1021 //__________________________________________________________
1022 void AliITSPlaneEffSPD::DeleteHistos() {
1023   if(fHisResX) {
1024     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResX[i];
1025     delete [] fHisResX; fHisResX=0;
1026   }
1027   if(fHisResZ) {
1028     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResZ[i];
1029     delete [] fHisResZ; fHisResZ=0;
1030   }
1031   if(fHisResXZ) {
1032     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResXZ[i];
1033     delete [] fHisResXZ; fHisResXZ=0;
1034   }
1035   if(fHisClusterSize) {
1036     for (Int_t i=0; i<kNHisto; i++ ) delete fHisClusterSize[i];
1037     delete [] fHisClusterSize; fHisClusterSize=0;
1038   }
1039   if(fHisResXclu) {
1040     for (Int_t i=0; i<kNHisto; i++ ) {
1041       for (Int_t clu=0; clu<kNclu; clu++) if (fHisResXclu[i][clu]) delete fHisResXclu[i][clu];
1042       delete [] fHisResXclu[i];
1043     }
1044     delete [] fHisResXclu;
1045     fHisResXclu = 0;
1046   }
1047   if(fHisResZclu) {
1048     for (Int_t i=0; i<kNHisto; i++ ) {
1049       for (Int_t clu=0; clu<kNclu; clu++) if (fHisResZclu[i][clu]) delete fHisResZclu[i][clu];
1050       delete [] fHisResZclu[i];
1051     }
1052     delete [] fHisResZclu;
1053     fHisResZclu = 0;
1054   }
1055   if(fHisResXchip) {
1056     for (Int_t i=0; i<kNHisto; i++ ) {
1057       for (Int_t chip=0; chip<kNChip; chip++) if (fHisResXchip[i][chip]) delete fHisResXchip[i][chip];
1058       delete [] fHisResXchip[i];
1059     }
1060     delete [] fHisResXchip;
1061     fHisResXchip = 0;
1062   }
1063   if(fHisResZchip) {
1064     for (Int_t i=0; i<kNHisto; i++ ) {
1065       for (Int_t chip=0; chip<kNChip; chip++) if (fHisResZchip[i][chip]) delete fHisResZchip[i][chip];
1066       delete [] fHisResZchip[i];
1067     }
1068     delete [] fHisResZchip;
1069     fHisResZchip = 0;
1070   }
1071   if(fHisTrackErrX) {
1072     for (Int_t i=0; i<kNHisto; i++ ) delete fHisTrackErrX[i];
1073     delete [] fHisTrackErrX; fHisTrackErrX=0;
1074   }
1075   if(fHisTrackErrZ) {
1076     for (Int_t i=0; i<kNHisto; i++ ) delete fHisTrackErrZ[i];
1077     delete [] fHisTrackErrZ; fHisTrackErrZ=0;
1078   }
1079   if(fHisClusErrX) {
1080     for (Int_t i=0; i<kNHisto; i++ ) delete fHisClusErrX[i];
1081     delete [] fHisClusErrX; fHisClusErrX=0;
1082   }
1083   if(fHisClusErrZ) {
1084     for (Int_t i=0; i<kNHisto; i++ ) delete fHisClusErrZ[i];
1085     delete [] fHisClusErrZ; fHisClusErrZ=0;
1086   }
1087   if(fProfResXvsPhi) {
1088     for (Int_t i=0; i<kNHisto; i++ ) delete fProfResXvsPhi[i];
1089     delete [] fProfResXvsPhi; fProfResXvsPhi=0;
1090   }
1091   if(fProfResZvsDip) {
1092     for (Int_t i=0; i<kNHisto; i++ ) delete fProfResZvsDip[i];
1093     delete [] fProfResZvsDip; fProfResZvsDip=0;
1094   }
1095   if(fProfResXvsPhiclu) {
1096     for (Int_t i=0; i<kNHisto; i++ ) {
1097       for (Int_t clu=0; clu<kNclu; clu++) if (fProfResXvsPhiclu[i][clu]) delete fProfResXvsPhiclu[i][clu];
1098       delete [] fProfResXvsPhiclu[i];
1099     }
1100     delete [] fProfResXvsPhiclu;
1101     fProfResXvsPhiclu = 0;
1102   }
1103   if(fProfResZvsDipclu) {
1104     for (Int_t i=0; i<kNHisto; i++ ) {
1105       for (Int_t clu=0; clu<kNclu; clu++) if (fProfResZvsDipclu[i][clu]) delete fProfResZvsDipclu[i][clu];
1106       delete [] fProfResZvsDipclu[i];
1107     }
1108     delete [] fProfResZvsDipclu;
1109     fProfResZvsDipclu = 0;
1110   }
1111   if(fHisTrackXFOtrue) {
1112     for (Int_t i=0; i<kNHisto; i++ ) {
1113       for (Int_t phas=0; phas<kNClockPhase; phas++) if (fHisTrackXFOtrue[i][phas]) delete fHisTrackXFOtrue[i][phas];
1114       delete [] fHisTrackXFOtrue[i];
1115     }
1116     delete [] fHisTrackXFOtrue;
1117     fHisTrackXFOtrue = 0;
1118   }
1119   if(fHisTrackZFOtrue) {
1120     for (Int_t i=0; i<kNHisto; i++ ) {
1121       for (Int_t phas=0; phas<kNClockPhase; phas++) if (fHisTrackZFOtrue[i][phas]) delete fHisTrackZFOtrue[i][phas];
1122       delete [] fHisTrackZFOtrue[i];
1123     }
1124     delete [] fHisTrackZFOtrue;
1125     fHisTrackZFOtrue = 0;
1126   }
1127   if(fHisTrackXFOfalse) {
1128     for (Int_t i=0; i<kNHisto; i++ ) {
1129       for (Int_t phas=0; phas<kNClockPhase; phas++) if (fHisTrackXFOfalse[i][phas]) delete fHisTrackXFOfalse[i][phas];
1130       delete [] fHisTrackXFOfalse[i];
1131     }
1132     delete [] fHisTrackXFOfalse;
1133     fHisTrackXFOfalse = 0;
1134   }
1135   if(fHisTrackZFOfalse) {
1136     for (Int_t i=0; i<kNHisto; i++ ) {
1137       for (Int_t phas=0; phas<kNClockPhase; phas++) if (fHisTrackZFOfalse[i][phas]) delete fHisTrackZFOfalse[i][phas];
1138       delete [] fHisTrackZFOfalse[i];
1139     }
1140     delete [] fHisTrackZFOfalse;
1141     fHisTrackZFOfalse = 0;
1142   }
1143 return;
1144 }
1145 //__________________________________________________________
1146 Bool_t AliITSPlaneEffSPD::FillHistos(UInt_t key, Bool_t found,
1147                                      Float_t *tr, Float_t *clu, Int_t *csize, Float_t *angtrkmod) {
1148 //
1149 // depending on the value of key this method
1150 // either call the standard one for clusters 
1151 // or the one for FO studies
1152 // if key <  1200 --> call FillHistosST
1153 // if key >= 1200 --> call FillHistosFO
1154 if(key>=kNModule*kNChip*(kNClockPhase+1))
1155   {AliError("GetChipFromKey: you asked for a non existing key"); return kFALSE;}
1156 if(key<kNModule*kNChip) return FillHistosStd(key,found,tr,clu,csize,angtrkmod);
1157 else return FillHistosFO(key,found,tr);
1158 return kFALSE;
1159 }
1160 //__________________________________________________________
1161 Bool_t AliITSPlaneEffSPD::FillHistosFO(UInt_t key, Bool_t found, Float_t *tr) {
1162 // this method fill the histograms for FastOr studies
1163 // input: - key: unique key of the basic block 
1164 //        - found: Boolean to asses whether a FastOr bit has been associated to the track or not 
1165 //        - tr[0],tr[1] local X and Z coordinates of the track prediction, respectively
1166 //        - tr[2],tr[3] error on local X and Z coordinates of the track prediction, respectively
1167 // output: kTRUE if filling was succesfull kFALSE otherwise
1168 // side effects: updating of the histograms.
1169   if (!fHis) {
1170     AliWarning("FillHistos: histograms do not exist! Call SetCreateHistos(kTRUE) first");
1171     return kFALSE;
1172   }
1173   if(key>=kNModule*kNChip*(kNClockPhase+1))
1174     {AliWarning("FillHistos: you asked for a non existing key"); return kFALSE;}
1175   if(key<kNModule*kNChip)
1176     {AliWarning("FillHistos: you asked for a key which is not for FO studies"); return kFALSE;}
1177   Int_t id=GetModFromKey(key);
1178   Int_t BCm4=GetBCm4FromKey(key);
1179   if(id>=kNHisto)
1180     {AliWarning("FillHistos: you want to fill a non-existing histos"); return kFALSE;}
1181   if(found) {
1182     fHisTrackXFOtrue[id][BCm4]->Fill(tr[0]);
1183     fHisTrackZFOtrue[id][BCm4]->Fill(tr[1]);
1184     fHisTrackXZFOtrue[id][BCm4]->Fill(tr[1],tr[0]);
1185   }
1186   else {
1187     fHisTrackXFOfalse[id][BCm4]->Fill(tr[0]);
1188     fHisTrackZFOfalse[id][BCm4]->Fill(tr[1]);
1189     fHisTrackXZFOfalse[id][BCm4]->Fill(tr[1],tr[0]);
1190   }
1191 return kTRUE;
1192 }
1193 //__________________________________________________________
1194 Bool_t AliITSPlaneEffSPD::FillHistosStd(UInt_t key, Bool_t found, 
1195                                      Float_t *tr, Float_t *clu, Int_t *csize, Float_t *angtrkmod) {
1196 // this method fill the histograms
1197 // input: - key: unique key of the basic block 
1198 //        - found: Boolean to asses whether a cluster has been associated to the track or not 
1199 //        - tr[0],tr[1] local X and Z coordinates of the track prediction, respectively
1200 //        - tr[2],tr[3] error on local X and Z coordinates of the track prediction, respectively
1201 //        - clu[0],clu[1] local X and Z coordinates of the cluster associated to the track, respectively
1202 //        - clu[2],clu[3] error on local X and Z coordinates of the cluster associated to the track, respectively
1203 //        - csize[0][1] cluster size in X and Z, respectively
1204 //        - angtrkmod[0],angtrkmod[1]  
1205 // output: kTRUE if filling was succesfull kFALSE otherwise
1206 // side effects: updating of the histograms. 
1207 //
1208   if (!fHis) {
1209     AliWarning("FillHistos: histograms do not exist! Call SetCreateHistos(kTRUE) first");
1210     return kFALSE;
1211   }
1212   if(key>=kNModule*kNChip)
1213     {AliWarning("FillHistos: you asked for a non existing key"); return kFALSE;}
1214   Int_t id=GetModFromKey(key);
1215   Int_t chip=GetChipFromKey(key);
1216   if(id>=kNHisto) 
1217     {AliWarning("FillHistos: you want to fill a non-existing histos"); return kFALSE;}
1218   if(found) {
1219     Float_t resx=tr[0]-clu[0];
1220     Float_t resz=tr[1]-clu[1];
1221     fHisResX[id]->Fill(resx);
1222     fHisResZ[id]->Fill(resz);
1223     fHisResXZ[id]->Fill(resx,resz);
1224     fHisClusterSize[id]->Fill((Double_t)csize[0],(Double_t)csize[1]);
1225     if(csize[0]>0 &&  csize[0]<=kNclu) fHisResXclu[id][csize[0]-1]->Fill(resx);
1226     if(csize[1]>0 &&  csize[1]<=kNclu) fHisResZclu[id][csize[1]-1]->Fill(resz);
1227     fHisResXchip[id][chip]->Fill(resx);
1228     fHisResZchip[id][chip]->Fill(resz);
1229     fProfResXvsPhi[id]->Fill(angtrkmod[0],resx);
1230     fProfResZvsDip[id]->Fill(angtrkmod[1],resz);
1231     if(csize[0]>0 &&  csize[0]<=kNclu) fProfResXvsPhiclu[id][csize[0]-1]->Fill(angtrkmod[0],resx);
1232     if(csize[1]>0 &&  csize[1]<=kNclu) fProfResZvsDipclu[id][csize[1]-1]->Fill(angtrkmod[1],resz);
1233   }
1234   fHisTrackErrX[id]->Fill(tr[2]);
1235   fHisTrackErrZ[id]->Fill(tr[3]);
1236   fHisClusErrX[id]->Fill(clu[2]);
1237   fHisClusErrZ[id]->Fill(clu[3]);
1238   return kTRUE;
1239 }
1240 //__________________________________________________________
1241 Bool_t AliITSPlaneEffSPD::WriteHistosToFile(TString filename, Option_t* option) {
1242   //
1243   // Saves the histograms into a tree and saves the trees into a file
1244   //
1245   if (!fHis) return kFALSE;
1246   if (filename.IsNull() || filename.IsWhitespace()) {
1247      AliWarning("WriteHistosToFile: null output filename!");
1248      return kFALSE;
1249   }
1250   char branchname[51];
1251   TFile *hFile=new TFile(filename.Data(),option,
1252                          "The File containing the TREEs with ITS PlaneEff Histos");
1253   TTree *SPDTree=new TTree("SPDTree","Tree whith Residuals and Cluster Type distributions for SPD");
1254   TH1F *histZ,*histX;
1255   TH2F *histXZ;
1256   TH2I *histClusterType;
1257   TH1F *histXclu[kNclu];
1258   TH1F *histZclu[kNclu];
1259   TH1F *histXchip[kNChip];
1260   TH1F *histZchip[kNChip];
1261   TH1F *histTrErrZ,*histTrErrX;
1262   TH1F *histClErrZ,*histClErrX;
1263   TProfile *profXvsPhi,*profZvsDip;
1264   TProfile *profXvsPhiclu[kNclu],*profZvsDipclu[kNclu];
1265   TH1F *histXtrkFOtrue[kNClockPhase];
1266   TH1F *histZtrkFOtrue[kNClockPhase];
1267   TH1F *histXtrkFOfalse[kNClockPhase];
1268   TH1F *histZtrkFOfalse[kNClockPhase];
1269   TH2F *histXZtrkFOtrue[kNClockPhase];
1270   TH2F *histXZtrkFOfalse[kNClockPhase];
1271
1272   histZ=new TH1F();
1273   histX=new TH1F();
1274   histXZ=new TH2F();
1275   histClusterType=new TH2I();
1276   for(Int_t clu=0;clu<kNclu;clu++) {
1277     histXclu[clu]=new TH1F();
1278     histZclu[clu]=new TH1F();
1279   }
1280   for(Int_t chip=0;chip<kNChip;chip++) {
1281     histXchip[chip]=new TH1F();
1282     histZchip[chip]=new TH1F();
1283   }
1284
1285   histTrErrX=new TH1F();
1286   histTrErrZ=new TH1F();
1287   histClErrX=new TH1F();
1288   histClErrZ=new TH1F();
1289   profXvsPhi=new TProfile();
1290   profZvsDip=new TProfile();
1291   for(Int_t clu=0;clu<kNclu;clu++) {
1292     profXvsPhiclu[clu]=new TProfile();
1293     profZvsDipclu[clu]=new TProfile();
1294   }
1295
1296   for(Int_t phas=0; phas<kNClockPhase;phas++){
1297     histXtrkFOtrue[phas]=new TH1F();
1298     histZtrkFOtrue[phas]=new TH1F();
1299     histXtrkFOfalse[phas]=new TH1F();
1300     histZtrkFOfalse[phas]=new TH1F();
1301     histXZtrkFOtrue[phas]=new TH2F();
1302     histXZtrkFOfalse[phas]=new TH2F();
1303   }
1304
1305   SPDTree->Branch("histX","TH1F",&histX,128000,0);
1306   SPDTree->Branch("histZ","TH1F",&histZ,128000,0);
1307   SPDTree->Branch("histXZ","TH2F",&histXZ,128000,0);
1308   SPDTree->Branch("histClusterType","TH2I",&histClusterType,128000,0);
1309   for(Int_t clu=0;clu<kNclu;clu++) {
1310     snprintf(branchname,50,"histXclu_%d",clu+1);
1311     SPDTree->Branch(branchname,"TH1F",&histXclu[clu],128000,0);
1312     snprintf(branchname,50,"histZclu_%d",clu+1);
1313     SPDTree->Branch(branchname,"TH1F",&histZclu[clu],128000,0);
1314   }
1315   for(Int_t chip=0;chip<kNChip;chip++) {
1316     snprintf(branchname,50,"histXchip_%d",chip);
1317     SPDTree->Branch(branchname,"TH1F",&histXchip[chip],128000,0);
1318     snprintf(branchname,50,"histZchip_%d",chip);
1319     SPDTree->Branch(branchname,"TH1F",&histZchip[chip],128000,0);
1320   }
1321   SPDTree->Branch("histTrErrX","TH1F",&histTrErrX,128000,0);
1322   SPDTree->Branch("histTrErrZ","TH1F",&histTrErrZ,128000,0);
1323   SPDTree->Branch("histClErrX","TH1F",&histClErrX,128000,0);
1324   SPDTree->Branch("histClErrZ","TH1F",&histClErrZ,128000,0);
1325   SPDTree->Branch("profXvsPhi","TProfile",&profXvsPhi,128000,0);
1326   SPDTree->Branch("profZvsDip","TProfile",&profZvsDip,128000,0);
1327   for(Int_t clu=0;clu<kNclu;clu++) {
1328     snprintf(branchname,50,"profXvsPhiclu_%d",clu+1);
1329     SPDTree->Branch(branchname,"TProfile",&profXvsPhiclu[clu],128000,0);
1330     snprintf(branchname,50,"profZvsDipclu_%d",clu+1);
1331     SPDTree->Branch(branchname,"TProfile",&profZvsDipclu[clu],128000,0);
1332   }
1333   for(Int_t phas=0; phas<kNClockPhase;phas++){
1334     snprintf(branchname,50,"histTrXFOokBCmod4_%d",phas);
1335     SPDTree->Branch(branchname,"TH1F",&histXtrkFOtrue[phas],128000,0);
1336     snprintf(branchname,50,"histTrZFOokBCmod4_%d",phas);
1337     SPDTree->Branch(branchname,"TH1F",&histZtrkFOtrue[phas],128000,0);
1338     snprintf(branchname,50,"histTrXFOkoBCmod4_%d",phas);
1339     SPDTree->Branch(branchname,"TH1F",&histXtrkFOfalse[phas],128000,0);
1340     snprintf(branchname,50,"histTrZFOkoBCmod4_%d",phas);
1341     SPDTree->Branch(branchname,"TH1F",&histZtrkFOfalse[phas],128000,0);
1342     snprintf(branchname,50,"histTrXZFOokBCmod4_%d",phas);
1343     SPDTree->Branch(branchname,"TH2F",&histXZtrkFOtrue[phas],128000,0);
1344     snprintf(branchname,50,"histTrXZFOkoBCmod4_%d",phas);
1345     SPDTree->Branch(branchname,"TH2F",&histXZtrkFOfalse[phas],128000,0);
1346   }
1347
1348   for(Int_t j=0;j<kNHisto;j++){
1349     histX=fHisResX[j];
1350     histZ=fHisResZ[j];
1351     histXZ=fHisResXZ[j];
1352     histClusterType=fHisClusterSize[j];
1353     for(Int_t clu=0;clu<kNclu;clu++) {
1354       histXclu[clu]=fHisResXclu[j][clu];
1355       histZclu[clu]=fHisResZclu[j][clu];
1356     }
1357     for(Int_t chip=0;chip<kNChip;chip++) {
1358       histXchip[chip]=fHisResXchip[j][chip];
1359       histZchip[chip]=fHisResZchip[j][chip];
1360     }
1361     histTrErrX=fHisTrackErrX[j];
1362     histTrErrZ=fHisTrackErrZ[j];
1363     histClErrX=fHisClusErrX[j];
1364     histClErrZ=fHisClusErrZ[j];
1365     profXvsPhi=fProfResXvsPhi[j];
1366     profZvsDip=fProfResZvsDip[j];
1367     for(Int_t clu=0;clu<kNclu;clu++) {
1368       profXvsPhiclu[clu]=fProfResXvsPhiclu[j][clu];
1369       profZvsDipclu[clu]=fProfResZvsDipclu[j][clu];
1370     }
1371     for(Int_t phas=0; phas<kNClockPhase;phas++){
1372       histXtrkFOtrue[phas]=fHisTrackXFOtrue[j][phas];
1373       histZtrkFOtrue[phas]=fHisTrackZFOtrue[j][phas];
1374       histXtrkFOfalse[phas]=fHisTrackXFOfalse[j][phas];
1375       histZtrkFOfalse[phas]=fHisTrackZFOfalse[j][phas];
1376       histXZtrkFOtrue[phas]=fHisTrackXZFOtrue[j][phas];
1377       histXZtrkFOfalse[phas]=fHisTrackXZFOfalse[j][phas];
1378     }
1379
1380     SPDTree->Fill();
1381   }
1382   hFile->Write();
1383   hFile->Close();
1384 return kTRUE;
1385 }
1386 //__________________________________________________________
1387 Bool_t AliITSPlaneEffSPD::ReadHistosFromFile(TString filename) {
1388   //
1389   // Read histograms from an already existing file 
1390   //
1391   if (!fHis) return kFALSE;
1392   if (filename.IsNull() || filename.IsWhitespace()) {
1393      AliWarning("ReadHistosFromFile: incorrect output filename!");
1394      return kFALSE;
1395   }
1396   char branchname[51];
1397
1398   TH1F *h  = 0;
1399   TH2F *h2 = 0;
1400   TH2I *h2i= 0;
1401   TProfile *p = 0;
1402
1403   TFile *file=TFile::Open(filename.Data(),"READONLY");
1404
1405   if (!file || file->IsZombie()) {
1406     AliWarning(Form("Can't open %s !",filename.Data()));
1407     delete file;
1408     return kFALSE;
1409   }
1410   TTree *tree = (TTree*) file->Get("SPDTree");
1411
1412   TBranch *histX = (TBranch*) tree->GetBranch("histX");
1413   TBranch *histZ = (TBranch*) tree->GetBranch("histZ");
1414   TBranch *histXZ = (TBranch*) tree->GetBranch("histXZ");
1415   TBranch *histClusterType = (TBranch*) tree->GetBranch("histClusterType");
1416    
1417   TBranch *histXclu[kNclu], *histZclu[kNclu];
1418   for(Int_t clu=0; clu<kNclu; clu++) {
1419     snprintf(branchname,50,"histXclu_%d",clu+1);
1420     histXclu[clu]= (TBranch*) tree->GetBranch(branchname);
1421     snprintf(branchname,50,"histZclu_%d",clu+1);
1422     histZclu[clu]= (TBranch*) tree->GetBranch(branchname);
1423   }
1424
1425   TBranch *histXchip[kNChip], *histZchip[kNChip];
1426   for(Int_t chip=0; chip<kNChip; chip++) {
1427     snprintf(branchname,50,"histXchip_%d",chip);
1428     histXchip[chip]= (TBranch*) tree->GetBranch(branchname);
1429     snprintf(branchname,50,"histZchip_%d",chip);
1430     histZchip[chip]= (TBranch*) tree->GetBranch(branchname);
1431   }
1432
1433   TBranch *histTrErrX = (TBranch*) tree->GetBranch("histTrErrX");
1434   TBranch *histTrErrZ = (TBranch*) tree->GetBranch("histTrErrZ");
1435   TBranch *histClErrX = (TBranch*) tree->GetBranch("histClErrX");
1436   TBranch *histClErrZ = (TBranch*) tree->GetBranch("histClErrZ");
1437   TBranch *profXvsPhi = (TBranch*) tree->GetBranch("profXvsPhi");
1438   TBranch *profZvsDip = (TBranch*) tree->GetBranch("profZvsDip");
1439
1440   TBranch *profXvsPhiclu[kNclu], *profZvsDipclu[kNclu];
1441   for(Int_t clu=0; clu<kNclu; clu++) {
1442     snprintf(branchname,50,"profXvsPhiclu_%d",clu+1);
1443     profXvsPhiclu[clu]= (TBranch*) tree->GetBranch(branchname);
1444     snprintf(branchname,50,"profZvsDipclu_%d",clu+1);
1445     profZvsDipclu[clu]= (TBranch*) tree->GetBranch(branchname);
1446   }
1447
1448   TBranch *histXtrkFOtrue[kNClockPhase], *histZtrkFOtrue[kNClockPhase],
1449           *histXtrkFOfalse[kNClockPhase], *histZtrkFOfalse[kNClockPhase],
1450           *histXZtrkFOtrue[kNClockPhase], *histXZtrkFOfalse[kNClockPhase];
1451   for(Int_t phas=0; phas<kNClockPhase;phas++){
1452     snprintf(branchname,50,"histTrXFOokBCmod4_%d",phas);
1453     histXtrkFOtrue[phas] = (TBranch*) tree->GetBranch(branchname);
1454     snprintf(branchname,50,"histTrZFOokBCmod4_%d",phas);
1455     histZtrkFOtrue[phas] = (TBranch*) tree->GetBranch(branchname);
1456     snprintf(branchname,50,"histTrXFOkoBCmod4_%d",phas);
1457     histXtrkFOfalse[phas] = (TBranch*) tree->GetBranch(branchname);
1458     snprintf(branchname,50,"histTrZFOkoBCmod4_%d",phas);
1459     histZtrkFOfalse[phas] = (TBranch*) tree->GetBranch(branchname);
1460     snprintf(branchname,50,"histTrXZFOokBCmod4_%d",phas);
1461     histXZtrkFOtrue[phas] = (TBranch*) tree->GetBranch(branchname);
1462     snprintf(branchname,50,"histTrXZFOkoBCmod4_%d",phas);
1463     histXZtrkFOfalse[phas] = (TBranch*) tree->GetBranch(branchname);
1464   }
1465
1466   gROOT->cd();
1467
1468   Int_t nevent = (Int_t)histX->GetEntries();
1469   if(nevent!=kNHisto) 
1470     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1471   histX->SetAddress(&h);
1472   for(Int_t j=0;j<kNHisto;j++){
1473     histX->GetEntry(j);
1474     fHisResX[j]->Add(h);
1475   }
1476
1477   nevent = (Int_t)histZ->GetEntries();
1478   if(nevent!=kNHisto) 
1479     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1480   histZ->SetAddress(&h);
1481   for(Int_t j=0;j<kNHisto;j++){
1482     histZ->GetEntry(j);
1483     fHisResZ[j]->Add(h);
1484   }
1485
1486   nevent = (Int_t)histXZ->GetEntries();
1487   if(nevent!=kNHisto) 
1488     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1489   histXZ->SetAddress(&h2);
1490   for(Int_t j=0;j<kNHisto;j++){
1491     histXZ->GetEntry(j);
1492     fHisResXZ[j]->Add(h2);
1493   }
1494
1495   nevent = (Int_t)histClusterType->GetEntries();
1496   if(nevent!=kNHisto) 
1497     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1498   histClusterType->SetAddress(&h2i);
1499   for(Int_t j=0;j<kNHisto;j++){
1500     histClusterType->GetEntry(j);
1501     fHisClusterSize[j]->Add(h2i);
1502   }
1503
1504   for(Int_t clu=0; clu<kNclu; clu++) {
1505
1506     nevent = (Int_t)histXclu[clu]->GetEntries();
1507     if(nevent!=kNHisto)
1508       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1509     histXclu[clu]->SetAddress(&h);
1510     for(Int_t j=0;j<kNHisto;j++){
1511       histXclu[clu]->GetEntry(j);
1512       fHisResXclu[j][clu]->Add(h);
1513     }
1514
1515    nevent = (Int_t)histZclu[clu]->GetEntries();
1516     if(nevent!=kNHisto)
1517       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1518     histZclu[clu]->SetAddress(&h);
1519     for(Int_t j=0;j<kNHisto;j++){
1520       histZclu[clu]->GetEntry(j);
1521       fHisResZclu[j][clu]->Add(h);
1522     }
1523   }
1524
1525
1526     for(Int_t chip=0; chip<kNChip; chip++) {
1527
1528     nevent = (Int_t)histXchip[chip]->GetEntries();
1529     if(nevent!=kNHisto)
1530       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1531     histXchip[chip]->SetAddress(&h);
1532     for(Int_t j=0;j<kNHisto;j++){
1533       histXchip[chip]->GetEntry(j);
1534       fHisResXchip[j][chip]->Add(h);
1535     }
1536
1537     nevent = (Int_t)histZchip[chip]->GetEntries();
1538     if(nevent!=kNHisto)
1539       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1540     histZchip[chip]->SetAddress(&h);
1541     for(Int_t j=0;j<kNHisto;j++){
1542       histZchip[chip]->GetEntry(j);
1543       fHisResZchip[j][chip]->Add(h);
1544     }
1545   }
1546
1547   nevent = (Int_t)histTrErrX->GetEntries(); 
1548   if(nevent!=kNHisto)
1549     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1550   histTrErrX->SetAddress(&h);
1551   for(Int_t j=0;j<kNHisto;j++){
1552     histTrErrX->GetEntry(j);
1553     fHisTrackErrX[j]->Add(h);
1554   }
1555
1556   nevent = (Int_t)histTrErrZ->GetEntries();
1557   if(nevent!=kNHisto)
1558     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1559   histTrErrZ->SetAddress(&h);
1560   for(Int_t j=0;j<kNHisto;j++){
1561     histTrErrZ->GetEntry(j);
1562     fHisTrackErrZ[j]->Add(h);
1563   }
1564
1565   nevent = (Int_t)histClErrX->GetEntries();
1566   if(nevent!=kNHisto)
1567     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1568   histClErrX->SetAddress(&h);
1569   for(Int_t j=0;j<kNHisto;j++){
1570     histClErrX->GetEntry(j);
1571     fHisClusErrX[j]->Add(h);
1572   }
1573
1574   nevent = (Int_t)histClErrZ->GetEntries();
1575   if(nevent!=kNHisto)
1576     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1577   histClErrZ->SetAddress(&h);
1578   for(Int_t j=0;j<kNHisto;j++){
1579     histClErrZ->GetEntry(j);
1580     fHisClusErrZ[j]->Add(h);
1581   }
1582   nevent = (Int_t)profXvsPhi->GetEntries();
1583   if(nevent!=kNHisto)
1584     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1585   profXvsPhi->SetAddress(&p);
1586   for(Int_t j=0;j<kNHisto;j++){
1587     profXvsPhi->GetEntry(j);
1588     fProfResXvsPhi[j]->Add(p);
1589   }
1590
1591   nevent = (Int_t)profZvsDip->GetEntries();
1592   if(nevent!=kNHisto)
1593     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1594   profZvsDip->SetAddress(&p);
1595   for(Int_t j=0;j<kNHisto;j++){ 
1596     profZvsDip->GetEntry(j);
1597     fProfResZvsDip[j]->Add(p);
1598   }
1599
1600     for(Int_t clu=0; clu<kNclu; clu++) {
1601
1602     nevent = (Int_t)profXvsPhiclu[clu]->GetEntries();
1603     if(nevent!=kNHisto)
1604       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1605     profXvsPhiclu[clu]->SetAddress(&p);
1606     for(Int_t j=0;j<kNHisto;j++){
1607       profXvsPhiclu[clu]->GetEntry(j);
1608       fProfResXvsPhiclu[j][clu]->Add(p);
1609     }
1610
1611     nevent = (Int_t)profZvsDipclu[clu]->GetEntries();
1612     if(nevent!=kNHisto)
1613       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1614     profZvsDipclu[clu]->SetAddress(&p);
1615     for(Int_t j=0;j<kNHisto;j++){ 
1616       profZvsDipclu[clu]->GetEntry(j);
1617       fProfResZvsDipclu[j][clu]->Add(p);
1618     }
1619   }
1620
1621     for(Int_t phas=0; phas<kNClockPhase;phas++){
1622
1623     nevent = (Int_t)histXtrkFOtrue[phas]->GetEntries();
1624     if(nevent!=kNHisto)
1625        {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1626     histXtrkFOtrue[phas]->SetAddress(&h);
1627     for(Int_t j=0;j<kNHisto;j++){
1628       histXtrkFOtrue[phas]->GetEntry(j);
1629       fHisTrackXFOtrue[j][phas]->Add(h);
1630     }
1631
1632     nevent = (Int_t)histZtrkFOtrue[phas]->GetEntries();
1633     if(nevent!=kNHisto)
1634        {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1635     histZtrkFOtrue[phas]->SetAddress(&h);
1636     for(Int_t j=0;j<kNHisto;j++){
1637       histZtrkFOtrue[phas]->GetEntry(j);
1638       fHisTrackZFOtrue[j][phas]->Add(h);
1639     }
1640
1641     nevent = (Int_t)histXtrkFOfalse[phas]->GetEntries();
1642     if(nevent!=kNHisto)
1643        {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1644     histXtrkFOfalse[phas]->SetAddress(&h);
1645     for(Int_t j=0;j<kNHisto;j++){
1646       histXtrkFOfalse[phas]->GetEntry(j);
1647       fHisTrackXFOfalse[j][phas]->Add(h);
1648     }
1649
1650     nevent = (Int_t)histZtrkFOfalse[phas]->GetEntries();
1651     if(nevent!=kNHisto)
1652        {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1653     histZtrkFOfalse[phas]->SetAddress(&h);
1654     for(Int_t j=0;j<kNHisto;j++){
1655       histZtrkFOfalse[phas]->GetEntry(j);
1656       fHisTrackZFOfalse[j][phas]->Add(h);
1657     }
1658
1659     nevent = (Int_t)histXZtrkFOtrue[phas]->GetEntries();
1660     if(nevent!=kNHisto)
1661        {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1662     histXZtrkFOtrue[phas]->SetAddress(&h2);
1663     for(Int_t j=0;j<kNHisto;j++){
1664       histXZtrkFOtrue[phas]->GetEntry(j);
1665       fHisTrackXZFOtrue[j][phas]->Add(h2);
1666     }
1667
1668     nevent = (Int_t)histXZtrkFOfalse[phas]->GetEntries();
1669     if(nevent!=kNHisto)
1670        {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1671     histXZtrkFOfalse[phas]->SetAddress(&h2);
1672     for(Int_t j=0;j<kNHisto;j++){
1673       histXZtrkFOfalse[phas]->GetEntry(j);
1674       fHisTrackXZFOfalse[j][phas]->Add(h2);
1675     }
1676
1677    }
1678
1679   delete h;   
1680   delete h2;  
1681   delete h2i; 
1682   delete p;   
1683
1684   if (file) {
1685     file->Close();
1686     delete file;
1687   }
1688 return kTRUE;
1689 }
1690