]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaAnalysis/AliTRDqaAT.cxx
29c3b94d47a8a16bcbf5dc0013bb39fc7d111a9f
[u/mrichter/AliRoot.git] / TRD / qaAnalysis / AliTRDqaAT.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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: AliTRDqaAT.cxx  $ */
17
18 //
19 // This class is a part of a package of high level QA monitoring for TRD.
20 // In this class provides a commonly used tools as static functions.
21 //
22 // S. Radomski
23 // radomski@physi.uni-heidelberg.de
24 // March 2008
25 //
26
27 #include "AliTRDqaAT.h"
28
29 #include "TMath.h"
30 #include "TH1D.h"
31 #include "AliExternalTrackParam.h"
32
33 //______________________________________________________________________________
34
35 AliTRDqaAT::AliTRDqaAT() {
36   //
37   // Dummy contructor
38   //
39
40 }
41
42 //______________________________________________________________________________
43 Int_t AliTRDqaAT::GetSector(const Double_t alpha) 
44 {
45   // Gets the sector number 
46
47   Double_t size = TMath::DegToRad() * 20.;
48   Int_t sector = (Int_t)((alpha + TMath::Pi())/size);
49   return sector;
50 }
51
52 //______________________________________________________________________________
53
54 Int_t AliTRDqaAT::GetStack(const AliExternalTrackParam *paramOut) 
55 {
56   //
57   // calculates the stack the track is in
58   //
59   
60   const Double_t l = -0.9;
61   const Double_t w = (2*l)/5;
62
63   Double_t tan = paramOut->GetZ() / paramOut->GetX();
64   Double_t pos = (tan - l) / w;
65   return (Int_t) pos;
66 }
67
68 //______________________________________________________________________________
69
70 void AliTRDqaAT::BuildRatio(TH1D *ratio, TH1D *histN, TH1D*histD) {
71   //
72   // Calculate the ratio of two histograms 
73   // error are calculated assuming the histos have the same counts
74   //
75
76   // calclate
77
78   Int_t nbins = histN->GetXaxis()->GetNbins();
79   for(Int_t i=1; i<nbins+2; i++) {
80     
81     Double_t valueN = histN->GetBinContent(i);
82     Double_t valueD = histD->GetBinContent(i);
83     
84     if (valueD < 1) {
85       ratio->SetBinContent(i, 0);
86       ratio->SetBinError(i, 0);
87       continue;
88     }
89
90     Double_t eps = (valueN < valueD-valueN)? valueN : valueD-valueN;
91     
92     ratio->SetBinContent(i, valueN/valueD);
93     ratio->SetBinError(i, TMath::Sqrt(eps)/valueD);
94   }
95
96   // style
97   ratio->SetMinimum(-0.1);
98   ratio->SetMaximum(1.1);
99   ratio->SetMarkerStyle(20);
100 }
101 //__________________________________________________________________________