BOSS 7.0.2
BESIII Offline Software System
Loading...
Searching...
No Matches
MdcxTrackFinder.cxx
Go to the documentation of this file.
1//--------------------------------------------------------------------------
2// File and Version Information:
3// $Id: MdcxTrackFinder.cxx,v 1.49 2012/07/20 05:48:16 zhangy Exp $
4//
5// Description:
6// Class MdcxTrackFinder
7//
8// Environment:
9// Software developed for the BaBar Detector at the SLAC B-Factory.
10//
11// Author List:
12// Steve Wagner Original Author
13// Zhang Yao([email protected])
14//
15// Copyright Information:
16// Copyright (C) 1997 BEPCII
17//
18// History:
19// Migration for BESIII MDC
20//
21//------------------------------------------------------------------------
22
23//-----------------------
24// This Class's Header --
25//-----------------------
26#include "MdcxReco/MdcxTrackFinder.h"
27
28//-------------
29// C Headers --
30//-------------
31#include <assert.h>
32#include <time.h>
33#include <math.h>
34
35//---------------
36// C++ Headers --
37//---------------
38#include <iostream>
39
40//-------------------------------
41// Collaborating Class Headers --
42//-------------------------------
43#include "GaudiKernel/MsgStream.h"
44#include "GaudiKernel/AlgFactory.h"
45#include "BField/BField.h"
46#include "MdcRawEvent/MdcDigi.h"
47#include "MdcGeom/MdcDetector.h"
48#include "EventModel/EventHeader.h"
49
50#include "MdcxReco/MdcxHit.h"
51#include "MdcxReco/MdcxHits.h"
52#include "MdcxReco/MdcxFindSegs.h"
53#include "MdcxReco/MdcxSeg.h"
54#include "MdcxReco/MdcxHel.h"
55#include "MdcxReco/MdcxFittedHel.h"
56#include "MdcxReco/MdcxFindTracks.h"
57#include "MdcxReco/Mdcxprobab.h"
58#include "CLHEP/Alist/AIterator.h"
59#include "TrkBase/TrkExchangePar.h"
60#include "TrkBase/TrkErrCode.h"
61#include "TrkBase/TrkRecoTrk.h"
62#include "TrkBase/TrkFit.h"
63#include "TrkBase/TrkHitList.h"
64#include "TrkFitter/TrkHelixMaker.h"
65#include "TrkFitter/TrkLineMaker.h"
66#include "MdcxReco/MdcxAddHits.h"
67#include "MdcxReco/MdcxMergeDups.h"
68#include "TrkBase/TrkFitStatus.h"
69#include "MdcRecEvent/RecMdcTrack.h"
70#include "MdcRecEvent/RecMdcHit.h"
71#include "ReconEvent/ReconEvent.h"
72#include "Identifier/MdcID.h"
73#include "Identifier/Identifier.h"
74#include "EvTimeEvent/RecEsTime.h"
75#include "McTruth/McParticle.h"
76#include "McTruth/MdcMcHit.h"
77#include "MdcCalibFunSvc/IMdcCalibFunSvc.h"
78#include "RawDataProviderSvc/IRawDataProviderSvc.h"
79#include "RawDataProviderSvc/MdcRawDataProvider.h"
80#include "RawEvent/RawDataUtil.h"
81#include "MdcRecoUtil/Pdt.h"
82#include "MdcTrkRecon/MdcTrack.h"
83//debug
84#include "TrkBase/TrkHotList.h"
85#include "GaudiKernel/INTupleSvc.h"
86#include "GaudiKernel/NTuple.h"
87#include "MdcxReco/MdcxHistItem.h"
88#include "MdcGeom/Constants.h"
89#include "MdcGeom/MdcTrkReconCut.h"
90#include "MdcxReco/MdcxSegPatterns.h"
91#include "TrigEvent/TrigData.h"
92
93#include "BesTimerSvc/IBesTimerSvc.h"
94#include "BesTimerSvc/BesTimerSvc.h"
95#include "BesTimerSvc/BesTimer.h"
96
97
98
99using std::cout;
100using std::endl;
101
102extern double MdcTrkReconCut_helix_fit[43];
103//----------------
104// Constructors --
105//----------------
106
107MdcxTrackFinder::MdcxTrackFinder(const std::string& name, ISvcLocator* pSvcLocator) :
108 Algorithm(name, pSvcLocator), m_mdcCalibFunSvc(0)
109{
110 //input
111 declareProperty("pdtFile", m_pdtFile = "pdt.table");
112 //debug control
113 declareProperty("debug", m_debug= 0);
114 declareProperty("hist", m_hist= 0);
115 declareProperty("mcHist", m_mcHist = false);
116 //cuts and control
117 declareProperty("cresol", m_cresol = 0.013);
118
119 declareProperty("getDigiFlag", m_getDigiFlag = 0);
120 declareProperty("maxMdcDigi", m_maxMdcDigi= 0);
121 declareProperty("keepBadTdc", m_keepBadTdc= 0);
122 declareProperty("dropHot", m_dropHot= 0);
123 declareProperty("keepUnmatch", m_keepUnmatch= 0);
124 declareProperty("salvageTrk", m_salvageTrk = false);
125 declareProperty("dropMultiHotInLayer",m_dropMultiHotInLayer = false);
126 declareProperty("dropTrkPt", m_dropTrkPt = -999.);
127 declareProperty("d0Cut", m_d0Cut = 999.);
128 declareProperty("z0Cut", m_z0Cut = 999.);
129
130 declareProperty("minMdcDigi", m_minMdcDigi = 0);
131 declareProperty("countPropTime",m_countPropTime = true);
132 declareProperty("addHitCut", m_addHitCut = 5.);
133 declareProperty("dropHitsSigma",m_dropHitsSigma);
134 declareProperty("helixFitCut", m_helixFitCut);
135 declareProperty("minTrkProb", m_minTrkProb = 0.01);
136 declareProperty("csmax4", m_csmax4 = 50.);
137 declareProperty("csmax3", m_csmax3 = 1.);
138 declareProperty("helixFitSigma",m_helixFitSigma= 5.);
139 declareProperty("maxRcsInAddSeg",m_maxRcsInAddSeg= 50.);
140 declareProperty("nSigAddHitTrk", m_nSigAddHitTrk= 5.);
141 declareProperty("maxProca", m_maxProca= 0.6);
142 declareProperty("doSag", m_doSag= false);
143 declareProperty("lineFit", m_lineFit = false);
144 //declareProperty("cosmicFit", m_cosmicFit= false);
145}
146
147//--------------
148// Destructor --
149//--------------
151 delete m_bfield;
152}
153
155 // Get Mdc Detector Geometry
156 m_gm = MdcDetector::instance(m_doSag);
157 if(NULL == m_gm) return StatusCode::FAILURE;
159
160 return StatusCode::SUCCESS;
161}
162//--------------
163// Operations --
164//--------------
166 MsgStream log(msgSvc(), name());
167 log << MSG::INFO << "in initialize()" << endreq;
168
169 t_nTkTot = 0;
170 for (int i=0;i<20;i++) t_nTkNum[i]=0;
171
172 //m_flags.readPar(m_paramFile);
173#ifdef MDCXTIMEDEBUG
174 StatusCode tsc = service( "BesTimerSvc", m_timersvc);
175 if( tsc.isFailure() ) {
176 log << MSG::WARNING << name() << ": Unable to locate BesTimer Service" << endreq;
177 return StatusCode::FAILURE;
178 }
179 m_timer[0] = m_timersvc->addItem("Execution");
180 m_timer[0]->propName("Execution");
181#endif
182
183 if (m_helixFitCut.size() == 43){
184 for(int i=0; i<43; i++){
185 //MdcTrkReconCut_helix_fit[i] = m_helixFitCut[i];
186 TrkHelixFitter::nSigmaCut[i] = m_helixFitCut[i];
187 }
188 }
189 MdcxParameters::debug = m_debug;
190 MdcxParameters::minTrkProb = m_minTrkProb;
191 MdcxParameters::csmax4 = m_csmax4;
192 MdcxParameters::csmax3 = m_csmax3;
193 MdcxParameters::helixFitSigma = m_helixFitSigma;
194 MdcxParameters::maxRcsInAddSeg= m_maxRcsInAddSeg;
195 MdcxParameters::nSigAddHitTrk = m_nSigAddHitTrk;
196 MdcxParameters::maxProca = m_maxProca;
197 TrkHelixFitter::m_debug = (m_debug>7);
198 Pdt::readMCppTable(m_pdtFile);
199 MdcxFittedHel::debug = m_debug;
200
201
202 // Get MdcCalibFunSvc
203 IMdcCalibFunSvc* imdcCalibSvc;
204 StatusCode sc = service ("MdcCalibFunSvc", imdcCalibSvc);
205 m_mdcCalibFunSvc = dynamic_cast<MdcCalibFunSvc*>(imdcCalibSvc);
206 if ( sc.isFailure() ){
207 log << MSG::FATAL << "Could not load MdcCalibFunSvc!" << endreq;
208 return StatusCode::FAILURE;
209 }
210 MdcxHit::setMdcCalibFunSvc(m_mdcCalibFunSvc);
211 MdcxHit::setCountPropTime(m_countPropTime);
212
213 // Get RawDataProviderSvc
214 IRawDataProviderSvc* iRawDataProvider;
215 sc = service ("RawDataProviderSvc", iRawDataProvider);
216 if ( sc.isFailure() ){
217 log << MSG::FATAL << "Could not load RawDataProviderSvc!" << endreq;
218 return StatusCode::FAILURE;
219 }
220 m_rawDataProviderSvc = dynamic_cast<RawDataProviderSvc*>(iRawDataProvider);
221
222
223 //Initailize magnetic filed
224 sc = service ("MagneticFieldSvc",m_pIMF);
225 if(sc!=StatusCode::SUCCESS) {
226 log << MSG::ERROR << "Unable to open Magnetic field service"<<endreq;
227 }
228 m_bfield = new BField(m_pIMF);
229 log << MSG::INFO << "field z = "<<m_bfield->bFieldNominal()<< endreq;
230 m_context = new TrkContextEv(m_bfield);
231
232 if (m_hist) {bookNTuple();}
233 if (m_dropHitsSigma.size()==43){
234 for (int ii=0;ii<43;ii++) {
235 MdcxParameters::dropHitsSigma[ii]=m_dropHitsSigma[ii];
236 }
237 }
238
239 return StatusCode::SUCCESS;
240}
241
243
244
245 b_saveEvent=false;
246 setFilterPassed(b_saveEvent);
247#ifdef MDCXTIMEDEBUG
248 m_timer[0]->start();
249#endif
250 MsgStream log(msgSvc(), name());
251 log << MSG::INFO << "in execute()" << endreq;
252 StatusCode sc;
253
254 nTk = 0; t_nTdsTk = 0; //yzhang for fill
255 //------------------------------------
256 // Get event No.
257 //------------------------------------
258 SmartDataPtr<Event::EventHeader> evtHead(eventSvc(),"/Event/EventHeader");
259 if (!evtHead) {
260 log << MSG::FATAL<< "Could not retrieve event header" << endreq;
261 return StatusCode::FAILURE;
262 }
263 m_eventNo = evtHead->eventNumber();
264 //std::cout << "x evt: "<<evtHead->runNumber()<<" " << m_eventNo<< std::endl;
265 long t_evtNo = m_eventNo;
266 g_eventNo = m_eventNo;
267 //if (t_evtNo % 1000 == 0) std::cout << "x evt: " << t_evtNo << std::endl;
268 IDataManagerSvc *dataManSvc;
269 DataObject *aTrackCol;
270 DataObject *aHitCol;
271 if (!m_salvageTrk) {
272 SmartIF<IDataManagerSvc> dataManSvc(eventSvc());
273 eventSvc()->findObject("/Event/Recon/RecMdcTrackCol",aTrackCol);
274 if(aTrackCol != NULL) {
275 dataManSvc->clearSubTree("/Event/Recon/RecMdcTrackCol");
276 eventSvc()->unregisterObject("/Event/Recon/RecMdcTrackCol");
277 }
278 eventSvc()->findObject("/Event/Recon/RecMdcHitCol",aHitCol);
279 if(aHitCol != NULL) {
280 dataManSvc->clearSubTree("/Event/Recon/RecMdcHitCol");
281 eventSvc()->unregisterObject("/Event/Recon/RecMdcHitCol");
282 }
283 }
284
285 //------------------------------------
286 // Initialize track collection in TDS
287 //------------------------------------
288 DataObject *aReconEvent;
289 eventSvc()->findObject("/Event/Recon",aReconEvent);
290 if (aReconEvent==NULL) {
291 aReconEvent = new ReconEvent();
292 sc = eventSvc()->registerObject("/Event/Recon",aReconEvent);
293 if(sc != StatusCode::SUCCESS) {
294 log << MSG::FATAL << "Could not register ReconEvent" <<endreq;
295 return StatusCode::FAILURE;
296 }
297 }
298 RecMdcTrackCol* trackList;
299 eventSvc()->findObject("/Event/Recon/RecMdcTrackCol",aTrackCol);
300 if (aTrackCol) {
301 trackList = dynamic_cast<RecMdcTrackCol*> (aTrackCol);
302 }else{
303 trackList = new RecMdcTrackCol;
304 sc = eventSvc()->registerObject(EventModel::Recon::RecMdcTrackCol, trackList);
305 if(!sc.isSuccess()) {
306 log << MSG::FATAL << " Could not register RecMdcTrack collection" <<endreq;
307 return StatusCode::FAILURE;
308 }
309 }
310 RecMdcHitCol* hitList;
311 eventSvc()->findObject("/Event/Recon/RecMdcHitCol",aHitCol);
312 if (aHitCol) {
313 hitList = dynamic_cast<RecMdcHitCol*> (aHitCol);
314 }else{
315 hitList = new RecMdcHitCol;
316 sc = eventSvc()->registerObject(EventModel::Recon::RecMdcHitCol, hitList);
317 if(!sc.isSuccess()) {
318 log << MSG::FATAL << " Could not register RecMdcHit collection" <<endreq;
319 return StatusCode::FAILURE;
320 }
321 }
322
323 //------------------------------------
324 // Initialize hit collection in TDS
325 //------------------------------------
326 DataObject *pnode = 0;
327 sc = eventSvc()->retrieveObject("/Event/Hit",pnode);
328 if(!sc.isSuccess()) {
329 pnode = new DataObject;
330 sc = eventSvc()->registerObject("/Event/Hit",pnode);
331 if(!sc.isSuccess()) {
332 log << MSG::FATAL << " Could not register /Event/Hit branch " <<endreq;
333 return StatusCode::FAILURE;
334 }
335 }
336 SmartDataPtr<MdcHitCol> m_hitCol(eventSvc(), "/Event/Hit/MdcHitCol");
337 if (!m_hitCol){
338 m_hitCol= new MdcHitCol;
339 sc = eventSvc()->registerObject("/Event/Hit/MdcHitCol",m_hitCol);
340 if(!sc.isSuccess()) {
341 log << MSG::FATAL << " Could not register hit collection" <<endreq;
342 return StatusCode::FAILURE;
343 }
344 }
345
346
347 //------------------------------------
348 // Get bunch time t0 (ns) and timing
349 //------------------------------------
350 m_bunchT0 = -999.;
351 SmartDataPtr<RecEsTimeCol> aevtimeCol(eventSvc(),"/Event/Recon/RecEsTimeCol");
352 if (!aevtimeCol || aevtimeCol->size()==0) {
353 log << MSG::WARNING<< "evt "<<m_eventNo<<" Could not find RecEsTimeCol"<< endreq;
354 return StatusCode::SUCCESS;
355 }
356
357 RecEsTimeCol::iterator iter_evt = aevtimeCol->begin();
358 for(; iter_evt!=aevtimeCol->end(); iter_evt++){
359 m_bunchT0 = (*iter_evt)->getTest();
360 m_t0Stat = (*iter_evt)->getStat();
361 if ((m_t0Stat==0) || (m_bunchT0 < 0.) || (m_bunchT0 > 9999.0) ){
362 log << MSG::WARNING << "Skip evt:"<<m_eventNo<< " by t0 = "<<m_bunchT0 << endreq;
363 //return StatusCode::SUCCESS;
364 }
365 }
366 if(m_debug>1) std::cout<<name()<<" t0 "<<m_bunchT0<<" t0Stat "<<m_t0Stat<<std::endl;
367 int trigtiming=-10;
368 SmartDataPtr<TrigData> trigData(eventSvc(),"/Event/Trig/TrigData");
369 if(trigData){
370 log << MSG::INFO <<"Trigger conditions 0--43:"<<endreq;
371 for(int i = 0; i < 48; i++) {
372 log << MSG::INFO << trigData->getTrigCondName(i)<<" ---- "<<trigData->getTrigCondition(i)<< endreq;
373 }
374 for(int i = 0; i < 16; i++) log << MSG::INFO << "Trigger channel "<< i << ": " << trigData->getTrigChannel(i) << endreq;
375 m_timing=trigData->getTimingType();
376 //cout<<"-----------------trigger timing type-----------------------: "<<trigtiming<<endl;
377 log << MSG::INFO <<"Tigger Timing type: "<< trigtiming << endreq;
378 }
379
380 //------------------------------------
381 // Initialize MdcxHits
382 //------------------------------------
383 m_mdcxHits.reset();
384 uint32_t getDigiFlag = 0;
385 getDigiFlag += m_maxMdcDigi;
386 if(m_dropHot||m_salvageTrk) getDigiFlag |= MdcRawDataProvider::b_dropHot;
387 if(m_keepBadTdc) getDigiFlag |= MdcRawDataProvider::b_keepBadTdc;
388 if(m_keepUnmatch) getDigiFlag |= MdcRawDataProvider::b_keepUnmatch;
389 mdcDigiVec = m_rawDataProviderSvc->getMdcDigiVec(getDigiFlag);
390 t_nDigi = mdcDigiVec.size();
391
392 if (0 == t_nDigi ){
393 log << MSG::WARNING << " No hits in MdcDigiVec" << endreq;
394 return StatusCode::SUCCESS;
395 }
396
397
398 // fill Mc truth
399 //if(m_hist) fillMcTruth();
400
401 //skip event by hit numbe
402 if (t_nDigi < m_minMdcDigi){
403 log << MSG::WARNING << " Skip this event for MdcDigiVec.size() < "<<m_minMdcDigi << endreq;
404 return StatusCode::SUCCESS;
405 }
406 m_mdcxHits.create(mdcDigiVec, m_bunchT0, m_cresol);
407 const HepAList<MdcxHit>& dchitlist = m_mdcxHits.GetMdcxHitList();
408
409 //--------------------------------------------
410 // Make segments (MdcxSeg's) out of MdcxHit's
411 //--------------------------------------------
412 MdcxFindSegs dcsegs(dchitlist,m_debug);
413 const HepAList<MdcxSeg>& seglist = dcsegs.GetMdcxSeglist();
414 if(m_debug > 1 ){ dumpMdcxSegs(seglist);}
415 //if(m_hist){ fillMdcxSegs(seglist);}
416
417 //--------------------------------------------
418 // Make tracks (MdcxFittedHel's) out of MdcxSeg's
419 //--------------------------------------------
420 MdcxFindTracks dctrks(seglist,m_debug);
422 if(m_debug>1) dumpTrackList(firsttrkl);
423
424 //if(m_hist){ fillTrkl(firsttrkl);}
425
426 //if(m_debug>1){
427 // std::cout << "dchitlist after find tracks before MergeDups, nhits=" << dchitlist.length() << std::endl;
428 // for (int ii = 0; ii < dchitlist.length(); ii++) {
429 // dchitlist[ii]->print(std::cout, ii);
430 // }
431 // std::cout<<std::endl;
432 //}
433 MdcxMergeDups dcmergeem(firsttrkl,m_debug);
435
436 //if (m_debug > 1 ){
437 // cout << "MdcxTrackFinder: after MergeDups, have "
438 // << trkl.length() << " track(s). nhits=" << dchitlist.length() << endl;
439 // for (int ii = 0; ii < dchitlist.length(); ii++) {
440 // dchitlist[ii]->print(std::cout, ii);
441 // }
442 // std::cout<<std::endl;
443 //}
444
445 //---------------------------------------------------------
446 // Put my tracks into official fitter and store to TDS
447 //----------------------------------------------------------
448
449
450 sc = FitMdcxTrack(trkl, dchitlist, m_hitCol, trackList, hitList);
451 if (!sc.isSuccess()) {return StatusCode::SUCCESS;}
452 t_nTdsTk = trackList->size();
453
454 t_nTkTot += trackList->size();
455 if(t_nTdsTk<20) t_nTkNum[t_nTdsTk]++;
456
457#ifdef MDCXTIMEDEBUG
458 m_timer[0]->stop();
459 m_xt4time = m_timer[0]->elapsed();
460 m_xtupleEvt->write();
461#endif
462
463 if(m_hist) fillEvent();
464 if (m_debug > 0) {
465 DataObject* pNode;
466 eventSvc()->retrieveObject("/Event/Recon/RecMdcTrackCol",pNode);
467 RecMdcTrackCol *tmpTrackCol = dynamic_cast<RecMdcTrackCol*> (pNode);
468 eventSvc()->retrieveObject("/Event/Recon/RecMdcHitCol",pNode);
469 int nTdsTk = 0;
470 if(tmpTrackCol) nTdsTk = tmpTrackCol->size();
471
472 //if (t_evtNo % 1000 == 0) {
473 std::cout<< "MdcxTrackFinder: evtNo "<< m_eventNo << " t0="<<m_bunchT0
474 <<" Found " <<trkl.length()
475 <<" keep "<< t_nTdsTk
476 <<" finialy keep "<< nTdsTk;
477
478 int ndelete =0; trkl.length() - trackList->size();
479 if( ndelete>0 ) std::cout <<" delete "<< ndelete;
480 std::cout <<" track(s)" <<endl;
481 //}
482
483 if(m_debug>1)dumpTdsTrack(tmpTrackCol);
484 if(m_debug>1)dumpTrack(tmpTrackCol);
485 //dumpTdsHits(tmpHitCol);
486
487 }
488 if((trackList->size()!=4) ) b_saveEvent = true;
489 setFilterPassed(b_saveEvent);
490 return StatusCode::SUCCESS;
491}
492
494 MsgStream log(msgSvc(), name());
495 log << MSG::INFO << "in finalize()" << endreq;
496
497 std::cout<< " --after " << name() << " keep " << t_nTkTot<< " tracks "<<std::endl;
498 for(int i=0; i<20; i++){
499 if(t_nTkNum[i]>0)std::cout<< " nTk="<< i << " "<< t_nTkNum[i]<< std::endl;
500 }
501
502 //tot evtNo, trkNum
503 return StatusCode::SUCCESS;
504}
505
506
507void MdcxTrackFinder::MdcxHitsToHots(MdcxHel& mdcxHelix, const HepAList<MdcxHit>& mdcxHits,
508 TrkHitList* m_trkHitList, MdcHitCol* mdcHitCol) {
509
510 if ( 0 == mdcxHits.length() ) return;
511
512 int ihits = 0;
513 while(mdcxHits[ihits]) {
514 int ambig = 0;//=mdcxHelix.Doca_Samb();//yzhang delete 2011-10-13
515 const MdcHit* newhit = mdcxHits[ihits]->getMdcHit();
516 if ( 0 == newhit ) {
517 const MdcDigi* theDigi = mdcxHits[ihits]->getDigi();
518 int layer = MdcID::layer(mdcxHits[ihits]->getDigi()->identify());
519 int wire = MdcID::wire(mdcxHits[ihits]->getDigi()->identify());
520 m_digiMap[layer][wire] = mdcxHits[ihits]->getDigi();
521 MdcHit *thehit = new MdcHit(theDigi, m_gm);
522 thehit->setCalibSvc(m_mdcCalibFunSvc);
523 thehit->setCountPropTime(m_countPropTime);
524 //thehit->setCosmicFit(m_cosmicFit);
525
526 mdcHitCol->push_back(thehit);
527 newhit = thehit;
528 }
529 MdcRecoHitOnTrack temp(*newhit, ambig, m_bunchT0);//m_bunchT0 nano second here
530 MdcHitOnTrack* newhot = &temp;
531 //double fltLen = mdcxHelix.Doca_FLen();
532 //std::cout<<" fltLen-- "<<fltLen<<" "<<std::endl;
533 //newhot->setFltLen(fltLen);
534 //newhot->print(std::cout);
535 //std::cout<< __FILE__ << " ambig " << ambig << " append "<<std::endl;
536
537 // Store MdcxHits to TrkHitList
538 m_trkHitList->appendHot(newhot); //yzhang TEMP FIXME
539 newhot->setUsability(true);//yzhang add 2011-10-12 TEMP
540 newhot->setActivity(true);//yzhang add 2011-10-12 TEMP
541
542 ihits++;
543 }
544}
545
546StatusCode MdcxTrackFinder::FitMdcxTrack(HepAList<MdcxFittedHel>& trkl,
547 const HepAList<MdcxHit>& dchitlist, MdcHitCol* m_hitCol,
548 RecMdcTrackCol* trackList, RecMdcHitCol* hitList){
549 StatusCode sc;
550 MsgStream log(msgSvc(), name());
551
552 //--Add Hit to MdcxFittedHel
553 //if(m_debug>1){
554 // std::cout << "dchitlist before addHits nhits=" << dchitlist.length() << std::endl;
555 // for (int ii = 0; ii < dchitlist.length(); ii++) {
556 // dchitlist[ii]->print(std::cout, ii);
557 // }
558 // std::cout<<std::endl;
559 //}
560 MdcxAddHits dcaddem(trkl, dchitlist, m_addHitCut);
561 //if (m_debug > 1){
562 // cout << "MdcxTrackFinder: after AddHits, have "
563 // << trkl.length() << " track(s). nhits=" << dchitlist.length() << endl;
564 // for (int ii = 0; ii < dchitlist.length(); ii++) {
565 // dchitlist[ii]->print(std::cout, ii);
566 // }
567 // std::cout<<std::endl;
568 //}
569
570 TrkLineMaker linefactory;
571 TrkHelixMaker helixfactory;
572
573
574 int tkLen = trkl.length();//FIXME
575 for (int kk=0; kk< tkLen; kk++){
576 const HepAList<MdcxHit>& xhits = trkl[kk]->XHitList();
577 if(m_debug>2){
578 std::cout<<__FILE__<<" FitMdcxTrack "<<kk<< std::endl;
579 for(int i=0; i<xhits.length(); i++){ xhits[i]->print(std::cout); }
580 std::cout<<std::endl;
581 }
582
583 if(m_debug>2) std::cout<<" before add hits nhits="<<xhits.length()<< std::endl;
584 HepAList<MdcxHit> xass = dcaddem.GetAssociates(kk);
585 if(m_debug>2){
586 std::cout<<" after,add "<<xass.length()<<" hits"<<std::endl;
587 }
588
589 MdcxFittedHel mdcxHelix = *trkl[kk];
590 double thechisq = mdcxHelix.Chisq();
591 TrkExchangePar tt(-mdcxHelix.D0(),mdcxHelix.Phi0(),mdcxHelix.Omega(),-mdcxHelix.Z0(),-mdcxHelix.Tanl());
592 TrkRecoTrk* aTrack;
593 int nparm;
594
595 if (m_lineFit) {
596 /// m_bunchT0 in "ns" here, but "second" in TrkLineMaker&TrkHelixMaker
597 aTrack = linefactory.makeTrack(tt, thechisq, *m_context, m_bunchT0*1.e-9);
598 nparm = 4;
599 linefactory.setFlipAndDrop(*aTrack, true, true);
600 } else {
601 aTrack = helixfactory.makeTrack(tt, thechisq, *m_context, m_bunchT0*1.e-9);
602 nparm = 5;
603 helixfactory.setFlipAndDrop(*aTrack, true, true);
604 }
605
606 TrkHitList* m_trkHitList = aTrack->hits();
607 if (0 == m_trkHitList) {
608 delete aTrack;
609 aTrack = NULL;
610 continue;
611 }
612
613
614 MdcxHitsToHots(mdcxHelix, xhits, m_trkHitList, m_hitCol);
615 //std::cout<<"xhits---------------------"<<std::endl;//debug
616 //m_trkHitList->hotList().printAll(cout);//debug
617 MdcxHitsToHots(mdcxHelix, xass, m_trkHitList, m_hitCol);
618 //std::cout<<"xass----------------------"<<std::endl;//debug
619 //m_trkHitList->hotList().printAll(cout);//debug
620 //std::cout<<"size "<<m_trkHitList->hotList().nHit()<<std::endl;
621 //int beforDrop = m_trkHitList->hotList().nHit();
622 if(m_dropMultiHotInLayer) dropMultiHotInLayer(m_trkHitList);//yzhang debug FIXME
623 //int afterDrop = m_trkHitList->hotList().nHit();
624 //std::cout<<"drop "<<beforDrop-afterDrop<<" keep:"<<afterDrop<<::endl;
625 if(m_debug>1){ std::cout<< "== put to official fitter "<<std::endl;}
626 TrkErrCode err = m_trkHitList->fit();
627
628 if(m_debug>1){
629 std::cout<< "== after official fitter "<<std::endl;
630 aTrack->printAll(std::cout);
631 }
632
633 const TrkFit* theFit = aTrack->fitResult();
634 float rcs = 10000.0;
635
636 if (theFit) {
637 int ndof = theFit->nActive()-nparm;
638 if (ndof > 0) rcs = theFit->chisq()/float(ndof);
639 if (m_debug>1) {
640 if (4 == nparm) cout << " TrkLineMaker";
641 else cout << " TrkHelixMaker";
642 cout << " trkNo. " << kk << " success " << err.success() << " rcs " << rcs
643 << " chi2 " << theFit->chisq() << " nactive " << theFit->nActive() << endl;
644 }
645 }
646
647 if ( (1 == err.success()) && (rcs < 20.0) ) {
648 if(m_debug>1) std::cout<<"== fit success "<<std::endl;//yzhang debug
649 if (4 == nparm) {
650 linefactory.setFlipAndDrop(*aTrack, false, false);
651 } else {
652 helixfactory.setFlipAndDrop(*aTrack, false, false);
653 }
654 //-------------Stick the found tracks into the list in RecMdcTrackCol--------
655 if (m_debug>1) { cout << "MdcxTrackFinder: accept a track " << endl; }
656 // update history
657 aTrack->status()->addHistory(err,name().c_str()); //yzhang FIXME
658 store(aTrack, trackList, hitList);//aTrack have been deleted
659 } else if ( (2 == err.success()) && (rcs < 150.0) ) { //////// zoujh
660 if(m_debug > 1) std::cout<<"== fit success = 2, refit now"<<std::endl;//yzhang debug
661 int nrefit = 0;
662 while (nrefit++ < 5) {
663 if (m_debug>1) std::cout << "refit time " << nrefit << std::endl;
664 err = m_trkHitList->fit();
665 if (err.success() == 1) break;
666 }
667 if (err.success() == 1) {
668 if (4 == nparm) {
669 linefactory.setFlipAndDrop(*aTrack, false, false);
670 } else {
671 helixfactory.setFlipAndDrop(*aTrack, false, false);
672 }
673 //-------------Stick the found tracks into the list in RecMdcTrackCol--------
674 //if (m_debug>1) { cout << "MdcxTrackFinder: accept a track and store to TDS" << endl; }
675 // update history
676 aTrack->status()->addHistory(err,name().c_str()); //yzhang FIXME
677 store(aTrack, trackList, hitList);//aTrack have been deleted
678 } //////////////////////////zoujh
679 } else {
680 if (m_debug >1) {
681 std::cout<<"== fit faild "<<std::endl;//yzhang debug
682 err.print(cout);
683 cout << endl;
684 }
685 delete aTrack;
686 aTrack = NULL;
687 //---------------------------------------
688 // Fit no good; try a better input helix
689 //---------------------------------------
690 if(m_debug>1) std::cout<<"== Fit no good; try a better input helix"<<std::endl;//yzhang debug
691 mdcxHelix.Grow(*trkl[kk],xass);
692 mdcxHelix.VaryRes();
693 mdcxHelix.SetChiDofBail(1500);//yzhang add 2009-11-03
694 int fail = mdcxHelix.ReFit();
695 if(m_debug>1)std::cout<<__FILE__<<" refit fail:"<<fail<< std::endl;
696 if (!mdcxHelix.Fail()) {
697 const HepAList<MdcxHit>& bxhits = mdcxHelix.XHitList();
698 thechisq = mdcxHelix.Chisq();
699 TrkExchangePar tb(mdcxHelix.D0(),mdcxHelix.Phi0(),mdcxHelix.Omega(),mdcxHelix.Z0(),mdcxHelix.Tanl());
700 TrkRecoTrk* bTrack;
701 if (4 == nparm){
702 bTrack = linefactory.makeTrack(tb,thechisq,*m_context,m_bunchT0*1.e-9);
703 linefactory.setFlipAndDrop(*bTrack, false, false);
704 }else{
705 bTrack = helixfactory.makeTrack(tb,thechisq,*m_context,m_bunchT0*1.e-9);
706 helixfactory.setFlipAndDrop(*bTrack, false, false);
707 }
708 TrkHitList* bhits = bTrack->hits();
709 if (0 == bhits){delete bTrack; bTrack = NULL; continue;}
710
711 MdcxHitsToHots(mdcxHelix, bxhits, bhits, m_hitCol);
712 TrkErrCode berr = bhits->fit();
713 const TrkFit* bFit = bTrack->fitResult();
714 rcs=10000.0;
715 if (bFit) {
716 int ndof = bFit->nActive() - nparm;
717 if (ndof > 0) rcs = bFit->chisq()/float(ndof);
718 if (m_debug >1) {
719 if (4 == nparm) cout << " TrkLineMaker";
720 else cout << " TrkHelixMaker";
721 cout << " success trkNo. " << kk << " status " << berr.success() << " rcs " << rcs
722 << " chi2 " << bFit->chisq() << " nactive "<< bFit->nActive() << endl;
723 }
724 }
725 if ( ( 1 == berr.success() ) && ( rcs < 50.0 ) ) {
726 // update history
727 bTrack->status()->addHistory(berr,name().c_str());//yzhang FIXME
728 if (m_debug>1) {
729 cout << "MdcxTrackFinder: accept b track and store to TDS" << endl;
730 bTrack->printAll(cout);
731 }
732 store(bTrack, trackList, hitList);//bTrack have been deleted
733 } else {
734 if (m_debug>1) {
735 cout<< " fit failed "<<endl;
736 berr.print(cout);
737 cout << endl;
738 }
739 if (bTrack!=NULL) { delete bTrack; bTrack = NULL; }
740 }
741 }else{
742 //cout<< " grow and refit failed "<<endl;
743 }
744 }
745 }
746 if(m_debug >1) dumpTdsTrack(trackList);
747 return StatusCode::SUCCESS;
748
749}// end of FitMdcxTrack
750
751void MdcxTrackFinder::store(TrkRecoTrk* aTrack, RecMdcTrackCol *trackList,
752 RecMdcHitCol *hitList) {
753 assert (aTrack != NULL);
754 nTk++;
755 int trackId = trackList->size();
756 TrkExchangePar helix = aTrack->fitResult()->helix(0.);
757 if(m_dropTrkPt>0. && (aTrack->fitResult()->pt()<m_dropTrkPt)) {
758 if(m_debug>1) std::cout<<"MdcxTrackFinder delete track by pt "
759 <<aTrack->fitResult()->pt()<<"<ptCut "<<m_dropTrkPt << std::endl;
760 return;
761 }
762
763 if( ( (fabs(helix.d0())>m_d0Cut) ||( fabs(helix.z0())>m_z0Cut) ) ){
764 if(m_debug>1) std::cout<< name()<<" delete track by d0 "<<helix.d0()<<">d0Cut "<<m_d0Cut
765 <<" or z0 "<<helix.z0()<<" >z0Cut "<<m_z0Cut << std::endl;
766 return;
767 }
768
769 if(m_hist) fillTrack(aTrack);
770 MdcTrack mdcTrack(aTrack);//aTrack have been deleted in ~MdcTrack()
771 //tkStat: 0,PatRec 1,MdcxReco 2,Tsf 3,CurlFinder -1,Combined cosmic
772 int tkStat = 1;
773 int nHitbefore = hitList->size();
774
775 mdcTrack.storeTrack(trackId, trackList, hitList, tkStat);
776 int nHitAfter = hitList->size();
777 if (nHitAfter - nHitbefore <10 ) setFilterPassed(true);
778}
779
780void MdcxTrackFinder::dumpTrack(RecMdcTrackCol* trackList){
781 RecMdcTrackCol::iterator i_tk = trackList->begin();
782 for (;i_tk != trackList->end(); i_tk++) {
783 printTrack(*(i_tk));
784 }
785}// dumpTrack
786
787void MdcxTrackFinder::printTrack(RecMdcTrack* tk){
788 //yzhang debug
789 std::cout<< " MdcTrack Id:"<<tk->trackId() <<" q:"<< tk->charge()<< std::endl;
790 std::cout<< "dr Fi0 Cpa Dz Tanl Chi2 Ndf nSt FiTerm poca" << std::endl;
791 std::cout<<"(" <<setw(5) << tk->helix(0)<<","<< setw(5) << tk->helix(1)<<"," << setw(5) << tk->helix(2) <<","
792 << setw(5) << tk->helix(3) << ","<< setw(5) << tk->helix(4) <<")"
793 << setw(5) << tk->chi2() << setw(4) << tk->ndof()
794 << setw(4) << tk->getNhits() << setw(4) << tk->nster()
795 << setw(5) << tk->getFiTerm() <<tk->poca()<<std::endl;
796 std::cout<< " ErrMat "<<tk->err() << std::endl;
797
798 std::cout<< "hitId tkId (l,w) fltLen lr dt ddl tdc "//<<chi2Add
799 <<"doca entr z tprop stat " << std::endl;
800
801 HitRefVec hl = tk->getVecHits();
802 HitRefVec::iterator it = hl.begin();
803 for (;it!=hl.end();++it){
804 RecMdcHit* h = *it;
805 int layer = MdcID::layer(h->getMdcId());
806 double _vprop = (layer<8) ? Constants::vpropInner : Constants::vpropOuter;
807 const MdcLayer* _layerPtr = m_gm->Layer(layer);
808 double _zlen = _layerPtr->zLength();
809 double z = h->getZhit();
810 double tprop;
811 if (0 == layer%2){
812 tprop = (0.5*_zlen + z)/_vprop; //odd
813 }else{
814 tprop = (0.5*_zlen - z)/_vprop; //even
815 }
816 // build the sense wires
817 //const MdcSWire* wire = m_gm->Wire(MdcID::layer(h->getMdcId()),MdcID::wire(h->getMdcId()));
818 //double z = wire->zForw();
819 //while(z<wire->zRear()){
820 //HepPoint3D pos;
821 //Hep3Vector dir;
822 //if(!(wire->getTraj()==NULL)){
823 //wire->getTraj()->getInfo(z,pos,dir);
824 //}
825 //std::cout<<"("<< wire->layer()->layNum()<<","<<wire->cell()<<" "<<wire->Id()<<")";
826 //std::cout<<" z, sag:"<<z
827 //<<", "<<wire->getTraj()->deltaY(z-wire->zForw())<<std::endl;
828 ////<<" pos:"<<pos <<" dir:"<<dir<<std::endl;
829 //z+=1.;
830 //}
831 std::cout<< setw(3) << h->getId() << setw(2) << h->getTrkId() <<
832 "(" << MdcID::layer(h->getMdcId()) <<"," << MdcID::wire(h->getMdcId()) <<")"<<
833 setw(10) << h->getFltLen() <<
834 setw(2) << h->getFlagLR() <<setw(10) << h->getDriftT() <<
835 setw(12) << h->getDriftDistLeft() <<//setw(8) << h->getErrDriftDistRight() <<
836 setw(8) << h->getTdc() <<//setw(10) << h->getChisqAdd() <<
837 setw(10) << h->getDoca() <<setw(10) << h->getEntra() <<
838 setw(10) << h->getZhit() << setw(10) << tprop<<
839 setw(2)<< h->getStat() << std::endl;
840 }
841}//print track
842
843void MdcxTrackFinder::bookNTuple(){
844 MsgStream log(msgSvc(), name());
845 StatusCode sc;
846 if(!m_hist) return;
847
848 g_poison = histoSvc()->book( "poison", "poison",43,0,42,288,0,288 );
849 g_csmax4 = histoSvc()->book( "csmax4", "csmax4",100,0,100 );
850 g_csmax3 = histoSvc()->book( "csmax3", "csmax3",50000,0,500000 );
851 g_omegag = histoSvc()->book( "omegag", "omegag",1000 ,0.,1.);
852 g_dPhiAU = histoSvc()->book( "dPhiAU", "dPhiAU",1000 ,0.,3.5);
853 g_dPhiAU_0 = histoSvc()->book( "dPhiAU_0", "dPhiAU_0",1000 ,0.,3.5);
854 g_dPhiAU_1 = histoSvc()->book( "dPhiAU_1", "dPhiAU_1",1000 ,0.,3.5);
855 g_dPhiAU_5 = histoSvc()->book( "dPhiAU_5", "dPhiAU_5",1000 ,0.,3.5);
856 g_dPhiAU_7 = histoSvc()->book( "dPhiAU_7", "dPhiAU_7",1000 ,0.,3.5);
857 g_dPhiAV = histoSvc()->book( "dPhiAV", "dPhiAV",1000 ,0.,3.5);
858 g_addSegPhi = histoSvc()->book( "addSegPhi", "addSegPhi",1000 ,0.,3.5);
859 g_dPhiAV_0 = histoSvc()->book( "dPhiAV_0", "dPhiAV_0",1000 ,0.,3.5);
860 g_dPhiAV_1 = histoSvc()->book( "dPhiAV_1", "dPhiAV_1",1000 ,0.,3.5);
861 g_dPhiAV_6 = histoSvc()->book( "dPhiAV_6", "dPhiAV_6",1000 ,0.,3.5);
862 g_dPhiAV_8 = histoSvc()->book( "dPhiAV_8", "dPhiAV_8",1000 ,0.,3.5);
863 g_trkllmk = histoSvc()->book( "trkllmk", "trkllmk",43,0.,42);
864 g_trklgood = histoSvc()->book( "trklgood", "trklgood",43,0.,42);
865 g_trklcircle = histoSvc()->book( "trklcircle", "trklcircle",43,0.,42);
866 g_trklhelix= histoSvc()->book( "trklhelix", "trklhelix",43,0.,42);
867 g_trkldrop1= histoSvc()->book( "trkldrop1", "trkldrop1",43,0.,42);
868 g_trkldrop2= histoSvc()->book( "trkldrop2", "trkldrop2",43,0.,42);
869 g_trklappend1= histoSvc()->book( "trklappend1", "trklappend1",43,0.,42);
870 g_trklappend2= histoSvc()->book( "trklappend2", "trklappend2",43,0.,42);
871 g_trklappend3= histoSvc()->book( "trklappend3", "trklappend3",43,0.,42);
872 //g_fitOmega= histoSvc()->book( "fitOmega", "fitOmega",200,0.,100.);
873 g_trklfirstProb= histoSvc()->book( "trklfirstProb", "trklfirstProb",200,0.,2.);
874 g_trkltemp= histoSvc()->book( "trkltemp", "trkltemp",200,-3.5,3.5);
875
876 g_trklproca= histoSvc()->book( "trklproca", "trklproca",100,0.,5.);
877 g_trkle= histoSvc()->book( "trkle", "trkle",200,0.,0.025);
878 g_trkld= histoSvc()->book( "trkld", "trkld",200,-1.2,1.2);
879 g_trklel= histoSvc()->book( "trklel", "trklel",200,0.,0.025,43,0,43);
880 g_trkldl= histoSvc()->book( "trkldl", "trkldl",200,-1.2,1.2,43,0,43);
881 g_trkldoca= histoSvc()->book( "trkldoca", "trkldoca",200,-1.2,1.2);
882 g_trkllayer= histoSvc()->book( "trkllayer", "trkllayer",43,0,43);
883 g_dropHitsSigma= histoSvc()->book( "dropHitsSigma", "dropHitsSigma",43,0,43,100,0,11);
884 g_addHitCut= histoSvc()->book( "addHitCut", "addHitCut",100,0,10);
885 g_addHitCut2d= histoSvc()->book( "addHitCut2d", "addHitCut2d",43,0,43,100,0,10);
886 //g_addSegPhiDiff = histoSvc()->book( "addSegPhiDiff", "addSegPhiDiff",100,-6.28,6.28);
887
888 NTuplePtr nt1(ntupleSvc(), "MdcxReco/rec");
889 if ( nt1 ) { m_xtuple1 = nt1;}
890 else {
891 m_xtuple1 = ntupleSvc()->book ("MdcxReco/rec", CLID_ColumnWiseTuple, "MdcxReco reconsturction results");
892 if ( m_xtuple1 ) {
893 sc = m_xtuple1->addItem ("t0", m_xt0);
894 sc = m_xtuple1->addItem ("timing", m_xtiming);
895 sc = m_xtuple1->addItem ("t0Stat", m_xt0Stat);
896 sc = m_xtuple1->addItem ("t0Truth", m_xt0Truth);
897
898 sc = m_xtuple1->addItem ("nSlay", m_xnSlay);
899 sc = m_xtuple1->addItem ("p", m_xp);
900 sc = m_xtuple1->addItem ("pt", m_xpt);
901 sc = m_xtuple1->addItem ("pz", m_xpz);
902 sc = m_xtuple1->addItem ("d0", m_xd0);
903 sc = m_xtuple1->addItem ("phi0", m_xphi0);
904 sc = m_xtuple1->addItem ("cpa", m_xcpa);
905 sc = m_xtuple1->addItem ("z0", m_xz0);
906 sc = m_xtuple1->addItem ("tanl", m_xtanl);
907 sc = m_xtuple1->addItem ("q", m_xq);
908 sc = m_xtuple1->addItem ("pocax", m_xpocax);
909 sc = m_xtuple1->addItem ("pocay", m_xpocay);
910 sc = m_xtuple1->addItem ("pocaz", m_xpocaz);
911
912 sc = m_xtuple1->addItem ("evtNo", m_xevtNo);
913 sc = m_xtuple1->addItem ("nSt", m_xnSt);
914 sc = m_xtuple1->addItem ("nAct", m_xnAct);
915 sc = m_xtuple1->addItem ("nDof", m_xnDof);
916 sc = m_xtuple1->addItem ("chi2", m_xchi2);
917 sc = m_xtuple1->addItem ("tkId", m_xtkId);
918 sc = m_xtuple1->addItem ("nHit", m_xnHit, 0, 7000);//# of hit/rec track
919 sc = m_xtuple1->addIndexedItem ("resid", m_xnHit, m_xresid);
920 sc = m_xtuple1->addIndexedItem ("sigma", m_xnHit, m_xsigma);
921 sc = m_xtuple1->addIndexedItem ("driftD", m_xnHit, m_xdriftD);
922 sc = m_xtuple1->addIndexedItem ("driftT", m_xnHit, m_xdriftT);
923 sc = m_xtuple1->addIndexedItem ("doca", m_xnHit, m_xdoca);
924 sc = m_xtuple1->addIndexedItem ("entra", m_xnHit, m_xentra);
925 sc = m_xtuple1->addIndexedItem ("ambig", m_xnHit, m_xambig);
926 sc = m_xtuple1->addIndexedItem ("fltLen", m_xnHit, m_xfltLen);
927 sc = m_xtuple1->addIndexedItem ("tof", m_xnHit, m_xtof);
928 sc = m_xtuple1->addIndexedItem ("act", m_xnHit, m_xact);
929 sc = m_xtuple1->addIndexedItem ("tdc", m_xnHit, m_xtdc);
930 sc = m_xtuple1->addIndexedItem ("adc", m_xnHit, m_xadc);
931 sc = m_xtuple1->addIndexedItem ("layer", m_xnHit, m_xlayer);
932 sc = m_xtuple1->addIndexedItem ("wire", m_xnHit, m_xwire);
933 sc = m_xtuple1->addIndexedItem ("x", m_xnHit, m_xx);
934 sc = m_xtuple1->addIndexedItem ("y", m_xnHit, m_xy);
935 sc = m_xtuple1->addIndexedItem ("z", m_xnHit, m_xz);
936 } else { // did not manage to book the N tuple....
937 log << MSG::ERROR << " Cannot book N-tuple: MdcxReco/rec" << endmsg;
938 //return;
939 }
940 }//end book of nt1
941 NTuplePtr nt4(ntupleSvc(), "MdcxReco/evt");
942 if ( nt4 ) { m_xtupleEvt = nt4;}
943 else {
944 m_xtupleEvt = ntupleSvc()->book ("MdcxReco/evt", CLID_ColumnWiseTuple, "MdcxReco event data");
945 if ( m_xtupleEvt ) {
946 sc = m_xtupleEvt->addItem ("evtNo", m_xt4EvtNo );
947 sc = m_xtupleEvt->addItem ("nRecTk", m_xt4nRecTk );
948 sc = m_xtupleEvt->addItem ("nTdsTk", m_xt4nTdsTk );
949 sc = m_xtupleEvt->addItem ("t0", m_xt4t0);
950 sc = m_xtupleEvt->addItem ("t0Stat", m_xt4t0Stat);
951 sc = m_xtupleEvt->addItem ("t0Truth", m_xt4t0Truth);
952 sc = m_xtupleEvt->addItem ("time", m_xt4time);
953 sc = m_xtupleEvt->addItem ("nDigi", m_xt4nDigi, 0, 7000);//# of hit/mc track
954 sc = m_xtupleEvt->addIndexedItem ("layer", m_xt4nDigi, m_xt4Layer);
955 sc = m_xtupleEvt->addIndexedItem ("rt", m_xt4nDigi, m_xt4Time);
956 sc = m_xtupleEvt->addIndexedItem ("rc", m_xt4nDigi, m_xt4Charge);
957 //sc = m_xtupleEvt->addIndexedItem ("rawHit", m_xt4nDigi, m_xt4rawHit);
958 //sc = m_xtupleEvt->addIndexedItem ("recHit", m_xt4nDigi, m_xt4recHit);
959 } else { // did not manage to book the N tuple....
960 log << MSG::ERROR << "Cannot book N-tuple: MdcxReco/evt" << endmsg;
961 //return;
962 }
963 }// end of book nt4
964
965 //book tuple of segment
966 NTuplePtr ntSeg(ntupleSvc(), "MdcxReco/seg");
967 if ( ntSeg ) { m_xtupleSeg = ntSeg;}
968 else {
969 m_xtupleSeg = ntupleSvc()->book ("MdcxReco/seg", CLID_ColumnWiseTuple, "MdcxTrackFinder segment data");
970 if ( m_xtupleSeg ) {
971 sc = m_xtupleSeg->addItem ("sl", m_xtsSl);
972 sc = m_xtupleSeg->addItem ("d0", m_xtsD0);
973 sc = m_xtupleSeg->addItem ("omega", m_xtsOmega);
974 sc = m_xtupleSeg->addItem ("phi0", m_xtsPhi0);
975 sc = m_xtupleSeg->addItem ("d0sl", m_xtsD0_sl_approx);
976 sc = m_xtupleSeg->addItem ("phi0sl", m_xtsPhi0_sl_approx);
977 sc = m_xtupleSeg->addItem ("xbbrrf", m_xtsXline_bbrrf);
978 sc = m_xtupleSeg->addItem ("ybbrrf", m_xtsYline_bbrrf);
979 sc = m_xtupleSeg->addItem ("xslope", m_xtsXline_slope);
980 sc = m_xtupleSeg->addItem ("yslope", m_xtsYline_slope);
981 sc = m_xtupleSeg->addItem ("pat", m_xtsPat);
982 sc = m_xtupleSeg->addItem ("chisq", m_xtsChisq);
983 sc = m_xtupleSeg->addItem ("nDigi", m_xtsNDigi, 0, 20);
984 sc = m_xtupleSeg->addIndexedItem ("layer", m_xtsNDigi, m_xtsLayer);
985 sc = m_xtupleSeg->addIndexedItem ("wire", m_xtsNDigi, m_xtsWire);
986 sc = m_xtupleSeg->addIndexedItem ("inSeg", m_xtsNDigi, m_xtsInSeg);
987 }else{
988 log << MSG::ERROR << "Cannot book N-tuple: MdcxReco/evt" << endmsg;
989 //return;
990 }
991 }
992 NTuplePtr nt5(ntupleSvc(), "MdcxReco/trkl");
993 if ( nt5 ) { m_xtupleTrkl= nt5;}
994 else {
995 m_xtupleTrkl= ntupleSvc()->book ("MdcxReco/trkl", CLID_RowWiseTuple, "MdcxReco track info");
996 if ( m_xtupleTrkl) {
997 sc = m_xtupleTrkl->addItem ("layer", m_xt5Layer);
998 sc = m_xtupleTrkl->addItem ("wire", m_xt5Wire);
999 }else{
1000 log << MSG::ERROR << "Cannot book N-tuple: MdcxReco/trkl" << endmsg;
1001 //return;
1002 }
1003 }
1004
1005 NTuplePtr ntCsmcSew(ntupleSvc(), "MdcxReco/csmc");
1006 if ( ntCsmcSew ) { m_xtupleCsmcSew = ntCsmcSew;}
1007 else {
1008 m_xtupleCsmcSew = ntupleSvc()->book ("MdcxReco/csmc", CLID_ColumnWiseTuple, "MdcxReco reconsturction results");
1009 if ( m_xtupleCsmcSew ) {
1010 sc = m_xtupleCsmcSew->addItem ("dD0", m_csmcD0);
1011 sc = m_xtupleCsmcSew->addItem ("dPhi0", m_csmcPhi0);
1012 sc = m_xtupleCsmcSew->addItem ("dZ0", m_csmcZ0);
1013 sc = m_xtupleCsmcSew->addItem ("dOmega", m_csmcOmega);
1014 sc = m_xtupleCsmcSew->addItem ("dPt", m_csmcPt);
1015 sc = m_xtupleCsmcSew->addItem ("dTanl", m_csmcTanl);
1016 }else{
1017 log << MSG::ERROR << "Cannot book N-tuple: MdcxReco/csmc" << endmsg;
1018 //return;
1019 }
1020 }
1021 NTuplePtr ntSegAdd1(ntupleSvc(), "MdcxReco/addSeg1");
1022 if ( ntSegAdd1 ) { m_xtupleAddSeg1 = ntSegAdd1;}
1023 else {
1024 m_xtupleAddSeg1 = ntupleSvc()->book ("MdcxReco/addSeg1", CLID_ColumnWiseTuple, "MdcxReco event data");
1025 if ( m_xtupleAddSeg1 ) {
1026 sc = m_xtupleAddSeg1->addItem ("same", m_addSegSame);
1027 sc = m_xtupleAddSeg1->addItem ("seedSl", m_addSegSeedSl );
1028 sc = m_xtupleAddSeg1->addItem ("seedPhi", m_addSegSeedPhi );
1029 sc = m_xtupleAddSeg1->addItem ("seedPhiLay",m_addSegSeedPhiLay );
1030 sc = m_xtupleAddSeg1->addItem ("seedLen", m_addSegSeedLen );
1031 sc = m_xtupleAddSeg1->addItem ("seedD0", m_addSegSeedD0 );
1032 sc = m_xtupleAddSeg1->addItem ("seedPhi0", m_addSegSeedPhi0 );
1033 sc = m_xtupleAddSeg1->addItem ("addSl", m_addSegAddSl );
1034 sc = m_xtupleAddSeg1->addItem ("addPhi", m_addSegAddPhi );
1035 sc = m_xtupleAddSeg1->addItem ("addPhiLay", m_addSegAddPhiLay );
1036 sc = m_xtupleAddSeg1->addItem ("addLen", m_addSegAddLen );
1037 sc = m_xtupleAddSeg1->addItem ("addD0", m_addSegAddD0 );
1038 sc = m_xtupleAddSeg1->addItem ("addPhi0", m_addSegAddPhi0 );
1039 }
1040 else{
1041 log << MSG::ERROR << "Cannot book N-tuple: MdcxReco/addSeg1" << endmsg;
1042 //return;
1043 }
1044 }
1045
1046 NTuplePtr ntSegComb(ntupleSvc(), "MdcxReco/segComb");
1047 if ( ntSegComb ) { m_xtupleSegComb = ntSegComb;}
1048 else {
1049 m_xtupleSegComb = ntupleSvc()->book ("MdcxReco/segComb", CLID_ColumnWiseTuple, "MdcxReco event data");
1050 if ( m_xtupleSegComb ) {
1051 sc = m_xtupleSegComb->addItem ("evtNo", m_segCombEvtNo);
1052 sc = m_xtupleSegComb->addItem ("omega", m_segCombOmega);
1053 sc = m_xtupleSegComb->addItem ("sameAU", m_segCombSameAU);
1054 sc = m_xtupleSegComb->addItem ("sameUV", m_segCombSameUV);
1055 sc = m_xtupleSegComb->addItem ("dlenAU", m_segCombDLenAU);
1056 sc = m_xtupleSegComb->addItem ("dlenUV", m_segCombDLenUV);
1057 sc = m_xtupleSegComb->addItem ("slA", m_segCombSlA);
1058 sc = m_xtupleSegComb->addItem ("slU", m_segCombSlU);
1059 sc = m_xtupleSegComb->addItem ("slV", m_segCombSlV);
1060 sc = m_xtupleSegComb->addItem ("phiA", m_segCombPhiA);
1061 sc = m_xtupleSegComb->addItem ("phiU", m_segCombPhiU);
1062 sc = m_xtupleSegComb->addItem ("phiV", m_segCombPhiV);
1063 }
1064 else{
1065 log << MSG::ERROR << "Cannot book N-tuple: MdcxReco/segComb" << endmsg;
1066 //return;
1067 }
1068 }
1069
1070 NTuplePtr ntDropHits(ntupleSvc(), "MdcxReco/dropHits");
1071 if ( ntDropHits ) { m_xtupleDropHits = ntDropHits;}
1072 else {
1073 m_xtupleDropHits = ntupleSvc()->book ("MdcxReco/dropHits", CLID_ColumnWiseTuple, "MdcxReco event data");
1074 if ( m_xtupleDropHits ) {
1075 sc = m_xtupleDropHits->addItem ("evtNo", m_segDropHitsEvtNo);
1076 sc = m_xtupleDropHits->addItem ("layer", m_segDropHitsLayer);
1077 sc = m_xtupleDropHits->addItem ("wire", m_segDropHitsWire);
1078 sc = m_xtupleDropHits->addItem ("pull", m_segDropHitsPull);
1079 sc = m_xtupleDropHits->addItem ("doca", m_segDropHitsDoca);
1080 sc = m_xtupleDropHits->addItem ("sigma", m_segDropHitsSigma);
1081 sc = m_xtupleDropHits->addItem ("drift", m_segDropHitsDrift);
1082 sc = m_xtupleDropHits->addItem ("mcTkId", m_segDropHitsMcTkId);
1083 } else{
1084 log << MSG::ERROR << "Cannot book N-tuple: MdcxReco/dropHits" << endmsg;
1085 //return;
1086 }
1087 }
1088 NTuplePtr ntAddSeg2(ntupleSvc(), "MdcxReco/addSeg2");
1089 if ( ntAddSeg2 ) { m_xtupleAddSeg2 = ntAddSeg2;}
1090 else {
1091 m_xtupleAddSeg2 = ntupleSvc()->book ("MdcxReco/addSeg2", CLID_ColumnWiseTuple, "MdcxReco event data");
1092 if ( m_xtupleAddSeg2 ) {
1093 sc = m_xtupleAddSeg2->addItem ("evtNo", m_addSegEvtNo);
1094 sc = m_xtupleAddSeg2->addItem ("slayer", m_addSegSlayer);
1095 sc = m_xtupleAddSeg2->addItem ("poca", m_addSegPoca);
1096 sc = m_xtupleAddSeg2->addItem ("len", m_addSegLen);
1097 } else{
1098 log << MSG::ERROR << "Cannot book N-tuple: MdcxReco/addSeg2" << endmsg;
1099 //return;
1100 }
1101 }
1102 //--------------end of book ntuple------------------
1103}
1104
1105
1106void MdcxTrackFinder::fillEvent(){
1107 //-----------get raw digi-----------------------
1108 if (t_nDigi<=0) return;
1109 if (m_xtupleEvt == NULL ) return;
1110 m_xt4EvtNo = m_eventNo;
1111 m_xt4t0 = m_bunchT0;
1112 m_xt4t0Stat = m_t0Stat;
1113 m_xt4t0Truth = m_t0Truth;
1114 m_xt4nRecTk = nTk;
1115 m_xt4nTdsTk = t_nTdsTk;
1116 m_xt4nDigi = t_nDigi;
1117 int iDigi=0;
1118 MdcDigiCol::iterator iter = mdcDigiVec.begin();
1119 for (;iDigi<t_nDigi; iter++ ) {
1120 int l = MdcID::layer((*iter)->identify());
1121 m_xt4Layer[iDigi] = l;
1122 //int w = MdcID::wire((*iter)->identify());
1123 m_xt4Time[iDigi] = RawDataUtil::MdcTime((*iter)->getTimeChannel());
1124 m_xt4Charge[iDigi] = RawDataUtil::MdcCharge((*iter)->getChargeChannel());
1125 //m_xt4rawHit[l]++;
1126 iDigi++;
1127 }//end for iter
1128 m_xtupleEvt->write();
1129}
1130
1131void MdcxTrackFinder::fillTrack(TrkRecoTrk* atrack){
1132
1133 if(!m_xtuple1) return;
1134 //-----------get MC truth data-------------------
1135 m_xevtNo = m_eventNo;
1136 int recHitMap[43][288]={0};
1137 //int nTk = (*m_tracks).nTrack();
1138 //for (int itrack = 0; itrack < nTk; itrack++) {
1139 if (atrack==NULL) return;
1140
1141 const TrkFit* fitResult = atrack->fitResult();
1142 if (fitResult == 0) return;//check the fit worked
1143
1144 //for fill ntuples
1145 int nSt=0; //int nSeg=0;
1146 int seg[11] = {0}; int segme;
1147 //-----------hit list-------------
1148 HitRefVec hitRefVec;
1149 const TrkHitList* hitlist = atrack->hits();
1150
1151 TrkHitList::hot_iterator hot = hitlist->begin();
1152 int layerCount[43]={0};
1153 int i=0;
1154 for (;hot!= hitlist->end();hot++){
1155
1156 const MdcRecoHitOnTrack* recoHot;
1157 recoHot = dynamic_cast<const MdcRecoHitOnTrack*> (&(*hot));
1158 int layer = recoHot->mdcHit()->layernumber();
1159 int wire = recoHot->mdcHit()->wirenumber();
1160 m_xlayer[i] = layer;
1161 //m_xt4recHit[layer]++;//fill rec hit for hit eff
1162 m_xwire[i] = wire;
1163 layerCount[layer]++;
1164 recHitMap[layer][wire]++;
1165 m_xambig[i] = recoHot->wireAmbig();// wire ambig
1166 //fltLen
1167 double fltLen = recoHot->fltLen();
1168 m_xfltLen[i] = fltLen;
1169 double tof = recoHot->getParentRep()->arrivalTime(fltLen);
1170 //position
1171 HepPoint3D pos; Hep3Vector dir;
1172 recoHot->getParentRep()->traj().getInfo(fltLen,pos,dir);
1173 m_xx[i] = pos.x();
1174 m_xy[i] = pos.y();
1175 m_xz[i] = pos.z();
1176 m_xdriftT[i] = recoHot->mdcHit()->driftTime(tof,pos.z());
1177 m_xtof[i] = tof;
1178 m_xdriftD[i] = recoHot->drift();
1179 m_xsigma[i] = recoHot->hitRms();
1180 m_xtdc[i] = recoHot->rawTime();
1181 m_xadc[i] = recoHot->mdcHit()->charge();
1182 m_xdoca[i] = recoHot->dcaToWire();//sign w.r.t. dirft() FIXME
1183 m_xentra[i] = recoHot->entranceAngle();
1184 //m_xentraHit[i] = recoHot->entranceAngleHit();
1185 m_xact[i] = recoHot->isActive();
1186 //resid
1187 double res=999.,rese=999.;
1188 if (recoHot->resid(res,rese,false)){
1189 }else{}
1190 m_xresid[i] = res;
1191 //for n seg
1192 segme=0;
1193 if ( layer >0 ) {segme=(layer-1)/4;}
1194 seg[segme]++;
1195 if (recoHot->layer()->view()) { ++nSt; }
1196 i++;
1197 }// end fill hit
1198
1199 int nSlay=0;
1200 for(int i=0;i<11;i++){
1201 if (seg[i]>0) nSlay++;
1202 }
1203 m_xnSlay = nSlay;
1204
1205 //------------fill track result-------------
1206 //m_xtkId = itrack;
1207 //track parameters at closest approach to beamline
1208 double fltLenPoca = 0.0;
1209 TrkExchangePar helix = fitResult->helix(fltLenPoca);
1210 m_xphi0 = helix.phi0();
1211 m_xtanl = helix.tanDip();
1212 m_xz0 = helix.z0();
1213 m_xd0 = helix.d0();
1214 if(fabs(m_xz0)>20||fabs(m_xd0)>2) {
1215 //b_saveEvent = true;
1216 if(m_debug>1) std::cout<<"evt:"<<m_eventNo<<" BigVtx:"
1217 <<" d0 "<<helix.d0()<<" z0 "<<helix.z0()<<std::endl;
1218 }
1219 m_xpt = fitResult->pt();
1220 if (m_xpt > 0.00001){
1221 m_xcpa = fitResult->charge()/fitResult->pt();
1222 }
1223 //momenta and position
1224 CLHEP::Hep3Vector p1 = fitResult->momentum(fltLenPoca);
1225 double px,py,pz,pxy;
1226 pxy = fitResult->pt();
1227 px = p1.x();
1228 py = p1.y();
1229 pz = p1.z();
1230 m_xp = p1.mag();
1231 m_xpz = pz;
1232 HepPoint3D poca = fitResult->position(fltLenPoca);
1233 m_xpocax = poca.x();
1234 m_xpocay = poca.y();
1235 m_xpocaz = poca.z();
1236
1237 m_xq = fitResult->charge();
1238 m_xnAct = fitResult->nActive();
1239 m_xnHit = hitlist->nHit();
1240 m_xnDof = fitResult->nDof();
1241 m_xnSt = nSt;
1242 m_xchi2 = fitResult->chisq();
1243 //for (int l=0;l<43;l++) m_xlayerCount[l] = layerCount[l];
1244 m_xt0 = m_bunchT0;
1245 m_xtiming = m_timing;
1246 m_xt0Stat = m_t0Stat;
1247 m_xt0Truth= m_t0Truth;
1248 m_xtuple1->write();
1249 //}//end of loop rec trk list
1250
1251}//end of
1252
1253void MdcxTrackFinder::dumpMdcxSegs(const HepAList<MdcxSeg>& segList) const {
1254 cout << name()<<" found " <<segList.length() << " segs :" << endl;
1255 for (int i =0; i< segList.length(); i++){
1256 std::cout<<i<<" ";
1257 segList[i]->printSeg();
1258 }
1259}
1260
1261void MdcxTrackFinder::fillMdcxSegs(const HepAList<MdcxSeg>& segList) const {
1262 if(!m_xtupleSeg) return;
1263
1264 int cellMax[43] ={
1265 40,44,48,56, 64,72,80,80, 76,76,88,88,
1266 100,100,112,112, 128,128,140,140, 160,160,160,160,
1267 176,176,176,176, 208,208,208,208, 240,240,240,240,
1268 256,256,256,256, 288,288,288 };
1269 // Fill hits of every layer after segment finding
1270 int hitInSegList[43][288];
1271 for (int ii=0;ii<43;ii++){
1272 for (int jj=0;jj<cellMax[ii];jj++){ hitInSegList[ii][jj] = 0; }
1273 }
1274 for (int i = 0; i < segList.length(); i++) {
1275 MdcxSeg* aSeg = segList[i];
1276 for (int ihit=0 ; ihit< aSeg->Nhits() ; ihit++){
1277 const MdcxHit* hit = aSeg->XHitList()[ihit];
1278 int layer = hit->Layer();
1279 int wire = hit->WireNo();
1280 hitInSegList[layer][wire]++;
1281 }
1282 }
1283 for (int i = 0; i < segList.length(); i++) {
1284 MdcxSeg* aSeg = segList[i];
1285 if(aSeg==NULL)continue;
1286 m_xtsSl = aSeg->SuperLayer();
1287 m_xtsD0 = aSeg->D0();
1288 m_xtsOmega = aSeg->Omega();
1289 m_xtsPhi0 = aSeg->Phi0();
1292 m_xtsXline_bbrrf = aSeg->Xline_bbrrf();
1293 m_xtsYline_bbrrf = aSeg->Yline_bbrrf();
1294 m_xtsXline_slope = aSeg->Xline_slope();
1295 m_xtsYline_slope = aSeg->Yline_slope();
1296 int patIndex = -1;
1297 for (int ii = 0;ii<14;ii++){
1298 if (aSeg->Pat() == (int) MdcxSegPatterns::patt4[ii]){
1299 patIndex = ii;
1300 break;
1301 }
1302 }
1303 for (int ii = 0;ii<20;ii++){
1304 if (aSeg->Pat() == (int) MdcxSegPatterns::patt3[ii]){
1305 patIndex = ii +15;
1306 break;
1307 }
1308 }
1309 m_xtsPat = patIndex;
1310 m_xtsChisq = aSeg->Chisq();
1311 m_xtsNDigi = aSeg->Nhits();
1312 for (int ihit=0 ; ihit< aSeg->Nhits() ; ihit++){
1313 const MdcxHit* hit = aSeg->XHitList()[ihit];
1314 int layer = hit->Layer();
1315 int wire = hit->WireNo();
1316 m_xtsLayer[ihit] = layer;
1317 m_xtsWire[ihit] = wire;
1318 m_xtsInSeg[ihit] = hitInSegList[layer][wire];
1319 //std::cout<< "hitInSegList "<<hitInSegList[layer][wire] << std::endl;//yzhang debug
1320 }
1321 m_xtupleSeg->write();
1322 }
1323
1324}//end of fillMdcxSegs
1325
1326void MdcxTrackFinder::dumpTrackList(const HepAList<MdcxFittedHel>& trackList) const {
1327 std::cout<< "dump track list " <<std::endl;
1328 for (int i=0;i<trackList.length();i++){
1329 std::cout<< "track "<<i<<std::endl;
1330 for (int ihit=0 ; ihit< trackList[i]->Nhits() ; ihit++){
1331 const MdcxHit* hit = trackList[i]->XHitList()[ihit];
1332 int layer = hit->Layer();
1333 int wire = hit->WireNo();
1334 std::cout<< " ("<<layer << ","<< wire<<") ";
1335 }
1336 std::cout<<std::endl;
1337 }
1338}// end of dumpTrackList
1339
1340void MdcxTrackFinder::fillTrkl(const HepAList<MdcxFittedHel>& firsttrkl) const {
1341 if(!m_xtupleTrkl) return;
1342 int nDigi = 0;
1343 int iDigi = 0;
1344 for (int i =0; i< firsttrkl.length(); i++){
1345 nDigi+= firsttrkl[i]->Nhits();
1346 }
1347 for (int i=0;i<firsttrkl.length();i++){
1348 for (int ihit=0 ; ihit< firsttrkl[i]->Nhits() ; ihit++){
1349 const MdcxHit* hit = firsttrkl[i]->XHitList()[ihit];
1350 int layer = hit->Layer();
1351 int wire = hit->WireNo();
1352 m_xt5Layer = layer;
1353 m_xt5Wire = wire;
1354 m_xtupleTrkl->write();
1355 iDigi++;
1356 }
1357 }
1358}//end of fillTrkl
1359
1360void MdcxTrackFinder::fillMcTruth(){
1361 MsgStream log(msgSvc(), name());
1362 StatusCode sc;
1363 //Initialize
1364 for (int ii=0;ii<43;ii++){
1365 for (int jj=0;jj<288;jj++){
1366 haveDigi[ii][jj] = -2;
1367 }
1368 }
1369 int nDigi = mdcDigiVec.size();
1370 for (int iDigi =0 ;iDigi<nDigi; iDigi++ ) {
1371 int l = MdcID::layer((mdcDigiVec[iDigi])->identify());
1372 int w = MdcID::wire((mdcDigiVec[iDigi])->identify());
1373 //haveDigi[l][w]=(mdcDigiVec[iDigi])->getTrackIndex();
1374 haveDigi[l][w]=1;
1375 }
1376
1377 if( m_mcHist ){
1378 //------------------get event start time truth-----------
1379 m_t0Truth = -10.;
1380 SmartDataPtr<Event::McParticleCol> mcParticleCol(eventSvc(),"/Event/MC/McParticleCol");
1381 if(!mcParticleCol){
1382 log << MSG::INFO << "Could not retrieve McParticelCol" << endreq;
1383 }else {
1384 Event::McParticleCol::iterator iter_mc = mcParticleCol->begin();
1385 for (;iter_mc != mcParticleCol->end(); iter_mc++){
1386 if ((*iter_mc)->primaryParticle()){
1387 m_t0Truth = (*iter_mc)->initialPosition().t();
1388 //px = (*iter_mc)->initialFourMomentum().x()/1000.;//GeV
1389 //py = (*iter_mc)->initialFourMomentum().y()/1000.;//GeV
1390 //pz = (*iter_mc)->initialFourMomentum().z()/1000.;//GeV
1391 }
1392 }
1393 }
1394 }
1395 //------------------Retrieve MC truth MdcMcHit------------
1396 /*
1397 SmartDataPtr<Event::MdcMcHitCol> mcMdcMcHitCol(eventSvc(),"/Event/MC/MdcMcHitCol");
1398 if (!mcMdcMcHitCol) {
1399 log << MSG::WARNING << "Could not find MdcMcHit" << endreq;
1400 }else{
1401 Event::MdcMcHitCol::iterator iter_mchit = mcMdcMcHitCol->begin();
1402 for (;iter_mchit != mcMdcMcHitCol->end(); iter_mchit++ ) {
1403 const Identifier id= (*iter_mchit)->identify();
1404 int layer = MdcID::layer(id);
1405 int wire = MdcID::wire(id);
1406 mcDrift[layer][wire] = (*iter_mchit)->getDriftDistance(); //drift in MC.
1407 mcLR[layer][wire] = (*iter_mchit)->getPositionFlag();
1408 mcX[layer][wire] = (*iter_mchit)->getPositionX();
1409 mcY[layer][wire] = (*iter_mchit)->getPositionY();
1410 mcZ[layer][wire] = (*iter_mchit)->getPositionZ();
1411 if (mcLR[layer][wire] == 0) mcLR[layer][wire] = -1;
1412 }
1413 }
1414 */
1415}
1416
1417void MdcxTrackFinder::dropMultiHotInLayer(TrkHitList* trkHitList){
1418 double tdr[43];
1419 double tdr_wire[43];
1420 for(int i=0; i<43; i++){tdr[i]=9999.;}
1421
1422 // make flag
1423 TrkHotList::hot_iterator hotIter = trkHitList->hotList().begin();
1424 while (hotIter!=trkHitList->hotList().end()) {
1425 MdcHitOnTrack* hot = const_cast<MdcHitOnTrack*> (&(*hotIter->mdcHitOnTrack()));
1426
1427 //driftTime(tof,z)
1428 double dt = hot->mdcHit()->driftTime(0.,0.);
1429 int layer = hot->mdcHit()->layernumber();
1430 int wire = hot->mdcHit()->wirenumber();
1431
1432 //std::cout<<__FILE__<<" "<<dt<< std::endl;
1433 if (dt < tdr[layer]) {
1434 tdr[layer] = dt;
1435 tdr_wire[layer] = wire;
1436 }
1437 hotIter++;
1438 }
1439
1440 //std::cout<<" tdr wire ";
1441 //for(int i=0;i<43;i++){
1442 //std::cout<<i<<","<<tdr_wire[i]<<" tdr="<<tdr[i]<<"\n";
1443 //}
1444 // inactive multi hit
1445 hotIter = trkHitList->hotList().begin();
1446 while (hotIter!=trkHitList->hotList().end()) {
1447 int layer = hotIter->mdcHitOnTrack()->mdcHit()->layernumber();
1448 int wire = hotIter->mdcHitOnTrack()->mdcHit()->wirenumber();
1449 //double dt = hotIter->mdcHitOnTrack()->mdcHit()->driftTime(0.,0.)-m_bunchT0;
1450
1451 if ((tdr[layer] <9998.) && (tdr_wire[layer]!=wire)){
1452 MdcHitOnTrack* hot = const_cast<MdcHitOnTrack*> (&(*hotIter->mdcHitOnTrack()));
1453 hot->setActivity(false);
1454 trkHitList->removeHit( hotIter->mdcHitOnTrack()->mdcHit() );
1455 //std::cout<<__FILE__<<" inactive "<< layer<<" "<<wire<<" dt "<<dt << std::endl;
1456 }else{
1457 hotIter++;
1458 }
1459 }
1460}
1461void MdcxTrackFinder::dumpTdsTrack(RecMdcTrackCol* trackList){
1462 std::cout<< "tksize = "<<trackList->size() << std::endl;//yzhang debug
1463 RecMdcTrackCol::iterator it = trackList->begin();
1464 for (;it!= trackList->end();it++){
1465 RecMdcTrack *tk = *it;
1466 std::cout<< "//====RecMdcTrack "<<tk->trackId()<<"====:" << std::endl;
1467 cout <<" d0 "<<tk->helix(0)
1468 <<" phi0 "<<tk->helix(1)
1469 <<" cpa "<<tk->helix(2)
1470 <<" z0 "<<tk->helix(3)
1471 <<" tanl "<<tk->helix(4)
1472 <<endl;
1473 std::cout<<" q "<<tk->charge()
1474 <<" theta "<<tk->theta()
1475 <<" phi "<<tk->phi()
1476 <<" x0 "<<tk->x()
1477 <<" y0 "<<tk->y()
1478 <<" z0 "<<tk->z()
1479 <<" r0 "<<tk->r()
1480 <<endl;
1481 std::cout <<" p "<<tk->p()
1482 <<" pt "<<tk->pxy()
1483 <<" px "<<tk->px()
1484 <<" py "<<tk->py()
1485 <<" pz "<<tk->pz()
1486 <<endl;
1487 std::cout<<" tkStat "<<tk->stat()
1488 <<" chi2 "<<tk->chi2()
1489 <<" ndof "<<tk->ndof()
1490 <<" nhit "<<tk->getNhits()
1491 <<" nst "<<tk->nster()
1492 <<endl;
1493
1494 int nhits = tk->getVecHits().size();
1495 std::cout<<nhits <<" Hits: " << std::endl;
1496 for(int ii=0; ii <nhits ; ii++){
1497 Identifier id(tk->getVecHits()[ii]->getMdcId());
1498 int layer = MdcID::layer(id);
1499 int wire = MdcID::wire(id);
1500 cout<<"("<< layer <<","<<wire<<","<<tk->getVecHits()[ii]->getStat()
1501 <<",lr:"<<tk->getVecHits()[ii]->getFlagLR()<<") ";
1502 }//end of hit list
1503 std::cout << " "<< std::endl;
1504 }//end of tk list
1505 std::cout << " "<< std::endl;
1506
1507
1508}//end of dumpTdsTrack
1509
1510void MdcxTrackFinder::dumpTdsHits(RecMdcHitCol* hitList){
1511 std::cout<<__FILE__<<" "<<__LINE__<<" All hits in TDS, nhit="<<hitList->size()<< std::endl;
1512 RecMdcHitCol::iterator it = hitList->begin();
1513 for(;it!= hitList->end(); it++){
1514 RecMdcHit* h = (*it);
1515 Identifier id(h->getMdcId());
1516 int layer = MdcID::layer(id);
1517 int wire = MdcID::wire(id);
1518 cout<<"("<< layer <<","<<wire<<") lr:"<<h->getFlagLR()<<" stat:"<<h->getStat()<<" tk:"<<h->getTrkId() <<" doca:"<<setw(10)<<h->getDoca()<<std::endl;
1519 }//end of hit list
1520}
double w
EvtStreamInputIterator< typename Generator::result_type > iter(Generator gen, int N=0)
ObjectVector< MdcHit > MdcHitCol
ObjectVector< RecMdcHit > RecMdcHitCol
ObjectVector< RecMdcTrack > RecMdcTrackCol
SmartRefVector< RecMdcHit > HitRefVec
NTuple::Array< double > m_xfltLen
NTuple::Item< double > m_addSegAddPhiLay
NTuple::Array< double > m_xwire
NTuple::Item< double > m_xtsXline_bbrrf
NTuple::Item< double > m_xt0Stat
NTuple::Array< double > m_xdriftT
AIDA::IHistogram1D * g_trkldrop1
NTuple::Array< double > m_xdriftD
NTuple::Item< double > m_xpocax
NTuple::Item< double > m_xtsChisq
NTuple::Item< double > m_xtsYline_slope
NTuple::Item< double > m_xtsYline_bbrrf
AIDA::IHistogram1D * g_trklappend2
NTuple::Item< double > m_addSegSeedPhi0
NTuple::Array< double > m_xambig
NTuple::Item< double > m_xtiming
NTuple::Item< double > m_segCombDLenUV
AIDA::IHistogram1D * g_trklappend1
NTuple::Item< long > m_addSegSlayer
NTuple::Item< double > m_xchi2
NTuple::Array< double > m_xresid
NTuple::Item< double > m_xpocaz
AIDA::IHistogram1D * g_trklproca
NTuple::Item< long > m_segDropHitsWire
NTuple::Item< double > m_addSegSeedSl
NTuple::Item< double > m_xnAct
NTuple::Array< double > m_xtof
NTuple::Item< double > m_addSegAddLen
AIDA::IHistogram1D * g_addSegPhi
NTuple::Item< double > m_xt4nTdsTk
NTuple::Item< double > m_addSegLen
NTuple::Item< double > m_xt0
NTuple::Item< double > m_segDropHitsDoca
NTuple::Item< double > m_xtsPhi0_sl_approx
NTuple::Item< double > m_csmcTanl
NTuple::Item< double > m_segDropHitsSigma
NTuple::Item< double > m_xtsD0
NTuple::Item< double > m_addSegSeedPhiLay
NTuple::Array< double > m_xact
NTuple::Item< double > m_csmcOmega
NTuple::Item< double > m_xpz
NTuple::Item< double > m_addSegSeedLen
NTuple::Item< double > m_xnSt
NTuple::Array< long > m_xtsInSeg
NTuple::Item< double > m_segCombSlV
NTuple::Item< long > m_segDropHitsLayer
NTuple::Item< double > m_addSegSeedD0
AIDA::IHistogram2D * g_dropHitsSigma
NTuple::Item< double > m_xtsPhi0
AIDA::IHistogram1D * g_trklappend3
NTuple::Item< double > m_xt4t0
AIDA::IHistogram1D * g_trklhelix
NTuple::Item< double > m_xt4t0Truth
NTuple::Item< double > m_segDropHitsMcTkId
NTuple::Array< double > m_xdoca
NTuple::Item< double > m_xevtNo
AIDA::IHistogram1D * g_addHitCut
NTuple::Item< double > m_xt4nRecTk
NTuple::Item< double > m_segCombPhiA
NTuple::Item< double > m_segCombSlA
NTuple::Array< double > m_xsigma
NTuple::Item< double > m_xt4time
NTuple::Item< double > m_segCombDLenAU
AIDA::IHistogram2D * g_addHitCut2d
AIDA::IHistogram1D * g_trkllayer
NTuple::Item< double > m_xphi0
NTuple::Item< double > m_segCombSlU
NTuple::Item< long > m_segCombEvtNo
NTuple::Item< double > m_segCombPhiU
NTuple::Item< double > m_csmcPhi0
NTuple::Array< double > m_xadc
NTuple::Array< double > m_xt4Time
NTuple::Item< double > m_segCombSameUV
NTuple::Item< double > m_addSegAddPhi
NTuple::Item< double > m_addSegSeedPhi
NTuple::Array< double > m_xx
NTuple::Item< double > m_xz0
NTuple::Item< double > m_csmcPt
AIDA::IHistogram1D * g_trkldrop2
NTuple::Array< double > m_xtdc
NTuple::Item< double > m_segCombPhiV
NTuple::Item< double > m_xtanl
NTuple::Item< double > m_addSegAddPhi0
NTuple::Item< double > m_xt0Truth
NTuple::Item< double > m_addSegAddSl
NTuple::Item< double > m_segCombOmega
NTuple::Array< long > m_xt4Layer
NTuple::Item< double > m_xtsD0_sl_approx
NTuple::Item< double > m_xtkId
NTuple::Item< double > m_addSegAddD0
NTuple::Item< long > m_segDropHitsEvtNo
NTuple::Item< double > m_segCombSameAU
NTuple::Item< double > m_xpt
NTuple::Item< double > m_xpocay
NTuple::Item< double > m_segDropHitsPull
NTuple::Item< double > m_xcpa
NTuple::Item< double > m_csmcD0
AIDA::IHistogram1D * g_trklfirstProb
NTuple::Array< long > m_xtsLayer
NTuple::Array< double > m_xentra
NTuple::Item< double > m_xnDof
NTuple::Array< double > m_xlayer
NTuple::Array< double > m_xy
NTuple::Array< double > m_xz
NTuple::Item< double > m_xtsOmega
NTuple::Item< double > m_xd0
NTuple::Item< long > m_addSegEvtNo
NTuple::Item< double > m_segDropHitsDrift
NTuple::Item< double > m_csmcZ0
NTuple::Item< double > m_addSegPoca
AIDA::IHistogram1D * g_trklcircle
NTuple::Item< double > m_xtsXline_slope
NTuple::Array< double > m_xt4Charge
TGraph2DErrors * dt
Definition: McCor.cxx:45
const HepSymMatrix err() const
const HepVector helix() const
......
const MdcLayer * Layer(unsigned id) const
static MdcDetector * instance()
Definition: MdcDetector.cxx:21
virtual const MdcHit * mdcHit() const
int wireAmbig() const
double rawTime() const
double dcaToWire() const
double entranceAngle() const
const MdcLayer * layer() const
void setCalibSvc(const MdcCalibFunSvc *calibSvc)
Definition: MdcHit.cxx:136
double driftTime(double tof, double z) const
Definition: MdcHit.cxx:142
void setCountPropTime(const bool count)
static int layer(const Identifier &id)
Values of different levels (failure returns 0)
Definition: MdcID.cxx:49
static int wire(const Identifier &id)
Definition: MdcID.cxx:54
const MdcHit * mdcHit() const
const HepAList< MdcxSeg > & GetMdcxSeglist()
const HepAList< MdcxFittedHel > & GetMdcxTrklist()
const HepAList< MdcxHit > & XHitList() const
int SuperLayer(int hitno=0) const
int Fail(float Probmin=0.0) const
MdcxFittedHel & Grow(const MdcxFittedHel &, const HepAList< MdcxHit > &)
static void setMdcDetector(const MdcDetector *gm)
Definition: MdcxHit.cxx:68
static void setMdcCalibFunSvc(const MdcCalibFunSvc *calibSvc)
Definition: MdcxHit.cxx:60
static void setCountPropTime(bool countPropTime)
Definition: MdcxHit.cxx:64
const HepAList< MdcxHit > & GetMdcxHitList()
void reset()
Definition: MdcxHits.cxx:55
void create(MdcDigiVec digiVec, float c0=0.0, float cresol=0.0180)
Definition: MdcxHits.cxx:59
const HepAList< MdcxFittedHel > & GetMergedTrklist()
virtual ~MdcxTrackFinder()
MdcxTrackFinder(const std::string &name, ISvcLocator *pSvcLocator)
StatusCode finalize()
StatusCode beginRun()
StatusCode execute()
StatusCode initialize()
static void readMCppTable(std::string filenm)
static double MdcCharge(int chargeChannel)
virtual void getInfo(double fltLen, HepPoint3D &pos, Hep3Vector &direction) const =0
virtual Hep3Vector momentum(double fltL=0.) const =0
virtual double pt(double fltL=0.) const =0
virtual double chisq() const =0
virtual int charge() const =0
virtual int nDof() const =0
virtual const TrkDifTraj & traj() const =0
virtual HepPoint3D position(double fltL) const =0
void print(std::ostream &ostr) const
virtual void addHistory(const TrkErrCode &status, const char *modulename)
virtual int nActive() const =0
virtual TrkExchangePar helix(double fltL) const =0
TrkErrCode fit()
Definition: TrkHitList.cxx:59
TrkHitOnTrk * appendHot(const TrkHitOnTrk *theHot)
Definition: TrkHitList.cxx:91
const TrkHotList & hotList() const
Definition: TrkHitList.cxx:51
bool removeHit(const TrkFundHit *theHit)
Definition: TrkHitList.cxx:69
void setActivity(bool turnOn)
Definition: TrkHitOnTrk.cxx:96
void setUsability(int usability)
double resid(bool exclude=false) const
const TrkFit * fitResult() const
Definition: TrkRecoTrk.cxx:387
virtual void printAll(std::ostream &) const
const TrkFitStatus * status() const
Definition: TrkRecoTrk.cxx:400
virtual double arrivalTime(double fltL) const
Definition: TrkRep.cxx:192
bool setFlipAndDrop(TrkRecoTrk &, bool allowFlips, bool allowDrops) const
TrkRecoTrk * makeTrack(const TrkExchangePar &helix, const double chi2, const TrkContext &, double trackT0) const
char * c_str(Index i)
Definition: EvtCyclic3.cc:252