CMS 3D CMS Logo

LowPtConversion.cc
Go to the documentation of this file.
2 
4 // Matched to any conversion (without selections)
5 //
6 bool LowPtConversion::wpOpen() const { return matched_; }
7 
9 // Nancy's baseline selections for conversions
10 // Based on: https://github.com/CMSBParking/BParkingNANO/blob/b2664ed/BParkingNano/plugins/ConversionSelector.cc#L253-L300
12  return (wpOpen() && ntracks_ == 2 && valid_ && quality_high_purity_ && chi2prob_ > 0.0005);
13 }
14 
16 // Nancy's selection for analysis of conversions
17 // Based on: slide 20 of https://indico.cern.ch/event/814402/contributions/3401312/
19  return (wpLoose() && sum_nhits_before_vtx_ <= 1 && l_xy_ > 0. && mass_from_conv_ > 0. && // sanity check
20  mass_from_conv_ < 0.05);
21 }
22 
24 // adds minimal set of flags to electron userData
26  ele.addUserInt("convOpen", matched_ ? 1 : 0);
27  ele.addUserInt("convLoose", wpLoose() ? 1 : 0);
28  ele.addUserInt("convTight", wpTight() ? 1 : 0);
29  ele.addUserInt("convLead", matched_lead_.isNonnull() ? 1 : 0);
30  ele.addUserInt("convTrail", matched_trail_.isNonnull() ? 1 : 0);
31  if (ele.hasUserInt("convExtra") == false) {
32  ele.addUserInt("convExtra", 0);
33  }
34 }
35 
37 // adds all variables to electron userData
39  // Flag that indicates if extra variables are added to electron userData
40  ele.addUserInt("convExtra", 1, true); // overwrite
41 
42  // quality
43  ele.addUserInt("convValid", valid_ ? 1 : 0);
44  ele.addUserFloat("convChi2Prob", chi2prob_);
45  ele.addUserInt("convQualityHighPurity", quality_high_purity_ ? 1 : 0);
46  ele.addUserInt("convQualityHighEff", quality_high_efficiency_ ? 1 : 0);
47 
48  // tracks
49  ele.addUserInt("convTracksN", ntracks_);
50  ele.addUserFloat("convMinTrkPt", min_trk_pt_);
51  ele.addUserInt("convLeadIdx", ilead_);
52  ele.addUserInt("convTrailIdx", itrail_);
53 
54  // displacement
55  ele.addUserFloat("convLxy", l_xy_);
56  ele.addUserFloat("convVtxRadius", vtx_radius_);
57 
58  // invariant mass
59  ele.addUserFloat("convMass", mass_from_conv_);
60  ele.addUserFloat("convMassFromPin", mass_from_Pin_);
61  ele.addUserFloat("convMassBeforeFit", mass_before_fit_);
62  ele.addUserFloat("convMassAfterFit", mass_after_fit_);
63 
64  // hits before vertex
65  ele.addUserInt("convLeadNHitsBeforeVtx", lead_nhits_before_vtx_);
66  ele.addUserInt("convTrailNHitsBeforeVtx", trail_nhits_before_vtx_);
67  ele.addUserInt("convMaxNHitsBeforeVtx", max_nhits_before_vtx_);
68  ele.addUserInt("convSumNHitsBeforeVtx", sum_nhits_before_vtx_);
69  ele.addUserInt("convDeltaExpectedNHitsInner", delta_expected_nhits_inner_);
70 
71  // opening angle
72  ele.addUserFloat("convDeltaCotFromPin", delta_cot_from_Pin_);
73 }
74 
76 //
79  const pat::Electron& ele) {
80  // Iterate through conversions and calculate quantities (requirement from Nancy)
81  for (const auto& conv : conversions) {
82  // Filter
83  if (conv.tracks().size() != 2) {
84  continue;
85  }
86 
87  // Quality
88  valid_ = conv.conversionVertex().isValid(); // (=true)
89  chi2prob_ = ChiSquaredProbability(conv.conversionVertex().chi2(), conv.conversionVertex().ndof()); // (<0.005)
92 
93  // Tracks
94  ntracks_ = conv.tracks().size(); // (=2)
95  min_trk_pt_ = -1.; // (>0.5)
96  for (const auto& trk : conv.tracks()) {
97  if (trk.isNonnull() && trk.isAvailable() && (min_trk_pt_ < 0. || trk->pt() < min_trk_pt_)) {
98  min_trk_pt_ = trk->pt();
99  }
100  }
101  ilead_ = -1;
102  itrail_ = -1;
103  if (conv.tracks().size() == 2) {
104  const edm::RefToBase<reco::Track>& trk1 = conv.tracks().front();
105  const edm::RefToBase<reco::Track>& trk2 = conv.tracks().back();
106  if (trk1.isNonnull() && trk1.isAvailable() && trk2.isNonnull() && trk2.isAvailable()) {
107  if (trk1->pt() > trk2->pt()) {
108  ilead_ = 0;
109  itrail_ = 1;
110  } else {
111  ilead_ = 1;
112  itrail_ = 0;
113  }
114  }
115  }
116 
117  // Transverse displacement (with respect to beamspot) and vertex radius
118  math::XYZVectorF p_refitted = conv.refittedPairMomentum();
119  float dx = conv.conversionVertex().x() - beamSpot.x0();
120  float dy = conv.conversionVertex().y() - beamSpot.y0();
121  l_xy_ = (p_refitted.x() * dx + p_refitted.y() * dy) / p_refitted.rho();
122  vtx_radius_ = sqrt(conv.conversionVertex().position().perp2()); // (1.5<r<4.)
123 
124  // invariant mass from track pair from conversion
125  mass_from_conv_ = conv.pairInvariantMass();
126 
127  // Invariant mass from Pin before fit to common vertex
128  if (conv.tracksPin().size() >= 2 && ilead_ > -1 && itrail_ > -1) {
129  math::XYZVectorF lead_Pin = conv.tracksPin().at(ilead_);
130  math::XYZVectorF trail_Pin = conv.tracksPin().at(itrail_);
131  mass_from_Pin_ = mee(lead_Pin.x(), lead_Pin.y(), lead_Pin.z(), trail_Pin.x(), trail_Pin.y(), trail_Pin.z());
132  // Opening angle
133  delta_cot_from_Pin_ = 1. / tan(trail_Pin.theta()) - 1. / tan(lead_Pin.theta());
134  }
135 
136  // Invariant mass before fit to common vertex
137  if (conv.tracks().size() >= 2 && ilead_ > -1 && itrail_ > -1) {
138  auto lead_before_vtx_fit = conv.tracks().at(ilead_)->momentum();
139  auto trail_before_vtx_fit = conv.tracks().at(itrail_)->momentum();
140  mass_before_fit_ = mee(lead_before_vtx_fit.x(),
141  lead_before_vtx_fit.y(),
142  lead_before_vtx_fit.z(),
143  trail_before_vtx_fit.x(),
144  trail_before_vtx_fit.y(),
145  trail_before_vtx_fit.z());
146  }
147 
148  // Invariant mass after the fit to common vertex
149  if (conv.conversionVertex().refittedTracks().size() >= 2 && ilead_ > -1 && itrail_ > -1) {
150  auto const& lead_after_vtx_fit = conv.conversionVertex().refittedTracks().at(ilead_);
151  auto const& trail_after_vtx_fit = conv.conversionVertex().refittedTracks().at(itrail_);
152  mass_after_fit_ = mee(lead_after_vtx_fit.px(),
153  lead_after_vtx_fit.py(),
154  lead_after_vtx_fit.pz(),
155  trail_after_vtx_fit.px(),
156  trail_after_vtx_fit.py(),
157  trail_after_vtx_fit.pz());
158  // Difference in expeted hits
160  lead_after_vtx_fit.hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS) -
161  trail_after_vtx_fit.hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
162  }
163 
164  // Hits prior to vertex
165  if (ilead_ > -1 && itrail_ > -1) {
166  auto const& nHits = conv.nHitsBeforeVtx();
167  bool enoughTracks = nHits.size() > 1;
168  lead_nhits_before_vtx_ = enoughTracks ? nHits.at(ilead_) : 0;
169  trail_nhits_before_vtx_ = enoughTracks ? nHits.at(itrail_) : 0;
170  max_nhits_before_vtx_ = enoughTracks ? std::max(nHits[0], nHits[1]) : 0;
171  sum_nhits_before_vtx_ = enoughTracks ? nHits[0] + nHits[1] : 0;
172  }
173 
174  // Attempt to match conversion track to electron
175  for (uint itrk = 0; itrk < conv.tracks().size(); ++itrk) {
176  const edm::RefToBase<reco::Track> trk = conv.tracks()[itrk];
177  if (trk.isNull()) {
178  continue;
179  }
180  reco::GsfTrackRef ref = ele.core()->gsfTrack();
182  if (gsf.isNull()) {
183  continue;
184  }
185  if (ref.id() == trk.id() && ref.key() == trk.key()) {
186  matched_ = true;
187  if (static_cast<int>(itrk) == ilead_) {
188  matched_lead_ = trk;
189  }
190  if (static_cast<int>(itrk) == itrail_) {
191  matched_trail_ = trk;
192  }
193  }
194  } // track loop
195  } // conversions loop
196 
197  return matched_;
198 }
199 
201 //
202 float LowPtConversion::mee(float px1, float py1, float pz1, float px2, float py2, float pz2) {
203  const float m = 0.000511;
204  const float px = px1 + px2;
205  const float py = py1 + py2;
206  const float pz = pz1 + pz2;
207  const float p1 = px1 * px1 + py1 * py1 + pz1 * pz1;
208  const float p2 = px2 * px2 + py2 * py2 + pz2 * pz2;
209  const float e = sqrt(p1 + m * m) + sqrt(p2 + m * m);
210  const float mass = (e * e - px * px - py * py - pz * pz);
211  return mass > 0. ? sqrt(mass) : -1.;
212 }
bool isAvailable() const
Definition: RefToBase.h:126
int delta_expected_nhits_inner_
bool hasUserInt(const std::string &key) const
Return true if there is a user-defined int with a given name.
Definition: PATObject.h:390
ProductID id() const
Accessor for product ID.
Definition: Ref.h:244
bool isNonnull() const
Checks for non-null.
Definition: RefToBase.h:290
bool wpOpen() const
bool quality_high_efficiency_
reco::GsfTrackRef gsfTrack() const override
override the reco::GsfElectron::gsfTrack method, to access the internal storage of the supercluster ...
void addUserFloat(const std::string &label, float data, const bool overwrite=false)
Set user-defined float.
Definition: PATObject.h:892
key_type key() const
Accessor for product key.
Definition: Ref.h:250
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
ProductID id() const
Definition: RefToBase.h:221
double pt() const
track transverse momentum
Definition: TrackBase.h:637
edm::RefToBase< reco::Track > matched_trail_
void addExtraUserVars(pat::Electron &ele) const
T sqrt(T t)
Definition: SSEVec.h:23
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
float ChiSquaredProbability(double chiSquared, double nrDOF)
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:16
bool isNull() const
Checks for null.
Definition: RefToBase.h:284
bool wpTight() const
void addUserInt(const std::string &label, int32_t data, const bool overwrite=false)
Set user-defined int.
Definition: PATObject.h:929
size_t key() const
Definition: RefToBase.h:226
edm::RefToBase< reco::Track > matched_lead_
bool wpLoose() const
Analysis-level electron class.
Definition: Electron.h:51
EPOS::IO_EPOS conv
bool match(const reco::BeamSpot &beamSpot, const reco::ConversionCollection &conversions, const pat::Electron &ele)
uint trail_nhits_before_vtx_
static float mee(float ipx1, float ipy1, float ipz1, float ipx2, float ipy2, float ipz2)
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
reco::GsfElectronCoreRef core() const override
override the virtual reco::GsfElectron::core method, so that the embedded core can be used by GsfElec...
void addUserVars(pat::Electron &ele) const