CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
PixelTripletHLTGenerator Class Reference

#include <PixelTripletHLTGenerator.h>

Inheritance diagram for PixelTripletHLTGenerator:
HitTripletGeneratorFromPairAndLayers HitTripletGenerator OrderedHitsGenerator

Public Member Functions

virtual void hitTriplets (const TrackingRegion &region, OrderedHitTriplets &trs, const edm::Event &ev, const edm::EventSetup &es)
 
virtual void init (const HitPairGenerator &pairs, const std::vector< ctfseeding::SeedingLayer > &layers, LayerCacheType *layerCache)
 
const HitPairGeneratorpairGenerator () const
 
 PixelTripletHLTGenerator (const edm::ParameterSet &cfg)
 
const std::vector
< ctfseeding::SeedingLayer > & 
thirdLayers () const
 
virtual ~PixelTripletHLTGenerator ()
 
- Public Member Functions inherited from HitTripletGeneratorFromPairAndLayers
virtual ~HitTripletGeneratorFromPairAndLayers ()
 
- Public Member Functions inherited from HitTripletGenerator
virtual void clear ()
 
 HitTripletGenerator (unsigned int size=500)
 
virtual void hitTriplets (const TrackingRegion &reg, OrderedHitTriplets &prs, const edm::EventSetup &es)
 
virtual const OrderedHitTripletsrun (const TrackingRegion &region, const edm::Event &ev, const edm::EventSetup &es)
 
virtual ~HitTripletGenerator ()
 
- Public Member Functions inherited from OrderedHitsGenerator
 OrderedHitsGenerator ()
 
virtual ~OrderedHitsGenerator ()
 

Private Types

typedef
CombinedHitTripletGenerator::LayerCacheType 
LayerCacheType
 

Private Member Functions

bool checkPhiInRange (float phi, float phi1, float phi2) const
 
std::pair< float, float > mergePhiRanges (const std::pair< float, float > &r1, const std::pair< float, float > &r2) const
 

Private Attributes

float dphi
 
float extraHitRPhitolerance
 
float extraHitRZtolerance
 
SeedComparitortheComparitor
 
LayerCacheTypetheLayerCache
 
std::vector
< ctfseeding::SeedingLayer
theLayers
 
HitPairGeneratorthePairGenerator
 
bool useBend
 
bool useFixedPreFiltering
 
bool useMScat
 

Additional Inherited Members

- Public Types inherited from HitTripletGeneratorFromPairAndLayers
typedef
CombinedHitTripletGenerator::LayerCacheType 
LayerCacheType
 
- Public Attributes inherited from OrderedHitsGenerator
unsigned int theMaxElement
 

Detailed Description

Definition at line 23 of file PixelTripletHLTGenerator.h.

Member Typedef Documentation

Definition at line 25 of file PixelTripletHLTGenerator.h.

Constructor & Destructor Documentation

PixelTripletHLTGenerator::PixelTripletHLTGenerator ( const edm::ParameterSet cfg)

Definition at line 26 of file PixelTripletHLTGenerator.cc.

References dphi, edm::ParameterSet::exists(), reco::get(), edm::ParameterSet::getParameter(), theComparitor, OrderedHitsGenerator::theMaxElement, and useFixedPreFiltering.

