CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DiJetVarProducer.cc
Go to the documentation of this file.
2 
4 
7 
10 
12 
14 
15 #include "TVector3.h"
16 #include "TLorentzVector.h"
17 
18 #include <memory>
19 #include <vector>
20 
21 //
22 // constructors and destructor
23 //
25  inputJetTag_ (iConfig.getParameter<edm::InputTag>("inputJetTag")),
26  wideJetDeltaR_ (iConfig.getParameter<double>("wideJetDeltaR"))
27 {
28  //register your products
29  //produces<std::vector<double> >("dijetvariables");
30  produces<std::vector<math::PtEtaPhiMLorentzVector> >("widejets");
31 
32  LogDebug("") << "Input Jet Tag: " << inputJetTag_.encode() << " ";
33  LogDebug("") << "Radius Parameter Wide Jet: " << wideJetDeltaR_ << ".";
34 }
35 
37 {
38 }
39 
40 // ------------ method called to produce the data ------------
41 void
43 {
44  using namespace std;
45  using namespace edm;
46  using namespace reco;
47 
48  // ## The output collections
49  //std::auto_ptr<std::vector<double> > dijetvariables(new std::vector<double>);
50  std::auto_ptr<std::vector<math::PtEtaPhiMLorentzVector> > widejets(new std::vector<math::PtEtaPhiMLorentzVector>);
51 
52  // ## Get jet collection
54  iEvent.getByLabel(inputJetTag_,calojets_handle);
55  //cout << "size: " << calojets_handle->size() << endl;
56 
57  // ## Wide Jet Algorithm
58  // At least two jets
59  if( calojets_handle->size() >=2 )
60  {
61  TLorentzVector wj1_tmp;
62  TLorentzVector wj2_tmp;
63  TLorentzVector wj1;
64  TLorentzVector wj2;
65  TLorentzVector wdijet;
66 
67  // // Loop over all the input jets
68  // for(reco::CaloJetCollection::const_iterator it = calojets_handle->begin(); it != calojets_handle->end(); ++it)
69  // {
70  // cout << "jet: " << it->pt() << " " << it->eta() << " " << it->phi() << endl;
71  // }
72 
73  // Find two leading jets
74  TLorentzVector jet1, jet2;
75 
76  reco::CaloJetCollection::const_iterator j1 = calojets_handle->begin();
77  reco::CaloJetCollection::const_iterator j2 = j1; ++j2;
78 
79  jet1.SetPtEtaPhiM(j1->pt(),j1->eta(),j1->phi(),j1->mass());
80  jet2.SetPtEtaPhiM(j2->pt(),j2->eta(),j2->phi(),j2->mass());
81 
82  //cout << "j1: " << jet1.Pt() << " " << jet1.Eta() << " " << jet1.Phi() << endl;
83  //cout << "j2: " << jet2.Pt() << " " << jet2.Eta() << " " << jet2.Phi() << endl;
84 
85  // Create wide jets (radiation recovery algorithm)
86  for(reco::CaloJetCollection::const_iterator it = calojets_handle->begin(); it != calojets_handle->end(); ++it)
87  {
88  TLorentzVector currentJet;
89  currentJet.SetPtEtaPhiM(it->pt(),it->eta(),it->phi(),it->mass());
90 
91  double DeltaR1 = currentJet.DeltaR(jet1);
92  double DeltaR2 = currentJet.DeltaR(jet2);
93 
94  if(DeltaR1 < DeltaR2 && DeltaR1 < wideJetDeltaR_) {
95  wj1_tmp += currentJet;
96  } else if(DeltaR2 < wideJetDeltaR_) {
97  wj2_tmp += currentJet;
98  }
99  }
100 
101  // Re-order the wide jets in pT
102  if( wj1_tmp.Pt() > wj2_tmp.Pt() )
103  {
104  wj1 = wj1_tmp;
105  wj2 = wj2_tmp;
106  }
107  else
108  {
109  wj1 = wj2_tmp;
110  wj2 = wj1_tmp;
111  }
112 
113  // Create dijet system
114  wdijet = wj1 + wj2;
115 
116  // cout << "j1 wide: " << wj1.Pt() << " " << wj1.Eta() << " " << wj1.Phi() << " " << wj1.M() << endl;
117  // cout << "j2 wide: " << wj2.Pt() << " " << wj2.Eta() << " " << wj2.Phi() << " " << wj2.M() << endl;
118  // cout << "MJJWide: " << wdijet.M() << endl;
119  // cout << "DeltaEtaJJWide: " << fabs(wj1.Eta()-wj2.Eta()) << endl;
120  // cout << "DeltaPhiJJWide: " << fabs(wj1.DeltaPhi(wj2)) << endl;
121 
122  // // Put variables in the container
123  // dijetvariables->push_back( wdijet.M() ); //0 = MJJWide
124  // dijetvariables->push_back( fabs(wj1.Eta()-wj2.Eta()) ); //1 = DeltaEtaJJWide
125  // dijetvariables->push_back( fabs(wj1.DeltaPhi(wj2)) ); //2 = DeltaPhiJJWide
126 
127  // Put widejets in the container
128  math::PtEtaPhiMLorentzVector wj1math(wj1.Pt(), wj1.Eta(), wj1.Phi(), wj1.M());
129  math::PtEtaPhiMLorentzVector wj2math(wj2.Pt(), wj2.Eta(), wj2.Phi(), wj2.M());
130  widejets->push_back( wj1math );
131  widejets->push_back( wj2math );
132  }
133  // else
134  // {
135  // // Put variables in the container
136  // dijetvariables->push_back( -1 ); //0 = MJJWide
137  // dijetvariables->push_back( -1 ); //1 = DeltaEtaJJWide
138  // dijetvariables->push_back( -1 ); //2 = DeltaPhiJJWide
139  // }
140 
141  // ## Put objects in the Event
142  //iEvent.put(dijetvariables, "dijetvariables");
143  iEvent.put(widejets, "widejets");
144 
145 }
146 
#define LogDebug(id)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::string encode() const
Definition: InputTag.cc:164
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:26
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
edm::InputTag inputJetTag_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
DiJetVarProducer(const edm::ParameterSet &)
virtual void produce(edm::Event &, const edm::EventSetup &)