CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
AreaSeededTrackingRegionsBuilder::Builder Class Reference

#include <AreaSeededTrackingRegionsBuilder.h>

Public Member Functions

 Builder (const AreaSeededTrackingRegionsBuilder *conf)
 
std::unique_ptr< TrackingRegionregion (const Origin &origin, const std::vector< Area > &areas) const
 
std::unique_ptr< TrackingRegionregion (const Origin &origin, const edm::VecArray< Area, 2 > &areas) const
 
std::vector< std::unique_ptr< TrackingRegion > > regions (const Origins &origins, const std::vector< Area > &areas) const
 
void setCandidates (const TrackingSeedCandidates::Objects cands)
 
void setMeasurementTracker (const MeasurementTrackerEvent *mte)
 
 ~Builder ()=default
 

Private Member Functions

template<typename T >
std::unique_ptr< TrackingRegionregionImpl (const Origin &origin, const T &areas) const
 

Private Attributes

TrackingSeedCandidates::Objects candidates
 
const AreaSeededTrackingRegionsBuilderm_conf = 0
 
const MeasurementTrackerEventm_measurementTracker = 0
 

Detailed Description

Definition at line 76 of file AreaSeededTrackingRegionsBuilder.h.

Constructor & Destructor Documentation

AreaSeededTrackingRegionsBuilder::Builder::Builder ( const AreaSeededTrackingRegionsBuilder conf)
inlineexplicit

Definition at line 78 of file AreaSeededTrackingRegionsBuilder.h.

78 : m_conf(conf) {}
const AreaSeededTrackingRegionsBuilder * m_conf
AreaSeededTrackingRegionsBuilder::Builder::~Builder ( )
default

Member Function Documentation

std::unique_ptr< TrackingRegion > AreaSeededTrackingRegionsBuilder::Builder::region ( const Origin origin,
const std::vector< Area > &  areas 
) const

Definition at line 79 of file AreaSeededTrackingRegionsBuilder.cc.

79  {
80  return regionImpl(origin, areas);
81 }
std::unique_ptr< TrackingRegion > regionImpl(const Origin &origin, const T &areas) const
std::unique_ptr< TrackingRegion > AreaSeededTrackingRegionsBuilder::Builder::region ( const Origin origin,
const edm::VecArray< Area, 2 > &  areas 
) const

Definition at line 82 of file AreaSeededTrackingRegionsBuilder.cc.

82  {
83  return regionImpl(origin, areas);
84 }
std::unique_ptr< TrackingRegion > regionImpl(const Origin &origin, const T &areas) const
template<typename T >
std::unique_ptr< TrackingRegion > AreaSeededTrackingRegionsBuilder::Builder::regionImpl ( const Origin origin,
const T areas 
) const
private

Definition at line 87 of file AreaSeededTrackingRegionsBuilder.cc.

References funct::abs(), jets_cff::area, objects.IsoTrackAnalyzer::candidates, funct::cos(), edmIntegrityCheck::d, reco::deltaPhi(), hiPixelPairStep_cff::deltaPhi, JetChargeProducer_cfi::exp, f, LogDebug, LogTrace, M_PI, SiStripPI::max, maxEta, trackingParticleSelector_cfi::maxPhi, min(), heepElectronID_HEEPV50_cff::minEta, trackingParticleSelector_cfi::minPhi, normalizedPhi(), perp2(), phimax, phimin, alignCSCRings::r, Validation_hcalonly_cfi::sign, funct::sin(), mathSSE::sqrt(), funct::tan(), tmp, x, y, and z.

