Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-02-14 13:34:07

0001 //
0002 // Package:    EgammaHFProducers
0003 // Class:      HFRecoEcalCandidateProducers
0004 //
0005 /**\class HFRecoEcalCandidateProducers.cc  
0006 */
0007 //
0008 // Original Author:  Kevin Klapoetke University of Minnesota
0009 //         Created:  Wed 26 Sept 2007
0010 // $Id:
0011 //
0012 //
0013 
0014 #include "DataFormats/Common/interface/Handle.h"
0015 #include "DataFormats/EgammaReco/interface/HFEMClusterShape.h"
0016 #include "DataFormats/EgammaReco/interface/HFEMClusterShapeAssociation.h"
0017 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0018 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0019 #include "FWCore/Framework/interface/Event.h"
0020 #include "FWCore/Framework/interface/EventSetup.h"
0021 #include "FWCore/Framework/interface/global/EDProducer.h"
0022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0023 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0024 #include "FWCore/Utilities/interface/Exception.h"
0025 
0026 #include "HFRecoEcalCandidateAlgo.h"
0027 #include "HFValueStruct.h"
0028 
0029 #include <vector>
0030 #include <memory>
0031 
0032 class HLTHFRecoEcalCandidateProducer : public edm::global::EDProducer<> {
0033 public:
0034   explicit HLTHFRecoEcalCandidateProducer(edm::ParameterSet const& conf);
0035   void produce(edm::StreamID, edm::Event&, edm::EventSetup const&) const override;
0036 
0037 private:
0038   const edm::InputTag hfclusters_, vertices_;
0039   const int HFDBversion_;
0040   const std::vector<double> HFDBvector_;
0041   const double Cut2D_;
0042   const double defaultSlope2D_;
0043   const reco::HFValueStruct hfvars_;
0044   const HFRecoEcalCandidateAlgo algo_;
0045 };
0046 
0047 HLTHFRecoEcalCandidateProducer::HLTHFRecoEcalCandidateProducer(edm::ParameterSet const& conf)
0048     : hfclusters_(conf.getParameter<edm::InputTag>("hfclusters")),
0049       HFDBversion_(conf.existsAs<bool>("HFDBversion") ? conf.getParameter<int>("HFDBversion") : 99),  //do nothing
0050       HFDBvector_(conf.existsAs<bool>("HFDBvector") ? conf.getParameter<std::vector<double> >("HFDBvector")
0051                                                     : std::vector<double>{}),
0052       Cut2D_(conf.getParameter<double>("intercept2DCut")),
0053       defaultSlope2D_(
0054           (Cut2D_ <= 0.83)
0055               ? (0.475)
0056               : ((Cut2D_ > 0.83 && Cut2D_ <= 0.9) ? (0.275) : (0.2))),  //fix for hlt unable to add slope variable now
0057       hfvars_(HFDBversion_, HFDBvector_),
0058       algo_(conf.existsAs<bool>("Correct") ? conf.getParameter<bool>("Correct") : true,
0059             conf.getParameter<double>("e9e25Cut"),
0060             conf.getParameter<double>("intercept2DCut"),
0061             conf.existsAs<bool>("intercept2DSlope") ? conf.getParameter<double>("intercept2DSlope") : defaultSlope2D_,
0062             conf.getParameter<std::vector<double> >("e1e9Cut"),
0063             conf.getParameter<std::vector<double> >("eCOREe9Cut"),
0064             conf.getParameter<std::vector<double> >("eSeLCut"),
0065             hfvars_) {
0066   produces<reco::RecoEcalCandidateCollection>();
0067 }
0068 
0069 void HLTHFRecoEcalCandidateProducer::produce(edm::StreamID sid, edm::Event& e, edm::EventSetup const& iSetup) const {
0070   edm::Handle<reco::SuperClusterCollection> super_clus;
0071   edm::Handle<reco::HFEMClusterShapeAssociationCollection> hf_assoc;
0072 
0073   e.getByLabel(hfclusters_, super_clus);
0074   e.getByLabel(hfclusters_, hf_assoc);
0075 
0076   int nvertex = 1;
0077 
0078   // create return data
0079   auto retdata1 = std::make_unique<reco::RecoEcalCandidateCollection>();
0080 
0081   algo_.produce(super_clus, *hf_assoc, *retdata1, nvertex);
0082 
0083   e.put(std::move(retdata1));
0084 }