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