CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
48  // general configurables
49  electronToken_ = consumes<edm::View<reco::GsfElectron> >(iConfig.getParameter<edm::InputTag>( "electronSource" ));
50  bsToken_ = consumes<reco::BeamSpot>(edm::InputTag("offlineBeamSpot"));
51  conversionsToken_ = consumes<reco::ConversionCollection>(edm::InputTag("allConversions"));
52 
53  // produces vector of muons
54  produces<std::vector<Conversion> >();
55 
56 }
57 
58 
60 }
61 
62 
64 
65  // Get the collection of electrons from the event
67  iEvent.getByToken(electronToken_, electrons);
68 
70  iEvent.getByToken(bsToken_, bsHandle);
71  const reco::BeamSpot &beamspot = *bsHandle.product();
72 
73  // for conversion veto selection
75  iEvent.getByToken(conversionsToken_, hConversions);
76 
77  std::vector<Conversion> * patConversions = new std::vector<Conversion>();
78 
79  for (reco::ConversionCollection::const_iterator conv = hConversions->begin(); conv!= hConversions->end(); ++conv) {
80 
81  reco::Vertex vtx = conv->conversionVertex();
82 
83  int index = 0;
84  for (edm::View<reco::GsfElectron>::const_iterator itElectron = electrons->begin(); itElectron != electrons->end(); ++itElectron) {
85 
86  //find matched conversions with electron and save those conversions with matched electron index
87  if (ConversionTools::matchesConversion(*itElectron, *conv)) {
88 
89  double vtxProb = TMath::Prob( vtx.chi2(), vtx.ndof());
90  math::XYZVector mom(conv->refittedPairMomentum());
91  double dbsx = vtx.x() - beamspot.position().x();
92  double dbsy = vtx.y() - beamspot.position().y();
93  double lxy = (mom.x()*dbsx + mom.y()*dbsy)/mom.rho();
94  int nHitsMax = 0;
95 
96  for (std::vector<uint8_t>::const_iterator it = conv->nHitsBeforeVtx().begin(); it!=conv->nHitsBeforeVtx().end(); ++it) {
97  if ((*it) > nHitsMax) nHitsMax = (*it);
98  }
99 
100  pat::Conversion anConversion( index );
101  anConversion.setVtxProb( vtxProb );
102  anConversion.setLxy( lxy );
103  anConversion.setNHitsMax( nHitsMax );
104 
105  patConversions->push_back(anConversion);
106  break;
107  }
108  index++;
109  }
110 
111  }
112 
113  // add the electrons to the event output
114  std::auto_ptr<std::vector<Conversion> > ptr(patConversions);
115  iEvent.put(ptr);
116 
117 }
118 
120 
T getParameter(std::string const &) const
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true, bool allowAmbiguousGsfMatch=false)
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
void setLxy(double lxy)
Definition: Conversion.cc:18
static HepMC::IO_HEPEVT conv
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
double y() const
y coordinate
Definition: Vertex.h:96
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
double chi2() const
chi-squares
Definition: Vertex.h:81
double ndof() const
Definition: Vertex.h:88
double x() const
x coordinate
Definition: Vertex.h:94
PATConversionProducer(const edm::ParameterSet &iConfig)
void setNHitsMax(int nHitsMax)
Definition: Conversion.cc:22
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
T const * product() const
Definition: Handle.h:81
const Point & position() const
position
Definition: BeamSpot.h:62
void setVtxProb(double vtxProb)
Definition: Conversion.cc:14