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 
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true, bool allowAmbiguousGsfMatch=false)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
const edm::EDGetTokenT< reco::ConversionCollection > conversionsToken_
void setLxy(double lxy)
Definition: Conversion.cc:10
static HepMC::IO_HEPEVT conv
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
double y() const
y coordinate
Definition: Vertex.h:117
const edm::EDGetTokenT< edm::View< reco::GsfElectron > > electronToken_
const edm::EDGetTokenT< reco::BeamSpot > bsToken_
Definition: HeavyIon.h:7
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double chi2() const
chi-squares
Definition: Vertex.h:102
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
double ndof() const
Definition: Vertex.h:109
std::vector< Conversion > ConversionCollection
Definition: Conversion.h:13
double x() const
x coordinate
Definition: Vertex.h:115
PATConversionProducer(const edm::ParameterSet &iConfig)
void setNHitsMax(int nHitsMax)
Definition: Conversion.cc:12
T const * product() const
Definition: Handle.h:69
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
const Point & position() const
position
Definition: BeamSpot.h:59
void setVtxProb(double vtxProb)
Definition: Conversion.cc:8
def move(src, dest)
Definition: eostools.py:511