Go to the documentation of this file.00001 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00002
00003 #include "DataFormats/Common/interface/Handle.h"
00004
00005 #include "DataFormats/JetReco/interface/CaloJet.h"
00006 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00007
00008 #include "FWCore/Framework/interface/MakerMacros.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010
00011 #include "FWCore/Utilities/interface/InputTag.h"
00012
00013 #include "DQM/DataScouting/interface/DiJetVarProducer.h"
00014
00015 #include "TVector3.h"
00016 #include "TLorentzVector.h"
00017
00018 #include <memory>
00019 #include <vector>
00020
00021
00022
00023
00024 DiJetVarProducer::DiJetVarProducer(const edm::ParameterSet& iConfig) :
00025 inputJetTag_ (iConfig.getParameter<edm::InputTag>("inputJetTag")),
00026 wideJetDeltaR_ (iConfig.getParameter<double>("wideJetDeltaR"))
00027 {
00028
00029
00030 produces<std::vector<math::PtEtaPhiMLorentzVector> >("widejets");
00031
00032 LogDebug("") << "Input Jet Tag: " << inputJetTag_.encode() << " ";
00033 LogDebug("") << "Radius Parameter Wide Jet: " << wideJetDeltaR_ << ".";
00034 }
00035
00036 DiJetVarProducer::~DiJetVarProducer()
00037 {
00038 }
00039
00040
00041 void
00042 DiJetVarProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00043 {
00044 using namespace std;
00045 using namespace edm;
00046 using namespace reco;
00047
00048
00049
00050 std::auto_ptr<std::vector<math::PtEtaPhiMLorentzVector> > widejets(new std::vector<math::PtEtaPhiMLorentzVector>);
00051
00052
00053 edm::Handle<reco::CaloJetCollection> calojets_handle;
00054 iEvent.getByLabel(inputJetTag_,calojets_handle);
00055
00056
00057
00058
00059 if( calojets_handle->size() >=2 )
00060 {
00061 TLorentzVector wj1_tmp;
00062 TLorentzVector wj2_tmp;
00063 TLorentzVector wj1;
00064 TLorentzVector wj2;
00065 TLorentzVector wdijet;
00066
00067
00068
00069
00070
00071
00072
00073
00074 TLorentzVector jet1, jet2;
00075
00076 reco::CaloJetCollection::const_iterator j1 = calojets_handle->begin();
00077 reco::CaloJetCollection::const_iterator j2 = j1; ++j2;
00078
00079 jet1.SetPtEtaPhiM(j1->pt(),j1->eta(),j1->phi(),j1->mass());
00080 jet2.SetPtEtaPhiM(j2->pt(),j2->eta(),j2->phi(),j2->mass());
00081
00082
00083
00084
00085
00086 for(reco::CaloJetCollection::const_iterator it = calojets_handle->begin(); it != calojets_handle->end(); ++it)
00087 {
00088 TLorentzVector currentJet;
00089 currentJet.SetPtEtaPhiM(it->pt(),it->eta(),it->phi(),it->mass());
00090
00091 double DeltaR1 = currentJet.DeltaR(jet1);
00092 double DeltaR2 = currentJet.DeltaR(jet2);
00093
00094 if(DeltaR1 < DeltaR2 && DeltaR1 < wideJetDeltaR_) {
00095 wj1_tmp += currentJet;
00096 } else if(DeltaR2 < wideJetDeltaR_) {
00097 wj2_tmp += currentJet;
00098 }
00099 }
00100
00101
00102 if( wj1_tmp.Pt() > wj2_tmp.Pt() )
00103 {
00104 wj1 = wj1_tmp;
00105 wj2 = wj2_tmp;
00106 }
00107 else
00108 {
00109 wj1 = wj2_tmp;
00110 wj2 = wj1_tmp;
00111 }
00112
00113
00114 wdijet = wj1 + wj2;
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128 math::PtEtaPhiMLorentzVector wj1math(wj1.Pt(), wj1.Eta(), wj1.Phi(), wj1.M());
00129 math::PtEtaPhiMLorentzVector wj2math(wj2.Pt(), wj2.Eta(), wj2.Phi(), wj2.M());
00130 widejets->push_back( wj1math );
00131 widejets->push_back( wj2math );
00132 }
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143 iEvent.put(widejets, "widejets");
00144
00145 }
00146
00147 DEFINE_FWK_MODULE(DiJetVarProducer);