CMS 3D CMS Logo

ElectronConversionRejectionVars.cc
Go to the documentation of this file.
1 //
2 
16 
19 
30 
32 public:
33  explicit ElectronConversionRejectionVars(const edm::ParameterSet& iConfig);
35 
36  void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
37 
38 private:
42 };
43 
45  : probesToken_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("probes"))),
46  tracksToken_(consumes<reco::TrackCollection>(edm::InputTag("generalTracks"))),
48  produces<edm::ValueMap<float>>("dist");
49  produces<edm::ValueMap<float>>("dcot");
50  produces<edm::ValueMap<float>>("convradius");
51  produces<edm::ValueMap<float>>("passConvRej");
52 }
53 
55 
57  using namespace edm;
58 
59  // read input
63 
64  iEvent.getByToken(probesToken_, probes);
65  iEvent.getByToken(tracksToken_, tracks_h);
66  iEvent.getByToken(gsfElecsToken_, elHandle);
67 
68  float evt_bField = 3.8;
69 
70  // prepare vector for output
71  std::vector<float> values;
72  std::vector<float> values2;
73  std::vector<float> values3;
74  std::vector<float> values4;
75 
76  // fill: use brute force
77  double dist = 0.0;
78  double dcot = 0.0;
79  double convradius = 0.0;
80  double passConvRej = 0.0;
81 
82  View<reco::Candidate>::const_iterator probe, endprobes = probes->end();
84  reco::GsfElectronCollection::const_iterator eleIt = electronCollection->begin();
85 
86  for (probe = probes->begin(); probe != endprobes; ++probe) {
87  for (eleIt = electronCollection->begin(); eleIt != electronCollection->end(); eleIt++) {
88  if (fabs(eleIt->et() - probe->et()) < 0.05 && fabs(eleIt->eta() - probe->eta()) < 0.01 &&
89  fabs(eleIt->phi() - probe->phi()) < 0.01) {
90  //we have a match
91  ConversionInfo convInfo = egammaTools::getConversionInfo(*eleIt, tracks_h, evt_bField);
92  dist = convInfo.dist;
93  dcot = convInfo.dcot;
94  convradius = convInfo.radiusOfConversion;
95  if (fabs(dist) > 0.02 && fabs(dcot) > 0.02)
96  passConvRej = 1.0;
97  break; //got our guy, so break
98  }
99  }
100  values.push_back(dist);
101  values2.push_back(dcot);
102  values3.push_back(convradius);
103  values4.push_back(passConvRej);
104  }
105 
106  // convert into ValueMap and store
107  auto valMap = std::make_unique<ValueMap<float>>();
109  filler.insert(probes, values.begin(), values.end());
110  filler.fill();
111  iEvent.put(std::move(valMap), "dist");
112 
113  // ---> same for dcot
114  auto valMap2 = std::make_unique<ValueMap<float>>();
115  ValueMap<float>::Filler filler2(*valMap2);
116  filler2.insert(probes, values2.begin(), values2.end());
117  filler2.fill();
118  iEvent.put(std::move(valMap2), "dcot");
119 
120  // ---> same for convradius
121  auto valMap3 = std::make_unique<ValueMap<float>>();
122  ValueMap<float>::Filler filler3(*valMap3);
123  filler3.insert(probes, values3.begin(), values3.end());
124  filler3.fill();
125  iEvent.put(std::move(valMap3), "convradius");
126 
127  // ---> same for passConvRej
128  auto valMap4 = std::make_unique<ValueMap<float>>();
129  ValueMap<float>::Filler filler4(*valMap4);
130  filler4.insert(probes, values4.begin(), values4.end());
131  filler4.fill();
132  iEvent.put(std::move(valMap4), "passConvRej");
133 }
134 
const double dist
Definition: ConversionInfo.h:9
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
edm::EDGetTokenT< edm::View< reco::Candidate > > probesToken_
ConversionInfo getConversionInfo(const reco::GsfElectronCore &, const edm::Handle< reco::TrackCollection > &ctftracks_h, const edm::Handle< reco::GsfTrackCollection > &gsftracks_h, const double bFieldAtOrigin, const double minFracSharedHits=0.45)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
edm::EDGetTokenT< reco::TrackCollection > tracksToken_
edm::EDGetTokenT< reco::GsfElectronCollection > gsfElecsToken_
const double radiusOfConversion
Store electron partner track conversion-rejection quantities ("dist" and "dcot") in the TP tree...
const double dcot
ElectronConversionRejectionVars(const edm::ParameterSet &iConfig)
T const * product() const
Definition: Handle.h:69
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
def move(src, dest)
Definition: eostools.py:511