]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSPlaneEffSPD.cxx
- New histograms added for residual distributions vs. incidence angle
[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   for (UInt_t i=0; i<kNModule*kNChip; i++){
62     fFound[i]=0;
63     fTried[i]=0;
64   }
65   // default constructor
66   AliDebug(1,Form("Calling default constructor"));
67 }
68 //______________________________________________________________________
69 AliITSPlaneEffSPD::~AliITSPlaneEffSPD(){
70     // destructor
71     // Inputs:
72     //    none.
73     // Outputs:
74     //    none.
75     // Return:
76     //     none.
77     DeleteHistos();
78 }
79 //______________________________________________________________________
80 AliITSPlaneEffSPD::AliITSPlaneEffSPD(const AliITSPlaneEffSPD &s) : AliITSPlaneEff(s), 
81 //fHis(s.fHis),
82 fHisResX(0),
83 fHisResZ(0),
84 fHisResXZ(0),
85 fHisClusterSize(0),
86 fHisResXclu(0),
87 fHisResZclu(0),
88 fHisResXchip(0),
89 fHisResZchip(0),
90 fProfResXvsPhi(0),
91 fProfResZvsDip(0),
92 fProfResXvsPhiclu(0),
93 fProfResZvsDipclu(0),
94 fHisTrackErrX(0),
95 fHisTrackErrZ(0),
96 fHisClusErrX(0),
97 fHisClusErrZ(0)
98 {
99     //     Copy Constructor
100     // Inputs:
101     //    AliITSPlaneEffSPD &s The original class for which
102     //                                this class is a copy of
103     // Outputs:
104     //    none.
105     // Return:
106
107  for (UInt_t i=0; i<kNModule*kNChip; i++){
108     fFound[i]=s.fFound[i];
109     fTried[i]=s.fTried[i];
110  }
111  if(fHis) { 
112    InitHistos();
113    for(Int_t i=0; i<kNHisto; i++) {
114       s.fHisResX[i]->Copy(*fHisResX[i]);
115       s.fHisResZ[i]->Copy(*fHisResZ[i]);
116       s.fHisResXZ[i]->Copy(*fHisResXZ[i]);
117       s.fHisClusterSize[i]->Copy(*fHisClusterSize[i]);
118       for(Int_t clu=0; clu<kNclu; clu++) {  // clu=0 --> cluster size 1
119         s.fHisResXclu[i][clu]->Copy(*fHisResXclu[i][clu]);
120         s.fHisResZclu[i][clu]->Copy(*fHisResZclu[i][clu]);
121         s.fProfResXvsPhiclu[i][clu]->Copy(*fProfResXvsPhiclu[i][clu]);
122         s.fProfResZvsDipclu[i][clu]->Copy(*fProfResZvsDipclu[i][clu]);
123       }
124       for(Int_t chip=0; chip<kNChip; chip++) { 
125         s.fHisResXchip[i][chip]->Copy(*fHisResXchip[i][chip]);
126         s.fHisResZchip[i][chip]->Copy(*fHisResZchip[i][chip]);
127       }
128       s.fProfResXvsPhi[i]->Copy(*fProfResXvsPhi[i]);
129       s.fProfResZvsDip[i]->Copy(*fProfResZvsDip[i]);
130       s.fHisTrackErrX[i]->Copy(*fHisTrackErrX[i]);
131       s.fHisTrackErrZ[i]->Copy(*fHisTrackErrZ[i]);
132       s.fHisClusErrX[i]->Copy(*fHisClusErrX[i]);
133       s.fHisClusErrZ[i]->Copy(*fHisClusErrZ[i]);
134    }
135  }
136 }
137 //_________________________________________________________________________
138 AliITSPlaneEffSPD& AliITSPlaneEffSPD::operator+=(const AliITSPlaneEffSPD &add){
139     //    Add-to-me operator
140     // Inputs:
141     //    const AliITSPlaneEffSPD &add  simulation class to be added
142     // Outputs:
143     //    none.
144     // Return:
145     //    none
146     for (UInt_t i=0; i<kNModule*kNChip; i++){
147       fFound[i] += add.fFound[i];
148       fTried[i] += add.fTried[i];
149     }
150     if(fHis && add.fHis) {
151       for(Int_t i=0; i<kNHisto; i++) {
152         fHisResX[i]->Add(add.fHisResX[i]); 
153         fHisResZ[i]->Add(add.fHisResZ[i]); 
154         fHisResXZ[i]->Add(add.fHisResXZ[i]); 
155         fHisClusterSize[i]->Add(add.fHisClusterSize[i]); 
156         for(Int_t clu=0; clu<kNclu; clu++) {  // clu=0 --> cluster size 1
157           fHisResXclu[i][clu]->Add(add.fHisResXclu[i][clu]); 
158           fHisResZclu[i][clu]->Add(add.fHisResZclu[i][clu]); 
159           fProfResXvsPhiclu[i][clu]->Add(add.fProfResXvsPhiclu[i][clu]);
160           fProfResZvsDipclu[i][clu]->Add(add.fProfResZvsDipclu[i][clu]);
161         }
162         for(Int_t chip=0; chip<kNChip; chip++) {  
163           fHisResXchip[i][chip]->Add(add.fHisResXchip[i][chip]); 
164           fHisResZchip[i][chip]->Add(add.fHisResZchip[i][chip]); 
165         }
166         fProfResXvsPhi[i]->Add(add.fProfResXvsPhi[i]);
167         fProfResZvsDip[i]->Add(add.fProfResZvsDip[i]);
168         fHisTrackErrX[i]->Add(add.fHisTrackErrX[i]);
169         fHisTrackErrZ[i]->Add(add.fHisTrackErrZ[i]);
170         fHisClusErrX[i]->Add(add.fHisClusErrX[i]);
171         fHisClusErrZ[i]->Add(add.fHisClusErrZ[i]);
172       }
173     }
174     return *this;
175 }
176 //______________________________________________________________________
177 AliITSPlaneEffSPD&  AliITSPlaneEffSPD::operator=(const
178                                            AliITSPlaneEffSPD &s){
179     //    Assignment operator
180     // Inputs:
181     //    AliITSPlaneEffSPD &s The original class for which
182     //                                this class is a copy of
183     // Outputs:
184     //    none.
185     // Return:
186  
187     if(this==&s) return *this;
188     s.Copy(*this);
189     return *this;
190 }
191 //______________________________________________________________________
192 void AliITSPlaneEffSPD::Copy(TObject &obj) const {
193   // protected method. copy this to obj
194   AliITSPlaneEff::Copy(obj);
195   AliITSPlaneEffSPD& target = (AliITSPlaneEffSPD &) obj;
196   for(Int_t i=0;i<kNModule*kNChip;i++) {
197       target.fFound[i] = fFound[i];
198       target.fTried[i] = fTried[i];
199   }
200   CopyHistos(target);
201   return;
202 }
203 //_______________________________________________________________________
204 void AliITSPlaneEffSPD::CopyHistos(AliITSPlaneEffSPD &target) const {
205   // protected method: copy histos from this to target
206   target.fHis  = fHis; // this is redundant only in some cases. Leave as it is.
207   if(fHis) {
208     target.fHisResX=new TH1F*[kNHisto];
209     target.fHisResZ=new TH1F*[kNHisto];
210     target.fHisResXZ=new TH2F*[kNHisto];
211     target.fHisClusterSize=new TH2I*[kNHisto];
212     target.fHisResXclu=new TH1F**[kNHisto];
213     target.fHisResZclu=new TH1F**[kNHisto];
214     target.fHisResXchip=new TH1F**[kNHisto];
215     target.fHisResZchip=new TH1F**[kNHisto];
216     target.fProfResXvsPhi=new TProfile*[kNHisto];
217     target.fProfResZvsDip=new TProfile*[kNHisto];
218     target.fProfResXvsPhiclu=new TProfile**[kNHisto];
219     target.fProfResZvsDipclu=new TProfile**[kNHisto];
220     target.fHisTrackErrX=new TH1F*[kNHisto];
221     target.fHisTrackErrZ=new TH1F*[kNHisto];
222     target.fHisClusErrX=new TH1F*[kNHisto];
223     target.fHisClusErrZ=new TH1F*[kNHisto];
224     for(Int_t i=0; i<kNHisto; i++) {
225       target.fHisResX[i] = new TH1F(*fHisResX[i]);
226       target.fHisResZ[i] = new TH1F(*fHisResZ[i]);
227       target.fHisResXZ[i] = new TH2F(*fHisResXZ[i]);
228       target.fHisClusterSize[i] = new TH2I(*fHisClusterSize[i]);
229       target.fHisResXclu[i]=new TH1F*[kNclu];
230       target.fHisResZclu[i]=new TH1F*[kNclu];
231       target.fProfResXvsPhiclu[i]=new TProfile*[kNclu];
232       target.fProfResZvsDipclu[i]=new TProfile*[kNclu];
233       for(Int_t clu=0; clu<kNclu; clu++) {  // clu=0 --> cluster size 1
234         target.fHisResXclu[i][clu] = new TH1F(*fHisResXclu[i][clu]);
235         target.fHisResZclu[i][clu] = new TH1F(*fHisResZclu[i][clu]);
236         target.fProfResXvsPhiclu[i][clu] = new TProfile(*fProfResXvsPhiclu[i][clu]);
237         target.fProfResZvsDipclu[i][clu] = new TProfile(*fProfResZvsDipclu[i][clu]);
238       }
239       target.fHisResXchip[i]=new TH1F*[kNChip];
240       target.fHisResZchip[i]=new TH1F*[kNChip];
241       for(Int_t chip=0; chip<kNChip; chip++) {  
242         target.fHisResXchip[i][chip] = new TH1F(*fHisResXchip[i][chip]);
243         target.fHisResZchip[i][chip] = new TH1F(*fHisResZchip[i][chip]);
244       }
245       target.fProfResXvsPhi[i] = new TProfile(*fProfResXvsPhi[i]);
246       target.fProfResZvsDip[i] = new TProfile(*fProfResZvsDip[i]);
247       target.fHisTrackErrX[i] = new TH1F(*fHisTrackErrX[i]);
248       target.fHisTrackErrZ[i] = new TH1F(*fHisTrackErrZ[i]);
249       target.fHisClusErrX[i] = new TH1F(*fHisClusErrX[i]);
250       target.fHisClusErrZ[i] = new TH1F(*fHisClusErrZ[i]);
251     }
252   }
253 return;
254 }
255 /*    commented out by M.Masera 8/3/08
256 //______________________________________________________________________
257 AliITSPlaneEff&  AliITSPlaneEffSPD::operator=(const
258                                            AliITSPlaneEff &s){
259     //    Assignment operator
260     // Inputs:
261     //    AliITSPlaneEffSPD &s The original class for which
262     //                                this class is a copy of
263     // Outputs:
264     //    none.
265     // Return:
266
267     if(&s == this) return *this;
268     AliError("operator=: Not allowed to make a =, use default creater instead");
269     return *this;
270 }
271 */
272 //_______________________________________________________________________
273 Int_t AliITSPlaneEffSPD::GetMissingTracksForGivenEff(Double_t eff, Double_t RelErr,
274           UInt_t im, UInt_t ic) const {
275    
276   //   Estimate the number of tracks still to be collected to attain a 
277   //   given efficiency eff, with relative error RelErr
278   //   Inputs:
279   //         eff    -> Expected efficiency (e.g. those from actual estimate)
280   //         RelErr -> tollerance [0,1] 
281   //         im     -> module number [0,249]
282   //         ic     -> chip number [0,4]
283   //   Outputs: none
284   //   Return: the estimated n. of tracks 
285   //
286 if (im>=kNModule || ic>=kNChip) 
287  {AliError("GetMissingTracksForGivenEff: you asked for a non existing chip");
288  return -1;}
289 else { 
290   UInt_t key=GetKey(im,ic);
291   if(key<kNModule*kNChip) return GetNTracksForGivenEff(eff,RelErr)-fTried[key];
292   else return -1;
293 }
294 }
295 //_________________________________________________________________________
296 Double_t  AliITSPlaneEffSPD::PlaneEff(const UInt_t im,const UInt_t ic) const {
297 // Compute the efficiency for a basic block, 
298 // Inputs:
299 //        im     -> module number [0,249]
300 //        ic     -> chip number [0,4] 
301 if (im>=kNModule || ic>=kNChip) 
302  {AliError("PlaneEff(Uint_t,Uint_t): you asked for a non existing chip"); return -1.;}
303 UInt_t key=GetKey(im,ic);
304 Int_t nf=-1;
305 Int_t nt=-1;
306 if(key<kNModule*kNChip) {
307   nf=fFound[key];
308   nt=fTried[key];
309 }
310 return AliITSPlaneEff::PlaneEff(nf,nt);
311 }
312 //_________________________________________________________________________
313 Double_t  AliITSPlaneEffSPD::ErrPlaneEff(const UInt_t im,const UInt_t ic) const {
314     // Compute the statistical error on efficiency for a basic block,
315     // using binomial statistics 
316     // Inputs:
317     //        im     -> module number [0,249]
318     //        ic     -> chip number [0,4] 
319 if (im>=kNModule || ic>=kNChip) 
320  {AliError("ErrPlaneEff(Uint_t,Uint_t): you asked for a non existing chip"); return -1.;}
321 UInt_t key=GetKey(im,ic);
322 Int_t nf=-1;
323 Int_t nt=-1;
324 if(key<kNModule*kNChip) {
325   nf=fFound[key];
326   nt=fTried[key];
327 }
328 return AliITSPlaneEff::ErrPlaneEff(nf,nt);
329
330 //_________________________________________________________________________
331 Bool_t AliITSPlaneEffSPD::UpDatePlaneEff(const Bool_t Kfound,
332                                          const UInt_t im, const UInt_t ic) {
333   // Update efficiency for a basic block
334 if (im>=kNModule || ic>=kNChip) 
335  {AliError("UpDatePlaneEff: you asked for a non existing chip"); return kFALSE;}
336  UInt_t key=GetKey(im,ic);
337  if(key<kNModule*kNChip) {
338    fTried[key]++;
339    if(Kfound) fFound[key]++;
340    return kTRUE;
341  }
342 }
343 //_________________________________________________________________________
344 UInt_t AliITSPlaneEffSPD::GetChipFromCol(const UInt_t col) const {
345   // get chip given the column
346 if(col>=kNCol*kNChip) 
347  {AliDebug(1,Form("GetChipFromCol: you asked for a non existing column %d",col)); return 10;}
348 return col/kNCol;
349 }
350 //__________________________________________________________________________
351 UInt_t AliITSPlaneEffSPD::GetKey(const UInt_t mod, const UInt_t chip) const {
352   // get key given a basic block
353 if(mod>=kNModule || chip>=kNChip)
354   {AliWarning("GetKey: you asked for a non existing block"); return 99999;}
355 return mod*kNChip+chip;
356 }
357 //__________________________________________________________________________
358 UInt_t AliITSPlaneEffSPD::GetModFromKey(const UInt_t key) const {
359   // get mod. from key
360 if(key>=kNModule*kNChip)
361   {AliError("GetModFromKey: you asked for a non existing key"); return 9999;}
362 return key/kNChip;
363 }
364 //__________________________________________________________________________
365 UInt_t AliITSPlaneEffSPD::GetChipFromKey(const UInt_t key) const {
366   // retrieves chip from key
367 if(key>=kNModule*kNChip)
368   {AliError("GetChipFromKey: you asked for a non existing key"); return 999;}
369 return (key%(kNModule*kNChip))%kNChip;
370 }
371 //__________________________________________________________________________
372 void AliITSPlaneEffSPD::GetModAndChipFromKey(const UInt_t key,UInt_t& mod,UInt_t& chip) const {
373   // get module and chip from a key
374 if(key>=kNModule*kNChip)
375   {AliError("GetModAndChipFromKey: you asked for a non existing key"); 
376   mod=9999;
377   chip=999;
378   return;}
379 mod=key/kNChip;
380 chip=(key%(kNModule*kNChip))%kNChip;
381 return;
382 }
383 //____________________________________________________________________________
384 Double_t AliITSPlaneEffSPD::LivePlaneEff(UInt_t key) const {
385   // returns plane efficieny after adding the fraction of sensor which is bad
386 if(key>=kNModule*kNChip)
387   {AliError("LivePlaneEff: you asked for a non existing key");
388    return -1.;}
389 Double_t leff=AliITSPlaneEff::LivePlaneEff(0); // this just for the Warning
390 leff=PlaneEff(key)+GetFracBad(key);
391 return leff>1?1:leff;
392 }
393 //____________________________________________________________________________
394 Double_t AliITSPlaneEffSPD::ErrLivePlaneEff(UInt_t key) const {
395   // returns error on live plane efficiency
396 if(key>=kNModule*kNChip)
397   {AliError("ErrLivePlaneEff: you asked for a non existing key");
398    return -1.;}
399 Int_t nf=fFound[key];
400 Double_t triedInLive=GetFracLive(key)*fTried[key];
401 Int_t nt=TMath::Max(nf,TMath::Nint(triedInLive));
402 return AliITSPlaneEff::ErrPlaneEff(nf,nt); // for the time being: to be checked
403 }
404 //_____________________________________________________________________________
405 Double_t AliITSPlaneEffSPD::GetFracLive(const UInt_t key) const {
406   // returns the fraction of the sensor which is OK
407 if(key>=kNModule*kNChip)
408   {AliError("GetFracLive: you asked for a non existing key");
409    return -1.;}
410     // Compute the fraction of bad (dead+noisy) detector 
411 UInt_t dead=0,noisy=0;
412 GetDeadAndNoisyInChip(key,dead,noisy);
413 Double_t live=dead+noisy;
414 live/=(kNRow*kNCol);
415 return 1.-live;
416 }
417 //_____________________________________________________________________________
418 void AliITSPlaneEffSPD::GetDeadAndNoisyInChip(const UInt_t key,
419       UInt_t& nrDeadInChip, UInt_t& nrNoisyInChip) const {
420   // returns the number of dead and noisy pixels
421 nrDeadInChip=0;
422 nrNoisyInChip=0;
423 if(key>=kNModule*kNChip)
424   {AliError("GetDeadAndNoisyInChip: you asked for a non existing key");
425    return;}
426     // Compute the number of bad (dead+noisy) pixel in a chip
427 //
428 if(!fInitCDBCalled) 
429   {AliError("GetDeadAndNoisyInChip: CDB not inizialized: call InitCDB first");
430    return;};
431 AliCDBManager* man = AliCDBManager::Instance();
432 // retrieve map of dead Pixel 
433 AliCDBEntry *cdbSPDDead = man->Get("ITS/Calib/SPDDead", fRunNumber);
434 TObjArray* spdDead;
435 if(cdbSPDDead) {
436   spdDead = (TObjArray*)cdbSPDDead->GetObject();
437   if(!spdDead) 
438   {AliError("GetDeadAndNoisyInChip: SPDDead not found in CDB");
439    return;}
440 } else {
441   AliError("GetDeadAndNoisyInChip: did not find Calib/SPDDead.");
442   return;
443 }
444 // retrieve map of noisy Pixel 
445 AliCDBEntry *cdbSPDNoisy = man->Get("ITS/Calib/SPDNoisy", fRunNumber);
446 TObjArray* spdNoisy;
447 if(cdbSPDNoisy) {
448   spdNoisy = (TObjArray*)cdbSPDNoisy->GetObject();
449   if(!spdNoisy) 
450   {AliError("GetDeadAndNoisyInChip: SPDNoisy not found in CDB");
451    return;}
452 } else {
453   AliError("GetDeadAndNoisyInChip: did not find Calib/SPDNoisy.");
454   return;
455 }
456 //
457 UInt_t mod=GetModFromKey(key);
458 UInt_t chip=GetChipFromKey(key);
459 // count number of dead
460 AliITSCalibrationSPD* calibSPD=(AliITSCalibrationSPD*) spdDead->At(mod);
461 UInt_t nrDead = calibSPD->GetNrBad();
462 for (UInt_t index=0; index<nrDead; index++) {
463   if(GetChipFromCol(calibSPD->GetBadColAt(index))==chip) nrDeadInChip++;
464 }
465 calibSPD=(AliITSCalibrationSPD*) spdNoisy->At(mod);
466 UInt_t nrNoisy = calibSPD->GetNrBad();
467 for (UInt_t index=0; index<nrNoisy; index++) {
468   if(GetChipFromCol(calibSPD->GetBadColAt(index))==chip) nrNoisyInChip++;
469 }
470 return;
471 }
472 //_____________________________________________________________________________
473 Double_t AliITSPlaneEffSPD::GetFracBad(const UInt_t key) const {
474   // returns 1-fractional live
475 if(key>=kNModule*kNChip)
476   {AliError("GetFracBad: you asked for a non existing key");
477    return -1.;}
478 return 1.-GetFracLive(key);
479 }
480 //_____________________________________________________________________________
481 Bool_t AliITSPlaneEffSPD::WriteIntoCDB() const {
482 // write onto CDB
483 if(!fInitCDBCalled)
484   {AliError("WriteIntoCDB: CDB not inizialized. Call InitCDB first");
485    return kFALSE;}
486 // to be written properly: now only for debugging 
487   AliCDBMetaData *md= new AliCDBMetaData(); // metaData describing the object
488   //md->SetObjectClassName("AliITSPlaneEff");
489   md->SetResponsible("Giuseppe Eugenio Bruno");
490   md->SetBeamPeriod(0);
491   md->SetAliRootVersion("head 19/11/07"); //root version
492   AliCDBId id("ITS/PlaneEff/PlaneEffSPD",0,AliCDBRunRange::Infinity()); 
493   AliITSPlaneEffSPD eff; 
494   eff=*this;
495   Bool_t r=AliCDBManager::Instance()->GetDefaultStorage()->Put(&eff,id,md);
496   delete md;
497   return r;
498 }
499 //_____________________________________________________________________________
500 Bool_t AliITSPlaneEffSPD::ReadFromCDB() {
501 // read from CDB
502 if(!fInitCDBCalled)
503   {AliError("ReadFromCDB: CDB not inizialized. Call InitCDB first");
504    return kFALSE;}
505 AliCDBEntry *cdbEntry = AliCDBManager::Instance()->Get("ITS/PlaneEff/PlaneEffSPD",fRunNumber);
506 if(!cdbEntry) return kFALSE;
507 AliITSPlaneEffSPD* eff= (AliITSPlaneEffSPD*)cdbEntry->GetObject();
508 if(this==eff) return kFALSE;
509 if(fHis) CopyHistos(*eff); // If histos already exist then copy them to eff
510 eff->Copy(*this);          // copy everything (statistics and histos) from eff to this
511 return kTRUE;
512 }
513 //_____________________________________________________________________________
514 Bool_t AliITSPlaneEffSPD::AddFromCDB(AliCDBId *cdbId) {
515 AliCDBEntry *cdbEntry=0;
516 if (!cdbId) {
517   if(!fInitCDBCalled)  
518     {AliError("ReadFromCDB: CDB not inizialized. Call InitCDB first"); return kFALSE;}
519   cdbEntry = AliCDBManager::Instance()->Get("ITS/PlaneEff/PlaneEffSPD",fRunNumber);
520 } else {
521   cdbEntry = AliCDBManager::Instance()->Get(*cdbId);
522 }
523 if(!cdbEntry) return kFALSE;
524 AliITSPlaneEffSPD* eff= (AliITSPlaneEffSPD*)cdbEntry->GetObject();
525 *this+=*eff;
526 return kTRUE;
527 }
528 //_____________________________________________________________________________
529 UInt_t AliITSPlaneEffSPD::GetKeyFromDetLocCoord(Int_t ilay, Int_t idet, 
530                                                 Float_t, Float_t locz) const {
531 // method to locate a basic block from Detector Local coordinate (to be used in tracking)
532 UInt_t key=999999;
533 if(ilay<0 || ilay>1) 
534   {AliError("GetKeyFromDetLocCoord: you asked for a non existing layer");
535    return key;}
536 if(ilay==0 && (idet<0 || idet>79))
537  {AliError("GetKeyFromDetLocCoord: you asked for a non existing detector");
538    return key;}
539 if(ilay==1 && (idet<0 || idet>159))
540  {AliError("GetKeyFromDetLocCoord: you asked for a non existing detector");
541    return key;}
542
543 UInt_t mod=idet;
544 if(ilay==1) mod+=80;
545 key=GetKey(mod,GetChipFromCol(GetColFromLocZ(locz)));
546 return key;
547 }
548 //_____________________________________________________________________________
549 UInt_t AliITSPlaneEffSPD::GetColFromLocZ(Float_t zloc) const {
550 // method to retrieve column number from the local z coordinate
551   UInt_t col=0;
552   AliITSsegmentationSPD spd;
553   Int_t ix,iz;
554   if(spd.LocalToDet(0,zloc,ix,iz)) col+=iz;
555   else {
556     AliDebug(1,Form("cannot compute column number from local z=%f",zloc));
557     col=99999;}
558   return col;
559 /*
560 const Float_t kconv = 1.0E-04; // converts microns to cm.
561 Float_t bz[160];
562 for(Int_t i=000;i<160;i++) bz[i] = 425.0; // most are 425 microns except below
563 bz[ 31] = bz[ 32] = 625.0; // first chip boundry
564 bz[ 63] = bz[ 64] = 625.0; // first chip boundry
565 bz[ 95] = bz[ 96] = 625.0; // first chip boundry
566 bz[127] = bz[128] = 625.0; // first chip boundry
567 //
568 Int_t j=-1;
569 Float_t dz=0;
570 for(Int_t i=000;i<160;i++) dz+=bz[i];
571 dz = -0.5*kconv*dz;
572 if(zloc<dz || zloc>-1*dz) { // outside z range
573   AliDebug(1,Form("GetColFromLocZ: cannot compute column number from local z=%f",zloc));
574   return 99999;}
575 for(j=0;j<160;j++){
576   dz += kconv*bz[j];
577   if(zloc<dz) break;
578 } // end for j
579 col+=j;
580 //
581 return col;
582 */
583 }
584 //________________________________________________________
585 Bool_t AliITSPlaneEffSPD::GetBlockBoundaries(const UInt_t key, Float_t& xmn,Float_t& xmx,
586                                              Float_t& zmn,Float_t& zmx) const {
587 //
588 //  This method return the geometrical boundaries of the active volume of a given 
589 //  basic block, in the detector reference system.
590 //  Input: unique key to locate a basic block.
591 //  
592 //  Output: Ymin, Ymax, Zmin, Zmax of a basic block (chip for SPD)
593 //  Return: kTRUE if computation was succesfully, kFALSE otherwise
594 //
595 if(key>=kNModule*kNChip)
596   {AliWarning("GetBlockBoundaries: you asked for a non existing key"); return kFALSE;}
597 UInt_t chip=GetChipFromKey(key);
598 zmn=GetLocZFromCol(chip*kNCol);
599 zmx=GetLocZFromCol((chip+1)*kNCol);
600 xmn=GetLocXFromRow(0);
601 xmx=GetLocXFromRow(kNRow);
602 //
603 Float_t tmp=zmn;
604 if(zmx<zmn) {zmn=zmx; zmx=tmp;}
605 tmp=xmn;
606 if(xmx<xmn) {xmn=xmx; xmx=tmp;}
607 return kTRUE;
608 }
609 //________________________________________________________
610 Float_t AliITSPlaneEffSPD::GetLocXFromRow(const UInt_t row) const {
611 // 
612 //  This method return the local (i.e. detector reference system) lower x coordinate 
613 //  of the row. To get the central value of a given row, you can do 
614 //  1/2*[LocXFromRow(row)+LocXFromRow(row+1)].
615 //
616 //  Input: row number in the range [0,kNRow] 
617 //  Output: lower local X coordinate of this row.
618 //
619 if(row>kNRow)  // not >= ! allow also computation of upper limit of the last row. 
620   {AliError("LocYFromRow: you asked for a non existing row"); return 9999999.;}
621 // Use only AliITSsegmentationSPD
622 AliITSsegmentationSPD spd;
623 Double_t dummy,x;
624 if(row==kNRow) spd.CellBoundries((Int_t)row-1,0,dummy,x,dummy,dummy);
625 else spd.CellBoundries((Int_t)row,0,x,dummy,dummy,dummy);
626 return (Float_t)x;
627
628 }
629 //________________________________________________________
630 Float_t AliITSPlaneEffSPD::GetLocZFromCol(const UInt_t col) const {
631 //
632 //  This method return the local (i.e. detector reference system) lower Z coordinate
633 //  of the column. To get the central value of a given column, you can do
634 //  1/2*[LocZFromCol(col)+LocZFromCol(col+1)].
635 //
636 //  Input: col number in the range [0,kNChip*kNCol]
637 //  Output: lower local Y coordinate of this row.
638 //
639 if(col>kNChip*kNCol) // not >= ! allow also computation of upper limit of the last column
640   {AliError("LocZFromCol: you asked for a non existing column"); return 9999999.;}
641 // Use only AliITSsegmentationSPD
642 AliITSsegmentationSPD spd;
643 Double_t dummy,y;
644 if(col==kNChip*kNCol) spd.CellBoundries(0,(Int_t)col-1,dummy,dummy,dummy,y);
645 else spd.CellBoundries(0,(Int_t)col,dummy,dummy,y,dummy);
646 return (Float_t)y;
647
648 }
649 //__________________________________________________________
650 void AliITSPlaneEffSPD::InitHistos() {
651   // for the moment let's create the histograms 
652   // module by  module
653   TString histnameResX="HistResX_mod_",aux;
654   TString histnameResZ="HistResZ_mod_";
655   TString histnameResXZ="HistResXZ_mod_";
656   TString histnameClusterType="HistClusterType_mod_";
657   TString histnameResXclu="HistResX_mod_";
658   TString histnameResZclu="HistResZ_mod_";
659   TString histnameResXchip="HistResX_mod_";
660   TString histnameResZchip="HistResZ_mod_";
661   TString profnameResXvsPhi="ProfResXvsPhi_mod_";
662   TString profnameResZvsDip="ProfResZvsDip_mod_";
663   TString profnameResXvsPhiclu="ProfResXvsPhi_mod_";
664   TString profnameResZvsDipclu="ProfResZvsDip_mod_";
665   TString histnameTrackErrX="HistTrackErrX_mod_";
666   TString histnameTrackErrZ="HistTrackErrZ_mod_";
667   TString histnameClusErrX="HistClusErrX_mod_";
668   TString histnameClusErrZ="HistClusErrZ_mod_";
669 //
670
671   TH1::AddDirectory(kFALSE);
672
673   fHisResX=new TH1F*[kNHisto];
674   fHisResZ=new TH1F*[kNHisto];
675   fHisResXZ=new TH2F*[kNHisto];
676   fHisClusterSize=new TH2I*[kNHisto];
677   fHisResXclu=new TH1F**[kNHisto];
678   fHisResZclu=new TH1F**[kNHisto];
679   fHisResXchip=new TH1F**[kNHisto];
680   fHisResZchip=new TH1F**[kNHisto];
681   fProfResXvsPhi=new TProfile*[kNHisto];
682   fProfResZvsDip=new TProfile*[kNHisto];
683   fProfResXvsPhiclu=new TProfile**[kNHisto];
684   fProfResZvsDipclu=new TProfile**[kNHisto];
685   fHisTrackErrX=new TH1F*[kNHisto];
686   fHisTrackErrZ=new TH1F*[kNHisto];
687   fHisClusErrX=new TH1F*[kNHisto];
688   fHisClusErrZ=new TH1F*[kNHisto];
689
690   for (Int_t nhist=0;nhist<kNHisto;nhist++){
691     aux=histnameResX;
692     aux+=nhist;
693     fHisResX[nhist]=new TH1F("histname","histname",1600,-0.32,0.32); // +- 3200 micron; 1 bin=4 micron
694     fHisResX[nhist]->SetName(aux.Data());
695     fHisResX[nhist]->SetTitle(aux.Data());
696
697     aux=histnameResZ;
698     aux+=nhist;
699     fHisResZ[nhist]=new TH1F("histname","histname",1200,-0.48,0.48); // +-4800 micron; 1 bin=8 micron
700     fHisResZ[nhist]->SetName(aux.Data());
701     fHisResZ[nhist]->SetTitle(aux.Data());
702
703     aux=histnameResXZ;
704     aux+=nhist;
705     fHisResXZ[nhist]=new TH2F("histname","histname",80,-0.16,0.16,80,-0.32,0.32); // binning:
706     fHisResXZ[nhist]->SetName(aux.Data());                                         // 40 micron in x;
707     fHisResXZ[nhist]->SetTitle(aux.Data());                                        // 80 micron in z;
708
709     aux=histnameClusterType;
710     aux+=nhist;
711     fHisClusterSize[nhist]=new TH2I("histname","histname",10,0.5,10.5,10,0.5,10.5);
712     fHisClusterSize[nhist]->SetName(aux.Data());
713     fHisClusterSize[nhist]->SetTitle(aux.Data());
714
715     fHisResXclu[nhist]=new TH1F*[kNclu];
716     fHisResZclu[nhist]=new TH1F*[kNclu];
717     for(Int_t clu=0; clu<kNclu; clu++) {  // clu=0 --> cluster size 1
718       aux=histnameResXclu;
719       aux+=nhist;
720       aux+="_clu_";
721       aux+=clu+1; // clu=0 --> cluster size 1
722       fHisResXclu[nhist][clu]=new TH1F("histname","histname",1600,-0.32,0.32); // +- 3200 micron; 1 bin=4 micron
723       fHisResXclu[nhist][clu]->SetName(aux.Data());
724       fHisResXclu[nhist][clu]->SetTitle(aux.Data());
725
726       aux=histnameResZclu;
727       aux+=nhist;
728       aux+="_clu_";
729       aux+=clu+1; // clu=0 --> cluster size 1
730       fHisResZclu[nhist][clu]=new TH1F("histname","histname",1200,-0.48,0.48); // +-4800 micron; 1 bin=8 micron
731       fHisResZclu[nhist][clu]->SetName(aux.Data());
732       fHisResZclu[nhist][clu]->SetTitle(aux.Data());
733     }
734
735     fHisResXchip[nhist]=new TH1F*[kNChip];
736     fHisResZchip[nhist]=new TH1F*[kNChip];
737     for(Int_t chip=0; chip<kNChip; chip++) { 
738       aux=histnameResXchip;
739       aux+=nhist;
740       aux+="_chip_";
741       aux+=chip; 
742       fHisResXchip[nhist][chip]=new TH1F("histname","histname",800,-0.32,0.32); // +- 3200 micron; 1 bin=8 micron
743       fHisResXchip[nhist][chip]->SetName(aux.Data());
744       fHisResXchip[nhist][chip]->SetTitle(aux.Data());
745
746       aux=histnameResZchip;
747       aux+=nhist;
748       aux+="_chip_";
749       aux+=chip;
750       fHisResZchip[nhist][chip]=new TH1F("histname","histname",300,-0.48,0.48); // +-4800 micron; 1 bin=32 micron
751       fHisResZchip[nhist][chip]->SetName(aux.Data());
752       fHisResZchip[nhist][chip]->SetTitle(aux.Data());
753     }
754
755     aux=histnameTrackErrX;
756     aux+=nhist;
757     fHisTrackErrX[nhist]=new TH1F("histname","histname",400,0.,0.32); // 0-3200 micron; 1 bin=8 micron
758     fHisTrackErrX[nhist]->SetName(aux.Data());
759     fHisTrackErrX[nhist]->SetTitle(aux.Data());
760
761     aux=histnameTrackErrZ;
762     aux+=nhist;
763     fHisTrackErrZ[nhist]=new TH1F("histname","histname",200,0.,0.32); // 0-3200 micron; 1 bin=16 micron
764     fHisTrackErrZ[nhist]->SetName(aux.Data());
765     fHisTrackErrZ[nhist]->SetTitle(aux.Data());
766
767     aux=histnameClusErrX;
768     aux+=nhist;
769     fHisClusErrX[nhist]=new TH1F("histname","histname",400,0.,0.08); //  0-800 micron; 1 bin=2 micron
770     fHisClusErrX[nhist]->SetName(aux.Data());
771     fHisClusErrX[nhist]->SetTitle(aux.Data());
772
773     aux=histnameClusErrZ;
774     aux+=nhist;
775     fHisClusErrZ[nhist]=new TH1F("histname","histname",400,0.,0.32); //  0-3200 micron; 1 bin=8 micron
776     fHisClusErrZ[nhist]->SetName(aux.Data());
777     fHisClusErrZ[nhist]->SetTitle(aux.Data());
778
779     aux=profnameResXvsPhi;
780     aux+=nhist;
781     fProfResXvsPhi[nhist]=new TProfile("histname","histname",40,-40.,40.0); // binning: range:  -40°- 40°
782     fProfResXvsPhi[nhist]->SetName(aux.Data());                             //          bin width: 2°
783     fProfResXvsPhi[nhist]->SetTitle(aux.Data());
784
785     aux=profnameResZvsDip;
786     aux+=nhist;
787     fProfResZvsDip[nhist]=new TProfile("histname","histname",48,-72.,72.0); // binning: range:  -70°-4°
788     fProfResZvsDip[nhist]->SetName(aux.Data());                             //          bin width: 3°
789     fProfResZvsDip[nhist]->SetTitle(aux.Data());
790
791     fProfResXvsPhiclu[nhist]=new TProfile*[kNclu];
792     fProfResZvsDipclu[nhist]=new TProfile*[kNclu];
793     for(Int_t clu=0; clu<kNclu; clu++) {  // clu=0 --> cluster size 1
794       aux=profnameResXvsPhiclu;
795       aux+=nhist;
796       aux+="_clu_";
797       aux+=clu+1; // clu=0 --> cluster size 1
798       fProfResXvsPhiclu[nhist][clu]=new TProfile("histname","histname",40,-40.,40.0); // binning: range:  -40°- 40
799       fProfResXvsPhiclu[nhist][clu]->SetName(aux.Data());                             //          bin width: 2°
800       fProfResXvsPhiclu[nhist][clu]->SetTitle(aux.Data());
801
802       aux=profnameResZvsDipclu;
803       aux+=nhist;
804       aux+="_clu_";
805       aux+=clu+1; // clu=0 --> cluster size 1
806       fProfResZvsDipclu[nhist][clu]= new TProfile("histname","histname",48,-72.,72.0); // binning: range:  -70°-7°
807       fProfResZvsDipclu[nhist][clu]->SetName(aux.Data());                              //      bin width: 3°
808       fProfResZvsDipclu[nhist][clu]->SetTitle(aux.Data());
809     }
810
811   } // end loop on module
812
813   TH1::AddDirectory(kTRUE);
814
815 return;
816 }
817 //__________________________________________________________
818 void AliITSPlaneEffSPD::DeleteHistos() {
819   if(fHisResX) {
820     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResX[i];
821     delete [] fHisResX; fHisResX=0;
822   }
823   if(fHisResZ) {
824     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResZ[i];
825     delete [] fHisResZ; fHisResZ=0;
826   }
827   if(fHisResXZ) {
828     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResXZ[i];
829     delete [] fHisResXZ; fHisResXZ=0;
830   }
831   if(fHisClusterSize) {
832     for (Int_t i=0; i<kNHisto; i++ ) delete fHisClusterSize[i];
833     delete [] fHisClusterSize; fHisClusterSize=0;
834   }
835   if(fHisResXclu) {
836     for (Int_t i=0; i<kNHisto; i++ ) {
837       for (Int_t clu=0; clu<kNclu; clu++) if (fHisResXclu[i][clu]) delete fHisResXclu[i][clu];
838       delete [] fHisResXclu[i];
839     }
840     delete [] fHisResXclu;
841     fHisResXclu = 0;
842   }
843   if(fHisResZclu) {
844     for (Int_t i=0; i<kNHisto; i++ ) {
845       for (Int_t clu=0; clu<kNclu; clu++) if (fHisResZclu[i][clu]) delete fHisResZclu[i][clu];
846       delete [] fHisResZclu[i];
847     }
848     delete [] fHisResZclu;
849     fHisResZclu = 0;
850   }
851   if(fHisResXchip) {
852     for (Int_t i=0; i<kNHisto; i++ ) {
853       for (Int_t chip=0; chip<kNChip; chip++) if (fHisResXchip[i][chip]) delete fHisResXchip[i][chip];
854       delete [] fHisResXchip[i];
855     }
856     delete [] fHisResXchip;
857     fHisResXchip = 0;
858   }
859   if(fHisResZchip) {
860     for (Int_t i=0; i<kNHisto; i++ ) {
861       for (Int_t chip=0; chip<kNChip; chip++) if (fHisResZchip[i][chip]) delete fHisResZchip[i][chip];
862       delete [] fHisResZchip[i];
863     }
864     delete [] fHisResZchip;
865     fHisResZchip = 0;
866   }
867   if(fHisTrackErrX) {
868     for (Int_t i=0; i<kNHisto; i++ ) delete fHisTrackErrX[i];
869     delete [] fHisTrackErrX; fHisTrackErrX=0;
870   }
871   if(fHisTrackErrZ) {
872     for (Int_t i=0; i<kNHisto; i++ ) delete fHisTrackErrZ[i];
873     delete [] fHisTrackErrZ; fHisTrackErrZ=0;
874   }
875   if(fHisClusErrX) {
876     for (Int_t i=0; i<kNHisto; i++ ) delete fHisClusErrX[i];
877     delete [] fHisClusErrX; fHisClusErrX=0;
878   }
879   if(fHisClusErrZ) {
880     for (Int_t i=0; i<kNHisto; i++ ) delete fHisClusErrZ[i];
881     delete [] fHisClusErrZ; fHisClusErrZ=0;
882   }
883   if(fProfResXvsPhi) {
884     for (Int_t i=0; i<kNHisto; i++ ) delete fProfResXvsPhi[i];
885     delete [] fProfResXvsPhi; fProfResXvsPhi=0;
886   }
887   if(fProfResZvsDip) {
888     for (Int_t i=0; i<kNHisto; i++ ) delete fProfResZvsDip[i];
889     delete [] fProfResZvsDip; fProfResZvsDip=0;
890   }
891   if(fProfResXvsPhiclu) {
892     for (Int_t i=0; i<kNHisto; i++ ) {
893       for (Int_t clu=0; clu<kNclu; clu++) if (fProfResXvsPhiclu[i][clu]) delete fProfResXvsPhiclu[i][clu];
894       delete [] fProfResXvsPhiclu[i];
895     }
896     delete [] fProfResXvsPhiclu;
897     fProfResXvsPhiclu = 0;
898   }
899   if(fProfResZvsDipclu) {
900     for (Int_t i=0; i<kNHisto; i++ ) {
901       for (Int_t clu=0; clu<kNclu; clu++) if (fProfResZvsDipclu[i][clu]) delete fProfResZvsDipclu[i][clu];
902       delete [] fProfResZvsDipclu[i];
903     }
904     delete [] fProfResZvsDipclu;
905     fProfResZvsDipclu = 0;
906   }
907
908 return;
909 }
910 //__________________________________________________________
911 Bool_t AliITSPlaneEffSPD::FillHistos(UInt_t key, Bool_t found, 
912                                      Float_t *tr, Float_t *clu, Int_t *csize, Float_t *angtrkmod) {
913 // this method fill the histograms
914 // input: - key: unique key of the basic block 
915 //        - found: Boolean to asses whether a cluster has been associated to the track or not 
916 //        - tr[0],tr[1] local X and Z coordinates of the track prediction, respectively
917 //        - tr[2],tr[3] error on local X and Z coordinates of the track prediction, respectively
918 //        - clu[0],clu[1] local X and Z coordinates of the cluster associated to the track, respectively
919 //        - clu[2],clu[3] error on local X and Z coordinates of the cluster associated to the track, respectively
920 //        - csize[0][1] cluster size in X and Z, respectively
921 //        - angtrkmod[0],angtrkmod[1]  
922 // output: kTRUE if filling was succesfull kFALSE otherwise
923 // side effects: updating of the histograms. 
924 //
925   if (!fHis) {
926     AliWarning("FillHistos: histograms do not exist! Call SetCreateHistos(kTRUE) first");
927     return kFALSE;
928   }
929   if(key>=kNModule*kNChip)
930     {AliWarning("FillHistos: you asked for a non existing key"); return kFALSE;}
931   Int_t id=GetModFromKey(key);
932   Int_t chip=GetChipFromKey(key);
933   if(id>=kNHisto) 
934     {AliWarning("FillHistos: you want to fill a non-existing histos"); return kFALSE;}
935   if(found) {
936     Float_t resx=tr[0]-clu[0];
937     Float_t resz=tr[1]-clu[1];
938     fHisResX[id]->Fill(resx);
939     fHisResZ[id]->Fill(resz);
940     fHisResXZ[id]->Fill(resx,resz);
941     fHisClusterSize[id]->Fill((Double_t)csize[0],(Double_t)csize[1]);
942     if(csize[0]>0 &&  csize[0]<=kNclu) fHisResXclu[id][csize[0]-1]->Fill(resx);
943     if(csize[1]>0 &&  csize[1]<=kNclu) fHisResZclu[id][csize[1]-1]->Fill(resz);
944     fHisResXchip[id][chip]->Fill(resx);
945     fHisResZchip[id][chip]->Fill(resz);
946     fProfResXvsPhi[id]->Fill(angtrkmod[0],resx);
947     fProfResZvsDip[id]->Fill(angtrkmod[1],resz);
948     if(csize[0]>0 &&  csize[0]<=kNclu) fProfResXvsPhiclu[id][csize[0]-1]->Fill(angtrkmod[0],resx);
949     if(csize[1]>0 &&  csize[1]<=kNclu) fProfResZvsDipclu[id][csize[1]-1]->Fill(angtrkmod[1],resz);
950   }
951   fHisTrackErrX[id]->Fill(tr[2]);
952   fHisTrackErrZ[id]->Fill(tr[3]);
953   fHisClusErrX[id]->Fill(clu[2]);
954   fHisClusErrZ[id]->Fill(clu[3]);
955   return kTRUE;
956 }
957 //__________________________________________________________
958 Bool_t AliITSPlaneEffSPD::WriteHistosToFile(TString filename, Option_t* option) {
959   //
960   // Saves the histograms into a tree and saves the trees into a file
961   //
962   if (!fHis) return kFALSE;
963   if (filename.IsNull() || filename.IsWhitespace()) {
964      AliWarning("WriteHistosToFile: null output filename!");
965      return kFALSE;
966   }
967   char branchname[30];
968   TFile *hFile=new TFile(filename.Data(),option,
969                          "The File containing the TREEs with ITS PlaneEff Histos");
970   TTree *SPDTree=new TTree("SPDTree","Tree whith Residuals and Cluster Type distributions for SPD");
971   TH1F *histZ,*histX;
972   TH2F *histXZ;
973   TH2I *histClusterType;
974   TH1F *histXclu[kNclu];
975   TH1F *histZclu[kNclu];
976   TH1F *histXchip[kNChip];
977   TH1F *histZchip[kNChip];
978   TH1F *histTrErrZ,*histTrErrX;
979   TH1F *histClErrZ,*histClErrX;
980   TProfile *profXvsPhi,*profZvsDip;
981   TProfile *profXvsPhiclu[kNclu],*profZvsDipclu[kNclu];
982
983   histZ=new TH1F();
984   histX=new TH1F();
985   histXZ=new TH2F();
986   histClusterType=new TH2I();
987   for(Int_t clu=0;clu<kNclu;clu++) {
988     histXclu[clu]=new TH1F();
989     histZclu[clu]=new TH1F();
990   }
991   for(Int_t chip=0;chip<kNChip;chip++) {
992     histXchip[chip]=new TH1F();
993     histZchip[chip]=new TH1F();
994   }
995   histTrErrX=new TH1F();
996   histTrErrZ=new TH1F();
997   histClErrX=new TH1F();
998   histClErrZ=new TH1F();
999   profXvsPhi=new TProfile();
1000   profZvsDip=new TProfile();
1001   for(Int_t clu=0;clu<kNclu;clu++) {
1002     profXvsPhiclu[clu]=new TProfile();
1003     profZvsDipclu[clu]=new TProfile();
1004   }
1005
1006
1007   SPDTree->Branch("histX","TH1F",&histX,128000,0);
1008   SPDTree->Branch("histZ","TH1F",&histZ,128000,0);
1009   SPDTree->Branch("histXZ","TH2F",&histXZ,128000,0);
1010   SPDTree->Branch("histClusterType","TH2I",&histClusterType,128000,0);
1011   for(Int_t clu=0;clu<kNclu;clu++) {
1012     sprintf(branchname,"histXclu_%d",clu+1);
1013     SPDTree->Branch(branchname,"TH1F",&histXclu[clu],128000,0);
1014     sprintf(branchname,"histZclu_%d",clu+1);
1015     SPDTree->Branch(branchname,"TH1F",&histZclu[clu],128000,0);
1016   }
1017   for(Int_t chip=0;chip<kNChip;chip++) {
1018     sprintf(branchname,"histXchip_%d",chip);
1019     SPDTree->Branch(branchname,"TH1F",&histXchip[chip],128000,0);
1020     sprintf(branchname,"histZchip_%d",chip);
1021     SPDTree->Branch(branchname,"TH1F",&histZchip[chip],128000,0);
1022   }
1023   SPDTree->Branch("histTrErrX","TH1F",&histTrErrX,128000,0);
1024   SPDTree->Branch("histTrErrZ","TH1F",&histTrErrZ,128000,0);
1025   SPDTree->Branch("histClErrX","TH1F",&histClErrX,128000,0);
1026   SPDTree->Branch("histClErrZ","TH1F",&histClErrZ,128000,0);
1027   SPDTree->Branch("profXvsPhi","TProfile",&profXvsPhi,128000,0);
1028   SPDTree->Branch("profZvsDip","TProfile",&profZvsDip,128000,0);
1029   for(Int_t clu=0;clu<kNclu;clu++) {
1030     sprintf(branchname,"profXvsPhiclu_%d",clu+1);
1031     SPDTree->Branch(branchname,"TProfile",&profXvsPhiclu[clu],128000,0);
1032     sprintf(branchname,"profZvsDipclu_%d",clu+1);
1033     SPDTree->Branch(branchname,"TProfile",&profZvsDipclu[clu],128000,0);
1034   }
1035
1036   for(Int_t j=0;j<kNHisto;j++){
1037     histX=fHisResX[j];
1038     histZ=fHisResZ[j];
1039     histXZ=fHisResXZ[j];
1040     histClusterType=fHisClusterSize[j];
1041     for(Int_t clu=0;clu<kNclu;clu++) {
1042       histXclu[clu]=fHisResXclu[j][clu];
1043       histZclu[clu]=fHisResZclu[j][clu];
1044     }
1045     for(Int_t chip=0;chip<kNChip;chip++) {
1046       histXchip[chip]=fHisResXchip[j][chip];
1047       histZchip[chip]=fHisResZchip[j][chip];
1048     }
1049     histTrErrX=fHisTrackErrX[j];
1050     histTrErrZ=fHisTrackErrZ[j];
1051     histClErrX=fHisClusErrX[j];
1052     histClErrZ=fHisClusErrZ[j];
1053     profXvsPhi=fProfResXvsPhi[j];
1054     profZvsDip=fProfResZvsDip[j];
1055     for(Int_t clu=0;clu<kNclu;clu++) {
1056       profXvsPhiclu[clu]=fProfResXvsPhiclu[j][clu];
1057       profZvsDipclu[clu]=fProfResZvsDipclu[j][clu];
1058     }
1059
1060     SPDTree->Fill();
1061   }
1062   hFile->Write();
1063   hFile->Close();
1064 return kTRUE;
1065 }
1066 //__________________________________________________________
1067 Bool_t AliITSPlaneEffSPD::ReadHistosFromFile(TString filename) {
1068   //
1069   // Read histograms from an already existing file 
1070   //
1071   if (!fHis) return kFALSE;
1072   if (filename.IsNull() || filename.IsWhitespace()) {
1073      AliWarning("ReadHistosFromFile: incorrect output filename!");
1074      return kFALSE;
1075   }
1076   char branchname[30];
1077
1078   TH1F *h  = 0;
1079   TH2F *h2 = 0;
1080   TH2I *h2i= 0;
1081   TProfile *p = 0;
1082
1083   TFile *file=TFile::Open(filename.Data(),"READONLY");
1084
1085   if (!file || file->IsZombie()) {
1086     AliWarning(Form("Can't open %s !",filename.Data()));
1087     delete file;
1088     return kFALSE;
1089   }
1090   TTree *tree = (TTree*) file->Get("SPDTree");
1091
1092   TBranch *histX = (TBranch*) tree->GetBranch("histX");
1093   TBranch *histZ = (TBranch*) tree->GetBranch("histZ");
1094   TBranch *histXZ = (TBranch*) tree->GetBranch("histXZ");
1095   TBranch *histClusterType = (TBranch*) tree->GetBranch("histClusterType");
1096    
1097   TBranch *histXclu[kNclu], *histZclu[kNclu];
1098   for(Int_t clu=0; clu<kNclu; clu++) {
1099     sprintf(branchname,"histXclu_%d",clu+1);
1100     histXclu[clu]= (TBranch*) tree->GetBranch(branchname);
1101     sprintf(branchname,"histZclu_%d",clu+1);
1102     histZclu[clu]= (TBranch*) tree->GetBranch(branchname);
1103   }
1104
1105   TBranch *histXchip[kNChip], *histZchip[kNChip];
1106   for(Int_t chip=0; chip<kNChip; chip++) {
1107     sprintf(branchname,"histXchip_%d",chip);
1108     histXchip[chip]= (TBranch*) tree->GetBranch(branchname);
1109     sprintf(branchname,"histZchip_%d",chip);
1110     histZchip[chip]= (TBranch*) tree->GetBranch(branchname);
1111   }
1112
1113   TBranch *histTrErrX = (TBranch*) tree->GetBranch("histTrErrX");
1114   TBranch *histTrErrZ = (TBranch*) tree->GetBranch("histTrErrZ");
1115   TBranch *histClErrX = (TBranch*) tree->GetBranch("histClErrX");
1116   TBranch *histClErrZ = (TBranch*) tree->GetBranch("histClErrZ");
1117   TBranch *profXvsPhi = (TBranch*) tree->GetBranch("profXvsPhi");
1118   TBranch *profZvsDip = (TBranch*) tree->GetBranch("profZvsDip");
1119
1120   TBranch *profXvsPhiclu[kNclu], *profZvsDipclu[kNclu];
1121   for(Int_t clu=0; clu<kNclu; clu++) {
1122     sprintf(branchname,"profXvsPhiclu_%d",clu+1);
1123     profXvsPhiclu[clu]= (TBranch*) tree->GetBranch(branchname);
1124     sprintf(branchname,"profZvsDipclu_%d",clu+1);
1125     profZvsDipclu[clu]= (TBranch*) tree->GetBranch(branchname);
1126   }
1127
1128   gROOT->cd();
1129
1130   Int_t nevent = (Int_t)histX->GetEntries();
1131   if(nevent!=kNHisto) 
1132     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1133   histX->SetAddress(&h);
1134   for(Int_t j=0;j<kNHisto;j++){
1135     delete h; h=0;
1136     histX->GetEntry(j);
1137     fHisResX[j]->Add(h);
1138   }
1139
1140   nevent = (Int_t)histZ->GetEntries();
1141   if(nevent!=kNHisto) 
1142     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1143   histZ->SetAddress(&h);
1144   for(Int_t j=0;j<kNHisto;j++){
1145     delete h; h=0;
1146     histZ->GetEntry(j);
1147     fHisResZ[j]->Add(h);
1148   }
1149
1150   nevent = (Int_t)histXZ->GetEntries();
1151   if(nevent!=kNHisto) 
1152     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1153   histXZ->SetAddress(&h2);
1154   for(Int_t j=0;j<kNHisto;j++){
1155     delete h2; h2=0;
1156     histXZ->GetEntry(j);
1157     fHisResXZ[j]->Add(h2);
1158   }
1159
1160   nevent = (Int_t)histClusterType->GetEntries();
1161   if(nevent!=kNHisto) 
1162     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1163   histClusterType->SetAddress(&h2i);
1164   for(Int_t j=0;j<kNHisto;j++){
1165     delete h2i; h2i=0;
1166     histClusterType->GetEntry(j);
1167     fHisClusterSize[j]->Add(h2i);
1168   }
1169
1170   for(Int_t clu=0; clu<kNclu; clu++) {
1171
1172     nevent = (Int_t)histXclu[clu]->GetEntries();
1173     if(nevent!=kNHisto)
1174       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1175     histXclu[clu]->SetAddress(&h);
1176     for(Int_t j=0;j<kNHisto;j++){
1177       delete h; h=0;
1178       histXclu[clu]->GetEntry(j);
1179       fHisResXclu[j][clu]->Add(h);
1180     }
1181
1182    nevent = (Int_t)histZclu[clu]->GetEntries();
1183     if(nevent!=kNHisto)
1184       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1185     histZclu[clu]->SetAddress(&h);
1186     for(Int_t j=0;j<kNHisto;j++){
1187       delete h; h=0;
1188       histZclu[clu]->GetEntry(j);
1189       fHisResZclu[j][clu]->Add(h);
1190     }
1191   }
1192
1193
1194     for(Int_t chip=0; chip<kNChip; chip++) {
1195
1196     nevent = (Int_t)histXchip[chip]->GetEntries();
1197     if(nevent!=kNHisto)
1198       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1199     histXchip[chip]->SetAddress(&h);
1200     for(Int_t j=0;j<kNHisto;j++){
1201       delete h; h=0;
1202       histXchip[chip]->GetEntry(j);
1203       fHisResXchip[j][chip]->Add(h);
1204     }
1205
1206     nevent = (Int_t)histZchip[chip]->GetEntries();
1207     if(nevent!=kNHisto)
1208       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1209     histZchip[chip]->SetAddress(&h);
1210     for(Int_t j=0;j<kNHisto;j++){
1211       delete h; h=0;
1212       histZchip[chip]->GetEntry(j);
1213       fHisResZchip[j][chip]->Add(h);
1214     }
1215   }
1216
1217   nevent = (Int_t)histTrErrX->GetEntries(); 
1218   if(nevent!=kNHisto)
1219     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1220   histTrErrX->SetAddress(&h);
1221   for(Int_t j=0;j<kNHisto;j++){
1222     delete h; h=0;
1223     histTrErrX->GetEntry(j);
1224     fHisTrackErrX[j]->Add(h);
1225   }
1226
1227   nevent = (Int_t)histTrErrZ->GetEntries();
1228   if(nevent!=kNHisto)
1229     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1230   histTrErrZ->SetAddress(&h);
1231   for(Int_t j=0;j<kNHisto;j++){
1232     delete h; h=0;
1233     histTrErrZ->GetEntry(j);
1234     fHisTrackErrZ[j]->Add(h);
1235   }
1236
1237   nevent = (Int_t)histClErrX->GetEntries();
1238   if(nevent!=kNHisto)
1239     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1240   histClErrX->SetAddress(&h);
1241   for(Int_t j=0;j<kNHisto;j++){
1242     delete h; h=0;
1243     histClErrX->GetEntry(j);
1244     fHisClusErrX[j]->Add(h);
1245   }
1246
1247   nevent = (Int_t)histClErrZ->GetEntries();
1248   if(nevent!=kNHisto)
1249     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1250   histClErrZ->SetAddress(&h);
1251   for(Int_t j=0;j<kNHisto;j++){
1252     delete h; h=0;
1253     histClErrZ->GetEntry(j);
1254     fHisClusErrZ[j]->Add(h);
1255   }
1256
1257   nevent = (Int_t)profXvsPhi->GetEntries();
1258   if(nevent!=kNHisto)
1259     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1260   profXvsPhi->SetAddress(&p);
1261   for(Int_t j=0;j<kNHisto;j++){
1262     delete p; p=0;
1263     profXvsPhi->GetEntry(j);
1264     fProfResXvsPhi[j]->Add(p);
1265   }
1266
1267   nevent = (Int_t)profZvsDip->GetEntries();
1268   if(nevent!=kNHisto)
1269     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1270   profZvsDip->SetAddress(&p);
1271   for(Int_t j=0;j<kNHisto;j++){
1272     delete p; p=0;
1273     profZvsDip->GetEntry(j);
1274     fProfResZvsDip[j]->Add(p);
1275   }
1276
1277     for(Int_t clu=0; clu<kNclu; clu++) {
1278
1279     nevent = (Int_t)profXvsPhiclu[clu]->GetEntries();
1280     if(nevent!=kNHisto)
1281       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1282     profXvsPhiclu[clu]->SetAddress(&p);
1283     for(Int_t j=0;j<kNHisto;j++){
1284       delete p; p=0;
1285       profXvsPhiclu[clu]->GetEntry(j);
1286       fProfResXvsPhiclu[j][clu]->Add(p);
1287     }
1288
1289     nevent = (Int_t)profZvsDipclu[clu]->GetEntries();
1290     if(nevent!=kNHisto)
1291       {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
1292     profZvsDipclu[clu]->SetAddress(&p);
1293     for(Int_t j=0;j<kNHisto;j++){
1294       delete p; p=0;
1295       profZvsDipclu[clu]->GetEntry(j);
1296       fProfResZvsDipclu[j][clu]->Add(p);
1297     }
1298   }
1299
1300
1301   delete h;   h=0;
1302   delete h2;  h2=0;
1303   delete h2i; h2i=0;
1304   delete p;   p=0;
1305
1306   if (file) {
1307     file->Close();
1308   }
1309 return kTRUE;
1310 }
1311