27  : thePairGenerator(0),
28  theLayerCache(0),
29  useFixedPreFiltering(cfg.getParameter<bool>("useFixedPreFiltering")),
30  extraHitRZtolerance(cfg.getParameter<double>("extraHitRZtolerance")),
31  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")),
32  useMScat(cfg.getParameter<bool>("useMultScattering")),
33  useBend(cfg.getParameter<bool>("useBending"))
34 {
35  theMaxElement=cfg.getParameter<unsigned int>("maxElement");
36  dphi = (useFixedPreFiltering) ? cfg.getParameter<double>("phiPreFiltering") : 0;
37 
38  if (cfg.exists("SeedComparitorPSet")){
39  edm::ParameterSet comparitorPSet =
40  cfg.getParameter<edm::ParameterSet>("SeedComparitorPSet");
41  std::string comparitorName = comparitorPSet.getParameter<std::string>("ComponentName");
42  theComparitor = (comparitorName == "none") ?
43  0 : SeedComparitorFactory::get()->create( comparitorName, comparitorPSet);
44  }
45  else
46  theComparitor=0;
47 }
T getParameter(std::string const &) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
T get(const Candidate &c)
Definition: component.h:56
PixelTripletHLTGenerator::~PixelTripletHLTGenerator ( )
virtual

Definition at line 49 of file PixelTripletHLTGenerator.cc.

References theComparitor, and thePairGenerator.

51 { delete thePairGenerator;
52  delete theComparitor;
53 }

Member Function Documentation

bool PixelTripletHLTGenerator::checkPhiInRange ( float  phi,
float  phi1,
float  phi2 
) const
private

Definition at line 223 of file PixelTripletHLTGenerator.cc.

References Geom::ftwoPi().

Referenced by hitTriplets().

224 {
225  while (phi > phi2) phi -= Geom::ftwoPi();
226  while (phi < phi1) phi += Geom::ftwoPi();
227  return ( (phi1 <= phi) && (phi <= phi2) );
228 }
float ftwoPi()
Definition: Pi.h:36
Definition: DDAxes.h:10
void PixelTripletHLTGenerator::hitTriplets ( const TrackingRegion region,
OrderedHitTriplets trs,
const edm::Event ev,
const edm::EventSetup es 
)
virtual

Implements HitTripletGenerator.

Definition at line 64 of file PixelTripletHLTGenerator.cc.

References GeomDetEnumerators::barrel, checkPhiInRange(), SeedComparitor::compatible(), PixelRecoUtilities::curvature(), dphi, extraHitRPhitolerance, extraHitRZtolerance, f, HitPairGenerator::hitPairs(), ThirdHitRZPredictionBase::initLayer(), ThirdHitRZPrediction< Propagator >::initPropagator(), PixelRecoRange< T >::intersection(), geometryCSVtoXML::line, DetLayer::location(), LogDebug, max(), PixelRecoRange< T >::max(), mergePhiRanges(), min, nSigmaPhi, nSigmaRZ, TrackingRegion::origin(), TrackingRegion::originRBound(), PV3DBase< T, PVType, FrameType >::perp(), point, TrackingRegion::ptMin(), CosmicsPD_Skims::radius, OrderedHitTriplets::size(), findQualityFiles::size, mathSSE::sqrt(), PixelRecoRange< T >::sum(), theComparitor, theLayers, OrderedHitsGenerator::theMaxElement, thePairGenerator, useBend, useFixedPreFiltering, useMScat, PV3DBase< T, PVType, FrameType >::x(), Hit::x, PV3DBase< T, PVType, FrameType >::y(), Hit::y, and PV3DBase< T, PVType, FrameType >::z().

