Go to the documentation of this file.00001 #include "GeneratorInterface/Pythia8Interface/interface/JetMatchingHook.h"
00002 #include "GeneratorInterface/PartonShowerVeto/interface/JetMatchingMadgraph.h"
00003 #include "FWCore/Utilities/interface/Exception.h"
00004
00005 #include "HepMC/HEPEVT_Wrapper.h"
00006
00007 extern "C" {
00008
00009
00010
00011 extern struct {
00012 int mint[400];
00013 double vint[400];
00014 } pyint1_;
00015
00016 }
00017
00018
00019 using namespace gen;
00020 using namespace Pythia8;
00021
00022
00023 JetMatchingHook::JetMatchingHook( const edm::ParameterSet& ps, Info* info )
00024 : UserHooks(), fRunBlock(0), fEventBlock(0), fEventNumber(0), fInfoPtr(info), fJetMatching(0)
00025 {
00026
00027 assert( fInfoPtr );
00028
00029 std::string scheme = ps.getParameter<std::string>("scheme");
00030
00031 if ( scheme == "Madgraph" )
00032 {
00033 fJetMatching = new JetMatchingMadgraph(ps);
00034 }
00035 else if ( scheme == "MLM" || scheme == "Alpgen" )
00036 {
00037 throw cms::Exception("JetMatching")
00038 << "Port of " << scheme << "scheme \"" << "\""
00039 " for parton-shower matching is still in progress."
00040 << std::endl;
00041 }
00042 else
00043 throw cms::Exception("InvalidJetMatching")
00044 << "Unknown scheme \"" << scheme << "\""
00045 " specified for parton-shower matching."
00046 << std::endl;
00047
00048 }
00049
00050 JetMatchingHook::~JetMatchingHook()
00051 {
00052 if ( fJetMatching ) delete fJetMatching;
00053 }
00054
00055 void JetMatchingHook::init ( lhef::LHERunInfo* runInfo )
00056 {
00057
00058 setLHERunInfo( runInfo );
00059 if ( !fRunBlock )
00060 {
00061 throw cms::Exception("JetMatching")
00062 << "Invalid RunInfo" << std::endl;
00063
00064 }
00065 fJetMatching->init( runInfo );
00066 return;
00067
00068 }
00069
00070 void JetMatchingHook::beforeHadronization( lhef::LHEEvent* lhee )
00071 {
00072
00073 setLHEEvent( lhee );
00074 fJetMatching->beforeHadronisation( lhee );
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085 fJetMatching->beforeHadronisationExec();
00086
00087 return;
00088
00089 }
00090
00091 bool JetMatchingHook::doVetoPartonLevel( const Event& event )
00092 {
00093
00094
00095
00096
00097 subEvent(event,true);
00098
00099 setHEPEVT( event );
00100
00101
00102
00103 if ( !hepeup_.nup || fJetMatching->isMatchingDone() )
00104 {
00105 return true;
00106 }
00107
00108
00109
00110
00111
00112
00113
00114
00115 bool jmtch = fJetMatching->match( 0, 0, true );
00116 if ( jmtch )
00117 {
00118 return true;
00119 }
00120
00121
00122
00123 return false;
00124
00125 }
00126
00127 void JetMatchingHook::setHEPEVT( const Event& event )
00128 {
00129
00130 HepMC::HEPEVT_Wrapper::zero_everything();
00131
00132
00133
00134 std::vector<int> Py8PartonIdx;
00135 Py8PartonIdx.clear();
00136 std::vector<int> HEPEVTPartonIdx;
00137 HEPEVTPartonIdx.clear();
00138
00139
00140
00141 int index = 0;
00142
00143 int Py8PartonCounter = 0;
00144 int HEPEVTPartonCounter = 0;
00145
00146
00147
00148 for ( int iprt=1; iprt<event.size(); iprt++ )
00149 {
00150 const Particle& part = event[iprt];
00151 if ( abs(part.status()) < 22 ) continue;
00152
00153
00154 Py8PartonCounter = iprt;
00155 break;
00156 }
00157
00158 const lhef::HEPEUP& hepeup = *fEventBlock->getHEPEUP();
00159
00160 for ( int iprt=2; iprt<hepeup.NUP; iprt++ )
00161 {
00162 index++;
00163 HepMC::HEPEVT_Wrapper::set_id( index, hepeup.IDUP[iprt] );
00164 HepMC::HEPEVT_Wrapper::set_status( index, 2 );
00165 HepMC::HEPEVT_Wrapper::set_momentum( index, hepeup.PUP[iprt][0], hepeup.PUP[iprt][1], hepeup.PUP[iprt][2], hepeup.PUP[iprt][4] );
00166 HepMC::HEPEVT_Wrapper::set_mass( index, hepeup.PUP[iprt][4] );
00167
00168 HepMC::HEPEVT_Wrapper::set_parents( index, 0, 0 );
00169 HepMC::HEPEVT_Wrapper::set_children( index, 0, 0 );
00170 if ( hepeup.MOTHUP[iprt].first > 2 && hepeup.MOTHUP[iprt].second > 2 )
00171 {
00172 HEPEVTPartonCounter++;
00173 continue;
00174 }
00175 Py8PartonIdx.push_back( Py8PartonCounter );
00176 Py8PartonCounter++;
00177 HEPEVTPartonIdx.push_back( HEPEVTPartonCounter);
00178 HEPEVTPartonCounter++;
00179 }
00180
00181 HepMC::HEPEVT_Wrapper::set_number_entries( index );
00182
00183
00184
00185
00186 for ( int iprt=1; iprt<workEvent.size(); iprt++ )
00187 {
00188
00189 const Particle& part = workEvent[iprt];
00190
00191
00192 if ( part.status() != 62 ) continue;
00193
00194 index++;
00195 HepMC::HEPEVT_Wrapper::set_id( index, part.id() );
00196
00197
00198 HepMC::HEPEVT_Wrapper::set_status( index, 1 );
00199 HepMC::HEPEVT_Wrapper::set_momentum( index, part.px(), part.py(), part.pz(), part.e() );
00200 HepMC::HEPEVT_Wrapper::set_mass( index, part.m() );
00201 HepMC::HEPEVT_Wrapper::set_position( index, part.xProd(), part.yProd(), part.zProd(), part.tProd() );
00202 HepMC::HEPEVT_Wrapper::set_parents( index, 0, 0 );
00203
00204 HepMC::HEPEVT_Wrapper::set_children( index, 0, 0 );
00205
00206
00207
00208 int parentId = getAncestor( part.daughter1(), event );
00209
00210 if ( parentId <= 0 ) continue;
00211
00212 for ( int idx=0; idx<(int)Py8PartonIdx.size(); idx++ )
00213 {
00214 if ( parentId == Py8PartonIdx[idx] )
00215 {
00216 int idx1 = HEPEVTPartonIdx[idx];
00217 HepMC::HEPEVT_Wrapper::set_parents( index, idx1+1, idx1+1 );
00218 break;
00219 }
00220 }
00221
00222 }
00223
00224 HepMC::HEPEVT_Wrapper::set_number_entries( index );
00225
00226 HepMC::HEPEVT_Wrapper::set_event_number( fEventNumber );
00227
00228
00229
00230 return;
00231
00232 }
00233
00234 int JetMatchingHook::getAncestor( int pos, const Event& fullEvent )
00235 {
00236
00237 int parentId = fullEvent[pos].mother1();
00238 int parentPrevId = 0;
00239 int counter = pos;
00240
00241 while ( parentId > 0 )
00242 {
00243 if ( parentId == fullEvent[counter].mother2() )
00244 {
00245 parentPrevId = parentId;
00246 counter = parentId;
00247 parentId = fullEvent[parentPrevId].mother1();
00248 continue;
00249 }
00250
00251
00252
00253
00254
00255 if ( (parentId < parentPrevId) || parentId < fullEvent[counter].mother2() )
00256 {
00257
00258
00259 if ( abs(fullEvent[counter].status()) == 22 || abs(fullEvent[counter].status()) == 23 )
00260 {
00261
00262
00263 parentId = counter;
00264 break;
00265 }
00266 else
00267 {
00268 parentPrevId = parentId;
00269 parentId = fullEvent[parentPrevId].mother1();
00270 }
00271 }
00272 else if ( parentId > parentPrevId || parentId > pos )
00273 {
00274 parentId = -1;
00275 break;
00276 }
00277
00278
00279
00280 if ( abs(fullEvent[parentId].status()) == 22 || abs(fullEvent[parentId].status())== 23 )
00281 {
00282 break;
00283 }
00284 if ( abs(fullEvent[parentId].status()) < 22 )
00285 {
00286 parentId = -1;
00287 break;
00288 }
00289 }
00290
00291 return parentId;
00292
00293 }