CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/PhysicsTools/JetMCAlgos/plugins/GenJetBCEnergyRatio.cc

Go to the documentation of this file.
00001 // 
00002 // Plugin to store B and C ratio for a GenJet in the event
00003 // Author: Attilio  
00004 // Date: 05.10.2007
00005 //
00006 
00007 //=======================================================================
00008 
00009 // user include files
00010 #include "FWCore/Framework/interface/EDProducer.h"
00011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSetfwd.h"
00013 #include "FWCore/Utilities/interface/InputTag.h"
00014  
00015 #include "FWCore/Framework/interface/Event.h"
00016 #include "FWCore/Framework/interface/EventSetup.h"
00017 #include "FWCore/Framework/interface/MakerMacros.h"
00018 #include "FWCore/Framework/interface/ESHandle.h"
00019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00020 
00021 #include "DataFormats/JetReco/interface/Jet.h"
00022 #include "DataFormats/JetReco/interface/GenJetCollection.h"
00023 #include "SimDataFormats/JetMatching/interface/JetFlavour.h"
00024 
00025 #include "DataFormats/Common/interface/Ref.h"
00026 #include "DataFormats/Candidate/interface/Candidate.h"
00027 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00028 #include "DataFormats/JetReco/interface/JetFloatAssociation.h"
00029 #include "DataFormats/Math/interface/Point3D.h"
00030 #include "DataFormats/Math/interface/LorentzVector.h"
00031 
00032 #include "DataFormats/JetReco/interface/GenJet.h"
00033 #include "PhysicsTools/JetMCUtils/interface/JetMCTag.h"
00034 #include "PhysicsTools/JetMCUtils/interface/CandMCTag.h"
00035 
00036 #include <memory>
00037 #include <string>
00038 #include <iostream>
00039 #include <vector>
00040 #include <Math/VectorUtil.h>
00041 #include <TMath.h>
00042 
00043 using namespace std;
00044 using namespace reco;
00045 using namespace edm;
00046 using namespace ROOT::Math::VectorUtil;
00047 using namespace JetMCTagUtils;
00048 using namespace CandMCTagUtils;
00049 
00050 class GenJetBCEnergyRatio : public edm::EDProducer 
00051 {
00052   public:
00053     GenJetBCEnergyRatio( const edm::ParameterSet & );
00054     ~GenJetBCEnergyRatio();
00055 
00056     typedef reco::JetFloatAssociation::Container JetBCEnergyRatioCollection;
00057 
00058   private:
00059     virtual void produce(edm::Event&, const edm::EventSetup& );
00060     Handle< View <Jet> > genjets;
00061     edm::InputTag m_genjetsSrc;
00062 
00063 };
00064 
00065 //=========================================================================
00066 
00067 GenJetBCEnergyRatio::GenJetBCEnergyRatio( const edm::ParameterSet& iConfig )
00068 {
00069     produces<JetBCEnergyRatioCollection>("bRatioCollection");
00070     produces<JetBCEnergyRatioCollection>("cRatioCollection");
00071     m_genjetsSrc = iConfig.getParameter<edm::InputTag>("genJets");
00072 }
00073 
00074 //=========================================================================
00075 
00076 GenJetBCEnergyRatio::~GenJetBCEnergyRatio() 
00077 {
00078 }
00079 
00080 // ------------ method called to produce the data  ------------
00081 
00082 void GenJetBCEnergyRatio::produce( Event& iEvent, const EventSetup& iEs ) 
00083 {
00084   iEvent.getByLabel(m_genjetsSrc, genjets);
00085 
00086   typedef edm::RefToBase<reco::Jet> JetRef;
00087 
00088   JetBCEnergyRatioCollection * jtc1;
00089   JetBCEnergyRatioCollection * jtc2;
00090 
00091   if (genjets.product()->size() > 0) {
00092     const JetRef jj = genjets->refAt(0);
00093     jtc1 = new JetBCEnergyRatioCollection(RefToBaseProd<Jet>(jj));
00094     jtc2 = new JetBCEnergyRatioCollection(RefToBaseProd<Jet>(jj));
00095   } else {
00096     jtc1 = new JetBCEnergyRatioCollection();
00097     jtc2 = new JetBCEnergyRatioCollection();
00098   }
00099 
00100   std::auto_ptr<JetBCEnergyRatioCollection> bRatioColl(jtc1);
00101   std::auto_ptr<JetBCEnergyRatioCollection> cRatioColl(jtc2);
00102 
00103   for( size_t j = 0; j != genjets->size(); ++j ) {
00104 
00105     float bRatio = EnergyRatioFromBHadrons( (*genjets)[j] );
00106     float cRatio = EnergyRatioFromCHadrons( (*genjets)[j] );
00107 
00108     const JetRef & aJet = genjets->refAt(j) ;
00109 
00110     JetFloatAssociation::setValue(*bRatioColl, aJet, bRatio);
00111     JetFloatAssociation::setValue(*cRatioColl, aJet, cRatio);
00112 
00113   }
00114 
00115 
00116   iEvent.put(bRatioColl, "bRatioCollection");
00117   iEvent.put(cRatioColl, "cRatioCollection");
00118 
00119 }
00120 
00121 //define this as a plug-in
00122 DEFINE_FWK_MODULE(GenJetBCEnergyRatio);
00123