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