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
00040 qCut = memain_.qcut;
00041 nQmatch = memain_.nqmatch;
00042 clFact = 1.;
00043
00044
00045 nJetMin = memain_.minjets;
00046 nJetMax = memain_.maxjets;
00047
00048 etaJetMax = memain_.etaclmax;
00049
00050 coneRadius = 1.0;
00051 jetAlgoPower = 1;
00052
00053
00054
00055 qCutSq = pow(qCut,2);
00056
00057 fExcLocal = true;
00058
00059
00060
00061
00062
00063
00064
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
00090
00091
00092
00093
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;
00101
00102
00103
00104 idx = 2;
00105 if ( hepeup.IDUP[i] == ID_GLUON || (fabs(hepeup.IDUP[i]) <= nQmatch) )
00106
00107 idx = 0;
00108 else if ( fabs(hepeup.IDUP[i]) > nQmatch && fabs(hepeup.IDUP[i]) <= ID_TOP)
00109 idx = 1;
00110
00111 typeIdx[idx].push_back(i);
00112 }
00113
00114
00115
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
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 )
00160 {
00161 if( ClusSeqNJets > NPartons ) return MORE_JETS;
00162 }
00163 else
00164 {
00165 fPtSortedJets = fastjet::sorted_by_pt( *jetInput );
00166 localQcutSq = std::max( qCutSq, fPtSortedJets[0].pt2() );
00167 fClusJets = ClusSequence.exclusive_jets( NPartons );
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
00193
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
00208
00209 int NBefore = MatchingInput.size();
00210
00211
00212
00213
00214
00215
00216
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
00225
00226
00227
00228
00229 if ( NClusJets < 1 )
00230 {
00231 return UNMATCHED_PARTON;
00232 }
00233
00234
00235
00236 if ( NClusJets >= NBefore )
00237 {
00238 return MORE_JETS;
00239 }
00240
00241
00242
00243
00244
00245
00246
00247
00248 bool matched = false;
00249 const std::vector<fastjet::ClusterSequence::history_element>& history = ClusSeq.history();
00250
00251
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
00261
00262
00263 int parent1 = history[hidx].parent1;
00264 int parent2 = history[hidx].parent2;
00265 if ( parent1 < 0 || parent2 < 0 ) break;
00266
00267
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
00284
00285
00286
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 }