CMS 3D CMS Logo

JetMatchingHook.cc
Go to the documentation of this file.
4 
6 
7 //#include "HepMC/HEPEVT_Wrapper.h"
8 #include <cassert>
9 
11 
12 extern "C" {
13 // this is patchup for Py6 common block because
14 // several elements of the VINT array are used in the matching process
15 
16 extern struct {
17  int mint[400];
18  double vint[400];
19 } pyint1_;
20 }
21 
22 using namespace gen;
23 using namespace Pythia8;
24 
26  : UserHooks(),
27  fRunBlock(nullptr),
28  fEventBlock(nullptr),
29  fEventNumber(0),
30 // fInfoPtr(info),
31  fJetMatching(nullptr),
32  fJetInputFill(nullptr),
33  fIsInitialized(false) {
34 // assert(fInfoPtr);
35 
37 
38  if (scheme == "Madgraph") {
39  fJetMatching = new JetMatchingMadgraph(ps);
40  fJetInputFill = new Py8toJetInputHEPEVT();
41  } else if (scheme == "MadgraphFastJet") {
42  fJetMatching = new JetMatchingMGFastJet(ps);
43  fJetInputFill = new Py8toJetInput();
44  } else if (scheme == "MLM" || scheme == "Alpgen") {
45  throw cms::Exception("JetMatching") << "Port of " << scheme << "scheme \""
46  << "\""
47  " for parton-shower matching is still in progress."
48  << std::endl;
49  } else
50  throw cms::Exception("InvalidJetMatching") << "Unknown scheme \"" << scheme
51  << "\""
52  " specified for parton-shower matching."
53  << std::endl;
54 }
55 
57  if (fJetMatching)
58  delete fJetMatching;
59 }
60 
63  if (!fRunBlock) {
64  throw cms::Exception("JetMatching") << "Invalid RunInfo" << std::endl;
65  }
67  double etaMax = fJetMatching->getJetEtaMax();
69  return;
70 }
71 
73  setLHEEvent(lhee);
75 
76  // here we'll have to adjust, if needed, for "massless" particles
77  // from earlier Madgraph version(s)
78  // also, we'll have to setup elements of the Py6 fortran array
79  // VINT(357), VINT(358), VINT(360) and VINT(390)
80  // if ( fJetMatching->getMatchingScheme() == "Madgraph" )
81  // {
82  //
83  // }
84 
86 
87  return;
88 }
89 
91 // JetMatchingHook::doVetoPartonLevelEarly( const Event& event )
92 {
93  // event.list();
94 
95  // extract "hardest" event - the output will go into workEvent,
96  // which is a data mamber of base class UserHooks
97  //
98  subEvent(event, true);
99 
100  if (!hepeup_.nup || fJetMatching->isMatchingDone()) {
101  return true;
102  }
103 
104  //
105  // bool jmtch = fJetMatching->match( 0, 0, true ); // true if veto-ed, false if accepted (not veto-ed)
106  std::vector<fastjet::PseudoJet> jetInput =
108  bool jmtch = fJetMatching->match(fEventBlock, &jetInput);
109  if (jmtch) {
110  return true;
111  }
112 
113  // Do not veto events that got this far
114  //
115  return false;
116 }
funct::false
false
Definition: Factorize.h:29
JetMatchingMadgraph.h
vint
double vint[400]
Definition: JetMatchingHook.cc:18
gen::JetMatching::getJetEtaMax
virtual double getJetEtaMax() const =0
mint
int mint[400]
Definition: JetMatchingHook.cc:17
JetMatchingHook::fJetMatching
gen::JetMatching * fJetMatching
Definition: JetMatchingHook.h:80
Py8toJetInput
Definition: Py8toJetInput.h:13
JetMatchingHook::beforeHadronization
virtual void beforeHadronization(lhef::LHEEvent *lhee)
Definition: JetMatchingHook.cc:72
JetMatchingHook.h
JetMatchingHook::fJetInputFill
Py8toJetInput * fJetInputFill
Definition: JetMatchingHook.h:81
Py8toJetInput::setJetEtaMax
void setJetEtaMax(double max)
Definition: Py8toJetInput.h:25
gen::JetMatchingMadgraph
Definition: JetMatchingMadgraph.h:8
info
static const TGPicture * info(bool iBackgroundIsBlack)
Definition: FWCollectionSummaryWidget.cc:153
JetMatchingHook::JetMatchingHook
JetMatchingHook()
Definition: JetMatchingHook.h:53
pyint1_
struct @760 pyint1_
JetMatchingHook::fEventBlock
lhef::LHEEvent * fEventBlock
Definition: JetMatchingHook.h:75
Pythia8
Definition: PythiaDecays.h:21
hepeup_
struct HEPEUP_ hepeup_
JetMatchingHook::~JetMatchingHook
~JetMatchingHook() override
Definition: JetMatchingHook.cc:56
gen::JetMatching::getPartonList
virtual const std::vector< int > * getPartonList()
Definition: JetMatching.h:76
Py8toJetInputHEPEVT
Definition: Py8toJetInput.h:39
gen
Definition: PythiaDecays.h:13
gen::JetMatching::isMatchingDone
bool isMatchingDone()
Definition: JetMatching.h:74
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
gen::JetMatching::match
virtual int match(const lhef::LHEEvent *partonLevel, const std::vector< fastjet::PseudoJet > *jetInput)=0
lhef::LHERunInfo
Definition: LHERunInfo.h:25
edm::ParameterSet
Definition: ParameterSet.h:47
Py8toJetInput::fillJetAlgoInput
virtual const std::vector< fastjet::PseudoJet > fillJetAlgoInput(const Event &, const Event &, const lhef::LHEEvent *lhee=nullptr, const std::vector< int > *partonList=nullptr)
Definition: Py8toJetInput.cc:7
lhef::LHEEvent
Definition: LHEEvent.h:23
JetMatchingHook::setLHEEvent
void setLHEEvent(lhef::LHEEvent *lhee)
Definition: JetMatchingHook.h:63
JetMatchingHook::init
virtual void init(lhef::LHERunInfo *runInfo)
Definition: JetMatchingHook.cc:61
generator_cfi.scheme
scheme
Definition: generator_cfi.py:22
JetMatchingHook::fRunBlock
lhef::LHERunInfo * fRunBlock
Definition: JetMatchingHook.h:74
gen::JetMatching::beforeHadronisationExec
virtual void beforeHadronisationExec()
Definition: JetMatching.cc:25
gen::JetMatchingMGFastJet
Definition: JetMatchingMGFastJet.h:31
Exception
Definition: hltDiff.cc:245
gen::JetMatching::beforeHadronisation
virtual void beforeHadronisation(const lhef::LHEEvent *event)
Definition: JetMatching.cc:23
ALCARECOTkAlBeamHalo_cff.etaMax
etaMax
Definition: ALCARECOTkAlBeamHalo_cff.py:33
HEPEUP_::nup
int nup
Definition: LHECommonBlocks.h:21
submitPVValidationJobs.runInfo
dictionary runInfo
Definition: submitPVValidationJobs.py:1013
JetMatchingHook::setLHERunInfo
void setLHERunInfo(lhef::LHERunInfo *lheri)
Definition: JetMatchingHook.h:55
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Exception.h
JetMatchingHook::doVetoPartonLevel
bool doVetoPartonLevel(const Pythia8::Event &event) override
Definition: JetMatchingHook.cc:90
Py8toJetInput.h
JetMatchingMGFastJet.h
event
Definition: event.py:1
gen::JetMatching::init
virtual void init(const lhef::LHERunInfo *runInfo)
Definition: JetMatching.cc:21