CMS 3D CMS Logo

lhef::JetMatchingMLM Class Reference

#include <GeneratorInterface/LHEInterface/interface/JetMatchingMLM.h>

Inheritance diagram for lhef::JetMatchingMLM:

lhef::JetMatching

List of all members.

Public Member Functions

 JetMatchingMLM (const edm::ParameterSet &params)
 ~JetMatchingMLM ()

Private Types

enum  MatchMode { kExclusive = 0, kInclusive }

Private Member Functions

std::set< std::string > capabilities () const
double match (const HepMC::GenEvent *partonLevel, const HepMC::GenEvent *finalState, bool showeredFinalState)

Private Attributes

std::auto_ptr< JetClusteringjetClustering
std::auto_ptr< JetInputjetInput
MatchMode matchMode
double matchPtFraction
const double maxDeltaR
double maxEta
const double minJetPt
std::auto_ptr< JetInputpartonInput
bool useEt


Detailed Description

Definition at line 19 of file JetMatchingMLM.h.


Member Enumeration Documentation

enum lhef::JetMatchingMLM::MatchMode [private]

Enumerator:
kExclusive 
kInclusive 

Definition at line 31 of file JetMatchingMLM.h.

00031                        {
00032                 kExclusive = 0,
00033                 kInclusive
00034         };


Constructor & Destructor Documentation

lhef::JetMatchingMLM::JetMatchingMLM ( const edm::ParameterSet params  ) 

Definition at line 65 of file JetMatchingMLM.cc.

References lat::endl(), Exception, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), jetClustering, kExclusive, kInclusive, matchMode, matchPtFraction, minJetPt, and partonInput.

00065                                                             :
00066         JetMatching(params),
00067         maxDeltaR(params.getParameter<double>("maxDeltaR")),
00068         minJetPt(params.getParameter<double>("jetPtMin")),
00069         maxEta(params.getParameter<double>("maxEta")),
00070         matchPtFraction(0.75),
00071         useEt(params.getParameter<bool>("useEt")),
00072         partonInput(new JetInput(params)),
00073         jetInput(new JetInput(*partonInput))
00074 {
00075         partonInput->setPartonicFinalState(false);
00076         partonInput->setHardProcessOnly(false);
00077 
00078         if (params.exists("matchPtFraction"))
00079                 matchPtFraction =
00080                         params.getParameter<double>("matchPtFraction");
00081 
00082         jetClustering.reset(
00083                 new JetClustering(params, minJetPt * matchPtFraction));
00084 
00085         std::string matchMode = params.getParameter<std::string>("matchMode");
00086         if (matchMode == "inclusive")
00087                 this->matchMode = kInclusive;
00088         else if (matchMode == "exclusive")
00089                 this->matchMode = kExclusive;
00090         else
00091                 throw cms::Exception("Configuration")
00092                         << "Invalid matching mode '" << matchMode
00093                         << "' specified." << std::endl;
00094 }

lhef::JetMatchingMLM::~JetMatchingMLM (  ) 

Definition at line 96 of file JetMatchingMLM.cc.

00097 {
00098 }


Member Function Documentation

std::set< std::string > lhef::JetMatchingMLM::capabilities (  )  const [private, virtual]

Reimplemented from lhef::JetMatching.

Definition at line 100 of file JetMatchingMLM.cc.

References lhef::JetMatching::capabilities(), and HLT_VtxMuL3::result.

00101 {
00102         std::set<std::string> result = JetMatching::capabilities();
00103         result.insert("matchSummary");
00104         return result;
00105 }

double lhef::JetMatchingMLM::match ( const HepMC::GenEvent *  partonLevel,
const HepMC::GenEvent *  finalState,
bool  showeredFinalState 
) [private, virtual]

Implements lhef::JetMatching.

Definition at line 110 of file JetMatchingMLM.cc.

References funct::abs(), c, lhef::convert(), GenMuonPlsPt100GeV_cfg::cout, lat::endl(), eta, i, iter, jetInput, pfTauBenchmarkGeneric_cfi::jets, kExclusive, kInclusive, matchMode, lhef::JetMatching::matchSummary, maxDeltaR, maxEta, minJetPt, lhef::JetClustering::Jet::pt(), python::multivaluedict::sort(), and useEt.

