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