CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Attributes
reco::Conversion Class Reference

#include <DataFormats/EgammaCandidates/interface/Conversion.h>

Public Types

enum  ConversionAlgorithm {
  undefined = 0, ecalSeeded = 1, trackerOnly = 2, mixed = 3,
  pflow = 4, algoSize = 5
}
 
enum  ConversionQuality {
  generalTracksOnly = 0, arbitratedEcalSeeded = 1, arbitratedMerged = 2, arbitratedMergedEcalGeneral = 3,
  gsfTracksOpenOnly = 4, highPurity = 8, highEfficiency = 9, ecalMatched1Track = 10,
  ecalMatched2Track = 11
}
 

Public Member Functions

ConversionAlgorithm algo () const
 
std::string algoName () const
 
const std::vector< reco::CaloClusterPtr > & bcMatchingWithTracks () const
 
reco::CaloClusterPtrVector caloCluster () const
 Pointer to CaloCluster (foe Egamma Conversions it points to a SuperCluster) More...
 
Conversionclone () const
 returns a clone of the candidate More...
 
 Conversion ()
 
 Conversion (const reco::CaloClusterPtrVector &clu, const std::vector< edm::RefToBase< reco::Track > > &tr, const std::vector< math::XYZPointF > &trackPositionAtEcal, const reco::Vertex &convVtx, const std::vector< reco::CaloClusterPtr > &matchingBC, const float DCA, const std::vector< math::XYZPointF > &innPoint, const std::vector< math::XYZVectorF > &trackPin, const std::vector< math::XYZVectorF > &trackPout, const std::vector< uint8_t > &nHitsBeforeVtx, const std::vector< Measurement1DFloat > &dlClosestHitToVtx, uint8_t nSharedHits, const float mva, ConversionAlgorithm=undefined)
 
 Conversion (const reco::CaloClusterPtrVector &clu, const std::vector< reco::TrackRef > &tr, const std::vector< math::XYZPointF > &trackPositionAtEcal, const reco::Vertex &convVtx, const std::vector< reco::CaloClusterPtr > &matchingBC, const float DCA, const std::vector< math::XYZPointF > &innPoint, const std::vector< math::XYZVectorF > &trackPin, const std::vector< math::XYZVectorF > &trackPout, const float mva, ConversionAlgorithm=undefined)
 
 Conversion (const reco::CaloClusterPtrVector &clu, const std::vector< reco::TrackRef > &tr, const reco::Vertex &convVtx, ConversionAlgorithm=undefined)
 
 Conversion (const reco::CaloClusterPtrVector &clu, const std::vector< edm::RefToBase< reco::Track > > &tr, const reco::Vertex &convVtx, ConversionAlgorithm=undefined)
 
const reco::VertexconversionVertex () const
 returns the reco conversion vertex More...
 
double dEtaTracksAtEcal () const
 
double distOfMinimumApproach () const
 
const std::vector< Measurement1DFloat > & dlClosestHitToVtx () const
 Vector of signed decay length with uncertainty from nearest hit on track to the conversion vtx positions. More...
 
double dPhiTracksAtEcal () const
 
double dPhiTracksAtVtx () const
 
