CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PATConversionProducer.cc
Go to the documentation of this file.
38 
39 #include <memory>
40 #include <string>
41 #include <vector>
42 
43 namespace pat {
44 
46  public:
48  ~PATConversionProducer() override;
49 
50  void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override;
51 
52  private:
53  // configurables
57  };
58 
59 } // namespace pat
60 
61 using namespace pat;
62 using namespace std;
63 
65  : electronToken_(consumes<edm::View<reco::GsfElectron> >(iConfig.getParameter<edm::InputTag>("electronSource"))),
66  bsToken_(consumes<reco::BeamSpot>(edm::InputTag("offlineBeamSpot"))),
67  conversionsToken_(consumes<reco::ConversionCollection>(edm::InputTag("allConversions"))) {
68  // produces vector of muons
69  produces<std::vector<Conversion> >();
70 }
71 
73 
75  // Get the collection of electrons from the event
77  iEvent.getByToken(electronToken_, electrons);
78 
80  iEvent.getByToken(bsToken_, bsHandle);
81  const reco::BeamSpot &beamspot = *bsHandle.product();
82 
83  // for conversion veto selection
85  iEvent.getByToken(conversionsToken_, hConversions);
86 
87  std::vector<Conversion> *patConversions = new std::vector<Conversion>();
88 
89  for (reco::ConversionCollection::const_iterator conv = hConversions->begin(); conv != hConversions->end(); ++conv) {
90  reco::Vertex vtx = conv->conversionVertex();
91 
92  int index = 0;
93  for (edm::View<reco::GsfElectron>::const_iterator itElectron = electrons->begin(); itElectron != electrons->end();
94  ++itElectron) {
95  //find matched conversions with electron and save those conversions with matched electron index
96  if (ConversionTools::matchesConversion(*itElectron, *conv)) {
97  double vtxProb = TMath::Prob(vtx.chi2(), vtx.ndof());
98  math::XYZVector mom(conv->refittedPairMomentum());
99  double dbsx = vtx.x() - beamspot.position().x();
100  double dbsy = vtx.y() - beamspot.position().y();
101  double lxy = (mom.x() * dbsx + mom.y() * dbsy) / mom.rho();
102  int nHitsMax = 0;
103 
104  for (std::vector<uint8_t>::const_iterator it = conv->nHitsBeforeVtx().begin();
105  it != conv->nHitsBeforeVtx().end();
106  ++it) {
107  if ((*it) > nHitsMax)
108  nHitsMax = (*it);
109  }
110 
111  pat::Conversion anConversion(index);
112  anConversion.setVtxProb(vtxProb);
113  anConversion.setLxy(lxy);
114  anConversion.setNHitsMax(nHitsMax);
115 
116  patConversions->push_back(anConversion);
117  break;
118  }
119  index++;
120  }
121  }
122 
123  // add the electrons to the event output
124  std::unique_ptr<std::vector<Conversion> > ptr(patConversions);
125  iEvent.put(std::move(ptr));
126 }
127 
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:133
const edm::EDGetTokenT< reco::ConversionCollection > conversionsToken_
void setLxy(double lxy)
Definition: Conversion.cc:10
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double y() const
y coordinate
Definition: Vertex.h:131
const edm::EDGetTokenT< edm::View< reco::GsfElectron > > electronToken_
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
const edm::EDGetTokenT< reco::BeamSpot > bsToken_
int iEvent
Definition: GenABIO.cc:224
def move
Definition: eostools.py:511
double chi2() const
chi-squares
Definition: Vertex.h:116
double ndof() const
Definition: Vertex.h:123
std::vector< Conversion > ConversionCollection
Definition: Conversion.h:13
double x() const
x coordinate
Definition: Vertex.h:129
PATConversionProducer(const edm::ParameterSet &iConfig)
void setNHitsMax(int nHitsMax)
Definition: Conversion.cc:12
T const * product() const
Definition: Handle.h:70
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
EPOS::IO_EPOS conv
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