CMS 3D CMS Logo

LowPtElectronsModifier.cc
Go to the documentation of this file.
7 
9 //
11 public:
13  ~LowPtElectronModifier() override = default;
14 
15  void setEvent(const edm::Event&) final;
16  void setEventContent(const edm::EventSetup&) final;
17 
18  void modifyObject(pat::Electron& ele) const final;
19 
20 private:
24  reco::BeamSpot const* beamSpot_ = nullptr;
27  bool extra_;
28 };
29 
31 //
33  : ModifyObjectValueBase(conf),
34  convT_(cc.consumes<reco::ConversionCollection>(conf.getParameter<edm::InputTag>("conversions"))),
35  conv_(),
36  beamSpotT_(cc.consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("beamSpot"))),
37  beamSpot_(),
38  verticesT_(cc.consumes<reco::VertexCollection>(conf.getParameter<edm::InputTag>("vertices"))),
39  vertices_(),
40  extra_(conf.getParameter<bool>("addExtraUserVars")) {
41  ;
42 }
43 
45 //
47  conv_ = &iEvent.get(convT_);
48  beamSpot_ = &iEvent.get(beamSpotT_);
49  vertices_ = &iEvent.get(verticesT_);
50 }
51 
53 //
55 
57 //
59  // Embed Conversion info
61  conv.match(*beamSpot_, *conv_, ele);
62  conv.addUserVars(ele);
63  if (extra_) {
64  conv.addExtraUserVars(ele);
65  }
66  // Set impact parameters
67  auto const& gsfTrack = *ele.gsfTrack();
68  if (!vertices_->empty()) {
69  const reco::Vertex& pv = vertices_->front();
70  ele.setDB(gsfTrack.dxy(pv.position()),
71  gsfTrack.dxyError(pv.position(), pv.covariance()),
72  pat::Electron::PV2D); // PV2D
73  ele.setDB(gsfTrack.dz(pv.position()), std::hypot(gsfTrack.dzError(), pv.zError()),
74  pat::Electron::PVDZ); // PVDZ
75  }
76  ele.setDB(gsfTrack.dxy(*beamSpot_), gsfTrack.dxyError(*beamSpot_),
77  pat::Electron::BS2D); // BS2D
78 }
79 
81 //
MPlex< T, D1, D2, N > hypot(const MPlex< T, D1, D2, N > &a, const MPlex< T, D1, D2, N > &b)
Definition: Matriplex.h:417
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
LowPtElectronModifier(const edm::ParameterSet &conf, edm::ConsumesCollector &)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
reco::GsfTrackRef gsfTrack() const override
override the reco::GsfElectron::gsfTrack method, to access the internal storage of the supercluster ...
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
void setDB(double dB, double edB, IPTYPE type)
Set impact parameter of a certain type and its uncertainty.
int iEvent
Definition: GenABIO.cc:224
reco::BeamSpot const * beamSpot_
const edm::EDGetTokenT< reco::BeamSpot > beamSpotT_
const edm::EDGetTokenT< reco::VertexCollection > verticesT_
const edm::EDGetTokenT< reco::ConversionCollection > convT_
~LowPtElectronModifier() override=default
Analysis-level electron class.
Definition: Electron.h:51
EPOS::IO_EPOS conv
fixed size matrix
HLT enums.
reco::ConversionCollection const * conv_
#define DEFINE_EDM_PLUGIN(factory, type, name)
void setEvent(const edm::Event &) final
void modifyObject(pat::Electron &ele) const final
void setEventContent(const edm::EventSetup &) final
reco::VertexCollection const * vertices_