69 {
70  OrderedHitPairs pairs; pairs.reserve(30000);
71  OrderedHitPairs::const_iterator ip;
72 
73  thePairGenerator->hitPairs(region,pairs,ev,es);
74 
75  if (pairs.empty()) return;
76 
77  int size = theLayers.size();
78 
79  typedef std::vector<ThirdHitRZPrediction<PixelRecoLineRZ> > Preds;
80  Preds preds(size);
81 
82  std::vector<const RecHitsSortedInPhi *> thirdHitMap(size);
84  vector<Hit> thirdHits;
85 
86  // fill the prediciton vetor
87  for (int il=0; il!=size; ++il) {
88  thirdHitMap[il] = &(*theLayerCache)(&theLayers[il], region, ev, es);
89  ThirdHitRZPrediction<PixelRecoLineRZ> & pred = preds[il];
90  pred.initLayer(theLayers[il].detLayer());
91  pred.initTolerance(extraHitRZtolerance);
92  }
93 
94 
95  double imppar = region.originRBound();
96  double curv = PixelRecoUtilities::curvature(1/region.ptMin(), es);
97 
98  for (ip = pairs.begin(); ip != pairs.end(); ip++) {
99 
100  GlobalPoint gp1tmp = (*ip).inner()->globalPosition();
101  GlobalPoint gp2tmp = (*ip).outer()->globalPosition();
102  GlobalPoint gp1(gp1tmp.x()-region.origin().x(), gp1tmp.y()-region.origin().y(), gp1tmp.z());
103  GlobalPoint gp2(gp2tmp.x()-region.origin().x(), gp2tmp.y()-region.origin().y(), gp2tmp.z());
104 
105  PixelRecoPointRZ point1(gp1.perp(), gp1.z());
106  PixelRecoPointRZ point2(gp2.perp(), gp2.z());
107  PixelRecoLineRZ line(point1, point2);
108  ThirdHitPredictionFromInvParabola predictionRPhi(gp1,gp2,imppar,curv,extraHitRPhitolerance);
109  ThirdHitPredictionFromInvParabola predictionRPhitmp(gp1tmp,gp2tmp,imppar+region.origin().perp(),curv,extraHitRPhitolerance);
110 
111 
112  for (int il=0; il!=size; ++il) {
113  const DetLayer * layer = theLayers[il].detLayer();
114 // bool pixelLayer = ( layer->subDetector() == GeomDetEnumerators::PixelBarrel
115 // || layer->subDetector() == GeomDetEnumerators::PixelEndcap);
116  bool barrelLayer = (layer->location() == GeomDetEnumerators::barrel);
117 
118  ThirdHitCorrection correction(es, region.ptMin(), layer, line, point2, useMScat, useBend);
119 
120  ThirdHitRZPrediction<PixelRecoLineRZ> & predictionRZ = preds[il];
121 
122  predictionRZ.initPropagator(&line);
123  Range rzRange = predictionRZ();
124 
125  correction.correctRZRange(rzRange);
126  Range phiRange;
127  if (useFixedPreFiltering) {
128  float phi0 = (*ip).outer()->globalPosition().phi();
129  phiRange = Range(phi0-dphi,phi0+dphi);
130  }
131  else {
132  Range radius;
133  if (barrelLayer) {
134  radius = predictionRZ.detRange();
135  } else {
136  radius = Range(
137  max(rzRange.min(), predictionRZ.detSize().min()),
138  min(rzRange.max(), predictionRZ.detSize().max()) );
139  }
140  if (radius.empty()) continue;
141  Range rPhi1m = predictionRPhitmp(radius.max(), -1);
142  Range rPhi1p = predictionRPhitmp(radius.max(), 1);
143  Range rPhi2m = predictionRPhitmp(radius.min(), -1);
144  Range rPhi2p = predictionRPhitmp(radius.min(), 1);
145  Range rPhi1 = rPhi1m.sum(rPhi1p);
146  Range rPhi2 = rPhi2m.sum(rPhi2p);
147  correction.correctRPhiRange(rPhi1);
148  correction.correctRPhiRange(rPhi2);
149  rPhi1.first /= radius.max();
150  rPhi1.second /= radius.max();
151  rPhi2.first /= radius.min();
152  rPhi2.second /= radius.min();
153  phiRange = mergePhiRanges(rPhi1,rPhi2);
154  }
155 
156 // LayerHitMapLoop thirdHits =
157 // pixelLayer ? thirdHitMap[il]->loop(phiRange, rzRange) :
158 // thirdHitMap[il]->loop();
159 
160  thirdHits.clear();
161  thirdHitMap[il]->hits(phiRange.min(),phiRange.max(), thirdHits);
162 
163  static float nSigmaRZ = std::sqrt(12.f);
164  static float nSigmaPhi = 3.f;
165 
166  typedef vector<Hit>::const_iterator IH;
167  for (IH th=thirdHits.begin(), eh=thirdHits.end(); th !=eh; ++th) {
168 
169  if (theMaxElement!=0 && result.size() >= theMaxElement){
170  result.clear();
171  edm::LogError("TooManyTriplets")<<" number of triples exceed maximum. no triplets produced.";
172  return;
173  }
174  const Hit& hit = (*th);
175  GlobalPoint point(hit->globalPosition().x()-region.origin().x(),
176  hit->globalPosition().y()-region.origin().y(),
177  hit->globalPosition().z() );
178  float p3_r = point.perp();
179  float p3_z = point.z();
180  float p3_phi = point.phi();
181 
182  if (barrelLayer) {
183  Range allowedZ = predictionRZ(p3_r);
184  correction.correctRZRange(allowedZ);
185 
186  float zErr = nSigmaRZ * std::sqrt(float(hit->globalPositionError().czz()));
187  Range hitRange(p3_z-zErr, p3_z+zErr);
188  Range crossingRange = allowedZ.intersection(hitRange);
189  if (crossingRange.empty()) continue;
190  } else {
191  Range allowedR = predictionRZ(p3_z);
192  correction.correctRZRange(allowedR);
193  float rErr = nSigmaRZ * std::sqrt(float(hit->globalPositionError().rerr( hit->globalPosition())));
194  Range hitRange(p3_r-rErr, p3_r+rErr);
195  Range crossingRange = allowedR.intersection(hitRange);
196  if (crossingRange.empty()) continue;
197  }
198 
199  float phiErr = nSigmaPhi*std::sqrt(float(hit->globalPositionError().phierr(hit->globalPosition())));
200  for (int icharge=-1; icharge <=1; icharge+=2) {
201  Range rangeRPhi = predictionRPhi(p3_r, icharge);
202  correction.correctRPhiRange(rangeRPhi);
203  if (checkPhiInRange(p3_phi, rangeRPhi.first/p3_r-phiErr, rangeRPhi.second/p3_r+phiErr)) {
204  // insert here check with comparitor
205  OrderedHitTriplet hittriplet( (*ip).inner(), (*ip).outer(), hit);
206  if(!theComparitor || theComparitor->compatible(hittriplet,es) ) {
207  result.push_back( hittriplet );
208  } else {
209  LogDebug("RejectedTriplet") << "rejected triplet from comparitor "
210  << hittriplet.outer()->globalPosition().x() << " "
211  << hittriplet.outer()->globalPosition().y() << " "
212  << hittriplet.outer()->globalPosition().z();
213  }
214  break;
215  }
216  }
217  }
218  }
219  }
220 
221 }
#define LogDebug(id)
std::vector< ctfseeding::SeedingLayer > theLayers
void initPropagator(const Propagator *propagator)
T perp() const
Definition: PV3DBase.h:66
T max() const
void initLayer(const DetLayer *layer)
virtual Location location() const =0
Which part of the detector (barrel, endcap)
static const double nSigmaPhi
PixelRecoRange< T > sum(const PixelRecoRange< T > &r) const
T y() const
Definition: PV3DBase.h:57
#define min(a, b)
Definition: mlp_lapack.h:161
std::pair< float, float > mergePhiRanges(const std::pair< float, float > &r1, const std::pair< float, float > &r2) const
T curvature(T InversePt, const edm::EventSetup &iSetup)
const T & max(const T &a, const T &b)
virtual void hitPairs(const TrackingRegion &reg, OrderedHitPairs &prs, const edm::EventSetup &es)
T sqrt(T t)
Definition: SSEVec.h:28
T z() const
Definition: PV3DBase.h:58
tuple result
Definition: query.py:137
double f[11][100]
static const double nSigmaRZ
PixelRecoRange< float > Range
PixelRecoRange< T > intersection(const PixelRecoRange< T > &r) const
bool checkPhiInRange(float phi, float phi1, float phi2) const
TransientTrackingRecHit::ConstRecHitPointer Hit
T x() const
Definition: PV3DBase.h:56
virtual bool compatible(const SeedingHitSet &hits, const edm::EventSetup &es)=0
tuple size
Write out results.
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
void PixelTripletHLTGenerator::init ( const HitPairGenerator pairs,
const std::vector< ctfseeding::SeedingLayer > &  layers,
LayerCacheType layerCache 
)
virtual

