00001 #include "RecoPixelVertexing/PixelTriplets/plugins/PixelTripletHLTGenerator.h"
00002
00003 #include "ThirdHitPredictionFromInvParabola.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 "ThirdHitCorrection.h"
00011 #include "RecoTracker/TkHitPairs/interface/RecHitsSortedInPhi.h"
00012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00013 #include <iostream>
00014
00015 #include "RecoTracker/TkSeedingLayers/interface/SeedComparitorFactory.h"
00016 #include "RecoTracker/TkSeedingLayers/interface/SeedComparitor.h"
00017
00018 #include "DataFormats/GeometryVector/interface/Pi.h"
00019
00020
00021 #include "RecoPixelVertexing/PixelTriplets/plugins/KDTreeLinkerAlgo.h"
00022 #include "RecoPixelVertexing/PixelTriplets/plugins/KDTreeLinkerTools.h"
00023
00024 #include<cstdio>
00025
00026 using pixelrecoutilities::LongitudinalBendingCorrection;
00027 typedef PixelRecoRange<float> Range;
00028
00029 using namespace std;
00030 using namespace ctfseeding;
00031
00032 PixelTripletHLTGenerator:: PixelTripletHLTGenerator(const edm::ParameterSet& cfg)
00033 : thePairGenerator(0),
00034 theLayerCache(0),
00035 useFixedPreFiltering(cfg.getParameter<bool>("useFixedPreFiltering")),
00036 extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
00037 extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
00038 useMScat(cfg.getParameter<bool>("useMultScattering")),
00039 useBend(cfg.getParameter<bool>("useBending"))
00040 {
00041 theMaxElement=cfg.getParameter<unsigned int>("maxElement");
00042 dphi = (useFixedPreFiltering) ? cfg.getParameter<double>("phiPreFiltering") : 0;
00043
00044 edm::ParameterSet comparitorPSet =
00045 cfg.getParameter<edm::ParameterSet>("SeedComparitorPSet");
00046 std::string comparitorName = comparitorPSet.getParameter<std::string>("ComponentName");
00047 theComparitor = (comparitorName == "none") ?
00048 0 : SeedComparitorFactory::get()->create( comparitorName, comparitorPSet);
00049 }
00050
00051 PixelTripletHLTGenerator::~PixelTripletHLTGenerator() {
00052 delete thePairGenerator;
00053 delete theComparitor;
00054 }
00055
00056 void PixelTripletHLTGenerator::init( const HitPairGenerator & pairs,
00057 const std::vector<SeedingLayer> & layers,
00058 LayerCacheType* layerCache)
00059 {
00060 thePairGenerator = pairs.clone();
00061 theLayers = layers;
00062 theLayerCache = layerCache;
00063 }
00064
00065 void PixelTripletHLTGenerator::hitTriplets(const TrackingRegion& region,
00066 OrderedHitTriplets & result,
00067 const edm::Event & ev,
00068 const edm::EventSetup& es)
00069 {
00070
00071 if (theComparitor) theComparitor->init(es);
00072
00073 auto const & doublets = thePairGenerator->doublets(region,ev,es);
00074
00075 if (doublets.empty()) return;
00076
00077 auto outSeq = doublets.detLayer(HitDoublets::outer)->seqNum();
00078
00079
00080
00081
00082 float regOffset = region.origin().perp();
00083 int size = theLayers.size();
00084
00085 ThirdHitRZPrediction<PixelRecoLineRZ> preds[size];
00086
00087 const RecHitsSortedInPhi * thirdHitMap[size];
00088 typedef RecHitsSortedInPhi::Hit Hit;
00089
00090 using NodeInfo = KDTreeNodeInfo<unsigned int>;
00091 std::vector<NodeInfo > layerTree;
00092
00093 KDTreeLinkerAlgo<unsigned int> hitTree[size];
00094 float rzError[size];
00095 float maxphi = Geom::ftwoPi(), minphi = -maxphi;
00096
00097
00098 for (int il=0; il!=size; ++il) {
00099 thirdHitMap[il] = &(*theLayerCache)(&theLayers[il], region, ev, es);
00100 auto const & hits = *thirdHitMap[il];
00101 ThirdHitRZPrediction<PixelRecoLineRZ> & pred = preds[il];
00102 pred.initLayer(theLayers[il].detLayer());
00103 pred.initTolerance(extraHitRZtolerance);
00104
00105 layerTree.clear();
00106 float minv=999999.0, maxv= -999999.0;
00107 float maxErr=0.0f;
00108 for (unsigned int i=0; i!=hits.size(); ++i) {
00109 auto angle = hits.phi(i);
00110 auto v = hits.gv(i);
00111
00112 minv = std::min(minv,v); maxv = std::max(maxv,v);
00113 float myerr = hits.dv[i];
00114 maxErr = std::max(maxErr,myerr);
00115 layerTree.emplace_back(i, angle, v);
00116 if (angle < 0)
00117 { layerTree.emplace_back(i, angle+Geom::ftwoPi(), v);}
00118 else
00119 { layerTree.emplace_back(i, angle-Geom::ftwoPi(), v);}
00120 }
00121 KDTreeBox phiZ(minphi, maxphi, minv-0.01f, maxv+0.01f);
00122
00123 hitTree[il].build(layerTree, phiZ);
00124 rzError[il] = maxErr;
00125
00126 }
00127
00128 float imppar = region.originRBound();
00129 float imppartmp = region.originRBound()+region.origin().perp();
00130 float curv = PixelRecoUtilities::curvature(1.f/region.ptMin(), es);
00131
00132 for (std::size_t ip =0; ip!=doublets.size(); ip++) {
00133 auto xi = doublets.x(ip,HitDoublets::inner);
00134 auto yi = doublets.y(ip,HitDoublets::inner);
00135 auto zi = doublets.z(ip,HitDoublets::inner);
00136 auto rvi = doublets.rv(ip,HitDoublets::inner);
00137 auto xo = doublets.x(ip,HitDoublets::outer);
00138 auto yo = doublets.y(ip,HitDoublets::outer);
00139 auto zo = doublets.z(ip,HitDoublets::outer);
00140 auto rvo = doublets.rv(ip,HitDoublets::outer);
00141
00142 PixelRecoPointRZ point1(rvi, zi);
00143 PixelRecoPointRZ point2(rvo, zo);
00144 PixelRecoLineRZ line(point1, point2);
00145 ThirdHitPredictionFromInvParabola predictionRPhi(xi-region.origin().x(),yi-region.origin().y(),
00146 xo-region.origin().x(),yo-region.origin().y(),
00147 imppar,curv,extraHitRPhitolerance);
00148 ThirdHitPredictionFromInvParabola predictionRPhitmp(xi,yi,xo,yo,imppartmp,curv,extraHitRPhitolerance);
00149
00150
00151
00152
00153
00154
00155 for (int il=0; il!=size; ++il) {
00156 if (hitTree[il].empty()) continue;
00157
00158 auto const & hits = *thirdHitMap[il];
00159
00160 const DetLayer * layer = theLayers[il].detLayer();
00161 auto barrelLayer = layer->isBarrel();
00162
00163 ThirdHitCorrection correction(es, region.ptMin(), layer, line, point2, outSeq, useMScat, useBend);
00164
00165 ThirdHitRZPrediction<PixelRecoLineRZ> & predictionRZ = preds[il];
00166
00167 predictionRZ.initPropagator(&line);
00168 Range rzRange = predictionRZ();
00169 correction.correctRZRange(rzRange);
00170
00171 Range phiRange;
00172 if (useFixedPreFiltering) {
00173 float phi0 = doublets.phi(ip,HitDoublets::outer);
00174 phiRange = Range(phi0-dphi,phi0+dphi);
00175 }
00176 else {
00177 Range radius;
00178 if (barrelLayer) {
00179 radius = predictionRZ.detRange();
00180 } else {
00181 radius = Range(max(rzRange.min(), predictionRZ.detSize().min()),
00182 min(rzRange.max(), predictionRZ.detSize().max()) );
00183 }
00184 if (radius.empty()) continue;
00185
00186
00187
00188 Range rPhi1m = predictionRPhitmp(radius.max(), -1);
00189 Range rPhi1p = predictionRPhitmp(radius.max(), 1);
00190 Range rPhi2m = predictionRPhitmp(radius.min(), -1);
00191 Range rPhi2p = predictionRPhitmp(radius.min(), 1);
00192 Range rPhi1 = rPhi1m.sum(rPhi1p);
00193 Range rPhi2 = rPhi2m.sum(rPhi2p);
00194 correction.correctRPhiRange(rPhi1);
00195 correction.correctRPhiRange(rPhi2);
00196 rPhi1.first /= radius.max();
00197 rPhi1.second /= radius.max();
00198 rPhi2.first /= radius.min();
00199 rPhi2.second /= radius.min();
00200 phiRange = mergePhiRanges(rPhi1,rPhi2);
00201 }
00202
00203 constexpr float nSigmaRZ = std::sqrt(12.f);
00204 constexpr float nSigmaPhi = 3.f;
00205
00206 layerTree.clear();
00207 float prmin=phiRange.min(), prmax=phiRange.max();
00208 if ((prmax-prmin) > Geom::ftwoPi())
00209 { prmax=Geom::fpi(); prmin = -Geom::fpi();}
00210 else
00211 { while (prmax>maxphi) { prmin -= Geom::ftwoPi(); prmax -= Geom::ftwoPi();}
00212 while (prmin<minphi) { prmin += Geom::ftwoPi(); prmax += Geom::ftwoPi();}
00213
00214 }
00215 if (barrelLayer)
00216 {
00217 Range regMax = predictionRZ.detRange();
00218 Range regMin = predictionRZ(regMax.min()-regOffset);
00219 regMax = predictionRZ(regMax.max()+regOffset);
00220 correction.correctRZRange(regMin);
00221 correction.correctRZRange(regMax);
00222 if (regMax.min() < regMin.min()) { swap(regMax, regMin);}
00223 KDTreeBox phiZ(prmin, prmax, regMin.min()-nSigmaRZ*rzError[il], regMax.max()+nSigmaRZ*rzError[il]);
00224 hitTree[il].search(phiZ, layerTree);
00225 }
00226 else
00227 {
00228 KDTreeBox phiZ(prmin, prmax,
00229 rzRange.min()-regOffset-nSigmaRZ*rzError[il],
00230 rzRange.max()+regOffset+nSigmaRZ*rzError[il]);
00231 hitTree[il].search(phiZ, layerTree);
00232 }
00233
00234
00235
00236
00237
00238 for (auto const & ih : layerTree) {
00239
00240 if (theMaxElement!=0 && result.size() >= theMaxElement){
00241 result.clear();
00242 edm::LogError("TooManyTriplets")<<" number of triples exceeds maximum. no triplets produced.";
00243 return;
00244 }
00245
00246 auto KDdata = ih.data;
00247 float p3_u = hits.u[KDdata];
00248 float p3_v = hits.v[KDdata];
00249 float p3_phi = hits.lphi[KDdata];
00250
00251
00252
00253
00254
00255 Range allowed = predictionRZ(p3_u);
00256 correction.correctRZRange(allowed);
00257 float vErr = nSigmaRZ *hits.dv[KDdata];
00258 Range hitRange(p3_v-vErr, p3_v+vErr);
00259 Range crossingRange = allowed.intersection(hitRange);
00260 if (crossingRange.empty()) continue;
00261
00262 float ir = 1.f/hits.rv(KDdata);
00263 float phiErr = nSigmaPhi * hits.drphi[KDdata]*ir;
00264 for (int icharge=-1; icharge <=1; icharge+=2) {
00265 Range rangeRPhi = predictionRPhi(hits.rv(KDdata), icharge);
00266 correction.correctRPhiRange(rangeRPhi);
00267 if (checkPhiInRange(p3_phi, rangeRPhi.first*ir-phiErr, rangeRPhi.second*ir+phiErr)) {
00268
00269 OrderedHitTriplet hittriplet( doublets.hit(ip,HitDoublets::inner), doublets.hit(ip,HitDoublets::outer), hits.theHits[KDdata].hit());
00270 if (!theComparitor || theComparitor->compatible(hittriplet,region) ) {
00271 result.push_back( hittriplet );
00272 } else {
00273 LogDebug("RejectedTriplet") << "rejected triplet from comparitor ";
00274 }
00275 break;
00276 }
00277 }
00278 }
00279 }
00280 }
00281
00282 }
00283
00284 bool PixelTripletHLTGenerator::checkPhiInRange(float phi, float phi1, float phi2) const
00285 {
00286 while (phi > phi2) phi -= Geom::ftwoPi();
00287 while (phi < phi1) phi += Geom::ftwoPi();
00288 return ( (phi1 <= phi) && (phi <= phi2) );
00289 }
00290
00291 std::pair<float,float> PixelTripletHLTGenerator::mergePhiRanges(const std::pair<float,float>& r1,
00292 const std::pair<float,float>& r2) const
00293 { float r2_min=r2.first;
00294 float r2_max=r2.second;
00295 while (r1.first-r2_min > Geom::fpi()) { r2_min += Geom::ftwoPi(); r2_max += Geom::ftwoPi();}
00296 while (r1.first-r2_min < -Geom::fpi()) { r2_min -= Geom::ftwoPi(); r2_max -= Geom::ftwoPi(); }
00297
00298 return std::make_pair(min(r1.first,r2_min),max(r1.second,r2_max));
00299 }