10 #include <boost/shared_ptr.hpp> 11 #include <boost/algorithm/string/trim.hpp> 26 std::memcpy(
value, str.c_str(), len);
34 extern struct UPPRIV {
40 extern struct MEMAIN {
41 double etcjet, rclmax, etaclmax, qcut, showerkt, clfact;
42 int maxjets, minjets, iexcfile, ktsche;
43 int mektsc,nexcres, excres[30];
44 int nqmatch,excproc,iexcproc[1000],iexcval[1000];
45 bool nosingrad,jetprocs;
48 extern struct OUTTREE {
58 std::istringstream ss(value);
68 if (!result.empty() && result[0] ==
'\'')
69 result = result.substr(1);
70 if (!result.empty() && result[result.length() - 1] ==
'\'')
71 result.resize(result.length() - 1);
80 value.begin(), (
int(*)(
int))std::toupper);
81 return value ==
"T" || value ==
"Y" || value==
"True" ||
82 value ==
"1" || value ==
".TRUE.";
88 const std::map<std::string, std::string> ¶ms,
91 std::map<std::string, std::string>::const_iterator
pos =
93 if (pos == params.end())
95 return parseParameter<T>(pos->second);
100 const T &defValue)
const 102 return getParameter(mgParams, var, defValue);
106 static std::map<std::string, std::string>
109 std::map<std::string, std::string> params;
111 for(std::vector<std::string>::const_iterator iter = header.begin();
112 iter != header.end(); ++iter) {
114 if (line.empty() || line[0] ==
'#')
118 if (pos != std::string::npos)
121 pos = line.find(
'=');
122 if (pos == std::string::npos)
126 boost::algorithm::trim_copy(line.substr(pos + 1));
128 boost::algorithm::trim_copy(line.substr(0, pos));
138 const std::map<std::string, std::string> ¶ms,
142 param = getParameter(params, name, param);
146 <<
"The MGParamCMS header does not specify the jet " 147 "matching parameter \"" << name <<
"\", but it " 148 "is requested by the CMSSW configuration." 155 runInitialized(
false),
160 if (mode ==
"inclusive") {
163 }
else if (mode ==
"exclusive") {
166 }
else if (mode ==
"auto")
170 <<
"MGFastJet jet matching scheme requires \"mode\" " 171 "parameter to be set to either \"inclusive\", " 172 "\"exclusive\" or \"auto\"." << std::endl;
186 std::vector<std::string> elems;
187 std::stringstream ss(list_excres);
190 while(std::getline(ss, item,
',')) {
191 elems.push_back(item);
256 <<
"Run not initialized in JetMatchingMGFastJet" 259 for (
int i = 0;
i < 3;
i++)
272 for (
int i=0;
i < hepeup.
NUP;
i++ )
274 if ( hepeup.
ISTUP[
i] < 0 )
continue;
275 if ( hepeup.
MOTHUP[
i].first != 1 && hepeup.
MOTHUP[
i].second !=2 )
continue;
293 int NPartons =
typeIdx[0].size();
310 std::map<std::string, std::string>
parameters;
315 <<
"In order to use MGFastJet jet matching, " 316 "the input file has to contain the corresponding " 317 "MadGraph headers." << std::endl;
323 std::vector<Param> params;
324 std::vector<Param>
values;
325 for(std::map<std::string, std::string>::const_iterator iter =
327 params.push_back(
" " + iter->first);
328 values.push_back(iter->second);
339 std::map<std::string, std::string> mgInfoCMS =
parseHeader(header);
341 for(std::map<std::string, std::string>::const_iterator iter =
342 mgInfoCMS.begin(); iter != mgInfoCMS.end(); ++iter) {
343 std::cout<<
"mgInfoCMS: "<<iter->first<<
" "<<iter->second<<std::endl;
365 const std::vector<fastjet::PseudoJet>* jetInput )
370 int NPartons =
typeIdx[0].size();
375 int ClusSeqNJets = 0;
377 fastjet::ClusterSequence ClusSequence( *jetInput, *
fJetFinder );
391 if ( ClusSeqNJets < NPartons )
return LESS_JETS;
393 double localQcutSq =
qCutSq;
397 if( ClusSeqNJets > NPartons )
return MORE_JETS;
403 fClusJets = ClusSequence.exclusive_jets( NPartons );
404 ClusSeqNJets = NPartons;
409 std::vector<fastjet::PseudoJet> MatchingInput;
411 std::vector<bool> jetAssigned;
418 while ( counter < NPartons )
421 MatchingInput.clear();
423 for (
int i=0;
i<ClusSeqNJets;
i++ )
425 if ( jetAssigned[
i] )
continue;
426 if ( i == NPartons )
break;
432 MatchingInput.push_back( fastjet::PseudoJet( exjet.px(), exjet.py(), exjet.pz(), exjet.e() ) );
433 MatchingInput.back().set_user_index(i);
437 MatchingInput.push_back( fastjet::PseudoJet( hepeup.
PUP[idx][0],
440 hepeup.
PUP[idx][3]) );
445 int NBefore = MatchingInput.size();
454 fastjet::ClusterSequence ClusSeq( MatchingInput, *
fJetFinder );
456 const std::vector<fastjet::PseudoJet>&
output = ClusSeq.jets();
457 int NClusJets = output.size() - NBefore;
472 if ( NClusJets >= NBefore )
485 const std::vector<fastjet::ClusterSequence::history_element>& history = ClusSeq.history();
488 for (
unsigned int i=NBefore;
i<output.size();
i++ )
491 int hidx = output[
i].cluster_sequence_history_index();
492 double dNext = history[hidx].dij;
493 if ( dNext < localQcutSq )
499 int parent1 = history[hidx].parent1;
500 int parent2 = history[hidx].parent2;
501 if ( parent1 < 0 || parent2 < 0 )
break;
505 int pidx = MatchingInput[parent1].user_index();
506 jetAssigned[pidx] =
true;
526 std::vector<double> MergingScale;
527 MergingScale.clear();
528 for (
int nj=0; nj<4; nj++ )
530 double dmscale2 = ClusSequence.exclusive_dmerge( nj );
531 double dmscale =
sqrt( dmscale2 );
532 MergingScale.push_back( dmscale );
534 fDJROutput.open(
"events.tree", std::ios_base::app );
535 double dNJets = (double)NPartons;
536 fDJROutput <<
" " << dNJets <<
" " << MergingScale[0] <<
" " 537 << MergingScale[1] <<
" " 538 << MergingScale[2] <<
" " 539 << MergingScale[3] << std::endl;
T getParameter(std::string const &) const
struct gen::OUTTREE outtree_
struct gen::MEMAIN memain_
fastjet::JetDefinition * fJetFinder
JetMatchingMGFastJet(const edm::ParameterSet ¶ms)
static T getParameter(const std::map< std::string, std::string > ¶ms, const std::string &var, const T &defValue=T())
static void updateOrDie(const std::map< std::string, std::string > ¶ms, T ¶m, const std::string &name)
const HEPEUP * getHEPEUP() const
void beforeHadronisation(const lhef::LHEEvent *) override
int match(const lhef::LHEEvent *partonLevel, const std::vector< fastjet::PseudoJet > *jetInput) override
std::vector< std::pair< int, int > > MOTHUP
std::vector< FiveVector > PUP
Abs< T >::type abs(const T &t)
bool initAfterBeams() override
static std::map< std::string, std::string > parseHeader(const std::vector< std::string > &header)
std::vector< fastjet::PseudoJet > fClusJets
static T parseParameter(const std::string &value)
std::vector< fastjet::PseudoJet > fPtSortedJets
std::map< std::string, std::string > mgParams
std::vector< int > typeIdx[3]
double getJetEtaMax() const override
void init(const lhef::LHERunInfo *runInfo) override
static std::atomic< unsigned int > counter
std::vector< std::string > findHeader(const std::string &tag) const
Power< A, B >::type pow(const A &a, const B &b)
struct gen::UPPRIV uppriv_