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