CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
L1JPTOffsetCorrectorImpl.cc
Go to the documentation of this file.
1 // Implementation of class L1JPTOffsetCorrector.
2 // L1JPTOffset jet corrector class.
3 
7 
20 
21 using namespace std;
22 
25  useOffset_(false)
26 {
27  auto const& offsetService = fConfig.getParameter<edm::InputTag>("offsetService");
28  if(not offsetService.label().empty()) {
29  useOffset_ =true;
30  offsetCorrectorToken_ = fCollector.consumes<reco::JetCorrector>(offsetService);
31  }
32 }
33 
34 std::unique_ptr<reco::JetCorrectorImpl>
36  reco::JetCorrector const* offset = nullptr;
37  if(useOffset_) {
39  fEvent.getByToken(offsetCorrectorToken_,hOffset);
40  offset = &(*hOffset);
41  }
42  auto calculator = getCalculator(fSetup,
43  [](std::string const& level)
44  {
45  if ( level != "L1JPTOffset") {
46  throw cms::Exception("L1OffsetCorrector")<<" correction level: "<<level<<" is not L1JPTOffset";
47  }
48  });
49  return std::unique_ptr<reco::JetCorrectorImpl>( new L1JPTOffsetCorrectorImpl(calculator,offset) );
50 }
51 
52 void
54 {
56  addToDescription(desc);
57  desc.add<edm::InputTag>("offsetService");
58  iDescriptions.addDefault(desc);
59 }
60 
61 
62 //------------------------------------------------------------------------
63 //--- L1OffsetCorrectorImpl constructor ------------------------------------------
64 //------------------------------------------------------------------------
65 L1JPTOffsetCorrectorImpl::L1JPTOffsetCorrectorImpl(std::shared_ptr<FactorizedJetCorrectorCalculator const> corrector,
66  const reco::JetCorrector* offsetService):
67  offsetService_(offsetService),
68  corrector_(corrector)
69 {
70 }
71 //------------------------------------------------------------------------
72 //--- L1OffsetCorrectorImpl destructor -------------------------------------------
73 //------------------------------------------------------------------------
74 
75 //------------------------------------------------------------------------
76 //--- Returns correction for a given 4-vector ----------------------------
77 //------------------------------------------------------------------------
79 {
80  throw cms::Exception("EventRequired")
81  <<"Wrong interface correction(LorentzVector), event required!";
82  return 1.0;
83 }
84 //------------------------------------------------------------------------
85 //--- Returns correction for a given jet ---------------------------------
86 //------------------------------------------------------------------------
88 {
89  double result = 1.;
90  const reco::JPTJet& jptjet = dynamic_cast <const reco::JPTJet&> (fJet);
91  edm::RefToBase<reco::Jet> jptjetRef = jptjet.getCaloJetRef();
92  reco::CaloJet const * rawcalojet = dynamic_cast<reco::CaloJet const *>( &* jptjetRef);
93  //------ access the offset correction service ----------------
94  double offset = 1.0;
95  if (offsetService_) {
96  offset = offsetService_->correction(*rawcalojet);
97  }
98  //------ calculate the correction for the JPT jet ------------
99  TLorentzVector JPTrawP4(rawcalojet->px(),rawcalojet->py(),rawcalojet->pz(),rawcalojet->energy());
101  values.setJPTrawP4(JPTrawP4);
102  values.setJPTrawOff(offset);
103  values.setJetE(fJet.energy());
104  result = corrector_->getCorrection(values);
105  return result;
106 }
107 
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
L1JPTOffsetCorrectorImplMaker(edm::ParameterSet const &, edm::ConsumesCollector)
Jets made from CaloTowers.
Definition: CaloJet.h:29
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
Base class for all types of Jets.
Definition: Jet.h:20
static void fillDescriptions(edm::ConfigurationDescriptions &iDescriptions)
reco::Particle::LorentzVector LorentzVector
double correction(const LorentzVector &fJet) const
get correction using Jet information only
Definition: JetCorrector.h:47
const edm::RefToBase< reco::Jet > & getCaloJetRef() const
Definition: JPTJet.h:130
const reco::JetCorrector * offsetService_
virtual double correction(const LorentzVector &fJet) const override
get correction using Jet information only
virtual double energy() const
energy
void addDefault(ParameterSetDescription const &psetDescription)
std::shared_ptr< FactorizedJetCorrectorCalculator const > corrector_
L1JPTOffsetCorrectorImpl(std::shared_ptr< FactorizedJetCorrectorCalculator const > corrector, const reco::JetCorrector *offsetService)
tuple corrector
Definition: mvaPFMET_cff.py:86
Jets made from CaloJets corrected for ZSP and tracks.
Definition: JPTJet.h:29
tuple result
Definition: query.py:137
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::unique_ptr< reco::JetCorrectorImpl > make(edm::Event const &, edm::EventSetup const &)
virtual double px() const
x coordinate of momentum vector
virtual double pz() const
z coordinate of momentum vector
edm::EDGetTokenT< reco::JetCorrector > offsetCorrectorToken_
static void addToDescription(edm::ParameterSetDescription &iDescription)
tuple level
Definition: testEve_cfg.py:34
volatile std::atomic< bool > shutdown_flag false
std::shared_ptr< FactorizedJetCorrectorCalculator const > getCalculator(edm::EventSetup const &, std::function< void(std::string const &)> levelCheck)
virtual double py() const
y coordinate of momentum vector