CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/RecoPixelVertexing/PixelTriplets/src/PixelTripletHLTGenerator.cc

Go to the documentation of this file.
00001 #include "RecoPixelVertexing/PixelTriplets/src/PixelTripletHLTGenerator.h"
00002 
00003 #include "RecoPixelVertexing/PixelTriplets/interface/ThirdHitPredictionFromInvParabola.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 "RecoPixelVertexing/PixelTriplets/src/ThirdHitCorrection.h"
00011 #include "RecoTracker/TkHitPairs/interface/RecHitsSortedInPhi.h"
00012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00013 #include <iostream>
00014 
00015 #include "DataFormats/GeometryVector/interface/Pi.h"
00016 
00017 using pixelrecoutilities::LongitudinalBendingCorrection;
00018 typedef PixelRecoRange<float> Range;
00019 
00020 using namespace std;
00021 using namespace ctfseeding;
00022 
00023 PixelTripletHLTGenerator:: PixelTripletHLTGenerator(const edm::ParameterSet& cfg)
00024     : thePairGenerator(0),
00025       theLayerCache(0),
00026       useFixedPreFiltering(cfg.getParameter<bool>("useFixedPreFiltering")),
00027       extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
00028       extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
00029       useMScat(cfg.getParameter<bool>("useMultScattering")),
00030       useBend(cfg.getParameter<bool>("useBending"))
00031 {
00032   theMaxElement=cfg.getParameter<unsigned int>("maxElement");
00033   dphi =  (useFixedPreFiltering) ?  cfg.getParameter<double>("phiPreFiltering") : 0;
00034 }
00035 
00036 
00037 void PixelTripletHLTGenerator::init( const HitPairGenerator & pairs,
00038       const std::vector<SeedingLayer> & layers,
00039       LayerCacheType* layerCache)
00040 {
00041   thePairGenerator = pairs.clone();
00042   theLayers = layers;
00043   theLayerCache = layerCache;
00044 }
00045 
00046 void PixelTripletHLTGenerator::hitTriplets( 
00047     const TrackingRegion& region, 
00048     OrderedHitTriplets & result,
00049     const edm::Event & ev,
00050     const edm::EventSetup& es)
00051 {
00052   OrderedHitPairs pairs; pairs.reserve(30000);
00053   OrderedHitPairs::const_iterator ip;
00054   
00055   thePairGenerator->hitPairs(region,pairs,ev,es);
00056 
00057   if (pairs.empty()) return;
00058 
00059   int size = theLayers.size();
00060 
00061   typedef std::vector<ThirdHitRZPrediction<PixelRecoLineRZ> >  Preds;
00062   Preds preds(size);
00063 
00064   std::vector<const RecHitsSortedInPhi *> thirdHitMap(size);
00065   typedef RecHitsSortedInPhi::Hit Hit;
00066   vector<Hit> thirdHits;
00067 
00068   // fill the prediciton vetor
00069   for (int il=0; il!=size; ++il) {
00070      thirdHitMap[il] = &(*theLayerCache)(&theLayers[il], region, ev, es);
00071      ThirdHitRZPrediction<PixelRecoLineRZ> & pred = preds[il];
00072      pred.initLayer(theLayers[il].detLayer());
00073      pred.initTolerance(extraHitRZtolerance);
00074   }
00075 
00076 
00077   double imppar = region.originRBound();
00078   double curv = PixelRecoUtilities::curvature(1/region.ptMin(), es);
00079 
00080   for (ip = pairs.begin(); ip != pairs.end(); ip++) {
00081   
00082     GlobalPoint gp1tmp = (*ip).inner()->globalPosition();
00083     GlobalPoint gp2tmp = (*ip).outer()->globalPosition();
00084     GlobalPoint gp1(gp1tmp.x()-region.origin().x(), gp1tmp.y()-region.origin().y(), gp1tmp.z());
00085     GlobalPoint gp2(gp2tmp.x()-region.origin().x(), gp2tmp.y()-region.origin().y(), gp2tmp.z());
00086 
00087     PixelRecoPointRZ point1(gp1.perp(), gp1.z());
00088     PixelRecoPointRZ point2(gp2.perp(), gp2.z());
00089     PixelRecoLineRZ  line(point1, point2);
00090     ThirdHitPredictionFromInvParabola predictionRPhi(gp1,gp2,imppar,curv,extraHitRPhitolerance);
00091     ThirdHitPredictionFromInvParabola predictionRPhitmp(gp1tmp,gp2tmp,imppar+region.origin().perp(),curv,extraHitRPhitolerance);
00092 
00093 
00094     for (int il=0; il!=size; ++il) {
00095       const DetLayer * layer = theLayers[il].detLayer();
00096 //      bool pixelLayer = (    layer->subDetector() == GeomDetEnumerators::PixelBarrel 
00097 //                          || layer->subDetector() == GeomDetEnumerators::PixelEndcap); 
00098       bool barrelLayer = (layer->location() == GeomDetEnumerators::barrel);
00099 
00100       ThirdHitCorrection correction(es, region.ptMin(), layer, line, point2, useMScat, useBend); 
00101       
00102       ThirdHitRZPrediction<PixelRecoLineRZ> & predictionRZ =  preds[il];
00103 
00104       predictionRZ.initPropagator(&line);
00105       Range rzRange = predictionRZ();
00106 
00107       correction.correctRZRange(rzRange);
00108       Range phiRange;
00109       if (useFixedPreFiltering) { 
00110         float phi0 = (*ip).outer()->globalPosition().phi();
00111         phiRange = Range(phi0-dphi,phi0+dphi);
00112       }
00113       else {
00114         Range radius;
00115         if (barrelLayer) {
00116           radius =  predictionRZ.detRange();
00117         } else {
00118           radius = Range(
00119               max(rzRange.min(), predictionRZ.detSize().min()),
00120               min(rzRange.max(), predictionRZ.detSize().max()) );
00121         }
00122         if (radius.empty()) continue;
00123         Range rPhi1m = predictionRPhitmp(radius.max(), -1);
00124         Range rPhi1p = predictionRPhitmp(radius.max(),  1);
00125         Range rPhi2m = predictionRPhitmp(radius.min(), -1);
00126         Range rPhi2p = predictionRPhitmp(radius.min(),  1);
00127         Range rPhi1 = rPhi1m.sum(rPhi1p);
00128         Range rPhi2 = rPhi2m.sum(rPhi2p);
00129         correction.correctRPhiRange(rPhi1);
00130         correction.correctRPhiRange(rPhi2);
00131         rPhi1.first  /= radius.max();
00132         rPhi1.second /= radius.max();
00133         rPhi2.first  /= radius.min();
00134         rPhi2.second /= radius.min();
00135         phiRange = mergePhiRanges(rPhi1,rPhi2);
00136       }
00137       
00138 //      LayerHitMapLoop thirdHits = 
00139 //          pixelLayer ? thirdHitMap[il]->loop(phiRange, rzRange) : 
00140 //          thirdHitMap[il]->loop();
00141 
00142       thirdHits.clear();
00143       thirdHitMap[il]->hits(phiRange.min(),phiRange.max(), thirdHits);
00144   
00145       static float nSigmaRZ = std::sqrt(12.f);
00146       static float nSigmaPhi = 3.f;
00147    
00148       typedef vector<Hit>::const_iterator IH;
00149       for (IH th=thirdHits.begin(), eh=thirdHits.end(); th !=eh; ++th) {
00150 
00151         if (theMaxElement!=0 && result.size() >= theMaxElement){
00152           result.clear();
00153           edm::LogError("TooManyTriplets")<<" number of triples exceed maximum. no triplets produced.";
00154           return;
00155         }
00156         const Hit& hit = (*th);
00157         GlobalPoint point(hit->globalPosition().x()-region.origin().x(),
00158                           hit->globalPosition().y()-region.origin().y(),
00159                           hit->globalPosition().z() ); 
00160         float p3_r = point.perp();
00161         float p3_z = point.z();
00162         float p3_phi = point.phi();
00163  
00164         if (barrelLayer) {
00165           Range allowedZ = predictionRZ(p3_r);
00166           correction.correctRZRange(allowedZ);
00167 
00168           float zErr = nSigmaRZ * std::sqrt(float(hit->globalPositionError().czz()));
00169           Range hitRange(p3_z-zErr, p3_z+zErr);
00170           Range crossingRange = allowedZ.intersection(hitRange);
00171           if (crossingRange.empty())  continue;
00172         } else {
00173           Range allowedR = predictionRZ(p3_z);
00174           correction.correctRZRange(allowedR); 
00175           float rErr = nSigmaRZ * std::sqrt(float(hit->globalPositionError().rerr( hit->globalPosition())));
00176           Range hitRange(p3_r-rErr, p3_r+rErr);
00177           Range crossingRange = allowedR.intersection(hitRange);
00178           if (crossingRange.empty())  continue;
00179         }
00180 
00181         float phiErr = nSigmaPhi*std::sqrt(float(hit->globalPositionError().phierr(hit->globalPosition())));
00182         for (int icharge=-1; icharge <=1; icharge+=2) {
00183           Range rangeRPhi = predictionRPhi(p3_r, icharge);
00184           correction.correctRPhiRange(rangeRPhi);
00185           if (checkPhiInRange(p3_phi, rangeRPhi.first/p3_r-phiErr, rangeRPhi.second/p3_r+phiErr)) {
00186             result.push_back( OrderedHitTriplet( (*ip).inner(), (*ip).outer(), hit)); 
00187             break;
00188           } 
00189         }
00190       } 
00191     }
00192   }
00193 
00194 }
00195 
00196 bool PixelTripletHLTGenerator::checkPhiInRange(float phi, float phi1, float phi2) const
00197 {
00198   while (phi > phi2) phi -=  Geom::ftwoPi();
00199   while (phi < phi1) phi +=  Geom::ftwoPi();
00200   return (  (phi1 <= phi) && (phi <= phi2) );
00201 }  
00202 
00203 std::pair<float,float> PixelTripletHLTGenerator::mergePhiRanges(
00204     const std::pair<float,float>& r1, const std::pair<float,float>& r2) const 
00205 {
00206   float r2_min=r2.first;
00207   float r2_max=r2.second;
00208   while (r1.first-r2_min > Geom::fpi()) { r2_min += Geom::ftwoPi(); r2_max += Geom::ftwoPi();}
00209   while (r1.first-r2_min < -Geom::fpi()) { r2_min -= Geom::ftwoPi();  r2_max -= Geom::ftwoPi(); }
00210   
00211   return std::make_pair(min(r1.first,r2_min),max(r1.second,r2_max));
00212 }