CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Classes | Typedefs | Functions | Variables
egamma::conv Namespace Reference

Classes

struct  ConversionInfo
 

Typedefs

using TrackRowView = TrackTable::const_iterator::value_type
 
using TrackTable = edm::soa::AddColumns< edm::soa::PtEtaPhiTable, TrackTableSpecificColumns >::type
 
using TrackTableSpecificColumns = std::tuple< edm::soa::col::Pz, edm::soa::col::PtError, edm::soa::col::MissingInnerHits, edm::soa::col::NumberOfValidHits, edm::soa::col::Charge, edm::soa::col::D0 >
 
using TrackTableView = edm::soa::ViewFromTable_t< TrackTable >
 

Functions

ConversionInfo const & arbitrateConversionPartnersbyR (const std::vector< ConversionInfo > &convCandidates)
 
ConversionInfo findBestConversionMatch (const std::vector< ConversionInfo > &v_convCandidates)
 
ConversionInfo findConversion (const reco::GsfElectronCore &gsfElectron, TrackTableView ctfTable, TrackTableView gsfTable, float bFieldAtOrigin, float minFracSharedHits=0.45f)
 
std::vector< ConversionInfofindConversions (const reco::GsfElectronCore &gsfElectron, TrackTableView ctfTable, TrackTableView gsfTable, float bFieldAtOrigin, float minFracSharedHits)
 
ConversionInfo getConversionInfo (reco::Track const &el_track, TrackRowView const &track, float bFieldAtOrigin)
 

Variables

constexpr float dR2Max = square(0.5f)
 
constexpr int maxDeltaMissingHits = 2
 
constexpr int maxDeltaMissingHitsForKFtoKF = 3
 
constexpr float maxDist2Dcot2 = square(0.05f)
 
constexpr float maxDistOrCotForKFtoKF = 0.02f
 
constexpr float maxRefPtErrorForKfConv = 0.05f
 
constexpr float maxRelGsfPtError = 0.5f
 
constexpr float maxRelGsfPtErrorForKfConv = 0.25f
 
constexpr float maxRelPtDiffForGsfConv = 0.25f
 
constexpr float maxRelPtDiffForKfConv = 0.2f
 
constexpr int minNumberOfValidHits = 5
 

Typedef Documentation

using egamma::conv::TrackRowView = typedef TrackTable::const_iterator::value_type

Definition at line 58 of file ConversionFinder.h.

Definition at line 56 of file ConversionFinder.h.

Definition at line 55 of file ConversionFinder.h.

Definition at line 57 of file ConversionFinder.h.

Function Documentation

ConversionInfo const& egamma::conv::arbitrateConversionPartnersbyR ( const std::vector< ConversionInfo > &  convCandidates)

Definition at line 209 of file ConversionFinder.cc.

References egamma::conv::ConversionInfo::dcot, egamma::conv::ConversionInfo::dist, funct::pow(), dttmaxenums::R, and groupFilesInBlocks::temp.

Referenced by findBestConversionMatch().

209  {
210  ConversionInfo const* closestConversion = &convCandidates.front();
211 
212  if (convCandidates.size() == 1)
213  return *closestConversion;
214 
215  float R = pow(closestConversion->dist, 2) + pow(closestConversion->dcot, 2);
216 
217  for (auto const& temp : convCandidates) {
218  float temp_R = pow(temp.dist, 2) + pow(temp.dcot, 2);
219  if (temp_R < R) {
220  R = temp_R;
221  closestConversion = &temp;
222  }
223  }
224 
225  return *closestConversion;
226  }
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
ConversionInfo egamma::conv::findBestConversionMatch ( const std::vector< ConversionInfo > &  v_convCandidates)

Definition at line 229 of file ConversionFinder.cc.

References funct::abs(), arbitrateConversionPartnersbyR(), maxDeltaMissingHits, maxDeltaMissingHitsForKFtoKF, maxDist2Dcot2, maxDistOrCotForKFtoKF, funct::pow(), square(), and groupFilesInBlocks::temp.

Referenced by findConversion().

