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