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