Implements HitTripletGeneratorFromPairAndLayers.

Definition at line 55 of file PixelTripletHLTGenerator.cc.

References HitPairGenerator::clone(), theLayerCache, theLayers, and thePairGenerator.

58 {
59  thePairGenerator = pairs.clone();
60  theLayers = layers;
61  theLayerCache = layerCache;
62 }
std::vector< ctfseeding::SeedingLayer > theLayers
virtual HitPairGenerator * clone() const =0
std::pair< float, float > PixelTripletHLTGenerator::mergePhiRanges ( const std::pair< float, float > &  r1,
const std::pair< float, float > &  r2 
) const
private

Definition at line 230 of file PixelTripletHLTGenerator.cc.

References Geom::fpi(), Geom::ftwoPi(), max(), and min.

Referenced by hitTriplets().

232 {
233  float r2_min=r2.first;
234  float r2_max=r2.second;
235  while (r1.first-r2_min > Geom::fpi()) { r2_min += Geom::ftwoPi(); r2_max += Geom::ftwoPi();}
236  while (r1.first-r2_min < -Geom::fpi()) { r2_min -= Geom::ftwoPi(); r2_max -= Geom::ftwoPi(); }
237 
238  return std::make_pair(min(r1.first,r2_min),max(r1.second,r2_max));
239 }
#define min(a, b)
Definition: mlp_lapack.h:161
float fpi()
Definition: Pi.h:35
const T & max(const T &a, const T &b)
float ftwoPi()
Definition: Pi.h:36
const HitPairGenerator& PixelTripletHLTGenerator::pairGenerator ( ) const
inline

