]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITSCorrMapSDD.cxx
Fix for coverity
[u/mrichter/AliRoot.git] / ITS / AliITSCorrMapSDD.cxx
CommitLineData
54af1add 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
29const Int_t AliITSCorrMapSDD::fgkNAnodePtsDefault = 1;
30const Int_t AliITSCorrMapSDD::fgkNDriftPtsDefault = 72;
31
32ClassImp(AliITSCorrMapSDD)
33//______________________________________________________________________
34AliITSCorrMapSDD::AliITSCorrMapSDD():
35TNamed("defaultmap",""),
36fNAnodePts(fgkNAnodePtsDefault),
37fNDriftPts(fgkNDriftPtsDefault)
38{
39 // default constructor
40}
41//______________________________________________________________________
42AliITSCorrMapSDD::AliITSCorrMapSDD(Char_t *mapname):
43TNamed(mapname,""),
44fNAnodePts(fgkNAnodePtsDefault),
45fNDriftPts(fgkNDriftPtsDefault)
46{
47 // standard constructor
48}
49//______________________________________________________________________
50Float_t AliITSCorrMapSDD::GetCorrection(Float_t z, Float_t x, AliITSsegmentationSDD *seg){
51 // returns correction in cm starting from local coordinates on the module
52 const Double_t kMicronTocm = 1.0e-4;
53 Int_t nAnodes=seg->Npz();
54 Int_t nAnodesHybrid=seg->NpzHalf();
55 Int_t bina =(Int_t) seg->GetAnodeFromLocal(x,z);
56 if(bina>nAnodes) AliError("Wrong anode anumber!");
57 if(bina>=nAnodesHybrid) bina-=nAnodesHybrid;
58 Float_t stept = seg->Dx()*kMicronTocm/(Float_t)fNDriftPts;
6c790b32 59 Float_t drLen= seg->Dx()*kMicronTocm-TMath::Abs(x);
60 Int_t bint = TMath::Abs((Int_t)(drLen/stept));
54af1add 61 if(bint==fNDriftPts) bint-=1;
62 if(bint>=fNDriftPts) AliError("Wrong bin number along drift direction!");
63 return kMicronTocm*GetCellContent(bina,bint);
64}
65//______________________________________________________________________
66TH2F* AliITSCorrMapSDD::GetMapHisto() const{
67 // Returns a TH2F histogram with map of residuals
2c4e6a6a 68 TString hname;
69 hname.Form("h%s",GetName());
70 TH2F* hmap=new TH2F(hname.Data(),"",fNAnodePts,-0.5,255.5,fNDriftPts,0.,35.);
54af1add 71 for(Int_t iAn=0;iAn<fNAnodePts; iAn++){
72 for(Int_t iDr=0;iDr<fNDriftPts; iDr++){
73 hmap->SetBinContent(iAn+1,iDr+1,GetCellContent(iAn,iDr));
74 }
75 }
76 return hmap;
77}
78//______________________________________________________________________
6c790b32 79TH1F* AliITSCorrMapSDD::GetMapProfile() const{
80 // Returns a TH1F with the projection of the map along drift coordinate
2c4e6a6a 81 TString hname;
82 hname.Form("p%s",GetName());
83 TH1F* hprof=new TH1F(hname.Data(),"",fNDriftPts,0.,35.);
6c790b32 84 for(Int_t iDr=0;iDr<fNDriftPts; iDr++){
85 Float_t meanval=0.;
86 for(Int_t iAn=0;iAn<fNAnodePts; iAn++){
87 meanval+=GetCellContent(iAn,iDr);
88 }
89 hprof->SetBinContent(iDr+1,meanval/fNAnodePts);
90 }
91 return hprof;
92
93}
94//______________________________________________________________________
54af1add 95TH1F* AliITSCorrMapSDD::GetResidualDistr(Float_t dmin, Float_t dmax) const{
96 // Returns a TH1F histogram with distribution of residual
2c4e6a6a 97 TString hname;
98 hname.Form("hd%s",GetName());
99 TH1F* hd=new TH1F(hname.Data(),"",100,dmin,dmax);
54af1add 100 for(Int_t iAn=0;iAn<fNAnodePts; iAn++){
101 for(Int_t iDr=0;iDr<fNDriftPts; iDr++){
102 hd->Fill(GetCellContent(iAn,iDr));
103 }
104 }
105 return hd;
106}