87  {
88  float minEta=std::numeric_limits<float>::max(), maxEta=std::numeric_limits<float>::lowest();
89  float minPhi=std::numeric_limits<float>::max(), maxPhi=std::numeric_limits<float>::lowest();
90 
91  const auto& orig = origin.first;
92 
93  LogDebug("AreaSeededTrackingRegionsProducer") << "Origin x,y,z " << orig.x() << "," << orig.y() << "," << orig.z();
94 
95  auto vecFromOrigPlusRadius = [&](const std::array<float, 2>& vec) {
96  const auto invlen = 1.f/std::sqrt(perp2(vec));
97  const auto tmp = (1.f - m_conf->m_originRadius*invlen);
98  return std::array<float, 2>{{vec[0]*tmp, vec[1]*tmp}};
99  };
100  auto tangentVec = [&](const std::array<float, 2>& vec, int sign) {
101  const auto d = std::sqrt(perp2(vec));
102  const auto r = m_conf->m_originRadius;
103  const auto tmp = r/std::sqrt(d*d - r*r);
104  // sign+ for counterclockwise, sign- for clockwise
105  const auto orthvec = sign > 0 ? std::array<float, 2>{{-vec[1]*tmp, vec[0]*tmp}}
106  : std::array<float, 2>{{ vec[1]*tmp, -vec[0]*tmp}};
107  return std::array<float, 2>{{vec[0]-orthvec[0], vec[1]-orthvec[1]}};
108  };
109 
110  for(const auto& area: areas) {
111  // straight line assumption is conservative, accounding for
112  // low-pT bending would only tighten the eta-phi window
113  LogTrace("AreaSeededTrackingRegionsBuilder") << " area x,y points "
114  << area.x_rmin_phimin() << "," << area.y_rmin_phimin() << " " // go around
115  << area.x_rmin_phimax() << "," << area.y_rmin_phimax() << " "
116  << area.x_rmax_phimax() << "," << area.y_rmax_phimax() << " "
117  << area.x_rmax_phimin() << "," << area.y_rmax_phimin() << " "
118  << "z " << area.zmin() << "," << area.zmax();
119 
120  // some common variables
121  const float x_rmin_phimin = area.x_rmin_phimin()-orig.x();
122  const float x_rmin_phimax = area.x_rmin_phimax()-orig.x();
123  const float y_rmin_phimin = area.y_rmin_phimin()-orig.y();
124  const float y_rmin_phimax = area.y_rmin_phimax()-orig.y();
125 
126  const std::array<float, 2> p_rmin_phimin = {{x_rmin_phimin, y_rmin_phimin}};
127  const std::array<float, 2> p_rmin_phimax = {{x_rmin_phimax, y_rmin_phimax}};
128  const std::array<float, 2> p_rmax_phimin = {{area.x_rmax_phimin() - orig.x(), area.y_rmax_phimin() - orig.y()}};
129  const std::array<float, 2> p_rmax_phimax = {{area.x_rmax_phimax() - orig.x(), area.y_rmax_phimax() - orig.y()}};
130 
131 
132  // eta
133  {
134  // find which of the two rmin points is closer to the origin
135  //
136  // closest point for p_rmin, farthest point for p_rmax
137  const std::array<float, 2> p_rmin = perp2(p_rmin_phimin) < perp2(p_rmin_phimax) ? p_rmin_phimin : p_rmin_phimax;
138  const std::array<float, 2> p_rmax = perp2(p_rmax_phimin) > perp2(p_rmax_phimax) ? p_rmax_phimin : p_rmax_phimax;
139 
140  // then calculate the xy vector from the point closest to p_rmin on
141  // the (origin,originRadius) circle to the p_rmin
142  // this will maximize the eta window
143  const auto p_min = vecFromOrigPlusRadius(p_rmin);
144  const auto p_max = vecFromOrigPlusRadius(p_rmax);
145 
146  // then we calculate the maximal eta window
147  const auto etamin = std::min(etaFromXYZ(p_min[0], p_min[1], area.zmin() - (orig.z()+origin.second)),
148  etaFromXYZ(p_max[0], p_max[1], area.zmin() - (orig.z()+origin.second)));
149  const auto etamax = std::max(etaFromXYZ(p_min[0], p_min[1], area.zmax() - (orig.z()-origin.second)),
150  etaFromXYZ(p_max[0], p_max[1], area.zmax() - (orig.z()-origin.second)));
151 
152  LogTrace("AreaSeededTrackingRegionBuilder") << " eta min,max " << etamin << "," << etamax;
153 
154  minEta = std::min(minEta, etamin);
155  maxEta = std::max(maxEta, etamax);
156  }
157 
158  // phi
159  {
160  // For phi we construct the tangent lines of (origin,
161  // originRadius) that go though each of the 4 points (in xy
162  // plane) of the area. Easiest is to make a vector orthogonal to
163  // origin->point vector which has a length of
164  //
165  // length = r*d/std::sqrt(d*d-r*r)
166  //
167  // where r is the originRadius and d is the distance from origin
168  // to the point (to derive draw the situation and start with
169  // definitions of sin/cos of one of the angles of the
170  // right-angled triangle.
171 
172  // But we start with a "reference phi" so that we can easily
173  // decide which phi is the largest/smallest without having too
174  // much of headache with the wrapping. The reference is in
175  // principle defined by the averages of y&x phimin/phimax at
176  // rmin, but the '0.5f*' factor is omitted here to reduce
177  // computations.
178  const auto phi_ref = std::atan2(y_rmin_phimin + y_rmin_phimax,
179  x_rmin_phimin + x_rmin_phimax);
180 
181  // for maximum phi we need the orthogonal vector to the left
182  const auto tan_rmin_phimax = tangentVec(p_rmin_phimax, +1);
183  const auto tan_rmax_phimax = tangentVec(p_rmax_phimax, +1);
184  const auto phi_rmin_phimax = std::atan2(tan_rmin_phimax[1], tan_rmin_phimax[0]);
185  const auto phi_rmax_phimax = std::atan2(tan_rmax_phimax[1], tan_rmax_phimax[0]);
186 
187  auto phimax = std::abs(reco::deltaPhi(phi_rmin_phimax, phi_ref)) > std::abs(reco::deltaPhi(phi_rmax_phimax, phi_ref)) ?
188  phi_rmin_phimax : phi_rmax_phimax;
189 
190  LogTrace("AreaSeededTrackingRegionBuilder") << " rmin_phimax vec " << p_rmin_phimax[0] << "," << p_rmin_phimax[1]
191  << " tangent " << tan_rmin_phimax[0] << "," << tan_rmin_phimax[1]
192  << " phi " << phi_rmin_phimax << "\n"
193  << " rmax_phimax vec " << p_rmax_phimax[0] << "," << p_rmax_phimax[1]
194  << " tangent " << tan_rmax_phimax[0] << "," << tan_rmax_phimax[1]
195  << " phi " << phi_rmax_phimax << "\n"
196  << " phimax " << phimax;
197 
198 
199 
200  // for minimum phi we need the orthogonal vector to the right
201  const auto tan_rmin_phimin = tangentVec(p_rmin_phimin, -1);
202  const auto tan_rmax_phimin = tangentVec(p_rmax_phimin, -1);
203  const auto phi_rmin_phimin = std::atan2(tan_rmin_phimin[1], tan_rmin_phimin[0]);
204  const auto phi_rmax_phimin = std::atan2(tan_rmax_phimin[1], tan_rmax_phimin[0]);
205 
206  auto phimin = std::abs(reco::deltaPhi(phi_rmin_phimin, phi_ref)) > std::abs(reco::deltaPhi(phi_rmax_phimin, phi_ref)) ?
207  phi_rmin_phimin : phi_rmax_phimin;
208 
209 
210  LogTrace("AreaSeededTrackingRegionBuilder") << " rmin_phimin vec " << p_rmin_phimin[0] << "," << p_rmin_phimin[1]
211  << " tangent " << tan_rmin_phimin[0] << "," << tan_rmin_phimin[1]
212  << " phi " << phi_rmin_phimin << "\n"
213  << " rmax_phimin vec " << p_rmax_phimin[0] << "," << p_rmax_phimin[1]
214  << " tangent " << tan_rmax_phimin[0] << "," << tan_rmax_phimin[1]
215  << " phi " << phi_rmax_phimin << "\n"
216  << " phimin " << phimin;
217 
218  // wrapped around, need to decide which one to wrap
219  if(phimax < phimin) {
220  if(phimax < 0) phimax += 2*M_PI;
221  else phimin -= 2*M_PI;
222  }
223 
224  LogTrace("AreaSeededTrackingRegionBuilder") << " phi min,max " << phimin << "," << phimax;
225 
226  minPhi = std::min(minPhi, phimin);
227  maxPhi = std::max(maxPhi, phimax);
228  }
229 
230  LogTrace("AreaSeededTrackingRegionBuilder") << " windows after this area eta " << minEta << "," << maxEta
231  << " phi " << minPhi << "," << maxPhi;
232  }
233 
234  const auto meanEta = (minEta+maxEta)/2.f;
235  const auto meanPhi = (minPhi+maxPhi)/2.f;
236  const auto dEta = maxEta-meanEta + m_conf->m_extraEta;
237  const auto dPhi = maxPhi-meanPhi + m_conf->m_extraPhi;
238 
239  auto useCandidates = false;
240  if (candidates.first) useCandidates = true;
241 
242  if (useCandidates){
243  // If we have objects used for seeding, loop over objects and find overlap with the found region. Return overlaps as tracking regions to use
244  for(const auto& object : *candidates.first) {
245  float dEta_Cand = candidates.second.first;
246  float dPhi_Cand = candidates.second.second;
247  float eta_Cand = object.eta();
248  float phi_Cand = object.phi();
249  float dEta_Cand_Point = std::abs(eta_Cand-meanEta);
250  float dPhi_Cand_Point = std::abs(deltaPhi(phi_Cand,meanPhi));
251 
252  if(dEta_Cand_Point > (dEta_Cand + dEta) || dPhi_Cand_Point > (dPhi_Cand + dPhi)) continue;
253 
254  float etaMin_RoI = std::max(eta_Cand-dEta_Cand,meanEta-dEta);
255  float etaMax_RoI = std::min(eta_Cand+dEta_Cand,meanEta+dEta);
256 
257  float phi_Cand_minus = normalizedPhi(phi_Cand-dPhi_Cand);
258  float phi_Point_minus = normalizedPhi(meanPhi-dPhi);
259  float phi_Cand_plus = normalizedPhi(phi_Cand+dPhi_Cand);
260  float phi_Point_plus = normalizedPhi(meanPhi+dPhi);
261 
262  float phiMin_RoI = deltaPhi(phi_Cand_minus,phi_Point_minus)>0. ? phi_Cand_minus : phi_Point_minus ;
263  float phiMax_RoI = deltaPhi(phi_Cand_plus,phi_Point_plus)<0. ? phi_Cand_plus : phi_Point_plus;
264 
265 
266  const auto meanEtaTemp = (etaMin_RoI+etaMax_RoI)/2.f;
267  auto meanPhiTemp = (phiMin_RoI+phiMax_RoI)/2.f;
268  if( phiMax_RoI < phiMin_RoI ) meanPhiTemp-=M_PI;
269  meanPhiTemp = normalizedPhi(meanPhiTemp);
270 
271  const auto dPhiTemp = deltaPhi(phiMax_RoI,meanPhiTemp);
272  const auto dEtaTemp = etaMax_RoI-meanEtaTemp;
273 
274  const auto x = std::cos(meanPhiTemp);
275  const auto y = std::sin(meanPhiTemp);
276  const auto z = 1./std::tan(2.f*std::atan(std::exp(-meanEtaTemp)));
277 
278  LogTrace("AreaSeededTrackingRegionsBuilder") << "Direction x,y,z " << x << "," << y << "," << z
279  << " eta,phi " << meanEtaTemp << "," << meanPhiTemp
280  << " window eta " << (meanEtaTemp-dEtaTemp) << "," << (meanEtaTemp+dEtaTemp)
281  << " phi " << (meanPhiTemp-dPhiTemp) << "," << (meanPhiTemp+dPhiTemp);
282 
283  return std::make_unique<RectangularEtaPhiTrackingRegion>(
284  GlobalVector(x,y,z),
285  origin.first, // GlobalPoint
286  m_conf->m_ptMin,
288  origin.second,
289  dEtaTemp,
290  dPhiTemp,
292  m_conf->m_precise,
295  );
296  }
297  // Have to retun nullptr here to ensure that we always return something
298  return nullptr;
299 
300  }
301  else{
302  const auto x = std::cos(meanPhi);
303  const auto y = std::sin(meanPhi);
304  const auto z = 1./std::tan(2.f*std::atan(std::exp(-meanEta)));
305 
306  LogTrace("AreaSeededTrackingRegionsBuilder") << "Direction x,y,z " << x << "," << y << "," << z
307  << " eta,phi " << meanEta << "," << meanPhi
308  << " window eta " << (meanEta-dEta) << "," << (meanEta+dEta)
309  << " phi " << (meanPhi-dPhi) << "," << (meanPhi+dPhi);
310 
311  return std::make_unique<RectangularEtaPhiTrackingRegion>(
312  GlobalVector(x,y,z),
313  origin.first, // GlobalPoint
314  m_conf->m_ptMin,
316  origin.second,
317  dEta,
318  dPhi,
320  m_conf->m_precise,
323  );
324  }
325 }
#define LogDebug(id)
const AreaSeededTrackingRegionsBuilder * m_conf
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double maxEta
T normalizedPhi(T phi)
Definition: normalizedPhi.h:9
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_whereToUseMeasurementTracker
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
T min(T a, T b)
Definition: MathUtil.h:58
#define LogTrace(id)
#define M_PI
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
T perp2() const
Squared magnitude of transverse component.
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
Global3DVector GlobalVector
Definition: GlobalVector.h:10
std::vector< std::unique_ptr< TrackingRegion > > AreaSeededTrackingRegionsBuilder::Builder::regions ( const Origins origins,
const std::vector< Area > &  areas 
) const