00113 {
00114         JetInput::ParticleVector partons = (*partonInput)(partonLevel);
00115         std::sort(partons.begin(), partons.end(), ParticlePtGreater());
00116 
00117         JetInput::ParticleVector jetInput =
00118                         showeredFinalState ? (*partonInput)(finalState)
00119                                            : (*this->jetInput)(finalState);
00120         std::sort(jetInput.begin(), jetInput.end(), ParticlePtGreater());
00121 
00122         std::vector<JetClustering::Jet> jets = (*jetClustering)(jetInput);
00123 
00124 #ifdef DEBUG
00125         std::cout << "===== Partons:" << std::endl;
00126         for(JetClustering::ParticleVector::const_iterator c = partons.begin();
00127             c != partons.end(); c++)
00128                 std::cout << "\tpid = " << (*c)->pdg_id()
00129                           << ", pt = " << (*c)->momentum().perp()
00130                           << ", eta = " << (*c)->momentum().eta()
00131                           << ", phi = " << (*c)->momentum().phi()
00132                           << std::endl;
00133         std::cout << "===== JetInput:" << std::endl;
00134         for(JetClustering::ParticleVector::const_iterator c = jetInput.begin();
00135             c != jetInput.end(); c++)
00136                 std::cout << "\tpid = " << (*c)->pdg_id()
00137                           << ", pt = " << (*c)->momentum().perp()
00138                           << ", eta = " << (*c)->momentum().eta()
00139                           << ", phi = " << (*c)->momentum().phi()
00140                           << std::endl;
00141         std::cout << "----- " << jets.size() << " jets:" << std::endl;
00142         for(std::vector<JetClustering::Jet>::const_iterator iter = jets.begin();
00143             iter != jets.end(); ++iter) {
00144                 std::cout << "* pt = " << iter->pt()
00145                           << ", eta = " << iter->eta()
00146                           << ", phi = " << iter->phi()
00147                           << std::endl;
00148                 for(JetClustering::ParticleVector::const_iterator c = iter->constituents().begin();
00149                     c != iter->constituents().end(); c++)
00150                         std::cout << "\tpid = " << (*c)->pdg_id()
00151                                   << ", pt = " << (*c)->momentum().perp()
00152                                   << ", eta = " << (*c)->momentum().eta()
00153                                   << ", phi = " << (*c)->momentum().phi()
00154                                   << std::endl;
00155         }
00156 
00157         using boost::bind;
00158         std::cout << partons.size() << " partons and "
00159                   << std::count_if(jets.begin(), jets.end(),
00160                         bind(std::greater<double>(),
00161                              bind(&JetClustering::Jet::pt, _1),
00162                              minJetPt)) << " jets." << std::endl;
00163 #endif
00164 
00165         Matching<double> matching(partons, jets,
00166                         DeltaRSeparation<JetInput::ParticleVector::value_type,
00167                                          JetClustering::Jet, double>());
00168 
00169         typedef Matching<double>::Match Match;
00170         std::vector<Match> matches =
00171                         matching.match(
00172                                 std::less<double>(),
00173                                 std::bind2nd(std::less<double>(), maxDeltaR));
00174 
00175 #ifdef DEBUG
00176         for(std::vector<Match>::const_iterator iter = matches.begin();
00177             iter != matches.end(); ++iter)
00178                 std::cout << "\tParton " << iter->index1 << " matches jet "
00179                           << iter->index2 << " with a Delta_R of "
00180                           << matching.delta(*iter) << std::endl;
00181 #endif
00182 
00183         unsigned int unmatchedPartons = 0;
00184         unsigned int unmatchedJets = 0;
00185 
00186         matchSummary.clear();
00187         for(std::vector<Match>::const_iterator iter = matches.begin();
00188             iter != matches.end(); ++iter) {
00189                 if ((useEt ? jets[iter->index2].et()
00190                            : jets[iter->index2].pt()) < minJetPt ||
00191                     std::abs(jets[iter->index2].eta()) > maxEta)
00192                         unmatchedPartons++;
00193                 matchSummary.push_back(
00194                         JetPartonMatch(partons[iter->index1]->momentum(),
00195                                        convert(jets[iter->index2]),
00196                                        matching.delta(*iter),
00197                                        partons[iter->index1]->pdg_id()));
00198         }
00199 
00200         for(Matching<double>::index_type i = 0; i < partons.size(); i++) {
00201                 if (!matching.isMatched1st(i)) {
00202                         unmatchedPartons++;
00203                         matchSummary.push_back(
00204                                 JetPartonMatch(partons[i]->momentum(),
00205                                                partons[i]->pdg_id()));
00206                 }
00207         }
00208 
00209         for(Matching<double>::index_type i = 0; i < jets.size(); i++) {
00210                 if (!matching.isMatched2nd(i)) {
00211                         if ((useEt ? jets[i].et()
00212                                    : jets[i].pt()) >= minJetPt &&
00213                             std::abs(jets[i].eta()) <= maxEta)
00214                                 unmatchedJets++;
00215                         matchSummary.push_back(
00216                                 JetPartonMatch(convert(jets[i])));
00217                 }
00218         }
00219 
00220         switch(matchMode) {
00221             case kExclusive:
00222                 if (!unmatchedJets && !unmatchedPartons)
00223                         return 1.0;
00224                 break;
00225             case kInclusive:
00226                 if (!unmatchedPartons)
00227                         return 1.0;
00228         }
00229 
00230         return 0.0;
00231 }


Member Data Documentation

std::auto_ptr<JetClustering> lhef::JetMatchingMLM::jetClustering [private]

Definition at line 45 of file JetMatchingMLM.h.

Referenced by JetMatchingMLM().

std::auto_ptr<JetInput> lhef::JetMatchingMLM::jetInput [private]

Definition at line 44 of file JetMatchingMLM.h.

Referenced by match().

MatchMode lhef::JetMatchingMLM::matchMode [private]

Definition at line 41 of file JetMatchingMLM.h.

Referenced by JetMatchingMLM(), and match().

double lhef::JetMatchingMLM::matchPtFraction [private]

Definition at line 39 of file JetMatchingMLM.h.

Referenced by JetMatchingMLM().

const double lhef::JetMatchingMLM::maxDeltaR [private]

Definition at line 36 of file JetMatchingMLM.h.

Referenced by match().

double lhef::JetMatchingMLM::maxEta [private]

Definition at line 38 of file JetMatchingMLM.h.

Referenced by match().

const double lhef::JetMatchingMLM::minJetPt [private]

Definition at line 37 of file JetMatchingMLM.h.

Referenced by JetMatchingMLM(), and match().

std::auto_ptr<JetInput> lhef::JetMatchingMLM::partonInput [private]

Definition at line 43 of file JetMatchingMLM.h.

Referenced by JetMatchingMLM().

bool lhef::JetMatchingMLM::useEt [private]

Definition at line 40 of file JetMatchingMLM.h.

Referenced by match().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:48:49 2009 for CMSSW by  doxygen 1.5.4