CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HFRecoEcalCandidateProducer.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:
34  explicit HFRecoEcalCandidateProducer(edm::ParameterSet const& conf);
35  void produce(edm::Event& e, edm::EventSetup const& iSetup) override;
36 
37 private:
38  std::vector<double> defaultDB_;
41  std::vector<double> HFDBvector_;
42  bool doPU_;
43  double Cut2D_;
47 };
48 
51 
53  : defaultDB_(std::vector<double>()),
54  hfclustersSC_(consumes<reco::SuperClusterCollection>(conf.getParameter<edm::InputTag>("hfclusters"))),
55  hfclustersHFEM_(
56  consumes<reco::HFEMClusterShapeAssociationCollection>(conf.getParameter<edm::InputTag>("hfclusters"))),
57 
58  HFDBversion_(conf.existsAs<int>("HFDBversion") ? conf.getParameter<int>("HFDBversion") : 99), //do nothing
59  HFDBvector_(conf.existsAs<std::vector<double> >("HFDBvector")
60  ? conf.getParameter<std::vector<double> >("HFDBvector")
61  : defaultDB_),
62  doPU_(false),
63  Cut2D_(conf.getParameter<double>("intercept2DCut")),
65  (Cut2D_ <= 0.83)
66  ? (0.475)
67  : ((Cut2D_ > 0.83 && Cut2D_ <= 0.9) ? (0.275) : (0.2))), //fix for hlt unable to add slope variable now
68  hfvars_(HFDBversion_, HFDBvector_),
69  algo_(conf.existsAs<bool>("Correct") ? conf.getParameter<bool>("Correct") : true,
70  conf.getParameter<double>("e9e25Cut"),
71  conf.getParameter<double>("intercept2DCut"),
72  conf.existsAs<double>("intercept2DSlope") ? conf.getParameter<double>("intercept2DSlope") : defaultSlope2D_,
73  conf.getParameter<std::vector<double> >("e1e9Cut"),
74  conf.getParameter<std::vector<double> >("eCOREe9Cut"),
75  conf.getParameter<std::vector<double> >("eSeLCut"),
76  hfvars_) {
77  if (conf.existsAs<edm::InputTag>("VertexCollection")) {
78  vertices_ = consumes<reco::VertexCollection>(conf.getParameter<edm::InputTag>("VertexCollection"));
79  } else
80  vertices_ = consumes<reco::VertexCollection>(edm::InputTag("offlinePrimaryVertices"));
81 
82  produces<reco::RecoEcalCandidateCollection>();
83 }
84 
88 
89  e.getByToken(hfclustersSC_, super_clus);
90  e.getByToken(hfclustersHFEM_, hf_assoc);
91 
92  int nvertex = 0;
93  if (HFDBversion_ != 99) {
95  e.getByToken(vertices_, pvHandle);
96  const reco::VertexCollection& vertices = *pvHandle.product();
97  static const int minNDOF = 4;
98  static const double maxAbsZ = 15.0;
99  static const double maxd0 = 2.0;
100 
101  //count verticies
102 
103  for (reco::VertexCollection::const_iterator vit = vertices.begin(); vit != vertices.end(); ++vit) {
104  if (vit->ndof() > minNDOF && ((maxAbsZ <= 0) || fabs(vit->z()) <= maxAbsZ) &&
105  ((maxd0 <= 0) || fabs(vit->position().rho()) <= maxd0))
106  nvertex++;
107  }
108  } else {
109  nvertex = 1;
110  }
111 
112  // create return data
113  auto retdata1 = std::make_unique<reco::RecoEcalCandidateCollection>();
114 
115  algo_.produce(super_clus, *hf_assoc, *retdata1, nvertex);
116 
117  e.put(std::move(retdata1));
118 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:171
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
HFRecoEcalCandidateProducer(edm::ParameterSet const &conf)
void produce(const edm::Handle< reco::SuperClusterCollection > &SuperClusters, const reco::HFEMClusterShapeAssociationCollection &AssocShapes, reco::RecoEcalCandidateCollection &RecoECand, int nvtx) const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
Cut2D_(conf.getParameter< double >("intercept2DCut"))
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
def move
Definition: eostools.py:511
defaultSlope2D_((Cut2D_<=0.83)?(0.475):((Cut2D_ > 0.83 &&Cut2D_<=0.9)?(0.275):(0.2)))
edm::AssociationMap< edm::OneToOne< SuperClusterCollection, HFEMClusterShapeCollection > > HFEMClusterShapeAssociationCollection
T const * product() const
Definition: Handle.h:70
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void produce(edm::Event &e, edm::EventSetup const &iSetup) override
hfvars_(HFDBversion_, HFDBvector_)