Definition at line 63 of file AreaSeededTrackingRegionsBuilder.cc.

References LogDebug, eostools::move(), and mps_fire::result.

Referenced by AreaSeededTrackingRegionsProducer::regions().

63  {
64  std::vector<std::unique_ptr<TrackingRegion> > result;
65 
66  // create tracking regions in directions of the points of interest
67  int n_regions = 0;
68  for(const auto& origin: origins) {
69  auto reg = region(origin, areas);
70  if (!reg) continue;
71  result.push_back(std::move(reg));
72  ++n_regions;
73  }
74  LogDebug("AreaSeededTrackingRegionsBuilder") << "produced "<<n_regions<<" regions";
75 
76  return result;
77 }
#define LogDebug(id)
std::unique_ptr< TrackingRegion > region(const Origin &origin, const std::vector< Area > &areas) const
def move(src, dest)
Definition: eostools.py:510
void AreaSeededTrackingRegionsBuilder::Builder::setCandidates ( const TrackingSeedCandidates::Objects  cands)
inline
void AreaSeededTrackingRegionsBuilder::Builder::setMeasurementTracker ( const MeasurementTrackerEvent mte)
inline

Definition at line 81 of file AreaSeededTrackingRegionsBuilder.h.

81 { m_measurementTracker = mte; }

Member Data Documentation

TrackingSeedCandidates::Objects AreaSeededTrackingRegionsBuilder::Builder::candidates
private

Definition at line 94 of file AreaSeededTrackingRegionsBuilder.h.

const AreaSeededTrackingRegionsBuilder* AreaSeededTrackingRegionsBuilder::Builder::m_conf = 0
private

Definition at line 92 of file AreaSeededTrackingRegionsBuilder.h.

const MeasurementTrackerEvent* AreaSeededTrackingRegionsBuilder::Builder::m_measurementTracker = 0
private

Definition at line 93 of file AreaSeededTrackingRegionsBuilder.h.