Line data Source code
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 : /* $Id$ */
17 :
18 : ///////////////////////////////////////////////////////////////////
19 : // //
20 : // Implementation of the base class for SDD map corrections //
21 : // Origin: F.Prino, Torino, prino@to.infn.it //
22 : // //
23 : ///////////////////////////////////////////////////////////////////
24 :
25 : #include "TH1F.h"
26 : #include "TH2F.h"
27 : #include "AliITSCorrMapSDD.h"
28 : #include "AliITSsegmentationSDD.h"
29 :
30 : const Int_t AliITSCorrMapSDD::fgkNAnodePtsDefault = 1;
31 : const Int_t AliITSCorrMapSDD::fgkNDriftPtsDefault = 72;
32 :
33 118 : ClassImp(AliITSCorrMapSDD)
34 : //______________________________________________________________________
35 1560 : AliITSCorrMapSDD::AliITSCorrMapSDD():TNamed("defaultmap",""),
36 1560 : fNAnodePts(fgkNAnodePtsDefault),
37 1560 : fNDriftPts(fgkNDriftPtsDefault),
38 1560 : fXt1(0.),
39 1560 : fXt2(0.),
40 1560 : fXm1(0.),
41 1560 : fXm2(0.),
42 1560 : fDrLen(0.)
43 4680 : {
44 : // default constructor
45 1560 : }
46 : //______________________________________________________________________
47 : AliITSCorrMapSDD::AliITSCorrMapSDD(Char_t *mapname):
48 0 : TNamed(mapname,""),
49 0 : fNAnodePts(fgkNAnodePtsDefault),
50 0 : fNDriftPts(fgkNDriftPtsDefault),
51 0 : fXt1(0.),
52 0 : fXt2(0.),
53 0 : fXm1(0.),
54 0 : fXm2(0.),
55 0 : fDrLen(0.)
56 0 : {
57 : // standard constructor
58 0 : }
59 : //______________________________________________________________________
60 : void AliITSCorrMapSDD::ComputeGridPoints(Float_t z, Float_t x, AliITSsegmentationSDD *seg, Bool_t isReco){
61 : // extracts points from the discrete grid with the correction map
62 :
63 : const Double_t kMicronTocm = 1.0e-4;
64 518 : Int_t nAnodes=seg->Npz();
65 518 : Int_t nAnodesHybrid=seg->NpzHalf();
66 518 : Int_t bina =(Int_t) seg->GetAnodeFromLocal(x,z);
67 518 : if(bina>nAnodes) AliError("Wrong anode anumber!");
68 776 : if(bina>=nAnodesHybrid) bina-=nAnodesHybrid;
69 518 : Float_t stept = seg->Dx()*kMicronTocm/(Float_t)fNDriftPts;
70 518 : fDrLen= seg->Dx()*kMicronTocm-TMath::Abs(x);
71 518 : if(fDrLen<0) fDrLen=0;
72 518 : Int_t bint = TMath::Abs((Int_t)(fDrLen/stept));
73 518 : if(bint==fNDriftPts) bint-=1;
74 518 : if(bint>=fNDriftPts){
75 0 : AliError("Wrong bin number along drift direction!");
76 0 : bint=fNDriftPts-1;
77 0 : }
78 518 : fXt1=stept*bint;
79 518 : fXm1=fXt1-GetCellContent(bina,bint)*kMicronTocm;
80 518 : if((bint+1)<fNDriftPts){
81 518 : fXt2=stept*(bint+1);
82 518 : fXm2=fXt2-GetCellContent(bina,bint+1)*kMicronTocm;
83 518 : }else{
84 0 : fXt2=stept*(bint-1);
85 0 : fXm2=fXt2-GetCellContent(bina,bint-1)*kMicronTocm;
86 : }
87 518 : if(isReco){
88 678 : if(fXm1<fDrLen && fXm2>fDrLen) return;
89 12 : if(bint==0 || bint==(fNDriftPts-1)) return;
90 6 : if(fXm1>fDrLen){
91 0 : for(Int_t itry=1; itry<=10; itry++){
92 0 : Float_t xmtest=(bint-itry)*stept-GetCellContent(bina,bint-itry)*kMicronTocm;
93 0 : if(xmtest<fDrLen){
94 0 : fXt1=stept*(bint-itry);
95 0 : fXt2=fXt1+stept;
96 0 : fXm1=fXt1-GetCellContent(bina,bint-itry)*kMicronTocm;
97 0 : fXm2=fXt2-GetCellContent(bina,bint+1-itry)*kMicronTocm;
98 0 : return;
99 : }
100 0 : }
101 : }
102 6 : if(fXm2<fDrLen){
103 12 : for(Int_t itry=1; itry<=10; itry++){
104 6 : Float_t xmtest=(bint+1+itry)*stept-GetCellContent(bina,bint+1+itry)*kMicronTocm;
105 6 : if(xmtest>fDrLen){
106 6 : fXt1=stept*(bint+itry);
107 6 : fXt2=fXt1+stept;
108 6 : fXm1=fXt1-GetCellContent(bina,bint+itry)*kMicronTocm;
109 6 : fXm2=fXt2-GetCellContent(bina,bint+1+itry)*kMicronTocm;
110 6 : return;
111 : }
112 0 : }
113 : }
114 : }
115 808 : }
116 : //______________________________________________________________________
117 : Float_t AliITSCorrMapSDD::GetCorrection(Float_t z, Float_t x, AliITSsegmentationSDD *seg){
118 : // returns correction in cm starting from local coordinates on the module
119 456 : ComputeGridPoints(z,x,seg,kTRUE);
120 228 : Float_t m=(fXt2-fXt1)/(fXm2-fXm1);
121 228 : Float_t q=fXt1-m*fXm1;
122 228 : Float_t xcorr=m*fDrLen+q;
123 : // fDrLen is the measured drift distance, xcorr is the corresponding true
124 684 : return GetInversionBit() ? fDrLen-xcorr : xcorr-fDrLen;
125 : }
126 : //______________________________________________________________________
127 : Float_t AliITSCorrMapSDD::GetShiftForSimulation(Float_t z, Float_t x, AliITSsegmentationSDD *seg){
128 : // returns shift to be appiled in digitizarion (in cm) starting from local coordinates on the module
129 580 : ComputeGridPoints(z,x,seg,kFALSE);
130 290 : Float_t m=(fXm2-fXm1)/(fXt2-fXt1);
131 290 : Float_t q=fXm1-m*fXt1;
132 290 : Float_t xshifted=m*fDrLen+q;
133 : // fDrLen is the true drift distance, xshifted is the one with map shift
134 870 : return GetInversionBit() ? xshifted-fDrLen : fDrLen-xshifted;
135 : }
136 : //______________________________________________________________________
137 : TH2F* AliITSCorrMapSDD::GetMapHisto() const{
138 : // Returns a TH2F histogram with map of residuals
139 0 : TString hname;
140 0 : hname.Form("h%s",GetName());
141 0 : TH2F* hmap=new TH2F(hname.Data(),"",fNAnodePts,-0.5,255.5,fNDriftPts,0.,35.);
142 0 : for(Int_t iAn=0;iAn<fNAnodePts; iAn++){
143 0 : for(Int_t iDr=0;iDr<fNDriftPts; iDr++){
144 0 : hmap->SetBinContent(iAn+1,iDr+1,GetCellContent(iAn,iDr));
145 : }
146 : }
147 : return hmap;
148 0 : }
149 : //______________________________________________________________________
150 : TH1F* AliITSCorrMapSDD::GetMapProfile() const{
151 : // Returns a TH1F with the projection of the map along drift coordinate
152 0 : TString hname;
153 0 : hname.Form("p%s",GetName());
154 0 : TH1F* hprof=new TH1F(hname.Data(),"",fNDriftPts,0.,35.);
155 0 : for(Int_t iDr=0;iDr<fNDriftPts; iDr++){
156 : Float_t meanval=0.;
157 0 : for(Int_t iAn=0;iAn<fNAnodePts; iAn++){
158 0 : meanval+=GetCellContent(iAn,iDr);
159 : }
160 0 : hprof->SetBinContent(iDr+1,meanval/fNAnodePts);
161 : }
162 : return hprof;
163 :
164 0 : }
165 : //______________________________________________________________________
166 : TH1F* AliITSCorrMapSDD::GetResidualDistr(Float_t dmin, Float_t dmax) const{
167 : // Returns a TH1F histogram with distribution of residual
168 0 : TString hname;
169 0 : hname.Form("hd%s",GetName());
170 0 : TH1F* hd=new TH1F(hname.Data(),"",100,dmin,dmax);
171 0 : for(Int_t iAn=0;iAn<fNAnodePts; iAn++){
172 0 : for(Int_t iDr=0;iDr<fNDriftPts; iDr++){
173 0 : hd->Fill(GetCellContent(iAn,iDr));
174 : }
175 : }
176 : return hd;
177 0 : }
178 :
|