CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/RecoPixelVertexing/PixelTriplets/plugins/PixelTripletLargeTipGenerator.cc

Go to the documentation of this file.
00001 #include "PixelTripletLargeTipGenerator.h"
00002 
00003 #include "RecoPixelVertexing/PixelTriplets/interface/ThirdHitPredictionFromCircle.h"
00004 #include "RecoPixelVertexing/PixelTriplets/interface/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 
00013 #include "MatchedHitRZCorrectionFromBending.h"
00014 
00015 #include <algorithm>
00016 #include <iostream>
00017 #include <vector>
00018 #include <cmath>
00019 #include <map>
00020 
00021 using namespace std;
00022 using namespace ctfseeding;
00023 
00024 typedef PixelRecoRange<float> Range;
00025 
00026 typedef ThirdHitPredictionFromCircle::HelixRZ HelixRZ;
00027 
00028 namespace {
00029   struct LayerRZPredictions {
00030     ThirdHitRZPrediction<PixelRecoLineRZ> line;
00031     ThirdHitRZPrediction<HelixRZ> helix1, helix2;
00032     MatchedHitRZCorrectionFromBending rzPositionFixup;
00033   };
00034 }
00035 
00036 static const double nSigmaRZ = 3.4641016151377544; // sqrt(12.)
00037 static const double nSigmaPhi = 3.;
00038 
00039 PixelTripletLargeTipGenerator::PixelTripletLargeTipGenerator(const edm::ParameterSet& cfg)
00040     : thePairGenerator(0),
00041       theLayerCache(0),
00042       useFixedPreFiltering(cfg.getParameter<bool>("useFixedPreFiltering")),
00043       extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
00044       extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
00045       useMScat(cfg.getParameter<bool>("useMultScattering")),
00046       useBend(cfg.getParameter<bool>("useBending"))
00047 {
00048   theMaxElement=cfg.getParameter<unsigned int>("maxElement");
00049   if (useFixedPreFiltering)
00050     dphi = cfg.getParameter<double>("phiPreFiltering");
00051 }
00052 
00053 void PixelTripletLargeTipGenerator::init(const HitPairGenerator & pairs,
00054       const std::vector<SeedingLayer> &layers, LayerCacheType *layerCache)
00055 {
00056   thePairGenerator = pairs.clone();
00057   theLayers = layers;
00058   theLayerCache = layerCache;
00059 }
00060 
00061 static bool intersect(Range &range, const Range &second)
00062 {
00063   if (range.first >= second.max() || range.second <= second.min())
00064     return false;
00065   if (range.first < second.min())
00066     range.first = second.min();
00067   if (range.second > second.max())
00068     range.second = second.max();
00069   return range.first < range.second;
00070 }
00071 
00072 void PixelTripletLargeTipGenerator::hitTriplets( 
00073     const TrackingRegion& region, 
00074     OrderedHitTriplets & result,
00075     const edm::Event & ev,
00076     const edm::EventSetup& es)
00077 {
00078   edm::ESHandle<TrackerGeometry> tracker;
00079   es.get<TrackerDigiGeometryRecord>().get(tracker);
00080 
00081   OrderedHitPairs pairs;
00082   pairs.reserve(30000);
00083   thePairGenerator->hitPairs(region,pairs,ev,es);
00084   if (pairs.empty())
00085     return;
00086 
00087   int size = theLayers.size();
00088 
00089   map<const DetLayer*, LayerRZPredictions> mapPred;
00090   const RecHitsSortedInPhi **thirdHitMap = new const RecHitsSortedInPhi*[size];
00091   for(int il = 0; il < size; il++) {
00092      thirdHitMap[il] = &(*theLayerCache)(&theLayers[il], region, ev, es);
00093      const DetLayer *layer = theLayers[il].detLayer();
00094      LayerRZPredictions &predRZ = mapPred[layer];
00095      predRZ.line.initLayer(layer);
00096      predRZ.helix1.initLayer(layer);
00097      predRZ.helix2.initLayer(layer);
00098      predRZ.line.initTolerance(extraHitRZtolerance);
00099      predRZ.helix1.initTolerance(extraHitRZtolerance);
00100      predRZ.helix2.initTolerance(extraHitRZtolerance);
00101      predRZ.rzPositionFixup = MatchedHitRZCorrectionFromBending(layer);
00102   }
00103 
00104   double curv = PixelRecoUtilities::curvature(1. / region.ptMin(), es);
00105 
00106   for(OrderedHitPairs::const_iterator ip = pairs.begin();
00107       ip != pairs.end(); ++ip) {
00108 
00109     GlobalPoint gp1 = ip->inner()->globalPosition();
00110     GlobalPoint gp2 = ip->outer()->globalPosition();
00111 
00112     PixelRecoLineRZ line(gp1, gp2);
00113     PixelRecoPointRZ point2(gp2.perp(), gp2.z());
00114     ThirdHitPredictionFromCircle predictionRPhi(gp1, gp2, extraHitRPhitolerance);
00115 
00116     Range generalCurvature = predictionRPhi.curvature(region.originRBound());
00117     if (!intersect(generalCurvature, Range(-curv, curv)))
00118       continue;
00119 
00120     for(int il = 0; il < size; il++) {
00121       const DetLayer *layer = theLayers[il].detLayer();
00122 //      bool pixelLayer = layer->subDetector() == GeomDetEnumerators::PixelBarrel ||
00123 //                        layer->subDetector() == GeomDetEnumerators::PixelEndcap;
00124       bool barrelLayer = layer->location() == GeomDetEnumerators::barrel;
00125 
00126       Range curvature = generalCurvature;
00127       ThirdHitCorrection correction(es, region.ptMin(), layer, line, point2, useMScat);
00128 
00129       LayerRZPredictions &predRZ = mapPred.find(layer)->second;
00130       predRZ.line.initPropagator(&line);
00131 
00132       HelixRZ helix;
00133       Range rzRange;
00134       if (useBend) {
00135         // For the barrel region:
00136         // swiping the helix passing through the two points across from
00137         // negative to positive bending, can give us a sort of U-shaped
00138         // projection onto the phi-z (barrel) or r-z plane (forward)
00139         // so we checking minimum/maximum of all three possible extrema
00140         // 
00141         // For the endcap region:
00142         // Checking minimum/maximum radius of the helix projection
00143         // onto an endcap plane, here we have to guard against
00144         // looping tracks, when phi(delta z) gets out of control.
00145         // HelixRZ::rAtZ should not follow looping tracks, but clamp
00146         // to the minimum reachable r with the next-best lower |curvature|.
00147         // So same procedure as for the barrel region can be applied.
00148         //
00149         // In order to avoid looking for potential looping tracks at all
00150         // we also clamp the allowed curvature range for this layer,
00151         // and potentially fail the layer entirely
00152 
00153         if (!barrelLayer) {
00154           Range z3s = predRZ.line.detRange();
00155           double z3 = z3s.first < 0 ? max(z3s.first, z3s.second)
00156                                     : min(z3s.first, z3s.second);
00157           double maxCurvature = HelixRZ::maxCurvature(&predictionRPhi,
00158                                                       gp1.z(), gp2.z(), z3);
00159           if (!intersect(curvature, Range(-maxCurvature, maxCurvature)))
00160             continue;
00161         }
00162 
00163         helix = HelixRZ(&predictionRPhi, gp1.z(), gp2.z(), curvature.first);
00164         HelixRZ helix2(&predictionRPhi, gp1.z(), gp2.z(), curvature.second);
00165 
00166         predRZ.helix1.initPropagator(&helix);
00167         predRZ.helix2.initPropagator(&helix2);
00168 
00169         Range rzRanges[2] = { predRZ.helix1(), predRZ.helix2() };
00170         rzRange.first = min(rzRanges[0].first, rzRanges[1].first);
00171         rzRange.second = max(rzRanges[0].second, rzRanges[1].second);
00172 
00173         // if the allowed curvatures include a straight line,
00174         // this can give us another extremum for allowed r/z
00175         if (curvature.first * curvature.second < 0.0) {
00176           Range rzLineRange = predRZ.line();
00177           rzRange.first = min(rzRange.first, rzLineRange.first);
00178           rzRange.second = max(rzRange.second, rzLineRange.second);
00179         }
00180       } else
00181         rzRange = predRZ.line();
00182 
00183       if (rzRange.first >= rzRange.second)
00184         continue;
00185 
00186       correction.correctRZRange(rzRange);
00187 
00188       Range phiRange;
00189       if (useFixedPreFiltering) { 
00190         float phi0 = ip->outer()->globalPosition().phi();
00191         phiRange = Range(phi0 - dphi, phi0 + dphi);
00192       } else {
00193         Range radius;
00194 
00195         if (barrelLayer) {
00196           radius = predRZ.line.detRange();
00197           if (!intersect(rzRange, predRZ.line.detSize()))
00198             continue;
00199         } else {
00200           radius = rzRange;
00201           if (!intersect(radius, predRZ.line.detSize()))
00202             continue;
00203         }
00204 
00205         Range rPhi1 = predictionRPhi(curvature, radius.first);
00206         Range rPhi2 = predictionRPhi(curvature, radius.second);
00207         correction.correctRPhiRange(rPhi1);
00208         correction.correctRPhiRange(rPhi2);
00209         rPhi1.first  /= radius.first;
00210         rPhi1.second /= radius.first;
00211         rPhi2.first  /= radius.second;
00212         rPhi2.second /= radius.second;
00213         phiRange = mergePhiRanges(rPhi1, rPhi2);
00214       }
00215 
00216       typedef RecHitsSortedInPhi::Hit Hit;
00217       vector<Hit> thirdHits = thirdHitMap[il]->hits(phiRange.min(), phiRange.max());
00218 
00219       MatchedHitRZCorrectionFromBending l2rzFixup(ip->outer()->det()->geographicalId());
00220       MatchedHitRZCorrectionFromBending l3rzFixup = predRZ.rzPositionFixup;
00221 
00222       typedef vector<Hit>::const_iterator IH;
00223       for (IH th=thirdHits.begin(), eh=thirdHits.end(); th < eh; ++th) {
00224          GlobalPoint p3 = (*th)->globalPosition();
00225          double p3_r = p3.perp();
00226          double p3_z = p3.z();
00227          double p3_phi = p3.phi();
00228 
00229          Range rangeRPhi = predictionRPhi(curvature, p3_r);
00230          correction.correctRPhiRange(rangeRPhi);
00231 
00232          double phiErr = nSigmaPhi * (*th)->errorGlobalRPhi()/p3_r;
00233          if (!checkPhiInRange(p3_phi, rangeRPhi.first/p3_r - phiErr, rangeRPhi.second/p3_r + phiErr))
00234            continue;
00235 
00236          const TransientTrackingRecHit::ConstRecHitPointer& hit = *th;
00237          Basic2DVector<double> thc(p3.x(), p3.y());
00238 
00239          double curv_ = predictionRPhi.curvature(thc);
00240          double p2_r = point2.r(), p2_z = point2.z();
00241 
00242          l2rzFixup(predictionRPhi, curv_, *ip->outer(), p2_r, p2_z);
00243          l3rzFixup(predictionRPhi, curv_, **th, p3_r, p3_z);
00244 
00245          Range rangeRZ;
00246          if (useBend) {
00247            HelixRZ updatedHelix(&predictionRPhi, gp1.z(), p2_z, curv_);
00248            rangeRZ = predRZ.helix1(barrelLayer ? p3_r : p3_z, updatedHelix);
00249          } else {
00250            float tIP = predictionRPhi.transverseIP(thc);
00251            PixelRecoPointRZ updatedPoint2(p2_r, p2_z);
00252            PixelRecoLineRZ updatedLine(line.origin(), point2, tIP);
00253            rangeRZ = predRZ.line(barrelLayer ? p3_r : p3_z, line);
00254          }
00255          correction.correctRZRange(rangeRZ);
00256 
00257          double err = nSigmaRZ * (barrelLayer
00258                                 ? hit->errorGlobalZ()
00259                                 : hit->errorGlobalR());
00260          rangeRZ.first -= err, rangeRZ.second += err;
00261 
00262          if (!rangeRZ.inside(barrelLayer ? p3_z : p3_r))
00263            continue;
00264          if (theMaxElement!=0 && result.size() >= theMaxElement){
00265            result.clear();
00266            edm::LogError("TooManyTriplets")<<" number of triples exceed maximum. no triplets produced.";
00267            delete[] thirdHitMap;
00268            return;
00269          }
00270          result.push_back(OrderedHitTriplet(ip->inner(), ip->outer(), *th)); 
00271       }
00272     }
00273   }
00274   delete[] thirdHitMap;
00275 }
00276 
00277 bool PixelTripletLargeTipGenerator::checkPhiInRange(float phi, float phi1, float phi2) const
00278 {
00279   while (phi > phi2) phi -= 2. * M_PI;
00280   while (phi < phi1) phi += 2. * M_PI;
00281   return phi <= phi2;
00282 }  
00283 
00284 std::pair<float, float> PixelTripletLargeTipGenerator::mergePhiRanges(
00285           const std::pair<float, float> &r1, const std::pair<float, float> &r2) const
00286 {
00287   float r2Min = r2.first;
00288   float r2Max = r2.second;
00289   while (r1.first - r2Min > +M_PI) r2Min += 2. * M_PI, r2Max += 2. * M_PI;
00290   while (r1.first - r2Min < -M_PI) r2Min -= 2. * M_PI, r2Max -= 2. * M_PI;
00291   return std::make_pair(min(r1.first, r2Min), max(r1.second, r2Max));
00292 }