Go to the documentation of this file.00001 #ifndef GeneratorInterface_PartonShowerVeto_JetMatching_h
00002 #define GeneratorInterface_PartonShowerVeto_JetMatching_h
00003
00004 #include <memory>
00005 #include <vector>
00006 #include <string>
00007 #include <set>
00008
00009 #include <boost/shared_ptr.hpp>
00010
00011 #include <HepMC/GenEvent.h>
00012 #include <HepMC/SimpleVector.h>
00013
00014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00015
00016 namespace lhef {
00017
00018 class LHERunInfo;
00019 class LHEEvent;
00020 class JetInput;
00021
00022 }
00023
00024 namespace gen {
00025
00026 class JetMatching {
00027 public:
00028 JetMatching(const edm::ParameterSet ¶ms);
00029 virtual ~JetMatching();
00030
00031 struct JetPartonMatch {
00032 JetPartonMatch(const HepMC::FourVector &parton,
00033 const HepMC::FourVector &jet,
00034 double delta,
00035 int pdgId) :
00036 parton(parton), jet(jet),
00037 delta(delta), pdgId(pdgId) {}
00038
00039 JetPartonMatch(const HepMC::FourVector &parton,
00040 int pdgId) :
00041 parton(parton), delta(-1.0), pdgId(pdgId) {}
00042
00043 JetPartonMatch(const HepMC::FourVector &jet) :
00044 jet(jet), delta(-1.0), pdgId(0) {}
00045
00046 inline bool isMatch() const { return delta >= 0 && pdgId; }
00047 inline bool hasParton() const { return pdgId; }
00048 inline bool hasJet() const { return delta >= 0 || !pdgId; }
00049
00050 HepMC::FourVector parton;
00051 HepMC::FourVector jet;
00052 double delta;
00053 int pdgId;
00054 };
00055
00056 virtual void init(const lhef::LHERunInfo* runInfo);
00057 virtual void beforeHadronisation(const lhef::LHEEvent* event);
00058 virtual void beforeHadronisationExec();
00059
00060 virtual int match(const HepMC::GenEvent *partonLevel,
00061 const HepMC::GenEvent *finalState,
00062 bool showeredFinalState = false) = 0;
00063
00064 virtual std::set<std::string> capabilities() const;
00065
00066 void resetMatchingStatus() { fMatchingStatus = false; }
00067 bool isMatchingDone() { return fMatchingStatus; }
00068
00069 const std::vector<JetPartonMatch> &getMatchSummary() const
00070 { return matchSummary; }
00071
00072 static std::auto_ptr<JetMatching> create(
00073 const edm::ParameterSet ¶ms);
00074
00075 protected:
00076 bool fMatchingStatus;
00077 std::vector<JetPartonMatch> matchSummary;
00078 };
00079
00080 }
00081
00082
00083 #endif // GeneratorCommon_PartonShowerVeto_JetMatching_h