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, const MagneticField *field, const MultipleScatteringParametrisationMaker *msmaker)
 
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 = nullptr
 
const MagneticFieldm_field = nullptr
 
const MeasurementTrackerEventm_measurementTracker = nullptr
 
const MultipleScatteringParametrisationMakerm_msmaker = nullptr
 

Detailed Description

Definition at line 79 of file AreaSeededTrackingRegionsBuilder.h.

Constructor & Destructor Documentation

◆ Builder()

AreaSeededTrackingRegionsBuilder::Builder::Builder ( const AreaSeededTrackingRegionsBuilder conf,
const MagneticField field,
const MultipleScatteringParametrisationMaker msmaker 
)
inlineexplicit

Definition at line 81 of file AreaSeededTrackingRegionsBuilder.h.

84  : m_conf(conf), m_field(field), m_msmaker(msmaker) {}
const AreaSeededTrackingRegionsBuilder * m_conf
const MultipleScatteringParametrisationMaker * m_msmaker

◆ ~Builder()

AreaSeededTrackingRegionsBuilder::Builder::~Builder ( )
default

Member Function Documentation

◆ region() [1/2]

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

Definition at line 89 of file AreaSeededTrackingRegionsBuilder.cc.

Referenced by regions().

90  {
91  return regionImpl(origin, areas);
92 }
std::unique_ptr< TrackingRegion > regionImpl(const Origin &origin, const T &areas) const

◆ region() [2/2]

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

Definition at line 93 of file AreaSeededTrackingRegionsBuilder.cc.

94  {
95  return regionImpl(origin, areas);
96 }
std::unique_ptr< TrackingRegion > regionImpl(const Origin &origin, const T &areas) const

◆ regionImpl()

template<typename T >
std::unique_ptr< TrackingRegion > AreaSeededTrackingRegionsBuilder::Builder::regionImpl ( const Origin origin,
const T areas 
) const
private

Definition at line 99 of file AreaSeededTrackingRegionsBuilder.cc.

References funct::abs(), custom_jme_cff::area, HLT_2022v15_cff::candidates, funct::cos(), ztail::d, SiPixelRawToDigiRegional_cfi::deltaPhi, reco::deltaPhi(), HLT_2022v15_cff::dEta, HLT_2022v15_cff::dPhi, muonTiming_cfi::etamax, muonTiming_cfi::etamin, JetChargeProducer_cfi::exp, f, LogDebug, LogTrace, M_PI, SiStripPI::max, razorScouting_cff::maxEta, HLT_2022v15_cff::maxPhi, SiStripPI::min, EgHLTOffEleSelection_cfi::minEta, HLT_2022v15_cff::minPhi, normalizedPhi(), perp2(), phimax, phimin, alignCSCRings::r, Validation_hcalonly_cfi::sign, funct::sin(), mathSSE::sqrt(), funct::tan(), createJobs::tmp, x, y, and z.

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

◆ regions()

std::vector< std::unique_ptr< TrackingRegion > > AreaSeededTrackingRegionsBuilder::Builder::regions ( const Origins origins,
const std::vector< Area > &  areas 
) const

Definition at line 71 of file AreaSeededTrackingRegionsBuilder.cc.

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

Referenced by AreaSeededTrackingRegionsProducer::regions().

72  {
73  std::vector<std::unique_ptr<TrackingRegion> > result;
74 
75  // create tracking regions in directions of the points of interest
76  int n_regions = 0;
77  for (const auto& origin : origins) {
78  auto reg = region(origin, areas);
79  if (!reg)
80  continue;
81  result.push_back(std::move(reg));
82  ++n_regions;
83  }
84  LogDebug("AreaSeededTrackingRegionsBuilder") << "produced " << n_regions << " regions";
85 
86  return result;
87 }
std::unique_ptr< TrackingRegion > region(const Origin &origin, const std::vector< Area > &areas) const
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)

◆ setCandidates()

void AreaSeededTrackingRegionsBuilder::Builder::setCandidates ( const TrackingSeedCandidates::Objects  cands)
inline

◆ setMeasurementTracker()

void AreaSeededTrackingRegionsBuilder::Builder::setMeasurementTracker ( const MeasurementTrackerEvent mte)
inline

Definition at line 87 of file AreaSeededTrackingRegionsBuilder.h.

References m_measurementTracker.

87 { m_measurementTracker = mte; }

Member Data Documentation

◆ candidates

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

Definition at line 102 of file AreaSeededTrackingRegionsBuilder.h.

Referenced by setCandidates().

◆ m_conf

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

Definition at line 98 of file AreaSeededTrackingRegionsBuilder.h.

◆ m_field

const MagneticField* AreaSeededTrackingRegionsBuilder::Builder::m_field = nullptr
private

Definition at line 100 of file AreaSeededTrackingRegionsBuilder.h.

◆ m_measurementTracker

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

Definition at line 99 of file AreaSeededTrackingRegionsBuilder.h.

Referenced by setMeasurementTracker().

◆ m_msmaker

const MultipleScatteringParametrisationMaker* AreaSeededTrackingRegionsBuilder::Builder::m_msmaker = nullptr
private

Definition at line 101 of file AreaSeededTrackingRegionsBuilder.h.