CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CkfDebugger.cc
Go to the documentation of this file.
1 #include "CkfDebugger.h"
12 #include "TSOSFromSimHitFactory.h"
17 
19 
25 // #include "RecoTracker/TkDetLayers/interface/PixelBarrelLayer.h"
26 
29 
30 #include <iostream>
31 #include <iomanip>
32 #include <sstream>
33 
34 using namespace std;
35 
36 inline DetId gluedId(const DetId& du) {
37  unsigned int mask = ~3; // mask the last two bits
38  return DetId(du.rawId() & mask);
39 }
40 
42  : trackerHitAssociatorConfig_(std::move(iC)), totSeeds(0) {
43  file = new TFile("out.root", "recreate");
44  hchi2seedAll = new TH1F("hchi2seedAll", "hchi2seedAll", 2000, 0, 200);
45  hchi2seedProb = new TH1F("hchi2seedProb", "hchi2seedProb", 2000, 0, 200);
46 
48  es.get<TrackerDigiGeometryRecord>().get(tracker);
49  theTrackerGeom = &(*tracker);
50 
52  es.get<IdealMagneticFieldRecord>().get(theField);
53  theMagField = &(*theField);
54 
55  //Retrieve tracker topology from geometry
57  es.get<IdealGeometryRecord>().get(tTopoHand);
58  theTopo = tTopoHand.product();
59 
61  es.get<NavigationSchoolRecord>().get("SimpleNavigationSchool", nav);
62  theNavSchool = nav.product();
63 
64  for (int i = 0; i != 17; i++) {
65  dump.push_back(0);
66  }
67 
68  std::stringstream title;
69  for (int i = 0; i != 6; i++)
70  for (int j = 0; j != 9; j++) {
71  if (i == 0 && j > 2)
72  break;
73  if (i == 1 && j > 1)
74  break;
75  if (i == 2 && j > 3)
76  break;
77  if (i == 3 && j > 2)
78  break;
79  if (i == 4 && j > 5)
80  break;
81  if (i == 5 && j > 8)
82  break;
83  dump2[pair<int, int>(i, j)] = 0;
84  dump3[pair<int, int>(i, j)] = 0;
85  dump4[pair<int, int>(i, j)] = 0;
86  dump5[pair<int, int>(i, j)] = 0;
87  dump6[pair<int, int>(i, j)] = 0;
88  title.str("");
89  title << "pullX_" << i + 1 << "-" << j + 1 << "_sh-rh";
90  hPullX_shrh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
91  title.str("");
92  title << "pullY_" << i + 1 << "-" << j + 1 << "_sh-rh";
93  hPullY_shrh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
94  title.str("");
95  title << "pullX_" << i + 1 << "-" << j + 1 << "_sh-st";
96  hPullX_shst[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
97  title.str("");
98  title << "pullY_" << i + 1 << "-" << j + 1 << "_sh-st";
99  hPullY_shst[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
100  title.str("");
101  title << "pullX_" << i + 1 << "-" << j + 1 << "_st-rh";
102  hPullX_strh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
103  title.str("");
104  title << "pullY_" << i + 1 << "-" << j + 1 << "_st-rh";
105  hPullY_strh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
106  title.str("");
107  title << "PullGP_X_" << i + 1 << "-" << j + 1 << "_sh-st";
108  hPullGP_X_shst[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
109  title.str("");
110  title << "PullGP_Y_" << i + 1 << "-" << j + 1 << "_sh-st";
111  hPullGP_Y_shst[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
112  title.str("");
113  title << "PullGP_Z_" << i + 1 << "-" << j + 1 << "_sh-st";
114  hPullGP_Z_shst[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
115  if (((i == 2 || i == 4) && (j == 0 || j == 1)) || (i == 3 || i == 5)) {
116  title.str("");
117  title << "pullM_" << i + 1 << "-" << j + 1 << "_sh-rh";
118  hPullM_shrh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
119  title.str("");
120  title << "pullS_" << i + 1 << "-" << j + 1 << "_sh-rh";
121  hPullS_shrh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
122  title.str("");
123  title << "pullM_" << i + 1 << "-" << j + 1 << "_sh-st";
124  hPullM_shst[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
125  title.str("");
126  title << "pullS_" << i + 1 << "-" << j + 1 << "_sh-st";
127  hPullS_shst[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
128  title.str("");
129  title << "pullM_" << i + 1 << "-" << j + 1 << "_st-rh";
130  hPullM_strh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
131  title.str("");
132  title << "pullS_" << i + 1 << "-" << j + 1 << "_st-rh";
133  hPullS_strh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
134  }
135  }
136 
137  hPullGPXvsGPX_shst = new TH2F("PullGPXvsGPX_shst", "PullGPXvsGPX_shst", 1000, -50, 50, 100, -50, 50);
138  hPullGPXvsGPY_shst = new TH2F("PullGPXvsGPY_shst", "PullGPXvsGPY_shst", 1000, -50, 50, 100, -50, 50);
139  hPullGPXvsGPZ_shst = new TH2F("PullGPXvsGPZ_shst", "PullGPXvsGPZ_shst", 1000, -50, 50, 200, -100, 100);
140  hPullGPXvsGPr_shst = new TH2F("PullGPXvsGPr_shst", "PullGPXvsGPr_shst", 1000, -50, 50, 300, -150, 150);
141  hPullGPXvsGPeta_shst = new TH2F("PullGPXvsGPeta_shst", "PullGPXvsGPeta_shst", 1000, -50, 50, 50, -2.5, 2.5);
142  hPullGPXvsGPphi_shst = new TH2F("PullGPXvsGPphi_shst", "PullGPXvsGPphi_shst", 1000, -50, 50, 63, 0, 6.3);
143 
144  seedWithDelta = 0;
145  problems = 0;
146  no_sim_hit = 0;
147  no_layer = 0;
148  layer_not_found = 0;
149  det_not_found = 0;
150  chi2gt30 = 0;
151  chi2gt30delta = 0;
152  chi2gt30deltaSeed = 0;
153  chi2ls30 = 0;
155  no_component = 0;
156  only_one_component = 0;
157  matched_not_found = 0;
161  propagation = 0;
162  other = 0;
163  totchi2gt30 = 0;
164 }
165 
167  edm::LogVerbatim("CkfDebugger") << "\nEVENT #" << iEvent.id();
168 
170  iEvent, trackerHitAssociatorConfig_); //delete deleteHitAssociator() in TrackCandMaker.cc
171 
172  std::map<unsigned int, std::vector<PSimHit> >& theHitsMap = hitAssociator->SimHitMap;
173  idHitsMap.clear();
174 
175  for (std::map<unsigned int, std::vector<PSimHit> >::iterator it = theHitsMap.begin(); it != theHitsMap.end(); it++) {
176  for (std::vector<PSimHit>::iterator isim = it->second.begin(); isim != it->second.end(); ++isim) {
177  idHitsMap[isim->trackId()].push_back(&*isim);
178  }
179  }
180 
181  for (std::map<unsigned int, std::vector<PSimHit*> >::iterator it = idHitsMap.begin(); it != idHitsMap.end(); it++) {
182  sort(it->second.begin(), it->second.end(), [](auto* a, auto* b) { return a->timeOfFlight() < b->timeOfFlight(); });
183  for (std::vector<PSimHit*>::iterator isim = it->second.begin(); isim != it->second.end(); ++isim) {
184  const GeomDetUnit* detUnit = theTrackerGeom->idToDetUnit(DetId((*isim)->detUnitId()));
185  dumpSimHit(SimHit((*isim), detUnit));
186  }
187  }
188 }
189 
190 void CkfDebugger::dumpSimHit(const SimHit& hit) const {
191  GlobalPoint pos = hit.globalPosition();
192  edm::LogVerbatim("CkfDebugger") << "SimHit pos" << pos << " r=" << pos.perp() << " phi=" << pos.phi()
193  << " trackId=" << hit.trackId() << " particleType=" << hit.particleType()
194  << " pabs=" << hit.pabs() << " processType=" << hit.processType();
195 }
196 
198  const std::vector<TrajectoryMeasurement>& meas,
199  const MeasurementTrackerEvent* aMeasurementTracker,
200  const Propagator* propagator,
202  const TransientTrackingRecHitBuilder* aTTRHBuilder) {
203  LogTrace("CkfDebugger") << "\nnow in analyseCompatibleMeasurements";
204  LogTrace("CkfDebugger") << "number of input hits:" << meas.size();
205  for (std::vector<TrajectoryMeasurement>::const_iterator tmpIt = meas.begin(); tmpIt != meas.end(); tmpIt++) {
206  if (tmpIt->recHit()->isValid())
207  LogTrace("CkfDebugger") << "valid hit at position:" << tmpIt->recHit()->globalPosition();
208  }
210  theChi2 = estimator;
211  theMeasurementTracker = aMeasurementTracker;
212  theGeomSearchTracker = theMeasurementTracker->geometricSearchTracker();
213  theTTRHBuilder = aTTRHBuilder;
214  unsigned int trajId = 0;
215  if (!correctTrajectory(traj, trajId)) {
216  LogTrace("CkfDebugger") << "trajectory not correct";
217  return true;
218  } // only correct trajectories analysed
219  LogTrace("CkfDebugger") << "correct trajectory";
220 
221  if (traj.measurements().size() == 2) {
222  if (testSeed(traj.firstMeasurement().recHit(),
223  traj.lastMeasurement().recHit(),
224  traj.lastMeasurement().updatedState()) == -1) {
225  LogTrace("CkfDebugger") << "Seed has delta";
226  seedWithDelta++;
227  return false; //true;//false?
228  }
229  }
230 
231  //const PSimHit* correctHit = nextCorrectHit(traj, trajId);
232  //if ( correctHit == 0) return true; // no more simhits on this track
233  std::vector<const PSimHit*> correctHits = nextCorrectHits(traj, trajId);
234  if (correctHits.empty())
235  return true; // no more simhits on this track
236 
237  for (std::vector<const PSimHit*>::iterator corHit = correctHits.begin(); corHit != correctHits.end(); corHit++) {
238  for (std::vector<TM>::const_iterator i = meas.begin(); i != meas.end(); i++) {
239  if (correctMeas(*i, *corHit)) {
240  LogTrace("CkfDebugger") << "Correct hit found at position " << i - meas.begin();
241  return true;
242  }
243  }
244  }
245 
246  //debug why the first hit in correctHits is not found
247  //FIXME should loop over all hits
248  const PSimHit* correctHit = *(correctHits.begin());
249 
250  // correct hit not found
251  edm::LogVerbatim("CkfDebugger") << std::endl
252  << "CkfDebugger: problem found: correct hit not found by findCompatibleMeasurements";
253  edm::LogVerbatim("CkfDebugger") << "The correct hit position is " << position(correctHit) << " lp "
254  << correctHit->localPosition();
255  edm::LogVerbatim("CkfDebugger") << "The size of the meas vector is " << meas.size();
256  dump[0]++;
257  problems++;
258 
259  for (std::vector<TM>::const_iterator i = meas.begin(); i != meas.end(); i++) {
260  edm::LogVerbatim("CkfDebugger") << "Is the hit valid? " << i->recHit()->isValid();
261  if (i->recHit()->isValid()) {
262  edm::LogVerbatim("CkfDebugger") << "RecHit at " << i->recHit()->globalPosition() << " layer "
263  << ((i->recHit()->det()->geographicalId().rawId() >> 16) & 0xF) << " subdet "
264  << i->recHit()->det()->geographicalId().subdetId() << " Chi2 " << i->estimate();
265  } else if (i->recHit()->det() == nullptr) {
266  edm::LogVerbatim("CkfDebugger") << "Invalid RecHit returned with zero Det pointer";
267  } else if (i->recHit()->det() == det(correctHit)) {
268  edm::LogVerbatim("CkfDebugger") << "Invalid hit returned in correct Det";
269  } else {
270  edm::LogVerbatim("CkfDebugger") << "Invalid hit returned in Det at gpos " << i->recHit()->det()->position()
271  << " correct Det is at " << det(correctHit)->position();
272  }
273  }
274 
275  //Look if the correct RecHit exists
276  std::pair<CTTRHp, double> correctRecHit = analyseRecHitExistance(*correctHit, traj.lastMeasurement().updatedState());
277  if (correctRecHit.first == nullptr) {
278  //the hit does not exist or is uncorrectly matched
279  if (fabs(correctRecHit.second - 0) < 0.01) {
280  dump[1]++;
281  } //other
282  if (fabs(correctRecHit.second + 1) < 0.01) {
283  dump[8]++;
284  } //propagation
285  if (fabs(correctRecHit.second + 2) < 0.01) {
286  dump[9]++;
287  } //No component is found
288  if (fabs(correctRecHit.second + 3) < 0.01) {
289  dump[10]++;
290  } //Partner measurementDet not found
291  if (fabs(correctRecHit.second + 4) < 0.01) {
292  dump[11]++;
293  } //glued MeasurementDet not found
294  if (fabs(correctRecHit.second + 5) < 0.01) {
295  dump[12]++;
296  } //matched not found
297  if (fabs(correctRecHit.second + 6) < 0.01) {
298  dump[13]++;
299  } //Matched not associated
300  if (fabs(correctRecHit.second + 7) < 0.01) {
301  dump[14]++;
302  } //Only one component is found
303  if (fabs(correctRecHit.second + 8) < 0.01) {
304  dump[15]++;
305  } //not found (is not a glued det)
306  } else {
307  //the hit exists: why wasn't it found?
308  int result = analyseRecHitNotFound(traj, correctRecHit.first);
309  if (result == 5) {
310  if (correctRecHit.second > 30) {
311  edm::LogVerbatim("CkfDebugger") << "Outling RecHit at pos=" << correctRecHit.first->globalPosition()
312  << " from SimHit at pos=" << position(correctHit)
313  << " det=" << correctHit->detUnitId()
314  << " process=" << correctHit->processType();
315  if (hasDelta(correctHit)) {
316  edm::LogVerbatim("CkfDebugger") << "there are deltas on this det";
317  chi2gt30delta++;
318  dump5[pair<int, int>((correctRecHit.first->det()->geographicalId().subdetId() - 1),
319  (layer(correctRecHit.first->det())) - 1)]++;
320  } else {
321  edm::LogVerbatim("CkfDebugger") << "no deltas on this det";
322  dump[5]++;
323  chi2gt30++;
324  dump3[pair<int, int>((correctRecHit.first->det()->geographicalId().subdetId() - 1),
325  (layer(correctRecHit.first->det())) - 1)]++;
326  CTTRHp h1 = traj.measurements()[0].recHit();
327  CTTRHp h2 = traj.measurements()[1].recHit();
328  TSOS t = traj.measurements()[1].updatedState();
329  double chi2 = testSeed(h1, h2, t);
330  if (chi2 == -1) {
331  edm::LogVerbatim("CkfDebugger") << "there were deltas in the seed";
333  } else {
334  hchi2seedProb->Fill(chi2);
335  edm::LogVerbatim("CkfDebugger") << "no deltas in the seed. What is wrong?";
336 
338  correctRecHit.first->det()->surface());
339  TSOS simDetState =
341 
342  if (true /*detState.globalMomentum().y()>0*/) {
343  int subdetId = correctRecHit.first->det()->geographicalId().subdetId();
344  int layerId = layer(correctRecHit.first->det());
345 
346  LogTrace("CkfDebugger") << "position(correctHit)=" << position(correctHit);
347  LogTrace("CkfDebugger") << "correctRecHit.first->globalPosition()="
348  << correctRecHit.first->globalPosition();
349  LogTrace("CkfDebugger") << "detState.globalPosition()=" << detState.globalPosition();
350  LogTrace("CkfDebugger") << "simDetState.globalPosition()=" << simDetState.globalPosition();
351 
352  LogTrace("CkfDebugger") << "correctHit->localPosition()=" << correctHit->localPosition();
353  LogTrace("CkfDebugger") << "correctRecHit.first->localPosition()="
354  << correctRecHit.first->localPosition();
355  LogTrace("CkfDebugger") << "correctRecHit.first->localPositionError()="
356  << correctRecHit.first->localPositionError();
357  LogTrace("CkfDebugger") << "detState.localPosition()=" << detState.localPosition();
358  LogTrace("CkfDebugger") << "detState.localError().positionError()="
359  << detState.localError().positionError();
360  LogTrace("CkfDebugger") << "simDetState.localPosition()=" << simDetState.localPosition();
361  LogTrace("CkfDebugger") << "simDetState.localError().positionError()="
362  << simDetState.localError().positionError();
363  double pullx_shrh = (correctHit->localPosition().x() - correctRecHit.first->localPosition().x()) /
364  sqrt(correctRecHit.first->localPositionError().xx());
365  double pully_shrh = 0;
366  if (correctRecHit.first->localPositionError().yy() != 0)
367  pully_shrh = (correctHit->localPosition().y() - correctRecHit.first->localPosition().y()) /
368  sqrt(correctRecHit.first->localPositionError().yy());
369  double pullx_shst = (correctHit->localPosition().x() - simDetState.localPosition().x()) /
370  sqrt(simDetState.localError().positionError().xx());
371  double pully_shst = (correctHit->localPosition().y() - simDetState.localPosition().y()) /
372  sqrt(simDetState.localError().positionError().yy());
373 
374  LogTrace("CkfDebugger") << "pullx(sh-rh)=" << pullx_shrh;
375  LogTrace("CkfDebugger") << "pully(sh-rh)=" << pully_shrh;
376  LogTrace("CkfDebugger") << "pullx(sh-st)=" << pullx_shst;
377  LogTrace("CkfDebugger") << "pully(sh-st)=" << pully_shst;
378 
379  LogTrace("CkfDebugger") << "pullx(st-rh)="
380  << (detState.localPosition().x() - correctRecHit.first->localPosition().x()) /
381  sqrt(correctRecHit.first->localPositionError().xx() +
382  detState.localError().positionError().xx());
383 
384  std::pair<double, double> pulls = computePulls(correctRecHit.first, detState);
385  if (subdetId > 0 && subdetId < 7 && layerId > 0 && layerId < 10) {
386  stringstream title;
387  title.str("");
388  title << "pullX_" << subdetId << "-" << layerId << "_sh-rh";
389  hPullX_shrh[title.str()]->Fill(pullx_shrh);
390  title.str("");
391  title << "pullY_" << subdetId << "-" << layerId << "_sh-rh";
392  hPullY_shrh[title.str()]->Fill(pully_shrh);
393  title.str("");
394  title << "pullX_" << subdetId << "-" << layerId << "_sh-st";
395  hPullX_shst[title.str()]->Fill(pullx_shst);
396  title.str("");
397  title << "pullY_" << subdetId << "-" << layerId << "_sh-st";
398  hPullY_shst[title.str()]->Fill(pully_shst);
399  title.str("");
400  title << "pullX_" << subdetId << "-" << layerId << "_st-rh";
401  hPullX_strh[title.str()]->Fill(pulls.first);
402  title.str("");
403  title << "pullY_" << subdetId << "-" << layerId << "_st-rh";
404  hPullY_strh[title.str()]->Fill(pulls.second);
405 
406  GlobalPoint shGPos = position(correctHit);
407  GlobalPoint stGPos = simDetState.globalPosition();
408  GlobalError stGPosErr = simDetState.cartesianError().position();
409  double pullGPx = (shGPos.x() - stGPos.x()) / sqrt(stGPosErr.cxx());
410  title.str("");
411  title << "PullGP_X_" << subdetId << "-" << layerId << "_sh-st";
412  hPullGP_X_shst[title.str()]->Fill(pullGPx);
413  title.str("");
414  title << "PullGP_Y_" << subdetId << "-" << layerId << "_sh-st";
415  hPullGP_Y_shst[title.str()]->Fill((shGPos.y() - stGPos.y()) / sqrt(stGPosErr.cyy()));
416  title.str("");
417  title << "PullGP_Z_" << subdetId << "-" << layerId << "_sh-st";
418  hPullGP_Z_shst[title.str()]->Fill((shGPos.z() - stGPos.z()) / sqrt(stGPosErr.czz()));
419 
420  if (subdetId == 3 && layerId == 1) {
421  hPullGPXvsGPX_shst->Fill(pullGPx, shGPos.x());
422  hPullGPXvsGPY_shst->Fill(pullGPx, shGPos.y());
423  hPullGPXvsGPZ_shst->Fill(pullGPx, shGPos.z());
424  hPullGPXvsGPr_shst->Fill(pullGPx, shGPos.mag());
425  hPullGPXvsGPeta_shst->Fill(pullGPx, shGPos.eta());
426  hPullGPXvsGPphi_shst->Fill(pullGPx, shGPos.phi());
427  }
428  if (dynamic_cast<const SiStripMatchedRecHit2D*>(correctRecHit.first->hit())) {
429  LogTrace("CkfDebugger") << "MONO HIT";
430  auto m = dynamic_cast<const SiStripMatchedRecHit2D*>(correctRecHit.first->hit())->monoHit();
431  CTTRHp tMonoHit = theTTRHBuilder->build(&m);
432  const PSimHit sMonoHit = *(hitAssociator->associateHit(*tMonoHit->hit()).begin());
434  tMonoHit->det()->surface());
435  double pullM_shrh = (sMonoHit.localPosition().x() - tMonoHit->localPosition().x()) /
436  sqrt(tMonoHit->localPositionError().xx());
437  double pullM_shst = (sMonoHit.localPosition().x() - monoState.localPosition().x()) /
438  sqrt(monoState.localError().positionError().xx());
439  std::pair<double, double> pullsMono = computePulls(tMonoHit, monoState);
440  title.str("");
441  title << "pullM_" << subdetId << "-" << layerId << "_sh-rh";
442  hPullM_shrh[title.str()]->Fill(pullM_shrh);
443  title.str("");
444  title << "pullM_" << subdetId << "-" << layerId << "_sh-st";
445  hPullM_shst[title.str()]->Fill(pullM_shst);
446  title.str("");
447  title << "pullM_" << subdetId << "-" << layerId << "_st-rh";
448  hPullM_strh[title.str()]->Fill(pullsMono.first);
449 
450  LogTrace("CkfDebugger") << "STEREO HIT";
451  auto s = dynamic_cast<const SiStripMatchedRecHit2D*>(correctRecHit.first->hit())->stereoHit();
452  CTTRHp tStereoHit = theTTRHBuilder->build(&s);
453  const PSimHit sStereoHit = *(hitAssociator->associateHit(*tStereoHit->hit()).begin());
455  tStereoHit->det()->surface());
456  double pullS_shrh = (sStereoHit.localPosition().x() - tStereoHit->localPosition().x()) /
457  sqrt(tStereoHit->localPositionError().xx());
458  double pullS_shst = (sStereoHit.localPosition().x() - stereoState.localPosition().x()) /
459  sqrt(stereoState.localError().positionError().xx());
460  std::pair<double, double> pullsStereo = computePulls(tStereoHit, stereoState);
461  title.str("");
462  title << "pullS_" << subdetId << "-" << layerId << "_sh-rh";
463  hPullS_shrh[title.str()]->Fill(pullS_shrh);
464  title.str("");
465  title << "pullS_" << subdetId << "-" << layerId << "_sh-st";
466  hPullS_shst[title.str()]->Fill(pullS_shst);
467  title.str("");
468  title << "pullS_" << subdetId << "-" << layerId << "_st-rh";
469  hPullS_strh[title.str()]->Fill(pullsStereo.first);
470  }
471  } else
472  edm::LogVerbatim("CkfDebugger")
473  << "unexpected result: wrong det or layer id " << subdetId << " " << layerId << " "
474  << correctRecHit.first->det()->geographicalId().rawId();
475  }
476  }
477  }
478  } else {
479  edm::LogVerbatim("CkfDebugger") << "unexpected result " << correctRecHit.second;
480  dump[6]++;
481  chi2ls30++;
482  }
483  } else
484  dump[result]++;
485  if (result == 3) {
486  dump2[pair<int, int>((correctRecHit.first->det()->geographicalId().subdetId() - 1),
487  (layer(correctRecHit.first->det())) - 1)]++;
488  }
489  if (result == 4) {
490  dump4[pair<int, int>((correctRecHit.first->det()->geographicalId().subdetId() - 1),
491  (layer(correctRecHit.first->det())) - 1)]++;
492  }
493  if (correctRecHit.second > 30) {
494  dump[7]++;
495  totchi2gt30++;
496  }
497  }
498  return false;
499 }
500 
501 bool CkfDebugger::correctTrajectory(const Trajectory& traj, unsigned int& trajId) const {
502  LogTrace("CkfDebugger") << "now in correctTrajectory";
503  Trajectory::RecHitContainer hits = traj.recHits();
504 
505  std::vector<SimHitIdpr> currentTrackId = hitAssociator->associateHitId(*hits.front()->hit());
506  if (currentTrackId.empty())
507  return false;
508 
509  for (Trajectory::RecHitContainer::const_iterator rh = hits.begin(); rh != hits.end(); ++rh) {
510  //if invalid hit exit
511  if (!(*rh)->hit()->isValid()) {
512  //LogTrace("CkfDebugger") << "invalid hit" ;
513  return false;
514  }
515 
516  //if hits from deltas exit
517  bool nogoodhit = true;
518  std::vector<PSimHit> assSimHits = hitAssociator->associateHit(*(*rh)->hit());
519  for (std::vector<PSimHit>::iterator shit = assSimHits.begin(); shit != assSimHits.end(); shit++) {
520  if (goodSimHit(*shit))
521  nogoodhit = false;
522  }
523  if (nogoodhit)
524  return false;
525 
526  //all hits must be associated to the same sim track
527  bool test = true;
528  std::vector<SimHitIdpr> nextTrackId = hitAssociator->associateHitId(*(*rh)->hit());
529  for (std::vector<SimHitIdpr>::iterator i = currentTrackId.begin(); i != currentTrackId.end(); i++) {
530  for (std::vector<SimHitIdpr>::iterator j = nextTrackId.begin(); j != nextTrackId.end(); j++) {
531  if (i->first == j->first)
532  test = false;
533  //LogTrace("CkfDebugger") << "valid " << *i << " " << *j ;
534  trajId = j->first;
535  }
536  }
537  if (test) { /*LogTrace("CkfDebugger") << "returning false" ;*/
538  return false;
539  }
540  // std::vector<PSimHit*> simTrackHits = idHitsMap[trajId];
541  // if (!goodSimHit(simTrackHits.))
542  }
543  //LogTrace("CkfDebugger") << "returning true" ;
544  return true;
545 }
546 
548  LogTrace("CkfDebugger") << "now in assocTrackId";
549 
550  if (!rechit->hit()->isValid()) {
551  return -1;
552  }
553 
554  std::vector<SimHitIdpr> ids = hitAssociator->associateHitId(*rechit->hit());
555  if (!ids.empty()) {
556  return ids[0].first; //FIXME if size>1!!
557  } else {
558  return -1;
559  }
560 }
561 
562 vector<const PSimHit*> CkfDebugger::nextCorrectHits(const Trajectory& traj, unsigned int& trajId) {
563  std::vector<const PSimHit*> result;
564  // find the component of the RecHit at largest distance from origin (FIXME: should depend on propagation direction)
565  LogTrace("CkfDebugger") << "now in nextCorrectHits";
567  TransientTrackingRecHit::RecHitContainer comp = lastRecHit->transientHits();
568  if (!comp.empty()) {
569  float maxR = 0;
570  for (TransientTrackingRecHit::RecHitContainer::const_iterator ch = comp.begin(); ch != comp.end(); ++ch) {
571  if ((*ch)->globalPosition().mag() > maxR)
572  lastRecHit = *ch;
573  maxR = (*ch)->globalPosition().mag();
574  }
575  }
576  edm::LogVerbatim("CkfDebugger") << "CkfDebugger: lastRecHit is at gpos " << lastRecHit->globalPosition() << " layer "
577  << layer((lastRecHit->det())) << " subdet "
578  << lastRecHit->det()->geographicalId().subdetId();
579 
580  //find the simHits associated to the recHit
581  const std::vector<PSimHit>& pSimHitVec = hitAssociator->associateHit(*lastRecHit->hit());
582  for (std::vector<PSimHit>::const_iterator shit = pSimHitVec.begin(); shit != pSimHitVec.end(); shit++) {
583  const GeomDetUnit* detUnit = theTrackerGeom->idToDetUnit(DetId(shit->detUnitId()));
584  LogTrace("CkfDebugger") << "from hitAssociator SimHits are at GP=" << detUnit->toGlobal(shit->localPosition())
585  << " traId=" << shit->trackId() << " particleType " << shit->particleType()
586  << " pabs=" << shit->pabs() << " detUnitId=" << shit->detUnitId() << " layer "
587  << layer((det(&*shit))) << " subdet " << det(&*shit)->geographicalId().subdetId();
588  }
589 
590  //choose the simHit from the same track that has the highest tof
591  const PSimHit* lastPSH = nullptr;
592  if (!pSimHitVec.empty()) {
593  float maxTOF = 0;
594  for (std::vector<PSimHit>::const_iterator ch = pSimHitVec.begin(); ch != pSimHitVec.end(); ++ch) {
595  if ((ch->trackId() == trajId) && (ch->timeOfFlight() > maxTOF) && (goodSimHit(*ch))) {
596  lastPSH = &*ch;
597  maxTOF = lastPSH->timeOfFlight();
598  }
599  }
600  } else
601  return result; //return empty vector: no more hits on the sim track
602  if (lastPSH == nullptr)
603  return result; //return empty vector: no more good hits on the sim track
604  edm::LogVerbatim("CkfDebugger") << "CkfDebugger: corresponding SimHit is at gpos " << position(&*lastPSH);
605 
606  //take the simHits on the simTrack that are in the nextLayer (could be > 1 if overlap or matched)
607  std::vector<PSimHit*> trackHits = idHitsMap[trajId];
608  if (fabs((double)(trackHits.back()->detUnitId() - lastPSH->detUnitId())) < 1)
609  return result; //end of sim track
610  std::vector<PSimHit*>::iterator currentIt = trackHits.end();
611  for (std::vector<PSimHit*>::iterator it = trackHits.begin(); it != trackHits.end(); it++) {
612  if (goodSimHit(**it) && //good hit
613  (lastPSH->timeOfFlight() < (*it)->timeOfFlight()) && //greater tof
614  //( fabs((double)((*it)->detUnitId()-(lastPSH->detUnitId()) ))>1) && //not components of the same matched hit
615  ((det(lastPSH)->geographicalId().subdetId() != det(*it)->geographicalId().subdetId()) ||
616  (layer(det(lastPSH)) != layer(det(*it)))) //change layer or detector(tib,tob,...)
617  ) {
618  edm::LogVerbatim("CkfDebugger") << "Next good PSimHit is at gpos " << position(*it);
619  result.push_back(*it);
620  currentIt = it;
621  break;
622  }
623  }
624  bool samelayer = true;
625  if (currentIt != (trackHits.end() - 1) && currentIt != trackHits.end()) {
626  for (std::vector<PSimHit*>::iterator nextIt = currentIt; (samelayer && nextIt != trackHits.end()); nextIt++) {
627  if (goodSimHit(**nextIt)) {
628  if ((det(*nextIt)->geographicalId().subdetId() == det(*currentIt)->geographicalId().subdetId()) &&
629  (layer(det(*nextIt)) == layer(det(*currentIt)))) {
630  result.push_back(*nextIt);
631  } else
632  samelayer = false;
633  }
634  }
635  }
636 
637  return result;
638 }
639 
640 bool CkfDebugger::goodSimHit(const PSimHit& sh) const {
641  if (sh.pabs() > 0.9)
642  return true; // GeV, reject delta rays from association
643  else
644  return false;
645 }
646 
647 bool CkfDebugger::associated(CTTRHp rechit, const PSimHit& pSimHit) const {
648  LogTrace("CkfDebugger") << "now in associated";
649 
650  if (!rechit->isValid())
651  return false;
652  // LogTrace("CkfDebugger") << "rec hit valid" ;
653  const std::vector<PSimHit>& pSimHitVec = hitAssociator->associateHit(*rechit->hit());
654  // LogTrace("CkfDebugger") << "size=" << pSimHitVec.size() ;
655  for (std::vector<PSimHit>::const_iterator shit = pSimHitVec.begin(); shit != pSimHitVec.end(); shit++) {
656  //const GeomDetUnit* detUnit = theTrackerGeom->idToDetUnit( DetId(shit->detUnitId()));
657  // LogTrace("CkfDebugger") << "pSimHit.timeOfFlight()=" << pSimHit.timeOfFlight()
658  // << " pSimHit.pabs()=" << pSimHit.pabs() << " GP=" << position(&pSimHit);
659  // LogTrace("CkfDebugger") << "(*shit).timeOfFlight()=" << (*shit).timeOfFlight()
660  // << " (*shit).pabs()=" << (*shit).pabs() << " GP=" << detUnit->toGlobal( shit->localPosition());
661  if ((fabs((*shit).timeOfFlight() - pSimHit.timeOfFlight()) < 1e-9) &&
662  (fabs((*shit).pabs() - pSimHit.pabs()) < 1e-9))
663  return true;
664  }
665  return false;
666 }
667 
668 bool CkfDebugger::correctMeas(const TM& tm, const PSimHit* correctHit) const {
669  LogTrace("CkfDebugger") << "now in correctMeas";
670  const CTTRHp& recHit = tm.recHit();
671  if (recHit->isValid())
672  LogTrace("CkfDebugger") << "hit at position:" << recHit->globalPosition();
673  TransientTrackingRecHit::RecHitContainer comp = recHit->transientHits();
674  if (comp.empty()) {
675  // LogTrace("CkfDebugger") << "comp.empty()==true" ;
676  return associated(recHit, *correctHit);
677  } else {
678  for (TransientTrackingRecHit::RecHitContainer::const_iterator ch = comp.begin(); ch != comp.end(); ++ch) {
679  if (associated(recHit, *correctHit)) {
680  // check if the other components are associated to the same trackId
681  for (TransientTrackingRecHit::RecHitContainer::const_iterator ch2 = comp.begin(); ch2 != comp.end(); ++ch2) {
682  if (ch2 == ch)
683  continue;
685  // LogTrace("CkfDebugger") << "correctHit->trackId()=" << correctHit->trackId() ;
686  bool test = true;
687  std::vector<SimHitIdpr> ids = hitAssociator->associateHitId(*(*ch2)->hit());
688  for (std::vector<SimHitIdpr>::iterator j = ids.begin(); j != ids.end(); j++) {
689  // LogTrace("CkfDebugger") << "j=" <<j->first;
690  if (correctHit->trackId() == j->first) {
691  test = false;
692  // LogTrace("CkfDebugger") << correctHit->trackId()<< " " <<j->first;
693  }
694  }
695  if (assocTrackId(*ch2) != ((int)(correctHit->trackId()))) {
696  LogTrace("CkfDebugger") << "returning false 1"; /*return false;*/
697  } //fixme
698  if (test) {
699  // LogTrace("CkfDebugger") << "returning false 2" ;
700  return false; // not all components from same simtrack
701  }
702  // if (assocTrackId( **ch2) != ((int)( correctHit->trackId())) ) {
703  // return false; // not all components from same simtrack
704  // }
705  }
706  return true; // if all components from same simtrack
707  }
708  }
709  return false;
710  }
711 }
712 
713 //this checks only if there is the rechit on the det where the sim hit is
714 pair<CTTRHp, double> CkfDebugger::analyseRecHitExistance(const PSimHit& sh, const TSOS& startingState) {
715  LogTrace("CkfDebugger") << "now in analyseRecHitExistance";
716 
717 #if 0
718  std::pair<CTTRHp, double> result;
719 
720  const MeasurementDet* simHitDet = theMeasurementTracker->idToDet( DetId( sh.detUnitId()));
721  TSOS simHitState = TSOSFromSimHitFactory()(sh, *det(&sh), *theMagField);
722  MeasurementDet::RecHitContainer recHits = simHitDet->recHits( simHitState);//take all hits from det
723 
724  //check if the hit is not present or is a problem of association
725  TSOS firstDetState = theForwardPropagator->propagate( startingState, det(&sh)->surface());
726  if (!firstDetState.isValid()) {
727  edm::LogVerbatim("CkfDebugger") << "CkfDebugger: propagation failed from state " << startingState << " to first det surface "
728  << position(&sh) ;
729  propagation++;
730  return std::pair<CTTRHp, double>((CTTRHp)(0),-1);
731  }
732 
733  bool found = false;
734  for ( MeasurementDet::RecHitContainer::const_iterator rh = recHits.begin(); rh != recHits.end(); rh++) {
735  if ( associated( *rh, sh)) {
736  found = true;
737  result = std::pair<CTTRHp, double>(*rh,theChi2->estimate( firstDetState, **rh).second);
738  edm::LogVerbatim("CkfDebugger") << "CkfDebugger: A RecHit associated to the correct Simhit exists at lpos "
739  << (**rh).localPosition()
740  << " gpos " << (**rh).globalPosition()
741  << " layer " << layer((**rh).det())
742  << " subdet " << (**rh).det()->geographicalId().subdetId()
743  << " Chi2 " << theChi2->estimate( firstDetState, **rh).second;
744  }
745  }
746  if (!found) {
747  edm::LogVerbatim("CkfDebugger") << "CkfDebugger: there is no RecHit associated to the correct SimHit." ;
748  edm::LogVerbatim("CkfDebugger") << " There are " << recHits.size() << " RecHits in the simHit DetUnit" ;
749  edm::LogVerbatim("CkfDebugger") << "SH GP=" << position(&sh) << " subdet=" << det(&sh)->geographicalId().subdetId()
750  << " layer=" << layer(det(&sh)) ;
751  int y=0;
752  for (MeasurementDet::RecHitContainer::const_iterator rh = recHits.begin(); rh != recHits.end(); rh++)
753  edm::LogVerbatim("CkfDebugger") << "RH#" << y++ << " GP=" << (**rh).globalPosition() << " subdet=" << (**rh).det()->geographicalId().subdetId()
754  << " layer=" << layer((**rh).det()) ;
755  for ( MeasurementDet::RecHitContainer::const_iterator rh = recHits.begin(); rh != recHits.end(); rh++) {
756  edm::LogVerbatim("CkfDebugger") << "Non-associated RecHit at pos " << (**rh).localPosition() ;
757  }
758  }
759 
760  bool found2 = false;
761  const PSimHit* sh2;
762  StripSubdetector subdet( det(&sh)->geographicalId());
763  if (!subdet.glued()) {
764  edm::LogVerbatim("CkfDebugger") << "The DetUnit is not part of a GluedDet" ;
765  if (found) {
766  if (result.second>30){
767  LogTrace("CkfDebugger") << "rh->parameters()=" << result.first->parameters() ;
768  LogTrace("CkfDebugger") << "rh->parametersError()=" << result.first->parametersError() ;
769  MeasurementExtractor me(firstDetState);
770  AlgebraicVector r(result.first->parameters() - me.measuredParameters(*result.first));
771  LogTrace("CkfDebugger") << "me.measuredParameters(**rh)=" << me.measuredParameters(*result.first) ;
772  LogTrace("CkfDebugger") << "me.measuredError(**rh)=" << me.measuredError(*result.first) ;
773  AlgebraicSymMatrix R(result.first->parametersError() + me.measuredError(*result.first));
774  LogTrace("CkfDebugger") << "r=" << r ;
775  LogTrace("CkfDebugger") << "R=" << R ;
776  int ierr;
777  R.invert(ierr);
778  LogTrace("CkfDebugger") << "R(-1)=" << R ;
779  LogTrace("CkfDebugger") << "chi2=" << R.similarity(r) ;
780  }
781  return result;
782  }
783  else {
785  return std::pair<CTTRHp, double>((CTTRHp)(0),-8);//not found (is not a glued det)
786  }
787  } else {
788  edm::LogVerbatim("CkfDebugger") << "The DetUnit is part of a GluedDet" ;
789  DetId partnerDetId = DetId( subdet.partnerDetId());
790 
791  sh2 = pSimHit( sh.trackId(), partnerDetId);
792  if (sh2 == 0) {
793  edm::LogVerbatim("CkfDebugger") << "Partner DetUnit does not have a SimHit from the same track" ;
794  if (found) {
795  //projected rec hit
796  TrackingRecHitProjector<ProjectedRecHit2D> proj;
797  DetId gid = gluedId( subdet);
798  const MeasurementDet* gluedDet = theMeasurementTracker->idToDet( gid);
799  TSOS gluedTSOS = theForwardPropagator->propagate(startingState, gluedDet->geomDet().surface());
800  CTTRHp projHit = proj.project( *result.first,gluedDet->geomDet(),gluedTSOS).get();
801  //LogTrace("CkfDebugger") << proj.project( *result.first,gluedDet->geomDet(),gluedTSOS)->parameters() ;
802  //LogTrace("CkfDebugger") << projHit->parametersError() ;
803  double chi2 = theChi2->estimate(gluedTSOS, *proj.project( *result.first,gluedDet->geomDet(),gluedTSOS)).second;
804  return std::pair<CTTRHp, double>(projHit,chi2);
805  }
806  }
807  else {
808  edm::LogVerbatim("CkfDebugger") << "Partner DetUnit has a good SimHit at gpos " << position(sh2)
809  << " lpos " << sh2->localPosition() ;
810  //}
811 
812  const MeasurementDet* partnerDet = theMeasurementTracker->idToDet( partnerDetId);
813  if (partnerDet == 0) {
814  edm::LogVerbatim("CkfDebugger") << "Partner measurementDet not found!!!" ;
816  return std::pair<CTTRHp, double>((CTTRHp)(0),-3);
817  }
818  TSOS simHitState2 = TSOSFromSimHitFactory()(*sh2, *det(sh2), *theMagField);
819  MeasurementDet::RecHitContainer recHits2 = partnerDet->recHits( simHitState2);
820 
821  TSOS secondDetState = theForwardPropagator->propagate( startingState, det(sh2)->surface());
822  if (!secondDetState.isValid()) {
823  edm::LogVerbatim("CkfDebugger") << "CkfDebugger: propagation failed from state " << startingState << " to second det surface "
824  << position(sh2) ;
825  propagation++;
826  return std::pair<CTTRHp, double>((CTTRHp)(0),-1);
827  }
828 
829  for ( MeasurementDet::RecHitContainer::const_iterator rh = recHits2.begin(); rh != recHits2.end(); rh++) {
830  if ( associated( *rh, *sh2)) {
831  found2 = true;
832  edm::LogVerbatim("CkfDebugger") << "CkfDebugger: A RecHit associated to the correct Simhit exists at lpos "
833  << (**rh).localPosition()
834  << " gpos " << (**rh).globalPosition()
835  << " Chi2 " << theChi2->estimate( secondDetState, **rh).second
836  ;
837  }
838  }
839  if (!found2) {
840  edm::LogVerbatim("CkfDebugger") << "CkfDebugger: there is no RecHit associated to the correct SimHit." ;
841  LogTrace("CkfDebugger") << " There are " << recHits.size() << " RecHits in the simHit DetUnit" ;
842  for ( MeasurementDet::RecHitContainer::const_iterator rh = recHits.begin(); rh != recHits.end(); rh++) {
843  LogTrace("CkfDebugger") << "Non-associated RecHit at pos " << (**rh).localPosition() ;
844  }
845  }
846  }
847  }
848 
850  if (found && found2) {
851  // look in the glued det
852  DetId gid = gluedId( subdet);
853  const MeasurementDet* gluedDet = theMeasurementTracker->idToDet( gid);
854  if ( gluedDet == 0) {
855  edm::LogVerbatim("CkfDebugger") << "CkfDebugger ERROR: glued MeasurementDet not found!" ;
857  return std::pair<CTTRHp, double>((CTTRHp)(0),-4);
858  }
859 
860  TSOS gluedDetState = theForwardPropagator->propagate( startingState, gluedDet->surface());
861  if (!gluedDetState.isValid()) {
862  edm::LogVerbatim("CkfDebugger") << "CkfDebugger: propagation failed from state " << startingState << " to det surface "
863  << gluedDet->position() ;
864  propagation++;
865  return std::pair<CTTRHp, double>((CTTRHp)(0),-1);
866  }
867 
868  gluedHits = gluedDet->recHits( gluedDetState);
869  edm::LogVerbatim("CkfDebugger") << "CkfDebugger: the GluedDet returned " << gluedHits.size() << " hits" ;
870  if (gluedHits.size()==0){
871  edm::LogVerbatim("CkfDebugger") << "Found and associated mono and stereo recHits but not matched!!!" ;
873  return std::pair<CTTRHp, double>((CTTRHp)(0),-5);
874  }
875  bool found3 = false;
876  for ( MeasurementDet::RecHitContainer::const_iterator rh = gluedHits.begin(); rh != gluedHits.end(); rh++) {
877  if ( associated( *rh, sh) && associated( *rh, *sh2)) {
878  double chi2 = theChi2->estimate(gluedDetState, **rh).second;
879  edm::LogVerbatim("CkfDebugger") << "Matched hit at lpos " << (**rh).localPosition()
880  << " gpos " << (**rh).globalPosition()
881  << " has Chi2 " << chi2
882  ;
883  result = std::pair<CTTRHp, double>(&**rh,chi2);
884  found3 = true;
885  if (chi2>30){
886  LogTrace("CkfDebugger") << "rh->parameters()=" << (*rh)->parameters() ;
887  LogTrace("CkfDebugger") << "rh->parametersError()=" << (*rh)->parametersError() ;
888  MeasurementExtractor me(gluedDetState);
889  AlgebraicVector r((*rh)->parameters() - me.measuredParameters(**rh));
890  LogTrace("CkfDebugger") << "me.measuredParameters(**rh)=" << me.measuredParameters(**rh) ;
891  LogTrace("CkfDebugger") << "me.measuredError(**rh)=" << me.measuredError(**rh) ;
892  AlgebraicSymMatrix R((*rh)->parametersError() + me.measuredError(**rh));
893  LogTrace("CkfDebugger") << "r=" << r ;
894  LogTrace("CkfDebugger") << "R=" << R ;
895  int ierr;
896  R.invert(ierr);
897  LogTrace("CkfDebugger") << "R(-1)=" << R ;
898  LogTrace("CkfDebugger") << "chi2=" << R.similarity(r) ;
899  }
900  break;
901  }
902  }
903  if (found3) return result;
904  else {
905  edm::LogVerbatim("CkfDebugger") << "Found and associated mono and stereo recHits. Matched found but not associated!!!" ;
907  return std::pair<CTTRHp, double>((CTTRHp)(0),-6);
908  }
909  }
910  else if ( (found && !found2) || (!found && found2) ) {
911  edm::LogVerbatim("CkfDebugger") << "Only one component is found" ;
913  return std::pair<CTTRHp, double>((CTTRHp)(0),-7);
914  }
915  else {
916  edm::LogVerbatim("CkfDebugger") << "No component is found" ;
917  no_component++;
918  return std::pair<CTTRHp, double>((CTTRHp)(0),-2);
919  }
920  other++;
921 #endif
922  return std::pair<CTTRHp, double>((CTTRHp)(nullptr), 0); //other
923 }
924 
925 const PSimHit* CkfDebugger::pSimHit(unsigned int tkId, DetId detId) {
926  for (std::vector<PSimHit*>::iterator shi = idHitsMap[tkId].begin(); shi != idHitsMap[tkId].end(); ++shi) {
927  if ((*shi)->detUnitId() == detId.rawId() &&
928  //(shi)->trackId() == tkId &&
929  goodSimHit(**shi)) {
930  return (*shi);
931  }
932  }
933  return nullptr;
934 }
935 
936 int CkfDebugger::analyseRecHitNotFound(const Trajectory& traj, CTTRHp correctRecHit) {
937  unsigned int correctDetId = correctRecHit->det()->geographicalId().rawId();
938  int correctLayId = layer(correctRecHit->det());
939  LogTrace("CkfDebugger") << "correct layer id=" << correctLayId;
940 
941  TSOS currentState(traj.lastMeasurement().updatedState());
942  std::vector<const DetLayer*> nl =
943  theNavSchool->nextLayers(*traj.lastLayer(), *currentState.freeState(), traj.direction());
944  if (nl.empty()) {
945  edm::LogVerbatim("CkfDebugger") << "no compatible layers";
946  no_layer++;
947  return 2;
948  }
949 
950  TkLayerLess lless; //FIXME - was lless(traj.direction())
951  const DetLayer* detLayer = nullptr;
952  bool navLayerAfter = false;
953  bool test = false;
954  for (std::vector<const DetLayer*>::iterator il = nl.begin(); il != nl.end(); il++) {
955  if (dynamic_cast<const BarrelDetLayer*>(*il)) {
956  const BarrelDetLayer* pbl = dynamic_cast<const BarrelDetLayer*>(*il);
957  LogTrace("CkfDebugger") << "pbl->specificSurface().bounds().length()="
958  << pbl->specificSurface().bounds().length();
959  LogTrace("CkfDebugger") << "pbl->specificSurface().bounds().width()=" << pbl->specificSurface().bounds().width();
960  }
961  int layId = layer(((*(*il)->basicComponents().begin())));
962  LogTrace("CkfDebugger") << " subdet=" << (*(*il)->basicComponents().begin())->geographicalId().subdetId()
963  << "layer id=" << layId;
964  if (layId == correctLayId) {
965  test = true;
966  detLayer = &**il;
967  break;
968  }
969  if (lless(*il, theGeomSearchTracker->detLayer(correctRecHit->det()->geographicalId())))
970  navLayerAfter = true; //it is enough that only one layer is after the correct one?
971  }
972 
973  if (test) {
974  edm::LogVerbatim("CkfDebugger") << "correct layer taken into account. layer id: " << correctLayId;
975  } else if (navLayerAfter) {
976  edm::LogVerbatim("CkfDebugger") << "SimHit layer after the layers returned by Navigation.";
977  edm::LogVerbatim("CkfDebugger") << "Probably a missing SimHit.";
978  edm::LogVerbatim("CkfDebugger") << "check: " << (correctRecHit->det()->geographicalId().subdetId()) << " "
979  << (layer(correctRecHit->det()));
980  dump6[pair<int, int>((correctRecHit->det()->geographicalId().subdetId() - 1), (layer(correctRecHit->det())) - 1)]++;
981  no_sim_hit++;
982  return 16;
983  } else {
984  edm::LogVerbatim("CkfDebugger") << "correct layer NOT taken into account. correct layer id: " << correctLayId;
985  layer_not_found++;
986  return 3;
987  }
988 
990  std::vector<DetWithState> compatDets = detLayer->compatibleDets(currentState, *theForwardPropagator, *theChi2);
991  // LogTrace("CkfDebugger") << "DEBUGGER" ;
992  // LogTrace("CkfDebugger") << "runned compatDets." ;
993  // LogTrace("CkfDebugger") << "started from the following TSOS:" ;
994  // LogTrace("CkfDebugger") << currentState ;
995  // LogTrace("CkfDebugger") << "number of dets found=" << compatDets.size() ;
996  // for (std::vector<DetWithState>::iterator det=compatDets.begin();det!=compatDets.end();det++){
997  // unsigned int detId = det->first->geographicalId().rawId();
998  // LogTrace("CkfDebugger") << "detId=" << detId ;
999  // }
1000  bool test2 = false;
1001  for (std::vector<DetWithState>::iterator det = compatDets.begin(); det != compatDets.end(); det++) {
1002  unsigned int detId = det->first->geographicalId().rawId();
1003  // LogTrace("CkfDebugger") << "detId=" << detId
1004  // << "\ncorrectRecHit->det()->geographicalId().rawId()=" << correctRecHit->det()->geographicalId().rawId()
1005  // << "\ngluedId(correctRecHit->det()->geographicalId()).rawId()=" << gluedId(correctRecHit->det()->geographicalId()).rawId()
1006  // ;
1007  if (detId == gluedId(correctRecHit->det()->geographicalId()).rawId()) {
1008  test2 = true;
1009  break;
1010  }
1011  }
1012 
1013  if (test2) {
1014  edm::LogVerbatim("CkfDebugger") << "correct det taken into account. correctDetId is: " << correctDetId
1015  << ". please check chi2.";
1016  return 5;
1017  } else {
1018  edm::LogVerbatim("CkfDebugger") << "correct det NOT taken into account. correctDetId: " << correctDetId;
1019  det_not_found++;
1020  return 4;
1021  }
1022 }
1023 
1024 double CkfDebugger::testSeed(CTTRHp recHit1, CTTRHp recHit2, TSOS state) {
1025  //edm::LogVerbatim("CkfDebugger") << "CkfDebugger::testSeed";
1026  //test Deltas
1027  const std::vector<PSimHit>& pSimHitVec1 = hitAssociator->associateHit(*recHit1->hit());
1028  const std::vector<PSimHit>& pSimHitVec2 = hitAssociator->associateHit(*recHit2->hit());
1029 
1030  if (pSimHitVec1.empty() || pSimHitVec2.empty() || hasDelta(&(*pSimHitVec1.begin())) ||
1031  hasDelta(&(*pSimHitVec2.begin()))) {
1032  edm::LogVerbatim("CkfDebugger") << "Seed has delta or problems";
1033  return -1;
1034  }
1035 
1036  // LogTrace("CkfDebugger") << "state=\n" << state ;
1037  // double stlp1 = state.localParameters().vector()[0];
1038  // double stlp2 = state.localParameters().vector()[1];
1039  // double stlp3 = state.localParameters().vector()[2];
1040  // double stlp4 = state.localParameters().vector()[3];
1041  // double stlp5 = state.localParameters().vector()[4];
1042 
1043  if (!pSimHitVec2.empty()) {
1044  const PSimHit& simHit = *pSimHitVec2.begin();
1045 
1046  double shlp1 = -1 / simHit.momentumAtEntry().mag();
1047  double shlp2 = simHit.momentumAtEntry().x() / simHit.momentumAtEntry().z();
1048  double shlp3 = simHit.momentumAtEntry().y() / simHit.momentumAtEntry().z();
1049  double shlp4 = simHit.localPosition().x();
1050  double shlp5 = simHit.localPosition().y();
1052  v[0] = shlp1;
1053  v[1] = shlp2;
1054  v[2] = shlp3;
1055  v[3] = shlp4;
1056  v[4] = shlp5;
1057 
1058  // LogTrace("CkfDebugger") << "simHit.localPosition()=" << simHit.localPosition() ;
1059  // LogTrace("CkfDebugger") << "simHit.momentumAtEntry()=" << simHit.momentumAtEntry() ;
1060  // LogTrace("CkfDebugger") << "recHit2->localPosition()=" << recHit2->localPosition() ;
1061  // LogTrace("CkfDebugger") << "recHit2->localPositionError()=" << recHit2->localPositionError() ;
1062  // LogTrace("CkfDebugger") << "state.localPosition()=" << state.localPosition() ;
1063  // LogTrace("CkfDebugger") << "state.localError().positionError()=" << state.localError().positionError() ;
1064 
1065  // LogTrace("CkfDebugger") << "pullx(sh-rh)=" << (simHit.localPosition().x()-recHit2->localPosition().x())/sqrt(recHit2->localPositionError().xx()) ;
1066  // LogTrace("CkfDebugger") << "pullx(sh-st)=" << (simHit.localPosition().x()-state.localPosition().x())/sqrt(state.localError().positionError().xx()) ;
1067  // LogTrace("CkfDebugger") << "pullx(st-rh)=" << (state.localPosition().x()-recHit2->localPosition().x())/
1068  // sqrt(recHit2->localPositionError().xx()+state.localError().positionError().xx()) ;
1069 
1070  // LogTrace("CkfDebugger") << "local parameters" ;
1071  // LogTrace("CkfDebugger") << left;
1072  // LogTrace("CkfDebugger") << setw(15) << stlp1 << setw(15) << shlp1 << setw(15) << sqrt(state.localError().matrix()[0][0])
1073  // << setw(15) << (stlp1-shlp1)/stlp1 << setw(15) << (stlp1-shlp1)/sqrt(state.localError().matrix()[0][0]) ;
1074  // LogTrace("CkfDebugger") << setw(15) << stlp2 << setw(15) << shlp2 << setw(15) << sqrt(state.localError().matrix()[1][1])
1075  // << setw(15) << (stlp2-shlp2)/stlp2 << setw(15) << (stlp2-shlp2)/sqrt(state.localError().matrix()[1][1]) ;
1076  // LogTrace("CkfDebugger") << setw(15) << stlp3 << setw(15) << shlp3 << setw(15) << sqrt(state.localError().matrix()[2][2])
1077  // << setw(15) << (stlp3-shlp3)/stlp3 << setw(15) << (stlp3-shlp3)/sqrt(state.localError().matrix()[2][2]) ;
1078  // LogTrace("CkfDebugger") << setw(15) << stlp4 << setw(15) << shlp4 << setw(15) << sqrt(state.localError().matrix()[3][3])
1079  // << setw(15) << (stlp4-shlp4)/stlp4 << setw(15) << (stlp4-shlp4)/sqrt(state.localError().matrix()[3][3]) ;
1080  // LogTrace("CkfDebugger") << setw(15) << stlp5 << setw(15) << shlp5 << setw(15) << sqrt(state.localError().matrix()[4][4]) <<
1081  // setw(15) << (stlp5-shlp5)/stlp5 << setw(15) << (stlp5-shlp5)/sqrt(state.localError().matrix()[4][4]) ;
1082 
1084  R.Invert();
1085  double chi2 = ROOT::Math::Similarity(v - state.localParameters().vector(), R);
1086  LogTrace("CkfDebugger") << "chi2=" << chi2;
1087  return chi2;
1088  }
1089 
1090  return 0; //fixme
1091 }
1092 
1094  for (int it = 0; it != ((int)(dump.size())); it++)
1095  edm::LogVerbatim("CkfDebugger") << "dump " << it << " " << dump[it];
1096 
1097  edm::LogVerbatim("CkfDebugger");
1098  edm::LogVerbatim("CkfDebugger") << "seedWithDelta=" << ((double)seedWithDelta / totSeeds);
1099  edm::LogVerbatim("CkfDebugger") << "problems=" << ((double)problems / totSeeds);
1100  edm::LogVerbatim("CkfDebugger") << "no_sim_hit=" << ((double)no_sim_hit / totSeeds);
1101  edm::LogVerbatim("CkfDebugger") << "no_layer=" << ((double)no_layer / totSeeds);
1102  edm::LogVerbatim("CkfDebugger") << "layer_not_found=" << ((double)layer_not_found / totSeeds);
1103  edm::LogVerbatim("CkfDebugger") << "det_not_found=" << ((double)det_not_found / totSeeds);
1104  edm::LogVerbatim("CkfDebugger") << "chi2gt30=" << ((double)chi2gt30 / totSeeds);
1105  edm::LogVerbatim("CkfDebugger") << "chi2gt30deltaSeed=" << ((double)chi2gt30deltaSeed / totSeeds);
1106  edm::LogVerbatim("CkfDebugger") << "chi2gt30delta=" << ((double)chi2gt30delta / totSeeds);
1107  edm::LogVerbatim("CkfDebugger") << "chi2ls30=" << ((double)chi2ls30 / totSeeds);
1108  edm::LogVerbatim("CkfDebugger") << "simple_hit_not_found=" << ((double)simple_hit_not_found / totSeeds);
1109  edm::LogVerbatim("CkfDebugger") << "no_component=" << ((double)no_component / totSeeds);
1110  edm::LogVerbatim("CkfDebugger") << "only_one_component=" << ((double)only_one_component / totSeeds);
1111  edm::LogVerbatim("CkfDebugger") << "matched_not_found=" << ((double)matched_not_found / totSeeds);
1112  edm::LogVerbatim("CkfDebugger") << "matched_not_associated=" << ((double)matched_not_associated / totSeeds);
1113  edm::LogVerbatim("CkfDebugger") << "partner_det_not_fuond=" << ((double)partner_det_not_fuond / totSeeds);
1114  edm::LogVerbatim("CkfDebugger") << "glued_det_not_fuond=" << ((double)glued_det_not_fuond / totSeeds);
1115  edm::LogVerbatim("CkfDebugger") << "propagation=" << ((double)propagation / totSeeds);
1116  edm::LogVerbatim("CkfDebugger") << "other=" << ((double)other / totSeeds);
1117  edm::LogVerbatim("CkfDebugger") << "totchi2gt30=" << ((double)totchi2gt30 / totSeeds);
1118  edm::LogVerbatim("CkfDebugger") << "totSeeds=" << totSeeds;
1119  edm::LogVerbatim("CkfDebugger");
1120 
1121  edm::LogVerbatim("CkfDebugger") << "layer navigation problems:";
1122  for (int i = 0; i != 6; i++)
1123  for (int j = 0; j != 9; j++) {
1124  if (i == 0 && j > 2)
1125  break;
1126  if (i == 1 && j > 1)
1127  break;
1128  if (i == 2 && j > 3)
1129  break;
1130  if (i == 3 && j > 2)
1131  break;
1132  if (i == 4 && j > 5)
1133  break;
1134  if (i == 5 && j > 8)
1135  break;
1136  edm::LogVerbatim("CkfDebugger") << "det=" << i + 1 << " lay=" << j + 1 << " " << dump2[pair<int, int>(i, j)];
1137  }
1138  edm::LogVerbatim("CkfDebugger") << "\nlayer with hit having chi2>30:";
1139  for (int i = 0; i != 6; i++)
1140  for (int j = 0; j != 9; j++) {
1141  if (i == 0 && j > 2)
1142  break;
1143  if (i == 1 && j > 1)
1144  break;
1145  if (i == 2 && j > 3)
1146  break;
1147  if (i == 3 && j > 2)
1148  break;
1149  if (i == 4 && j > 5)
1150  break;
1151  if (i == 5 && j > 8)
1152  break;
1153  edm::LogVerbatim("CkfDebugger") << "det=" << i + 1 << " lay=" << j + 1 << " " << dump3[pair<int, int>(i, j)];
1154  }
1155  edm::LogVerbatim("CkfDebugger") << "\nlayer with hit having chi2>30 for delta rays:";
1156  for (int i = 0; i != 6; i++)
1157  for (int j = 0; j != 9; j++) {
1158  if (i == 0 && j > 2)
1159  break;
1160  if (i == 1 && j > 1)
1161  break;
1162  if (i == 2 && j > 3)
1163  break;
1164  if (i == 3 && j > 2)
1165  break;
1166  if (i == 4 && j > 5)
1167  break;
1168  if (i == 5 && j > 8)
1169  break;
1170  edm::LogVerbatim("CkfDebugger") << "det=" << i + 1 << " lay=" << j + 1 << " " << dump5[pair<int, int>(i, j)];
1171  }
1172  edm::LogVerbatim("CkfDebugger") << "\nlayer with det not found:";
1173  for (int i = 0; i != 6; i++)
1174  for (int j = 0; j != 9; j++) {
1175  if (i == 0 && j > 2)
1176  break;
1177  if (i == 1 && j > 1)
1178  break;
1179  if (i == 2 && j > 3)
1180  break;
1181  if (i == 3 && j > 2)
1182  break;
1183  if (i == 4 && j > 5)
1184  break;
1185  if (i == 5 && j > 8)
1186  break;
1187  edm::LogVerbatim("CkfDebugger") << "det=" << i + 1 << " lay=" << j + 1 << " " << dump4[pair<int, int>(i, j)];
1188  }
1189  edm::LogVerbatim("CkfDebugger") << "\nlayer with correct RecHit after missing Sim Hit:";
1190  for (int i = 0; i != 6; i++)
1191  for (int j = 0; j != 9; j++) {
1192  if (i == 0 && j > 2)
1193  break;
1194  if (i == 1 && j > 1)
1195  break;
1196  if (i == 2 && j > 3)
1197  break;
1198  if (i == 3 && j > 2)
1199  break;
1200  if (i == 4 && j > 5)
1201  break;
1202  if (i == 5 && j > 8)
1203  break;
1204  edm::LogVerbatim("CkfDebugger") << "det=" << i + 1 << " lay=" << j + 1 << " " << dump6[pair<int, int>(i, j)];
1205  }
1206  hchi2seedAll->Write();
1207  hchi2seedProb->Write();
1208  std::stringstream title;
1209  for (int i = 0; i != 6; i++)
1210  for (int j = 0; j != 9; j++) {
1211  if (i == 0 && j > 2)
1212  break;
1213  if (i == 1 && j > 1)
1214  break;
1215  if (i == 2 && j > 3)
1216  break;
1217  if (i == 3 && j > 2)
1218  break;
1219  if (i == 4 && j > 5)
1220  break;
1221  if (i == 5 && j > 8)
1222  break;
1223  title.str("");
1224  title << "pullX_" << i + 1 << "-" << j + 1 << "_sh-rh";
1225  hPullX_shrh[title.str()]->Write();
1226  title.str("");
1227  title << "pullY_" << i + 1 << "-" << j + 1 << "_sh-rh";
1228  hPullY_shrh[title.str()]->Write();
1229  title.str("");
1230  title << "pullX_" << i + 1 << "-" << j + 1 << "_sh-st";
1231  hPullX_shst[title.str()]->Write();
1232  title.str("");
1233  title << "pullY_" << i + 1 << "-" << j + 1 << "_sh-st";
1234  hPullY_shst[title.str()]->Write();
1235  title.str("");
1236  title << "pullX_" << i + 1 << "-" << j + 1 << "_st-rh";
1237  hPullX_strh[title.str()]->Write();
1238  title.str("");
1239  title << "pullY_" << i + 1 << "-" << j + 1 << "_st-rh";
1240  hPullY_strh[title.str()]->Write();
1241  title.str("");
1242  title << "PullGP_X_" << i + 1 << "-" << j + 1 << "_sh-st";
1243  hPullGP_X_shst[title.str()]->Write();
1244  title.str("");
1245  title << "PullGP_Y_" << i + 1 << "-" << j + 1 << "_sh-st";
1246  hPullGP_Y_shst[title.str()]->Write();
1247  title.str("");
1248  title << "PullGP_Z_" << i + 1 << "-" << j + 1 << "_sh-st";
1249  hPullGP_Z_shst[title.str()]->Write();
1250  if (((i == 2 || i == 4) && (j == 0 || j == 1)) || (i == 3 || i == 5)) {
1251  title.str("");
1252  title << "pullM_" << i + 1 << "-" << j + 1 << "_sh-rh";
1253  hPullM_shrh[title.str()]->Write();
1254  title.str("");
1255  title << "pullS_" << i + 1 << "-" << j + 1 << "_sh-rh";
1256  hPullS_shrh[title.str()]->Write();
1257  title.str("");
1258  title << "pullM_" << i + 1 << "-" << j + 1 << "_sh-st";
1259  hPullM_shst[title.str()]->Write();
1260  title.str("");
1261  title << "pullS_" << i + 1 << "-" << j + 1 << "_sh-st";
1262  hPullS_shst[title.str()]->Write();
1263  title.str("");
1264  title << "pullM_" << i + 1 << "-" << j + 1 << "_st-rh";
1265  hPullM_strh[title.str()]->Write();
1266  title.str("");
1267  title << "pullS_" << i + 1 << "-" << j + 1 << "_st-rh";
1268  hPullS_strh[title.str()]->Write();
1269  }
1270  }
1271  hPullGPXvsGPX_shst->Write();
1272  hPullGPXvsGPY_shst->Write();
1273  hPullGPXvsGPZ_shst->Write();
1274  hPullGPXvsGPr_shst->Write();
1275  hPullGPXvsGPeta_shst->Write();
1276  hPullGPXvsGPphi_shst->Write();
1277 
1278  //file->Write();
1279  file->Close();
1280 }
Log< level::Info, true > LogVerbatim
unsigned short processType() const
Definition: CkfDebugger.h:86
std::map< std::string, TH1F * > hPullM_shst
Definition: CkfDebugger.h:216
GlobalPoint globalPosition() const
Definition: CkfDebugger.h:75
tuple propagator
float xx() const
Definition: LocalError.h:22
const TrackerTopology * theTopo
Definition: CkfDebugger.h:103
bool hasDelta(const PSimHit *correctHit)
Definition: CkfDebugger.h:131
std::pair< CTTRHp, double > analyseRecHitExistance(const PSimHit &sh, const TSOS &startingState)
Definition: CkfDebugger.cc:714
bool correctMeas(const TM &tm, const PSimHit *correctHit) const
Definition: CkfDebugger.cc:668
void printSimHits(const edm::Event &iEvent)
Definition: CkfDebugger.cc:166
T perp() const
Definition: PV3DBase.h:69
ConstRecHitPointer const & recHit() const
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
std::map< std::string, TH1F * > hPullS_strh
Definition: CkfDebugger.h:219
TrackerHitAssociator::Config trackerHitAssociatorConfig_
Definition: CkfDebugger.h:99
std::map< std::string, TH1F * > hPullX_shrh
Definition: CkfDebugger.h:207
const LocalTrajectoryParameters & localParameters() const
TH2F * hPullGPXvsGPeta_shst
Definition: CkfDebugger.h:229
void dumpSimHit(const SimHit &hit) const
Definition: CkfDebugger.cc:190
std::pair< double, double > computePulls(CTTRHp recHit, TSOS startingState)
Definition: CkfDebugger.h:155
LocalVector momentumAtEntry() const
The momentum of the track that produced the hit, at entry point.
Definition: PSimHit.h:55
const MagneticField * theMagField
Definition: CkfDebugger.h:95
const TransientTrackingRecHitBuilder * theTTRHBuilder
Definition: CkfDebugger.h:102
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
std::map< std::pair< int, int >, int > dump6
Definition: CkfDebugger.h:202
const GeometricSearchTracker * theGeomSearchTracker
Definition: CkfDebugger.h:96
AlgebraicVector measuredParameters(const TrackingRecHit &)
const MeasurementTrackerEvent * theMeasurementTracker
Definition: CkfDebugger.h:101
virtual const GeomDet & geomDet() const
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
const CartesianTrajectoryError cartesianError() const
T y() const
Definition: PV3DBase.h:60
std::map< std::pair< int, int >, int > dump3
Definition: CkfDebugger.h:199
std::vector< ConstRecHitPointer > RecHitContainer
GlobalPoint globalPosition() const
ConstRecHitContainer recHits() const
Definition: Trajectory.h:186
const Surface & surface() const
int partner_det_not_fuond
Definition: CkfDebugger.h:247
std::vector< const PSimHit * > nextCorrectHits(const Trajectory &, unsigned int &)
Definition: CkfDebugger.cc:562
int det_not_found
Definition: CkfDebugger.h:237
unsigned int partnerDetId() const
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
std::map< std::pair< int, int >, int > dump2
Definition: CkfDebugger.h:198
int seedWithDelta
Definition: CkfDebugger.h:232
LocalError positionError() const
virtual std::vector< DetWithState > compatibleDets(const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
const MeasurementEstimator * theChi2
Definition: CkfDebugger.h:97
int matched_not_found
Definition: CkfDebugger.h:245
#define LogTrace(id)
unsigned int trackId() const
Definition: CkfDebugger.h:77
tuple result
Definition: mps_fire.py:311
PropagationDirection const & direction() const
Definition: Trajectory.cc:133
DataContainer const & measurements() const
Definition: Trajectory.h:178
AlgebraicVector5 vector() const
int analyseRecHitNotFound(const Trajectory &, CTTRHp)
Definition: CkfDebugger.cc:936
int iEvent
Definition: GenABIO.cc:224
const SurfaceType & surface() const
unsigned int glued() const
glued
T mag() const
Definition: PV3DBase.h:64
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:43
float timeOfFlight() const
Definition: PSimHit.h:73
double testSeed(CTTRHp, CTTRHp, TrajectoryStateOnSurface)
std::map< std::pair< int, int >, int > dump4
Definition: CkfDebugger.h:200
int particleType() const
Definition: CkfDebugger.h:84
int matched_not_associated
Definition: CkfDebugger.h:246
std::pair< const GeomDet *, TrajectoryStateOnSurface > DetWithState
float yy() const
Definition: LocalError.h:24
std::map< std::string, TH1F * > hPullS_shrh
Definition: CkfDebugger.h:215
Local3DPoint localPosition() const
Definition: PSimHit.h:52
bool associated(CTTRHp rechit, const PSimHit &sh) const
Definition: CkfDebugger.cc:647
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
int only_one_component
Definition: CkfDebugger.h:244
T sqrt(T t)
Definition: SSEVec.h:19
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
AlgebraicSymMatrix measuredError(const TrackingRecHit &)
TrajectoryMeasurement const & lastMeasurement() const
Definition: Trajectory.h:150
T z() const
Definition: PV3DBase.h:61
int assocTrackId(CTTRHp rechit) const
Definition: CkfDebugger.cc:547
def move
Definition: eostools.py:511
TH2F * hPullGPXvsGPr_shst
Definition: CkfDebugger.h:228
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
int glued_det_not_fuond
Definition: CkfDebugger.h:248
TH1F * hchi2seedAll
Definition: CkfDebugger.h:205
ROOT::Math::SVector< double, 5 > AlgebraicVector5
const DetLayer * detLayer(const DetId &id) const
obsolete method. Use idToLayer() instead.
TrackerHitAssociator * hitAssociator
Definition: CkfDebugger.h:100
TH2F * hPullGPXvsGPY_shst
Definition: CkfDebugger.h:226
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
CkfDebugger(edm::EventSetup const &es, edm::ConsumesCollector &&iC)
Definition: CkfDebugger.cc:41
TransientTrackingRecHit::ConstRecHitPointer CTTRHp
Definition: CkfDebugger.h:39
const AlgebraicSymMatrix55 & matrix() const
const LocalTrajectoryError & localError() const
std::map< std::string, TH1F * > hPullY_shrh
Definition: CkfDebugger.h:208
virtual HitReturnType estimate(const TrajectoryStateOnSurface &ts, const TrackingRecHit &hit) const =0
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:67
const Surface::PositionType & position() const
virtual RecHitContainer recHits(const TrajectoryStateOnSurface &, const MeasurementTrackerEvent &) const =0
DetId gluedId(const DetId &du)
Definition: CkfDebugger.cc:36
const Propagator * theForwardPropagator
Definition: CkfDebugger.h:98
bool correctTrajectory(const Trajectory &, unsigned int &) const
Definition: CkfDebugger.cc:501
int chi2gt30deltaSeed
Definition: CkfDebugger.h:240
Definition: DetId.h:17
CLHEP::HepVector AlgebraicVector
TrajectoryMeasurement const & firstMeasurement() const
Definition: Trajectory.h:166
ConstRecHitContainer RecHitContainer
Definition: Trajectory.h:42
float pabs() const
Definition: CkfDebugger.h:81
std::map< std::string, TH1F * > hPullX_shst
Definition: CkfDebugger.h:209
int simple_hit_not_found
Definition: CkfDebugger.h:242
const GlobalError position() const
Position error submatrix.
std::vector< const DetLayer * > nextLayers(const DetLayer &detLayer, Args &&...args) const
std::map< std::string, TH1F * > hPullGP_Y_shst
Definition: CkfDebugger.h:222
T const * product() const
Definition: ESHandle.h:86
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
virtual const BoundCylinder & specificSurface() const final
Extension of the interface.
double b
Definition: hdecay.h:118
TrajectoryStateOnSurface propagate(STA const &state, SUR const &surface) const
Definition: Propagator.h:50
unsigned short processType() const
Definition: PSimHit.h:120
TFile * file
Definition: CkfDebugger.h:204
T eta() const
Definition: PV3DBase.h:73
const TrackerGeometry * theTrackerGeom
Definition: CkfDebugger.h:94
std::map< std::string, TH1F * > hPullM_shrh
Definition: CkfDebugger.h:214
TrackingRecHit::ConstRecHitContainer RecHitContainer
const GeomDetUnit * det(const PSimHit *sh) const
Definition: CkfDebugger.h:147
edm::EventID id() const
Definition: EventBase.h:59
std::vector< SimHitIdpr > associateHitId(const TrackingRecHit &thit) const
NavigationSchool const * theNavSchool
Definition: CkfDebugger.h:104
bool analyseCompatibleMeasurements(const Trajectory &, const std::vector< TrajectoryMeasurement > &, const MeasurementTrackerEvent *, const Propagator *, const Chi2MeasurementEstimatorBase *, const TransientTrackingRecHitBuilder *)
Definition: CkfDebugger.cc:197
bool goodSimHit(const PSimHit &sh) const
Definition: CkfDebugger.cc:640
std::map< std::pair< int, int >, int > dump5
Definition: CkfDebugger.h:201
double a
Definition: hdecay.h:119
T get() const
Definition: EventSetup.h:82
TH2F * hPullGPXvsGPX_shst
Definition: CkfDebugger.h:225
std::map< std::string, TH1F * > hPullGP_Z_shst
Definition: CkfDebugger.h:223
CLHEP::HepSymMatrix AlgebraicSymMatrix
std::map< std::string, TH1F * > hPullX_strh
Definition: CkfDebugger.h:211
unsigned int trackId() const
Definition: PSimHit.h:106
const DetLayer * lastLayer() const
Redundant method, returns the layer of lastMeasurement() .
Definition: Trajectory.h:288
std::vector< PSimHit > associateHit(const TrackingRecHit &thit) const
int layer_not_found
Definition: CkfDebugger.h:236
tuple MeasurementTrackerEvent
int no_component
Definition: CkfDebugger.h:243
TrajectoryStateOnSurface const & updatedState() const
TH2F * hPullGPXvsGPZ_shst
Definition: CkfDebugger.h:227
std::vector< int > dump
Definition: CkfDebugger.h:197
std::map< std::string, TH1F * > hPullS_shst
Definition: CkfDebugger.h:217
const PSimHit * pSimHit(unsigned int tkId, DetId detId)
Definition: CkfDebugger.cc:925
int layer(const GeomDet *det)
Definition: CkfDebugger.h:149
std::map< std::string, TH1F * > hPullY_shst
Definition: CkfDebugger.h:210
std::pair< const GeomDet *, TrajectoryStateOnSurface > DetWithState
Global3DPoint position(const PSimHit *sh) const
Definition: CkfDebugger.h:144
std::map< std::string, TH1F * > hPullM_strh
Definition: CkfDebugger.h:218
std::map< unsigned int, std::vector< PSimHit * > > idHitsMap
Definition: CkfDebugger.h:106
T x() const
Definition: PV3DBase.h:59
TH1F * hchi2seedProb
Definition: CkfDebugger.h:205
std::map< std::string, TH1F * > hPullY_strh
Definition: CkfDebugger.h:212
unsigned int detUnitId() const
Definition: PSimHit.h:97
TH2F * hPullGPXvsGPphi_shst
Definition: CkfDebugger.h:230
int chi2gt30delta
Definition: CkfDebugger.h:239
std::map< std::string, TH1F * > hPullGP_X_shst
Definition: CkfDebugger.h:221