00001 #include "RecoPixelVertexing/PixelTriplets/plugins/PixelTripletLargeTipGenerator.h"
00002
00003 #include "ThirdHitPredictionFromCircle.h"
00004 #include "ThirdHitRZPrediction.h"
00005 #include "RecoTracker/TkMSParametrization/interface/PixelRecoUtilities.h"
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007
00008 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00009 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00010 #include "RecoPixelVertexing/PixelTriplets/plugins/ThirdHitCorrection.h"
00011 #include "RecoTracker/TkHitPairs/interface/RecHitsSortedInPhi.h"
00012
00013 #include "MatchedHitRZCorrectionFromBending.h"
00014
00015
00016 #include "RecoPixelVertexing/PixelTriplets/plugins/KDTreeLinkerAlgo.h"
00017 #include "RecoPixelVertexing/PixelTriplets/plugins/KDTreeLinkerTools.h"
00018
00019 #include <algorithm>
00020 #include <iostream>
00021 #include <vector>
00022 #include <cmath>
00023 #include <map>
00024
00025 using namespace std;
00026 using namespace ctfseeding;
00027
00028 typedef PixelRecoRange<float> Range;
00029
00030 typedef ThirdHitPredictionFromCircle::HelixRZ HelixRZ;
00031
00032 namespace {
00033 struct LayerRZPredictions {
00034 ThirdHitRZPrediction<PixelRecoLineRZ> line;
00035 ThirdHitRZPrediction<HelixRZ> helix1, helix2;
00036 MatchedHitRZCorrectionFromBending rzPositionFixup;
00037 };
00038 }
00039
00040 constexpr double nSigmaRZ = 3.4641016151377544;
00041 constexpr double nSigmaPhi = 3.;
00042 static float fnSigmaRZ = std::sqrt(12.f);
00043
00044 PixelTripletLargeTipGenerator::PixelTripletLargeTipGenerator(const edm::ParameterSet& cfg)
00045 : thePairGenerator(0),
00046 theLayerCache(0),
00047 useFixedPreFiltering(cfg.getParameter<bool>("useFixedPreFiltering")),
00048 extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
00049 extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
00050 useMScat(cfg.getParameter<bool>("useMultScattering")),
00051 useBend(cfg.getParameter<bool>("useBending"))
00052 { theMaxElement=cfg.getParameter<unsigned int>("maxElement");
00053 if (useFixedPreFiltering)
00054 dphi = cfg.getParameter<double>("phiPreFiltering");
00055 }
00056
00057 void PixelTripletLargeTipGenerator::init(const HitPairGenerator & pairs,
00058 const std::vector<SeedingLayer> &layers,
00059 LayerCacheType *layerCache)
00060 {
00061 thePairGenerator = pairs.clone();
00062 theLayers = layers;
00063 theLayerCache = layerCache;
00064 }
00065
00066
00067 namespace {
00068 inline
00069 bool intersect(Range &range, const Range &second)
00070 {
00071 if (range.first > second.max() || range.second < second.min())
00072 return false;
00073 if (range.first < second.min())
00074 range.first = second.min();
00075 if (range.second > second.max())
00076 range.second = second.max();
00077 return range.first < range.second;
00078 }
00079 }
00080
00081 void PixelTripletLargeTipGenerator::hitTriplets(const TrackingRegion& region,
00082 OrderedHitTriplets & result,
00083 const edm::Event & ev,
00084 const edm::EventSetup& es)
00085 {
00086 edm::ESHandle<TrackerGeometry> tracker;
00087 es.get<TrackerDigiGeometryRecord>().get(tracker);
00088
00089
00090 edm::ESHandle<TrackerTopology> tTopoHand;
00091 es.get<IdealGeometryRecord>().get(tTopoHand);
00092 const TrackerTopology *tTopo=tTopoHand.product();
00093
00094 auto const & doublets = thePairGenerator->doublets(region,ev,es);
00095
00096 if (doublets.empty()) return;
00097
00098 auto outSeq = doublets.detLayer(HitDoublets::outer)->seqNum();
00099
00100
00101 int size = theLayers.size();
00102
00103
00104 using NodeInfo = KDTreeNodeInfo<unsigned int>;
00105 std::vector<NodeInfo > layerTree;
00106 KDTreeLinkerAlgo<unsigned int> hitTree[size];
00107
00108 float rzError[size];
00109 float maxphi = Geom::ftwoPi(), minphi = -maxphi;
00110
00111 LayerRZPredictions mapPred[size];
00112
00113 const RecHitsSortedInPhi * thirdHitMap[size];
00114
00115 for(int il = 0; il < size; il++) {
00116 thirdHitMap[il] = &(*theLayerCache)(&theLayers[il], region, ev, es);
00117 auto const & hits = *thirdHitMap[il];
00118
00119 const DetLayer *layer = theLayers[il].detLayer();
00120 LayerRZPredictions &predRZ = mapPred[il];
00121 predRZ.line.initLayer(layer);
00122 predRZ.helix1.initLayer(layer);
00123 predRZ.helix2.initLayer(layer);
00124 predRZ.line.initTolerance(extraHitRZtolerance);
00125 predRZ.helix1.initTolerance(extraHitRZtolerance);
00126 predRZ.helix2.initTolerance(extraHitRZtolerance);
00127 predRZ.rzPositionFixup = MatchedHitRZCorrectionFromBending(layer,tTopo);
00128
00129 layerTree.clear();
00130 float minv=999999.0; float maxv = -999999.0;
00131 float maxErr=0.0f;
00132 for (unsigned int i=0; i!=hits.size(); ++i) {
00133 auto angle = hits.phi(i);
00134 auto v = hits.gv(i);
00135
00136 minv = std::min(minv,v); maxv = std::max(maxv,v);
00137 float myerr = hits.dv[i];
00138 maxErr = std::max(maxErr,myerr);
00139 layerTree.emplace_back(i, angle, v);
00140 if (angle < 0)
00141 { layerTree.emplace_back(i, angle+Geom::ftwoPi(), v);}
00142 else
00143 { layerTree.emplace_back(i, angle-Geom::ftwoPi(), v);}
00144 }
00145 KDTreeBox phiZ(minphi, maxphi, minv-0.01f, maxv+0.01f);
00146
00147 hitTree[il].build(layerTree, phiZ);
00148 rzError[il] = maxErr;
00149 }
00150
00151 double curv = PixelRecoUtilities::curvature(1. / region.ptMin(), es);
00152
00153 for (std::size_t ip =0; ip!=doublets.size(); ip++) {
00154 auto xi = doublets.x(ip,HitDoublets::inner);
00155 auto yi = doublets.y(ip,HitDoublets::inner);
00156 auto zi = doublets.z(ip,HitDoublets::inner);
00157
00158 auto xo = doublets.x(ip,HitDoublets::outer);
00159 auto yo = doublets.y(ip,HitDoublets::outer);
00160 auto zo = doublets.z(ip,HitDoublets::outer);
00161
00162 GlobalPoint gp1(xi,yi,zi);
00163 GlobalPoint gp2(xo,yo,zo);
00164
00165 PixelRecoLineRZ line(gp1, gp2);
00166 PixelRecoPointRZ point2(gp2.perp(), zo);
00167 ThirdHitPredictionFromCircle predictionRPhi(gp1, gp2, extraHitRPhitolerance);
00168
00169 Range generalCurvature = predictionRPhi.curvature(region.originRBound());
00170 if (!intersect(generalCurvature, Range(-curv, curv))) continue;
00171
00172 for(int il = 0; il < size; il++) {
00173 if (hitTree[il].empty()) continue;
00174 const DetLayer *layer = theLayers[il].detLayer();
00175 bool barrelLayer = layer->isBarrel();
00176
00177 Range curvature = generalCurvature;
00178 ThirdHitCorrection correction(es, region.ptMin(), layer, line, point2, outSeq, useMScat);
00179
00180 LayerRZPredictions &predRZ = mapPred[il];
00181 predRZ.line.initPropagator(&line);
00182
00183 Range rzRange;
00184 if (useBend) {
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203 if (!barrelLayer) {
00204 Range z3s = predRZ.line.detRange();
00205 double z3 = z3s.first < 0 ? std::max(z3s.first, z3s.second)
00206 : std::min(z3s.first, z3s.second);
00207 double maxCurvature = HelixRZ::maxCurvature(&predictionRPhi,
00208 gp1.z(), gp2.z(), z3);
00209 if (!intersect(curvature, Range(-maxCurvature, maxCurvature)))
00210 continue;
00211 }
00212
00213 HelixRZ helix1(&predictionRPhi, gp1.z(), gp2.z(), curvature.first);
00214 HelixRZ helix2(&predictionRPhi, gp1.z(), gp2.z(), curvature.second);
00215
00216 predRZ.helix1.initPropagator(&helix1);
00217 predRZ.helix2.initPropagator(&helix2);
00218
00219 Range rzRanges[2] = { predRZ.helix1(), predRZ.helix2() };
00220 predRZ.helix1.initPropagator(nullptr);
00221 predRZ.helix2.initPropagator(nullptr);
00222
00223 rzRange.first = std::min(rzRanges[0].first, rzRanges[1].first);
00224 rzRange.second = std::max(rzRanges[0].second, rzRanges[1].second);
00225
00226
00227
00228 if (curvature.first * curvature.second < 0.0) {
00229 Range rzLineRange = predRZ.line();
00230 rzRange.first = std::min(rzRange.first, rzLineRange.first);
00231 rzRange.second = std::max(rzRange.second, rzLineRange.second);
00232 }
00233 } else {
00234 rzRange = predRZ.line();
00235 }
00236
00237 if (rzRange.first >= rzRange.second)
00238 continue;
00239
00240 correction.correctRZRange(rzRange);
00241
00242 Range phiRange;
00243 if (useFixedPreFiltering) {
00244 float phi0 = doublets.phi(ip,HitDoublets::outer);
00245 phiRange = Range(phi0 - dphi, phi0 + dphi);
00246 } else {
00247 Range radius;
00248
00249 if (barrelLayer) {
00250 radius = predRZ.line.detRange();
00251 if (!intersect(rzRange, predRZ.line.detSize()))
00252 continue;
00253 } else {
00254 radius = rzRange;
00255 if (!intersect(radius, predRZ.line.detSize()))
00256 continue;
00257 }
00258
00259 Range rPhi1 = predictionRPhi(curvature, radius.first);
00260 Range rPhi2 = predictionRPhi(curvature, radius.second);
00261 correction.correctRPhiRange(rPhi1);
00262 correction.correctRPhiRange(rPhi2);
00263 rPhi1.first /= radius.first;
00264 rPhi1.second /= radius.first;
00265 rPhi2.first /= radius.second;
00266 rPhi2.second /= radius.second;
00267 phiRange = mergePhiRanges(rPhi1, rPhi2);
00268 }
00269
00270 layerTree.clear();
00271 float prmin=phiRange.min(), prmax=phiRange.max();
00272 if ((prmax-prmin) > Geom::ftwoPi())
00273 { prmax=Geom::fpi(); prmin = -Geom::fpi();}
00274 else
00275 { while (prmax>maxphi) { prmin -= Geom::ftwoPi(); prmax -= Geom::ftwoPi();}
00276 while (prmin<minphi) { prmin += Geom::ftwoPi(); prmax += Geom::ftwoPi();}
00277
00278 }
00279 if (barrelLayer) {
00280 Range regMax = predRZ.line.detRange();
00281 Range regMin = predRZ.line(regMax.min());
00282 regMax = predRZ.line(regMax.max());
00283 correction.correctRZRange(regMin);
00284 correction.correctRZRange(regMax);
00285 if (regMax.min() < regMin.min()) { swap(regMax, regMin);}
00286 KDTreeBox phiZ(prmin, prmax,
00287 regMin.min()-fnSigmaRZ*rzError[il],
00288 regMax.max()+fnSigmaRZ*rzError[il]);
00289 hitTree[il].search(phiZ, layerTree);
00290 }
00291 else {
00292 KDTreeBox phiZ(prmin, prmax,
00293 rzRange.min()-fnSigmaRZ*rzError[il],
00294 rzRange.max()+fnSigmaRZ*rzError[il]);
00295 hitTree[il].search(phiZ, layerTree);
00296 }
00297
00298 MatchedHitRZCorrectionFromBending l2rzFixup(doublets.hit(ip,HitDoublets::outer)->det()->geographicalId(), tTopo);
00299 MatchedHitRZCorrectionFromBending l3rzFixup = predRZ.rzPositionFixup;
00300
00301 thirdHitMap[il] = &(*theLayerCache)(&theLayers[il], region, ev, es);
00302 auto const & hits = *thirdHitMap[il];
00303 for (auto const & ih : layerTree) {
00304 auto KDdata = ih.data;
00305 GlobalPoint p3 = hits.gp(KDdata);
00306 double p3_r = p3.perp();
00307 double p3_z = p3.z();
00308 float p3_phi = hits.phi(KDdata);
00309
00310 Range rangeRPhi = predictionRPhi(curvature, p3_r);
00311 correction.correctRPhiRange(rangeRPhi);
00312
00313 float ir = 1.f/p3_r;
00314 float phiErr = nSigmaPhi * hits.drphi[KDdata]*ir;
00315 if (!checkPhiInRange(p3_phi, rangeRPhi.first*ir-phiErr, rangeRPhi.second*ir+phiErr))
00316 continue;
00317
00318 Basic2DVector<double> thc(p3.x(), p3.y());
00319
00320 auto curv_ = predictionRPhi.curvature(thc);
00321 double p2_r = point2.r(); double p2_z = point2.z();
00322
00323 l2rzFixup(predictionRPhi, curv_, *doublets.hit(ip,HitDoublets::outer), p2_r, p2_z, tTopo);
00324 l3rzFixup(predictionRPhi, curv_, *hits.theHits[KDdata].hit(), p3_r, p3_z, tTopo);
00325
00326 Range rangeRZ;
00327 if (useBend) {
00328 HelixRZ updatedHelix(&predictionRPhi, gp1.z(), p2_z, curv_);
00329 rangeRZ = predRZ.helix1(barrelLayer ? p3_r : p3_z, updatedHelix);
00330 } else {
00331 float tIP = predictionRPhi.transverseIP(thc);
00332 PixelRecoPointRZ updatedPoint2(p2_r, p2_z);
00333 PixelRecoLineRZ updatedLine(line.origin(), point2, tIP);
00334 rangeRZ = predRZ.line(barrelLayer ? p3_r : p3_z, line);
00335 }
00336 correction.correctRZRange(rangeRZ);
00337
00338 double err = nSigmaRZ * hits.dv[KDdata];
00339
00340 rangeRZ.first -= err, rangeRZ.second += err;
00341
00342 if (!rangeRZ.inside(barrelLayer ? p3_z : p3_r)) continue;
00343
00344 if (theMaxElement!=0 && result.size() >= theMaxElement) {
00345 result.clear();
00346 edm::LogError("TooManyTriplets")<<" number of triples exceed maximum. no triplets produced.";
00347 return;
00348 }
00349 result.emplace_back( doublets.hit(ip,HitDoublets::inner), doublets.hit(ip,HitDoublets::outer), hits.theHits[KDdata].hit());
00350 }
00351 }
00352 }
00353
00354 }
00355
00356 bool PixelTripletLargeTipGenerator::checkPhiInRange(float phi, float phi1, float phi2) const
00357 { while (phi > phi2) phi -= 2. * M_PI;
00358 while (phi < phi1) phi += 2. * M_PI;
00359 return phi <= phi2;
00360 }
00361
00362 std::pair<float, float>
00363 PixelTripletLargeTipGenerator::mergePhiRanges(const std::pair<float, float> &r1,
00364 const std::pair<float, float> &r2) const
00365 { float r2Min = r2.first;
00366 float r2Max = r2.second;
00367 while (r1.first - r2Min > +M_PI) r2Min += 2. * M_PI, r2Max += 2. * M_PI;
00368 while (r1.first - r2Min < -M_PI) r2Min -= 2. * M_PI, r2Max -= 2. * M_PI;
00369 return std::make_pair(min(r1.first, r2Min), max(r1.second, r2Max));
00370 }