BOSS 7.0.8
BESIII Offline Software System
Loading...
Searching...
No Matches
ResiAlign.cxx
Go to the documentation of this file.
4
5#include "GaudiKernel/INTupleSvc.h"
6#include "GaudiKernel/IService.h"
7#include "GaudiKernel/IMessageSvc.h"
8#include "GaudiKernel/StatusCode.h"
9#include "GaudiKernel/ISvcLocator.h"
10#include "GaudiKernel/Bootstrap.h"
11
12#include "EventModel/Event.h"
14#include "MdcRawEvent/MdcDigi.h"
16#include "Identifier/MdcID.h"
17
18#include "TF1.h"
19#include "TCanvas.h"
20
21#include <iostream>
22#include <fstream>
23#include <iomanip>
24#include <string>
25#include <cstring>
26
27using namespace std;
28
30 //m_ndiv = 6;
31 m_ndiv = 12;
32 m_resiCut = 1.0;
33 for(int i=0; i<NEP; i++) m_npoint[i] = 0;
34
35 for(int lay=0; lay<LAYERNMAX; lay++){
36 if(lay < 8){
37 m_docaMin[lay] = 0.6; // mm
38 m_docaMax[lay] = 4.8; // mm
39 } else{
40 m_docaMin[lay] = 0.8; // mm
41 m_docaMax[lay] = 6.4; // mm
42 }
43 }
44 for(int lay=0; lay<LAYERNMAX; lay++){
45 if((0==lay) || (7==lay) || (8==lay) || (19==lay) || (20==lay) ||
46 (35==lay) || (36==lay) || (42==lay) ) m_layBound[lay] = true;
47 else m_layBound[lay] = false;
48 }
49
50 m_ncut1 = 0;
51 m_ncut2 = 0;
52 m_ncut3 = 0;
53 m_ncut4 = 0;
54 m_ncut5 = 0;
55 m_ncut6 = 0;
56 m_ncut7 = 0;
57 m_ncut8 = 0;
58 m_ncut9 = 0;
59 m_ncut10 = 0;
60 m_ncut11 = 0;
61 m_ncut12 = 0;
62 m_ncut13 = 0;
63}
64
66}
67
69 delete m_hresAll;
70 delete m_hresInn;
71 delete m_hresStp;
72 delete m_hresOut;
73 for(int lay=0; lay<LAYERNMAX; lay++) delete m_hresLay[lay];
74 for(int i=0; i<NEP; i++) delete m_gr[i];
75}
76
77void ResiAlign::initialize(TObjArray* hlist, IMdcGeomSvc* mdcGeomSvc,
78 IMdcCalibFunSvc* mdcFunSvc){
79 IMessageSvc* msgSvc;
80 Gaudi::svcLocator() -> service("MessageSvc", msgSvc);
81 MsgStream log(msgSvc, "ResiAlign");
82 log << MSG::INFO << "ResiAlign::initialize()" << endreq;
83
84 m_hlist = hlist;
85 m_mdcGeomSvc = mdcGeomSvc;
86 m_mdcFunSvc = mdcFunSvc;
87
88 double zeast;
89 for(int lay=0; lay<43; lay++){
90 zeast = m_mdcGeomSvc->Wire(lay, 0)->Backward().z();
91 m_zrange[lay][1] = 2.0 * fabs(zeast) / (double)m_ndiv;
92 m_zrange[lay][0] = -1.0 * m_zrange[lay][1];
93
94 m_radii[lay] = m_mdcGeomSvc->Layer(lay)->Radius();
95 }
96
97 for(int wir=0; wir<WIRENMAX; wir++){
98 m_xe[wir] = m_mdcGeomSvc->Wire(wir)->Backward().x();
99 m_ye[wir] = m_mdcGeomSvc->Wire(wir)->Backward().y();
100 m_ze[wir] = m_mdcGeomSvc->Wire(wir)->Backward().z();
101 m_xw[wir] = m_mdcGeomSvc->Wire(wir)->Forward().x();
102 m_yw[wir] = m_mdcGeomSvc->Wire(wir)->Forward().y();
103 m_zw[wir] = m_mdcGeomSvc->Wire(wir)->Forward().z();
104 }
105
106 char hname[200];
107 int iEP;
108
109 INTupleSvc* ntupleSvc;
110 Gaudi::svcLocator() -> service("NTupleSvc", ntupleSvc);
111 for(iEP=0; iEP<=NEP; iEP++){
112 if(iEP < NEP) sprintf(hname, "FILE137/align%02d", iEP);
113 else sprintf(hname, "FILE137/alignAll");
114
115 NTuplePtr nt(ntupleSvc, hname);
116 if( nt ) m_tuple[iEP] = nt;
117 else{
118 m_tuple[iEP] = ntupleSvc->book(hname, CLID_ColumnWiseTuple,"align");
119 if (m_tuple[iEP]) {
120 m_tuple[iEP]->addItem ("run", m_iRun[iEP]);
121 m_tuple[iEP]->addItem ("evt", m_iEvt[iEP]);
122 m_tuple[iEP]->addItem ("resi", m_resi[iEP]);
123 m_tuple[iEP]->addItem ("p", m_p[iEP]);
124 m_tuple[iEP]->addItem ("pt", m_pt[iEP]);
125 m_tuple[iEP]->addItem ("phi", m_phi[iEP]);
126 m_tuple[iEP]->addItem ("lay", m_lay[iEP]);
127 m_tuple[iEP]->addItem ("lr", m_lr[iEP]);
128 m_tuple[iEP]->addItem ("cel", m_cel[iEP]);
129 }
130 else {
131 log << MSG::FATAL << "Cannot book N-tuple:"
132 << long(m_tuple[iEP]) << endmsg;
133 }
134 }
135 }
136
137 m_hnTrk = new TH1F("HNtrack", "", 10, -0.5, 9.5);
138 m_hlist->Add(m_hnTrk);
139
140 m_hnHit = new TH1F("HNhit", "", 100, -0.5, 99.5);
141 m_hlist->Add(m_hnHit);
142
143 m_hlayHitmap = new TH1F("Hitmap", "", 43, -0.5, 42.5);
144 m_hlist->Add(m_hlayHitmap);
145
146 m_hresAll = new TH1F("HResAllInc", "", 200, -1.0, 1.0);
147 m_hlist->Add(m_hresAll);
148
149 m_hresInn = new TH1F("HResInnInc", "", 200, -1.0, 1.0);
150 m_hlist->Add(m_hresInn);
151
152 m_hresStp = new TH1F("HResStpInc", "", 200, -1.0, 1.0);
153 m_hlist->Add(m_hresStp);
154
155 m_hresOut = new TH1F("HResOutInc", "", 200, -1.0, 1.0);
156 m_hlist->Add(m_hresOut);
157
158 int lay;
159 for(lay=0; lay<LAYERNMAX; lay++){
160 sprintf(hname, "Res_Layer%02d", lay);
161 m_hresLay[lay] = new TH1F(hname, "", 200, -1.0, 1.0);
162 m_hlist->Add(m_hresLay[lay]);
163 }
164
165 for(iEP=0; iEP<NEP; iEP++){
166 m_gr[iEP] = new TGraph();
167 sprintf(hname, "grResi%02d", iEP);
168 m_gr[iEP]->SetName(hname);
169 m_hlist->Add(m_gr[iEP]);
170 }
171 m_fevt.open("evt.txt");
172}
173
175 IMessageSvc* msgSvc;
176 Gaudi::svcLocator() -> service("MessageSvc", msgSvc);
177 MsgStream log(msgSvc, "ResiAlign");
178 log << MSG::DEBUG << "ResiAlign::fillHist()" << endreq;
179
180 bool esCutFg = event->getEsCutFlag();
181 if( ! esCutFg ){
182 m_ncut1++;
183 return true;
184 }
185
186 int i = 0;
187 int k;
188
189 int trkStat;
190 double dr;
191 double phi0;
192 double kappa;
193 double dz;
194 double tanl;
195 double chisq;
196 double p;
197 double pt;
198
199 int nhits;
200 int lay;
201 int cel;
202 int wir;
203 int lr;
204 int iEnd;
205 int iEP;
206
207 double doca;
208 double resi;
209 double zhit;
210 double wphi;
211 double dphi;
212 double hitPhi;
213 double xx;
214 double yy;
215 double rr;
216 int stat;
217 MdcAliRecTrk* rectrk;
218 MdcAliRecHit* rechit;
219 int nhitlay;
220 bool fgHitLay[LAYERNMAX];
221
222 IDataProviderSvc* eventSvc = NULL;
223 Gaudi::svcLocator()->service("EventDataSvc", eventSvc);
224 SmartDataPtr<Event::EventHeader> eventHeader(eventSvc,"/Event/EventHeader");
225 if (!eventHeader) {
226 log << MSG::FATAL << "Could not find Event Header" << endreq;
227 m_ncut3++;
228 return( StatusCode::FAILURE);
229 }
230 int iEvt = eventHeader->eventNumber();
231 int iRun = eventHeader->runNumber();
232
233 int nTrk = event -> getNTrk();
234 m_hnTrk->Fill(nTrk);
235 m_fevt << setw(10) << iRun << setw(10) << iEvt << setw(10) << nTrk << endl;
236 if((nTrk < m_param.nTrkCut[0]) || (nTrk > m_param.nTrkCut[1])){
237 m_ncut2++;
238 return true;
239 }
240
241 for(i=0; i<nTrk; i++){
242 rectrk = event->getRecTrk(i);
243 nhits = rectrk->getNHits();
244 trkStat = rectrk->getStat();
245// if (0 != trkStat) continue;
246
247 // dr cut
248 dr = rectrk -> getDr();
249 if(fabs(dr) > m_param.drCut){ m_ncut4++; continue; }
250
251 phi0 = rectrk -> getPhi0();
252 kappa = rectrk -> getKappa();
253
254 // dz cut
255 dz = rectrk -> getDz();
256 if(fabs(dz) > m_param.dzCut){ m_ncut5++; continue; }
257
258 for(lay=0; lay<LAYERNMAX; lay++){
259 fgHitLay[lay] = false;
260 }
261
262 m_hnHit->Fill(nhits);
263 for(k=0; k<nhits; k++){
264 rechit = rectrk -> getRecHit(k);
265 lay = rechit -> getLayid();
266 fgHitLay[lay] = true;
267 }
268
269 nhitlay = 0;
270 for(lay=0; lay<LAYERNMAX; lay++){
271 if(fgHitLay[lay]) nhitlay++;
272 }
273 if(nhitlay < m_param.nHitLayCut){ m_ncut6++; continue; }
274
275 tanl = rectrk -> getTanLamda();
276 chisq = rectrk -> getChisq();
277 p = rectrk -> getP();
278 pt = rectrk -> getPt();
279
280 if((fabs(pt)<m_param.ptCut[0]) || (fabs(pt)>m_param.ptCut[1])){ m_ncut7++; continue;}
281
282 for(k=0; k<nhits; k++){
283 rechit = rectrk->getRecHit(k);
284 lay = rechit->getLayid();
285 cel = rechit->getCellid();
286 lr = rechit->getLR();
287 doca = rechit -> getDocaInc();
288 zhit = rechit->getZhit();
289
290 stat = rechit -> getStat();
291 if((1 == m_param.hitStatCut) && (1 != stat)){ m_ncut8++; continue; }
292
293 if (1 == m_param.resiType) {
294 resi = rechit->getResiExcLR();
295 } else {
296 resi = rechit->getResiIncLR();
297 }
298 resi *= -1.0;
299 if( (1==isnan(resi)) || (fabs(resi) > m_resiCut) ||
300 (fabs(doca) > m_docaMax[lay]) || (fabs(doca) < m_docaMin[lay]) ){
301 m_ncut9++;
302 continue;
303 }
304
305 if(m_param.fgAdjacLayerCut){
306 if(0 == lay){
307 if( ! fgHitLay[1] ){ m_ncut10++; continue; }
308 } else if(42 == lay){
309 if( ! fgHitLay[41] ){ m_ncut11++; continue; }
310 } else{
311 if( (!fgHitLay[lay-1]) && (!fgHitLay[lay+1]) ){ m_ncut12++; continue; }
312 // for boundary layers
313 if( m_param.fgBoundLayerCut && m_layBound[lay] &&
314 ((!fgHitLay[lay-1]) || (!fgHitLay[lay+1])) ){
315 m_ncut13++;
316 continue;
317 }
318 }
319 }
320
321 m_hlayHitmap->Fill(lay);
322
323 // fill alignment trees
324 if((zhit < m_zrange[lay][0]) || (zhit > m_zrange[lay][1])){
325 wir = m_mdcGeomSvc -> Wire(lay, cel) -> Id();
326 xx = (zhit - m_zw[wir]) * (m_xe[wir] - m_xw[wir]) /
327 (m_ze[wir] - m_zw[wir]) + m_xw[wir];
328 yy = (zhit - m_zw[wir]) * (m_ye[wir] - m_yw[wir]) /
329 (m_ze[wir] - m_zw[wir]) + m_yw[wir];
330 rr = sqrt( (xx * xx) + (yy * yy) );
331 dphi = fabs(doca) / m_radii[lay];
332
333 if( yy >= 0 ) wphi = acos(xx / rr);
334 else wphi = PI2 - acos(xx / rr);
335 if(1 == lr) hitPhi = wphi + dphi; // mention
336 else hitPhi = wphi - dphi;
337 if(hitPhi < 0) hitPhi += PI2;
338 else if(hitPhi > PI2) hitPhi -= PI2;
339
340 if(zhit < m_zrange[lay][0]) iEnd = 0; // west
341 else iEnd = 1; // east
342 iEP = Alignment::getEpId(lay, iEnd);
343
344 m_iRun[iEP] = iRun;
345 m_iEvt[iEP] = iEvt;
346 m_resi[iEP] = resi;
347 m_p[iEP] = p;
348 m_pt[iEP] = pt;
349 m_phi[iEP] = hitPhi;
350 m_lay[iEP] = lay;
351 m_lr[iEP] = lr;
352 m_cel[iEP] = cel;
353 m_tuple[iEP]->write();
354
355 m_resi[NEP] = resi;
356 m_p[NEP] = p;
357 m_pt[NEP] = pt;
358 m_phi[NEP] = hitPhi;
359 m_lay[NEP] = lay;
360 m_lr[NEP] = lr;
361 m_cel[NEP] = cel;
362 m_tuple[NEP]->write();
363
364 m_hresAll->Fill(resi);
365 if(lay < 8) m_hresInn->Fill(resi);
366 else if(lay < 20) m_hresStp->Fill(resi);
367 else m_hresOut->Fill(resi);
368 m_hresLay[lay]->Fill(resi);
369
370 m_gr[iEP]->SetPoint(m_npoint[iEP], hitPhi, resi);
371 m_npoint[iEP]++;
372 }
373 }
374 }
375
376 return true;
377}
378
380 IMessageSvc* msgSvc;
381 Gaudi::svcLocator() -> service("MessageSvc", msgSvc);
382 MsgStream log(msgSvc, "ResiAlign");
383 log << MSG::INFO << "ResiAlign::updateConst()" << endreq;
384 m_fevt.close();
385
386 int iEP;
387 double par[3];
388 double err[3];
389 double dx;
390 double dy;
391 double rz;
392 double rLayer[] = { 120.225, 205.0, 237.55, 270.175,
393 302.625, 334.775, 366.65, 500.0,
394 120.225, 205.0, 237.55, 270.175,
395 302.625, 334.775, 366.65, 500.0 };
396
397 TCanvas c1("c1", "c1", 10, 10, 700, 500);
398
399 TF1* fResPhi = new TF1("fResPhi", funResi, 0, PI2, 3);
400 fResPhi->SetParameter(0, 0.0);
401 fResPhi->SetParameter(1, 0.0);
402 fResPhi->SetParameter(2, 0.0);
403
404 for(iEP=0; iEP<NEP; iEP++){
405 if((m_gr[iEP]->GetN()) > 500){
406 m_gr[iEP]->Fit("fResPhi", "V");
407 par[0] = fResPhi->GetParameter(0);
408 par[1] = fResPhi->GetParameter(1);
409 par[2] = fResPhi->GetParameter(2);
410
411 err[0] = fResPhi->GetParError(0);
412 err[1] = fResPhi->GetParError(1);
413 err[2] = fResPhi->GetParError(2);
414
415 dx = -1.0 * par[1];
416 dy = par[2];
417 rz = par[0] / rLayer[iEP];
418
419 // assume the shift of the outer section is 0
420 if (7==iEP || 15==iEP) {
421 dx = 0.0;
422 dy = 0.0;
423 rz = 0.0;
424 par[0] = 0.0;
425 par[1] = 0.0;
426 par[2] = 0.0;
427 }
428 alignPar->setDelDx(iEP, dx);
429 alignPar->setDelDy(iEP, dy);
430 alignPar->setDelRz(iEP, rz);
431
432 alignPar->setErrDx(iEP, err[1]);
433 alignPar->setErrDy(iEP, err[2]);
434 alignPar->setErrRz(iEP, err[0]/rLayer[iEP]);
435 }
436 }
437
438 cout << "TrackCut: cut1: " << m_ncut1 << ", cut2: " << m_ncut2 << ", cut3: " << m_ncut3
439 << ", cut4: " << m_ncut4 << ", cut5: " << m_ncut5 << ", cut6: " << m_ncut6
440 << ", cut7: " << m_ncut7 << endl;
441 cout << "HitCut: cut8: " << m_ncut8 << ", cut9: " << m_ncut9 << ", cut10: " << m_ncut10
442 << ", cut11: " << m_ncut11 << ", cut12: " << m_ncut12 << ", cut13: " << m_ncut13 << endl;
443
444 delete fResPhi;
445}
446
447Double_t ResiAlign::funResi(Double_t* x, Double_t* par){
448 Double_t val;
449 val = par[0] + par[1]*sin(x[0]) + par[2]*cos(x[0]);
450 return val;
451}
452
double sin(const BesAngle a)
Definition: BesAngle.h:210
double cos(const BesAngle a)
Definition: BesAngle.h:213
sprintf(cut,"kal_costheta0_em>-0.93&&kal_costheta0_em<0.93&&kal_pxy0_em>=0.05+%d*0.1&&kal_pxy0_em<0.15+%d*0.1&&NGch>=2", j, j)
INTupleSvc * ntupleSvc()
IMessageSvc * msgSvc()
#define NULL
virtual const MdcGeoLayer *const Layer(unsigned id)=0
virtual const MdcGeoWire *const Wire(unsigned id)=0
double dzCut
Definition: MdcAliParams.h:29
double drCut
Definition: MdcAliParams.h:28
int nTrkCut[2]
Definition: MdcAliParams.h:21
int fgAdjacLayerCut
Definition: MdcAliParams.h:19
double ptCut[2]
Definition: MdcAliParams.h:26
int fgBoundLayerCut
Definition: MdcAliParams.h:20
int getCellid() const
Definition: MdcAliRecHit.h:21
int getLR() const
Definition: MdcAliRecHit.h:22
double getResiExcLR() const
Definition: MdcAliRecHit.h:34
int getLayid() const
Definition: MdcAliRecHit.h:20
double getZhit() const
Definition: MdcAliRecHit.h:39
double getResiIncLR() const
Definition: MdcAliRecHit.h:31
MdcAliRecHit * getRecHit(int index) const
Definition: MdcAliRecTrk.h:38
int getStat() const
Definition: MdcAliRecTrk.h:26
int getNHits() const
Definition: MdcAliRecTrk.h:37
void setErrRz(int iEP, double val)
void setDelDy(int iEP, double val)
void setDelRz(int iEP, double val)
void setErrDx(int iEP, double val)
void setDelDx(int iEP, double val)
void setErrDy(int iEP, double val)
double Radius(void) const
Definition: MdcGeoLayer.h:160
HepPoint3D Forward(void) const
Definition: MdcGeoWire.h:129
HepPoint3D Backward(void) const
Definition: MdcGeoWire.h:128
static Double_t funResi(double *x, double *par)
void initialize(TObjArray *hlist, IMdcGeomSvc *mdcGeomSvc, IMdcCalibFunSvc *mdcFunSvc)
Definition: ResiAlign.cxx:77
bool fillHist(MdcAliEvent *event)
Definition: ResiAlign.cxx:174
void clear()
Definition: ResiAlign.cxx:68
void updateConst(MdcAlignPar *alignPar)
Definition: ResiAlign.cxx:379
const double PI2
Definition: Alignment.h:41
const int LAYERNMAX
Definition: Alignment.h:45
const int NEP
Definition: Alignment.h:48
int getEpId(int lay, int iEnd)
Definition: Alignment.cxx:18
const int WIRENMAX
Definition: Alignment.h:44
int nhits