]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITSPlaneEffSPD.cxx
Possibility to use alien or xrd zipped files
[u/mrichter/AliRoot.git] / ITS / AliITSPlaneEffSPD.cxx
CommitLineData
4a66240a 1/**************************************************************************
2 * Copyright(c) 2007-2009, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15///////////////////////////////////////////////////////////////////////////
16// Plane Efficiency class for ITS
17// It is used for chip by chip efficiency of the SPD,
18// evaluated by tracks
19// (Inherits from AliITSPlaneEff)
20// Author: G.E. Bruno
21// giuseppe.bruno@ba.infn.it
22//
23///////////////////////////////////////////////////////////////////////////
24
25/* $Id$ */
26
6344adcc 27#include <TMath.h>
5fbd4fd6 28#include <TH1F.h>
29#include <TFile.h>
30#include <TTree.h>
31#include <TROOT.h>
4a66240a 32#include "AliITSPlaneEffSPD.h"
33#include "AliLog.h"
34#include "AliCDBStorage.h"
35#include "AliCDBEntry.h"
36#include "AliCDBManager.h"
37//#include "AliCDBRunRange.h"
7167ae53 38//#include "AliITSsegmentationSPD.h"
4a66240a 39#include "AliITSCalibrationSPD.h"
40
41ClassImp(AliITSPlaneEffSPD)
42//______________________________________________________________________
43AliITSPlaneEffSPD::AliITSPlaneEffSPD():
5fbd4fd6 44 AliITSPlaneEff(),
45 //fHis(kFALSE),
46 fHisResX(0),
47 fHisResZ(0),
48 fHisResXZ(0),
49 fHisClusterSize(0),
50 fHisResXclu(0),
51 fHisResZclu(0){
4a66240a 52 for (UInt_t i=0; i<kNModule*kNChip; i++){
53 fFound[i]=0;
54 fTried[i]=0;
55 }
56 // default constructor
57 AliDebug(1,Form("Calling default constructor"));
58}
59//______________________________________________________________________
60AliITSPlaneEffSPD::~AliITSPlaneEffSPD(){
61 // destructor
62 // Inputs:
63 // none.
64 // Outputs:
65 // none.
66 // Return:
67 // none.
5fbd4fd6 68 DeleteHistos();
4a66240a 69}
70//______________________________________________________________________
5fbd4fd6 71AliITSPlaneEffSPD::AliITSPlaneEffSPD(const AliITSPlaneEffSPD &s) : AliITSPlaneEff(s),
4a66240a 72//fHis(s.fHis),
5fbd4fd6 73fHisResX(0),
74fHisResZ(0),
75fHisResXZ(0),
76fHisClusterSize(0),
77fHisResXclu(0),
78fHisResZclu(0)
4a66240a 79{
80 // Copy Constructor
81 // Inputs:
82 // AliITSPlaneEffSPD &s The original class for which
83 // this class is a copy of
84 // Outputs:
85 // none.
86 // Return:
87
5fbd4fd6 88 for (UInt_t i=0; i<kNModule*kNChip; i++){
89 fFound[i]=s.fFound[i];
90 fTried[i]=s.fTried[i];
91 }
92 if(fHis) {
93 InitHistos();
94 for(Int_t i=0; i<kNHisto; i++) {
95 s.fHisResX[i]->Copy(*fHisResX[i]);
96 s.fHisResZ[i]->Copy(*fHisResZ[i]);
97 s.fHisResXZ[i]->Copy(*fHisResXZ[i]);
98 s.fHisClusterSize[i]->Copy(*fHisClusterSize[i]);
99 for(Int_t clu=0; clu<kNclu; clu++) { // clu=0 --> cluster size 1
100 s.fHisResXclu[i][clu]->Copy(*fHisResXclu[i][clu]);
101 s.fHisResZclu[i][clu]->Copy(*fHisResZclu[i][clu]);
102 }
103 }
104 }
4a66240a 105}
106//_________________________________________________________________________
107AliITSPlaneEffSPD& AliITSPlaneEffSPD::operator+=(const AliITSPlaneEffSPD &add){
108 // Add-to-me operator
109 // Inputs:
110 // const AliITSPlaneEffSPD &add simulation class to be added
111 // Outputs:
112 // none.
113 // Return:
114 // none
115 for (UInt_t i=0; i<kNModule*kNChip; i++){
116 fFound[i] += add.fFound[i];
117 fTried[i] += add.fTried[i];
118 }
5fbd4fd6 119 if(fHis && add.fHis) {
120 for(Int_t i=0; i<kNHisto; i++) {
121 fHisResX[i]->Add(add.fHisResX[i]);
122 fHisResZ[i]->Add(add.fHisResZ[i]);
123 fHisResXZ[i]->Add(add.fHisResXZ[i]);
124 fHisClusterSize[i]->Add(add.fHisClusterSize[i]);
125 for(Int_t clu=0; clu<kNclu; clu++) { // clu=0 --> cluster size 1
126 fHisResXclu[i][clu]->Add(add.fHisResXclu[i][clu]);
127 fHisResZclu[i][clu]->Add(add.fHisResZclu[i][clu]);
128 }
129 }
130 }
4a66240a 131 return *this;
132}
133//______________________________________________________________________
134AliITSPlaneEffSPD& AliITSPlaneEffSPD::operator=(const
135 AliITSPlaneEffSPD &s){
136 // Assignment operator
137 // Inputs:
138 // AliITSPlaneEffSPD &s The original class for which
139 // this class is a copy of
140 // Outputs:
141 // none.
142 // Return:
143
144 if(this==&s) return *this;
145 s.Copy(*this);
4a66240a 146 return *this;
147}
148//______________________________________________________________________
149void AliITSPlaneEffSPD::Copy(TObject &obj) const {
150 // protected method. copy this to obj
151 AliITSPlaneEff::Copy(obj);
5fbd4fd6 152 AliITSPlaneEffSPD& target = (AliITSPlaneEffSPD &) obj;
4a66240a 153 for(Int_t i=0;i<kNModule*kNChip;i++) {
5fbd4fd6 154 target.fFound[i] = fFound[i];
155 target.fTried[i] = fTried[i];
156 }
157 CopyHistos(target);
158 return;
159}
160//_______________________________________________________________________
161void AliITSPlaneEffSPD::CopyHistos(AliITSPlaneEffSPD &target) const {
162 // protected method: copy histos from this to target
163 target.fHis = fHis; // this is redundant only in some cases. Leave as it is.
164 if(fHis) {
165 target.fHisResX=new TH1F*[kNHisto];
166 target.fHisResZ=new TH1F*[kNHisto];
167 target.fHisResXZ=new TH2F*[kNHisto];
168 target.fHisClusterSize=new TH2I*[kNHisto];
169 target.fHisResXclu=new TH1F**[kNHisto];
170 target.fHisResZclu=new TH1F**[kNHisto];
171 for(Int_t i=0; i<kNHisto; i++) {
172 target.fHisResX[i] = new TH1F(*fHisResX[i]);
173 target.fHisResZ[i] = new TH1F(*fHisResZ[i]);
174 target.fHisResXZ[i] = new TH2F(*fHisResXZ[i]);
175 target.fHisClusterSize[i] = new TH2I(*fHisClusterSize[i]);
176 target.fHisResXclu[i]=new TH1F*[kNclu];
177 target.fHisResZclu[i]=new TH1F*[kNclu];
178 for(Int_t clu=0; clu<kNclu; clu++) { // clu=0 --> cluster size 1
179 target.fHisResXclu[i][clu] = new TH1F(*fHisResXclu[i][clu]);
180 target.fHisResZclu[i][clu] = new TH1F(*fHisResZclu[i][clu]);
181 }
182 }
4a66240a 183 }
5fbd4fd6 184return;
4a66240a 185}
186//______________________________________________________________________
187AliITSPlaneEff& AliITSPlaneEffSPD::operator=(const
188 AliITSPlaneEff &s){
189 // Assignment operator
190 // Inputs:
191 // AliITSPlaneEffSPD &s The original class for which
192 // this class is a copy of
193 // Outputs:
194 // none.
195 // Return:
196
197 if(&s == this) return *this;
7167ae53 198 AliError("operator=: Not allowed to make a =, use default creater instead");
4a66240a 199 return *this;
200}
201//_______________________________________________________________________
202Int_t AliITSPlaneEffSPD::GetMissingTracksForGivenEff(Double_t eff, Double_t RelErr,
203 UInt_t im, UInt_t ic) const {
204
205 // Estimate the number of tracks still to be collected to attain a
206 // given efficiency eff, with relative error RelErr
207 // Inputs:
208 // eff -> Expected efficiency (e.g. those from actual estimate)
209 // RelErr -> tollerance [0,1]
210 // im -> module number [0,249]
6344adcc 211 // ic -> chip number [0,4]
4a66240a 212 // Outputs: none
213 // Return: the estimated n. of tracks
214 //
215if (im>=kNModule || ic>=kNChip)
7167ae53 216 {AliError("GetMissingTracksForGivenEff: you asked for a non existing chip");
4a66240a 217 return -1;}
218else return GetNTracksForGivenEff(eff,RelErr)-fTried[GetKey(im,ic)];
219}
220//_________________________________________________________________________
221Double_t AliITSPlaneEffSPD::PlaneEff(const UInt_t im,const UInt_t ic) const {
222// Compute the efficiency for a basic block,
223// Inputs:
6344adcc 224// im -> module number [0,249]
225// ic -> chip number [0,4]
4a66240a 226if (im>=kNModule || ic>=kNChip)
7167ae53 227 {AliError("PlaneEff(Uint_t,Uint_t): you asked for a non existing chip"); return -1.;}
4a66240a 228 Int_t nf=fFound[GetKey(im,ic)];
229 Int_t nt=fTried[GetKey(im,ic)];
230return AliITSPlaneEff::PlaneEff(nf,nt);
231}
232//_________________________________________________________________________
233Double_t AliITSPlaneEffSPD::ErrPlaneEff(const UInt_t im,const UInt_t ic) const {
234 // Compute the statistical error on efficiency for a basic block,
235 // using binomial statistics
236 // Inputs:
6344adcc 237 // im -> module number [0,249]
238 // ic -> chip number [0,4]
4a66240a 239if (im>=kNModule || ic>=kNChip)
7167ae53 240 {AliError("ErrPlaneEff(Uint_t,Uint_t): you asked for a non existing chip"); return -1.;}
4a66240a 241Int_t nf=fFound[GetKey(im,ic)];
242Int_t nt=fTried[GetKey(im,ic)];
243return AliITSPlaneEff::ErrPlaneEff(nf,nt);
244}
245//_________________________________________________________________________
246Bool_t AliITSPlaneEffSPD::UpDatePlaneEff(const Bool_t Kfound,
247 const UInt_t im, const UInt_t ic) {
248 // Update efficiency for a basic block
249if (im>=kNModule || ic>=kNChip)
7167ae53 250 {AliError("UpDatePlaneEff: you asked for a non existing chip"); return kFALSE;}
4a66240a 251 fTried[GetKey(im,ic)]++;
252 if(Kfound) fFound[GetKey(im,ic)]++;
253 return kTRUE;
254}
255//_________________________________________________________________________
6344adcc 256UInt_t AliITSPlaneEffSPD::GetChipFromCol(const UInt_t col) const {
4a66240a 257 // get chip given the column
258if(col>=kNCol*kNChip)
7167ae53 259 {AliDebug(1,Form("GetChipFromCol: you asked for a non existing column %d",col)); return 10;}
4a66240a 260return col/kNCol;
261}
262//__________________________________________________________________________
263UInt_t AliITSPlaneEffSPD::GetKey(const UInt_t mod, const UInt_t chip) const {
264 // get key given a basic block
265if(mod>=kNModule || chip>=kNChip)
7167ae53 266 {AliWarning("GetKey: you asked for a non existing block"); return 99999;}
4a66240a 267return mod*kNChip+chip;
268}
269//__________________________________________________________________________
270UInt_t AliITSPlaneEffSPD::GetModFromKey(const UInt_t key) const {
271 // get mod. from key
272if(key>=kNModule*kNChip)
7167ae53 273 {AliError("GetModFromKey: you asked for a non existing key"); return 9999;}
4a66240a 274return key/kNChip;
275}
276//__________________________________________________________________________
277UInt_t AliITSPlaneEffSPD::GetChipFromKey(const UInt_t key) const {
278 // retrieves chip from key
279if(key>=kNModule*kNChip)
7167ae53 280 {AliError("GetChipFromKey: you asked for a non existing key"); return 999;}
4a66240a 281return (key%(kNModule*kNChip))%kNChip;
282}
283//__________________________________________________________________________
284void AliITSPlaneEffSPD::GetModAndChipFromKey(const UInt_t key,UInt_t& mod,UInt_t& chip) const {
285 // get module and chip from a key
286if(key>=kNModule*kNChip)
7167ae53 287 {AliError("GetModAndChipFromKey: you asked for a non existing key");
4a66240a 288 mod=9999;
289 chip=999;
290 return;}
291mod=key/kNChip;
292chip=(key%(kNModule*kNChip))%kNChip;
293return;
294}
295//____________________________________________________________________________
296Double_t AliITSPlaneEffSPD::LivePlaneEff(UInt_t key) const {
6344adcc 297 // returns plane efficieny after adding the fraction of sensor which is bad
4a66240a 298if(key>=kNModule*kNChip)
7167ae53 299 {AliError("LivePlaneEff: you asked for a non existing key");
4a66240a 300 return -1.;}
6344adcc 301Double_t leff=AliITSPlaneEff::LivePlaneEff(0); // this just for the Warning
302leff=PlaneEff(key)+GetFracBad(key);
303return leff>1?1:leff;
4a66240a 304}
305//____________________________________________________________________________
306Double_t AliITSPlaneEffSPD::ErrLivePlaneEff(UInt_t key) const {
307 // returns error on live plane efficiency
308if(key>=kNModule*kNChip)
7167ae53 309 {AliError("ErrLivePlaneEff: you asked for a non existing key");
4a66240a 310 return -1.;}
6344adcc 311Int_t nf=fFound[key];
312Double_t triedInLive=GetFracLive(key)*fTried[key];
313Int_t nt=TMath::Max(nf,TMath::Nint(triedInLive));
314return AliITSPlaneEff::ErrPlaneEff(nf,nt); // for the time being: to be checked
4a66240a 315}
316//_____________________________________________________________________________
317Double_t AliITSPlaneEffSPD::GetFracLive(const UInt_t key) const {
318 // returns the fraction of the sensor which is OK
319if(key>=kNModule*kNChip)
7167ae53 320 {AliError("GetFracLive: you asked for a non existing key");
4a66240a 321 return -1.;}
322 // Compute the fraction of bad (dead+noisy) detector
323UInt_t dead=0,noisy=0;
324GetDeadAndNoisyInChip(key,dead,noisy);
325Double_t live=dead+noisy;
326live/=(kNRow*kNCol);
327return 1.-live;
328}
329//_____________________________________________________________________________
330void AliITSPlaneEffSPD::GetDeadAndNoisyInChip(const UInt_t key,
331 UInt_t& nrDeadInChip, UInt_t& nrNoisyInChip) const {
332 // returns the number of dead and noisy pixels
333nrDeadInChip=0;
334nrNoisyInChip=0;
335if(key>=kNModule*kNChip)
7167ae53 336 {AliError("GetDeadAndNoisyInChip: you asked for a non existing key");
4a66240a 337 return;}
338 // Compute the number of bad (dead+noisy) pixel in a chip
339//
340if(!fInitCDBCalled)
7167ae53 341 {AliError("GetDeadAndNoisyInChip: CDB not inizialized: call InitCDB first");
4a66240a 342 return;};
343AliCDBManager* man = AliCDBManager::Instance();
344// retrieve map of dead Pixel
345AliCDBEntry *cdbSPDDead = man->Get("ITS/Calib/SPDDead", fRunNumber);
346TObjArray* spdDead;
347if(cdbSPDDead) {
348 spdDead = (TObjArray*)cdbSPDDead->GetObject();
349 if(!spdDead)
7167ae53 350 {AliError("GetDeadAndNoisyInChip: SPDDead not found in CDB");
4a66240a 351 return;}
352} else {
7167ae53 353 AliError("GetDeadAndNoisyInChip: did not find Calib/SPDDead.");
4a66240a 354 return;
355}
356// retrieve map of noisy Pixel
357AliCDBEntry *cdbSPDNoisy = man->Get("ITS/Calib/SPDNoisy", fRunNumber);
358TObjArray* spdNoisy;
359if(cdbSPDNoisy) {
360 spdNoisy = (TObjArray*)cdbSPDNoisy->GetObject();
361 if(!spdNoisy)
7167ae53 362 {AliError("GetDeadAndNoisyInChip: SPDNoisy not found in CDB");
4a66240a 363 return;}
364} else {
7167ae53 365 AliError("GetDeadAndNoisyInChip: did not find Calib/SPDNoisy.");
4a66240a 366 return;
367}
368//
369UInt_t mod=GetModFromKey(key);
370UInt_t chip=GetChipFromKey(key);
371// count number of dead
372AliITSCalibrationSPD* calibSPD=(AliITSCalibrationSPD*) spdDead->At(mod);
373UInt_t nrDead = calibSPD->GetNrBad();
374for (UInt_t index=0; index<nrDead; index++) {
6344adcc 375 if(GetChipFromCol(calibSPD->GetBadColAt(index))==chip) nrDeadInChip++;
4a66240a 376}
377calibSPD=(AliITSCalibrationSPD*) spdNoisy->At(mod);
378UInt_t nrNoisy = calibSPD->GetNrBad();
379for (UInt_t index=0; index<nrNoisy; index++) {
6344adcc 380 if(GetChipFromCol(calibSPD->GetBadColAt(index))==chip) nrNoisyInChip++;
4a66240a 381}
382return;
383}
384//_____________________________________________________________________________
385Double_t AliITSPlaneEffSPD::GetFracBad(const UInt_t key) const {
386 // returns 1-fractional live
387if(key>=kNModule*kNChip)
7167ae53 388 {AliError("GetFracBad: you asked for a non existing key");
4a66240a 389 return -1.;}
390return 1.-GetFracLive(key);
391}
392//_____________________________________________________________________________
393Bool_t AliITSPlaneEffSPD::WriteIntoCDB() const {
394// write onto CDB
395if(!fInitCDBCalled)
7167ae53 396 {AliError("WriteIntoCDB: CDB not inizialized. Call InitCDB first");
4a66240a 397 return kFALSE;}
398// to be written properly: now only for debugging
399 AliCDBMetaData *md= new AliCDBMetaData(); // metaData describing the object
400 md->SetObjectClassName("AliITSPlaneEff");
401 md->SetResponsible("Giuseppe Eugenio Bruno");
402 md->SetBeamPeriod(0);
403 md->SetAliRootVersion("head 19/11/07"); //root version
404 AliCDBId id("ITS/PlaneEff/PlaneEffSPD",0,AliCDBRunRange::Infinity());
405 AliITSPlaneEffSPD eff;
406 eff=*this;
407 Bool_t r=AliCDBManager::Instance()->GetDefaultStorage()->Put(&eff,id,md);
408 delete md;
409 return r;
410}
411//_____________________________________________________________________________
412Bool_t AliITSPlaneEffSPD::ReadFromCDB() {
413// read from CDB
414if(!fInitCDBCalled)
7167ae53 415 {AliError("ReadFromCDB: CDB not inizialized. Call InitCDB first");
4a66240a 416 return kFALSE;}
417//if(!AliCDBManager::Instance()->IsDefaultStorageSet()) {
418// AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT");
419// }
420AliCDBEntry *cdbEntry = AliCDBManager::Instance()->Get("ITS/PlaneEff/PlaneEffSPD",fRunNumber);
421AliITSPlaneEffSPD* eff= (AliITSPlaneEffSPD*)cdbEntry->GetObject();
422if(this==eff) return kFALSE;
5fbd4fd6 423if(fHis) CopyHistos(*eff); // If histos already exist then copy them to eff
424eff->Copy(*this); // copy everything (statistics and histos) from eff to this
4a66240a 425return kTRUE;
426}
7167ae53 427//_____________________________________________________________________________
428UInt_t AliITSPlaneEffSPD::GetKeyFromDetLocCoord(Int_t ilay, Int_t idet,
429 Float_t, Float_t locz) const {
430// method to locate a basic block from Detector Local coordinate (to be used in tracking)
431UInt_t key=999999;
432if(ilay<0 || ilay>1)
433 {AliError("GetKeyFromDetLocCoord: you asked for a non existing layer");
434 return key;}
435if(ilay==0 && (idet<0 || idet>79))
436 {AliError("GetKeyFromDetLocCoord: you asked for a non existing detector");
437 return key;}
438if(ilay==1 && (idet<0 || idet>159))
439 {AliError("GetKeyFromDetLocCoord: you asked for a non existing detector");
440 return key;}
441
442UInt_t mod=idet;
443if(ilay==1) mod+=80;
444key=GetKey(mod,GetChipFromCol(GetColFromLocZ(locz)));
445return key;
446}
447//_____________________________________________________________________________
448UInt_t AliITSPlaneEffSPD::GetColFromLocZ(Float_t zloc) const {
449UInt_t col=0;
450/* note: as it is now, the AliITSsegmentationSPD::Init() does not properly initialize (6 chips !!!)
451AliITSsegmentationSPD spd;
452spd.Init();
453Int_t ix,iz;
454if(spd.LocalToDet(0,zloc,ix,iz)) col+=iz;
455else {
456 AliError("GetColFromLocZ: cannot compute column number from local z");
457 col=99999;}
458return col;
459*/
460const Float_t kconv = 1.0E-04; // converts microns to cm.
461Float_t bz[160];
462for(Int_t i=000;i<160;i++) bz[i] = 425.0; // most are 425 microns except below
463bz[ 31] = bz[ 32] = 625.0; // first chip boundry
464bz[ 63] = bz[ 64] = 625.0; // first chip boundry
465bz[ 95] = bz[ 96] = 625.0; // first chip boundry
466bz[127] = bz[128] = 625.0; // first chip boundry
467//
468Int_t j=-1;
469Float_t dz=0;
470for(Int_t i=000;i<160;i++) dz+=bz[i];
471dz = -0.5*kconv*dz;
472if(zloc<dz || zloc>-1*dz) { // outside z range
473 AliDebug(1,Form("GetColFromLocZ: cannot compute column number from local z=%f",zloc));
474 return 99999;}
475for(j=0;j<160;j++){
476 dz += kconv*bz[j];
477 if(zloc<dz) break;
478} // end for j
479col+=j;
480//
481return col;
482}
483//________________________________________________________
aa0de373 484Bool_t AliITSPlaneEffSPD::GetBlockBoundaries(const UInt_t key, Float_t& xmn,Float_t& xmx,
485 Float_t& zmn,Float_t& zmx) const {
486//
487// This method return the geometrical boundaries of the active volume of a given
488// basic block, in the detector reference system.
489// Input: unique key to locate a basic block.
490//
491// Output: Ymin, Ymax, Zmin, Zmax of a basic block (chip for SPD)
492// Return: kTRUE if computation was succesfully, kFALSE otherwise
493//
494if(key>=kNModule*kNChip)
495 {AliWarning("GetBlockBoundaries: you asked for a non existing key"); return kFALSE;}
496UInt_t chip=GetChipFromKey(key);
497zmn=GetLocZFromCol(chip*kNCol);
498zmx=GetLocZFromCol((chip+1)*kNCol);
499xmn=GetLocXFromRow(0);
500xmx=GetLocXFromRow(kNRow);
501Float_t tmp=zmn;
502if(zmx<zmn) {zmn=zmx; zmx=tmp;}
503tmp=xmn;
504if(xmx<xmn) {xmn=xmx; xmx=tmp;}
505return kTRUE;
506}
507//________________________________________________________
508Float_t AliITSPlaneEffSPD::GetLocXFromRow(const UInt_t row) const {
509//
510// This method return the local (i.e. detector reference system) lower x coordinate
511// of the row. To get the central value of a given row, you can do
512// 1/2*[LocXFromRow(row)+LocXFromRow(row+1)].
513//
514// Input: row number in the range [0,kNRow]
515// Output: lower local X coordinate of this row.
516//
517if(row>kNRow) // not >= ! allow also computation of upper limit of the last row.
518 {AliError("LocYFromRow: you asked for a non existing row"); return 9999999.;}
519const Float_t kconv = 1.0E-04; // converts microns to cm.
520Float_t bx[256];
521for(Int_t i=000;i<256;i++) bx[i] = 50.0; // in x all are 50 microns.
522//
523Float_t dx=0;
524for(Int_t i=000;i<256;i++) dx+=bx[i];
525dx = -0.5*kconv*dx;
526for(UInt_t j=0;j<row;j++){
527 dx += kconv*bx[j];
528} // end for j
529return dx;
530}
531//________________________________________________________
532Float_t AliITSPlaneEffSPD::GetLocZFromCol(const UInt_t col) const {
533//
534// This method return the local (i.e. detector reference system) lower Z coordinate
535// of the column. To get the central value of a given column, you can do
536// 1/2*[LocZFromCol(col)+LocZFromCol(col+1)].
537//
538// Input: col number in the range [0,kNChip*kNCol]
539// Output: lower local Y coordinate of this row.
540//
541if(col>kNChip*kNCol) // not >= ! allow also computation of upper limit of the last column
542 {AliError("LocZFromCol: you asked for a non existing column"); return 9999999.;}
543const Float_t kconv = 1.0E-04; // converts microns to cm.
544Float_t bz[160];
545for(Int_t i=000;i<160;i++) bz[i] = 425.0; // most are 425 microns except below
546bz[ 31] = bz[ 32] = 625.0; // first chip boundry
547bz[ 63] = bz[ 64] = 625.0; // first chip boundry
548bz[ 95] = bz[ 96] = 625.0; // first chip boundry
549bz[127] = bz[128] = 625.0; // first chip boundry
550//
551Float_t dz=0;
552for(Int_t i=000;i<160;i++) dz+=bz[i];
553dz = -0.5*kconv*dz;
554for(UInt_t j=0;j<col;j++){
555 dz += kconv*bz[j];
556} // end for j
557return dz;
558}
5fbd4fd6 559//__________________________________________________________
560void AliITSPlaneEffSPD::InitHistos() {
561 // for the moment let's create the histograms
562 // module by module
563 TString histnameResX="HistResX_mod_",aux;
564 TString histnameResZ="HistResZ_mod_";
565 TString histnameResXZ="HistResXZ_mod_";
566 TString histnameClusterType="HistClusterType_mod_";
567 TString histnameResXclu="HistResX_mod_";
568 TString histnameResZclu="HistResZ_mod_";
569//
570 fHisResX=new TH1F*[kNHisto];
571 fHisResZ=new TH1F*[kNHisto];
572 fHisResXZ=new TH2F*[kNHisto];
573 fHisClusterSize=new TH2I*[kNHisto];
574 fHisResXclu=new TH1F**[kNHisto];
575 fHisResZclu=new TH1F**[kNHisto];
576
577 for (Int_t nhist=0;nhist<kNHisto;nhist++){
578 aux=histnameResX;
579 aux+=nhist;
580 fHisResX[nhist]=new TH1F("histname","histname",800,-0.04,0.04); // +- 400 micron; 1 bin=1 micron
581 fHisResX[nhist]->SetName(aux.Data());
582 fHisResX[nhist]->SetTitle(aux.Data());
583
584 aux=histnameResZ;
585 aux+=nhist;
586 fHisResZ[nhist]=new TH1F("histname","histname",800,-0.32,0.32); // +-3200 micron; 1 bin=8 micron
587 fHisResZ[nhist]->SetName(aux.Data());
588 fHisResZ[nhist]->SetTitle(aux.Data());
589
590 aux=histnameResXZ;
591 aux+=nhist;
592 fHisResXZ[nhist]=new TH2F("histname","histname",40,-0.02,0.02,40,-0.16,0.16); // binning:
593 // 10 micron in x;
594 // 80 micron in z;
595 fHisResXZ[nhist]->SetName(aux.Data());
596 fHisResXZ[nhist]->SetTitle(aux.Data());
597
598 aux=histnameClusterType;
599 aux+=nhist;
600 fHisClusterSize[nhist]=new TH2I("histname","histname",10,0.5,10.5,10,0.5,10.5);
601 fHisClusterSize[nhist]->SetName(aux.Data());
602 fHisClusterSize[nhist]->SetTitle(aux.Data());
603
604 fHisResXclu[nhist]=new TH1F*[kNclu];
605 fHisResZclu[nhist]=new TH1F*[kNclu];
606 for(Int_t clu=0; clu<kNclu; clu++) { // clu=0 --> cluster size 1
607 aux=histnameResXclu;
608 aux+=nhist;
609 aux+="_clu_";
610 aux+=clu+1; // clu=0 --> cluster size 1
611 fHisResXclu[nhist][clu]=new TH1F("histname","histname",800,-0.04,0.04); // +- 400 micron; 1 bin=1 micron
612 fHisResXclu[nhist][clu]->SetName(aux.Data());
613 fHisResXclu[nhist][clu]->SetTitle(aux.Data());
614
615 aux=histnameResZclu;
616 aux+=nhist;
617 aux+="_clu_";
618 aux+=clu+1; // clu=0 --> cluster size 1
619 fHisResZclu[nhist][clu]=new TH1F("histname","histname",800,-0.32,0.32); // +-3200 micron; 1 bin=8 micron
620 fHisResZclu[nhist][clu]->SetName(aux.Data());
621 fHisResZclu[nhist][clu]->SetTitle(aux.Data());
622 }
623
624 }
625return;
626}
627//__________________________________________________________
628void AliITSPlaneEffSPD::DeleteHistos() {
629 if(fHisResX) {
630 for (Int_t i=0; i<kNHisto; i++ ) delete fHisResX[i];
631 delete [] fHisResX;
632 }
633 if(fHisResZ) {
634 for (Int_t i=0; i<kNHisto; i++ ) delete fHisResZ[i];
635 delete [] fHisResZ;
636 }
637 if(fHisResXZ) {
638 for (Int_t i=0; i<kNHisto; i++ ) delete fHisResXZ[i];
639 delete [] fHisResXZ;
640 }
641 if(fHisClusterSize) {
642 for (Int_t i=0; i<kNHisto; i++ ) delete fHisClusterSize[i];
643 delete [] fHisClusterSize;
644 }
645 if(fHisResXclu) {
646 for (Int_t i=0; i<kNHisto; i++ ) {
647 for (Int_t clu=0; clu<kNclu; clu++) if (fHisResXclu[i][clu]) delete fHisResXclu[i][clu];
648 delete [] fHisResXclu[i];
649 }
650 delete [] fHisResXclu;
651 fHisResXclu = 0;
652 }
653 if(fHisResZclu) {
654 for (Int_t i=0; i<kNHisto; i++ ) {
655 for (Int_t clu=0; clu<kNclu; clu++) if (fHisResZclu[i][clu]) delete fHisResZclu[i][clu];
656 delete [] fHisResZclu[i];
657 }
658 delete [] fHisResZclu;
659 fHisResZclu = 0;
660 }
661
662return;
663}
664//__________________________________________________________
665Bool_t AliITSPlaneEffSPD::FillHistos(UInt_t key, Bool_t found,
666 Float_t tXZ[2], Float_t cXZ[2], Int_t ctXZ[2]) {
667// this method fill the histograms
668// input: - key: unique key of the basic block
669// - found: Boolean to asses whether a cluster has been associated to the track or not
670// - tXZ[2] local X and Z coordinates of the track prediction
671// - cXZ[2] local X and Z coordinates of the cluster associated to the track
672// output: kTRUE if filling was succesfull kFALSE otherwise
673// side effects: updating of the histograms.
674//
675 if (!fHis) {
676 AliWarning("FillHistos: histograms do not exist! Call SetCreateHistos(kTRUE) first");
677 return kFALSE;
678 }
679 if(key>=kNModule*kNChip)
680 {AliWarning("FillHistos: you asked for a non existing key"); return kFALSE;}
681 Int_t id=GetModFromKey(key);
682 if(id>=kNHisto)
683 {AliWarning("FillHistos: you want to fill a non-existing histos"); return kFALSE;}
684 if(found) {
685 Float_t resx=tXZ[0]-cXZ[0];
686 Float_t resz=tXZ[1]-cXZ[1];
687 fHisResX[id]->Fill(resx);
688 fHisResZ[id]->Fill(resz);
689 fHisResXZ[id]->Fill(resx,resz);
690 fHisClusterSize[id]->Fill((Double_t)ctXZ[0],(Double_t)ctXZ[1]);
691 if(ctXZ[0]>0 && ctXZ[0]<=kNclu) fHisResXclu[id][ctXZ[0]-1]->Fill(resx);
692 if(ctXZ[1]>0 && ctXZ[1]<=kNclu) fHisResZclu[id][ctXZ[1]-1]->Fill(resx);
693 }
694 return kTRUE;
695}
696//__________________________________________________________
697Bool_t AliITSPlaneEffSPD::WriteHistosToFile(TString filename, Option_t* option) {
698 //
699 // Saves the histograms into a tree and saves the trees into a file
700 //
701 if (!fHis) return kFALSE;
702 if (filename.Data()=="") {
703 AliWarning("WriteHistosToFile: null output filename!");
704 return kFALSE;
705 }
706 char branchname[30];
707 TFile *hFile=new TFile(filename.Data(),option,
708 "The File containing the TREEs with ITS PlaneEff Histos");
709 TTree *SPDTree=new TTree("SPDTree","Tree whith Residuals and Cluster Type distributions for SPD");
710 TH1F *histZ,*histX;
711 TH2F *histXZ;
712 TH2I *histClusterType;
713 TH1F *histXclu[kNclu];
714 TH1F *histZclu[kNclu];
715
716 histZ=new TH1F();
717 histX=new TH1F();
718 histXZ=new TH2F();
719 histClusterType=new TH2I();
720 for(Int_t clu=0;clu<kNclu;clu++) {
721 histXclu[clu]=new TH1F();
722 histZclu[clu]=new TH1F();
723 }
724
725 SPDTree->Branch("histX","TH1F",&histX,128000,0);
726 SPDTree->Branch("histZ","TH1F",&histZ,128000,0);
727 SPDTree->Branch("histXZ","TH2F",&histXZ,128000,0);
728 SPDTree->Branch("histClusterType","TH2I",&histClusterType,128000,0);
729 for(Int_t clu=0;clu<kNclu;clu++) {
730 sprintf(branchname,"histXclu_%d",clu+1);
731 SPDTree->Branch(branchname,"TH1F",&histXclu[clu],128000,0);
732 sprintf(branchname,"histZclu_%d",clu+1);
733 SPDTree->Branch(branchname,"TH1F",&histZclu[clu],128000,0);
734 }
735
736 for(Int_t j=0;j<kNHisto;j++){
737 histX=fHisResX[j];
738 histZ=fHisResZ[j];
739 histXZ=fHisResXZ[j];
740 histClusterType=fHisClusterSize[j];
741 for(Int_t clu=0;clu<kNclu;clu++) {
742 histXclu[clu]=fHisResXclu[j][clu];
743 histZclu[clu]=fHisResZclu[j][clu];
744 }
745 SPDTree->Fill();
746 }
747 hFile->Write();
748 hFile->Close();
749return kTRUE;
750}
751//__________________________________________________________
752Bool_t AliITSPlaneEffSPD::ReadHistosFromFile(TString filename) {
753 //
754 // Read histograms from an already existing file
755 //
756 if (!fHis) return kFALSE;
757 if (filename.Data()=="") {
758 AliWarning("ReadHistosFromFile: incorrect output filename!");
759 return kFALSE;
760 }
761 char branchname[30];
762
763 TH1F *h = 0;
764 TH2F *h2 = 0;
765 TH2I *h2i= 0;
766
767 TFile *file=TFile::Open(filename.Data(),"READONLY");
768
769 if (!file || file->IsZombie()) {
770 AliWarning(Form("Can't open %s !",filename.Data()));
771 delete file;
772 return kFALSE;
773 }
774 TTree *tree = (TTree*) file->Get("SPDTree");
775
776 TBranch *histX = (TBranch*) tree->GetBranch("histX");
777 TBranch *histZ = (TBranch*) tree->GetBranch("histZ");
778 TBranch *histXZ = (TBranch*) tree->GetBranch("histXZ");
779 TBranch *histClusterType = (TBranch*) tree->GetBranch("histClusterType");
780
781 TBranch *histXclu[kNclu], *histZclu[kNclu];
782 for(Int_t clu=0; clu<kNclu; clu++) {
783 sprintf(branchname,"histXclu_%d",clu+1);
784 histXclu[clu]= (TBranch*) tree->GetBranch(branchname);
785 sprintf(branchname,"histZclu_%d",clu+1);
786 histZclu[clu]= (TBranch*) tree->GetBranch(branchname);
787 }
788
789 gROOT->cd();
790
791 Int_t nevent = (Int_t)histX->GetEntries();
792 if(nevent!=kNHisto)
793 {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
794 histX->SetAddress(&h);
795 for(Int_t j=0;j<kNHisto;j++){
796 delete h; h=0;
797 histX->GetEntry(j);
798 fHisResX[j]->Add(h);
799 }
800
801 nevent = (Int_t)histZ->GetEntries();
802 if(nevent!=kNHisto)
803 {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
804 histZ->SetAddress(&h);
805 for(Int_t j=0;j<kNHisto;j++){
806 delete h; h=0;
807 histZ->GetEntry(j);
808 fHisResZ[j]->Add(h);
809 }
810
811 nevent = (Int_t)histXZ->GetEntries();
812 if(nevent!=kNHisto)
813 {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
814 histXZ->SetAddress(&h2);
815 for(Int_t j=0;j<kNHisto;j++){
816 delete h2; h2=0;
817 histXZ->GetEntry(j);
818 fHisResXZ[j]->Add(h2);
819 }
820
821 nevent = (Int_t)histClusterType->GetEntries();
822 if(nevent!=kNHisto)
823 {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
824 histClusterType->SetAddress(&h2i);
825 for(Int_t j=0;j<kNHisto;j++){
826 delete h2i; h2i=0;
827 histClusterType->GetEntry(j);
828 fHisClusterSize[j]->Add(h2i);
829 }
830
831 for(Int_t clu=0; clu<kNclu; clu++) {
832
833 nevent = (Int_t)histXclu[clu]->GetEntries();
834 if(nevent!=kNHisto)
835 {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
836 histXclu[clu]->SetAddress(&h);
837 for(Int_t j=0;j<kNHisto;j++){
838 delete h; h=0;
839 histXclu[clu]->GetEntry(j);
840 fHisResXclu[j][clu]->Add(h);
841 }
842
843 nevent = (Int_t)histZclu[clu]->GetEntries();
844 if(nevent!=kNHisto)
845 {AliWarning("ReadHistosFromFile: trying to read too many or too few histos!"); return kFALSE;}
846 histZclu[clu]->SetAddress(&h);
847 for(Int_t j=0;j<kNHisto;j++){
848 delete h; h=0;
849 histZclu[clu]->GetEntry(j);
850 fHisResZclu[j][clu]->Add(h);
851 }
852 }
853
854 delete h; h=0;
855 delete h2; h2=0;
856 delete h2i; h2i=0;
857
858 if (file) {
859 file->Close();
860 }
861
862 fHisResZclu[0][0]->Draw();
863return kTRUE;
864}