CMS 3D CMS Logo

ElectronSeedProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: ElectronProducers
4 // Class: ElectronSeedProducer
5 //
13 //
14 // Original Author: Ursula Berthon, Claude Charlot
15 // Created: Mon Mar 27 13:22:06 CEST 2006
16 //
17 //
18 
32 
34 public:
35  explicit ElectronSeedProducer(const edm::ParameterSet&);
36 
37  void produce(edm::Event&, const edm::EventSetup&) final;
38 
39  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
40 
41 private:
44 
46  std::vector<edm::EDGetTokenT<TrajectorySeedCollection>> initialSeeds_;
48 
49  std::unique_ptr<ElectronSeedGenerator> matcher_;
50 
51  bool applyHOverECut_ = true;
52  std::unique_ptr<ElectronHcalHelper> hcalHelper_ = nullptr;
55  double SCEtCut_;
56 
58  std::unique_ptr<hgcal::ClusterTools> hgcClusterTools_;
59 };
60 
61 using namespace reco;
62 
64  :
65 
66  initialSeeds_{
67  edm::vector_transform(conf.getParameter<std::vector<edm::InputTag>>("initialSeedsVector"),
68  [this](edm::InputTag const& tag) { return consumes<TrajectorySeedCollection>(tag); })}
69 
70 {
71  SCEtCut_ = conf.getParameter<double>("SCEtCut");
72  auto theconsumes = consumesCollector();
73 
74  // new beamSpot tag
75  beamSpotTag_ = consumes(conf.getParameter<edm::InputTag>("beamSpot"));
76 
77  // for H/E
78  applyHOverECut_ = conf.getParameter<bool>("applyHOverECut");
79  if (applyHOverECut_) {
81  hcalCfg.hOverEConeSize = conf.getParameter<double>("hOverEConeSize");
82  if (hcalCfg.hOverEConeSize > 0) {
83  hcalCfg.onlyBehindCluster = false;
84  hcalCfg.checkHcalStatus = false;
85 
86  hcalCfg.hbheRecHits = consumes<HBHERecHitCollection>(conf.getParameter<edm::InputTag>("hbheRecHits"));
87 
88  hcalCfg.eThresHB = conf.getParameter<EgammaHcalIsolation::arrayHB>("recHitEThresholdHB");
89  hcalCfg.maxSeverityHB = conf.getParameter<int>("maxHcalRecHitSeverity");
90  hcalCfg.eThresHE = conf.getParameter<EgammaHcalIsolation::arrayHE>("recHitEThresholdHE");
91  hcalCfg.maxSeverityHE = hcalCfg.maxSeverityHB;
92  }
93  hcalHelper_ = std::make_unique<ElectronHcalHelper>(hcalCfg, consumesCollector());
94 
95  allowHGCal_ = conf.getParameter<bool>("allowHGCal");
96  if (allowHGCal_) {
97  const edm::ParameterSet& hgcCfg = conf.getParameterSet("HGCalConfig");
98  hgcClusterTools_ = std::make_unique<hgcal::ClusterTools>(hgcCfg, theconsumes);
99  }
100 
101  maxHOverEBarrel_ = conf.getParameter<double>("maxHOverEBarrel");
102  maxHOverEEndcaps_ = conf.getParameter<double>("maxHOverEEndcaps");
103  }
104 
106  esg_tokens.token_bs = beamSpotTag_;
107  esg_tokens.token_vtx = mayConsume<reco::VertexCollection>(conf.getParameter<edm::InputTag>("vertices"));
108 
109  matcher_ = std::make_unique<ElectronSeedGenerator>(conf, esg_tokens, consumesCollector());
110 
111  superClusters_[0] = consumes(conf.getParameter<edm::InputTag>("barrelSuperClusters"));
112  superClusters_[1] = consumes(conf.getParameter<edm::InputTag>("endcapSuperClusters"));
113 
114  //register your products
115  produces<ElectronSeedCollection>();
116 }
117 
119  LogDebug("ElectronSeedProducer") << "[ElectronSeedProducer::produce] entering ";
120 
121  std::vector<TrajectorySeedCollection const*> initialSeedCollections;
122  std::unique_ptr<TrajectorySeedCollection> initialSeedCollectionPtr = nullptr; //created on the fly
123 
124  if (hcalHelper_) {
125  hcalHelper_->beginEvent(e, iSetup);
126  if (allowHGCal_) {
127  hgcClusterTools_->getEventSetup(iSetup);
128  hgcClusterTools_->getEvent(e);
129  }
130  }
131 
132  matcher_->setupES(iSetup);
133 
134  // get initial TrajectorySeeds
135  initialSeedCollections.clear();
136  for (auto const& seeds : initialSeeds_) {
137  initialSeedCollections.push_back(&e.get(seeds));
138  }
139 
140  auto seeds = std::make_unique<ElectronSeedCollection>();
141  auto const& beamSportPosition = e.get(beamSpotTag_).position();
142 
143  // loop over barrel + endcap
144  for (unsigned int i = 0; i < 2; i++) {
145  auto clusterRefs = filterClusters(beamSportPosition, e.getHandle(superClusters_[i]));
146  matcher_->run(e, clusterRefs, initialSeedCollections, *seeds);
147  }
148 
149  // store the accumulated result
150 #ifdef EDM_ML_DEBUG
151  for (auto const& seed : *seeds) {
152  SuperClusterRef superCluster = seed.caloCluster().castTo<SuperClusterRef>();
153  LogDebug("ElectronSeedProducer") << "new seed with " << seed.nHits() << " hits"
154  << ", charge " << seed.getCharge() << " and cluster energy "
155  << superCluster->energy() << " PID " << superCluster.id();
156  }
157 #endif
158  e.put(std::move(seeds));
159 }
160 
161 //===============================
162 // Filter the superclusters
163 // - with EtCut
164 // - with HoE using calo cone
165 //===============================
166 
168  math::XYZPoint const& beamSpotPosition, const edm::Handle<reco::SuperClusterCollection>& superClusters) const {
169  SuperClusterRefVector sclRefs;
170 
171  for (unsigned int i = 0; i < superClusters->size(); ++i) {
172  auto const& scl = (*superClusters)[i];
173  double sclEta = EleRelPoint(scl.position(), beamSpotPosition).eta();
174  if (scl.energy() / cosh(sclEta) > SCEtCut_) {
175  if (applyHOverECut_) {
176  bool hoeVeto = false;
177  double had = hcalHelper_->hcalESum(scl, 0);
178  double scle = scl.energy();
179  int det_group = scl.seed()->hitsAndFractions()[0].first.det();
180  int detector = scl.seed()->hitsAndFractions()[0].first.subdetId();
181  if (detector == EcalBarrel && had / scle < maxHOverEBarrel_)
182  hoeVeto = true;
183  else if (!allowHGCal_ && detector == EcalEndcap && had / scle < maxHOverEEndcaps_)
184  hoeVeto = true;
185  else if (allowHGCal_ && EcalTools::isHGCalDet((DetId::Detector)det_group)) {
186  float had_fraction = hgcClusterTools_->getClusterHadronFraction(*(scl.seed()));
187  hoeVeto = (had_fraction >= 0.f && had_fraction < maxHOverEEndcaps_);
188  }
189  if (hoeVeto) {
190  sclRefs.push_back({superClusters, i});
191  }
192  } else {
193  sclRefs.push_back({superClusters, i});
194  }
195  }
196  }
197  LogDebug("ElectronSeedProducer") << "Filtered out " << sclRefs.size() << " superclusters from "
198  << superClusters->size();
199 
200  return sclRefs;
201 }
202 
205 
206  // steering
207  desc.add<std::vector<edm::InputTag>>("initialSeedsVector", {});
208  desc.add<bool>("useRecoVertex", false);
209  desc.add<edm::InputTag>("vertices", {"offlinePrimaryVerticesWithBS"});
210  desc.add<edm::InputTag>("beamSpot", {"offlineBeamSpot"});
211  desc.add<bool>("dynamicPhiRoad", true);
212 
213  // SC filtering
214  desc.add<double>("SCEtCut", 0.0);
215 
216  // H/E
217  desc.add<bool>("applyHOverECut", true);
218  desc.add<double>("hOverEConeSize", 0.15);
219  desc.add<double>("maxHOverEBarrel", 0.15);
220  desc.add<double>("maxHOverEEndcaps", 0.15);
221  desc.add<edm::InputTag>("hbheRecHits", {"hbhereco"});
222  desc.add<std::vector<double>>("recHitEThresholdHB", {0., 0., 0., 0.});
223  desc.add<std::vector<double>>("recHitEThresholdHE", {0., 0., 0., 0., 0., 0., 0.});
224  desc.add<int>("maxHcalRecHitSeverity", 999999);
225 
226  // H/E equivalent for HGCal
227  desc.add<bool>("allowHGCal", false);
229  psd4.add<edm::InputTag>("HGCEEInput", {"HGCalRecHit", "HGCEERecHits"});
230  psd4.add<edm::InputTag>("HGCFHInput", {"HGCalRecHit", "HGCHEFRecHits"});
231  psd4.add<edm::InputTag>("HGCBHInput", {"HGCalRecHit", "HGCHEBRecHits"});
232  desc.add<edm::ParameterSetDescription>("HGCalConfig", psd4);
233 
234  // r/z windows
235  desc.add<double>("nSigmasDeltaZ1", 5.0); // in case beam spot is used for the matching
236  desc.add<double>("deltaZ1WithVertex", 25.0); // in case reco vertex is used for the matching
237  desc.add<double>("z2MaxB", 0.09);
238  desc.add<double>("r2MaxF", 0.15);
239  desc.add<double>("rMaxI", 0.2); // intermediate region SC in EB and 2nd hits in PXF
240 
241  // phi windows (dynamic)
242  desc.add<double>("LowPtThreshold", 5.0);
243  desc.add<double>("HighPtThreshold", 35.0);
244  desc.add<double>("SizeWindowENeg", 0.675);
245  desc.add<double>("DeltaPhi1Low", 0.23);
246  desc.add<double>("DeltaPhi1High", 0.08);
247  desc.add<double>("DeltaPhi2B", 0.008);
248  desc.add<double>("DeltaPhi2F", 0.012);
249 
250  // phi windows (non dynamic, overwritten in case dynamic is selected)
251  desc.add<double>("ePhiMin1", -0.125);
252  desc.add<double>("ePhiMax1", 0.075);
253  desc.add<double>("PhiMax2B", 0.002);
254  desc.add<double>("PhiMax2F", 0.003);
255 
256  desc.add<edm::InputTag>("barrelSuperClusters",
257  {"particleFlowSuperClusterECAL", "particleFlowSuperClusterECALBarrel"});
258  desc.add<edm::InputTag>("endcapSuperClusters",
259  {"particleFlowSuperClusterECAL", "particleFlowSuperClusterECALEndcapWithPreshower"});
260 
261  descriptions.add("ecalDrivenElectronSeedsDefault", desc);
262 }
263 
static bool isHGCalDet(DetId::Detector thedet)
identify HGCal cells
Definition: EcalTools.h:49
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
ProductID id() const
Accessor for product ID.
Definition: Ref.h:244
edm::EDGetTokenT< reco::BeamSpot > token_bs
ParameterSet const & getParameterSet(std::string const &) const
std::unique_ptr< ElectronHcalHelper > hcalHelper_
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
ElectronSeedProducer(const edm::ParameterSet &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::unique_ptr< ElectronSeedGenerator > matcher_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
Detector
Definition: DetId.h:24
reco::SuperClusterRefVector filterClusters(math::XYZPoint const &beamSpotPosition, const edm::Handle< reco::SuperClusterCollection > &superClusters) const
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< reco::SuperClusterCollection > superClusters_[2]
fixed size matrix
void produce(edm::Event &, const edm::EventSetup &) final
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
edm::EDGetTokenT< reco::BeamSpot > beamSpotTag_
std::array< double, 4 > arrayHB
edm::EDGetTokenT< std::vector< reco::Vertex > > token_vtx
def move(src, dest)
Definition: eostools.py:511
std::unique_ptr< hgcal::ClusterTools > hgcClusterTools_
std::vector< edm::EDGetTokenT< TrajectorySeedCollection > > initialSeeds_
#define LogDebug(id)
std::array< double, 7 > arrayHE