229  {
230  using namespace std;
231 
232  if (v_convCandidates.empty())
233  return {};
234 
235  if (v_convCandidates.size() == 1)
236  return v_convCandidates.at(0);
237 
238  vector<ConversionInfo> v_0;
239  vector<ConversionInfo> v_1;
240  vector<ConversionInfo> v_2;
241  vector<ConversionInfo> v_3;
242  //loop over the candidates
243 
244  for (auto const& temp : v_convCandidates) {
245  if (temp.radiusOfConversion <= -2)
246  continue;
247 
248  if (temp.flag == 0) {
250  temp.deltaMissingHits < maxDeltaMissingHitsForKFtoKF) ||
251  (pow(temp.dist, 2) + pow(temp.dcot, 2) < maxDist2Dcot2 && temp.deltaMissingHits < maxDeltaMissingHits))
252  v_0.push_back(temp);
253  }
254 
255  if (temp.flag == 1) {
256  if (square(temp.dist) + square(temp.dcot) < maxDist2Dcot2 && temp.deltaMissingHits < maxDeltaMissingHits)
257  v_1.push_back(temp);
258  }
259  if (temp.flag == 2) {
260  if (square(temp.dist) + square(temp.dcot) < maxDist2Dcot2 && temp.deltaMissingHits < maxDeltaMissingHits)
261  v_2.push_back(temp);
262  }
263  if (temp.flag == 3) {
264  if (square(temp.dist) + square(temp.dcot) < maxDist2Dcot2 && temp.deltaMissingHits < maxDeltaMissingHits)
265  v_3.push_back(temp);
266  }
267 
268  } //candidate conversion loop
269 
270  //now do some arbitration
271 
272  //give preference to conversion partners found in the CTF collection
273  //using the electron's CTF track
274  if (!v_0.empty())
275  return arbitrateConversionPartnersbyR(v_0);
276 
277  if (!v_1.empty())
278  return arbitrateConversionPartnersbyR(v_1);
279 
280  if (!v_2.empty())
281  return arbitrateConversionPartnersbyR(v_2);
282 
283  if (!v_3.empty())
284  return arbitrateConversionPartnersbyR(v_3);
285 
286  //if we get here, we didn't find a candidate conversion partner that
287  //satisfied even the loose selections
288  //return the the closest partner by R
289  return arbitrateConversionPartnersbyR(v_convCandidates);
290  }
constexpr int maxDeltaMissingHitsForKFtoKF
constexpr float maxDistOrCotForKFtoKF
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ConversionInfo const & arbitrateConversionPartnersbyR(const std::vector< ConversionInfo > &convCandidates)
static double square(double x)
constexpr float maxDist2Dcot2
constexpr int maxDeltaMissingHits
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
ConversionInfo egamma::conv::findConversion ( const reco::GsfElectronCore gsfElectron,
TrackTableView  ctfTable,
TrackTableView  gsfTable,
float  bFieldAtOrigin,
float  minFracSharedHits = 0.45f 
)
inline

Definition at line 71 of file ConversionFinder.h.

References findBestConversionMatch(), and findConversions().

Referenced by GsfElectronAlgo::createElectron().

75  {
76  return findBestConversionMatch(findConversions(gsfElectron, ctfTable, gsfTable, bFieldAtOrigin, minFracSharedHits));
77  }
ConversionInfo findBestConversionMatch(const std::vector< ConversionInfo > &v_convCandidates)
std::vector< ConversionInfo > findConversions(const reco::GsfElectronCore &gsfElectron, TrackTableView ctfTable, TrackTableView gsfTable, float bFieldAtOrigin, float minFracSharedHits)
std::vector< ConversionInfo > egamma::conv::findConversions ( const reco::GsfElectronCore gsfElectron,
TrackTableView  ctfTable,
TrackTableView  gsfTable,
float  bFieldAtOrigin,
float  minFracSharedHits 
)

Definition at line 29 of file ConversionFinder.cc.

References funct::abs(), cuy::col, reco::GsfElectronCore::ctfGsfOverlap(), reco::GsfElectronCore::ctfTrack(), egamma::conv::ConversionInfo::dcot, reco::deltaR2(), egamma::conv::ConversionInfo::dist, dR2Max, getConversionInfo(), reco::GsfElectronCore::gsfTrack(), edm::Ref< C, T, F >::isNonnull(), maxRefPtErrorForKfConv, maxRelGsfPtError, maxRelGsfPtErrorForKfConv, maxRelPtDiffForGsfConv, maxRelPtDiffForKfConv, minNumberOfValidHits, egamma::conv::ConversionInfo::radiusOfConversion, and dt_dqm_sourceclient_common_cff::reco.

