]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDdEdxUtils.h
Fix
[u/mrichter/AliRoot.git] / TRD / AliTRDdEdxUtils.h
CommitLineData
fb3d369e 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//
17//
18// Xianguo Lu
19// lu@physi.uni-heidelberg.de
20// Xianguo.Lu@cern.ch
21//
22/*
23grep " AliTRDdEdxUtils::" AliTRDdEdxUtils.cxx | grep "=" -v | grep -v "[6]" | grep -v printf |wc
24grep "(" AliTRDdEdxUtils.h | grep ";" | grep -v grep | grep -v ClassDef | grep -v "{" | grep -v typedef | wc
25*/
26
27
28#ifndef ALITRDDEDXUTILS_H
29#define ALITRDDEDXUTILS_H
30
31#ifndef TVECTORD_H
32#include "TVectorD.h"
33#endif
34
35#ifndef THNSPARSE_H
36#include "THnSparse.h"
37#endif
38
39#ifndef TTREESTREAM_H
40#include "TTreeStream.h"
41#endif
42
43class TH1D;
44class TH2D;
45class TObjArray;
46
47class AliESDEvent;
48class AliESDtrack;
49class AliTRDcluster;
50class AliTRDtrackV1;
51class AliTRDseedV1;
52
53class AliTRDdEdxUtils
54{
55 public:
56
57 //===================================================================================
58 // Math and Histogram
59 //===================================================================================
60 static void GetCDFCuts(const TH1D *hh, const Int_t ncut, Double_t cuts[], const Double_t cdfs[], const Double_t thres);
61 static Double_t GetMeanRMS(const Double_t nn, const Double_t sum, const Double_t w2s, Double_t * grms=0x0, Double_t * gerr=0x0);
62 static Double_t TruncatedMean(const Int_t nx, const Double_t xdata[], const Double_t lowfrac, const Double_t highfrac, Double_t * grms=0x0, Double_t * gerr=0x0, Double_t *wws=0x0);
63 static Double_t TruncatedMean(const TH1 *hh, const Double_t lowfrac, const Double_t highfrac, Double_t * grms=0x0, Double_t * gerr=0x0);
64 static void FitSlicesY(const TH2D *hh, TH1D *&hnor, TH1D *&hmpv, TH1D *&hwid, TH1D *&hres, const Double_t thres, const Double_t lowfrac, const Double_t highfrac);
65
66 //===================================================================================
67 // TRD Analysis Fast Tool
68 //===================================================================================
69 static Int_t GetNtracklet(const AliESDEvent *esd);
70 static AliTRDtrackV1 * GetTRDtrackV1(const AliESDtrack * esdtrack);
71 static Bool_t IsInSameStack(const AliTRDtrackV1 *trdtrack);
72 static Bool_t GetFirstSectorStackMomentum(const AliTRDtrackV1 *trdtrack, Int_t & isec, Int_t & istk, Double_t & mom);
73
74 //===================================================================================
75 // Calibration
76 //===================================================================================
77 static void SetCalibFile(const TString file){fgCalibFile = file;}
78 static void DeleteCalibObj();
79 static void IniCalibObj();
80
81 static void SetObjPHQ(TObjArray * obj){fgObjPHQ = obj;}
82 static Bool_t GenerateDefaultPHQOCDB(const TString path="local://./");
83
84 static Double_t GetCalibTPCscale(const Int_t tpcncls, const Double_t tpcsig);
85
86 static void DeleteCalibHist();
87 static void IniCalibHist(TList *list, const Bool_t kPHQonly=kFALSE);
88 static Bool_t ReadCalibHist(const TString filename, const TString listname);
89
26733ed4 90 static TObjArray * GetObjPHQ();
fb3d369e 91 static TObjArray * GetHistPHQ(){return fgHistPHQ;}
92 static TObjArray * GetObjPHQ(const Bool_t kinvq, const Double_t mag, const Int_t charge);
93 static THnSparseD * GetHistPHQ(const Bool_t kinvq, const Double_t mag, const Int_t charge);
94
95 static void FillCalibHist(const Int_t ncls, const TVectorD *arrayQ, const TVectorD *arrayX, THnSparseD * hcalib, const Double_t scale);
96 static void FillCalibHist(const AliTRDtrackV1 *trdv1, const Bool_t kinvq, const Double_t mag, const Int_t charge, const Double_t scale);
97
98 static void CalibOutput(const TList *l0, Int_t run);
99 static TObjArray* GetCalibObj(const THnSparseD *hh, Int_t run=-999, TList *lout=0x0, TTreeSRedirector *calibStream=0x0);
100
101 static THnSparseD * GetHistGain(){return fgHistGain;}
102 static THnSparseD * GetHistT0(){return fgHistT0;}
103 static THnSparseD * GetHistVd(){return fgHistVd;}
104
105 //===================================================================================
106 // dEdx calculation
107 //===================================================================================
108 static Double_t GetSignal(const Int_t nch, const Int_t ncls, const Double_t qq);
109 static Int_t GetNch(const Double_t signal);
110 static Int_t GetNcls(const Double_t signal);
111 static Double_t GetQ(const Double_t signal);
112
113 static Double_t ToyCook(const Bool_t kinvq, Int_t &ncluster, TVectorD *arrayQ, TVectorD *arrayX, const TObjArray *cobj=0x0);
114 static Double_t CombineddEdx(const Bool_t kinvq, Int_t &concls, TVectorD *coarrayQ, TVectorD *coarrayX, const Int_t tpcncls, const TVectorD *tpcarrayQ, const TVectorD *tpcarrayX, const Int_t trdncls, const TVectorD *trdarrayQ, const TVectorD *trdarrayX);
115
116 //===================================================================================
117 // dEdx Getter and Setter
118 //===================================================================================
119 static Int_t GetArrayClusterQ(const Bool_t kinvq, TVectorD *arrayQ, TVectorD *arrayX, const AliTRDtrackV1 *trdtrack, Int_t timeBin0=-1, Int_t timeBin1=1000, Int_t tstep=1);
120 static Int_t UpdateArrayX(const Int_t ncls, TVectorD* arrayX);
121 static void SetChamberQT(const AliTRDtrackV1 *trdtrack, const Int_t kcalib, THnSparseD * hgain=0x0, THnSparseD * ht0=0x0, THnSparseD * hvd=0x0);
122
123 static Int_t GetNchamber() {return fgNchamber;}
124 static Double_t GetChamberQ(const Int_t ich) {return fgChamberQ[ich];}
125 static Double_t GetChamberTmean(const Int_t ich) {return fgChamberTmean[ich];}
126 static Double_t GetTrackTmean() {return fgTrackTmean;}
127
128 //===================================================================================
129 // dEdx Parameterization
130 //===================================================================================
131
132 static Double_t MeandEdx(const Double_t * xx, const Double_t * par);
133 static Double_t MeanTR(const Double_t * xx, const Double_t * par);
134 static Double_t MeandEdxTR(const Double_t * xx, const Double_t * par);
135
136 static Double_t QMeanTPC(const Double_t bg);
137 static Double_t Q0MeanTRDpp(const Double_t bg);
138 static Double_t Q1MeanTRDpp(const Double_t bg);
139 static Double_t Q0MeanTRDPbPb(const Double_t bg);
140 static Double_t Q1MeanTRDPbPb(const Double_t bg);
141
142 typedef Double_t (*FFunc)(const Double_t *xx, const Double_t *par);
143
144 static Double_t MeandEdxLogx(const Double_t * xx, const Double_t * par){return ToLogx(MeandEdx, xx, par);}
145 static Double_t MeanTRLogx(const Double_t * xx, const Double_t * par){return ToLogx(MeanTR, xx, par);}
146 static Double_t MeandEdxTRLogx(const Double_t * xx, const Double_t * par){return ToLogx(MeandEdxTR, xx, par);}
147
148 //===================================================================================
149 // Detector, Data and Control Constant
150 //===================================================================================
151 static Int_t NTRDchamber(){return 18*5*6;} //540
152 static Int_t NTRDtimebin(){return NTRDchamber()*31;} //16740
153 static Int_t ToDetector(const Int_t gtb);
154 static Int_t ToTimeBin(const Int_t gtb);
155 static Int_t ToSector(const Int_t gtb);
156 static Int_t ToStack(const Int_t gtb);
157 static Int_t ToLayer(const Int_t gtb);
158
159 static TString GetRunType(const Int_t run);
160
161 static Bool_t IsExBOn(){return kTRUE;}
162
163 static void PrintControl();
164 //===================================================================================
165 //===================================================================================
166
167 private:
168
169 //dEdx Getter and Setter
170 static Double_t GetAngularCorrection(const AliTRDseedV1 *seed);
171 static Double_t GetClusterQ(const Bool_t kinvq, const AliTRDseedV1 * seed, const Int_t itb);
172
173 //dEdx Parameterization
174 static Double_t ToLogx(FFunc func, const Double_t * xx, const Double_t * par);
175
176 //Calibration
177 static Int_t ApplyCalib(const Int_t nc0, TVectorD *arrayQ, TVectorD *arrayX, const TObjArray *cobj);
178 static void GetPHCountMeanRMS(const TH1D *hnor, TH1D *&hmean);
179 static Int_t GetPHQIterator(const Bool_t kinvq, const Double_t mag, const Int_t charge);
180 static TString GetPHQName(const Bool_t kobj, const Int_t iter);
181
182 //Detector, Data and Control Constant
183 static Double_t Q0Frac(){return 0.3;}
184 static Double_t Q1Frac(){return 0.5;}
185 static Double_t TimeBinCountCut(){return 0.0;}
186 static Int_t CalibTPCnclsCut(){return 70;}
187
188 static THnSparseD *fgHistGain;//PH hist
189 static THnSparseD *fgHistT0;//PH hist
190 static THnSparseD *fgHistVd;//PH hist
191 static TObjArray * fgHistPHQ;//array containing 8 THnSparseD!
192
193 static TString fgCalibFile; //private path for calibration object
194
195 static TObjArray * fgObjGain;//chamber gain obj
196 static TObjArray * fgObjT0;//t0 obj
197 static TObjArray * fgObjVd;//Vd obj
198 static TObjArray * fgObjPHQ;//array containing 8 TObjArray!
199
200 static Int_t fgNchamber; //number of chamber used in cookdedx
201 static Double_t fgChamberQ[6]; //dqdl in chamber [i]
202 static Double_t fgChamberTmean[6]; //Q-weighted timebin \sum Q*T / \sum Q
203
204 static Double_t fgTrackTmean; //mean timebin over track
205};
206
207#endif