double dxy (const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
 
double dz (const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
 
const std::vector< math::XYZPointF > & ecalImpactPosition () const
 
double EoverP () const
 
double EoverPrefittedTracks () const
 
bool isConverted () const
 Bool flagging objects having track size >0. More...
 
double lxy (const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
 
double lz (const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
 
double MVAout () const
 get the value of the TMVA output More...
 
const std::vector< uint8_t > & nHitsBeforeVtx () const
 Vector of the number of hits before the vertex along each track trajector. More...
 
uint8_t nSharedHits () const
 number of shared hits btw the two track More...
 
unsigned int nTracks () const
 Number of tracks= 0,1,2. More...
 
std::vector< float > const oneLegMVA ()
 get the MVS output from PF for one leg conversions More...
 
double pairCotThetaSeparation () const
 Delta cot(Theta) where Theta is the angle in the (y,z) plane between the two tracks. Original tracks are used. More...
 
double pairInvariantMass () const
 if nTracks=2 returns the pair invariant mass. Original tracks are used here More...
 
math::XYZVectorF pairMomentum () const
 Conversion tracks momentum from the tracks inner momentum. More...
 
bool quality (ConversionQuality q) const
 
math::XYZTLorentzVectorF refittedPair4Momentum () const
 Conversion track pair 4-momentum from the tracks refitted with vertex constraint. More...
 
math::XYZVectorF refittedPairMomentum () const
 Conversion tracks momentum from the tracks refitted with vertex constraint. More...
 
void setConversionAlgorithm (const ConversionAlgorithm a, bool set=true)
 Conversion Track algorithm/provenance. More...
 
void setMatchingSuperCluster (const reco::CaloClusterPtrVector &sc)
 
void setMVAout (const float &mva)
 set the value of the TMVA output More...
 
void setOneLegMVA (const std::vector< float > &mva)
 set the MVS output from PF for one leg conversions More...
 
void setQuality (ConversionQuality q, bool b)
 
std::vector< edm::RefToBase< reco::Track > > const & tracks () const
 vector of track to base references More...
 
const std::vector< math::XYZPointF > & tracksInnerPosition () const
 Vector containing the position of the innermost hit of each track. More...
 
const std::vector< math::XYZVectorF > & tracksPin () const
 Vector of track momentum measured at the innermost hit. More...
 
const std::vector< math::XYZVectorF > & tracksPout () const
 Vector of track momentum measured at the outermost hit. More...
 
std::vector< double > tracksSigned_d0 () const
 signed transverse impact parameter for each track More...
 
double zOfPrimaryVertexFromTracks (const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
 

Static Public Member Functions

static ConversionAlgorithm algoByName (const std::string &name)
 
static std::string algoName (ConversionAlgorithm)
 

Static Public Attributes

static const std::string algoNames [] = {"undefined", "ecalSeeded", "trackerOnly", "mixed", "pflow"}
 

Private Attributes

uint8_t algorithm_
 conversion algorithm/provenance More...
 
reco::CaloClusterPtrVector caloCluster_
 vector pointer to a/multiple seed CaloCluster(s) More...
 
std::vector< Measurement1DFloatdlClosestHitToVtx_
 signed decay length and uncertainty from nearest hit on track to conversion vertex More...
 
std::vector< uint8_t > nHitsBeforeVtx_
 number of hits before the vertex on each trackerOnly More...
 
uint8_t nSharedHits_
 number of shared hits between tracks More...
 
uint16_t qualityMask_
 
reco::Vertex theConversionVertex_
 Fitted Kalman conversion vertex. More...
 
std::vector< reco::CaloClusterPtrtheMatchingBCs_
 Clusters mathing the tracks (these are not the seeds) More...
 
float theMinDistOfApproach_
 Distance of min approach of the two tracks. More...
 
float theMVAout_
 TMVA output. More...
 
std::vector< float > theOneLegMVA_
 vectors of TMVA outputs from pflow for one leg conversions More...
 
std::vector< math::XYZPointFthePositionAtEcal_
 position at the ECAl surface of the track extrapolation More...
 
std::vector< math::XYZPointFtheTrackInnerPosition_
 P_in of tracks. More...
 
std::vector< math::XYZVectorFtheTrackPin_
 P_in of tracks. More...
 
std::vector< math::XYZVectorFtheTrackPout_
 P_out of tracks. More...
 
std::vector< edm::RefToBase< reco::Track > > trackToBaseRefs_
 vector Track RefToBase More...
 

Detailed Description

Author
N.Marinelli University of Notre Dame, US

Definition at line 23 of file Conversion.h.

Member Enumeration Documentation

◆ ConversionAlgorithm

Enumerator
undefined 
ecalSeeded 
trackerOnly 
mixed 
pflow 
algoSize 

Definition at line 25 of file Conversion.h.

◆ ConversionQuality

Enumerator
generalTracksOnly 
arbitratedEcalSeeded 
arbitratedMerged 
arbitratedMergedEcalGeneral 
gsfTracksOpenOnly 
highPurity 
highEfficiency 
ecalMatched1Track 
ecalMatched2Track 

Definition at line 27 of file Conversion.h.

Constructor & Destructor Documentation

◆ Conversion() [1/5]

Conversion::Conversion ( )

Definition at line 121 of file Conversion.cc.

References algorithm_, nSharedHits_, qualityMask_, theMinDistOfApproach_, theMVAout_, thePositionAtEcal_, theTrackInnerPosition_, theTrackPin_, and theTrackPout_.

Referenced by clone().

121  {
122  algorithm_ = 0;
123  qualityMask_ = 0;
124  theMinDistOfApproach_ = 9999.;
125  nSharedHits_ = 0;
126  theMVAout_ = 9999.;
127  thePositionAtEcal_.push_back(math::XYZPointF(0., 0., 0.));
128  thePositionAtEcal_.push_back(math::XYZPointF(0., 0., 0.));
129  theTrackInnerPosition_.push_back(math::XYZPointF(0., 0., 0.));
130  theTrackInnerPosition_.push_back(math::XYZPointF(0., 0., 0.));
131  theTrackPin_.push_back(math::XYZVectorF(0., 0., 0.));
132  theTrackPin_.push_back(math::XYZVectorF(0., 0., 0.));
133  theTrackPout_.push_back(math::XYZVectorF(0., 0., 0.));
134  theTrackPout_.push_back(math::XYZVectorF(0., 0., 0.));
135 }
std::vector< math::XYZVectorF > theTrackPout_
P_out of tracks.
Definition: Conversion.h:197
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:10
uint8_t nSharedHits_
number of shared hits between tracks
Definition: Conversion.h:210
uint16_t qualityMask_
Definition: Conversion.h:208
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:16
std::vector< math::XYZPointF > theTrackInnerPosition_
P_in of tracks.
Definition: Conversion.h:193
uint8_t algorithm_
conversion algorithm/provenance
Definition: Conversion.h:212
std::vector< math::XYZVectorF > theTrackPin_
P_in of tracks.
Definition: Conversion.h:195
float theMVAout_
TMVA output.
Definition: Conversion.h:207
float theMinDistOfApproach_
Distance of min approach of the two tracks.
Definition: Conversion.h:205
std::vector< math::XYZPointF > thePositionAtEcal_
position at the ECAl surface of the track extrapolation
Definition: Conversion.h:187

◆ Conversion() [2/5]

Conversion::Conversion ( const reco::CaloClusterPtrVector clu,
const std::vector< edm::RefToBase< reco::Track > > &  tr,
const std::vector< math::XYZPointF > &  trackPositionAtEcal,
const reco::Vertex convVtx,
const std::vector< reco::CaloClusterPtr > &  matchingBC,
const float  DCA,
const std::vector< math::XYZPointF > &  innPoint,
const std::vector< math::XYZVectorF > &  trackPin,
const std::vector< math::XYZVectorF > &  trackPout,
const std::vector< uint8_t > &  nHitsBeforeVtx,
const std::vector< Measurement1DFloat > &  dlClosestHitToVtx,
uint8_t  nSharedHits,
const float  mva,
ConversionAlgorithm  algo = undefined 
)

Definition at line 40 of file Conversion.cc.

54  :
55 
56  caloCluster_(sc),
57  trackToBaseRefs_(tr),
58  thePositionAtEcal_(trackPositionAtEcal),
59  theConversionVertex_(convVtx),
60  theMatchingBCs_(matchingBC),
61  theTrackInnerPosition_(innPoint),
62  theTrackPin_(trackPin),
63  theTrackPout_(trackPout),
67  theMVAout_(mva),
68  qualityMask_(0),
70  algorithm_(algo) {}
const std::vector< Measurement1DFloat > & dlClosestHitToVtx() const
Vector of signed decay length with uncertainty from nearest hit on track to the conversion vtx positi...
Definition: Conversion.h:156
std::vector< math::XYZVectorF > theTrackPout_
P_out of tracks.
Definition: Conversion.h:197
reco::CaloClusterPtrVector caloCluster_
vector pointer to a/multiple seed CaloCluster(s)
Definition: Conversion.h:183
uint8_t nSharedHits() const
number of shared hits btw the two track
Definition: Conversion.h:158
ConversionAlgorithm algo() const
Definition: Conversion.h:215
uint8_t nSharedHits_
number of shared hits between tracks
Definition: Conversion.h:210
uint16_t qualityMask_
Definition: Conversion.h:208
std::vector< math::XYZPointF > theTrackInnerPosition_
P_in of tracks.
Definition: Conversion.h:193
uint8_t algorithm_
conversion algorithm/provenance
Definition: Conversion.h:212
std::vector< math::XYZVectorF > theTrackPin_
P_in of tracks.
Definition: Conversion.h:195
float theMVAout_
TMVA output.
Definition: Conversion.h:207
std::vector< reco::CaloClusterPtr > theMatchingBCs_
Clusters mathing the tracks (these are not the seeds)
Definition: Conversion.h:191
float theMinDistOfApproach_
Distance of min approach of the two tracks.
Definition: Conversion.h:205
std::vector< math::XYZPointF > thePositionAtEcal_
position at the ECAl surface of the track extrapolation
Definition: Conversion.h:187
const std::vector< uint8_t > & nHitsBeforeVtx() const
Vector of the number of hits before the vertex along each track trajector.
Definition: Conversion.h:154
std::vector< edm::RefToBase< reco::Track > > trackToBaseRefs_
vector Track RefToBase
Definition: Conversion.h:185
std::vector< Measurement1DFloat > dlClosestHitToVtx_
signed decay length and uncertainty from nearest hit on track to conversion vertex ...
Definition: Conversion.h:201
reco::Vertex theConversionVertex_
Fitted Kalman conversion vertex.
Definition: Conversion.h:189
std::vector< uint8_t > nHitsBeforeVtx_
number of hits before the vertex on each trackerOnly
Definition: Conversion.h:199

◆ Conversion() [3/5]

Conversion::Conversion ( const reco::CaloClusterPtrVector clu,
const std::vector< reco::TrackRef > &  tr,
const std::vector< math::XYZPointF > &  trackPositionAtEcal,
const reco::Vertex convVtx,
const std::vector< reco::CaloClusterPtr > &  matchingBC,
const float  DCA,
const std::vector< math::XYZPointF > &  innPoint,
const std::vector< math::XYZVectorF > &  trackPin,
const std::vector< math::XYZVectorF > &  trackPout,
const float  mva,
ConversionAlgorithm  algo = undefined 
)

Definition at line 8 of file Conversion.cc.

References HLT_2024v13_cff::track, and trackToBaseRefs_.

19  :
20 
21  caloCluster_(sc),
23  thePositionAtEcal_(trackPositionAtEcal),
24  theConversionVertex_(convVtx),
25  theMatchingBCs_(matchingBC),
26  theTrackInnerPosition_(innPoint),
27  theTrackPin_(trackPin),
28  theTrackPout_(trackPout),
30  theMVAout_(mva),
31  qualityMask_(0),
32  nSharedHits_(0),
33  algorithm_(algo) {
34  trackToBaseRefs_.reserve(tr.size());
35  for (auto const& track : tr) {
36  trackToBaseRefs_.emplace_back(track);
37  }
38 }
std::vector< math::XYZVectorF > theTrackPout_
P_out of tracks.
Definition: Conversion.h:197
reco::CaloClusterPtrVector caloCluster_
vector pointer to a/multiple seed CaloCluster(s)
Definition: Conversion.h:183
ConversionAlgorithm algo() const
Definition: Conversion.h:215
uint8_t nSharedHits_
number of shared hits between tracks
Definition: Conversion.h:210
uint16_t qualityMask_
Definition: Conversion.h:208
std::vector< math::XYZPointF > theTrackInnerPosition_
P_in of tracks.
Definition: Conversion.h:193
uint8_t algorithm_
conversion algorithm/provenance
Definition: Conversion.h:212
std::vector< math::XYZVectorF > theTrackPin_
P_in of tracks.
Definition: Conversion.h:195
float theMVAout_
TMVA output.
Definition: Conversion.h:207
std::vector< reco::CaloClusterPtr > theMatchingBCs_
Clusters mathing the tracks (these are not the seeds)
Definition: Conversion.h:191
float theMinDistOfApproach_
Distance of min approach of the two tracks.
Definition: Conversion.h:205
std::vector< math::XYZPointF > thePositionAtEcal_
position at the ECAl surface of the track extrapolation
Definition: Conversion.h:187
std::vector< edm::RefToBase< reco::Track > > trackToBaseRefs_
vector Track RefToBase
Definition: Conversion.h:185
reco::Vertex theConversionVertex_
Fitted Kalman conversion vertex.
Definition: Conversion.h:189

◆ Conversion() [4/5]

Conversion::Conversion ( const reco::CaloClusterPtrVector clu,
const std::vector< reco::TrackRef > &  tr,
const reco::Vertex convVtx,
ConversionAlgorithm  algo = undefined 
)

Definition at line 72 of file Conversion.cc.

References theMinDistOfApproach_, theMVAout_, thePositionAtEcal_, theTrackInnerPosition_, theTrackPin_, theTrackPout_, HLT_2024v13_cff::track, and trackToBaseRefs_.

76  : caloCluster_(sc),
78  theConversionVertex_(convVtx),
79  qualityMask_(0),
80  nSharedHits_(0),
81  algorithm_(algo) {
82  trackToBaseRefs_.reserve(tr.size());
83  for (auto const& track : tr) {
84  trackToBaseRefs_.emplace_back(track);
85  }
86 
87  theMinDistOfApproach_ = 9999.;
88  theMVAout_ = 9999.;
89  thePositionAtEcal_.push_back(math::XYZPointF(0., 0., 0.));
90  thePositionAtEcal_.push_back(math::XYZPointF(0., 0., 0.));
91  theTrackInnerPosition_.push_back(math::XYZPointF(0., 0., 0.));
92  theTrackInnerPosition_.push_back(math::XYZPointF(0., 0., 0.));
93  theTrackPin_.push_back(math::XYZVectorF(0., 0., 0.));
94  theTrackPin_.push_back(math::XYZVectorF(0., 0., 0.));
95  theTrackPout_.push_back(math::XYZVectorF(0., 0., 0.));
96  theTrackPout_.push_back(math::XYZVectorF(0., 0., 0.));
97 }
std::vector< math::XYZVectorF > theTrackPout_
P_out of tracks.
Definition: Conversion.h:197
reco::CaloClusterPtrVector caloCluster_
vector pointer to a/multiple seed CaloCluster(s)
Definition: Conversion.h:183
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:10
ConversionAlgorithm algo() const
Definition: Conversion.h:215
uint8_t nSharedHits_
number of shared hits between tracks
Definition: Conversion.h:210
uint16_t qualityMask_
Definition: Conversion.h:208
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:16
std::vector< math::XYZPointF > theTrackInnerPosition_
P_in of tracks.
Definition: Conversion.h:193
uint8_t algorithm_
conversion algorithm/provenance
Definition: Conversion.h:212
std::vector< math::XYZVectorF > theTrackPin_
P_in of tracks.
Definition: Conversion.h:195
float theMVAout_
TMVA output.
Definition: Conversion.h:207
float theMinDistOfApproach_
Distance of min approach of the two tracks.
Definition: Conversion.h:205
std::vector< math::XYZPointF > thePositionAtEcal_
position at the ECAl surface of the track extrapolation
Definition: Conversion.h:187
std::vector< edm::RefToBase< reco::Track > > trackToBaseRefs_
vector Track RefToBase
Definition: Conversion.h:185
reco::Vertex theConversionVertex_
Fitted Kalman conversion vertex.
Definition: Conversion.h:189

◆ Conversion() [5/5]

Conversion::Conversion ( const reco::CaloClusterPtrVector clu,
const std::vector< edm::RefToBase< reco::Track > > &  tr,
const reco::Vertex convVtx,
ConversionAlgorithm  algo = undefined 
)

Definition at line 99 of file Conversion.cc.

References theMinDistOfApproach_, theMVAout_, thePositionAtEcal_, theTrackInnerPosition_, theTrackPin_, and theTrackPout_.

103  : caloCluster_(sc),
104  trackToBaseRefs_(tr),
105  theConversionVertex_(convVtx),
106  qualityMask_(0),
107  nSharedHits_(0),
108  algorithm_(algo) {
109  theMinDistOfApproach_ = 9999.;
110  theMVAout_ = 9999.;
111  thePositionAtEcal_.push_back(math::XYZPointF(0., 0., 0.));
112  thePositionAtEcal_.push_back(math::XYZPointF(0., 0., 0.));
113  theTrackInnerPosition_.push_back(math::XYZPointF(0., 0., 0.));
114  theTrackInnerPosition_.push_back(math::XYZPointF(0., 0., 0.));
115  theTrackPin_.push_back(math::XYZVectorF(0., 0., 0.));
116  theTrackPin_.push_back(math::XYZVectorF(0., 0., 0.));
117  theTrackPout_.push_back(math::XYZVectorF(0., 0., 0.));
118  theTrackPout_.push_back(math::XYZVectorF(0., 0., 0.));
119 }
std::vector< math::XYZVectorF > theTrackPout_
P_out of tracks.
Definition: Conversion.h:197
reco::CaloClusterPtrVector caloCluster_
vector pointer to a/multiple seed CaloCluster(s)
Definition: Conversion.h:183
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:10
ConversionAlgorithm algo() const
Definition: Conversion.h:215
uint8_t nSharedHits_
number of shared hits between tracks
Definition: Conversion.h:210
uint16_t qualityMask_
Definition: Conversion.h:208
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:16
std::vector< math::XYZPointF > theTrackInnerPosition_
P_in of tracks.
Definition: Conversion.h:193
uint8_t algorithm_
conversion algorithm/provenance
Definition: Conversion.h:212
std::vector< math::XYZVectorF > theTrackPin_
P_in of tracks.
Definition: Conversion.h:195
float theMVAout_
TMVA output.
Definition: Conversion.h:207
float theMinDistOfApproach_
Distance of min approach of the two tracks.
Definition: Conversion.h:205
std::vector< math::XYZPointF > thePositionAtEcal_
position at the ECAl surface of the track extrapolation
Definition: Conversion.h:187
std::vector< edm::RefToBase< reco::Track > > trackToBaseRefs_
vector Track RefToBase
Definition: Conversion.h:185
reco::Vertex theConversionVertex_
Fitted Kalman conversion vertex.
Definition: Conversion.h:189

Member Function Documentation

◆ algo()

Conversion::ConversionAlgorithm reco::Conversion::algo ( ) const
inline

Definition at line 215 of file Conversion.h.

References algorithm_.

215 { return (ConversionAlgorithm)algorithm_; }
uint8_t algorithm_
conversion algorithm/provenance
Definition: Conversion.h:212

◆ algoByName()

Conversion::ConversionAlgorithm Conversion::algoByName ( const std::string &  name)
static

Definition at line 139 of file Conversion.cc.

References algoNames, algoSize, spr::find(), Skims_PA_cff::name, findQualityFiles::size, and undefined.

Referenced by ConversionProducer::buildCollection(), and ConvertedPhotonProducer::buildCollections().

139  {
142  if (index == size)
143  return undefined;
144 
145  return ConversionAlgorithm(index);
146 }
size
Write out results.
static const std::string algoNames[]
Definition: Conversion.h:39
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19

◆ algoName() [1/2]

std::string reco::Conversion::algoName ( ) const
inline

Definition at line 217 of file Conversion.h.

References algorithm_, ecalSeeded, mixed, pflow, trackerOnly, and undefined.

217  {
218  switch (algorithm_) {
219  case undefined:
220  return "undefined";
221  case ecalSeeded:
222  return "ecalSeeded";
223  case trackerOnly:
224  return "trackerOnly";
225  case mixed:
226  return "mixed";
227  case pflow:
228  return "pflow";
229  }
230  return "undefined";
231  }
uint8_t algorithm_
conversion algorithm/provenance
Definition: Conversion.h:212

◆ algoName() [2/2]

std::string reco::Conversion::algoName ( ConversionAlgorithm  a)
inlinestatic

Definition at line 233 of file Conversion.h.

References a, algoNames, algoSize, and createfilelist::int.

233  {
234  if (int(a) < int(algoSize) && int(a) > 0)
235  return algoNames[int(a)];
236  return "undefined";
237  }
static const std::string algoNames[]
Definition: Conversion.h:39
double a
Definition: hdecay.h:121

◆ bcMatchingWithTracks()

const std::vector<reco::CaloClusterPtr>& reco::Conversion::bcMatchingWithTracks ( ) const
inline

Definition at line 144 of file Conversion.h.

References theMatchingBCs_.

Referenced by dEtaTracksAtEcal(), and dPhiTracksAtEcal().

144 { return theMatchingBCs_; }
std::vector< reco::CaloClusterPtr > theMatchingBCs_
Clusters mathing the tracks (these are not the seeds)
Definition: Conversion.h:191

◆ caloCluster()

reco::CaloClusterPtrVector reco::Conversion::caloCluster ( ) const
inline

Pointer to CaloCluster (foe Egamma Conversions it points to a SuperCluster)

Definition at line 84 of file Conversion.h.

References caloCluster_.

Referenced by EoverP(), and EoverPrefittedTracks().

84 { return caloCluster_; }
reco::CaloClusterPtrVector caloCluster_
vector pointer to a/multiple seed CaloCluster(s)
Definition: Conversion.h:183

◆ clone()

Conversion * Conversion::clone ( void  ) const

returns a clone of the candidate

Definition at line 148 of file Conversion.cc.

References Conversion().

Referenced by ConvertedPhotonProducer::cleanCollections().

148 { return new Conversion(*this); }

◆ conversionVertex()

const reco::Vertex& reco::Conversion::conversionVertex ( ) const
inline

returns the reco conversion vertex

Definition at line 88 of file Conversion.h.

References theConversionVertex_.

Referenced by TkConvValidator::analyze(), ConversionLessByChi2(), dxy(), dz(), lxy(), lz(), refittedPair4Momentum(), and refittedPairMomentum().

88 { return theConversionVertex_; }
reco::Vertex theConversionVertex_
Fitted Kalman conversion vertex.
Definition: Conversion.h:189

◆ dEtaTracksAtEcal()

double Conversion::dEtaTracksAtEcal ( ) const

Definition at line 303 of file Conversion.cc.

References bcMatchingWithTracks(), ecalImpactPosition(), nTracks(), and mps_fire::result.

303  {
304  double result = -99.;
305 
306  if (nTracks() == 2 && bcMatchingWithTracks()[0].isNonnull() && bcMatchingWithTracks()[1].isNonnull()) {
307  result = ecalImpactPosition()[0].eta() - ecalImpactPosition()[1].eta();
308  }
309 
310  return result;
311 }
const std::vector< math::XYZPointF > & ecalImpactPosition() const
Definition: Conversion.h:142
const std::vector< reco::CaloClusterPtr > & bcMatchingWithTracks() const
Definition: Conversion.h:144
unsigned int nTracks() const
Number of tracks= 0,1,2.
Definition: Conversion.h:92

◆ distOfMinimumApproach()

double reco::Conversion::distOfMinimumApproach ( ) const
inline

Definition at line 116 of file Conversion.h.

References theMinDistOfApproach_.

Referenced by TkConvValidator::analyze().

116 { return theMinDistOfApproach_; }
float theMinDistOfApproach_
Distance of min approach of the two tracks.
Definition: Conversion.h:205

◆ dlClosestHitToVtx()

const std::vector<Measurement1DFloat>& reco::Conversion::dlClosestHitToVtx ( ) const
inline

Vector of signed decay length with uncertainty from nearest hit on track to the conversion vtx positions.

Definition at line 156 of file Conversion.h.

References dlClosestHitToVtx_.

Referenced by TkConvValidator::analyze().

156 { return dlClosestHitToVtx_; }
std::vector< Measurement1DFloat > dlClosestHitToVtx_
signed decay length and uncertainty from nearest hit on track to conversion vertex ...
Definition: Conversion.h:201

◆ dPhiTracksAtEcal()

double Conversion::dPhiTracksAtEcal ( ) const

Definition at line 270 of file Conversion.cc.

References bcMatchingWithTracks(), ecalImpactPosition(), nTracks(), pi, and mps_fire::result.

270  {
271  double result = -99.;
272 
273  if (nTracks() == 2 && bcMatchingWithTracks()[0].isNonnull() && bcMatchingWithTracks()[1].isNonnull()) {
274  float recoPhi1 = ecalImpactPosition()[0].phi();
275  if (recoPhi1 > pi) {
276  recoPhi1 = recoPhi1 - twopi;
277  }
278  if (recoPhi1 < -pi) {
279  recoPhi1 = recoPhi1 + twopi;
280  }
281 
282  float recoPhi2 = ecalImpactPosition()[1].phi();
283  if (recoPhi2 > pi) {
284  recoPhi2 = recoPhi2 - twopi;
285  }
286  if (recoPhi2 < -pi) {
287  recoPhi2 = recoPhi2 + twopi;
288  }
289 
290  result = recoPhi1 - recoPhi2;
291 
292  if (result > pi) {
293  result = result - twopi;
294  }
295  if (result < -pi) {
296  result = result + twopi;
297  }
298  }
299 
300  return result;
301 }
const std::vector< math::XYZPointF > & ecalImpactPosition() const
Definition: Conversion.h:142
const Double_t pi
const std::vector< reco::CaloClusterPtr > & bcMatchingWithTracks() const
Definition: Conversion.h:144
unsigned int nTracks() const
Number of tracks= 0,1,2.
Definition: Conversion.h:92

◆ dPhiTracksAtVtx()

double Conversion::dPhiTracksAtVtx ( ) const

Definition at line 255 of file Conversion.cc.

References nTracks(), pi, mps_fire::result, and tracksPin().

255  {
256  double result = -99;
257  if (nTracks() == 2) {
258  result = tracksPin()[0].phi() - tracksPin()[1].phi();
259  if (result > pi) {
260  result = result - twopi;
261  }
262  if (result < -pi) {
263  result = result + twopi;
264  }
265  }
266 
267  return result;
268 }
const Double_t pi
const std::vector< math::XYZVectorF > & tracksPin() const
Vector of track momentum measured at the innermost hit.
Definition: Conversion.h:152
unsigned int nTracks() const
Number of tracks= 0,1,2.
Definition: Conversion.h:92

◆ dxy()

double Conversion::dxy ( const math::XYZPoint myBeamSpot = math::XYZPoint()) const

Definition at line 313 of file Conversion.cc.

References conversionVertex(), refittedPairMomentum(), and L1BJetProducer_cff::vtx.

Referenced by Electron.Electron::cutBasedId(), and ntupleDataFormat.Track::dxyPull().

313  {
314  const reco::Vertex& vtx = conversionVertex();
315  if (!vtx.isValid())
316  return -9999.;
317 
319 
320  double dxy = (-(vtx.x() - myBeamSpot.x()) * mom.y() + (vtx.y() - myBeamSpot.y()) * mom.x()) / mom.rho();
321  return dxy;
322 }
double dxy(const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
Definition: Conversion.cc:313
math::XYZVectorF refittedPairMomentum() const
Conversion tracks momentum from the tracks refitted with vertex constraint.
Definition: Conversion.cc:206
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:16
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:88

◆ dz()

double Conversion::dz ( const math::XYZPoint myBeamSpot = math::XYZPoint()) const

Definition at line 324 of file Conversion.cc.

References conversionVertex(), refittedPairMomentum(), and L1BJetProducer_cff::vtx.

Referenced by Electron.Electron::cutBasedId(), ntupleDataFormat.Track::dzPull(), and zOfPrimaryVertexFromTracks().

324  {
325  const reco::Vertex& vtx = conversionVertex();
326  if (!vtx.isValid())
327  return -9999.;
328 
330 
331  double dz =
332  (vtx.z() - myBeamSpot.z()) -
333  ((vtx.x() - myBeamSpot.x()) * mom.x() + (vtx.y() - myBeamSpot.y()) * mom.y()) / mom.rho() * mom.z() / mom.rho();
334  return dz;
335 }
math::XYZVectorF refittedPairMomentum() const
Conversion tracks momentum from the tracks refitted with vertex constraint.
Definition: Conversion.cc:206
double dz(const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
Definition: Conversion.cc:324
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:16
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:88

◆ ecalImpactPosition()

const std::vector<math::XYZPointF>& reco::Conversion::ecalImpactPosition ( ) const
inline

The following are variables provided per each track positions of the track extrapolation at the ECAL front face

Definition at line 142 of file Conversion.h.

References thePositionAtEcal_.

Referenced by dEtaTracksAtEcal(), and dPhiTracksAtEcal().

142 { return thePositionAtEcal_; }
std::vector< math::XYZPointF > thePositionAtEcal_
position at the ECAl surface of the track extrapolation
Definition: Conversion.h:187

◆ EoverP()

double Conversion::EoverP ( ) const

Super Cluster energy divided by track pair momentum if Standard seeding method. If a pointer to two (or more clusters) is stored in the conversion, this method returns the energy sum of clusters divided by the track pair momentum Track innermost momentum is used here

Definition at line 213 of file Conversion.cc.

References caloCluster(), SiStripBadComponentsDQMServiceTemplate_cfg::ep, mps_fire::i, nTracks(), pairMomentum(), edm::PtrVectorBase::size(), findQualityFiles::size, and mathSSE::sqrt().

213  {
214  double ep = -99.;
215 
216  if (nTracks() > 0) {
217  unsigned int size = this->caloCluster().size();
218  float etot = 0.;
219  for (unsigned int i = 0; i < size; i++) {
220  etot += caloCluster()[i]->energy();
221  }
222  if (this->pairMomentum().Mag2() != 0)
223  ep = etot / sqrt(this->pairMomentum().Mag2());
224  }
225 
226  return ep;
227 }
size
Write out results.
size_type size() const
Size of the RefVector.
Definition: PtrVectorBase.h:75
math::XYZVectorF pairMomentum() const
Conversion tracks momentum from the tracks inner momentum.
Definition: Conversion.cc:191
T sqrt(T t)
Definition: SSEVec.h:19
unsigned int nTracks() const
Number of tracks= 0,1,2.
Definition: Conversion.h:92
reco::CaloClusterPtrVector caloCluster() const
Pointer to CaloCluster (foe Egamma Conversions it points to a SuperCluster)
Definition: Conversion.h:84

◆ EoverPrefittedTracks()

double Conversion::EoverPrefittedTracks ( ) const

Super Cluster energy divided by track pair momentum if Standard seeing method. If a pointer to two (or more clusters) is stored in the conversion, this method returns the energy sum of clusters divided by the track pair momentum Track momentum refitted with vertex constraint is used

Definition at line 229 of file Conversion.cc.

References caloCluster(), SiStripBadComponentsDQMServiceTemplate_cfg::ep, mps_fire::i, nTracks(), refittedPairMomentum(), edm::PtrVectorBase::size(), findQualityFiles::size, and mathSSE::sqrt().

229  {
230  double ep = -99.;
231 
232  if (nTracks() > 0) {
233  unsigned int size = this->caloCluster().size();
234  float etot = 0.;
235  for (unsigned int i = 0; i < size; i++) {
236  etot += caloCluster()[i]->energy();
237  }
238  if (this->refittedPairMomentum().Mag2() != 0)
239  ep = etot / sqrt(this->refittedPairMomentum().Mag2());
240  }
241 
242  return ep;
243 }
size
Write out results.
size_type size() const
Size of the RefVector.
Definition: PtrVectorBase.h:75
math::XYZVectorF refittedPairMomentum() const
Conversion tracks momentum from the tracks refitted with vertex constraint.
Definition: Conversion.cc:206
T sqrt(T t)
Definition: SSEVec.h:19
unsigned int nTracks() const
Number of tracks= 0,1,2.
Definition: Conversion.h:92
reco::CaloClusterPtrVector caloCluster() const
Pointer to CaloCluster (foe Egamma Conversions it points to a SuperCluster)
Definition: Conversion.h:84

◆ isConverted()

bool Conversion::isConverted ( ) const

Bool flagging objects having track size >0.

Definition at line 152 of file Conversion.cc.

References nTracks().

152  {
153  if (this->nTracks() == 2)
154  return true;
155  else
156  return false;
157 }
unsigned int nTracks() const
Number of tracks= 0,1,2.
Definition: Conversion.h:92

◆ lxy()

double Conversion::lxy ( const math::XYZPoint myBeamSpot = math::XYZPoint()) const

Definition at line 337 of file Conversion.cc.

References conversionVertex(), refittedPairMomentum(), and L1BJetProducer_cff::vtx.

337  {
338  const reco::Vertex& vtx = conversionVertex();
339  if (!vtx.isValid())
340  return -9999.;
341 
343 
344  double dbsx = vtx.x() - myBeamSpot.x();
345  double dbsy = vtx.y() - myBeamSpot.y();
346  double lxy = (mom.x() * dbsx + mom.y() * dbsy) / mom.rho();
347  return lxy;
348 }
math::XYZVectorF refittedPairMomentum() const
Conversion tracks momentum from the tracks refitted with vertex constraint.
Definition: Conversion.cc:206
double lxy(const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
Definition: Conversion.cc:337
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:16
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:88

◆ lz()

double Conversion::lz ( const math::XYZPoint myBeamSpot = math::XYZPoint()) const

Definition at line 350 of file Conversion.cc.

References funct::abs(), conversionVertex(), refittedPairMomentum(), and L1BJetProducer_cff::vtx.

350  {
351  const reco::Vertex& vtx = conversionVertex();
352  if (!vtx.isValid())
353  return -9999.;
354 
356 
357  double lz = (vtx.z() - myBeamSpot.z()) * mom.z() / std::abs(mom.z());
358  return lz;
359 }
math::XYZVectorF refittedPairMomentum() const
Conversion tracks momentum from the tracks refitted with vertex constraint.
Definition: Conversion.cc:206
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:16
double lz(const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
Definition: Conversion.cc:350
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:88

◆ MVAout()

double reco::Conversion::MVAout ( ) const
inline

get the value of the TMVA output

Definition at line 94 of file Conversion.h.

References theMVAout_.

94 { return theMVAout_; }
float theMVAout_
TMVA output.
Definition: Conversion.h:207

◆ nHitsBeforeVtx()

const std::vector<uint8_t>& reco::Conversion::nHitsBeforeVtx ( ) const
inline

Vector of the number of hits before the vertex along each track trajector.

Definition at line 154 of file Conversion.h.

References nHitsBeforeVtx_.

Referenced by TkConvValidator::analyze().

154 { return nHitsBeforeVtx_; }
std::vector< uint8_t > nHitsBeforeVtx_
number of hits before the vertex on each trackerOnly
Definition: Conversion.h:199

◆ nSharedHits()

uint8_t reco::Conversion::nSharedHits ( ) const
inline

number of shared hits btw the two track

Definition at line 158 of file Conversion.h.

References nSharedHits_.

Referenced by TkConvValidator::analyze().

158 { return nSharedHits_; }
uint8_t nSharedHits_
number of shared hits between tracks
Definition: Conversion.h:210

◆ nTracks()

unsigned int reco::Conversion::nTracks ( ) const
inline

Number of tracks= 0,1,2.

Definition at line 92 of file Conversion.h.

References tracks().

Referenced by dEtaTracksAtEcal(), dPhiTracksAtEcal(), dPhiTracksAtVtx(), EoverP(), EoverPrefittedTracks(), isConverted(), pairCotThetaSeparation(), pairInvariantMass(), and pairMomentum().

92 { return tracks().size(); }
std::vector< edm::RefToBase< reco::Track > > const & tracks() const
vector of track to base references
Definition: Conversion.cc:150

◆ oneLegMVA()

std::vector<float> const reco::Conversion::oneLegMVA ( )
inline

get the MVS output from PF for one leg conversions

Definition at line 96 of file Conversion.h.

References theOneLegMVA_.

96 { return theOneLegMVA_; }
std::vector< float > theOneLegMVA_
vectors of TMVA outputs from pflow for one leg conversions
Definition: Conversion.h:203

◆ pairCotThetaSeparation()

double Conversion::pairCotThetaSeparation ( ) const

Delta cot(Theta) where Theta is the angle in the (y,z) plane between the two tracks. Original tracks are used.

Definition at line 179 of file Conversion.cc.

References nTracks(), funct::tan(), and tracksPin().

Referenced by TkConvValidator::analyze().

179  {
180  double dCotTheta = -99.;
181 
182  if (nTracks() == 2) {
183  double theta1 = this->tracksPin()[0].Theta();
184  double theta2 = this->tracksPin()[1].Theta();
185  dCotTheta = 1. / tan(theta1) - 1. / tan(theta2);
186  }
187 
188  return dCotTheta;
189 }
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
const std::vector< math::XYZVectorF > & tracksPin() const
Vector of track momentum measured at the innermost hit.
Definition: Conversion.h:152
unsigned int nTracks() const
Number of tracks= 0,1,2.
Definition: Conversion.h:92

◆ pairInvariantMass()

double Conversion::pairInvariantMass ( ) const

if nTracks=2 returns the pair invariant mass. Original tracks are used here

Definition at line 159 of file Conversion.cc.

References MillePedeFileConverter_cfg::e, nTracks(), multPhiCorr_741_25nsDY_cfi::px, multPhiCorr_741_25nsDY_cfi::py, mathSSE::sqrt(), and tracksPin().

Referenced by TkConvValidator::analyze().

159  {
160  double invMass = -99.;
161  const float mElec = 0.000511;
162  if (nTracks() == 2) {
163  double px = tracksPin()[0].x() + tracksPin()[1].x();
164  double py = tracksPin()[0].y() + tracksPin()[1].y();
165  double pz = tracksPin()[0].z() + tracksPin()[1].z();
166  double mom1 = tracksPin()[0].Mag2();
167  double mom2 = tracksPin()[1].Mag2();
168  double e = sqrt(mom1 + mElec * mElec) + sqrt(mom2 + mElec * mElec);
169  invMass = (e * e - px * px - py * py - pz * pz);
170  if (invMass > 0)
171  invMass = sqrt(invMass);
172  else
173  invMass = -1;
174  }
175 
176  return invMass;
177 }
T sqrt(T t)
Definition: SSEVec.h:19
const std::vector< math::XYZVectorF > & tracksPin() const
Vector of track momentum measured at the innermost hit.
Definition: Conversion.h:152
unsigned int nTracks() const
Number of tracks= 0,1,2.
Definition: Conversion.h:92

◆ pairMomentum()

math::XYZVectorF Conversion::pairMomentum ( ) const

Conversion tracks momentum from the tracks inner momentum.

Definition at line 191 of file Conversion.cc.

References nTracks(), and tracksPin().

Referenced by EoverP().

191  {
192  if (nTracks() == 2) {
193  return this->tracksPin()[0] + this->tracksPin()[1];
194  }
195  return math::XYZVectorF(0., 0., 0.);
196 }
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:16
const std::vector< math::XYZVectorF > & tracksPin() const
Vector of track momentum measured at the innermost hit.
Definition: Conversion.h:152
unsigned int nTracks() const
Number of tracks= 0,1,2.
Definition: Conversion.h:92

◆ quality()

bool reco::Conversion::quality ( ConversionQuality  q) const
inline

Definition at line 178 of file Conversion.h.

References submitPVResolutionJobs::q, and qualityMask_.

Referenced by TkConvValidator::analyze().

178 { return (qualityMask_ & (1 << q)) >> q; }
uint16_t qualityMask_
Definition: Conversion.h:208

◆ refittedPair4Momentum()

math::XYZTLorentzVectorF Conversion::refittedPair4Momentum ( ) const

Conversion track pair 4-momentum from the tracks refitted with vertex constraint.

Definition at line 198 of file Conversion.cc.

References conversionVertex(), sistrip::SpyUtilities::isValid(), and reco::Vertex::p4().

Referenced by refittedPairMomentum().

198  {
200  if (this->conversionVertex().isValid())
201  p4 = this->conversionVertex().p4(0.000511, 0.5);
202 
203  return p4;
204 }
const bool isValid(const Frame &aFrame, const FrameQuality &aQuality, const uint16_t aExpectedPos)
math::XYZTLorentzVectorD p4(float mass=0.13957018, float minWeight=0.5) const
Returns the four momentum of the sum of the tracks, assuming the given mass for the decay products...
Definition: Vertex.cc:103
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:88
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:22

◆ refittedPairMomentum()

math::XYZVectorF Conversion::refittedPairMomentum ( ) const

Conversion tracks momentum from the tracks refitted with vertex constraint.

Definition at line 206 of file Conversion.cc.

References conversionVertex(), sistrip::SpyUtilities::isValid(), and refittedPair4Momentum().

Referenced by TkConvValidator::analyze(), dxy(), dz(), EoverPrefittedTracks(), lxy(), lz(), and ConversionProducer::matchingSC().

206  {
207  if (this->conversionVertex().isValid()) {
208  return this->refittedPair4Momentum().Vect();
209  }
210  return math::XYZVectorF(0., 0., 0.);
211 }
const bool isValid(const Frame &aFrame, const FrameQuality &aQuality, const uint16_t aExpectedPos)
math::XYZTLorentzVectorF refittedPair4Momentum() const
Conversion track pair 4-momentum from the tracks refitted with vertex constraint. ...
Definition: Conversion.cc:198
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:16
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:88

◆ setConversionAlgorithm()

void reco::Conversion::setConversionAlgorithm ( const ConversionAlgorithm  a,
bool  set = true 
)
inline

Conversion Track algorithm/provenance.

Definition at line 167 of file Conversion.h.

References a, algorithm_, and undefined.

167  {
168  if (set)
169  algorithm_ = a;
170  else
172  }
uint8_t algorithm_
conversion algorithm/provenance
Definition: Conversion.h:212
double a
Definition: hdecay.h:121

◆ setMatchingSuperCluster()

void reco::Conversion::setMatchingSuperCluster ( const reco::CaloClusterPtrVector sc)
inline

Definition at line 165 of file Conversion.h.

References caloCluster_.

Referenced by ConversionProducer::buildCollection().

165 { caloCluster_ = sc; }
reco::CaloClusterPtrVector caloCluster_
vector pointer to a/multiple seed CaloCluster(s)
Definition: Conversion.h:183

◆ setMVAout()

void reco::Conversion::setMVAout ( const float &  mva)
inline

set the value of the TMVA output

Definition at line 161 of file Conversion.h.

References beam_dqm_sourceclient-live_cfg::mva, and theMVAout_.

Referenced by ConvertedPhotonProducer::buildCollections().

◆ setOneLegMVA()

void reco::Conversion::setOneLegMVA ( const std::vector< float > &  mva)
inline

set the MVS output from PF for one leg conversions

Definition at line 163 of file Conversion.h.

References beam_dqm_sourceclient-live_cfg::mva, and theOneLegMVA_.

Referenced by PFEGammaProducer::createSingleLegConversions().

163 { theOneLegMVA_ = mva; }
std::vector< float > theOneLegMVA_
vectors of TMVA outputs from pflow for one leg conversions
Definition: Conversion.h:203

◆ setQuality()

void reco::Conversion::setQuality ( ConversionQuality  q,
bool  b 
)
inline

Definition at line 239 of file Conversion.h.

References b, submitPVResolutionJobs::q, and qualityMask_.

Referenced by ConversionProducer::buildCollection().

239  {
240  if (b) //regular OR if setting value to true
241  qualityMask_ |= (1 << q);
242  else // doing "half-XOR" if unsetting value
243  qualityMask_ &= (~(1 << q));
244  }
uint16_t qualityMask_
Definition: Conversion.h:208
double b
Definition: hdecay.h:120

◆ tracks()

std::vector< edm::RefToBase< reco::Track > > const & Conversion::tracks ( void  ) const

vector of track to base references

Definition at line 150 of file Conversion.cc.

References trackToBaseRefs_.

Referenced by TkConvValidator::analyze(), ConversionEqualByTrack(), nTracks(), reco::operator==(), and tracksSigned_d0().

150 { return trackToBaseRefs_; }
std::vector< edm::RefToBase< reco::Track > > trackToBaseRefs_
vector Track RefToBase
Definition: Conversion.h:185

◆ tracksInnerPosition()

const std::vector<math::XYZPointF>& reco::Conversion::tracksInnerPosition ( ) const
inline

Vector containing the position of the innermost hit of each track.

Definition at line 148 of file Conversion.h.

References theTrackInnerPosition_.

148 { return theTrackInnerPosition_; }
std::vector< math::XYZPointF > theTrackInnerPosition_
P_in of tracks.
Definition: Conversion.h:193

◆ tracksPin()

const std::vector<math::XYZVectorF>& reco::Conversion::tracksPin ( ) const
inline

Vector of track momentum measured at the innermost hit.

Definition at line 152 of file Conversion.h.

References theTrackPin_.

Referenced by dPhiTracksAtVtx(), pairCotThetaSeparation(), pairInvariantMass(), and pairMomentum().

152 { return theTrackPin_; }
std::vector< math::XYZVectorF > theTrackPin_
P_in of tracks.
Definition: Conversion.h:195

◆ tracksPout()

const std::vector<math::XYZVectorF>& reco::Conversion::tracksPout ( ) const
inline

Vector of track momentum measured at the outermost hit.

Definition at line 150 of file Conversion.h.

References theTrackPout_.

150 { return theTrackPout_; }
std::vector< math::XYZVectorF > theTrackPout_
P_out of tracks.
Definition: Conversion.h:197

◆ tracksSigned_d0()

std::vector< double > Conversion::tracksSigned_d0 ( ) const

signed transverse impact parameter for each track

Definition at line 245 of file Conversion.cc.

References mps_fire::result, findQualityFiles::size, HLT_2024v13_cff::track, and tracks().

245  {
246  std::vector<double> result;
247  result.reserve(tracks().size());
248 
249  for (auto const& track : tracks()) {
250  result.emplace_back(track->d0() * track->charge());
251  }
252  return result;
253 }
size
Write out results.
std::vector< edm::RefToBase< reco::Track > > const & tracks() const
vector of track to base references
Definition: Conversion.cc:150

◆ zOfPrimaryVertexFromTracks()

double reco::Conversion::zOfPrimaryVertexFromTracks ( const math::XYZPoint myBeamSpot = math::XYZPoint()) const
inline

Definition at line 136 of file Conversion.h.

References dz().

Referenced by TkConvValidator::analyze(), and ConversionProducer::matchingSC().

136  {
137  return dz(myBeamSpot) + myBeamSpot.z();
138  }
double dz(const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
Definition: Conversion.cc:324

Member Data Documentation

◆ algoNames

std::string const Conversion::algoNames = {"undefined", "ecalSeeded", "trackerOnly", "mixed", "pflow"}
static

Definition at line 39 of file Conversion.h.

Referenced by algoByName(), and algoName().

◆ algorithm_

uint8_t reco::Conversion::algorithm_
private

conversion algorithm/provenance

Definition at line 212 of file Conversion.h.

Referenced by algo(), algoName(), Conversion(), and setConversionAlgorithm().

◆ caloCluster_

reco::CaloClusterPtrVector reco::Conversion::caloCluster_
private

vector pointer to a/multiple seed CaloCluster(s)

Definition at line 183 of file Conversion.h.

Referenced by caloCluster(), and setMatchingSuperCluster().

◆ dlClosestHitToVtx_

std::vector<Measurement1DFloat> reco::Conversion::dlClosestHitToVtx_
private

signed decay length and uncertainty from nearest hit on track to conversion vertex

Definition at line 201 of file Conversion.h.

Referenced by dlClosestHitToVtx().

◆ nHitsBeforeVtx_

std::vector<uint8_t> reco::Conversion::nHitsBeforeVtx_
private

number of hits before the vertex on each trackerOnly

Definition at line 199 of file Conversion.h.

Referenced by nHitsBeforeVtx().

◆ nSharedHits_

uint8_t reco::Conversion::nSharedHits_
private

number of shared hits between tracks

Definition at line 210 of file Conversion.h.

Referenced by Conversion(), and nSharedHits().

◆ qualityMask_

uint16_t reco::Conversion::qualityMask_
private

Definition at line 208 of file Conversion.h.

Referenced by Conversion(), quality(), and setQuality().

◆ theConversionVertex_

reco::Vertex reco::Conversion::theConversionVertex_
private

Fitted Kalman conversion vertex.

Definition at line 189 of file Conversion.h.

Referenced by conversionVertex().

◆ theMatchingBCs_

std::vector<reco::CaloClusterPtr> reco::Conversion::theMatchingBCs_
private

Clusters mathing the tracks (these are not the seeds)

Definition at line 191 of file Conversion.h.

Referenced by bcMatchingWithTracks().

◆ theMinDistOfApproach_

float reco::Conversion::theMinDistOfApproach_
private

Distance of min approach of the two tracks.

Definition at line 205 of file Conversion.h.

Referenced by Conversion(), and distOfMinimumApproach().

◆ theMVAout_

float reco::Conversion::theMVAout_
private

TMVA output.

Definition at line 207 of file Conversion.h.

Referenced by Conversion(), MVAout(), and setMVAout().

◆ theOneLegMVA_

std::vector<float> reco::Conversion::theOneLegMVA_
private

vectors of TMVA outputs from pflow for one leg conversions

Definition at line 203 of file Conversion.h.

Referenced by oneLegMVA(), and setOneLegMVA().

◆ thePositionAtEcal_

std::vector<math::XYZPointF> reco::Conversion::thePositionAtEcal_
private

position at the ECAl surface of the track extrapolation

Definition at line 187 of file Conversion.h.

Referenced by Conversion(), and ecalImpactPosition().

◆ theTrackInnerPosition_

std::vector<math::XYZPointF> reco::Conversion::theTrackInnerPosition_
private

P_in of tracks.

Definition at line 193 of file Conversion.h.

Referenced by Conversion(), and tracksInnerPosition().

◆ theTrackPin_

std::vector<math::XYZVectorF> reco::Conversion::theTrackPin_
private

P_in of tracks.

Definition at line 195 of file Conversion.h.

Referenced by Conversion(), and tracksPin().

◆ theTrackPout_

std::vector<math::XYZVectorF> reco::Conversion::theTrackPout_
private

P_out of tracks.

Definition at line 197 of file Conversion.h.

Referenced by Conversion(), and tracksPout().

◆ trackToBaseRefs_

std::vector<edm::RefToBase<reco::Track> > reco::Conversion::trackToBaseRefs_
private

vector Track RefToBase

Definition at line 185 of file Conversion.h.

Referenced by Conversion(), and tracks().