CMS 3D CMS Logo

PATConversionProducer.cc
Go to the documentation of this file.
1 //
2 //
4 
7 
12 
15 
18 
21 
24 
29 
33 
37 
38 #include <vector>
39 #include <memory>
40 #include "TMath.h"
41 
42 using namespace pat;
43 using namespace std;
44 
46  : electronToken_(consumes<edm::View<reco::GsfElectron> >(iConfig.getParameter<edm::InputTag>("electronSource"))),
47  bsToken_(consumes<reco::BeamSpot>(edm::InputTag("offlineBeamSpot"))),
48  conversionsToken_(consumes<reco::ConversionCollection>(edm::InputTag("allConversions"))) {
49  // produces vector of muons
50  produces<std::vector<Conversion> >();
51 }
52 
54 
56  // Get the collection of electrons from the event
58  iEvent.getByToken(electronToken_, electrons);
59 
61  iEvent.getByToken(bsToken_, bsHandle);
62  const reco::BeamSpot &beamspot = *bsHandle.product();
63 
64  // for conversion veto selection
66  iEvent.getByToken(conversionsToken_, hConversions);
67 
68  std::vector<Conversion> *patConversions = new std::vector<Conversion>();
69 
70  for (reco::ConversionCollection::const_iterator conv = hConversions->begin(); conv != hConversions->end(); ++conv) {
71  reco::Vertex vtx = conv->conversionVertex();
72 
73  int index = 0;
74  for (edm::View<reco::GsfElectron>::const_iterator itElectron = electrons->begin(); itElectron != electrons->end();
75  ++itElectron) {
76  //find matched conversions with electron and save those conversions with matched electron index
77  if (ConversionTools::matchesConversion(*itElectron, *conv)) {
78  double vtxProb = TMath::Prob(vtx.chi2(), vtx.ndof());
79  math::XYZVector mom(conv->refittedPairMomentum());
80  double dbsx = vtx.x() - beamspot.position().x();
81  double dbsy = vtx.y() - beamspot.position().y();
82  double lxy = (mom.x() * dbsx + mom.y() * dbsy) / mom.rho();
83  int nHitsMax = 0;
84 
85  for (std::vector<uint8_t>::const_iterator it = conv->nHitsBeforeVtx().begin();
86  it != conv->nHitsBeforeVtx().end();
87  ++it) {
88  if ((*it) > nHitsMax)
89  nHitsMax = (*it);
90  }
91 
92  pat::Conversion anConversion(index);
93  anConversion.setVtxProb(vtxProb);
94  anConversion.setLxy(lxy);
95  anConversion.setNHitsMax(nHitsMax);
96 
97  patConversions->push_back(anConversion);
98  break;
99  }
100  index++;
101  }
102  }
103 
104  // add the electrons to the event output
105  std::unique_ptr<std::vector<Conversion> > ptr(patConversions);
106  iEvent.put(std::move(ptr));
107 }
108 
110 
ConfigurationDescriptions.h
edm::StreamID
Definition: StreamID.h:30
MessageLogger.h
pat::PATConversionProducer::conversionsToken_
const edm::EDGetTokenT< reco::ConversionCollection > conversionsToken_
Definition: PATConversionProducer.h:46
edm::Handle::product
T const * product() const
Definition: Handle.h:70
align::BeamSpot
Definition: StructureType.h:89
sistrip::View
View
Definition: ConstantsForView.h:26
pat::PATConversionProducer
Definition: PATConversionProducer.h:35
PFCandidate.h
conv
static HepMC::IO_HEPEVT conv
Definition: BeamHaloProducer.cc:48
pat::ConversionCollection
std::vector< Conversion > ConversionCollection
Definition: Conversion.h:13
pat::PATConversionProducer::bsToken_
const edm::EDGetTokenT< reco::BeamSpot > bsToken_
Definition: PATConversionProducer.h:45
edm
HLT enums.
Definition: AlignableModifier.h:19
TrackerIsolationPt.h
pat::Conversion::setLxy
void setLxy(double lxy)
Definition: Conversion.cc:10
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89287
Conversion.h
PATConversionProducer.h
EcalClusterLazyTools.h
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
TransientTrack.h
Association.h
edm::Handle
Definition: AssociativeIterator.h:50
FileInPath.h
GenParticle.h
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
BeamSpot.h
GenParticleFwd.h
reco::BeamSpot
Definition: BeamSpot.h:21
beamspot
Definition: BeamSpotWrite2Txt.h:8
CaloIsolationEnergy.h
pat::Conversion
Definition: Conversion.h:21
pat::PATConversionProducer::produce
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
Definition: PATConversionProducer.cc:55
ParameterSetDescription.h
Vertex.h
pat::Conversion::setNHitsMax
void setNHitsMax(int nHitsMax)
Definition: Conversion.cc:12
TransientTrackBuilder.h
edm::ParameterSet
Definition: ParameterSet.h:47
ConversionTools::matchesConversion
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true, bool allowAmbiguousGsfMatch=false)
Definition: ConversionTools.cc:53
math::XYZVector
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
LorentzVector.h
iEvent
int iEvent
Definition: GenABIO.cc:224
GsfTrack.h
edm::EventSetup
Definition: EventSetup.h:57
pat
Definition: HeavyIon.h:7
TransientTrackRecord.h
pat::PATConversionProducer::PATConversionProducer
PATConversionProducer(const edm::ParameterSet &iConfig)
Definition: PATConversionProducer.cc:45
ValueMap.h
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
IPTools.h
GsfTrackFwd.h
extraflags_cff.vtx
vtx
Definition: extraflags_cff.py:18
pwdgSkimBPark_cfi.electrons
electrons
Definition: pwdgSkimBPark_cfi.py:6
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
edm::View::const_iterator
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
edm::Event
Definition: Event.h:73
pat::PATConversionProducer::~PATConversionProducer
~PATConversionProducer() override
Definition: PATConversionProducer.cc:53
pat::PATConversionProducer::electronToken_
const edm::EDGetTokenT< edm::View< reco::GsfElectron > > electronToken_
Definition: PATConversionProducer.h:44
pat::Conversion::setVtxProb
void setVtxProb(double vtxProb)
Definition: Conversion.cc:8
reco::Vertex
Definition: Vertex.h:35
PFCandidateFwd.h
ConversionTools.h