]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITSMapSDD.cxx
Added option TGeo in the constructor in order to initialize some parameters directly...
[u/mrichter/AliRoot.git] / ITS / AliITSMapSDD.cxx
CommitLineData
efb451bf 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 "AliITSMapSDD.h"
28
29ClassImp(AliITSMapSDD)
30//______________________________________________________________________
31AliITSMapSDD::AliITSMapSDD():TNamed("defaultmap","")
32{
33 // default constructor
34 for(Int_t iAn=0;iAn<fgkNAnodPts; iAn++){
35 for(Int_t iDr=0;iDr<fgkNDrifPts; iDr++){
36 fMap[iAn][iDr]=0;
37 }
38 }
39}
40//______________________________________________________________________
41AliITSMapSDD::AliITSMapSDD(Char_t *mapname):TNamed(mapname,"")
42{
43 // standard constructor
44 for(Int_t iAn=0;iAn<fgkNAnodPts; iAn++){
45 for(Int_t iDr=0;iDr<fgkNDrifPts; iDr++){
46 fMap[iAn][iDr]=0;
47 }
48 }
49}
50
51//______________________________________________________________________
52void AliITSMapSDD::SetMap(TH2F* hmap){
53 // Fill map staring from 2D histo
54 // with anodes on x axis and drift dist. on y axis
55 for(Int_t iAn=0;iAn<fgkNAnodPts; iAn++){
56 for(Int_t iDr=0;iDr<fgkNDrifPts; iDr++){
57 fMap[iAn][iDr]=hmap->GetBinContent(iAn+1,iDr+1);
58 }
59 }
60}
61//______________________________________________________________________
c9a38d3d 62Float_t AliITSMapSDD::GetCorrection(Float_t z, Float_t x, AliITSsegmentationSDD *seg){
63 // returns correction in cm starting from local coordinates on the module
64 const Double_t kMicronTocm = 1.0e-4;
65 Int_t nAnodes=seg->Npz();
66 Int_t nAnodesHybrid=seg->NpzHalf();
67 Int_t bina =(Int_t) seg->GetAnodeFromLocal(x,z);
68 if(bina>nAnodes) AliError("Wrong anode anumber!");
69 if(bina>=nAnodesHybrid) bina-=nAnodesHybrid;
70 Float_t stept = seg->Dx()*kMicronTocm/(Float_t)fgkNDrifPts;
71 Int_t bint = TMath::Abs((Int_t)(x/stept));
72 if(bint==fgkNDrifPts) bint-=1;
73 if(bint>=fgkNDrifPts) AliError("Wrong bin number along drift direction!");
74 return kMicronTocm*GetCellContent(bina,bint);
75}
76//______________________________________________________________________
efb451bf 77TH2F* AliITSMapSDD::GetMapHisto() const{
78 // Returns a TH2F histogram with map of residuals
79 Char_t hname[50];
80 sprintf(hname,"h%s",GetName());
81 TH2F* hmap=new TH2F(hname,"",fgkNAnodPts,-0.5,255.5,fgkNDrifPts,0.,35.);
82 for(Int_t iAn=0;iAn<fgkNAnodPts; iAn++){
83 for(Int_t iDr=0;iDr<fgkNDrifPts; iDr++){
84 hmap->SetBinContent(iAn+1,iDr+1,fMap[iAn][iDr]);
85 }
86 }
87 return hmap;
88}
89//______________________________________________________________________
90TH1F* AliITSMapSDD::GetResidualDistr(Float_t dmin, Float_t dmax) const{
91 // Returns a TH1F histogram with distribution of residual
92 Char_t hname[50];
93 sprintf(hname,"hd%s",GetName());
94 TH1F* hd=new TH1F(hname,"",100,dmin,dmax);
95 for(Int_t iAn=0;iAn<fgkNAnodPts; iAn++){
96 for(Int_t iDr=0;iDr<fgkNDrifPts; iDr++){
97 hd->Fill(fMap[iAn][iDr]);
98 }
99 }
100 return hd;
101}