Referenced by findConversion().

33  {
34  using namespace reco;
35  using namespace std;
36  using namespace edm;
37  using namespace edm::soa::col;
38 
39  //get the references to the gsf and ctf tracks that are made
40  //by the electron
41  const reco::TrackRef eleCtfTk = gsfElectron.ctfTrack();
42  const reco::GsfTrackRef& eleGsfTk = gsfElectron.gsfTrack();
43 
44  float eleGsfPt = eleGsfTk->pt();
45  float eleGsfEta = eleGsfTk->eta();
46  float eleGsfPhi = eleGsfTk->phi();
47 
48  const bool useEleCtfTrack = eleCtfTk.isNonnull() && gsfElectron.ctfGsfOverlap() > minFracSharedHits;
49 
50  std::optional<float> eleCtfPt = std::nullopt;
51  std::optional<float> eleCtfEta = std::nullopt;
52  std::optional<float> eleCtfPhi = std::nullopt;
53 
54  //the electron's CTF track must share at least 45% of the inner hits
55  //with the electron's GSF track
56  std::optional<int> ctfidx = std::nullopt;
57  int gsfidx = static_cast<int>(eleGsfTk.key());
58 
59  if (useEleCtfTrack) {
60  eleCtfPt = eleCtfTk->pt();
61  eleCtfEta = eleCtfTk->eta();
62  eleCtfPhi = eleCtfTk->phi();
63  ctfidx = static_cast<int>(eleCtfTk.key());
64  }
65 
66  //these vectors are for those candidate partner tracks that pass our cuts
67  std::vector<ConversionInfo> v_candidatePartners;
68  //track indices required to make references
69  int ctftk_i = 0;
70  int gsftk_i = 0;
71 
72  //loop over the CTF tracks and try to find the partner track
73  for (auto ctftkItr = ctfTable.begin(); ctftkItr != ctfTable.end(); ++ctftkItr, ctftk_i++) {
74  if (useEleCtfTrack && ctftk_i == ctfidx.value())
75  continue;
76 
77  auto ctftk = *ctftkItr;
78 
79  //apply quality cuts to remove bad tracks
80  if (ctftk.get<PtError>() > maxRefPtErrorForKfConv * ctftk.get<Pt>() ||
82  continue;
83 
84  if (useEleCtfTrack && std::abs(ctftk.get<Pt>() - eleCtfPt.value()) < maxRelPtDiffForKfConv * eleCtfPt.value())
85  continue;
86 
87  //use the electron's CTF track, if not null, to search for the partner track
88  //look only in a cone of 0.5 to save time, and require that the track is opp. sign
89  if (useEleCtfTrack && eleCtfTk->charge() + ctftk.get<Charge>() == 0 &&
90  deltaR2(eleCtfEta.value(), eleCtfPhi.value(), ctftk.get<Eta>(), ctftk.get<Phi>()) < dR2Max) {
91  ConversionInfo convInfo = getConversionInfo(*eleCtfTk, ctftk, bFieldAtOrigin);
92 
93  //need to add the track reference information for completeness
94  //because the overloaded fnc above does not make a trackRef
95  int deltaMissingHits = ctftk.get<MissingInnerHits>() - eleCtfTk->missingInnerHits();
96 
97  v_candidatePartners.push_back(
98  {convInfo.dist, convInfo.dcot, convInfo.radiusOfConversion, ctftk_i, std::nullopt, deltaMissingHits, 0});
99 
100  } //using the electron's CTF track
101 
102  //now we check using the electron's gsf track
103  if (eleGsfTk->charge() + ctftk.get<Charge>() == 0 &&
104  deltaR2(eleGsfEta, eleGsfPhi, ctftk.get<Eta>(), ctftk.get<Phi>()) < dR2Max &&
105  eleGsfTk->ptError() < maxRelGsfPtErrorForKfConv * eleGsfPt) {
106  int deltaMissingHits = ctftk.get<MissingInnerHits>() - eleGsfTk->missingInnerHits();
107 
108  ConversionInfo convInfo = getConversionInfo(*eleGsfTk, ctftk, bFieldAtOrigin);
109 
110  v_candidatePartners.push_back(
111  {convInfo.dist, convInfo.dcot, convInfo.radiusOfConversion, ctftk_i, std::nullopt, deltaMissingHits, 1});
112  } //using the electron's GSF track
113 
114  } //loop over the CTF track collection
115 
116  //------------------------------------------------------ Loop over GSF collection ----------------------------------//
117  for (auto gsftkItr = gsfTable.begin(); gsftkItr != gsfTable.end(); ++gsftkItr, gsftk_i++) {
118  //reject the electron's own gsfTrack
119  if (gsfidx == gsftk_i)
120  continue;
121 
122  auto gsftk = *gsftkItr;
123 
124  //apply quality cuts to remove bad tracks
125  if (gsftk.get<PtError>() > maxRelGsfPtError * gsftk.get<Pt>() ||
127  continue;
128 
129  if (std::abs(gsftk.get<Pt>() - eleGsfPt) < maxRelPtDiffForGsfConv * eleGsfPt)
130  continue;
131 
132  //try using the electron's CTF track first if it exists
133  //look only in a cone of 0.5 around the electron's track
134  //require opposite sign
135  if (useEleCtfTrack && eleCtfTk->charge() + gsftk.get<Charge>() == 0 &&
136  deltaR2(eleCtfEta.value(), eleCtfPhi.value(), gsftk.get<Eta>(), gsftk.get<Phi>()) < dR2Max) {
137  int deltaMissingHits = gsftk.get<MissingInnerHits>() - eleCtfTk->missingInnerHits();
138 
139  ConversionInfo convInfo = getConversionInfo(*eleCtfTk, gsftk, bFieldAtOrigin);
140  //fill the Ref info
141  v_candidatePartners.push_back(
142  {convInfo.dist, convInfo.dcot, convInfo.radiusOfConversion, std::nullopt, gsftk_i, deltaMissingHits, 2});
143  }
144 
145  //use the electron's gsf track
146  if (eleGsfTk->charge() + gsftk.get<Charge>() == 0 &&
147  deltaR2(eleGsfEta, eleGsfPhi, gsftk.get<Eta>(), gsftk.get<Phi>()) < dR2Max &&
148  (eleGsfTk->ptError() < maxRelGsfPtError * eleGsfPt)) {
149  ConversionInfo convInfo = getConversionInfo(*eleGsfTk, gsftk, bFieldAtOrigin);
150  //fill the Ref info
151 
152  int deltaMissingHits = gsftk.get<MissingInnerHits>() - eleGsfTk->missingInnerHits();
153 
154  v_candidatePartners.push_back(
155  {convInfo.dist, convInfo.dcot, convInfo.radiusOfConversion, std::nullopt, gsftk_i, deltaMissingHits, 3});
156  }
157  } //loop over the gsf track collection
158 
159  return v_candidatePartners;
160  }
const double dR2Max
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
constexpr float maxRelGsfPtError
float ctfGsfOverlap() const
const GsfTrackRef & gsfTrack() const
ConversionInfo getConversionInfo(reco::Track const &el_track, TrackRowView const &track, float bFieldAtOrigin)
constexpr float maxRelGsfPtErrorForKfConv
constexpr float maxRefPtErrorForKfConv
constexpr float maxRelPtDiffForGsfConv
constexpr float maxRelPtDiffForKfConv
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
TrackRef ctfTrack() const
constexpr int minNumberOfValidHits
int col
Definition: cuy.py:1009
ConversionInfo egamma::conv::getConversionInfo ( reco::Track const &  el_track,
TrackRowView const &  track,
float  bFieldAtOrigin 
)

Definition at line 163 of file ConversionFinder.cc.

References funct::abs(), reco::TrackBase::charge(), cuy::col, funct::cos(), ztail::d, reco::TrackBase::d0(), validate-o2o-wbm::f, reco::TrackBase::phi(), funct::pow(), reco::TrackBase::pt(), reco::TrackBase::px(), reco::TrackBase::py(), reco::TrackBase::pz(), funct::sin(), and mathSSE::sqrt().

Referenced by findConversions().

163  {
164  using namespace edm::soa::col;
165 
166  //now calculate the conversion related information
167  float rEl = 100.f * ele.pt() / (-0.3f * bFieldAtOrigin * ele.charge());
168  float xEl = -1.f * (rEl - ele.d0()) * sin(ele.phi());
169  float yEl = (rEl - ele.d0()) * cos(ele.phi());
170  rEl = std::abs(rEl);
171 
172  float rCand = 100.f * track.get<Pt>() / (-0.3f * bFieldAtOrigin * track.get<Charge>());
173  float xCand = -1.f * (rCand - track.get<D0>()) * sin(track.get<Phi>());
174  float yCand = (rCand - track.get<D0>()) * cos(track.get<Phi>());
175  rCand = std::abs(rCand);
176 
177  float d = sqrt(pow(xEl - xCand, 2) + pow(yEl - yCand, 2));
178  float dist = d - (rEl + rCand);
179 
180  // this is equivalent to `1/tan(theta_1) - 1/tan(theta_2)` but requires less trigonometry
181  float dcot = ele.pz() / ele.pt() - track.get<Pz>() / track.get<Pt>();
182 
183  //get the point of conversion
184  float xa1 = xEl + (xCand - xEl) * rEl / d;
185  float xa2 = xCand + (xEl - xCand) * rCand / d;
186  float ya1 = yEl + (yCand - yEl) * rEl / d;
187  float ya2 = yCand + (yEl - yCand) * rCand / d;
188 
189  float x = .5f * (xa1 + xa2);
190  float y = .5f * (ya1 + ya2);
191  float rconv = sqrt(pow(x, 2) + pow(y, 2));
192  // The z-position of the conversion is unused, but here is how it could be computed if needed:
193  // float z = ele.dz() + rEl * ele.pz() * std::acos(1 - pow(rconv, 2) / (2. * pow(rEl, 2))) / ele.pt();
194 
195  //now assign a sign to the radius of conversion
196  float tempsign = ele.px() * x + ele.py() * y;
197  tempsign = tempsign / std::abs(tempsign);
198  rconv = tempsign * rconv;
199 
200  //return an instance of ConversionInfo, but with a NULL track refs
201  return {dist, dcot, rconv};
202  }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
tuple d
Definition: ztail.py:151
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int col
Definition: cuy.py:1009
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

Variable Documentation

constexpr float egamma::conv::dR2Max = square(0.5f)

Definition at line 18 of file ConversionFinder.cc.

Referenced by findConversions().

constexpr int egamma::conv::maxDeltaMissingHits = 2

Definition at line 21 of file ConversionFinder.cc.

Referenced by findBestConversionMatch().

constexpr int egamma::conv::maxDeltaMissingHitsForKFtoKF = 3

Definition at line 23 of file ConversionFinder.cc.

Referenced by findBestConversionMatch().

constexpr float egamma::conv::maxDist2Dcot2 = square(0.05f)

Definition at line 20 of file ConversionFinder.cc.

Referenced by findBestConversionMatch().

constexpr float egamma::conv::maxDistOrCotForKFtoKF = 0.02f

Definition at line 24 of file ConversionFinder.cc.

Referenced by findBestConversionMatch().

constexpr float egamma::conv::maxRefPtErrorForKfConv = 0.05f

Definition at line 14 of file ConversionFinder.cc.

Referenced by findConversions().

constexpr float egamma::conv::maxRelGsfPtError = 0.5f

Definition at line 22 of file ConversionFinder.cc.

Referenced by findConversions().

constexpr float egamma::conv::maxRelGsfPtErrorForKfConv = 0.25f

Definition at line 15 of file ConversionFinder.cc.

Referenced by findConversions().

constexpr float egamma::conv::maxRelPtDiffForGsfConv = 0.25f

Definition at line 17 of file ConversionFinder.cc.

Referenced by findConversions().

constexpr float egamma::conv::maxRelPtDiffForKfConv = 0.2f

Definition at line 16 of file ConversionFinder.cc.

Referenced by findConversions().

constexpr int egamma::conv::minNumberOfValidHits = 5

Definition at line 19 of file ConversionFinder.cc.

Referenced by findConversions().