Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Pages
RecoEgamma
EgammaTools
interface
ConversionTools.h
Go to the documentation of this file.
1
//--------------------------------------------------------------------------------------------------
2
// $Id: ConversionTools.h,v 1.2 2011/05/10 19:27:22 eulisse Exp $
3
//
4
// ConversionTools
5
//
6
// Utility to match electrons/photons/tracks to conversions and perform various conversion
7
// selection criteria.
8
//
9
// Matching to photons is by geometrical match with the SuperCluster (defined by angle between
10
// conversion momentum and vector joining conversion vertex to SuperCluster position)
11
//
12
// Matching to tracks and electrons is by reference.
13
//
14
// Also implemented here is a "conversion-safe electron veto" for photons through the
15
// matchedPromptElectron and hasMatchedPromptElectron functions
16
//
17
//
18
// Authors: J.Bendavid
19
//--------------------------------------------------------------------------------------------------
20
21
#ifndef EgammaTools_ConversionTools_h
22
#define EgammaTools_ConversionTools_h
23
24
#include "
DataFormats/TrackReco/interface/Track.h
"
25
#include "
DataFormats/Math/interface/Point3D.h
"
26
#include "
DataFormats/EgammaCandidates/interface/GsfElectron.h
"
27
#include "
DataFormats/EgammaCandidates/interface/GsfElectronFwd.h
"
28
#include "
DataFormats/EgammaCandidates/interface/Photon.h
"
29
#include "
DataFormats/EgammaCandidates/interface/ConversionFwd.h
"
30
#include "
DataFormats/TrackReco/interface/Track.h
"
31
#include "
DataFormats/TrackReco/interface/TrackExtra.h
"
32
#include "
DataFormats/TrackReco/interface/TrackFwd.h
"
33
#include "
DataFormats/Math/interface/Point3D.h
"
34
#include "
MagneticField/Engine/interface/MagneticField.h
"
35
#include "
DataFormats/GsfTrackReco/interface/GsfTrack.h
"
36
#include "
DataFormats/GsfTrackReco/interface/GsfTrackFwd.h
"
37
class
ConversionTools
38
{
39
public
:
40
ConversionTools
() {}
41
42
43
static
bool
isGoodConversion
(
const
reco::Conversion
&
conv
,
const
math::XYZPoint
&beamspot,
float
lxyMin=2.0,
float
probMin=1
e
-6,
unsigned
int
nHitsBeforeVtxMax=1);
44
45
static
bool
matchesConversion
(
const
reco::GsfElectron
&ele,
const
reco::Conversion
&
conv
,
bool
allowCkfMatch=
true
);
46
static
bool
matchesConversion
(
const
reco::SuperCluster
&sc,
const
reco::Conversion
&
conv
,
float
dRMax = 0.1,
float
dEtaMax = 999.,
float
dPhiMax = 999.);
47
static
bool
matchesConversion
(
const
edm::RefToBase<reco::Track>
&trk,
const
reco::Conversion
&
conv
);
48
static
bool
matchesConversion
(
const
reco::TrackRef
&trk,
const
reco::Conversion
&
conv
);
49
static
bool
matchesConversion
(
const
reco::GsfTrackRef
&trk,
const
reco::Conversion
&
conv
);
50
51
52
static
bool
hasMatchedConversion
(
const
reco::GsfElectron
&ele,
53
const
edm::Handle<reco::ConversionCollection>
&convCol,
const
math::XYZPoint
&beamspot,
bool
allowCkfMatch=
true
,
float
lxyMin=2.0,
float
probMin=1
e
-6,
unsigned
int
nHitsBeforeVtxMax=0);
54
55
static
bool
hasMatchedConversion
(
const
reco::TrackRef
&trk,
56
const
edm::Handle<reco::ConversionCollection>
&convCol,
const
math::XYZPoint
&beamspot,
float
lxyMin=2.0,
float
probMin=1
e
-6,
unsigned
int
nHitsBeforeVtxMax=1);
57
58
static
bool
hasMatchedConversion
(
const
reco::SuperCluster
&sc,
59
const
edm::Handle<reco::ConversionCollection>
&convCol,
const
math::XYZPoint
&beamspot,
float
dRMax = 0.1,
float
dEtaMax = 999.,
float
dPhiMax = 999.,
float
lxyMin=2.0,
float
probMin=1
e
-6,
unsigned
int
nHitsBeforeVtxMax=1);
60
61
62
static
reco::ConversionRef
matchedConversion
(
const
reco::GsfElectron
&ele,
63
const
edm::Handle<reco::ConversionCollection>
&convCol,
const
math::XYZPoint
&beamspot,
bool
allowCkfMatch=
true
,
float
lxyMin=2.0,
float
probMin=1
e
-6,
unsigned
int
nHitsBeforeVtxMax=0);
64
65
static
reco::ConversionRef
matchedConversion
(
const
reco::TrackRef
&trk,
66
const
edm::Handle<reco::ConversionCollection>
&convCol,
const
math::XYZPoint
&beamspot,
float
lxyMin=2.0,
float
probMin=1
e
-6,
unsigned
int
nHitsBeforeVtxMax=1);
67
68
static
reco::ConversionRef
matchedConversion
(
const
reco::SuperCluster
&sc,
69
const
edm::Handle<reco::ConversionCollection>
&convCol,
const
math::XYZPoint
&beamspot,
float
dRMax = 0.1,
float
dEtaMax = 999.,
float
dPhiMax = 999.,
float
lxyMin=2.0,
float
probMin=1
e
-6,
unsigned
int
nHitsBeforeVtxMax=1);
70
71
static
bool
hasMatchedPromptElectron
(
const
reco::SuperClusterRef
&sc,
const
edm::Handle<reco::GsfElectronCollection>
&eleCol,
72
const
edm::Handle<reco::ConversionCollection>
&convCol,
const
math::XYZPoint
&beamspot,
float
lxyMin=2.0,
float
probMin=1
e
-6,
unsigned
int
nHitsBeforeVtxMax=0);
73
74
75
static
reco::GsfElectronRef
matchedPromptElectron
(
const
reco::SuperClusterRef
&sc,
const
edm::Handle<reco::GsfElectronCollection>
&eleCol,
76
const
edm::Handle<reco::ConversionCollection>
&convCol,
const
math::XYZPoint
&beamspot,
float
lxyMin=2.0,
float
probMin=1
e
-6,
unsigned
int
nHitsBeforeVtxMax=0);
77
78
};
79
#endif
reco::Conversion
Definition:
Conversion.h:25
conv
static HepMC::IO_HEPEVT conv
Definition:
BeamHaloProducer.cc:51
reco::GsfElectron
Definition:
GsfElectron.h:37
Photon.h
ConversionTools::hasMatchedPromptElectron
static bool hasMatchedPromptElectron(const reco::SuperClusterRef &sc, const edm::Handle< reco::GsfElectronCollection > &eleCol, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
Definition:
ConversionTools.cc:292
TrackFwd.h
edm::RefToBase< reco::Track >
edm::Handle< reco::ConversionCollection >
Point3D.h
ConversionTools::ConversionTools
ConversionTools()
Definition:
ConversionTools.h:40
MagneticField.h
ConversionTools
Definition:
ConversionTools.h:37
reco::SuperCluster
Definition:
SuperCluster.h:20
ConversionTools::hasMatchedConversion
static bool hasMatchedConversion(const reco::GsfElectron &ele, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
Definition:
ConversionTools.cc:145
ConversionTools::isGoodConversion
static bool isGoodConversion(const reco::Conversion &conv, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1)
Definition:
ConversionTools.cc:19
GsfElectron.h
ConversionTools::matchedPromptElectron
static reco::GsfElectronRef matchedPromptElectron(const reco::SuperClusterRef &sc, const edm::Handle< reco::GsfElectronCollection > &eleCol, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
Definition:
ConversionTools.cc:322
GsfTrack.h
GsfElectronFwd.h
TrackExtra.h
ConversionTools::matchesConversion
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true)
Definition:
ConversionTools.cc:51
math::XYZPoint
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition:
Point3D.h:13
alignCSCRings.e
list e
Definition:
alignCSCRings.py:90
ConversionFwd.h
Track.h
GsfTrackFwd.h
ConversionTools::matchedConversion
static reco::ConversionRef matchedConversion(const reco::GsfElectron &ele, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
Definition:
ConversionTools.cc:207
edm::Ref
Definition:
AssociativeIterator.h:52
Generated for CMSSW Reference Manual by
1.8.5