1#include "MdcAlignAlg/ResiAlign.h"
2#include "MdcAlignAlg/MdcAlignAlg.h"
3#include "MdcAlignAlg/Alignment.h"
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"
12#include "EventModel/Event.h"
13#include "EventModel/EventHeader.h"
14#include "MdcRawEvent/MdcDigi.h"
15#include "Identifier/Identifier.h"
16#include "Identifier/MdcID.h"
33 for(
int i=0; i<
NEP; i++) m_npoint[i] = 0;
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;
73 for(
int lay=0; lay<
LAYERNMAX; lay++)
delete m_hresLay[lay];
74 for(
int i=0; i<
NEP; i++)
delete m_gr[i];
80 Gaudi::svcLocator() -> service(
"MessageSvc",
msgSvc);
81 MsgStream log(
msgSvc,
"ResiAlign");
82 log << MSG::INFO <<
"ResiAlign::initialize()" << endreq;
85 m_mdcGeomSvc = mdcGeomSvc;
86 m_mdcFunSvc = mdcFunSvc;
89 for(
int lay=0; lay<43; lay++){
91 m_zrange[lay][1] = 2.0 * fabs(zeast) / (double)m_ndiv;
92 m_zrange[lay][0] = -1.0 * m_zrange[lay][1];
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();
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");
116 if( nt ) m_tuple[iEP] = nt;
118 m_tuple[iEP] =
ntupleSvc->book(hname, CLID_ColumnWiseTuple,
"align");
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]);
131 log << MSG::FATAL <<
"Cannot book N-tuple:"
132 << long(m_tuple[iEP]) << endmsg;
137 m_hnTrk =
new TH1F(
"HNtrack",
"", 10, -0.5, 9.5);
138 m_hlist->Add(m_hnTrk);
140 m_hnHit =
new TH1F(
"HNhit",
"", 100, -0.5, 99.5);
141 m_hlist->Add(m_hnHit);
143 m_hlayHitmap =
new TH1F(
"Hitmap",
"", 43, -0.5, 42.5);
144 m_hlist->Add(m_hlayHitmap);
146 m_hresAll =
new TH1F(
"HResAllInc",
"", 200, -1.0, 1.0);
147 m_hlist->Add(m_hresAll);
149 m_hresInn =
new TH1F(
"HResInnInc",
"", 200, -1.0, 1.0);
150 m_hlist->Add(m_hresInn);
152 m_hresStp =
new TH1F(
"HResStpInc",
"", 200, -1.0, 1.0);
153 m_hlist->Add(m_hresStp);
155 m_hresOut =
new TH1F(
"HResOutInc",
"", 200, -1.0, 1.0);
156 m_hlist->Add(m_hresOut);
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]);
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]);
171 m_fevt.open(
"evt.txt");
176 Gaudi::svcLocator() -> service(
"MessageSvc",
msgSvc);
177 MsgStream log(
msgSvc,
"ResiAlign");
178 log << MSG::DEBUG <<
"ResiAlign::fillHist()" << endreq;
180 bool esCutFg =
event->getEsCutFlag();
222 IDataProviderSvc* eventSvc = NULL;
223 Gaudi::svcLocator()->service(
"EventDataSvc", eventSvc);
224 SmartDataPtr<Event::EventHeader> eventHeader(eventSvc,
"/Event/EventHeader");
226 log << MSG::FATAL <<
"Could not find Event Header" << endreq;
228 return( StatusCode::FAILURE);
230 int iEvt = eventHeader->eventNumber();
231 int iRun = eventHeader->runNumber();
233 int nTrk =
event -> getNTrk();
235 m_fevt << setw(10) << iRun << setw(10) << iEvt << setw(10) << nTrk << endl;
241 for(i=0; i<nTrk; i++){
242 rectrk =
event->getRecTrk(i);
248 dr = rectrk -> getDr();
249 if(fabs(dr) > m_param.
drCut){ m_ncut4++;
continue; }
251 phi0 = rectrk -> getPhi0();
252 kappa = rectrk -> getKappa();
255 dz = rectrk -> getDz();
256 if(fabs(dz) > m_param.
dzCut){ m_ncut5++;
continue; }
259 fgHitLay[lay] =
false;
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;
271 if(fgHitLay[lay]) nhitlay++;
273 if(nhitlay < m_param.
nHitLayCut){ m_ncut6++;
continue; }
275 tanl = rectrk -> getTanLamda();
276 chisq = rectrk -> getChisq();
277 p = rectrk -> getP();
278 pt = rectrk -> getPt();
280 if((fabs(pt)<m_param.
ptCut[0]) || (fabs(pt)>m_param.
ptCut[1])){ m_ncut7++;
continue;}
282 for(k=0; k<nhits; k++){
286 lr = rechit->
getLR();
287 doca = rechit -> getDocaInc();
290 stat = rechit -> getStat();
291 if((1 == m_param.
hitStatCut) && (1 != stat)){ m_ncut8++;
continue; }
299 if( (1==isnan(resi)) || (fabs(resi) > m_resiCut) ||
300 (fabs(doca) > m_docaMax[lay]) || (fabs(doca) < m_docaMin[lay]) ){
307 if( ! fgHitLay[1] ){ m_ncut10++;
continue; }
308 }
else if(42 == lay){
309 if( ! fgHitLay[41] ){ m_ncut11++;
continue; }
311 if( (!fgHitLay[lay-1]) && (!fgHitLay[lay+1]) ){ m_ncut12++;
continue; }
314 ((!fgHitLay[lay-1]) || (!fgHitLay[lay+1])) ){
321 m_hlayHitmap->Fill(lay);
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];
333 if( yy >= 0 ) wphi = acos(xx / rr);
334 else wphi =
PI2 - acos(xx / rr);
335 if(1 == lr) hitPhi = wphi + dphi;
336 else hitPhi = wphi - dphi;
337 if(hitPhi < 0) hitPhi +=
PI2;
338 else if(hitPhi >
PI2) hitPhi -=
PI2;
340 if(zhit < m_zrange[lay][0]) iEnd = 0;
353 m_tuple[iEP]->write();
362 m_tuple[
NEP]->write();
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);
370 m_gr[iEP]->SetPoint(m_npoint[iEP], hitPhi, resi);
381 Gaudi::svcLocator() -> service(
"MessageSvc",
msgSvc);
382 MsgStream log(
msgSvc,
"ResiAlign");
383 log << MSG::INFO <<
"ResiAlign::updateConst()" << endreq;
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 };
397 TCanvas c1(
"c1",
"c1", 10, 10, 700, 500);
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);
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);
411 err[0] = fResPhi->GetParError(0);
412 err[1] = fResPhi->GetParError(1);
413 err[2] = fResPhi->GetParError(2);
417 rz = par[0] / rLayer[iEP];
420 if (7==iEP || 15==iEP) {
434 alignPar->
setErrRz(iEP, err[0]/rLayer[iEP]);
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;
449 val = par[0] + par[1]*
sin(x[0]) + par[2]*
cos(x[0]);
double sin(const BesAngle a)
double cos(const BesAngle a)
virtual const MdcGeoLayer *const Layer(unsigned id)=0
virtual const MdcGeoWire *const Wire(unsigned id)=0
double getResiExcLR() const
double getResiIncLR() const
MdcAliRecHit * getRecHit(int index) const
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
HepPoint3D Forward(void) const
HepPoint3D Backward(void) const
static Double_t funResi(double *x, double *par)
void initialize(TObjArray *hlist, IMdcGeomSvc *mdcGeomSvc, IMdcCalibFunSvc *mdcFunSvc)
bool fillHist(MdcAliEvent *event)
void updateConst(MdcAlignPar *alignPar)
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)
int getEpId(int lay, int iEnd)