CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/GeneratorInterface/Pythia8Interface/src/JetMatchingHook.cc

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 #include <cassert>
00007 
00008 extern "C" {
00009 // this is patchup for Py6 common block because 
00010 // several elements of the VINT array are used in the matching process
00011 
00012     extern struct {
00013         int mint[400];
00014         double vint[400];
00015     } pyint1_;
00016 
00017 }
00018 
00019 
00020 using namespace gen;
00021 using namespace Pythia8;
00022 
00023 
00024 JetMatchingHook::JetMatchingHook( const edm::ParameterSet& ps, Info* info )
00025    : UserHooks(), fRunBlock(0), fEventBlock(0), fEventNumber(0), fInfoPtr(info), fJetMatching(0)
00026 {
00027 
00028    assert( fInfoPtr );
00029    
00030    std::string scheme = ps.getParameter<std::string>("scheme");
00031 
00032    if ( scheme == "Madgraph" )
00033    {
00034       fJetMatching = new JetMatchingMadgraph(ps);
00035    }
00036    else if ( scheme == "MLM" || scheme == "Alpgen" )
00037    {
00038         throw cms::Exception("JetMatching")
00039         << "Port of " << scheme << "scheme \"" << "\""
00040         " for parton-shower matching is still in progress."
00041         << std::endl;
00042    }
00043    else
00044       throw cms::Exception("InvalidJetMatching")
00045       << "Unknown scheme \"" << scheme << "\""
00046       " specified for parton-shower matching."
00047       << std::endl;
00048  
00049 }
00050 
00051 JetMatchingHook::~JetMatchingHook()
00052 {
00053    if ( fJetMatching ) delete fJetMatching;
00054 }
00055 
00056 void JetMatchingHook::init ( lhef::LHERunInfo* runInfo )
00057 {
00058 
00059    setLHERunInfo( runInfo );
00060    if ( !fRunBlock )
00061    {
00062       throw cms::Exception("JetMatching")
00063       << "Invalid RunInfo" << std::endl;
00064    
00065    }
00066    fJetMatching->init( runInfo );
00067    return; 
00068 
00069 }
00070 
00071 void JetMatchingHook::beforeHadronization( lhef::LHEEvent* lhee )
00072 {
00073 
00074    setLHEEvent( lhee );   
00075    fJetMatching->beforeHadronisation( lhee );
00076       
00077    // here we'll have to adjust, if needed, for "massless" particles
00078    // from earlier Madgraph version(s)
00079    // also, we'll have to setup elements of the Py6 fortran array 
00080    // VINT(357), VINT(358), VINT(360) and VINT(390)
00081    // if ( fJetMatching->getMatchingScheme() == "Madgraph" )
00082    // {
00083    //    
00084    // }
00085    
00086    fJetMatching->beforeHadronisationExec();  
00087    
00088    return;   
00089 
00090 }
00091 
00092 bool JetMatchingHook::doVetoPartonLevel( const Event& event )
00093 {
00094                   
00095    // extract "hardest" event - the output will go into workEvent, 
00096    // which is a data mamber of base class UserHooks
00097    //
00098    subEvent(event,true);
00099 
00100    setHEPEVT( event ); // here we pass in "full" event, not workEvent !
00101 
00102    // std::cout << " NPartons= " << hepeup_.nup << std::endl;
00103    
00104    if ( !hepeup_.nup || fJetMatching->isMatchingDone() )
00105    {
00106       return true;
00107    }
00108    
00109    // Note from JY:
00110    // match(...)input agrs here are reserved for future development and are irrelevat at this point
00111    // just for future references, they're: 
00112    // const HepMC::GenEvent* partonLevel,
00113    // const HepMC::GenEvent *finalState,
00114    // bool showeredFinalState
00115    //
00116    bool jmtch = fJetMatching->match( 0, 0, true ); // true if veto-ed, false if accepted (not veto-ed)
00117    if ( jmtch )
00118    {
00119       return true;
00120    }
00121          
00122    // Do not veto events that got this far
00123    //
00124    return false;
00125 
00126 }
00127 
00128 void JetMatchingHook::setHEPEVT( const Event& event )
00129 {
00130            
00131    HepMC::HEPEVT_Wrapper::zero_everything();   
00132       
00133    // service container for further mother-daughters links
00134    //
00135    std::vector<int> Py8PartonIdx; // position of original (LHE) partons in Py8::Event
00136    Py8PartonIdx.clear(); 
00137    std::vector<int> HEPEVTPartonIdx; // position of LHE partons in HEPEVT (incl. ME-generated decays)
00138    HEPEVTPartonIdx.clear(); 
00139 
00140    // general counter
00141    //
00142    int index = 0;
00143 
00144    int Py8PartonCounter = 0;
00145    int HEPEVTPartonCounter = 0;
00146    
00147    // find the fisrt parton that comes from LHE (ME-generated)
00148    // skip the incoming particles/partons
00149    for ( int iprt=1; iprt<event.size(); iprt++ )
00150    {
00151       const Particle& part = event[iprt];
00152       if ( abs(part.status()) < 22 ) continue; // below 10 is "service"
00153                                                // 11-19 are beam particles; below 10 is "service"
00154                                                // 21 is incoming partons      
00155       Py8PartonCounter = iprt;
00156       break;
00157    }
00158 
00159    const lhef::HEPEUP& hepeup = *fEventBlock->getHEPEUP();
00160    // start the counter from 2, because we don't want the incoming particles/oartons !
00161    for ( int iprt=2; iprt<hepeup.NUP; iprt++ )
00162    {
00163       index++;
00164       HepMC::HEPEVT_Wrapper::set_id( index, hepeup.IDUP[iprt] );
00165       HepMC::HEPEVT_Wrapper::set_status( index, 2 );
00166       HepMC::HEPEVT_Wrapper::set_momentum( index, hepeup.PUP[iprt][0], hepeup.PUP[iprt][1], hepeup.PUP[iprt][2], hepeup.PUP[iprt][4] );
00167       HepMC::HEPEVT_Wrapper::set_mass( index, hepeup.PUP[iprt][4] );
00168       // --> FIXME HepMC::HEPEVT_Wrapper::set_position( index, part.xProd(), part.yProd(), part.zProd(), part.tProd() );
00169       HepMC::HEPEVT_Wrapper::set_parents( index, 0, 0 ); // NO, not anymore to the "system particle"
00170       HepMC::HEPEVT_Wrapper::set_children( index, 0, 0 ); 
00171       if (  hepeup.MOTHUP[iprt].first > 2 && hepeup.MOTHUP[iprt].second > 2 ) // decay from LHE, will NOT show at the start of Py8 event !!!
00172       {
00173          HEPEVTPartonCounter++;
00174          continue;
00175       }
00176       Py8PartonIdx.push_back( Py8PartonCounter );
00177       Py8PartonCounter++;
00178       HEPEVTPartonIdx.push_back( HEPEVTPartonCounter);
00179       HEPEVTPartonCounter++;   
00180    }
00181       
00182    HepMC::HEPEVT_Wrapper::set_number_entries( index );   
00183          
00184    // now that the initial partons are in, attach parton-level from Pythia8
00185    // do NOT reset index as we need to *add* more particles sequentially
00186    //
00187    for ( int iprt=1; iprt<workEvent.size(); iprt++ ) // from 0-entry (system) or from 1st entry ???
00188    {
00189    
00190       const Particle& part = workEvent[iprt];
00191       
00192 
00193       if ( part.status() != 62 ) continue;
00194                   
00195       index++;
00196       HepMC::HEPEVT_Wrapper::set_id( index, part.id() );
00197       
00198       // HepMC::HEPEVT_Wrapper::set_status( index, event.statusHepMC(iprt) ); 
00199       HepMC::HEPEVT_Wrapper::set_status( index, 1 );      
00200       HepMC::HEPEVT_Wrapper::set_momentum( index, part.px(), part.py(), part.pz(), part.e() );
00201       HepMC::HEPEVT_Wrapper::set_mass( index, part.m() );
00202       HepMC::HEPEVT_Wrapper::set_position( index, part.xProd(), part.yProd(), part.zProd(), part.tProd() );
00203       HepMC::HEPEVT_Wrapper::set_parents( index, 0, 0 ); // just set to 0 like in Py6...
00204                                                          // although for some, mother will need to be re-set properly !
00205       HepMC::HEPEVT_Wrapper::set_children( index, 0, 0 );
00206 
00207       // now refine mother-daughters links, where applicable
00208       
00209       int parentId = getAncestor( part.daughter1(), event );
00210       
00211       if ( parentId <= 0 ) continue;
00212 
00213       for ( int idx=0; idx<(int)Py8PartonIdx.size(); idx++ )
00214       {
00215          if ( parentId == Py8PartonIdx[idx] )
00216          {
00217             int idx1 = HEPEVTPartonIdx[idx]; 
00218             HepMC::HEPEVT_Wrapper::set_parents( index, idx1+1, idx1+1 ); 
00219             break;
00220          }
00221       }
00222 
00223    }  
00224      
00225    HepMC::HEPEVT_Wrapper::set_number_entries( index );
00226    
00227    HepMC::HEPEVT_Wrapper::set_event_number( fEventNumber ); // well, if you know it... well, it's one of the counters...
00228    
00229 //   HepMC::HEPEVT_Wrapper::print_hepevt();
00230    
00231    return;
00232 
00233 }
00234 
00235 int JetMatchingHook::getAncestor( int pos, const Event& fullEvent )
00236 {
00237 
00238    int parentId = fullEvent[pos].mother1();
00239    int parentPrevId = 0;
00240    int counter = pos;
00241    
00242    while ( parentId > 0 )
00243    {               
00244          if ( parentId == fullEvent[counter].mother2() ) // carbon copy, keep walking up
00245          {
00246             parentPrevId = parentId;
00247             counter = parentId;
00248             parentId = fullEvent[parentPrevId].mother1();
00249             continue;
00250          }
00251          
00252          // we get here if not a carbon copy
00253          
00254          // let's check if it's a normal process, etc.
00255          //
00256          if ( (parentId < parentPrevId) || parentId < fullEvent[counter].mother2() ) // normal process
00257          {
00258             
00259             // first of all, check if hard block
00260             if ( abs(fullEvent[counter].status()) == 22 || abs(fullEvent[counter].status()) == 23 )
00261             {
00262                // yes, it's the hard block
00263                // we got what we want, and can exit now !
00264                parentId = counter;
00265                break;
00266             }
00267             else
00268             {
00269                parentPrevId = parentId;
00270                parentId = fullEvent[parentPrevId].mother1();
00271             }
00272          }
00273          else if ( parentId > parentPrevId || parentId > pos ) // "circular"/"forward-pointing parent" - intermediate process
00274          {
00275             parentId = -1;
00276             break;
00277          }
00278 
00279          // additional checks... although we shouldn't be geeting here all that much...
00280          //      
00281          if ( abs(fullEvent[parentId].status()) == 22 || abs(fullEvent[parentId].status())== 23 ) // hard block
00282          {
00283             break;
00284          }       
00285          if ( abs(fullEvent[parentId].status()) < 22 ) // incoming
00286          {
00287             parentId = -1;
00288             break;
00289          } 
00290    }
00291    
00292    return parentId;
00293 
00294 }