Definition at line 38 of file PixelTripletHLTGenerator.h.

References thePairGenerator.

38 { return *thePairGenerator; }
const std::vector<ctfseeding::SeedingLayer>& PixelTripletHLTGenerator::thirdLayers ( ) const
inline

Definition at line 39 of file PixelTripletHLTGenerator.h.

References theLayers.

39 { return theLayers; }
std::vector< ctfseeding::SeedingLayer > theLayers

Member Data Documentation

float PixelTripletHLTGenerator::dphi
private

Definition at line 56 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets(), and PixelTripletHLTGenerator().

float PixelTripletHLTGenerator::extraHitRPhitolerance
private

Definition at line 53 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

float PixelTripletHLTGenerator::extraHitRZtolerance
private

Definition at line 52 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

SeedComparitor* PixelTripletHLTGenerator::theComparitor
private
LayerCacheType* PixelTripletHLTGenerator::theLayerCache
private

Definition at line 49 of file PixelTripletHLTGenerator.h.

Referenced by init().

std::vector<ctfseeding::SeedingLayer> PixelTripletHLTGenerator::theLayers
private

Definition at line 48 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets(), init(), and thirdLayers().

HitPairGenerator* PixelTripletHLTGenerator::thePairGenerator
private
bool PixelTripletHLTGenerator::useBend
private

Definition at line 55 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().

bool PixelTripletHLTGenerator::useFixedPreFiltering
private

Definition at line 51 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets(), and PixelTripletHLTGenerator().

bool PixelTripletHLTGenerator::useMScat
private

Definition at line 54 of file PixelTripletHLTGenerator.h.

Referenced by hitTriplets().