BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
MdcCalib.h
Go to the documentation of this file.
1#ifndef MDCCALIB_H
2#define MDCCALIB_H
3
4#include "MdcGeomSvc/IMdcGeomSvc.h"
5#include "MdcGeomSvc/MdcGeoLayer.h"
6#include "MdcGeomSvc/MdcGeoWire.h"
7#include "MdcGeomSvc/MdcGeomSvc.h"
8
9#include "MdcCalibFunSvc/IMdcCalibFunSvc.h"
10#include "MdcCalibFunSvc/MdcCalibFunSvc.h"
11
15#include "MdcUtilitySvc/IMdcUtilitySvc.h"
16#include "MdcUtilitySvc/MdcUtilitySvc.h"
17
18#include "GaudiKernel/INTupleSvc.h"
19#include "GaudiKernel/NTuple.h"
20
21#include "TF1.h"
22#include "TFile.h"
23#include "TFolder.h"
24#include "TGraph.h"
25#include "TH1F.h"
26#include "TObjArray.h"
27#include "TROOT.h"
28#include "TTree.h"
29
30#include <fstream>
31#include <map>
32#include <vector>
33
34// const double pi = 3.14159265;
35class MdcCalib {
36public:
37 MdcCalib();
38 virtual ~MdcCalib();
39 virtual void initialize( TObjArray* hlist, IMdcGeomSvc* mdcGeomSvc,
40 IMdcCalibFunSvc* mdcFunSvc, IMdcUtilitySvc* mdcUtilitySvc ) = 0;
41
42 virtual void setParam( MdcCalParams& param ) = 0;
43 virtual int fillHist( MdcCalEvent* event ) = 0;
44 virtual int updateConst( MdcCalibConst* calconst ) = 0;
45 virtual void printCut() const = 0;
46 virtual void clear() = 0;
47
48private:
49 int getHresId( int lay, int entr, int lr, int bin ) const;
50 int calDetEffi();
51 bool getCellTrkPass( MdcCalEvent* event, int iTrk, int cellTrkPass[] );
52
53 MdcCalParams m_param;
54
55 TObjArray* m_hlist;
56 IMdcGeomSvc* m_mdcGeomSvc;
57 IMdcCalibFunSvc* m_mdcFunSvc;
58 IMdcUtilitySvc* m_mdcUtilitySvc;
59
60 int m_nEvt;
61 int m_cut1;
62 int m_nTrkAfTrkCut;
63 int m_cut2;
64 int m_cut3;
65 int m_cut4;
66 int m_cut5;
67 int m_cut6;
68 int m_cut7;
69 int m_nTrkCal;
70
71 int m_nlayer;
72 int m_nEntr[43];
73 int m_nBin[MdcCalNLayer];
74 bool fgReadWireEff;
75
76 bool m_layBound[MdcCalNLayer];
77
78 /* for calculating efficiency */
79 /* int m_nTrk[43]; */
80 /* int m_nGoodHit[43]; */
81 bool m_fgGoodWire[43][288];
82 TH1F* m_effNtrk;
83 TH1F* m_effNtrkRecHit;
84
85 /* calculating efficiency without the impact of track fitting */
86 int m_hitNum[43][2];
87
88 /* for track cut */
89 bool m_fgIni;
90 double m_docaMax[MdcCalNLayer];
91
92 /* for calculating hitphi */
93 double m_radii[MdcCalNLayer];
94 double m_xe[MdcCalTotCell];
95 double m_ye[MdcCalTotCell];
96 double m_ze[MdcCalTotCell];
97 double m_xw[MdcCalTotCell];
98 double m_yw[MdcCalTotCell];
99 double m_zw[MdcCalTotCell];
100
101 TFolder* m_fdcom;
102
103 TH1F* m_hresAllInc;
104 TH1F* m_hresAllExc;
105 TH1F* m_hresAllAve;
106
107 TH1F* m_hresInnInc;
108 TH1F* m_hresInnExc;
109
110 TH1F* m_hresStpInc;
111 TH1F* m_hresStpExc;
112
113 TH1F* m_hresOutInc;
114 TH1F* m_hresOutExc;
115
116 TH1F* m_hresAllExcUp;
117 TH1F* m_hresAllExcDw;
118 TH1F* m_hresInnExcUp;
119 TH1F* m_hresInnExcDw;
120 TH1F* m_hresStpExcUp;
121 TH1F* m_hresStpExcDw;
122 TH1F* m_hresOutExcUp;
123 TH1F* m_hresOutExcDw;
124
125 TFolder* m_fdResQ;
126 TH1F* m_hresAveAllQ[14];
127 TH1F* m_hresAveOutQ[14];
128 TH1F* m_hresAveLayQ[43][14];
129
130 TH1F* m_hratio;
131
132 TH1F* m_hTes[10];
133 TH1F* m_hbbTrkFlg;
134 TH1F* m_hTesAll;
135 TH1F* m_hTesGood;
136 TH1F* m_hTesAllFlag;
137 TH1F* m_hTesRec;
138 TH1F* m_hTesCalFlag;
139 TH1F* m_hTesCalUse;
140
141 TH1F* m_hnRawHit;
142
143 TH1F* m_hpt;
144 TH1F* m_hpMax;
145 TH1F* m_hpMaxCms;
146 TH1F* m_hptPos;
147 TH1F* m_hptNeg;
148
149 TH1F* m_hp;
150 TH1F* m_hp_cms;
151 TH1F* m_hpPos;
152 TH1F* m_hpNeg;
153 TH1F* m_hpPoscms;
154 TH1F* m_hpNegcms;
155
156 // results after four momentum cut
157 TH1F* m_hp_cut;
158
159 TH1F* m_hnTrk;
160 TH1F* m_hnTrkCal;
161 TH1F* m_hnhitsRec;
162 TH1F* m_hnhitsRecInn;
163 TH1F* m_hnhitsRecStp;
164 TH1F* m_hnhitsRecOut;
165 TH1F* m_hnhitsCal;
166 TH1F* m_hnhitsCalInn;
167 TH1F* m_hnhitsCalStp;
168 TH1F* m_hnhitsCalOut;
169 TH1F* m_wirehitmap;
170 TH1F* m_layerhitmap;
171 TH1F* m_hchisq;
172
173 // histogram for distribution of noise
174 TH1F* m_hnoisephi;
175 TH1F* m_hnoiselay;
176 TH1F* m_hnoisenhits;
177
178 TH1F* m_hdr;
179 TH1F* m_hphi0;
180 TH1F* m_hkap;
181 TH1F* m_hdz;
182 TH1F* m_htanl;
183 TH1F* m_hcosthe;
184 TH1F* m_hcostheNeg;
185 TH1F* m_hcosthePos;
186
187 TH1F* m_hx0;
188 TH1F* m_hy0;
189 TH1F* m_hdelZ0;
190 TGraph* m_grX0Y0;
191 int m_nGrPoint;
192
193 TH1F* m_hitEffAll;
194 TH1F* m_hitEffRaw;
195 TH1F* m_hitEffRec;
196
197 /* X-T ntuple */
198 int m_nEvtNtuple;
199 NTuple::Tuple* m_xtTuple[MdcCalNLayer];
200 NTuple::Item<long> m_cel[MdcCalNLayer];
201 NTuple::Item<long> m_lr[MdcCalNLayer];
202 NTuple::Item<long> m_run[MdcCalNLayer];
203 NTuple::Item<long> m_evt[MdcCalNLayer];
204 NTuple::Item<double> m_doca[MdcCalNLayer];
205 NTuple::Item<double> m_dm[MdcCalNLayer];
206 NTuple::Item<double> m_tdr[MdcCalNLayer];
207 NTuple::Item<double> m_tdc[MdcCalNLayer];
208 NTuple::Item<double> m_entr[MdcCalNLayer];
209 NTuple::Item<double> m_zhit[MdcCalNLayer];
210 NTuple::Item<double> m_qhit[MdcCalNLayer];
211 NTuple::Item<double> m_p[MdcCalNLayer];
212 NTuple::Item<double> m_pt[MdcCalNLayer];
213 NTuple::Item<double> m_phi0[MdcCalNLayer];
214 NTuple::Item<double> m_tanl[MdcCalNLayer];
215 NTuple::Item<double> m_hitphi[MdcCalNLayer];
216
217 /* histograms drift time */
218 TFolder* m_fdTime;
219 TH1F* m_htraw[MdcCalNLayer];
220 TH1F* m_htdr[MdcCalNLayer];
221 TH1F* m_htdrlr[MdcCalNLayer][2];
222
223 /* histograms of adc */
224 TFolder* m_fdAdc;
225 TH1F* m_hadc[MdcCalNLayer];
226
227 /* histograms for spatial resolution */
228 TFolder* m_fdres;
229 TH1F* m_hresInc[MdcCalNLayer];
230 TH1F* m_hreslrInc[MdcCalNLayer][2];
231 TH1F* m_hresExc[MdcCalNLayer];
232 TH1F* m_hreslrExc[MdcCalNLayer][2];
233 TH1F* m_hresphi[MdcCalNLayer][20];
234
235 TFolder* m_fdresAve;
236 TH1F* m_hresAve[MdcCalNLayer];
237 TH1F* m_hreslrAve[MdcCalNLayer][2];
238
239 /* histograms for momentum vs phi */
240 static const int NPhiBin = 20;
241 static const int NThetaBin = 9;
242 double m_phiWid;
243 double m_theWid;
244 TFolder* m_fdmomPhi;
245 // in experimental reference frame
246 TH1F* m_ppPhi[NPhiBin];
247 TH1F* m_pnPhi[NPhiBin];
248 TH1F* m_ppThe[NThetaBin];
249 TH1F* m_pnThe[NThetaBin];
250 TH1F* m_ppThePhi[NThetaBin][NPhiBin];
251 TH1F* m_pnThePhi[NThetaBin][NPhiBin];
252
253 // in CMS
254 TH1F* m_ppPhiCms[NPhiBin];
255 TH1F* m_pnPhiCms[NPhiBin];
256 TH1F* m_ppTheCms[NThetaBin];
257 TH1F* m_pnTheCms[NThetaBin];
258 TH1F* m_ppThePhiCms[NThetaBin][NPhiBin];
259 TH1F* m_pnThePhiCms[NThetaBin][NPhiBin];
260
261 /* ntuple for cosmic-ray */
262 NTuple::Tuple* m_cosmic;
263 NTuple::Item<double> m_pUpcos;
264 NTuple::Item<double> m_pDwcos;
265 NTuple::Item<double> m_ptUpcos;
266 NTuple::Item<double> m_ptDwcos;
267 NTuple::Item<double> m_phiUpcos;
268 NTuple::Item<double> m_phiDwcos;
269 NTuple::Item<double> m_drUpcos;
270 NTuple::Item<double> m_drDwcos;
271 NTuple::Item<double> m_dzUpcos;
272 NTuple::Item<double> m_dzDwcos;
273 NTuple::Item<double> m_ctheUpcos;
274 NTuple::Item<double> m_ctheDwcos;
275 NTuple::Item<long> m_nhitUpcos;
276 NTuple::Item<long> m_nhitDwcos;
277 NTuple::Item<long> m_chargecos;
278 NTuple::Item<long> m_tesFlagcos;
279 NTuple::Item<double> m_tescos;
280
281 /* histograms for spatial resolution vs distance */
282 double m_dwid;
283 TFolder* m_fdres2d;
284 std::map<int, int> m_mapr2d;
285 std::vector<TH1F*> m_hr2dInc;
286 std::vector<TH1F*> m_hr2dExc;
287
288 TFolder* m_fdres2t;
289 TH1F* m_hr2t[MdcCalNLayer][MdcCalNENTRXT][2][45];
290
291 /* for the index of resolution histograms */
292 static const int HRESLAYER_INDEX = 10;
293 static const int HRESLAYER_MASK = 0xFC00;
294
295 static const int HRESENTRA_INDEX = 7;
296 static const int HRESENTRA_MASK = 0x380;
297
298 static const int HRESLR_INDEX = 5;
299 static const int HRESLR_MASK = 0x60;
300
301 static const int HRESBIN_INDEX = 0;
302 static const int HRESBIN_MASK = 0x1F;
303};
304
305inline void MdcCalib::setParam( MdcCalParams& param ) { m_param = param; }
306
307#endif /* MDCCALIB_H */
*******INTEGER m_nBinMax INTEGER m_NdiMax !No of bins in histogram for cell exploration division $ !Last vertex $ !Last active cell $ !Last cell in buffer $ !No of sampling when dividing cell $ !No of function total $ !Flag for random ceel for $ !Flag for type of for WtMax $ !Flag which decides whether vertices are included in the sampling $ entire domain is hyp !Maximum effective eevents per bin
Definition FoamA.h:85
const int MdcCalNLayer
Definition MdcCalParams.h:6
const int MdcCalNENTRXT
const int MdcCalTotCell
Definition MdcCalParams.h:9
virtual void printCut() const =0
virtual void clear()=0
Definition MdcCalib.cxx:83
virtual ~MdcCalib()
Definition MdcCalib.cxx:81
virtual void initialize(TObjArray *hlist, IMdcGeomSvc *mdcGeomSvc, IMdcCalibFunSvc *mdcFunSvc, IMdcUtilitySvc *mdcUtilitySvc)=0
Definition MdcCalib.cxx:225
virtual int updateConst(MdcCalibConst *calconst)=0
virtual int fillHist(MdcCalEvent *event)=0
Definition MdcCalib.cxx:767
virtual void setParam(MdcCalParams &param)=0
Definition MdcCalib.h:305