Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:23:33

0001 //
0002 // Plugin to store B and C ratio for a GenJet in the event
0003 // Author: Attilio
0004 // Date: 05.10.2007
0005 //
0006 
0007 //=======================================================================
0008 
0009 // user include files
0010 #include "FWCore/Framework/interface/global/EDProducer.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSetfwd.h"
0013 #include "FWCore/Utilities/interface/InputTag.h"
0014 
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "FWCore/Framework/interface/EventSetup.h"
0017 #include "FWCore/Framework/interface/MakerMacros.h"
0018 #include "FWCore/Framework/interface/ESHandle.h"
0019 #include "FWCore/Framework/interface/makeRefToBaseProdFrom.h"
0020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0021 
0022 #include "DataFormats/JetReco/interface/Jet.h"
0023 #include "DataFormats/JetReco/interface/GenJetCollection.h"
0024 #include "DataFormats/JetMatching/interface/JetFlavour.h"
0025 
0026 #include "DataFormats/Common/interface/Ref.h"
0027 #include "DataFormats/Candidate/interface/Candidate.h"
0028 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0029 #include "DataFormats/JetReco/interface/JetFloatAssociation.h"
0030 #include "DataFormats/Math/interface/Point3D.h"
0031 #include "DataFormats/Math/interface/LorentzVector.h"
0032 
0033 #include "DataFormats/JetReco/interface/GenJet.h"
0034 #include "PhysicsTools/JetMCUtils/interface/JetMCTag.h"
0035 #include "PhysicsTools/JetMCUtils/interface/CandMCTag.h"
0036 
0037 #include <memory>
0038 #include <string>
0039 #include <iostream>
0040 #include <vector>
0041 #include <Math/VectorUtil.h>
0042 #include <TMath.h>
0043 
0044 using namespace std;
0045 using namespace reco;
0046 using namespace edm;
0047 using namespace ROOT::Math::VectorUtil;
0048 using namespace JetMCTagUtils;
0049 using namespace CandMCTagUtils;
0050 
0051 class GenJetBCEnergyRatio : public edm::global::EDProducer<> {
0052 public:
0053   GenJetBCEnergyRatio(const edm::ParameterSet&);
0054   ~GenJetBCEnergyRatio() override;
0055 
0056   typedef reco::JetFloatAssociation::Container JetBCEnergyRatioCollection;
0057 
0058 private:
0059   void produce(StreamID, edm::Event&, const edm::EventSetup&) const override;
0060   edm::EDGetTokenT<View<Jet> > m_genjetsSrcToken;
0061 };
0062 
0063 //=========================================================================
0064 
0065 GenJetBCEnergyRatio::GenJetBCEnergyRatio(const edm::ParameterSet& iConfig) {
0066   produces<JetBCEnergyRatioCollection>("bRatioCollection");
0067   produces<JetBCEnergyRatioCollection>("cRatioCollection");
0068   m_genjetsSrcToken = consumes<View<Jet> >(iConfig.getParameter<edm::InputTag>("genJets"));
0069 }
0070 
0071 //=========================================================================
0072 
0073 GenJetBCEnergyRatio::~GenJetBCEnergyRatio() {}
0074 
0075 // ------------ method called to produce the data  ------------
0076 
0077 void GenJetBCEnergyRatio::produce(StreamID, Event& iEvent, const EventSetup& iEs) const {
0078   Handle<View<Jet> > genjets;
0079   iEvent.getByToken(m_genjetsSrcToken, genjets);
0080 
0081   typedef edm::RefToBase<reco::Jet> JetRef;
0082 
0083   JetBCEnergyRatioCollection* jtc1;
0084   JetBCEnergyRatioCollection* jtc2;
0085 
0086   if (!genjets.product()->empty()) {
0087     const JetRef jj = genjets->refAt(0);
0088     jtc1 = new JetBCEnergyRatioCollection(edm::makeRefToBaseProdFrom(jj, iEvent));
0089     jtc2 = new JetBCEnergyRatioCollection(edm::makeRefToBaseProdFrom(jj, iEvent));
0090   } else {
0091     jtc1 = new JetBCEnergyRatioCollection();
0092     jtc2 = new JetBCEnergyRatioCollection();
0093   }
0094 
0095   std::unique_ptr<JetBCEnergyRatioCollection> bRatioColl(jtc1);
0096   std::unique_ptr<JetBCEnergyRatioCollection> cRatioColl(jtc2);
0097 
0098   for (size_t j = 0; j != genjets->size(); ++j) {
0099     float bRatio = EnergyRatioFromBHadrons((*genjets)[j]);
0100     float cRatio = EnergyRatioFromCHadrons((*genjets)[j]);
0101 
0102     const JetRef& aJet = genjets->refAt(j);
0103 
0104     JetFloatAssociation::setValue(*bRatioColl, aJet, bRatio);
0105     JetFloatAssociation::setValue(*cRatioColl, aJet, cRatio);
0106   }
0107 
0108   iEvent.put(std::move(bRatioColl), "bRatioCollection");
0109   iEvent.put(std::move(cRatioColl), "cRatioCollection");
0110 }
0111 
0112 //define this as a plug-in
0113 DEFINE_FWK_MODULE(GenJetBCEnergyRatio);