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 
45 
47  electronToken_(consumes<edm::View<reco::GsfElectron> >(iConfig.getParameter<edm::InputTag>( "electronSource" ))),
48  bsToken_(consumes<reco::BeamSpot>(edm::InputTag("offlineBeamSpot"))),
49  conversionsToken_(consumes<reco::ConversionCollection>(edm::InputTag("allConversions"))) {
50  // produces vector of muons
51  produces<std::vector<Conversion> >();
52 }
53 
54 
56 }
57 
58 
60 
61  // Get the collection of electrons from the event
63  iEvent.getByToken(electronToken_, electrons);
64 
66  iEvent.getByToken(bsToken_, bsHandle);
67  const reco::BeamSpot &beamspot = *bsHandle.product();
68 
69  // for conversion veto selection
71  iEvent.getByToken(conversionsToken_, hConversions);
72 
73  std::vector<Conversion> * patConversions = new std::vector<Conversion>();
74 
75  for (reco::ConversionCollection::const_iterator conv = hConversions->begin(); conv!= hConversions->end(); ++conv) {
76 
77  reco::Vertex vtx = conv->conversionVertex();
78 
79  int index = 0;
80  for (edm::View<reco::GsfElectron>::const_iterator itElectron = electrons->begin(); itElectron != electrons->end(); ++itElectron) {
81 
82  //find matched conversions with electron and save those conversions with matched electron index
83  if (ConversionTools::matchesConversion(*itElectron, *conv)) {
84 
85  double vtxProb = TMath::Prob( vtx.chi2(), vtx.ndof());
86  math::XYZVector mom(conv->refittedPairMomentum());
87  double dbsx = vtx.x() - beamspot.position().x();
88  double dbsy = vtx.y() - beamspot.position().y();
89  double lxy = (mom.x()*dbsx + mom.y()*dbsy)/mom.rho();
90  int nHitsMax = 0;
91 
92  for (std::vector<uint8_t>::const_iterator it = conv->nHitsBeforeVtx().begin(); it!=conv->nHitsBeforeVtx().end(); ++it) {
93  if ((*it) > nHitsMax) nHitsMax = (*it);
94  }
95 
96  pat::Conversion anConversion( index );
97  anConversion.setVtxProb( vtxProb );
98  anConversion.setLxy( lxy );
99  anConversion.setNHitsMax( nHitsMax );
100 
101  patConversions->push_back(anConversion);
102  break;
103  }
104  index++;
105  }
106 
107  }
108 
109  // add the electrons to the event output
110  std::unique_ptr<std::vector<Conversion> > ptr(patConversions);
111  iEvent.put(std::move(ptr));
112 
113 }
114 
116 
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:137
const edm::EDGetTokenT< reco::ConversionCollection > conversionsToken_
void setLxy(double lxy)
Definition: Conversion.cc:18
static HepMC::IO_HEPEVT conv
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
double y() const
y coordinate
Definition: Vertex.h:113
const edm::EDGetTokenT< edm::View< reco::GsfElectron > > electronToken_
const edm::EDGetTokenT< reco::BeamSpot > bsToken_
Definition: HeavyIon.h:7
int iEvent
Definition: GenABIO.cc:230
double chi2() const
chi-squares
Definition: Vertex.h:98
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
double ndof() const
Definition: Vertex.h:105
std::vector< Conversion > ConversionCollection
Definition: Conversion.h:13
double x() const
x coordinate
Definition: Vertex.h:111
PATConversionProducer(const edm::ParameterSet &iConfig)
void setNHitsMax(int nHitsMax)
Definition: Conversion.cc:22
T const * product() const
Definition: Handle.h:81
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
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:62
void setVtxProb(double vtxProb)
Definition: Conversion.cc:14
def move(src, dest)
Definition: eostools.py:510