CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/RecoPixelVertexing/PixelTriplets/plugins/PixelTripletNoTipGenerator.cc

Go to the documentation of this file.
00001 #include "PixelTripletNoTipGenerator.h"
00002 #include "FWCore/Framework/interface/EventSetup.h"
00003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00004 
00005 #include "RecoTracker/TkHitPairs/interface/RecHitsSortedInPhi.h"
00006 #include "FWCore/Framework/interface/Event.h"
00007 #include "RecoTracker/TkMSParametrization/interface/PixelRecoUtilities.h"
00008 #include "RecoPixelVertexing/PixelTriplets/interface/ThirdHitPredictionFromInvLine.h"
00009 #include "ThirdHitZPrediction.h"
00010 
00011 #include "RecoTracker/TkMSParametrization/interface/MultipleScatteringParametrisation.h"
00012 #include "RecoTracker/TkHitPairs/interface/HitPairGeneratorFromLayerPair.h"
00013 
00014 #include "TrackingTools/DetLayers/interface/BarrelDetLayer.h"
00015 #include "TrackingTools/DetLayers/interface/ForwardDetLayer.h"
00016 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00017 #include "DataFormats/GeometrySurface/interface/SimpleDiskBounds.h"
00018 
00019 typedef PixelRecoRange<float> Range;
00020 template<class T> T sqr(T t) { return t * t;}
00021 
00022 using namespace std;
00023 using namespace ctfseeding;
00024 
00025 PixelTripletNoTipGenerator:: PixelTripletNoTipGenerator(const edm::ParameterSet& cfg)
00026     : thePairGenerator(0),
00027       theLayerCache(0),
00028       extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
00029       extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
00030       extraHitPhiToleranceForPreFiltering(cfg.getParameter<double>("extraHitPhiToleranceForPreFiltering")), 
00031       theNSigma(cfg.getParameter<double>("nSigma")),
00032       theChi2Cut(cfg.getParameter<double>("chi2Cut"))
00033 { }
00034 
00035 void PixelTripletNoTipGenerator::init( const HitPairGenerator & pairs,
00036       const std::vector<SeedingLayer> & layers,
00037       LayerCacheType* layerCache)
00038 {
00039   thePairGenerator = pairs.clone();
00040   theLayers = layers;
00041   theLayerCache = layerCache;
00042 }
00043 
00044 void PixelTripletNoTipGenerator::hitTriplets(
00045     const TrackingRegion& region,
00046     OrderedHitTriplets & result,
00047     const edm::Event & ev,
00048     const edm::EventSetup& es)
00049 {
00050 
00051 //
00052 //edm::Handle<reco::BeamSpot> bsHandle;
00053 //ev.getByLabel( theBeamSpotTag, bsHandle);
00054 //if(!bsHandle.isValid()) return;
00055 //const reco::BeamSpot & bs = *bsHandle;
00056 //double errorXY = sqrt( sqr(bs.BeamWidthX()) + sqr(bs.BeamWidthY()) );
00057 //
00058 
00059   GlobalPoint bsPoint = region.origin();
00060   double errorXY = region.originRBound();
00061   GlobalVector shift =   bsPoint - GlobalPoint(0.,0.,0.);
00062 
00063   OrderedHitPairs pairs; pairs.reserve(30000);
00064   OrderedHitPairs::const_iterator ip;
00065   thePairGenerator->hitPairs(region,pairs,ev,es);
00066 
00067   if (pairs.size() ==0) return;
00068 
00069   int size = theLayers.size();
00070 
00071   const RecHitsSortedInPhi **thirdHitMap = new const RecHitsSortedInPhi*[size];
00072   for (int il=0; il <=size-1; il++) {
00073      thirdHitMap[il] = &(*theLayerCache)(&theLayers[il], region, ev, es);
00074   }
00075 
00076   const HitPairGeneratorFromLayerPair * pairGen = dynamic_cast<const HitPairGeneratorFromLayerPair *>(thePairGenerator);
00077   const DetLayer * firstLayer = pairGen->innerLayer().detLayer();
00078   const DetLayer * secondLayer = pairGen->outerLayer().detLayer();
00079   if (!firstLayer || !secondLayer) return;
00080 
00081   MultipleScatteringParametrisation sigma1RPhi( firstLayer, es);
00082   MultipleScatteringParametrisation sigma2RPhi( secondLayer, es);
00083 
00084   typedef RecHitsSortedInPhi::Hit Hit;
00085   for (ip = pairs.begin(); ip != pairs.end(); ip++) {
00086 
00087     GlobalPoint p1((*ip).inner()->globalPosition()-shift);
00088     GlobalPoint p2((*ip).outer()->globalPosition()-shift);
00089 
00090     ThirdHitPredictionFromInvLine  predictionRPhiTMP(p1, p2 );
00091     double pt_p1p2 = 1./PixelRecoUtilities::inversePt(predictionRPhiTMP.curvature(),es);
00092 
00093     PixelRecoPointRZ point1(p1.perp(), p1.z());
00094     PixelRecoPointRZ point2(p2.perp(), p2.z());
00095 
00096     PixelRecoLineRZ  line(point1, point2);
00097     double msRPhi1 = sigma1RPhi(pt_p1p2, line.cotLine());
00098     double msRPhi2 = sigma2RPhi(pt_p1p2,  line.cotLine(),point1);
00099     double sinTheta = 1/sqrt(1+sqr(line.cotLine()));
00100     double cosTheta = fabs(line.cotLine())/sqrt(1+sqr(line.cotLine()));
00101 
00102     double p1_errorRPhi = sqrt(sqr((*ip).inner()->errorGlobalRPhi())+sqr(msRPhi1) +sqr(errorXY));
00103     double p2_errorRPhi = sqrt(sqr((*ip).outer()->errorGlobalRPhi())+sqr(msRPhi2) +sqr(errorXY));
00104 
00105     ThirdHitPredictionFromInvLine  predictionRPhi(p1, p2, p1_errorRPhi, p2_errorRPhi );
00106 
00107     for (int il=0; il <=size-1; il++) {
00108 
00109       const DetLayer * layer = theLayers[il].detLayer();
00110       bool barrelLayer = (layer->location() == GeomDetEnumerators::barrel);
00111       MultipleScatteringParametrisation sigma3RPhi( layer, es);
00112       double msRPhi3 = sigma3RPhi(pt_p1p2, line.cotLine(),point2);
00113 
00114       Range rRange;
00115       if (barrelLayer) {
00116          const BarrelDetLayer& bl = dynamic_cast<const BarrelDetLayer&>(*layer);
00117          float halfThickness  = bl.surface().bounds().thickness()/2;
00118          float radius = bl.specificSurface().radius();
00119          rRange = Range(radius-halfThickness, radius+halfThickness);
00120       } else {
00121         const ForwardDetLayer& fl = dynamic_cast<const ForwardDetLayer&>(*layer);
00122         float halfThickness  = fl.surface().bounds().thickness()/2;
00123         float zLayer = fl.position().z() ;
00124         float zMin = zLayer-halfThickness;
00125         float zMax = zLayer+halfThickness;
00126         GlobalVector dr = p2-p1;
00127         GlobalPoint p3_a = p2+dr*(zMin-p2.z())/dr.z();
00128         GlobalPoint p3_b = p2+dr*(zMax-p2.z())/dr.z();
00129         if (zLayer * p3_a.z() < 0) continue;
00130         rRange = Range(p3_a.perp(), p3_b.perp());
00131         rRange.sort();
00132       }
00133       double displacment = shift.perp();
00134       GlobalPoint crossing1 = predictionRPhi.crossing(rRange.min()-displacment)+shift;
00135       GlobalPoint crossing2 = predictionRPhi.crossing(rRange.max()+displacment)+shift;
00136       float c1_phi= crossing1.phi(); 
00137       float c2_phi= crossing2.phi(); 
00138       if (c2_phi < c1_phi) swap(c1_phi,c2_phi); 
00139       if (c2_phi-c1_phi > M_PI) { c2_phi -= 2*M_PI;  swap(c1_phi,c2_phi); }
00140       double extraAngle = (displacment+theNSigma*msRPhi3)/rRange.min()+extraHitPhiToleranceForPreFiltering;
00141       c1_phi -= extraAngle; 
00142       c2_phi += extraAngle;
00143       vector<Hit> thirdHits = thirdHitMap[il]->hits(c1_phi, c2_phi) ;
00144 
00145       typedef vector<Hit>::const_iterator IH;
00146       for (IH th=thirdHits.begin(), eh=thirdHits.end(); th < eh; ++th) {
00147         GlobalPoint p3((*th)->globalPosition()-shift);
00148         double p3_errorRPhi = sqrt(sqr((*th)->errorGlobalRPhi()) +sqr(msRPhi3) + sqr(errorXY));
00149 
00150         predictionRPhi.add(p3,p3_errorRPhi);
00151 
00152         double curvature = predictionRPhi.curvature();
00153            
00154         ThirdHitZPrediction zPrediction( (*ip).inner()->globalPosition(), sqrt(sqr((*ip).inner()->errorGlobalR())+sqr(msRPhi1/cosTheta)), sqrt( sqr((*ip).inner()->errorGlobalZ())+ sqr(msRPhi1/sinTheta)), 
00155                                          (*ip).outer()->globalPosition(), sqrt(sqr((*ip).outer()->errorGlobalR())+sqr(msRPhi2/cosTheta)), sqrt( sqr((*ip).outer()->errorGlobalZ())+sqr(msRPhi2/sinTheta)), 
00156                                          curvature, theNSigma);
00157          ThirdHitZPrediction::Range zRange = zPrediction((*th)->globalPosition(), sqrt(sqr((*th)->errorGlobalR()))+sqr(msRPhi3/cosTheta));
00158          
00159          double z3Hit = (*th)->globalPosition().z(); 
00160          double z3HitError = theNSigma*(sqrt(sqr((*th)->errorGlobalZ()) + sqr(msRPhi3/sinTheta) ))+extraHitRZtolerance; 
00161          Range hitZRange(z3Hit-z3HitError, z3Hit+z3HitError); 
00162          bool inside = hitZRange.hasIntersection(zRange); 
00163 
00164          double curvatureMS = PixelRecoUtilities::curvature(1./region.ptMin(),es);
00165          bool ptCut = (predictionRPhi.curvature()-theNSigma*predictionRPhi.errorCurvature() < curvatureMS); 
00166          bool chi2Cut = (predictionRPhi.chi2() < theChi2Cut);
00167          if (inside && ptCut && chi2Cut) { 
00168            if (theMaxElement!=0 && result.size() >= theMaxElement){
00169              result.clear();
00170              edm::LogError("TooManyTriplets")<<" number of triples exceed maximum. no triplets produced.";
00171              delete [] thirdHitMap;
00172              return;
00173            }
00174            result.push_back( OrderedHitTriplet( (*ip).inner(), (*ip).outer(), *th)); 
00175          }
00176          predictionRPhi.remove(p3,p3_errorRPhi);
00177       } 
00178     }
00179   }
00180   delete [] thirdHitMap;
00181 }