CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HFRecoEcalCandidateProducer.cc
Go to the documentation of this file.
1 
9 #include <iostream>
10 #include <vector>
11 #include <memory>
12 
13 // Framework
20 //
26 
30 
32  defaultDB_(std::vector<double>()),
33  hfclustersSC_(consumes<reco::SuperClusterCollection>(conf.getParameter<edm::InputTag>("hfclusters"))),
34  hfclustersHFEM_(consumes<reco::HFEMClusterShapeAssociationCollection>(conf.getParameter<edm::InputTag>("hfclusters"))),
35 
36  HFDBversion_(conf.existsAs<int>("HFDBversion") ? conf.getParameter<int>("HFDBversion"):99),//do nothing
37  HFDBvector_(conf.existsAs<std::vector<double> >("HFDBvector") ? conf.getParameter<std::vector<double> >("HFDBvector"):defaultDB_),
38  doPU_(false),
39  Cut2D_(conf.getParameter<double>("intercept2DCut")),
40  defaultSlope2D_((Cut2D_<=0.83)?(0.475):((Cut2D_>0.83 && Cut2D_<=0.9)?(0.275):(0.2))),//fix for hlt unable to add slope variable now
41  hfvars_(HFDBversion_,HFDBvector_),
42  algo_(conf.existsAs<bool>("Correct") ? conf.getParameter<bool>("Correct") : true,
43  conf.getParameter<double>("e9e25Cut"),
44  conf.getParameter<double>("intercept2DCut"),
45  conf.existsAs<double>("intercept2DSlope") ? conf.getParameter<double>("intercept2DSlope") : defaultSlope2D_,
46  conf.getParameter<std::vector<double> >("e1e9Cut"),
47  conf.getParameter<std::vector<double> >("eCOREe9Cut"),
48  conf.getParameter<std::vector<double> >("eSeLCut"),
49  hfvars_)
50 {
51  if(conf.existsAs<edm::InputTag>("VertexCollection")){
52  vertices_ = consumes<reco::VertexCollection>(conf.getParameter<edm::InputTag>("VertexCollection"));
53  }else vertices_ = consumes<reco::VertexCollection>(edm::InputTag("offlinePrimaryVertices"));
54 
55  produces<reco::RecoEcalCandidateCollection>();
56 
57 }
58 
60 
61 
64 
65  e.getByToken(hfclustersSC_,super_clus);
66  e.getByToken(hfclustersHFEM_,hf_assoc);
67 
68  int nvertex = 0;
69  if(HFDBversion_!=99){
71  e.getByToken(vertices_, pvHandle);
72  const reco::VertexCollection & vertices = *pvHandle.product();
73  static const int minNDOF = 4;
74  static const double maxAbsZ = 15.0;
75  static const double maxd0 = 2.0;
76 
77  //count verticies
78 
79  for(reco::VertexCollection::const_iterator vit = vertices.begin(); vit != vertices.end(); ++vit){
80  if(vit->ndof() > minNDOF && ((maxAbsZ <= 0) || fabs(vit->z()) <= maxAbsZ) && ((maxd0 <= 0) || fabs(vit->position().rho()) <= maxd0))
81  nvertex++;
82  }
83  }else{
84  nvertex = 1;
85  }
86 
87 
88 
89  // create return data
90  std::auto_ptr<reco::RecoEcalCandidateCollection> retdata1(new reco::RecoEcalCandidateCollection());
91 
92 
93  algo_.produce(super_clus,*hf_assoc,*retdata1,nvertex);
94 
95  e.put(retdata1);
96 
97 }
98 
99 
100 
101 
102 
103 
104 
105 
106 
107 
108 
109 
T getParameter(std::string const &) const
algo_(conf.existsAs< bool >("Correct")?conf.getParameter< bool >("Correct"):true, conf.getParameter< double >("e9e25Cut"), conf.getParameter< double >("intercept2DCut"), conf.existsAs< bool >("intercept2DSlope")?conf.getParameter< double >("intercept2DSlope"):defaultSlope2D_, conf.getParameter< std::vector< double > >("e1e9Cut"), conf.getParameter< std::vector< double > >("eCOREe9Cut"), conf.getParameter< std::vector< double > >("eSeLCut"), hfvars_)
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:184
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
HFRecoEcalCandidateProducer(edm::ParameterSet const &conf)
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
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
tuple conf
Definition: dbtoconf.py:185
edm::AssociationMap< edm::OneToOne< SuperClusterCollection, HFEMClusterShapeCollection > > HFEMClusterShapeAssociationCollection
void produce(const edm::Handle< reco::SuperClusterCollection > &SuperClusters, const reco::HFEMClusterShapeAssociationCollection &AssocShapes, reco::RecoEcalCandidateCollection &RecoECand, int nvtx)
T const * product() const
Definition: Handle.h:81
std::vector< RecoEcalCandidate > RecoEcalCandidateCollection
collectin of RecoEcalCandidate objects
volatile std::atomic< bool > shutdown_flag false
defaultSlope2D_((Cut2D_<=0.83)?(0.475):((Cut2D_ >0.83 &&Cut2D_<=0.9)?(0.275):(0.2)))
virtual void produce(edm::Event &e, edm::EventSetup const &iSetup)
hfvars_(HFDBversion_, HFDBvector_)