CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/GeneratorInterface/PartonShowerVeto/src/JetMatchingMGFastJet.cc

Go to the documentation of this file.
00001 
00002 #include "GeneratorInterface/PartonShowerVeto/interface/JetMatchingMGFastJet.h"
00003 #include "FWCore/Utilities/interface/Exception.h"
00004 #include <iomanip>
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00007 #include "SimDataFormats/GeneratorProducts/interface/LHECommonBlocks.h"
00008 #include "GeneratorInterface/LHEInterface/interface/LHERunInfo.h"
00009 #include "GeneratorInterface/LHEInterface/interface/LHEEvent.h"
00010 
00011 namespace gen
00012 { 
00013 
00014 extern "C" {
00015 
00016         extern struct UPPRIV {
00017                 int     lnhin, lnhout;
00018                 int     mscal, ievnt;
00019                 int     ickkw, iscale;
00020         } uppriv_;
00021 
00022         extern struct MEMAIN {
00023                 double  etcjet, rclmax, etaclmax, qcut, showerkt, clfact;
00024                 int     maxjets, minjets, iexcfile, ktsche;
00025                 int     mektsc,nexcres, excres[30];
00026                 int     nqmatch,excproc,iexcproc[1000],iexcval[1000];
00027                 bool    nosingrad,jetprocs;
00028         } memain_;
00029 
00030 }
00031 
00032 bool JetMatchingMGFastJet::initAfterBeams()
00033 {
00034   
00035   if ( fIsInit ) return true;
00036   
00037   //
00038   doMerge = uppriv_.ickkw; 
00039   // doMerge = true;
00040   qCut    = memain_.qcut; // 
00041   nQmatch        = memain_.nqmatch;
00042   clFact         = 1.; // Default value 
00043                       // NOTE: ME2pythia seems to default to 1.5 - need to check !!!
00044                       // in general, needs to read key ALPSFACT from LHE file - fix CMSSW code !!!
00045   nJetMin        = memain_.minjets;
00046   nJetMax        = memain_.maxjets;
00047 
00048   etaJetMax      = memain_.etaclmax; 
00049 
00050   coneRadius     = 1.0; 
00051   jetAlgoPower   = 1; //   this is the kT algorithm !!! 
00052 
00053   // Matching procedure
00054   //
00055   qCutSq         = pow(qCut,2);  
00056   // this should be something like memaev_.iexc
00057   fExcLocal = true;
00058 
00059 
00060    // If not merging, then done (?????)
00061    //
00062    // if (!doMerge) return true;
00063   
00064    // Initialise chosen jet algorithm. 
00065    //
00066    fJetFinder = new fastjet::JetDefinition( fastjet::kt_algorithm, coneRadius );
00067    fClusJets.clear();
00068    fPtSortedJets.clear();  
00069    
00070    fIsInit = true;
00071    
00072    return true;
00073 
00074 }
00075 
00076 void JetMatchingMGFastJet::beforeHadronisation( const lhef::LHEEvent* lhee )
00077 {
00078 
00079    if (!runInitialized)
00080                 throw cms::Exception("Generator|PartonShowerVeto")
00081                         << "Run not initialized in JetMatchingMadgraph"
00082                         << std::endl;
00083    
00084    for (int i = 0; i < 3; i++) 
00085    {
00086     typeIdx[i].clear();
00087    }
00088 
00089    // Sort original process final state into light/heavy jets and 'other'.
00090    // Criteria:
00091    //   1 <= ID <= 5 and massless, or ID == 21 --> light jet (typeIdx[0])
00092    //   4 <= ID <= 6 and massive               --> heavy jet (typeIdx[1])
00093    //   All else                               --> other     (typeIdx[2])
00094    //
00095    const lhef::HEPEUP& hepeup = *lhee->getHEPEUP(); 
00096    int idx = 2; 
00097    for ( int i=0; i < hepeup.NUP; i++ )
00098    {
00099       if ( hepeup.ISTUP[i] < 0 ) continue;
00100       if ( hepeup.MOTHUP[i].first != 1 && hepeup.MOTHUP[i].second !=2 ) continue; // this way we skip entries that come
00101                                                                                   // from resonance decays;
00102                                                                                   // we only take those that descent
00103                                                                                   // directly from "incoming partons"
00104       idx = 2;
00105       if ( hepeup.IDUP[i] == ID_GLUON || (fabs(hepeup.IDUP[i]) <= nQmatch) ) // light jet
00106          // light jet
00107          idx = 0;
00108       else if ( fabs(hepeup.IDUP[i]) > nQmatch && fabs(hepeup.IDUP[i]) <= ID_TOP) // heavy jet
00109          idx = 1;
00110       // Store
00111       typeIdx[idx].push_back(i); 
00112    } 
00113    
00114    // NOTE: In principle, I should use exclusive, inclusive, or soup !!!   
00115    // should be like this:
00116    if ( soup )
00117    {
00118       int NPartons = typeIdx[0].size(); 
00119       fExcLocal = ( NPartons < nJetMax );
00120    }
00121    else
00122       fExcLocal = exclusive;
00123         
00124    return;
00125 
00126 }
00127 
00128 int JetMatchingMGFastJet::match( const lhef::LHEEvent* partonLevel, 
00129                                  const std::vector<fastjet::PseudoJet>* jetInput )
00130 {
00131          
00132    // Number of hard partons
00133    //
00134    int NPartons = typeIdx[0].size();
00135    
00136    fClusJets.clear();
00137    fPtSortedJets.clear();
00138    
00139    int ClusSeqNJets = 0;
00140 
00141    fastjet::ClusterSequence ClusSequence( *jetInput, *fJetFinder );
00142    
00143    if ( fExcLocal )
00144    {
00145       fClusJets = ClusSequence.exclusive_jets( qCutSq );
00146    }
00147    else
00148    {
00149       fClusJets = ClusSequence.inclusive_jets( qCut ); 
00150    }
00151    
00152    ClusSeqNJets = fClusJets.size(); 
00153    
00154    
00155    if ( ClusSeqNJets < NPartons ) return LESS_JETS;
00156    
00157    double localQcutSq = qCutSq;
00158    
00159    if ( fExcLocal ) // exclusive
00160    {
00161       if( ClusSeqNJets > NPartons ) return MORE_JETS;
00162    }
00163    else // inclusive 
00164    {
00165       fPtSortedJets = fastjet::sorted_by_pt( *jetInput );
00166       localQcutSq = std::max( qCutSq, fPtSortedJets[0].pt2() );
00167       fClusJets = ClusSequence.exclusive_jets( NPartons ); // override
00168       ClusSeqNJets = NPartons;
00169    }
00170    
00171   if( clFact != 0 ) localQcutSq *= pow(clFact,2);
00172   
00173   std::vector<fastjet::PseudoJet> MatchingInput;
00174   
00175   std::vector<bool> jetAssigned;
00176   jetAssigned.assign( fClusJets.size(), false );
00177 
00178   int counter = 0;
00179   
00180   const lhef::HEPEUP& hepeup = *partonLevel->getHEPEUP();
00181     
00182   while ( counter < NPartons )
00183   {
00184 
00185      MatchingInput.clear();
00186 
00187      for ( int i=0; i<ClusSeqNJets; i++ )
00188      {
00189         if ( jetAssigned[i] ) continue;
00190         if ( i == NPartons )  break;
00191         //
00192         // this looks "awkward" but this way we do NOT pass in cluster_hist_index
00193         // which would actually indicate position in "history" of the "master" ClusSeq
00194         //
00195         fastjet::PseudoJet exjet = fClusJets[i];
00196         MatchingInput.push_back( fastjet::PseudoJet( exjet.px(), exjet.py(), exjet.pz(), exjet.e() ) );
00197         MatchingInput.back().set_user_index(i);
00198      }
00199 
00200      int idx = typeIdx[0][counter];
00201      MatchingInput.push_back( fastjet::PseudoJet( hepeup.PUP[idx][0],
00202                                                   hepeup.PUP[idx][1],
00203                                                   hepeup.PUP[idx][2],
00204                                                   hepeup.PUP[idx][3]) );
00205 
00206      //
00207      // in principle, one can use ClusterSequence::n_particles() 
00208      //
00209      int NBefore = MatchingInput.size();
00210 
00211      // Create new clustering object - which includes the 1st clustering run !!!
00212      // NOTE-1: it better be a new object for each try, or history will be "too crowded".
00213      // NOTE-2: when created, the it ALWAYS makes at least 1 clustering step, so in most
00214      //         cases at least 1 new jet object will be added; although in some cases 
00215      //         no jet is added, if the system doesn't cluster to anything (for example,
00216      //         input jet(s) and parton(s) are going opposite directions)
00217      //
00218      fastjet::ClusterSequence ClusSeq( MatchingInput, *fJetFinder );
00219      
00220      const std::vector<fastjet::PseudoJet>&        output = ClusSeq.jets();
00221      int NClusJets = output.size() - NBefore;
00222      
00223      // 
00224      // JVY - I think this is the right idea:
00225      // at least 1 (one) new jet needs to be added
00226      // however, need to double check details and refine, especially for inclusive mode
00227      //
00228      // 
00229      if ( NClusJets < 1 )
00230      {
00231         return UNMATCHED_PARTON;
00232      }
00233      //
00234      // very unlikely case but let's do it just to be safe
00235      //
00236      if ( NClusJets >= NBefore )
00237      {
00238         return MORE_JETS;
00239      }
00240       
00241      // Now browse history and see how close the clustering distance
00242      //
00243      // NOTE: Remember, there maybe more than one new jet in the list (for example,
00244      //       for process=2,3,...);
00245      //       in this case we take the ** first ** one that satisfies the distance/cut,
00246      //       which is ** typically ** also the best one 
00247      //
00248      bool matched = false;
00249      const std::vector<fastjet::ClusterSequence::history_element>&  history = ClusSeq.history();
00250               
00251      // for ( unsigned int i=nBefore; i<history.size(); i++ )
00252      for ( unsigned int i=NBefore; i<output.size(); i++ )
00253      {
00254 
00255         int hidx = output[i].cluster_sequence_history_index();
00256         double dNext = history[hidx].dij;
00257         if ( dNext < localQcutSq )
00258         {
00259            //
00260            // the way we form input, parent1 is always jet, and parent2 can be any,
00261            // but if it's a good match/cluster, it should point at the parton
00262            //
00263            int parent1 = history[hidx].parent1;
00264            int parent2 = history[hidx].parent2;
00265            if ( parent1 < 0 || parent2 < 0 ) break; // bad cluster, no match
00266            //
00267            // pull up jet's "global" index 
00268            //
00269            int pidx = MatchingInput[parent1].user_index();
00270            jetAssigned[pidx] = true; 
00271            matched = true;
00272            break;
00273         }
00274      }
00275      if ( !matched )
00276      {
00277         return UNMATCHED_PARTON;
00278      }
00279      
00280      counter++;   
00281   }
00282 
00283    // Now save some info for DJR analysis (if requested).
00284    // This mimics what is done in ME2pythia.f
00285    // Basically, NJets and matching scale for these 4 cases:
00286    // 1->0, 2->1, 3->2, and 4->3
00287    //
00288    if ( fDJROutFlag > 0 )
00289    {   
00290       std::vector<double> MergingScale;
00291       MergingScale.clear();
00292       for ( int nj=0; nj<4; nj++ )
00293       {
00294          double dmscale2 = ClusSequence.exclusive_dmerge( nj );
00295          double dmscale = sqrt( dmscale2 );
00296          MergingScale.push_back( dmscale );
00297       }
00298       fDJROutput.open( "events.tree", std::ios_base::app );  
00299       double dNJets = (double)NPartons;
00300       fDJROutput << " " << dNJets << " " << MergingScale[0] << " " 
00301                                          << MergingScale[1] << " "
00302                                          << MergingScale[2] << " "
00303                                          << MergingScale[3] << std::endl;  
00304       fDJROutput.close();
00305    }
00306  
00307   return NONE;
00308 
00309 }
00310 
00311 } // end namespace