CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetMatchingMLM.cc
Go to the documentation of this file.
1 #include <functional>
2 #include <algorithm>
3 #include <iostream>
4 #include <vector>
5 #include <memory>
6 #include <string>
7 
8 #include <boost/bind.hpp>
9 
10 #include "Math/GenVector/Cartesian3D.h"
11 #include "Math/GenVector/VectorUtil.h"
12 
13 #include "HepMC/GenEvent.h"
14 #include "HepMC/SimpleVector.h"
15 
18 
22 
23 #include "Matching.h"
24 
25 // #define DEBUG
26 
27 namespace lhef {
28 
29 namespace {
30  template<typename T1, typename T2, typename R>
31  struct DeltaRSeparation : public std::binary_function<T1, T2, R> {
32  double operator () (const T1 &v1_, const T2 &v2_) const
33  {
34  using namespace ROOT::Math;
35  Cartesian3D<R> v1(v1_.px(), v1_.py(), v1_.pz());
36  Cartesian3D<R> v2(v2_.px(), v2_.py(), v2_.pz());
37  return VectorUtil::DeltaR(v1, v2);
38  }
39  };
40 
41  // stupid pointer indirection... ugly specialization
42  template<typename T2, typename R>
43  struct DeltaRSeparation<const HepMC::GenParticle*, T2, R> {
44  double operator () (const HepMC::GenParticle *v1_,
45  const T2 &v2_) const
46  {
47  using namespace ROOT::Math;
48  Cartesian3D<R> v1(v1_->momentum().px(),
49  v1_->momentum().py(),
50  v1_->momentum().pz());
51  Cartesian3D<R> v2(v2_.px(), v2_.py(), v2_.pz());
52  return VectorUtil::DeltaR(v1, v2);
53  }
54  };
55 
56  struct ParticlePtGreater {
57  double operator () (const HepMC::GenParticle *v1,
58  const HepMC::GenParticle *v2)
59  { return v1->momentum().perp() > v2->momentum().perp(); }
60  };
61 
62  inline HepMC::FourVector convert(const JetClustering::Jet &jet)
63  { return HepMC::FourVector(jet.px(), jet.py(), jet.pz(), jet.e()); }
64 } // anonymous namespace
65 
67  JetMatching(params),
68  maxDeltaR(params.getParameter<double>("maxDeltaR")),
69  minJetPt(params.getParameter<double>("jetPtMin")),
70  maxEta(params.getParameter<double>("maxEta")),
71  matchPtFraction(0.75),
72  useEt(params.getParameter<bool>("useEt")),
73  partonInput(new JetInput(params)),
74  jetInput(new JetInput(*partonInput))
75 {
76  partonInput->setPartonicFinalState(false);
77  partonInput->setHardProcessOnly(false);
78 
79  if (params.exists("matchPtFraction"))
81  params.getParameter<double>("matchPtFraction");
82 
83  jetClustering.reset(
84  new JetClustering(params, minJetPt * matchPtFraction));
85 
86  std::string matchMode = params.getParameter<std::string>("matchMode");
87  if (matchMode == "inclusive")
88  this->matchMode = kInclusive;
89  else if (matchMode == "exclusive")
90  this->matchMode = kExclusive;
91  else
92  throw cms::Exception("Configuration")
93  << "Invalid matching mode '" << matchMode
94  << "' specified." << std::endl;
95 }
96 
98 {
99 }
100 
101 std::set<std::string> JetMatchingMLM::capabilities() const
102 {
103  std::set<std::string> result = JetMatching::capabilities();
104  result.insert("matchSummary");
105  return result;
106 }
107 
108 // implements the MLM method - simply reject or accept
109 // use polymorphic JetMatching subclasses when modularizing
110 
111 double JetMatchingMLM::match(const HepMC::GenEvent *partonLevel,
112  const HepMC::GenEvent *finalState,
113  bool showeredFinalState)
114 {
115  JetInput::ParticleVector partons = (*partonInput)(partonLevel);
116  std::sort(partons.begin(), partons.end(), ParticlePtGreater());
117 
119  showeredFinalState ? (*partonInput)(finalState)
120  : (*this->jetInput)(finalState);
121  std::sort(jetInput.begin(), jetInput.end(), ParticlePtGreater());
122 
123  std::vector<JetClustering::Jet> jets = (*jetClustering)(jetInput);
124 
125 #ifdef DEBUG
126  std::cout << "===== Partons:" << std::endl;
127  for(JetClustering::ParticleVector::const_iterator c = partons.begin();
128  c != partons.end(); c++)
129  std::cout << "\tpid = " << (*c)->pdg_id()
130  << ", pt = " << (*c)->momentum().perp()
131  << ", eta = " << (*c)->momentum().eta()
132  << ", phi = " << (*c)->momentum().phi()
133  << std::endl;
134  std::cout << "===== JetInput:" << std::endl;
135  for(JetClustering::ParticleVector::const_iterator c = jetInput.begin();
136  c != jetInput.end(); c++)
137  std::cout << "\tpid = " << (*c)->pdg_id()
138  << ", pt = " << (*c)->momentum().perp()
139  << ", eta = " << (*c)->momentum().eta()
140  << ", phi = " << (*c)->momentum().phi()
141  << std::endl;
142  std::cout << "----- " << jets.size() << " jets:" << std::endl;
143  for(std::vector<JetClustering::Jet>::const_iterator iter = jets.begin();
144  iter != jets.end(); ++iter) {
145  std::cout << "* pt = " << iter->pt()
146  << ", eta = " << iter->eta()
147  << ", phi = " << iter->phi()
148  << std::endl;
149  for(JetClustering::ParticleVector::const_iterator c = iter->constituents().begin();
150  c != iter->constituents().end(); c++)
151  std::cout << "\tpid = " << (*c)->pdg_id()
152  << ", pt = " << (*c)->momentum().perp()
153  << ", eta = " << (*c)->momentum().eta()
154  << ", phi = " << (*c)->momentum().phi()
155  << std::endl;
156  }
157 
158  using boost::bind;
159  std::cout << partons.size() << " partons and "
160  << std::count_if(jets.begin(), jets.end(),
161  bind(std::greater<double>(),
162  bind(&JetClustering::Jet::pt, _1),
163  minJetPt)) << " jets." << std::endl;
164 #endif
165 
166  Matching<double> matching(partons, jets,
167  DeltaRSeparation<JetInput::ParticleVector::value_type,
168  JetClustering::Jet, double>());
169 
171  std::vector<Match> matches =
172  matching.match(
173  std::less<double>(),
174  std::bind2nd(std::less<double>(), maxDeltaR));
175 
176 #ifdef DEBUG
177  for(std::vector<Match>::const_iterator iter = matches.begin();
178  iter != matches.end(); ++iter)
179  std::cout << "\tParton " << iter->index1 << " matches jet "
180  << iter->index2 << " with a Delta_R of "
181  << matching.delta(*iter) << std::endl;
182 #endif
183 
184  unsigned int unmatchedPartons = 0;
185  unsigned int unmatchedJets = 0;
186 
187  matchSummary.clear();
188  for(std::vector<Match>::const_iterator iter = matches.begin();
189  iter != matches.end(); ++iter) {
190  if ((useEt ? jets[iter->index2].et()
191  : jets[iter->index2].pt()) < minJetPt ||
192  std::abs(jets[iter->index2].eta()) > maxEta)
193  unmatchedPartons++;
194  matchSummary.push_back(
195  JetPartonMatch(partons[iter->index1]->momentum(),
196  convert(jets[iter->index2]),
197  matching.delta(*iter),
198  partons[iter->index1]->pdg_id()));
199  }
200 
201  for(Matching<double>::index_type i = 0; i < partons.size(); i++) {
202  if (!matching.isMatched1st(i)) {
203  unmatchedPartons++;
204  matchSummary.push_back(
205  JetPartonMatch(partons[i]->momentum(),
206  partons[i]->pdg_id()));
207  }
208  }
209 
210  for(Matching<double>::index_type i = 0; i < jets.size(); i++) {
211  if (!matching.isMatched2nd(i)) {
212  if ((useEt ? jets[i].et()
213  : jets[i].pt()) >= minJetPt &&
214  std::abs(jets[i].eta()) <= maxEta)
215  unmatchedJets++;
216  matchSummary.push_back(
217  JetPartonMatch(convert(jets[i])));
218  }
219  }
220 
221  switch(matchMode) {
222  case kExclusive:
223  if (!unmatchedJets && !unmatchedPartons)
224  return 1.0;
225  break;
226  case kInclusive:
227  if (!unmatchedPartons)
228  return 1.0;
229  }
230 
231  return 0.0;
232 }
233 
234 } // namespace lhef
const double maxDeltaR
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::set< std::string > capabilities() const
bool exists(std::string const &parameterName) const
checks if a parameter exists
double maxEta
std::auto_ptr< JetInput > partonInput
double match(const HepMC::GenEvent *partonLevel, const HepMC::GenEvent *finalState, bool showeredFinalState)
std::auto_ptr< JetInput > jetInput
std::auto_ptr< JetClustering > jetClustering
vector< PseudoJet > jets
tuple result
Definition: query.py:137
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Container::value_type value_type
SimpleMatrix< Delta >::size_type index_type
Definition: Matching.h:16
const double minJetPt
string const
Definition: compareJSON.py:14
virtual std::set< std::string > capabilities() const
Definition: JetMatching.cc:38
std::vector< const HepMC::GenParticle * > ParticleVector
Definition: JetInput.h:17
tuple cout
Definition: gather_cfg.py:121
std::vector< JetPartonMatch > matchSummary
Definition: JetMatching.h:74
JetMatchingMLM(const edm::ParameterSet &params)