]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSPlaneEffSSD.cxx
Fixes for Coverity warnings
[u/mrichter/AliRoot.git] / ITS / AliITSPlaneEffSSD.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 module by module efficiency of the SSD,        
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 "AliITSPlaneEffSSD.h"
33 #include "AliLog.h"
34 #include "AliCDBStorage.h"
35 #include "AliCDBEntry.h"
36 #include "AliCDBManager.h"
37 //#include "AliCDBRunRange.h"
38 #include "AliITSCalibrationSSD.h"
39
40 ClassImp(AliITSPlaneEffSSD)     
41 //______________________________________________________________________
42 AliITSPlaneEffSSD::AliITSPlaneEffSSD():
43   AliITSPlaneEff(),
44   fHisResX(0),
45   fHisResZ(0),
46   fHisResXZ(0),
47   fHisClusterSize(0),
48   fHisTrackErrX(0),
49   fHisTrackErrZ(0),
50   fHisClusErrX(0),
51   fHisClusErrZ(0){
52   for (UInt_t i=0; i<kNModule; i++){
53     fFound[i]=0;
54     fTried[i]=0;
55   }
56   // default constructor
57   AliDebug(1,Form("Calling default constructor"));
58 }
59 //______________________________________________________________________
60 AliITSPlaneEffSSD::~AliITSPlaneEffSSD(){
61     // destructor
62     // Inputs:
63     //    none.
64     // Outputs:
65     //    none.
66     // Return:
67     //     none.
68     DeleteHistos();
69 }
70 //______________________________________________________________________
71 AliITSPlaneEffSSD::AliITSPlaneEffSSD(const AliITSPlaneEffSSD &s) : AliITSPlaneEff(s),
72 fHisResX(0),
73 fHisResZ(0),
74 fHisResXZ(0),
75 fHisClusterSize(0),
76 fHisTrackErrX(0),
77 fHisTrackErrZ(0),
78 fHisClusErrX(0),
79 fHisClusErrZ(0)
80 {
81     //     Copy Constructor
82     // Inputs:
83     //    AliITSPlaneEffSSD &s The original class for which
84     //                                this class is a copy of
85     // Outputs:
86     //    none.
87     // Return:
88
89  for (UInt_t i=0; i<kNModule; i++){
90     fFound[i]=s.fFound[i];
91     fTried[i]=s.fTried[i];
92  }
93  if(fHis) {
94    InitHistos();
95    for(Int_t i=0; i<kNHisto; i++) {
96       s.fHisResX[i]->Copy(*fHisResX[i]);
97       s.fHisResZ[i]->Copy(*fHisResZ[i]);
98       s.fHisResXZ[i]->Copy(*fHisResXZ[i]);
99       s.fHisClusterSize[i]->Copy(*fHisClusterSize[i]);
100       s.fHisTrackErrX[i]->Copy(*fHisTrackErrX[i]);
101       s.fHisTrackErrZ[i]->Copy(*fHisTrackErrZ[i]);
102       s.fHisClusErrX[i]->Copy(*fHisTrackErrZ[i]);
103       s.fHisClusErrZ[i]->Copy(*fHisClusErrZ[i]);
104    }
105  }
106 }
107 //_________________________________________________________________________
108 AliITSPlaneEffSSD& AliITSPlaneEffSSD::operator+=(const AliITSPlaneEffSSD &add){
109     //    Add-to-me operator
110     // Inputs:
111     //    const AliITSPlaneEffSSD &add  simulation class to be added
112     // Outputs:
113     //    none.
114     // Return:
115     //    none
116     for (UInt_t i=0; i<kNModule; i++){
117       fFound[i] += add.fFound[i];
118       fTried[i] += add.fTried[i];
119     }
120     if(fHis && add.fHis) {
121       for(Int_t i=0; i<kNHisto; i++) {
122         fHisResX[i]->Add(add.fHisResX[i]);
123         fHisResZ[i]->Add(add.fHisResZ[i]);
124         fHisResXZ[i]->Add(add.fHisResXZ[i]);
125         fHisClusterSize[i]->Add(add.fHisClusterSize[i]);
126         fHisTrackErrX[i]->Add(add.fHisTrackErrX[i]);
127         fHisTrackErrZ[i]->Add(add.fHisTrackErrZ[i]);
128         fHisClusErrX[i]->Add(add.fHisTrackErrZ[i]);
129         fHisClusErrZ[i]->Add(add.fHisClusErrZ[i]);
130       }
131     }
132     return *this;
133 }
134 //______________________________________________________________________
135 AliITSPlaneEffSSD&  AliITSPlaneEffSSD::operator=(const
136                                            AliITSPlaneEffSSD &s){
137     //    Assignment operator
138     // Inputs:
139     //    AliITSPlaneEffSSD &s The original class for which
140     //                                this class is a copy of
141     // Outputs:
142     //    none.
143     // Return:
144  
145     if(this==&s) return *this;
146     s.Copy(*this);
147     return *this;
148 }
149 //______________________________________________________________________
150 void AliITSPlaneEffSSD::Copy(TObject &obj) const {
151   // protected method. copy this to obj
152   AliITSPlaneEff::Copy(obj);
153   AliITSPlaneEffSSD& target = (AliITSPlaneEffSSD &) obj;
154   for(Int_t i=0;i<kNModule;i++) {
155       target.fFound[i] = fFound[i];
156       target.fTried[i] = fTried[i];
157   }
158   CopyHistos(target);
159   return;
160 }
161 //_______________________________________________________________________
162 void AliITSPlaneEffSSD::CopyHistos(AliITSPlaneEffSSD &target) const {
163   // protected method: copy histos from this to target
164   target.fHis  = fHis; // this is redundant only in some cases. Leave as it is.
165   if(fHis) {
166     target.fHisResX=new TH1F*[kNHisto];
167     target.fHisResZ=new TH1F*[kNHisto];
168     target.fHisResXZ=new TH2F*[kNHisto];
169     target.fHisClusterSize=new TH2I*[kNHisto];
170     target.fHisTrackErrX=new TH1F*[kNHisto];
171     target.fHisTrackErrZ=new TH1F*[kNHisto];
172     target.fHisClusErrX=new TH1F*[kNHisto];
173     target.fHisClusErrZ=new TH1F*[kNHisto];
174     for(Int_t i=0; i<kNHisto; i++) {
175       target.fHisResX[i] = new TH1F(*fHisResX[i]);
176       target.fHisResZ[i] = new TH1F(*fHisResZ[i]);
177       target.fHisResXZ[i] = new TH2F(*fHisResXZ[i]);
178       target.fHisClusterSize[i] = new TH2I(*fHisClusterSize[i]);
179       target.fHisTrackErrX[i] = new TH1F(*fHisTrackErrX[i]);
180       target.fHisTrackErrZ[i] = new TH1F(*fHisTrackErrZ[i]);
181       target.fHisClusErrX[i]  = new TH1F(*fHisClusErrX[i]);
182       target.fHisClusErrZ[i]  = new TH1F(*fHisClusErrZ[i]);
183     }
184   }
185 return;
186 }
187
188 //_______________________________________________________________________
189 Int_t AliITSPlaneEffSSD::GetMissingTracksForGivenEff(Double_t eff, Double_t RelErr,
190           UInt_t im) const {
191    
192   //   Estimate the number of tracks still to be collected to attain a 
193   //   given efficiency eff, with relative error RelErr
194   //   Inputs:
195   //         eff    -> Expected efficiency (e.g. those from actual estimate)
196   //         RelErr -> tollerance [0,1] 
197   //         im     -> module number [0,1697]
198   //   Outputs: none
199   //   Return: the estimated n. of tracks 
200   //
201 if (im>=kNModule) 
202  {AliError("GetMissingTracksForGivenEff: you asked for a non existing module");
203  return -1;}
204 else return GetNTracksForGivenEff(eff,RelErr)-fTried[GetKey(im)];
205 }
206 //_________________________________________________________________________
207 Double_t  AliITSPlaneEffSSD::PlaneEff(const UInt_t im) const {
208 // Compute the efficiency for a basic block, 
209 // Inputs:
210 //        im     -> module number [0,1697]
211 if (im>=kNModule) 
212  {AliError("PlaneEff(UInt_t): you asked for a non existing module"); return -1.;}
213  Int_t nf=fFound[GetKey(im)];
214  Int_t nt=fTried[GetKey(im)];
215 return AliITSPlaneEff::PlaneEff(nf,nt);
216 }
217 //_________________________________________________________________________
218 Double_t  AliITSPlaneEffSSD::ErrPlaneEff(const UInt_t im) const {
219     // Compute the statistical error on efficiency for a basic block,
220     // using binomial statistics 
221     // Inputs:
222     //        im     -> module number [0,1697]
223 if (im>=kNModule) 
224  {AliError("ErrPlaneEff(UInt_t): you asked for a non existing module"); return -1.;}
225 Int_t nf=fFound[GetKey(im)];
226 Int_t nt=fTried[GetKey(im)];
227 return AliITSPlaneEff::ErrPlaneEff(nf,nt);
228
229 //_________________________________________________________________________
230 Bool_t AliITSPlaneEffSSD::UpDatePlaneEff(const Bool_t Kfound, const UInt_t im) {
231   // Update efficiency for a basic block
232 if (im>=kNModule) 
233  {AliError("UpDatePlaneEff: you asked for a non existing module"); return kFALSE;}
234  fTried[GetKey(im)]++;
235  if(Kfound) fFound[GetKey(im)]++;
236  return kTRUE;
237 }
238 //_________________________________________________________________________
239 UInt_t AliITSPlaneEffSSD::GetKey(const UInt_t mod) const {
240   // get key given a basic block
241 if(mod>=kNModule)
242   {AliError("GetKey: you asked for a non existing block"); return 99999;}
243 return mod;
244 }
245 //__________________________________________________________________________
246 UInt_t AliITSPlaneEffSSD::GetModFromKey(const UInt_t key) const {
247   // get mod. from key
248 if(key>=kNModule)
249   {AliError("GetModFromKey: you asked for a non existing key"); return 9999;}
250 return key;
251 }
252 //__________________________________________________________________________
253 Double_t AliITSPlaneEffSSD::LivePlaneEff(UInt_t key) const {
254   // returns plane efficieny after adding the fraction of sensor which is bad
255 if(key>=kNModule)
256   {AliError("LivePlaneEff: you asked for a non existing key");
257    return -1.;}
258 Double_t leff=AliITSPlaneEff::LivePlaneEff(0); // this just for the Warning
259 leff=PlaneEff(key)+GetFracBad(key);
260 return leff>1?1:leff;
261 }
262 //____________________________________________________________________________
263 Double_t AliITSPlaneEffSSD::ErrLivePlaneEff(UInt_t key) const {
264   // returns error on live plane efficiency
265 if(key>=kNModule)
266   {AliError("ErrLivePlaneEff: you asked for a non existing key");
267    return -1.;}
268 Int_t nf=fFound[key];
269 Double_t triedInLive=GetFracLive(key)*fTried[key];
270 Int_t nt=TMath::Max(nf,TMath::Nint(triedInLive));
271 return AliITSPlaneEff::ErrPlaneEff(nf,nt); // for the time being: to be checked
272 }
273 //_____________________________________________________________________________
274 Double_t AliITSPlaneEffSSD::GetFracLive(const UInt_t key) const {
275   // returns the fraction of the sensor area which is OK (neither noisy nor dead)
276   // As for now, it computes only the fraction of good strips / total strips. 
277   // If this fraction is large, then the computation is a good approximation. 
278   // In any case, it is a lower limit of the fraction of the live area. 
279   // The next upgrades would be to add the fraction of area of superoposition 
280   // between bad N-side strips and bad P-side strips.  
281 if(key>=kNModule)
282   {AliError("GetFracLive: you asked for a non existing key");
283    return -1.;}
284 AliInfo("GetFracLive: it computes only the fraction of working strips (N+P side) / total strips");
285 UInt_t bad=0;
286 GetBadInModule(key,bad);
287 Double_t live=bad;
288 live/=(kNChip*kNSide*kNStrip);
289 return 1.-live;
290 }
291 //_____________________________________________________________________________
292 void AliITSPlaneEffSSD::GetBadInModule(const UInt_t key, UInt_t& nrBadInMod) const {
293   // returns the number of dead and noisy strips (sum of P and N sides).
294 nrBadInMod=0;
295 if(key>=kNModule)
296   {AliError("GetBadInModule: you asked for a non existing key");
297    return;}
298     // Compute the number of bad (dead+noisy) pixel in a module
299 //
300 if(!fInitCDBCalled) 
301   {AliError("GetBadInModule: CDB not inizialized: call InitCDB first");
302    return;};
303 AliCDBManager* man = AliCDBManager::Instance();
304 // retrieve map of dead Pixel 
305 AliCDBEntry *cdbSSD = man->Get("ITS/Calib/BadChannelsSSD", fRunNumber);
306 TObjArray* ssdEntry;
307 if(cdbSSD) {
308   ssdEntry = (TObjArray*)cdbSSD->GetObject();
309   if(!ssdEntry) 
310   {AliError("GetBadInChip: SSDEntry not found in CDB");
311    return;}
312 } else {
313   AliError("GetBadInChip: did not find Calib/BadChannelsSSD");
314   return;
315 }
316 //
317 //UInt_t mod=GetModFromKey(key);
318 //
319 //AliITSBadChannelsSSD* badchannels=(AliITSBadChannelsSSD*) ssdEntry->At(mod);
320 // count the  number of bad channels on the p side
321 //nrBadInMod += (badchannels->GetBadPChannelsList()).GetSize();
322 // add the  number of bad channels on the s side
323 //nrBadInMod += (badchannels->GetBadNChannelsList()).GetSize();
324 return;
325 }
326 //_____________________________________________________________________________
327 Double_t AliITSPlaneEffSSD::GetFracBad(const UInt_t key) const {
328   // returns 1-fractional live
329 if(key>=kNModule)
330   {AliError("GetFracBad: you asked for a non existing key");
331    return -1.;}
332 return 1.-GetFracLive(key);
333 }
334 //_____________________________________________________________________________
335 Bool_t AliITSPlaneEffSSD::WriteIntoCDB() const {
336 // write onto CDB
337 if(!fInitCDBCalled)
338   {AliError("WriteIntoCDB: CDB not inizialized. Call InitCDB first");
339    return kFALSE;}
340 // to be written properly: now only for debugging 
341   AliCDBMetaData *md= new AliCDBMetaData(); // metaData describing the object
342   md->SetObjectClassName("AliITSPlaneEff");
343   md->SetResponsible("Giuseppe Eugenio Bruno");
344   md->SetBeamPeriod(0);
345   md->SetAliRootVersion("head 02/01/08"); //root version
346   AliCDBId id("ITS/PlaneEff/PlaneEffSSD",0,AliCDBRunRange::Infinity()); 
347   AliITSPlaneEffSSD eff; 
348   eff=*this;
349   Bool_t r=AliCDBManager::Instance()->GetDefaultStorage()->Put(&eff,id,md);
350   delete md;
351   return r;
352 }
353 //_____________________________________________________________________________
354 Bool_t AliITSPlaneEffSSD::ReadFromCDB() {
355 // read from CDB
356 if(!fInitCDBCalled)
357   {AliError("ReadFromCDB: CDB not inizialized. Call InitCDB first");
358    return kFALSE;}
359 AliCDBEntry *cdbEntry = AliCDBManager::Instance()->Get("ITS/PlaneEff/PlaneEffSSD",fRunNumber);
360 if(!cdbEntry) return kFALSE;
361 AliITSPlaneEffSSD* eff= (AliITSPlaneEffSSD*)cdbEntry->GetObject();
362 if(this==eff) return kFALSE;
363 if(fHis) CopyHistos(*eff); // If histos already exist then copy them to eff
364 eff->Copy(*this);          // copy everything (statistics and histos) from eff to this
365 return kTRUE;
366 }
367 //_____________________________________________________________________________
368 Bool_t AliITSPlaneEffSSD::AddFromCDB(AliCDBId *cdbId) {
369 AliCDBEntry *cdbEntry=0;
370 if (!cdbId) {
371   if(!fInitCDBCalled)
372     {AliError("ReadFromCDB: CDB not inizialized. Call InitCDB first"); return kFALSE;}
373   cdbEntry = AliCDBManager::Instance()->Get("ITS/PlaneEff/PlaneEffSSD",fRunNumber);
374 } else {
375   cdbEntry = AliCDBManager::Instance()->Get(*cdbId);
376 }
377 if(!cdbEntry) return kFALSE;
378 AliITSPlaneEffSSD* eff= (AliITSPlaneEffSSD*)cdbEntry->GetObject();
379 *this+=*eff;
380 return kTRUE;
381 }
382 //_____________________________________________________________________________
383 UInt_t AliITSPlaneEffSSD::GetKeyFromDetLocCoord(Int_t ilay, Int_t idet,
384                                                 Float_t, Float_t) const {
385 // method to locate a basic block from Detector Local coordinate (to be used in tracking)
386 UInt_t key=999999;
387 if(ilay<4 || ilay>5)
388   {AliError("GetKeyFromDetLocCoord: you asked for a non existing layer");
389    return key;}
390 if(ilay==4 && (idet<0 || idet>747))
391  {AliError("GetKeyFromDetLocCoord: you asked for a non existing detector");
392    return key;}
393 if(ilay==5 && (idet<0 || idet>949))
394  {AliError("GetKeyFromDetLocCoord: you asked for a non existing detector");
395    return key;}
396
397 UInt_t mod=idet;
398 if(ilay==5) mod+=748;
399 key=GetKey(mod);
400 return key;
401 }
402 //__________________________________________________________
403 void AliITSPlaneEffSSD::InitHistos() {
404   // for the moment let's create the histograms
405   // module by  module
406   TString histnameResX="HistResX_mod_",aux;
407   TString histnameResZ="HistResZ_mod_";
408   TString histnameResXZ="HistResXZ_mod_";
409   TString histnameClusterType="HistClusterType_mod_";
410   TString histnameTrackErrX="HistTrackErrX_mod_";
411   TString histnameTrackErrZ="HistTrackErrZ_mod_";
412   TString histnameClusErrX="HistClusErrX_mod_";
413   TString histnameClusErrZ="HistClusErrZ_mod_";
414 //
415
416   TH1::AddDirectory(kFALSE);
417
418   fHisResX=new TH1F*[kNHisto];
419   fHisResZ=new TH1F*[kNHisto];
420   fHisResXZ=new TH2F*[kNHisto];
421   fHisClusterSize=new TH2I*[kNHisto];
422   fHisTrackErrX=new TH1F*[kNHisto]; 
423   fHisTrackErrZ=new TH1F*[kNHisto]; 
424   fHisClusErrX=new TH1F*[kNHisto]; 
425   fHisClusErrZ=new TH1F*[kNHisto];
426
427   for (Int_t nhist=0;nhist<kNHisto;nhist++){
428     aux=histnameResX;
429     aux+=nhist;
430     fHisResX[nhist]=new TH1F("histname","histname",500,-0.10,0.10); // +- 1000 micron; 1 bin=4 micron
431     fHisResX[nhist]->SetName(aux.Data());
432     fHisResX[nhist]->SetTitle(aux.Data());
433
434     aux=histnameResZ;
435     aux+=nhist;
436     fHisResZ[nhist]=new TH1F("histname","histname",750,-0.75,0.75); // +-5000 micron; 1 bin=20 micron
437     fHisResZ[nhist]->SetName(aux.Data());
438     fHisResZ[nhist]->SetTitle(aux.Data());
439
440     aux=histnameResXZ;
441     aux+=nhist;
442     fHisResXZ[nhist]=new TH2F("histname","histname",40,-0.04,0.04,40,-0.32,0.32); // binning:
443                                                                                    // 20 micron in x;
444                                                                                    // 160 micron in z;
445     fHisResXZ[nhist]->SetName(aux.Data());
446     fHisResXZ[nhist]->SetTitle(aux.Data());
447
448     aux=histnameClusterType;
449     aux+=nhist;
450     fHisClusterSize[nhist]=new TH2I("histname","histname",6,0.5,6.5,6,0.5,6.5);
451     fHisClusterSize[nhist]->SetName(aux.Data());
452     fHisClusterSize[nhist]->SetTitle(aux.Data());
453
454     aux=histnameTrackErrX;
455     aux+=nhist;
456     fHisTrackErrX[nhist]=new TH1F("histname","histname",300,0.,0.24); // 0-2400 micron; 1 bin=8 micron
457     fHisTrackErrX[nhist]->SetName(aux.Data());
458     fHisTrackErrX[nhist]->SetTitle(aux.Data());
459
460     aux=histnameTrackErrZ;
461     aux+=nhist;
462     fHisTrackErrZ[nhist]=new TH1F("histname","histname",300,0.,0.48); // 0-4800 micron; 1 bin=16 micron
463     fHisTrackErrZ[nhist]->SetName(aux.Data());
464     fHisTrackErrZ[nhist]->SetTitle(aux.Data());
465
466     aux=histnameClusErrX;
467     aux+=nhist;
468     fHisClusErrX[nhist]=new TH1F("histname","histname",300,0.,0.24); //  0-2400 micron; 1 bin=8 micron
469     fHisClusErrX[nhist]->SetName(aux.Data());
470     fHisClusErrX[nhist]->SetTitle(aux.Data());
471
472     aux=histnameClusErrZ;
473     aux+=nhist;
474     fHisClusErrZ[nhist]=new TH1F("histname","histname",200,0.,0.32); //  0-1600 micron; 1 bin=16 micron
475     fHisClusErrZ[nhist]->SetName(aux.Data());
476     fHisClusErrZ[nhist]->SetTitle(aux.Data());
477
478   }
479
480   TH1::AddDirectory(kTRUE);
481
482 return;
483 }
484 //__________________________________________________________
485 void AliITSPlaneEffSSD::DeleteHistos() {
486   if(fHisResX) {
487     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResX[i];
488     delete [] fHisResX; fHisResX=0;
489   }
490   if(fHisResZ) {
491     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResZ[i];
492     delete [] fHisResZ; fHisResZ=0;
493   }
494   if(fHisResXZ) {
495     for (Int_t i=0; i<kNHisto; i++ ) delete fHisResXZ[i];
496     delete [] fHisResXZ; fHisResXZ=0;
497   }
498   if(fHisClusterSize) {
499     for (Int_t i=0; i<kNHisto; i++ ) delete fHisClusterSize[i];
500     delete [] fHisClusterSize; fHisClusterSize=0;
501   }
502   if(fHisTrackErrX) {
503     for (Int_t i=0; i<kNHisto; i++ ) delete fHisTrackErrX[i];
504     delete [] fHisTrackErrX; fHisTrackErrX=0;
505   }
506   if(fHisTrackErrZ) {
507     for (Int_t i=0; i<kNHisto; i++ ) delete fHisTrackErrZ[i];
508     delete [] fHisTrackErrZ; fHisTrackErrZ=0;
509   }
510   if(fHisClusErrX) {
511     for (Int_t i=0; i<kNHisto; i++ ) delete fHisClusErrX[i];
512     delete [] fHisClusErrX; fHisClusErrX=0;
513   }
514   if(fHisClusErrZ) {
515     for (Int_t i=0; i<kNHisto; i++ ) delete fHisClusErrZ[i];
516     delete [] fHisClusErrZ; fHisClusErrZ=0;
517   }
518
519 return;
520 }
521 //__________________________________________________________
522 Bool_t AliITSPlaneEffSSD::FillHistos(UInt_t key, Bool_t found,
523                                      Float_t *tr, Float_t *clu, Int_t *csize, Float_t*) {
524 // this method fill the histograms
525 // input: - key: unique key of the basic block
526 //        - found: Boolean to asses whether a cluster has been associated to the track or not
527 //        - tr[0],tr[1] local X and Z coordinates of the track prediction, respectively
528 //        - tr[2],tr[3] error on local X and Z coordinates of the track prediction, respectively
529 //        - clu[0],clu[1] local X and Z coordinates of the cluster associated to the track, respectively
530 //        - clu[2],clu[3] error on local X and Z coordinates of the cluster associated to the track, respectively
531 //        - csize[0][1] cluster size in X and Z, respectively
532 // output: kTRUE if filling was succesfull kFALSE otherwise
533 // side effects: updating of the histograms.
534 //
535   if (!fHis) {
536     AliWarning("FillHistos: histograms do not exist! Call SetCreateHistos(kTRUE) first");
537     return kFALSE;
538   }
539   if(key>=kNModule)
540     {AliWarning("FillHistos: you asked for a non existing key"); return kFALSE;}
541   Int_t id=GetModFromKey(key);
542   if(id>=kNHisto)
543     {AliWarning("FillHistos: you want to fill a non-existing histos"); return kFALSE;}
544   if(found) {
545     Float_t resx=tr[0]-clu[0];
546     Float_t resz=tr[1]-clu[1];
547     fHisResX[id]->Fill(resx);
548     fHisResZ[id]->Fill(resz);
549     fHisResXZ[id]->Fill(resx,resz);
550     fHisClusterSize[id]->Fill((Double_t)csize[0],(Double_t)csize[1]);
551   }
552   fHisTrackErrX[id]->Fill(tr[2]);
553   fHisTrackErrZ[id]->Fill(tr[3]);
554   fHisClusErrX[id]->Fill(clu[2]);
555   fHisClusErrZ[id]->Fill(clu[3]);
556   return kTRUE;
557 }
558 //__________________________________________________________
559 Bool_t AliITSPlaneEffSSD::WriteHistosToFile(TString filename, Option_t* option) {
560   //
561   // Saves the histograms into a tree and saves the trees into a file
562   //
563   if (!fHis) return kFALSE;
564   if (filename.IsNull() || filename.IsWhitespace()) {
565      AliWarning("WriteHistosToFile: null output filename!");
566      return kFALSE;
567   }
568
569   TFile *hFile=new TFile(filename.Data(),option,
570                          "The File containing the TREEs with ITS PlaneEff Histos");
571   TTree *SSDTree=new TTree("SSDTree","Tree whith Residuals and Cluster Type distributions for SSD");
572   TH1F *histZ,*histX;
573   TH2F *histXZ;
574   TH2I *histClusterType;
575   TH1F *histTrErrZ,*histTrErrX;
576   TH1F *histClErrZ,*histClErrX;
577
578   histZ=new TH1F();
579   histX=new TH1F();
580   histXZ=new TH2F();
581   histClusterType=new TH2I();
582   histTrErrX=new TH1F();
583   histTrErrZ=new TH1F();
584   histClErrX=new TH1F();
585   histClErrZ=new TH1F();
586
587   SSDTree->Branch("histX","TH1F",&histX,128000,0);
588   SSDTree->Branch("histZ","TH1F",&histZ,128000,0);
589   SSDTree->Branch("histXZ","TH2F",&histXZ,128000,0);
590   SSDTree->Branch("histClusterType","TH2I",&histClusterType,128000,0);
591   SSDTree->Branch("histTrErrX","TH1F",&histTrErrX,128000,0);
592   SSDTree->Branch("histTrErrZ","TH1F",&histTrErrZ,128000,0);
593   SSDTree->Branch("histClErrX","TH1F",&histClErrX,128000,0);
594   SSDTree->Branch("histClErrZ","TH1F",&histClErrZ,128000,0);
595
596   for(Int_t j=0;j<kNHisto;j++){
597     histX=fHisResX[j];
598     histZ=fHisResZ[j];
599     histXZ=fHisResXZ[j];
600     histClusterType=fHisClusterSize[j];
601     histTrErrX=fHisTrackErrX[j];
602     histTrErrZ=fHisTrackErrZ[j];
603     histClErrX=fHisClusErrX[j];
604     histClErrZ=fHisClusErrZ[j];
605
606     SSDTree->Fill();
607   }
608   hFile->Write();
609   hFile->Close();
610 return kTRUE;
611 }
612 //__________________________________________________________
613 Bool_t AliITSPlaneEffSSD::ReadHistosFromFile(TString filename) {
614   //
615   // Read histograms from an already existing file
616   //
617   if (!fHis) return kFALSE;
618   if (filename.IsNull() || filename.IsWhitespace()) {
619      AliWarning("ReadHistosFromFile: incorrect output filename!");
620      return kFALSE;
621   }
622
623   TH1F *h  = 0;
624   TH2F *h2 = 0;
625   TH2I *h2i= 0;
626
627   TFile *file=TFile::Open(filename.Data(),"READONLY");
628
629   if (!file || file->IsZombie()) {
630     AliWarning(Form("Can't open %s !",filename.Data()));
631     delete file;
632     return kFALSE;
633   }
634   TTree *tree = (TTree*) file->Get("SSDTree");
635
636   TBranch *histX = (TBranch*) tree->GetBranch("histX");
637   TBranch *histZ = (TBranch*) tree->GetBranch("histZ");
638   TBranch *histXZ = (TBranch*) tree->GetBranch("histXZ");
639   TBranch *histClusterType = (TBranch*) tree->GetBranch("histClusterType");
640   TBranch *histTrErrX = (TBranch*) tree->GetBranch("histTrErrX");
641   TBranch *histTrErrZ = (TBranch*) tree->GetBranch("histTrErrZ");
642   TBranch *histClErrX = (TBranch*) tree->GetBranch("histClErrX");
643   TBranch *histClErrZ = (TBranch*) tree->GetBranch("histClErrZ");
644
645   gROOT->cd();
646
647   Int_t nevent = (Int_t)histX->GetEntries();
648   if(nevent!=kNHisto)
649     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
650   histX->SetAddress(&h);
651   for(Int_t j=0;j<kNHisto;j++){
652     histX->GetEntry(j);
653     fHisResX[j]->Add(h);
654   }
655
656   nevent = (Int_t)histZ->GetEntries();
657   if(nevent!=kNHisto)
658     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
659   histZ->SetAddress(&h);
660   for(Int_t j=0;j<kNHisto;j++){
661     histZ->GetEntry(j);
662     fHisResZ[j]->Add(h);
663   }
664
665   nevent = (Int_t)histXZ->GetEntries();
666   if(nevent!=kNHisto)
667     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
668   histXZ->SetAddress(&h2);
669   for(Int_t j=0;j<kNHisto;j++){
670     histXZ->GetEntry(j);
671     fHisResXZ[j]->Add(h2);
672   }
673
674   nevent = (Int_t)histClusterType->GetEntries();
675   if(nevent!=kNHisto)
676     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
677   histClusterType->SetAddress(&h2i);
678   for(Int_t j=0;j<kNHisto;j++){
679     histClusterType->GetEntry(j);
680     fHisClusterSize[j]->Add(h2i);
681   }
682
683   nevent = (Int_t)histTrErrX->GetEntries();
684   if(nevent!=kNHisto)
685     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
686   histTrErrX->SetAddress(&h);
687   for(Int_t j=0;j<kNHisto;j++){
688     histTrErrX->GetEntry(j);
689     fHisTrackErrX[j]->Add(h);
690   }
691
692   nevent = (Int_t)histTrErrZ->GetEntries();
693   if(nevent!=kNHisto)
694     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
695   histTrErrZ->SetAddress(&h);
696   for(Int_t j=0;j<kNHisto;j++){
697     histTrErrZ->GetEntry(j);
698     fHisTrackErrZ[j]->Add(h);
699   }
700
701   nevent = (Int_t)histClErrX->GetEntries();
702   if(nevent!=kNHisto)
703     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
704   histClErrX->SetAddress(&h);
705   for(Int_t j=0;j<kNHisto;j++){
706     //delete h; h=0;
707     histClErrX->GetEntry(j);
708     fHisClusErrX[j]->Add(h);
709   }
710
711   nevent = (Int_t)histClErrZ->GetEntries();
712   if(nevent!=kNHisto)
713     {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
714   histClErrZ->SetAddress(&h);
715   for(Int_t j=0;j<kNHisto;j++){
716     //delete h; h=0; 
717     histClErrZ->GetEntry(j);
718     fHisClusErrZ[j]->Add(h);
719   }
720
721   delete h;   
722   delete h2;  
723   delete h2i; 
724
725   if (file) {
726     file->Close();
727     delete file;
728   }
729 return kTRUE;
730 }
731