CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetMatchingHook.cc
Go to the documentation of this file.
4 
5 #include "HepMC/HEPEVT_Wrapper.h"
6 
7 extern "C" {
8 // this is patchup for Py6 common block because
9 // several elements of the VINT array are used in the matching process
10 
11  extern struct {
12  int mint[400];
13  double vint[400];
14  } pyint1_;
15 
16 }
17 
18 
19 using namespace gen;
20 using namespace Pythia8;
21 
22 
24  : UserHooks(), fRunBlock(0), fEventBlock(0), fEventNumber(0), fInfoPtr(info), fJetMatching(0)
25 {
26 
27  assert( fInfoPtr );
28 
29  std::string scheme = ps.getParameter<std::string>("scheme");
30 
31  if ( scheme == "Madgraph" )
32  {
34  }
35  else if ( scheme == "MLM" || scheme == "Alpgen" )
36  {
37  throw cms::Exception("JetMatching")
38  << "Port of " << scheme << "scheme \"" << "\""
39  " for parton-shower matching is still in progress."
40  << std::endl;
41  }
42  else
43  throw cms::Exception("InvalidJetMatching")
44  << "Unknown scheme \"" << scheme << "\""
45  " specified for parton-shower matching."
46  << std::endl;
47 
48 }
49 
51 {
52  if ( fJetMatching ) delete fJetMatching;
53 }
54 
56 {
57 
58  setLHERunInfo( runInfo );
59  if ( !fRunBlock )
60  {
61  throw cms::Exception("JetMatching")
62  << "Invalid RunInfo" << std::endl;
63 
64  }
65  fJetMatching->init( runInfo );
66  return;
67 
68 }
69 
71 {
72 
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 }
90 
92 {
93 
94  // extract "hardest" event - the output will go into workEvent,
95  // which is a data mamber of base class UserHooks
96  //
97  subEvent(event,true);
98 
99  setHEPEVT( event ); // here we pass in "full" event, not workEvent !
100 
101  // std::cout << " NPartons= " << hepeup_.nup << std::endl;
102 
103  if ( !hepeup_.nup || fJetMatching->isMatchingDone() )
104  {
105  return true;
106  }
107 
108  // Note from JY:
109  // match(...)input agrs here are reserved for future development and are irrelevat at this point
110  // just for future references, they're:
111  // const HepMC::GenEvent* partonLevel,
112  // const HepMC::GenEvent *finalState,
113  // bool showeredFinalState
114  //
115  bool jmtch = fJetMatching->match( 0, 0, true ); // true if veto-ed, false if accepted (not veto-ed)
116  if ( jmtch )
117  {
118  return true;
119  }
120 
121  // Do not veto events that got this far
122  //
123  return false;
124 
125 }
126 
127 void JetMatchingHook::setHEPEVT( const Event& event )
128 {
129 
130  HepMC::HEPEVT_Wrapper::zero_everything();
131 
132  // service container for further mother-daughters links
133  //
134  std::vector<int> Py8PartonIdx; // position of original (LHE) partons in Py8::Event
135  Py8PartonIdx.clear();
136  std::vector<int> HEPEVTPartonIdx; // position of LHE partons in HEPEVT (incl. ME-generated decays)
137  HEPEVTPartonIdx.clear();
138 
139  // general counter
140  //
141  int index = 0;
142 
143  int Py8PartonCounter = 0;
144  int HEPEVTPartonCounter = 0;
145 
146  // find the fisrt parton that comes from LHE (ME-generated)
147  // skip the incoming particles/partons
148  for ( int iprt=1; iprt<event.size(); iprt++ )
149  {
150  const Particle& part = event[iprt];
151  if ( abs(part.status()) < 22 ) continue; // below 10 is "service"
152  // 11-19 are beam particles; below 10 is "service"
153  // 21 is incoming partons
154  Py8PartonCounter = iprt;
155  break;
156  }
157 
158  const lhef::HEPEUP& hepeup = *fEventBlock->getHEPEUP();
159  // start the counter from 2, because we don't want the incoming particles/oartons !
160  for ( int iprt=2; iprt<hepeup.NUP; iprt++ )
161  {
162  index++;
163  HepMC::HEPEVT_Wrapper::set_id( index, hepeup.IDUP[iprt] );
164  HepMC::HEPEVT_Wrapper::set_status( index, 2 );
165  HepMC::HEPEVT_Wrapper::set_momentum( index, hepeup.PUP[iprt][0], hepeup.PUP[iprt][1], hepeup.PUP[iprt][2], hepeup.PUP[iprt][4] );
166  HepMC::HEPEVT_Wrapper::set_mass( index, hepeup.PUP[iprt][4] );
167  // --> FIXME HepMC::HEPEVT_Wrapper::set_position( index, part.xProd(), part.yProd(), part.zProd(), part.tProd() );
168  HepMC::HEPEVT_Wrapper::set_parents( index, 0, 0 ); // NO, not anymore to the "system particle"
169  HepMC::HEPEVT_Wrapper::set_children( index, 0, 0 );
170  if ( hepeup.MOTHUP[iprt].first > 2 && hepeup.MOTHUP[iprt].second > 2 ) // decay from LHE, will NOT show at the start of Py8 event !!!
171  {
172  HEPEVTPartonCounter++;
173  continue;
174  }
175  Py8PartonIdx.push_back( Py8PartonCounter );
176  Py8PartonCounter++;
177  HEPEVTPartonIdx.push_back( HEPEVTPartonCounter);
178  HEPEVTPartonCounter++;
179  }
180 
181  HepMC::HEPEVT_Wrapper::set_number_entries( index );
182 
183  // now that the initial partons are in, attach parton-level from Pythia8
184  // do NOT reset index as we need to *add* more particles sequentially
185  //
186  for ( int iprt=1; iprt<workEvent.size(); iprt++ ) // from 0-entry (system) or from 1st entry ???
187  {
188 
189  const Particle& part = workEvent[iprt];
190 
191 
192  if ( part.status() != 62 ) continue;
193 
194  index++;
195  HepMC::HEPEVT_Wrapper::set_id( index, part.id() );
196 
197  // HepMC::HEPEVT_Wrapper::set_status( index, event.statusHepMC(iprt) );
198  HepMC::HEPEVT_Wrapper::set_status( index, 1 );
199  HepMC::HEPEVT_Wrapper::set_momentum( index, part.px(), part.py(), part.pz(), part.e() );
200  HepMC::HEPEVT_Wrapper::set_mass( index, part.m() );
201  HepMC::HEPEVT_Wrapper::set_position( index, part.xProd(), part.yProd(), part.zProd(), part.tProd() );
202  HepMC::HEPEVT_Wrapper::set_parents( index, 0, 0 ); // just set to 0 like in Py6...
203  // although for some, mother will need to be re-set properly !
204  HepMC::HEPEVT_Wrapper::set_children( index, 0, 0 );
205 
206  // now refine mother-daughters links, where applicable
207 
208  int parentId = getAncestor( part.daughter1(), event );
209 
210  if ( parentId <= 0 ) continue;
211 
212  for ( int idx=0; idx<(int)Py8PartonIdx.size(); idx++ )
213  {
214  if ( parentId == Py8PartonIdx[idx] )
215  {
216  int idx1 = HEPEVTPartonIdx[idx];
217  HepMC::HEPEVT_Wrapper::set_parents( index, idx1+1, idx1+1 );
218  break;
219  }
220  }
221 
222  }
223 
224  HepMC::HEPEVT_Wrapper::set_number_entries( index );
225 
226  HepMC::HEPEVT_Wrapper::set_event_number( fEventNumber ); // well, if you know it... well, it's one of the counters...
227 
228 // HepMC::HEPEVT_Wrapper::print_hepevt();
229 
230  return;
231 
232 }
233 
234 int JetMatchingHook::getAncestor( int pos, const Event& fullEvent )
235 {
236 
237  int parentId = fullEvent[pos].mother1();
238  int parentPrevId = 0;
239  int counter = pos;
240 
241  while ( parentId > 0 )
242  {
243  if ( parentId == fullEvent[counter].mother2() ) // carbon copy, keep walking up
244  {
245  parentPrevId = parentId;
246  counter = parentId;
247  parentId = fullEvent[parentPrevId].mother1();
248  continue;
249  }
250 
251  // we get here if not a carbon copy
252 
253  // let's check if it's a normal process, etc.
254  //
255  if ( (parentId < parentPrevId) || parentId < fullEvent[counter].mother2() ) // normal process
256  {
257 
258  // first of all, check if hard block
259  if ( abs(fullEvent[counter].status()) == 22 || abs(fullEvent[counter].status()) == 23 )
260  {
261  // yes, it's the hard block
262  // we got what we want, and can exit now !
263  parentId = counter;
264  break;
265  }
266  else
267  {
268  parentPrevId = parentId;
269  parentId = fullEvent[parentPrevId].mother1();
270  }
271  }
272  else if ( parentId > parentPrevId || parentId > pos ) // "circular"/"forward-pointing parent" - intermediate process
273  {
274  parentId = -1;
275  break;
276  }
277 
278  // additional checks... although we shouldn't be geeting here all that much...
279  //
280  if ( abs(fullEvent[parentId].status()) == 22 || abs(fullEvent[parentId].status())== 23 ) // hard block
281  {
282  break;
283  }
284  if ( abs(fullEvent[parentId].status()) < 22 ) // incoming
285  {
286  parentId = -1;
287  break;
288  }
289  }
290 
291  return parentId;
292 
293 }
T getParameter(std::string const &) const
lhef::LHEEvent * fEventBlock
Pythia8::Info * fInfoPtr
struct HEPEUP_ hepeup_
virtual ~JetMatchingHook()
void setLHEEvent(lhef::LHEEvent *lhee)
#define abs(x)
Definition: mlp_lapack.h:159
int mint[400]
virtual void beforeHadronisation(const lhef::LHEEvent *event)
Definition: JetMatching.cc:31
const HEPEUP * getHEPEUP() const
Definition: LHEEvent.h:43
void beforeHadronization(lhef::LHEEvent *lhee)
double vint[400]
std::vector< std::pair< int, int > > MOTHUP
Definition: LesHouches.h:236
JetMatchingHook(const edm::ParameterSet &, Pythia8::Info *)
std::vector< FiveVector > PUP
Definition: LesHouches.h:248
virtual bool doVetoPartonLevel(const Pythia8::Event &event)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
std::vector< int > IDUP
Definition: LesHouches.h:225
void setLHERunInfo(lhef::LHERunInfo *lheri)
void init(lhef::LHERunInfo *runInfo)
part
Definition: HCALResponse.h:21
virtual void beforeHadronisationExec()
Definition: JetMatching.cc:35
bool isMatchingDone()
Definition: JetMatching.h:67
void setHEPEVT(const Pythia8::Event &)
gen::JetMatching * fJetMatching
virtual int match(const HepMC::GenEvent *partonLevel, const HepMC::GenEvent *finalState, bool showeredFinalState=false)=0
int getAncestor(int, const Pythia8::Event &)
struct @324 pyint1_
tuple status
Definition: ntuplemaker.py:245
virtual void init(const lhef::LHERunInfo *runInfo)
Definition: JetMatching.cc:27
lhef::LHERunInfo * fRunBlock