]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSCorrMapSDD.cxx
Coverity fix
[u/mrichter/AliRoot.git] / ITS / AliITSCorrMapSDD.cxx
1 /**************************************************************************
2  * Copyright(c) 2007-2009, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $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
29 const Int_t AliITSCorrMapSDD::fgkNAnodePtsDefault = 1;
30 const Int_t AliITSCorrMapSDD::fgkNDriftPtsDefault = 72;
31
32 ClassImp(AliITSCorrMapSDD)
33 //______________________________________________________________________
34 AliITSCorrMapSDD::AliITSCorrMapSDD():TNamed("defaultmap",""),
35   fNAnodePts(fgkNAnodePtsDefault),
36   fNDriftPts(fgkNDriftPtsDefault),
37   fXt1(0.),
38   fXt2(0.),
39   fXm1(0.),
40   fXm2(0.),
41   fDrLen(0.)
42 {
43   // default constructor  
44 }
45 //______________________________________________________________________
46 AliITSCorrMapSDD::AliITSCorrMapSDD(Char_t *mapname):
47   TNamed(mapname,""),
48   fNAnodePts(fgkNAnodePtsDefault),
49   fNDriftPts(fgkNDriftPtsDefault),
50   fXt1(0.),
51   fXt2(0.),
52   fXm1(0.),
53   fXm2(0.),
54   fDrLen(0.)
55 {
56   // standard constructor
57 }
58 //______________________________________________________________________
59 void AliITSCorrMapSDD::ComputeGridPoints(Float_t z, Float_t x, AliITSsegmentationSDD *seg, Bool_t isReco){
60   // extracts points from the discrete grid with the correction map
61
62   const Double_t kMicronTocm = 1.0e-4; 
63   Int_t nAnodes=seg->Npz();
64   Int_t nAnodesHybrid=seg->NpzHalf();
65   Int_t bina =(Int_t) seg->GetAnodeFromLocal(x,z);
66   if(bina>nAnodes)  AliError("Wrong anode anumber!");
67   if(bina>=nAnodesHybrid) bina-=nAnodesHybrid;
68   Float_t stept = seg->Dx()*kMicronTocm/(Float_t)fNDriftPts;
69   fDrLen= seg->Dx()*kMicronTocm-TMath::Abs(x);
70   if(fDrLen<0) fDrLen=0;
71   Int_t bint = TMath::Abs((Int_t)(fDrLen/stept));
72   if(bint==fNDriftPts) bint-=1;
73   if(bint>=fNDriftPts){
74     AliError("Wrong bin number along drift direction!");
75     bint=fNDriftPts-1;
76   }
77   fXt1=stept*bint;
78   fXm1=fXt1-GetCellContent(bina,bint)*kMicronTocm;
79   if((bint+1)<fNDriftPts){
80     fXt2=stept*(bint+1);
81     fXm2=fXt2-GetCellContent(bina,bint+1)*kMicronTocm;
82   }else{
83     fXt2=stept*(bint-1);
84     fXm2=fXt2-GetCellContent(bina,bint-1)*kMicronTocm;
85   }
86   if(isReco){
87     if(fXm1<fDrLen && fXm2>fDrLen) return;
88     if(bint==0 || bint==(fNDriftPts-1)) return;
89     if(fXm1>fDrLen){
90       for(Int_t itry=1; itry<=10; itry++){
91         Float_t xmtest=(bint-itry)*stept-GetCellContent(bina,bint-itry)*kMicronTocm;
92         if(xmtest<fDrLen){
93           fXt1=stept*(bint-itry);
94           fXt2=fXt1+stept;
95           fXm1=fXt1-GetCellContent(bina,bint-itry)*kMicronTocm;
96           fXm2=fXt2-GetCellContent(bina,bint+1-itry)*kMicronTocm;
97           return;
98         }
99       }
100     }
101     if(fXm2<fDrLen){
102       for(Int_t itry=1; itry<=10; itry++){
103         Float_t xmtest=(bint+1+itry)*stept-GetCellContent(bina,bint+1+itry)*kMicronTocm;
104         if(xmtest>fDrLen){
105           fXt1=stept*(bint+itry);
106           fXt2=fXt1+stept;
107           fXm1=fXt1-GetCellContent(bina,bint+itry)*kMicronTocm;
108           fXm2=fXt2-GetCellContent(bina,bint+1+itry)*kMicronTocm;
109           return;
110         }
111       }
112     }
113   }
114 }
115 //______________________________________________________________________
116 Float_t AliITSCorrMapSDD::GetCorrection(Float_t z, Float_t x, AliITSsegmentationSDD *seg){
117   // returns correction in cm starting from local coordinates on the module
118   ComputeGridPoints(z,x,seg,kTRUE);
119   Float_t m=(fXt2-fXt1)/(fXm2-fXm1);
120   Float_t q=fXt1-m*fXm1;
121   Float_t xcorr=m*fDrLen+q;
122   // fDrLen is the measured drift distance, xcorr is the corresponding true
123   return GetInversionBit() ? fDrLen-xcorr : xcorr-fDrLen; 
124 }
125 //______________________________________________________________________
126 Float_t AliITSCorrMapSDD::GetShiftForSimulation(Float_t z, Float_t x, AliITSsegmentationSDD *seg){
127   // returns shift to be appiled in digitizarion (in cm) starting from local coordinates on the module
128   ComputeGridPoints(z,x,seg,kFALSE);
129   Float_t m=(fXm2-fXm1)/(fXt2-fXt1);
130   Float_t q=fXm1-m*fXt1;
131   Float_t xshifted=m*fDrLen+q;
132   // fDrLen is the true drift distance, xshifted is the one with map shift
133   return GetInversionBit() ? xshifted-fDrLen : fDrLen-xshifted;
134 }
135 //______________________________________________________________________
136 TH2F* AliITSCorrMapSDD::GetMapHisto() const{
137   // Returns a TH2F histogram with map of residuals
138   TString hname;
139   hname.Form("h%s",GetName());
140   TH2F* hmap=new TH2F(hname.Data(),"",fNAnodePts,-0.5,255.5,fNDriftPts,0.,35.);
141   for(Int_t iAn=0;iAn<fNAnodePts; iAn++){
142     for(Int_t iDr=0;iDr<fNDriftPts; iDr++){
143       hmap->SetBinContent(iAn+1,iDr+1,GetCellContent(iAn,iDr));
144     }
145   }
146   return hmap;
147 }
148 //______________________________________________________________________
149 TH1F* AliITSCorrMapSDD::GetMapProfile() const{
150   // Returns a TH1F with the projection of the map along drift coordinate
151   TString hname;
152   hname.Form("p%s",GetName());
153   TH1F* hprof=new TH1F(hname.Data(),"",fNDriftPts,0.,35.);
154   for(Int_t iDr=0;iDr<fNDriftPts; iDr++){
155     Float_t meanval=0.;
156     for(Int_t iAn=0;iAn<fNAnodePts; iAn++){
157       meanval+=GetCellContent(iAn,iDr);
158     }
159     hprof->SetBinContent(iDr+1,meanval/fNAnodePts);
160   }
161   return hprof;
162   
163 }
164 //______________________________________________________________________
165 TH1F* AliITSCorrMapSDD::GetResidualDistr(Float_t dmin, Float_t dmax) const{
166   // Returns a TH1F histogram with distribution of residual
167   TString hname;
168   hname.Form("hd%s",GetName());
169   TH1F* hd=new TH1F(hname.Data(),"",100,dmin,dmax);
170   for(Int_t iAn=0;iAn<fNAnodePts; iAn++){
171     for(Int_t iDr=0;iDr<fNDriftPts; iDr++){
172       hd->Fill(GetCellContent(iAn,iDr));
173     }
174   }
175   return hd;
176 }
177