00001 #include <functional>
00002 #include <algorithm>
00003 #include <iostream>
00004 #include <vector>
00005 #include <memory>
00006 #include <string>
00007
00008 #include <boost/bind.hpp>
00009
00010 #include <Math/GenVector/Cartesian3D.h>
00011 #include <Math/GenVector/VectorUtil.h>
00012
00013 #include <HepMC/GenEvent.h>
00014 #include <HepMC/SimpleVector.h>
00015
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 #include "FWCore/Utilities/interface/Exception.h"
00018
00019 #include "GeneratorInterface/LHEInterface/interface/JetInput.h"
00020 #include "GeneratorInterface/LHEInterface/interface/JetClustering.h"
00021 #include "GeneratorInterface/LHEInterface/interface/JetMatchingMLM.h"
00022
00023 #include "Matching.h"
00024
00025
00026
00027 namespace lhef {
00028
00029 namespace {
00030 template<typename T1, typename T2, typename R>
00031 struct DeltaRSeparation : public std::binary_function<T1, T2, R> {
00032 double operator () (const T1 &v1_, const T2 &v2_) const
00033 {
00034 using namespace ROOT::Math;
00035 Cartesian3D<R> v1(v1_.px(), v1_.py(), v1_.pz());
00036 Cartesian3D<R> v2(v2_.px(), v2_.py(), v2_.pz());
00037 return VectorUtil::DeltaR(v1, v2);
00038 }
00039 };
00040
00041
00042 template<typename T2, typename R>
00043 struct DeltaRSeparation<const HepMC::GenParticle*, T2, R> {
00044 double operator () (const HepMC::GenParticle *v1_,
00045 const T2 &v2_) const
00046 {
00047 using namespace ROOT::Math;
00048 Cartesian3D<R> v1(v1_->momentum().px(),
00049 v1_->momentum().py(),
00050 v1_->momentum().pz());
00051 Cartesian3D<R> v2(v2_.px(), v2_.py(), v2_.pz());
00052 return VectorUtil::DeltaR(v1, v2);
00053 }
00054 };
00055
00056 struct ParticlePtGreater {
00057 double operator () (const HepMC::GenParticle *v1,
00058 const HepMC::GenParticle *v2)
00059 { return v1->momentum().perp() > v2->momentum().perp(); }
00060 };
00061
00062 inline HepMC::FourVector convert(const JetClustering::Jet &jet)
00063 { return HepMC::FourVector(jet.px(), jet.py(), jet.pz(), jet.e()); }
00064 }
00065
00066 JetMatchingMLM::JetMatchingMLM(const edm::ParameterSet ¶ms) :
00067 JetMatching(params),
00068 maxDeltaR(params.getParameter<double>("maxDeltaR")),
00069 minJetPt(params.getParameter<double>("jetPtMin")),
00070 maxEta(params.getParameter<double>("maxEta")),
00071 matchPtFraction(0.75),
00072 useEt(params.getParameter<bool>("useEt")),
00073 partonInput(new JetInput(params)),
00074 jetInput(new JetInput(*partonInput))
00075 {
00076 partonInput->setPartonicFinalState(false);
00077 partonInput->setHardProcessOnly(false);
00078
00079 if (params.exists("matchPtFraction"))
00080 matchPtFraction =
00081 params.getParameter<double>("matchPtFraction");
00082
00083 jetClustering.reset(
00084 new JetClustering(params, minJetPt * matchPtFraction));
00085
00086 std::string matchMode = params.getParameter<std::string>("matchMode");
00087 if (matchMode == "inclusive")
00088 this->matchMode = kInclusive;
00089 else if (matchMode == "exclusive")
00090 this->matchMode = kExclusive;
00091 else
00092 throw cms::Exception("Configuration")
00093 << "Invalid matching mode '" << matchMode
00094 << "' specified." << std::endl;
00095 }
00096
00097 JetMatchingMLM::~JetMatchingMLM()
00098 {
00099 }
00100
00101 std::set<std::string> JetMatchingMLM::capabilities() const
00102 {
00103 std::set<std::string> result = JetMatching::capabilities();
00104 result.insert("matchSummary");
00105 return result;
00106 }
00107
00108
00109
00110
00111 double JetMatchingMLM::match(const HepMC::GenEvent *partonLevel,
00112 const HepMC::GenEvent *finalState,
00113 bool showeredFinalState)
00114 {
00115 JetInput::ParticleVector partons = (*partonInput)(partonLevel);
00116 std::sort(partons.begin(), partons.end(), ParticlePtGreater());
00117
00118 JetInput::ParticleVector jetInput =
00119 showeredFinalState ? (*partonInput)(finalState)
00120 : (*this->jetInput)(finalState);
00121 std::sort(jetInput.begin(), jetInput.end(), ParticlePtGreater());
00122
00123 std::vector<JetClustering::Jet> jets = (*jetClustering)(jetInput);
00124
00125 #ifdef DEBUG
00126 std::cout << "===== Partons:" << std::endl;
00127 for(JetClustering::ParticleVector::const_iterator c = partons.begin();
00128 c != partons.end(); c++)
00129 std::cout << "\tpid = " << (*c)->pdg_id()
00130 << ", pt = " << (*c)->momentum().perp()
00131 << ", eta = " << (*c)->momentum().eta()
00132 << ", phi = " << (*c)->momentum().phi()
00133 << std::endl;
00134 std::cout << "===== JetInput:" << std::endl;
00135 for(JetClustering::ParticleVector::const_iterator c = jetInput.begin();
00136 c != jetInput.end(); c++)
00137 std::cout << "\tpid = " << (*c)->pdg_id()
00138 << ", pt = " << (*c)->momentum().perp()
00139 << ", eta = " << (*c)->momentum().eta()
00140 << ", phi = " << (*c)->momentum().phi()
00141 << std::endl;
00142 std::cout << "----- " << jets.size() << " jets:" << std::endl;
00143 for(std::vector<JetClustering::Jet>::const_iterator iter = jets.begin();
00144 iter != jets.end(); ++iter) {
00145 std::cout << "* pt = " << iter->pt()
00146 << ", eta = " << iter->eta()
00147 << ", phi = " << iter->phi()
00148 << std::endl;
00149 for(JetClustering::ParticleVector::const_iterator c = iter->constituents().begin();
00150 c != iter->constituents().end(); c++)
00151 std::cout << "\tpid = " << (*c)->pdg_id()
00152 << ", pt = " << (*c)->momentum().perp()
00153 << ", eta = " << (*c)->momentum().eta()
00154 << ", phi = " << (*c)->momentum().phi()
00155 << std::endl;
00156 }
00157
00158 using boost::bind;
00159 std::cout << partons.size() << " partons and "
00160 << std::count_if(jets.begin(), jets.end(),
00161 bind(std::greater<double>(),
00162 bind(&JetClustering::Jet::pt, _1),
00163 minJetPt)) << " jets." << std::endl;
00164 #endif
00165
00166 Matching<double> matching(partons, jets,
00167 DeltaRSeparation<JetInput::ParticleVector::value_type,
00168 JetClustering::Jet, double>());
00169
00170 typedef Matching<double>::Match Match;
00171 std::vector<Match> matches =
00172 matching.match(
00173 std::less<double>(),
00174 std::bind2nd(std::less<double>(), maxDeltaR));
00175
00176 #ifdef DEBUG
00177 for(std::vector<Match>::const_iterator iter = matches.begin();
00178 iter != matches.end(); ++iter)
00179 std::cout << "\tParton " << iter->index1 << " matches jet "
00180 << iter->index2 << " with a Delta_R of "
00181 << matching.delta(*iter) << std::endl;
00182 #endif
00183
00184 unsigned int unmatchedPartons = 0;
00185 unsigned int unmatchedJets = 0;
00186
00187 matchSummary.clear();
00188 for(std::vector<Match>::const_iterator iter = matches.begin();
00189 iter != matches.end(); ++iter) {
00190 if ((useEt ? jets[iter->index2].et()
00191 : jets[iter->index2].pt()) < minJetPt ||
00192 std::abs(jets[iter->index2].eta()) > maxEta)
00193 unmatchedPartons++;
00194 matchSummary.push_back(
00195 JetPartonMatch(partons[iter->index1]->momentum(),
00196 convert(jets[iter->index2]),
00197 matching.delta(*iter),
00198 partons[iter->index1]->pdg_id()));
00199 }
00200
00201 for(Matching<double>::index_type i = 0; i < partons.size(); i++) {
00202 if (!matching.isMatched1st(i)) {
00203 unmatchedPartons++;
00204 matchSummary.push_back(
00205 JetPartonMatch(partons[i]->momentum(),
00206 partons[i]->pdg_id()));
00207 }
00208 }
00209
00210 for(Matching<double>::index_type i = 0; i < jets.size(); i++) {
00211 if (!matching.isMatched2nd(i)) {
00212 if ((useEt ? jets[i].et()
00213 : jets[i].pt()) >= minJetPt &&
00214 std::abs(jets[i].eta()) <= maxEta)
00215 unmatchedJets++;
00216 matchSummary.push_back(
00217 JetPartonMatch(convert(jets[i])));
00218 }
00219 }
00220
00221 switch(matchMode) {
00222 case kExclusive:
00223 if (!unmatchedJets && !unmatchedPartons)
00224 return 1.0;
00225 break;
00226 case kInclusive:
00227 if (!unmatchedPartons)
00228 return 1.0;
00229 }
00230
00231 return 0.0;
00232 }
00233
00234 }