CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MultiHitGeneratorFromChi2.cc
Go to the documentation of this file.
2 
7 
12 
15 
22 
28 
30 
32 
33 #include <algorithm>
34 #include <iostream>
35 #include <vector>
36 #include <cmath>
37 #include <map>
38 #include <limits>
39 
40 // Having this macro reduces the need to pollute the code with
41 // #ifdefs. The idea is that the condition is checked only if
42 // debugging is enabled. That way the condition expression may use
43 // variables that are declared only if EDM_ML_DEBUG is enabled. If it
44 // is disabled, rely on the fact that LogTrace should compile to
45 // no-op.
46 #ifdef EDM_ML_DEBUG
47 #define IfLogTrace(cond, cat) if(cond) LogTrace(cat)
48 #else
49 #define IfLogTrace(cond, cat) LogTrace(cat)
50 #endif
51 
52 using namespace std;
53 
55 
56 namespace {
57  struct LayerRZPredictions {
59  };
60 }
61 
64  useFixedPreFiltering(cfg.getParameter<bool>("useFixedPreFiltering")),
65  extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),//extra window in ThirdHitRZPrediction range
66  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),//extra window in ThirdHitPredictionFromCircle range (divide by R to get phi)
67  extraZKDBox(cfg.getParameter<double>("extraZKDBox")),//extra windown in Z when building the KDTree box (used in barrel)
68  extraRKDBox(cfg.getParameter<double>("extraRKDBox")),//extra windown in R when building the KDTree box (used in endcap)
69  extraPhiKDBox(cfg.getParameter<double>("extraPhiKDBox")),//extra windown in Phi when building the KDTree box
70  fnSigmaRZ(cfg.getParameter<double>("fnSigmaRZ")),//this multiplies the max hit error on the layer when building the KDTree box
71  chi2VsPtCut(cfg.getParameter<bool>("chi2VsPtCut")),
72  maxChi2(cfg.getParameter<double>("maxChi2")),
73  refitHits(cfg.getParameter<bool>("refitHits")),
74  filterName_(cfg.getParameter<std::string>("ClusterShapeHitFilterName")),
75  builderName_(cfg.existsAs<std::string>("TTRHBuilder") ? cfg.getParameter<std::string>("TTRHBuilder") : std::string("WithTrackAngle")),
76  useSimpleMF_(false),
77  mfName_("")
78 {
80  dphi = cfg.getParameter<double>("phiPreFiltering");
81  if (chi2VsPtCut) {
82  pt_interv = cfg.getParameter<std::vector<double> >("pt_interv");
83  chi2_cuts = cfg.getParameter<std::vector<double> >("chi2_cuts");
84  }
85 #ifdef EDM_ML_DEBUG
86  detIdsToDebug = cfg.getParameter<std::vector<int> >("detIdsToDebug");
87  //if (detIdsToDebug.size()<3) //fixme
88 #else
89  detIdsToDebug.push_back(0);
90  detIdsToDebug.push_back(0);
91  detIdsToDebug.push_back(0);
92 #endif
93  // 2014/02/11 mia:
94  // we should get rid of the boolean parameter useSimpleMF,
95  // and use only a string magneticField [instead of SimpleMagneticField]
96  // or better an edm::ESInputTag (at the moment HLT does not handle ESInputTag)
97  if (cfg.exists("SimpleMagneticField")) {
98  useSimpleMF_ = true;
99  mfName_ = cfg.getParameter<std::string>("SimpleMagneticField");
100  }
101  filter = 0;
102  bfield = 0;
103  nomField = -1.;
104 }
105 
108 {
109 
111  if (useSimpleMF_) es.get<IdealMagneticFieldRecord>().get(mfName_, bfield_h);
112  else es.get<IdealMagneticFieldRecord>().get(bfield_h);
113  bfield = bfield_h.product();
115 
117  es.get<CkfComponentsRecord>().get(filterName_, filterHandle_);
118  filter = filterHandle_.product();
119 
121  es.get<TransientRecHitRecord>().get(builderName_, builderH);
122  builder = (TkTransientTrackingRecHitBuilder const *)(builderH.product());
123  cloner = (*builder).cloner();
124 }
125 
126 
127 namespace {
128  inline
129  bool intersect(Range &range, const Range &second)
130  {
131  if (range.first > second.max() || range.second < second.min())
132  return false;
133  if (range.first < second.min())
134  range.first = second.min();
135  if (range.second > second.max())
136  range.second = second.max();
137  return range.first < range.second;
138  }
139 }
140 
143  const edm::Event & ev,
144  const edm::EventSetup& es,
146  std::vector<SeedingLayerSetsHits::SeedingLayer> thirdLayers)
147 {
148  unsigned int debug_Id0 = detIdsToDebug[0];
149  unsigned int debug_Id1 = detIdsToDebug[1];
150  unsigned int debug_Id2 = detIdsToDebug[2];
151 
152  LogDebug("MultiHitGeneratorFromChi2") << "pair: " << thePairGenerator->innerLayer(pairLayers).name() << "+" << thePairGenerator->outerLayer(pairLayers).name() << " 3rd lay size: " << thirdLayers.size();
153 
154  //gc: first get the pairs
155  OrderedHitPairs pairs;
156  pairs.reserve(30000);
157  thePairGenerator->hitPairs(region,pairs,ev,es, pairLayers);
158  LogTrace("MultiHitGeneratorFromChi2") << "";
159  if (pairs.empty()) {
160  // LogDebug("MultiHitGeneratorFromChi2") << "empy pairs";
161  return;
162  }
163 
164  //gc: these are all the layers compatible with the layer pairs (as defined in the config file)
165  int size = thirdLayers.size();
166 
167 
168  //gc: initialize a KDTree per each 3rd layer
169  std::vector<KDTreeNodeInfo<RecHitsSortedInPhi::HitIter> > layerTree; // re-used throughout
170  std::vector<RecHitsSortedInPhi::HitIter> foundNodes; // re-used thoughout
171  foundNodes.reserve(100);
172  #ifdef __clang__
173  std::vector<KDTreeLinkerAlgo<RecHitsSortedInPhi::HitIter>> hitTree(size);
174  #else
176  #endif
177  float rzError[size]; //save maximum errors
178  double maxphi = Geom::twoPi(), minphi = -maxphi; //increase to cater for any range
179 
180  map<std::string, LayerRZPredictions> mapPred;//need to use the name as map key since we may have more than one SeedingLayer per DetLayer (e.g. TID and MTID)
181  const RecHitsSortedInPhi * thirdHitMap[size];//gc: this comes from theLayerCache
182 
183  //gc: loop over each layer
184  for(int il = 0; il < size; il++) {
185  thirdHitMap[il] = &(*theLayerCache)(thirdLayers[il], region, ev, es);
186  LogTrace("MultiHitGeneratorFromChi2") << "considering third layer: " << thirdLayers[il].name() << " with hits: " << thirdHitMap[il]->all().second-thirdHitMap[il]->all().first;
187  const DetLayer *layer = thirdLayers[il].detLayer();
188  LayerRZPredictions &predRZ = mapPred[thirdLayers[il].name()];
189  predRZ.line.initLayer(layer);
190  predRZ.line.initTolerance(extraHitRZtolerance);
191 
192  //gc: now we take all hits in the layer and fill the KDTree
193  RecHitsSortedInPhi::Range hitRange = thirdHitMap[il]->all(); // Get iterators
194  layerTree.clear();
195  double minz=999999.0, maxz= -999999.0; // Initialise to extreme values in case no hits
196  float maxErr=0.0f;
197  bool barrelLayer = (thirdLayers[il].detLayer()->location() == GeomDetEnumerators::barrel);
198  if (hitRange.first != hitRange.second)
199  { minz = barrelLayer? hitRange.first->hit()->globalPosition().z() : hitRange.first->hit()->globalPosition().perp();
200  maxz = minz; //In case there's only one hit on the layer
201  for (RecHitsSortedInPhi::HitIter hi=hitRange.first; hi != hitRange.second; ++hi)
202  { double angle = hi->phi();
203  double myz = barrelLayer? hi->hit()->globalPosition().z() : hi->hit()->globalPosition().perp();
204 
205  IfLogTrace(hi->hit()->rawId()==debug_Id2, "MultiHitGeneratorFromChi2") << "filling KDTree with hit in id=" << debug_Id2
206  << " with pos: " << hi->hit()->globalPosition()
207  << " phi=" << hi->hit()->globalPosition().phi()
208  << " z=" << hi->hit()->globalPosition().z()
209  << " r=" << hi->hit()->globalPosition().perp();
210  //use (phi,r) for endcaps rather than (phi,z)
211  if (myz < minz) { minz = myz;} else { if (myz > maxz) {maxz = myz;}}
212  float myerr = barrelLayer? hi->hit()->errorGlobalZ(): hi->hit()->errorGlobalR();
213  if (myerr > maxErr) { maxErr = myerr;}
214  layerTree.push_back(KDTreeNodeInfo<RecHitsSortedInPhi::HitIter>(hi, angle, myz)); // save it
215  if (angle < 0) // wrap all points in phi
216  { layerTree.push_back(KDTreeNodeInfo<RecHitsSortedInPhi::HitIter>(hi, angle+Geom::twoPi(), myz));}
217  else
218  { layerTree.push_back(KDTreeNodeInfo<RecHitsSortedInPhi::HitIter>(hi, angle-Geom::twoPi(), myz));}
219  }
220  }
221  KDTreeBox phiZ(minphi, maxphi, minz-0.01, maxz+0.01); // declare our bounds
222  //add fudge factors in case only one hit and also for floating-point inaccuracy
223  hitTree[il].build(layerTree, phiZ); // make KDtree
224  rzError[il] = maxErr; //save error
225  }
226  //gc: now we have initialized the KDTrees and we are out of the layer loop
227 
228  //gc: this sets the minPt of the triplet
229  double curv = PixelRecoUtilities::curvature(1. / region.ptMin(), es);
230 
231  LogTrace("MultiHitGeneratorFromChi2") << "pair size=" << pairs.size() << std::endl;
232 
233  //gc: now we loop over all pairs
234  for (OrderedHitPairs::const_iterator ip = pairs.begin(); ip != pairs.end(); ++ip) {
235 
236  int foundTripletsFromPair = 0;
237  bool usePair = false;
238  cacheHitPointer bestH2;
240 
241  SeedingHitSet::ConstRecHitPointer oriHit0 = ip->inner();
242  SeedingHitSet::ConstRecHitPointer oriHit1 = ip->outer();
243 
244  HitOwnPtr hit0(*oriHit0);
245  HitOwnPtr hit1(*oriHit1);
246  GlobalPoint gp0 = hit0->globalPosition();
247  GlobalPoint gp1 = hit1->globalPosition();
248 
249 #ifdef EDM_ML_DEBUG
250  bool debugPair = ip->inner()->rawId()==debug_Id0 && ip->outer()->rawId()==debug_Id1;
251 #endif
252  IfLogTrace(debugPair, "MultiHitGeneratorFromChi2") << endl << endl
253  << "found new pair with ids "<<debug_Id0<<" "<<debug_Id1<<" with pos: " << gp0 << " " << gp1;
254 
255  if (refitHits) {
256 
257  TrajectoryStateOnSurface tsos0, tsos1;
258  assert(!hit0.isOwn()); assert(!hit1.isOwn());
259 #ifdef EDM_ML_DEBUG
260  refit2Hits(hit0,hit1,tsos0,tsos1,region,nomField,debugPair);
261 #else
262  refit2Hits(hit0,hit1,tsos0,tsos1,region,nomField,false);
263 #endif
264  assert(hit0.isOwn()); assert(hit1.isOwn());
265 
266  //fixme add pixels
267  bool passFilterHit0 = true;
268  if (//hit0->geographicalId().subdetId() > 2
269  hit0->geographicalId().subdetId()==SiStripDetId::TIB
270  || hit0->geographicalId().subdetId()==SiStripDetId::TID
271  //|| hit0->geographicalId().subdetId()==SiStripDetId::TOB
272  //|| hit0->geographicalId().subdetId()==SiStripDetId::TEC
273  ) {
274  const std::type_info &tid = typeid(*hit0->hit());
275  if (tid == typeid(SiStripMatchedRecHit2D)) {
276  const SiStripMatchedRecHit2D* matchedHit = dynamic_cast<const SiStripMatchedRecHit2D *>(hit0->hit());
277  if (filter->isCompatible(DetId(matchedHit->monoId()), matchedHit->monoCluster(), tsos0.localMomentum())==0 ||
278  filter->isCompatible(DetId(matchedHit->stereoId()), matchedHit->stereoCluster(), tsos0.localMomentum())==0) passFilterHit0 = false;
279  } else if (tid == typeid(SiStripRecHit2D)) {
280  const SiStripRecHit2D* recHit = dynamic_cast<const SiStripRecHit2D *>(hit0->hit());
281  if (filter->isCompatible(*recHit, tsos0.localMomentum())==0) passFilterHit0 = false;
282  } else if (tid == typeid(ProjectedSiStripRecHit2D)) {
283  const ProjectedSiStripRecHit2D* precHit = dynamic_cast<const ProjectedSiStripRecHit2D *>(hit0->hit());
284  if (filter->isCompatible(precHit->originalHit(), tsos0.localMomentum())==0) passFilterHit0 = false; //FIXME
285  }
286  }
287  IfLogTrace(debugPair&&!passFilterHit0, "MultiHitGeneratorFromChi2") << "hit0 did not pass cluster shape filter";
288  if (!passFilterHit0) continue;
289  bool passFilterHit1 = true;
290  if (//hit1->geographicalId().subdetId() > 2
291  hit1->geographicalId().subdetId()==SiStripDetId::TIB
292  || hit1->geographicalId().subdetId()==SiStripDetId::TID
293  //|| hit1->geographicalId().subdetId()==SiStripDetId::TOB
294  //|| hit1->geographicalId().subdetId()==SiStripDetId::TEC
295  ) {
296  const std::type_info &tid = typeid(*hit1->hit());
297  if (tid == typeid(SiStripMatchedRecHit2D)) {
298  const SiStripMatchedRecHit2D* matchedHit = dynamic_cast<const SiStripMatchedRecHit2D *>(hit1->hit());
299  if (filter->isCompatible(DetId(matchedHit->monoId()), matchedHit->monoCluster(), tsos1.localMomentum())==0 ||
300  filter->isCompatible(DetId(matchedHit->stereoId()), matchedHit->stereoCluster(), tsos1.localMomentum())==0) passFilterHit1 = false;
301  } else if (tid == typeid(SiStripRecHit2D)) {
302  const SiStripRecHit2D* recHit = dynamic_cast<const SiStripRecHit2D *>(hit1->hit());
303  if (filter->isCompatible(*recHit, tsos1.localMomentum())==0) passFilterHit1 = false;
304  } else if (tid == typeid(ProjectedSiStripRecHit2D)) {
305  const ProjectedSiStripRecHit2D* precHit = dynamic_cast<const ProjectedSiStripRecHit2D *>(hit1->hit());
306  if (filter->isCompatible(precHit->originalHit(), tsos1.localMomentum())==0) passFilterHit1 = false; //FIXME
307  }
308  }
309  IfLogTrace(debugPair&&!passFilterHit1, "MultiHitGeneratorFromChi2") << "hit1 did not pass cluster shape filter";
310  if (!passFilterHit1) continue;
311 
312  } else {
313  // not refit clone anyhow
314  hit0.reset((BaseTrackerRecHit *)hit0->clone());
315  hit1.reset((BaseTrackerRecHit *)hit1->clone());
316  }
317 
318  //gc: create the RZ line for the pair
319  SimpleLineRZ line(PixelRecoPointRZ(gp0.perp(),gp0.z()), PixelRecoPointRZ(gp1.perp(),gp1.z()));
320  ThirdHitPredictionFromCircle predictionRPhi(gp0, gp1, extraHitRPhitolerance);
321 
322  //gc: this is the curvature of the two hits assuming the region
323  Range pairCurvature = predictionRPhi.curvature(region.originRBound());
324  //gc: intersect not only returns a bool but may change pairCurvature to intersection with curv
325  if (!intersect(pairCurvature, Range(-curv, curv))) {
326  IfLogTrace(debugPair, "MultiHitGeneratorFromChi2") << "curvature cut: curv=" << curv
327  << " gc=(" << pairCurvature.first << ", " << pairCurvature.second << ")";
328  continue;
329  }
330 
331  //gc: loop over all third layers compatible with the pair
332  for(int il = 0; (il < size) & (!usePair); il++) {
333 
334  IfLogTrace(debugPair, "MultiHitGeneratorFromChi2") << "cosider layer: " << thirdLayers[il].name() << " for this pair. Location: " << thirdLayers[il].detLayer()->location();
335 
336  if (hitTree[il].empty()) {
337  IfLogTrace(debugPair, "MultiHitGeneratorFromChi2") << "empty hitTree";
338  continue; // Don't bother if no hits
339  }
340 
341  cacheHitPointer bestL2;
342  float chi2FromThisLayer = std::numeric_limits<float>::max();
343 
344  const DetLayer *layer = thirdLayers[il].detLayer();
345  bool barrelLayer = layer->location() == GeomDetEnumerators::barrel;
346 
347  LayerRZPredictions &predRZ = mapPred.find(thirdLayers[il].name())->second;
348  predRZ.line.initPropagator(&line);
349 
350  //gc: this takes the z at R-thick/2 and R+thick/2 according to
351  // the line from the two points and the adds the extra tolerance
352  Range rzRange = predRZ.line();
353 
354  if (rzRange.first >= rzRange.second) {
355  IfLogTrace(debugPair, "MultiHitGeneratorFromChi2") << "rzRange empty";
356  continue;
357  }
358  //gc: check that rzRange is compatible with detector bounds
359  // note that intersect may change rzRange to intersection with bounds
360  if (!intersect(rzRange, predRZ.line.detSize())) {// theDetSize = Range(-maxZ, maxZ);
361  IfLogTrace(debugPair, "MultiHitGeneratorFromChi2") << "rzRange and detector do not intersect";
362  continue;
363  }
364  Range radius = barrelLayer ? predRZ.line.detRange() : rzRange;
365 
366  //gc: define the phi range of the hits
367  Range phiRange;
368  if (useFixedPreFiltering) {
369  //gc: in this case it takes as range the phi of the outer
370  // hit +/- the phiPreFiltering value from cfg
371  float phi0 = ip->outer()->globalPosition().phi();
372  phiRange = Range(phi0 - dphi, phi0 + dphi);
373  } else {
374  //gc: predictionRPhi uses the cosine rule to find the phi of the 3rd point at radius, assuming the pairCurvature range [-c,+c]
375  if (pairCurvature.first<0. && pairCurvature.second<0.) {
376  float phi12 = predictionRPhi.phi(pairCurvature.first,radius.second);
377  float phi21 = predictionRPhi.phi(pairCurvature.second,radius.first);
378  while(unlikely(phi12 < phi21)) phi12 += float(2. * M_PI);
379  phiRange = Range(phi21,phi12);
380  } else if (pairCurvature.first>=0. && pairCurvature.second>=0.) {
381  float phi11 = predictionRPhi.phi(pairCurvature.first,radius.first);
382  float phi22 = predictionRPhi.phi(pairCurvature.second,radius.second);
383  while(unlikely(phi11 < phi22)) phi11 += float(2. * M_PI);
384  phiRange = Range(phi22,phi11);
385  } else {
386  float phi12 = predictionRPhi.phi(pairCurvature.first,radius.second);
387  float phi22 = predictionRPhi.phi(pairCurvature.second,radius.second);
388  while(unlikely(phi12 < phi22)) phi12 += float(2. * M_PI);
389  phiRange = Range(phi22,phi12);
390  }
391  }
392 
393  //gc: this is the place where hits in the compatible region are put in the foundNodes
395  foundNodes.clear(); // Now recover hits in bounding box...
396  // This needs range -twoPi to +twoPi to work
397  float prmin=phiRange.min(), prmax=phiRange.max(); //get contiguous range
398  if ((prmax-prmin) > Geom::twoPi()) {
399  prmax=Geom::pi(); prmin = -Geom::pi();
400  } else {
401  while (prmax>maxphi) { prmin -= Geom::twoPi(); prmax -= Geom::twoPi();}
402  while (prmin<minphi) { prmin += Geom::twoPi(); prmax += Geom::twoPi();}
403  }
404 
405  IfLogTrace(debugPair, "MultiHitGeneratorFromChi2") << "defining kd tree box";
406 
407  if (barrelLayer) {
408  KDTreeBox phiZ(prmin-extraPhiKDBox, prmax+extraPhiKDBox,
409  rzRange.min()-fnSigmaRZ*rzError[il]-extraZKDBox,
410  rzRange.max()+fnSigmaRZ*rzError[il]+extraZKDBox);
411  hitTree[il].search(phiZ, foundNodes);
412 
413  IfLogTrace(debugPair, "MultiHitGeneratorFromChi2") << "kd tree box bounds, phi: " << prmin-extraPhiKDBox <<","<< prmax+extraPhiKDBox
414  << " z: "<< rzRange.min()-fnSigmaRZ*rzError[il]-extraZKDBox <<","<<rzRange.max()+fnSigmaRZ*rzError[il]+extraZKDBox
415  << " rzRange: " << rzRange.min() <<","<<rzRange.max();
416 
417  } else {
418  KDTreeBox phiR(prmin-extraPhiKDBox, prmax+extraPhiKDBox,
419  rzRange.min()-fnSigmaRZ*rzError[il]-extraRKDBox,
420  rzRange.max()+fnSigmaRZ*rzError[il]+extraRKDBox);
421  hitTree[il].search(phiR, foundNodes);
422 
423  IfLogTrace(debugPair, "MultiHitGeneratorFromChi2") << "kd tree box bounds, phi: " << prmin-extraPhiKDBox <<","<< prmax+extraPhiKDBox
424  << " r: "<< rzRange.min()-fnSigmaRZ*rzError[il]-extraRKDBox <<","<<rzRange.max()+fnSigmaRZ*rzError[il]+extraRKDBox
425  << " rzRange: " << rzRange.min() <<","<<rzRange.max();
426  }
427 
428  IfLogTrace(debugPair, "MultiHitGeneratorFromChi2") << "kd tree box size: " << foundNodes.size();
429 
430 
431  //gc: now we loop over the hits in the box for this layer
432  for (std::vector<RecHitsSortedInPhi::HitIter>::iterator ih = foundNodes.begin();
433  ih !=foundNodes.end() && !usePair; ++ih) {
434 
435  IfLogTrace(debugPair, "MultiHitGeneratorFromChi2") << endl << "triplet candidate";
436 
437  const RecHitsSortedInPhi::HitIter KDdata = *ih;
438 
439  SeedingHitSet::ConstRecHitPointer oriHit2 = KDdata->hit();
440  cacheHitPointer hit2;
441 
442  if (refitHits) {//fixme
443 
444  //fitting all 3 hits takes too much time... do it quickly only for 3rd hit
445  GlobalVector initMomentum(oriHit2->globalPosition() - gp1);
446  initMomentum *= (1./initMomentum.perp()); //set pT=1
447  GlobalTrajectoryParameters kine = GlobalTrajectoryParameters(oriHit2->globalPosition(), initMomentum, 1, &*bfield);
448  TrajectoryStateOnSurface state(kine,*oriHit2->surface());
449  hit2.reset((SeedingHitSet::RecHitPointer)(cloner(*oriHit2,state)));
450 
451  //fixme add pixels
452  bool passFilterHit2 = true;
453  if (hit2->geographicalId().subdetId()==SiStripDetId::TIB
454  || hit2->geographicalId().subdetId()==SiStripDetId::TID
455  // || hit2->geographicalId().subdetId()==SiStripDetId::TOB
456  // || hit2->geographicalId().subdetId()==SiStripDetId::TEC
457  ) {
458  const std::type_info &tid = typeid(*hit2->hit());
459  if (tid == typeid(SiStripMatchedRecHit2D)) {
460  const SiStripMatchedRecHit2D* matchedHit = dynamic_cast<const SiStripMatchedRecHit2D *>(hit2->hit());
461  if (filter->isCompatible(DetId(matchedHit->monoId()), matchedHit->monoCluster(), initMomentum)==0 ||
462  filter->isCompatible(DetId(matchedHit->stereoId()), matchedHit->stereoCluster(), initMomentum)==0) passFilterHit2 = false;
463  } else if (tid == typeid(SiStripRecHit2D)) {
464  const SiStripRecHit2D* recHit = dynamic_cast<const SiStripRecHit2D *>(hit2->hit());
465  if (filter->isCompatible(*recHit, initMomentum)==0) passFilterHit2 = false;
466  } else if (tid == typeid(ProjectedSiStripRecHit2D)) {
467  const ProjectedSiStripRecHit2D* precHit = dynamic_cast<const ProjectedSiStripRecHit2D *>(hit2->hit());
468  if (filter->isCompatible(precHit->originalHit(), initMomentum)==0) passFilterHit2 = false;
469  }
470  }
471  IfLogTrace(debugPair&&!passFilterHit2, "MultiHitGeneratorFromChi2") << "hit2 did not pass cluster shape filter";
472  if (!passFilterHit2) continue;
473 
474  // fitting all 3 hits takes too much time :-(
475  // TrajectoryStateOnSurface tsos0, tsos1, tsos2;
476  // refit3Hits(hit0,hit1,hit2,tsos0,tsos1,tsos2,nomField,debugPair);
477  // if (hit2->geographicalId().subdetId()==SiStripDetId::TIB
478  // // || hit2->geographicalId().subdetId()==SiStripDetId::TOB
479  // // || hit2->geographicalId().subdetId()==SiStripDetId::TID
480  // // || hit2->geographicalId().subdetId()==SiStripDetId::TEC
481  // ) {
482  // const std::type_info &tid = typeid(*hit2->hit());
483  // if (tid == typeid(SiStripMatchedRecHit2D)) {
484  // const SiStripMatchedRecHit2D* matchedHit = dynamic_cast<const SiStripMatchedRecHit2D *>(hit2->hit());
485  // if (filter->isCompatible(DetId(matchedHit->monoId()), matchedHit->monoCluster(), tsos2.localMomentum())==0 ||
486  // filter->isCompatible(DetId(matchedHit->stereoId()), matchedHit->stereoCluster(), tsos2.localMomentum())==0) continue;
487  // } else if (tid == typeid(SiStripRecHit2D)) {
488  // const SiStripRecHit2D* recHit = dynamic_cast<const SiStripRecHit2D *>(hit2->hit());
489  // if (filter->isCompatible(*recHit, tsos2.localMomentum())==0) continue;
490  // } else if (tid == typeid(ProjectedSiStripRecHit2D)) {
491  // const ProjectedSiStripRecHit2D* precHit = dynamic_cast<const ProjectedSiStripRecHit2D *>(hit2->hit());
492  // if (filter->isCompatible(precHit->originalHit(), tsos2.localMomentum())==0) continue;;
493  // }
494  // }
495 
496  } else {
497  // not refit clone anyhow
498  hit2.reset((BaseTrackerRecHit*)oriHit2->clone());
499  }
500 
501  //gc: add the chi2 cut
502  vector<GlobalPoint> gp(3);
503  vector<GlobalError> ge(3);
504  vector<bool> bl(3);
505  gp[0] = hit0->globalPosition();
506  ge[0] = hit0->globalPositionError();
507  int subid0 = hit0->geographicalId().subdetId();
508  bl[0] = (subid0 == StripSubdetector::TIB || subid0 == StripSubdetector::TOB || subid0 == (int) PixelSubdetector::PixelBarrel);
509  gp[1] = hit1->globalPosition();
510  ge[1] = hit1->globalPositionError();
511  int subid1 = hit1->geographicalId().subdetId();
512  bl[1] = (subid1 == StripSubdetector::TIB || subid1 == StripSubdetector::TOB || subid1 == (int) PixelSubdetector::PixelBarrel);
513  gp[2] = hit2->globalPosition();
514  ge[2] = hit2->globalPositionError();
515  int subid2 = hit2->geographicalId().subdetId();
516  bl[2] = (subid2 == StripSubdetector::TIB || subid2 == StripSubdetector::TOB || subid2 == (int) PixelSubdetector::PixelBarrel);
517  RZLine rzLine(gp,ge,bl);
518  float cottheta, intercept, covss, covii, covsi;
519  rzLine.fit(cottheta, intercept, covss, covii, covsi);
520  float chi2 = rzLine.chi2(cottheta, intercept);
521 
522 #ifdef EDM_ML_DEBUG
523  bool debugTriplet = debugPair && hit2->rawId()==debug_Id2;
524 #endif
525  IfLogTrace(debugTriplet, "MultiHitGeneratorFromChi2") << endl << "triplet candidate in debug id" << std::endl
526  << "hit in id="<<hit2->rawId()<<" (from KDTree) with pos: " << KDdata->hit()->globalPosition()
527  << " refitted: " << hit2->globalPosition()
528  << " chi2: " << chi2;
529  // should fix nan
530  if ( (chi2 > maxChi2) | edm::isNotFinite(chi2) ) continue;
531 
532  if (chi2VsPtCut) {
533 
534  FastCircle theCircle(hit2->globalPosition(),hit1->globalPosition(),hit0->globalPosition());
535  float tesla0 = 0.1*nomField;
536  float rho = theCircle.rho();
537  float cm2GeV = 0.01 * 0.3*tesla0;
538  float pt = cm2GeV * rho;
539  IfLogTrace(debugTriplet, "MultiHitGeneratorFromChi2") << "triplet pT=" << pt;
540  if (pt<region.ptMin()) continue;
541 
542  if (chi2_cuts.size()==4) {
543  if ( ( pt>pt_interv[0] && pt<=pt_interv[1] && chi2 > chi2_cuts[0] ) ||
544  ( pt>pt_interv[1] && pt<=pt_interv[2] && chi2 > chi2_cuts[1] ) ||
545  ( pt>pt_interv[2] && pt<=pt_interv[3] && chi2 > chi2_cuts[2] ) ||
546  ( pt>pt_interv[3] && chi2 > chi2_cuts[3] ) ) continue;
547  }
548 
549  // apparently this takes too much time...
550  // if (chi2_cuts.size()>1) {
551  // int ncuts = chi2_cuts.size();
552  // if ( pt<=pt_interv[0] && chi2 > chi2_cuts[0] ) continue;
553  // bool pass = true;
554  // for (int icut=1; icut<ncuts-1; icut++){
555  // if ( pt>pt_interv[icut-1] && pt<=pt_interv[icut] && chi2 > chi2_cuts[icut] ) pass=false;
556  // }
557  // if (!pass) continue;
558  // if ( pt>pt_interv[ncuts-2] && chi2 > chi2_cuts[ncuts-1] ) continue;
559  // if (hit0->rawId()==debug_Id0 && hit1->rawId()==debug_Id1 && hit2->rawId()==debug_Id2) {
560  // LogTrace("MultiHitGeneratorFromChi2") << "triplet passed chi2 vs pt cut" << std::endl;
561  // }
562  // }
563 
564  }
565 
566  if (theMaxElement!=0 && result.size() >= theMaxElement) {
567  result.clear();
568  edm::LogError("TooManyTriplets")<<" number of triples exceed maximum. no triplets produced.";
569  return;
570  }
571  IfLogTrace(debugPair, "MultiHitGeneratorFromChi2") << "triplet made";
572  //result.push_back(SeedingHitSet(hit0, hit1, hit2));
573  /* no refit so keep only hit2
574  assert(tripletFromThisLayer.empty());
575  assert(hit0.isOwn()); assert(hit1.isOwn());assert(hit2.isOwn());
576  tripletFromThisLayer.emplace_back(std::move(hit0));
577  tripletFromThisLayer.emplace_back(std::move(hit1));
578  tripletFromThisLayer.emplace_back(std::move(hit2));
579  assert(hit0.isEmpty()); assert(hit1.isEmpty());assert(hit2.isEmpty());
580  */
581  bestL2 = std::move(hit2);
582  chi2FromThisLayer = chi2;
583  foundTripletsFromPair++;
584  if (foundTripletsFromPair>=2) {
585  usePair=true;
586  IfLogTrace(debugPair, "MultiHitGeneratorFromChi2") << "using pair";
587  break;
588  }
589  }//loop over hits in KDTree
590 
591  if (usePair) break;
592  else {
593  //if there is one triplet in more than one layer, try picking the one with best chi2
594  if (chi2FromThisLayer<minChi2) {
595  bestH2 = std::move(bestL2);
596  minChi2 = chi2FromThisLayer;
597  }
598  /*
599  else {
600  if (!bestH2 && foundTripletsFromPair>0)
601  LogTrace("MultiHitGeneratorFromChi2") << "what?? " << minChi2 << ' ' << chi2FromThisLayer;
602  }
603  */
604  }
605 
606  }//loop over layers
607 
608  if (foundTripletsFromPair==0) continue;
609 
610  //push back only (max) once per pair
611  IfLogTrace(debugPair, "MultiHitGeneratorFromChi2") << "Done seed #" << result.size();
612  if (usePair) result.push_back(SeedingHitSet(ip->inner(), ip->outer()));
613  else {
614  assert(1==foundTripletsFromPair);
615  assert(bestH2);
616  result.emplace_back(&*hit0,&*hit1,&*bestH2);
617  assert(hit0.isOwn()); assert(hit1.isOwn());
618  cache.emplace_back(const_cast<BaseTrackerRecHit*>(hit0.release()));
619  cache.emplace_back(const_cast<BaseTrackerRecHit*>(hit1.release()));
620  cache.emplace_back(std::move(bestH2));
621  assert(hit0.empty()); assert(hit1.empty());assert(!bestH2);
622  }
623  // LogTrace("MultiHitGeneratorFromChi2") << (usePair ? "pair " : "triplet ") << minChi2 <<' ' << cache.size();
624 
625 
626  }//loop over pairs
627  LogTrace("MultiHitGeneratorFromChi2") << "triplet size=" << result.size();
628 }
629 
630 bool MultiHitGeneratorFromChi2::checkPhiInRange(float phi, float phi1, float phi2) const
631 { while (phi > phi2) phi -= 2. * M_PI;
632  while (phi < phi1) phi += 2. * M_PI;
633  return phi <= phi2;
634 }
635 
636 std::pair<float, float>
637 MultiHitGeneratorFromChi2::mergePhiRanges(const std::pair<float, float> &r1,
638  const std::pair<float, float> &r2) const
639 { float r2Min = r2.first;
640  float r2Max = r2.second;
641  while (r1.first - r2Min > +M_PI) r2Min += 2. * M_PI, r2Max += 2. * M_PI;
642  while (r1.first - r2Min < -M_PI) r2Min -= 2. * M_PI, r2Max -= 2. * M_PI;
643  //LogTrace("MultiHitGeneratorFromChi2") << "mergePhiRanges " << fabs(r1.first-r2Min) << " " << fabs(r1.second-r2Max);
644  return std::make_pair(min(r1.first, r2Min), max(r1.second, r2Max));
645 }
646 
648  HitOwnPtr & hit2,
649  TrajectoryStateOnSurface& state1,
650  TrajectoryStateOnSurface& state2,
651  const TrackingRegion& region, float nomField, bool isDebug) {
652 
653  //these need to be sorted in R
654  GlobalPoint gp0 = region.origin();
655  GlobalPoint gp1 = hit1->globalPosition();
656  GlobalPoint gp2 = hit2->globalPosition();
657 
658  IfLogTrace(isDebug, "MultiHitGeneratorFromChi2") << "positions before refitting: " << hit1->globalPosition() << " " << hit2->globalPosition();
659 
660  FastCircle theCircle(gp2,gp1,gp0);
661  GlobalPoint cc(theCircle.x0(),theCircle.y0(),0);
662  float tesla0 = 0.1*nomField;
663  float rho = theCircle.rho();
664  float cm2GeV = 0.01 * 0.3*tesla0;
665  float pt = cm2GeV * rho;
666 
667  GlobalVector vec20 = gp2-gp0;
668  //if (isDebug) { cout << "vec20.eta=" << vec20.eta() << endl; }
669 
670  GlobalVector p0( gp0.y()-cc.y(), -gp0.x()+cc.x(), 0. );
671  p0 = p0*pt/p0.perp();
672  GlobalVector p1( gp1.y()-cc.y(), -gp1.x()+cc.x(), 0. );
673  p1 = p1*pt/p1.perp();
674  GlobalVector p2( gp2.y()-cc.y(), -gp2.x()+cc.x(), 0. );
675  p2 = p2*pt/p2.perp();
676 
677  //check sign according to scalar product
678  if ( (p0.x()*(gp1.x()-gp0.x())+p0.y()*(gp1.y()-gp0.y()) ) < 0 ) {
679  p0*=-1.;
680  p1*=-1.;
681  p2*=-1.;
682  }
683 
684  //now set z component
685  p0 = GlobalVector(p0.x(),p0.y(),p0.perp()/tan(vec20.theta()));
686  p1 = GlobalVector(p1.x(),p1.y(),p1.perp()/tan(vec20.theta()));
687  p2 = GlobalVector(p2.x(),p2.y(),p2.perp()/tan(vec20.theta()));
688 
689  //get charge from vectorial product
690  TrackCharge q = 1;
691  if ((gp1-cc).x()*p1.y() - (gp1-cc).y()*p1.x() > 0) q =-q;
692 
694  state1 = TrajectoryStateOnSurface(kine1,*hit1->surface());
695  hit1.reset((SeedingHitSet::RecHitPointer)(cloner(*hit1,state1)));
696 
698  state2 = TrajectoryStateOnSurface(kine2,*hit2->surface());
699  hit2.reset((SeedingHitSet::RecHitPointer)(cloner(*hit2,state2)));
700 
701  IfLogTrace(isDebug, "MultiHitGeneratorFromChi2") << "charge=" << q << std::endl
702  << "state1 pt=" << state1.globalMomentum().perp() << " eta=" << state1.globalMomentum().eta() << " phi=" << state1.globalMomentum().phi() << std::endl
703  << "state2 pt=" << state2.globalMomentum().perp() << " eta=" << state2.globalMomentum().eta() << " phi=" << state2.globalMomentum().phi() << std::endl
704  << "positions after refitting: " << hit1->globalPosition() << " " << hit2->globalPosition();
705 }
706 
707 /*
708 void MultiHitGeneratorFromChi2::refit3Hits(HitOwnPtr & hit0,
709  HitOwnPtr & hit1,
710  HitOwnPtr & hit2,
711  TrajectoryStateOnSurface& state0,
712  TrajectoryStateOnSurface& state1,
713  TrajectoryStateOnSurface& state2,
714  float nomField, bool isDebug) {
715 
716  //these need to be sorted in R
717  GlobalPoint gp0 = hit0->globalPosition();
718  GlobalPoint gp1 = hit1->globalPosition();
719  GlobalPoint gp2 = hit2->globalPosition();
720 
721  IfLogTrace(isDebug, "MultiHitGeneratorFromChi2") << "positions before refitting: " << hit0->globalPosition() << " " << hit1->globalPosition() << " " << hit2->globalPosition();
722  }
723 
724  FastCircle theCircle(gp2,gp1,gp0);
725  GlobalPoint cc(theCircle.x0(),theCircle.y0(),0);
726  float tesla0 = 0.1*nomField;
727  float rho = theCircle.rho();
728  float cm2GeV = 0.01 * 0.3*tesla0;
729  float pt = cm2GeV * rho;
730 
731  GlobalVector vec20 = gp2-gp0;
732  //IfLogTrace(isDebug, "MultiHitGeneratorFromChi2") << "vec20.eta=" << vec20.eta();
733 
734  GlobalVector p0( gp0.y()-cc.y(), -gp0.x()+cc.x(), 0. );
735  p0 = p0*pt/p0.perp();
736  GlobalVector p1( gp1.y()-cc.y(), -gp1.x()+cc.x(), 0. );
737  p1 = p1*pt/p1.perp();
738  GlobalVector p2( gp2.y()-cc.y(), -gp2.x()+cc.x(), 0. );
739  p2 = p2*pt/p2.perp();
740 
741  //check sign according to scalar product
742  if ( (p0.x()*(gp1.x()-gp0.x())+p0.y()*(gp1.y()-gp0.y()) ) < 0 ) {
743  p0*=-1.;
744  p1*=-1.;
745  p2*=-1.;
746  }
747 
748  //now set z component
749  p0 = GlobalVector(p0.x(),p0.y(),p0.perp()/tan(vec20.theta()));
750  p1 = GlobalVector(p1.x(),p1.y(),p1.perp()/tan(vec20.theta()));
751  p2 = GlobalVector(p2.x(),p2.y(),p2.perp()/tan(vec20.theta()));
752 
753  //get charge from vectorial product
754  TrackCharge q = 1;
755  if ((gp1-cc).x()*p1.y() - (gp1-cc).y()*p1.x() > 0) q =-q;
756 
757  GlobalTrajectoryParameters kine0 = GlobalTrajectoryParameters(gp0, p0, q, &*bfield);
758  state0 = TrajectoryStateOnSurface(kine0,*hit0->surface());
759  hit0 = hit0->clone(state0);
760 
761  GlobalTrajectoryParameters kine1 = GlobalTrajectoryParameters(gp1, p1, q, &*bfield);
762  state1 = TrajectoryStateOnSurface(kine1,*hit1->surface());
763  hit1 = hit1->clone(state1);
764 
765  GlobalTrajectoryParameters kine2 = GlobalTrajectoryParameters(gp2, p2, q, &*bfield);
766  state2 = TrajectoryStateOnSurface(kine2,*hit2->surface());
767  hit2 = hit2->clone(state2);
768 
769  IfLogTrace(isDebug, "MultiHitGeneratorFromChi2") << "charge=" << q << std::endl
770  << "state0 pt=" << state0.globalMomentum().perp() << " eta=" << state0.globalMomentum().eta() << " phi=" << state0.globalMomentum().phi() << std::endl
771  << "state1 pt=" << state1.globalMomentum().perp() << " eta=" << state1.globalMomentum().eta() << " phi=" << state1.globalMomentum().phi() << std::endl
772  << "state2 pt=" << state2.globalMomentum().perp() << " eta=" << state2.globalMomentum().eta() << " phi=" << state2.globalMomentum().phi() << std::endl
773  << "positions after refitting: " << hit0->globalPosition() << " " << hit1->globalPosition() << " " << hit2->globalPosition();
774 }
775 */
776 
#define LogDebug(id)
float originRBound() const
bounds the particle vertex in the transverse plane
T getParameter(std::string const &) const
std::pair< float, float > mergePhiRanges(const std::pair< float, float > &r1, const std::pair< float, float > &r2) const
void reset()
Definition: ProxyBase11.h:64
unsigned int stereoId() const
void build(std::vector< KDTreeNodeInfo > &eltList, const KDTreeBox &region)
tuple cfg
Definition: looper.py:293
bool empty() const
Definition: mayown_ptr.h:59
T perp() const
Definition: PV3DBase.h:72
std::unique_ptr< BaseTrackerRecHit > cacheHitPointer
GlobalPoint const & origin() const
virtual Location location() const =0
Which part of the detector (barrel, endcap)
int nominalValue() const
The nominal field value for this map in kGauss.
Definition: MagneticField.h:56
double x0() const
Definition: FastCircle.h:50
TkTransientTrackingRecHitBuilder const * builder
bool isOwn() const
Definition: mayown_ptr.h:26
SiStripCluster const & monoCluster() const
const ClusterShapeHitFilter * filter
assert(m_qm.get())
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
bool exists(std::string const &parameterName) const
checks if a parameter exists
#define IfLogTrace(cond, cat)
bool ev
virtual unsigned int size() const
void search(const KDTreeBox &searchBox, std::vector< KDTreeNodeInfo > &resRecHitList)
double rho() const
Definition: FastCircle.h:54
void reset()
Definition: mayown_ptr.h:56
T const * release()
Definition: mayown_ptr.h:55
Range curvature(double transverseIP) const
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
LocalVector localMomentum() const
#define unlikely(x)
U second(std::pair< T, U > const &p)
virtual unsigned int size() const
std::pair< HitIter, HitIter > Range
int TrackCharge
Definition: TrackCharge.h:4
BaseTrackerRecHit const * ConstRecHitPointer
Definition: SeedingHitSet.h:11
bool isNotFinite(T x)
Definition: isFinite.h:10
virtual void hitSets(const TrackingRegion &region, OrderedMultiHits &trs, const edm::Event &ev, const edm::EventSetup &es, SeedingLayerSetsHits::SeedingLayerSet pairLayers, std::vector< SeedingLayerSetsHits::SeedingLayer > thirdLayers)
T curvature(T InversePt, const edm::EventSetup &iSetup)
std::vector< HitWithPhi >::const_iterator HitIter
T z() const
Definition: PV3DBase.h:64
tuple result
Definition: query.py:137
def move
Definition: eostools.py:510
void initES(const edm::EventSetup &es) override
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
MultiHitGeneratorFromChi2(const edm::ParameterSet &cfg)
T min(T a, T b)
Definition: MathUtil.h:58
double p2[4]
Definition: TauolaWrapper.h:90
#define LogTrace(id)
#define M_PI
bool isCompatible(const SiPixelRecHit &recHit, const LocalVector &ldir, const SiPixelClusterShapeCache &clusterShapeCache, PixelData const *pd=0) const
SeedingHitSet::ConstRecHitPointer Hit
BaseTrackerRecHit const * Hit
Definition: RZLine.h:8
SiStripRecHit2D originalHit() const
Definition: DetId.h:18
PixelRecoRange< float > Range
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
float ptMin() const
minimal pt of interest
Geom::Phi< T > phi() const
T eta() const
Definition: PV3DBase.h:76
double y0() const
Definition: FastCircle.h:52
SiStripCluster const & stereoCluster() const
double p1[4]
Definition: TauolaWrapper.h:89
GlobalVector globalMomentum() const
double pi()
Definition: Pi.h:31
double twoPi()
Definition: Pi.h:32
std::unique_ptr< HitPairGeneratorFromLayerPair > thePairGenerator
void fit(float &cotTheta, float &intercept, float &covss, float &covii, float &covsi) const
Definition: RZLine.cc:44
float chi2(float cotTheta, float intercept) const
Definition: RZLine.cc:50
unsigned int monoId() const
volatile std::atomic< bool > shutdown_flag false
T x() const
Definition: PV3DBase.h:62
float phi(float curvature, float radius) const
static const float fnSigmaRZ
tuple size
Write out results.
bool checkPhiInRange(float phi, float phi1, float phi2) const
Global3DVector GlobalVector
Definition: GlobalVector.h:10
void refit2Hits(HitOwnPtr &hit0, HitOwnPtr &hit1, TrajectoryStateOnSurface &tsos0, TrajectoryStateOnSurface &tsos1, const TrackingRegion &region, float nomField, bool isDebug)
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11