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 // $Id: PATConversionProducer.cc,v 1.1 2012/04/14 02:12:39 tjkim Exp $
3 //
5 
8 
13 
16 
19 
22 
25 
30 
34 
38 
39 #include <vector>
40 #include <memory>
41 #include "TMath.h"
42 
43 using namespace pat;
44 using namespace std;
45 
46 
48 
49  // general configurables
50  electronSrc_ = iConfig.getParameter<edm::InputTag>( "electronSource" );
51 
52  // produces vector of muons
53  produces<std::vector<Conversion> >();
54 
55 }
56 
57 
59 }
60 
61 
63 
64  // Get the collection of electrons from the event
66  iEvent.getByLabel(electronSrc_, electrons);
67 
69  iEvent.getByLabel("offlineBeamSpot", bsHandle);
70  const reco::BeamSpot &beamspot = *bsHandle.product();
71 
72  // for conversion veto selection
74  iEvent.getByLabel("allConversions", hConversions);
75 
76  std::vector<Conversion> * patConversions = new std::vector<Conversion>();
77 
78  for (reco::ConversionCollection::const_iterator conv = hConversions->begin(); conv!= hConversions->end(); ++conv) {
79 
80  reco::Vertex vtx = conv->conversionVertex();
81 
82  int index = 0;
83  for (edm::View<reco::GsfElectron>::const_iterator itElectron = electrons->begin(); itElectron != electrons->end(); ++itElectron) {
84 
85  //find matched conversions with electron and save those conversions with matched electron index
86  if (ConversionTools::matchesConversion(*itElectron, *conv)) {
87 
88  double vtxProb = TMath::Prob( vtx.chi2(), vtx.ndof());
89  math::XYZVector mom(conv->refittedPairMomentum());
90  double dbsx = vtx.x() - beamspot.position().x();
91  double dbsy = vtx.y() - beamspot.position().y();
92  double lxy = (mom.x()*dbsx + mom.y()*dbsy)/mom.rho();
93  int nHitsMax = 0;
94 
95  for (std::vector<uint8_t>::const_iterator it = conv->nHitsBeforeVtx().begin(); it!=conv->nHitsBeforeVtx().end(); ++it) {
96  if ((*it) > nHitsMax) nHitsMax = (*it);
97  }
98 
99  pat::Conversion anConversion( index );
100  anConversion.setVtxProb( vtxProb );
101  anConversion.setLxy( lxy );
102  anConversion.setNHitsMax( nHitsMax );
103 
104  patConversions->push_back(anConversion);
105  break;
106  }
107  index++;
108  }
109 
110  }
111 
112  // add the electrons to the event output
113  std::auto_ptr<std::vector<Conversion> > ptr(patConversions);
114  iEvent.put(ptr);
115 
116 }
117 
119 
T getParameter(std::string const &) const
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
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
double y() const
y coordinate
Definition: Vertex.h:97
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
double chi2() const
chi-squares
Definition: Vertex.h:82
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
double ndof() const
Definition: Vertex.h:89
double x() const
x coordinate
Definition: Vertex.h:95
PATConversionProducer(const edm::ParameterSet &iConfig)
void setNHitsMax(int nHitsMax)
Definition: Conversion.cc:22
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
T const * product() const
Definition: Handle.h:74
const Point & position() const
position
Definition: BeamSpot.h:63
void setVtxProb(double vtxProb)
Definition: Conversion.cc:14