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