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