Main Page
Namespaces
Classes
Package Documentation
GeneratorInterface
PartonShowerVeto
interface
JetMatching.h
Go to the documentation of this file.
1
#ifndef GeneratorInterface_PartonShowerVeto_JetMatching_h
2
#define GeneratorInterface_PartonShowerVeto_JetMatching_h
3
4
#include <memory>
5
#include <vector>
6
#include <string>
7
#include <set>
8
9
#include <boost/shared_ptr.hpp>
10
11
// #include <HepMC/GenEvent.h>
12
// #include <HepMC/SimpleVector.h>
13
14
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
15
16
#include "fastjet/ClusterSequence.hh"
// gives both PseudoJet & JetDefinition
17
// #include "fastjet/Selector.hh"
18
19
namespace
lhef
{
20
21
class
LHERunInfo
;
22
class
LHEEvent;
23
class
JetInput;
24
// class JetClustering;
25
}
26
27
namespace
gen
{
28
29
class
JetMatching
{
30
public
:
31
JetMatching
(
const
edm::ParameterSet
¶ms);
32
virtual
~
JetMatching
();
33
34
/*
35
struct JetPartonMatch {
36
JetPartonMatch(const HepMC::FourVector &parton,
37
const HepMC::FourVector &jet,
38
double delta,
39
int pdgId) :
40
parton(parton), jet(jet),
41
delta(delta), pdgId(pdgId) {}
42
43
JetPartonMatch(const HepMC::FourVector &parton,
44
int pdgId) :
45
parton(parton), delta(-1.0), pdgId(pdgId) {}
46
47
JetPartonMatch(const HepMC::FourVector &jet) :
48
jet(jet), delta(-1.0), pdgId(0) {}
49
50
inline bool isMatch() const { return delta >= 0 && pdgId; }
51
inline bool hasParton() const { return pdgId; }
52
inline bool hasJet() const { return delta >= 0 || !pdgId; }
53
54
HepMC::FourVector parton;
55
HepMC::FourVector jet;
56
double delta;
57
int pdgId;
58
};
59
*/
60
virtual
void
init
(
const
lhef::LHERunInfo
* runInfo);
61
virtual
bool
initAfterBeams
() {
return
true
; }
62
virtual
void
beforeHadronisation(
const
lhef::LHEEvent
*
event
);
63
virtual
void
beforeHadronisationExec();
64
65
// void setJetInput( const std::vector<fastjet::PseudoJet> input ) { fJetInput=input; return; }
66
67
virtual
int
match
(
const
lhef::LHEEvent
* partonLevel,
const
std::vector<fastjet::PseudoJet>* jetInput ) = 0;
68
/*
69
virtual int match(const HepMC::GenEvent *partonLevel,
70
const HepMC::GenEvent *finalState,
71
bool showeredFinalState = false) = 0;
72
*/
73
virtual
std::set<std::string> capabilities()
const
;
74
75
void
resetMatchingStatus
() { fMatchingStatus =
false
; }
76
bool
isMatchingDone
() {
return
fMatchingStatus; }
77
78
virtual
const
std::vector<int>*
getPartonList
() {
return
nullptr
; }
79
virtual
double
getJetEtaMax()
const
= 0;
80
81
/*
82
const std::vector<JetPartonMatch> &getMatchSummary() const
83
{ return matchSummary; }
84
*/
85
static
std::unique_ptr<JetMatching>
create
(
86
const
edm::ParameterSet
¶ms);
87
88
protected
:
89
bool
fMatchingStatus
;
90
/* std::vector<JetPartonMatch> matchSummary; */
91
// std::vector<fastjet::PseudoJet> fJetInput;
92
93
};
94
95
}
// namespace gen
96
97
98
#endif // GeneratorCommon_PartonShowerVeto_JetMatching_h
lhef::LHEEvent
Definition:
LHEEvent.h:23
beamerCreator.create
def create(alignables, pedeDump, additionalData, outputFile, config)
Definition:
beamerCreator.py:44
init
int init
Definition:
HydjetWrapper.h:67
particlelevel_cff.LHERunInfo
LHERunInfo
Definition:
particlelevel_cff.py:56
lhef
Definition:
ExhumeHadronizer.h:14
ParameterSet.h
gen::JetMatching::fMatchingStatus
bool fMatchingStatus
Definition:
JetMatching.h:89
gen::JetMatching::resetMatchingStatus
void resetMatchingStatus()
Definition:
JetMatching.h:75
gen
Definition:
PythiaDecays.h:13
lhef::LHERunInfo
Definition:
LHERunInfo.h:25
gen::JetMatching::getPartonList
virtual const std::vector< int > * getPartonList()
Definition:
JetMatching.h:78
gen::JetMatching::initAfterBeams
virtual bool initAfterBeams()
Definition:
JetMatching.h:61
gen::JetMatching::isMatchingDone
bool isMatchingDone()
Definition:
JetMatching.h:76
edm::ParameterSet
Definition:
ParameterSet.h:36
match
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition:
Utils.h:10
event
Definition:
event.py:1
gen::JetMatching
Definition:
JetMatching.h:29
Generated for CMSSW Reference Manual by
1.8.11