CMS 3D CMS Logo

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