CMS 3D CMS Logo

ElectronConversionRejectionVars.cc
Go to the documentation of this file.
1 //
2 
17 
20 
31 
32 
34  public:
35  explicit ElectronConversionRejectionVars(const edm::ParameterSet & iConfig);
37 
38  void produce(edm::Event & iEvent, const edm::EventSetup & iSetup) override;
39 
40  private:
44 };
45 
47  probesToken_(consumes<edm::View<reco::Candidate> >(iConfig.getParameter<edm::InputTag>("probes"))),
48  tracksToken_(consumes<reco::TrackCollection>(edm::InputTag("generalTracks"))) ,
49  gsfElecsToken_(consumes<reco::GsfElectronCollection>(edm::InputTag("gsfElectrons")))
50 {
51  produces<edm::ValueMap<float> >("dist");
52  produces<edm::ValueMap<float> >("dcot");
53  produces<edm::ValueMap<float> >("convradius");
54  produces<edm::ValueMap<float> >("passConvRej");
55 }
56 
57 
59 {
60 }
61 
62 void
64  using namespace edm;
65 
66  // read input
70 
71  iEvent.getByToken(probesToken_, probes);
72  iEvent.getByToken(tracksToken_, tracks_h );
73  iEvent.getByToken(gsfElecsToken_, elHandle);
74 
75  float evt_bField = 3.8;
76 
77 
78  // prepare vector for output
79  std::vector<float> values;
80  std::vector<float> values2;
81  std::vector<float> values3;
82  std::vector<float> values4;
83 
84  // fill: use brute force
85  double dist = 0.0;
86  double dcot = 0.0;
87  double convradius = 0.0;
88  double passConvRej = 0.0;
89 
90  View<reco::Candidate>::const_iterator probe, endprobes = probes->end();
92  reco::GsfElectronCollection::const_iterator eleIt = electronCollection->begin();
93 
94  for (probe = probes->begin(); probe != endprobes; ++probe) {
95  for (eleIt=electronCollection->begin(); eleIt!=electronCollection->end(); eleIt++) {
96  if( fabs(eleIt->et() - probe->et() ) < 0.05 && fabs(eleIt->eta() - probe->eta() ) < 0.01
97  && fabs(eleIt->phi() - probe->phi() ) < 0.01 ){
98  //we have a match
99  ConversionInfo convInfo = egammaTools::getConversionInfo(*eleIt, tracks_h, evt_bField);
100  dist = convInfo.dist;
101  dcot = convInfo.dcot;
102  convradius = convInfo.radiusOfConversion;
103  if( fabs(dist)>0.02 && fabs(dcot)>0.02) passConvRej = 1.0;
104  break; //got our guy, so break
105  }
106  }
107  values.push_back(dist);
108  values2.push_back(dcot);
109  values3.push_back(convradius);
110  values4.push_back(passConvRej);
111  }
112 
113 
114  // convert into ValueMap and store
115  auto valMap = std::make_unique<ValueMap<float>>();
117  filler.insert(probes, values.begin(), values.end());
118  filler.fill();
119  iEvent.put(std::move(valMap), "dist");
120 
121 
122  // ---> same for dcot
123  auto valMap2 = std::make_unique<ValueMap<float>>();
124  ValueMap<float>::Filler filler2(*valMap2);
125  filler2.insert(probes, values2.begin(), values2.end());
126  filler2.fill();
127  iEvent.put(std::move(valMap2), "dcot");
128 
129  // ---> same for convradius
130  auto valMap3 = std::make_unique<ValueMap<float>>();
131  ValueMap<float>::Filler filler3(*valMap3);
132  filler3.insert(probes, values3.begin(), values3.end());
133  filler3.fill();
134  iEvent.put(std::move(valMap3), "convradius");
135 
136 
137  // ---> same for passConvRej
138  auto valMap4 = std::make_unique<ValueMap<float>>();
139  ValueMap<float>::Filler filler4(*valMap4);
140  filler4.insert(probes, values4.begin(), values4.end());
141  filler4.fill();
142  iEvent.put(std::move(valMap4), "passConvRej");
143 }
144 
145 
const double dist
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
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:517
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:15
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:74
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