CMS 3D CMS Logo

HLTHFRecoEcalCandidateProducer.cc
Go to the documentation of this file.
1 //
2 // Package: EgammaHFProducers
3 // Class: HFRecoEcalCandidateProducers
4 //
7 //
8 // Original Author: Kevin Klapoetke University of Minnesota
9 // Created: Wed 26 Sept 2007
10 // $Id:
11 //
12 //
13 
25 
27 #include "HFValueStruct.h"
28 
29 #include <vector>
30 #include <memory>
31 
33 public:
35  void produce(edm::StreamID, edm::Event&, edm::EventSetup const&) const override;
36 
37 private:
39  const int HFDBversion_;
40  const std::vector<double> HFDBvector_;
41  const double Cut2D_;
42  const double defaultSlope2D_;
45 };
46 
48  : hfclusters_(conf.getParameter<edm::InputTag>("hfclusters")),
49  HFDBversion_(conf.existsAs<bool>("HFDBversion") ? conf.getParameter<int>("HFDBversion") : 99), //do nothing
50  HFDBvector_(conf.existsAs<bool>("HFDBvector") ? conf.getParameter<std::vector<double> >("HFDBvector")
51  : std::vector<double>{}),
52  Cut2D_(conf.getParameter<double>("intercept2DCut")),
53  defaultSlope2D_(
54  (Cut2D_ <= 0.83)
55  ? (0.475)
56  : ((Cut2D_ > 0.83 && Cut2D_ <= 0.9) ? (0.275) : (0.2))), //fix for hlt unable to add slope variable now
57  hfvars_(HFDBversion_, HFDBvector_),
58  algo_(conf.existsAs<bool>("Correct") ? conf.getParameter<bool>("Correct") : true,
59  conf.getParameter<double>("e9e25Cut"),
60  conf.getParameter<double>("intercept2DCut"),
61  conf.existsAs<bool>("intercept2DSlope") ? conf.getParameter<double>("intercept2DSlope") : defaultSlope2D_,
62  conf.getParameter<std::vector<double> >("e1e9Cut"),
63  conf.getParameter<std::vector<double> >("eCOREe9Cut"),
64  conf.getParameter<std::vector<double> >("eSeLCut"),
65  hfvars_) {
66  produces<reco::RecoEcalCandidateCollection>();
67 }
68 
72 
73  e.getByLabel(hfclusters_, super_clus);
74  e.getByLabel(hfclusters_, hf_assoc);
75 
76  int nvertex = 1;
77 
78  // create return data
79  auto retdata1 = std::make_unique<reco::RecoEcalCandidateCollection>();
80 
81  algo_.produce(super_clus, *hf_assoc, *retdata1, nvertex);
82 
83  e.put(std::move(retdata1));
84 }
HLTHFRecoEcalCandidateProducer::HFDBversion_
const int HFDBversion_
Definition: HLTHFRecoEcalCandidateProducer.cc:39
edm::StreamID
Definition: StreamID.h:30
Handle.h
electrons_cff.bool
bool
Definition: electrons_cff.py:366
MessageLogger.h
HFEMClusterShapeAssociation.h
HFRecoEcalCandidateAlgo::produce
void produce(const edm::Handle< reco::SuperClusterCollection > &SuperClusters, const reco::HFEMClusterShapeAssociationCollection &AssocShapes, reco::RecoEcalCandidateCollection &RecoECand, int nvtx) const
Definition: HFRecoEcalCandidateAlgo.cc:105
edm
HLT enums.
Definition: AlignableModifier.h:19
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89285
HLTHFRecoEcalCandidateProducer::produce
void produce(edm::StreamID, edm::Event &, edm::EventSetup const &) const override
Definition: HLTHFRecoEcalCandidateProducer.cc:69
HLTHFRecoEcalCandidateProducer::Cut2D_
const double Cut2D_
Definition: HLTHFRecoEcalCandidateProducer.cc:41
reco::HFValueStruct
Definition: HFValueStruct.h:8
edm::Handle
Definition: AssociativeIterator.h:50
HFRecoEcalCandidateAlgo
Definition: HFRecoEcalCandidateAlgo.h:23
HLTHFRecoEcalCandidateProducer
Definition: HLTHFRecoEcalCandidateProducer.cc:32
HLTHFRecoEcalCandidateProducer::vertices_
const edm::InputTag vertices_
Definition: HLTHFRecoEcalCandidateProducer.cc:38
edm::global::EDProducer
Definition: EDProducer.h:32
HFRecoEcalCandidateAlgo.h
HLTHFRecoEcalCandidateProducer::HLTHFRecoEcalCandidateProducer
HLTHFRecoEcalCandidateProducer(edm::ParameterSet const &conf)
Definition: HLTHFRecoEcalCandidateProducer.cc:47
edm::ParameterSet
Definition: ParameterSet.h:47
HFValueStruct.h
Event.h
createfilelist.int
int
Definition: createfilelist.py:10
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
HLTHFRecoEcalCandidateProducer::algo_
const HFRecoEcalCandidateAlgo algo_
Definition: HLTHFRecoEcalCandidateProducer.cc:44
edm::EventSetup
Definition: EventSetup.h:58
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
RecoEcalCandidate.h
HLTHFRecoEcalCandidateProducer::defaultSlope2D_
const double defaultSlope2D_
Definition: HLTHFRecoEcalCandidateProducer.cc:42
SuperCluster.h
HLTHFRecoEcalCandidateProducer::hfvars_
const reco::HFValueStruct hfvars_
Definition: HLTHFRecoEcalCandidateProducer.cc:43
HLTHFRecoEcalCandidateProducer::hfclusters_
const edm::InputTag hfclusters_
Definition: HLTHFRecoEcalCandidateProducer.cc:38
EventSetup.h
Exception.h
ParameterSet.h
EDProducer.h
edm::Event
Definition: Event.h:73
HFEMClusterShape.h
edm::InputTag
Definition: InputTag.h:15
HLTHFRecoEcalCandidateProducer::HFDBvector_
const std::vector< double > HFDBvector_
Definition: HLTHFRecoEcalCandidateProducer.cc:40
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37