CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Member Functions | Private Types | Private Member Functions | Static Private Member Functions | Private Attributes
gen::JetMatchingMGFastJet Class Reference

#include <JetMatchingMGFastJet.h>

Inheritance diagram for gen::JetMatchingMGFastJet:
gen::JetMatching

Public Member Functions

const std::vector< int > * getPartonList ()
 
 JetMatchingMGFastJet (const edm::ParameterSet &params)
 
template<>
std::string parseParameter (const std::string &value)
 
template<>
bool parseParameter (const std::string &value_)
 
 ~JetMatchingMGFastJet ()
 
- Public Member Functions inherited from gen::JetMatching
virtual std::set< std::string > capabilities () const
 
bool isMatchingDone ()
 
 JetMatching (const edm::ParameterSet &params)
 
void resetMatchingStatus ()
 
virtual ~JetMatching ()
 

Protected Member Functions

virtual void beforeHadronisation (const lhef::LHEEvent *)
 
virtual void beforeHadronisationExec ()
 
virtual double getJetEtaMax () const
 
virtual void init (const lhef::LHERunInfo *runInfo)
 
virtual bool initAfterBeams ()
 
virtual int match (const lhef::LHEEvent *partonLevel, const std::vector< fastjet::PseudoJet > *jetInput)
 

Private Types

enum  partonTypes { ID_TOP =6, ID_GLUON =21, ID_PHOTON =22 }
 
enum  vetoStatus {
  NONE, LESS_JETS, MORE_JETS, HARD_JET,
  UNMATCHED_PARTON
}
 

Private Member Functions

template<typename T >
T getParameter (const std::string &var, const T &defValue=T()) const
 

Static Private Member Functions

template<typename T >
static T getParameter (const std::map< std::string, std::string > &params, const std::string &var, const T &defValue=T())
 
template<typename T >
static T parseParameter (const std::string &value)
 
template<typename T >
static void updateOrDie (const std::map< std::string, std::string > &params, T &param, const std::string &name)
 

Private Attributes

double clFact
 
double coneRadius
 
bool doMerge
 
double etaJetMax
 
bool exclusive
 
std::vector< fastjet::PseudoJet > fClusJets
 
int fDJROutFlag
 
std::ofstream fDJROutput
 
bool fExcLocal
 
bool fIsInit
 
fastjet::JetDefinition * fJetFinder
 
std::vector< fastjet::PseudoJet > fPtSortedJets
 
int jetAlgoPower
 
std::map< std::string,
std::string > 
mgParams
 
int nJetMax
 
int nJetMin
 
int nQmatch
 
double qCut
 
double qCutSq
 
bool runInitialized
 
bool soup
 
std::vector< int > typeIdx [3]
 

Additional Inherited Members

- Static Public Member Functions inherited from gen::JetMatching
static std::auto_ptr< JetMatchingcreate (const edm::ParameterSet &params)
 
- Protected Attributes inherited from gen::JetMatching
bool fMatchingStatus
 

Detailed Description

Definition at line 32 of file JetMatchingMGFastJet.h.

Member Enumeration Documentation

Constructor & Destructor Documentation

gen::JetMatchingMGFastJet::JetMatchingMGFastJet ( const edm::ParameterSet params)

Definition at line 153 of file JetMatchingMGFastJet.cc.

References gen::MEMAIN::clfact, gen::MEMAIN::etaclmax, gen::MEMAIN::etcjet, Exception, exclusive, gen::MEMAIN::excres, fDJROutFlag, gen::OUTTREE::flag, edm::ParameterSet::getParameter(), cmsHarvester::index, gen::MEMAIN::ktsche, gen::MEMAIN::maxjets, gen::memain_, gen::MEMAIN::minjets, alignBH_cfg::mode, gen::MEMAIN::nexcres, gen::MEMAIN::nqmatch, gen::outtree_, gen::MEMAIN::qcut, gen::MEMAIN::rclmax, gen::MEMAIN::showerkt, soup, contentValuesCheck::ss, and AlCaHLTBitMon_QueryRunRegistry::string.

153  :
154  JetMatching(params),
155  runInitialized(false),
156  fJetFinder(0),
157  fIsInit(false)
158 {
159  std::string mode = params.getParameter<std::string>("mode");
160  if (mode == "inclusive") {
161  soup = false;
162  exclusive = false;
163  } else if (mode == "exclusive") {
164  soup = false;
165  exclusive = true;
166  } else if (mode == "auto")
167  soup = true;
168  else
169  throw cms::Exception("Generator|LHEInterface")
170  << "MGFastJet jet matching scheme requires \"mode\" "
171  "parameter to be set to either \"inclusive\", "
172  "\"exclusive\" or \"auto\"." << std::endl;
173 
174  memain_.etcjet = 0.;
175  memain_.rclmax = 0.0;
176  memain_.clfact = 0.0;
177  memain_.ktsche = 0.0;
178  memain_.etaclmax = params.getParameter<double>("MEMAIN_etaclmax");
179  memain_.qcut = params.getParameter<double>("MEMAIN_qcut");
180  memain_.minjets = params.getParameter<int>("MEMAIN_minjets");
181  memain_.maxjets = params.getParameter<int>("MEMAIN_maxjets");
182  memain_.showerkt = params.getParameter<double>("MEMAIN_showerkt");
183  memain_.nqmatch = params.getParameter<int>("MEMAIN_nqmatch");
184  outtree_.flag = params.getParameter<int>("outTree_flag");
185  std::string list_excres = params.getParameter<std::string>("MEMAIN_excres");
186  std::vector<std::string> elems;
187  std::stringstream ss(list_excres);
188  std::string item;
189  int index=0;
190  while(std::getline(ss, item, ',')) {
191  elems.push_back(item);
192  memain_.excres[index]=std::atoi(item.c_str());
193  index++;
194  }
196 
197  fDJROutFlag = params.getParameter<int>("outTree_flag");
198 }
T getParameter(std::string const &) const
struct gen::OUTTREE outtree_
struct gen::MEMAIN memain_
fastjet::JetDefinition * fJetFinder
JetMatching(const edm::ParameterSet &params)
Definition: JetMatching.cc:18
gen::JetMatchingMGFastJet::~JetMatchingMGFastJet ( )
inline

Definition at line 39 of file JetMatchingMGFastJet.h.

References fJetFinder.

39 { if (fJetFinder) delete fJetFinder; }
fastjet::JetDefinition * fJetFinder

Member Function Documentation

void gen::JetMatchingMGFastJet::beforeHadronisation ( const lhef::LHEEvent lhee)
protectedvirtual

Reimplemented from gen::JetMatching.

Definition at line 251 of file JetMatchingMGFastJet.cc.

References Exception, exclusive, fExcLocal, lhef::LHEEvent::getHEPEUP(), i, ID_GLUON, ID_TOP, lhef::HEPEUP::IDUP, customizeTrackingMonitorSeedNumber::idx, lhef::HEPEUP::ISTUP, lhef::HEPEUP::MOTHUP, nJetMax, nQmatch, lhef::HEPEUP::NUP, runInitialized, soup, and typeIdx.

252 {
253 
254  if (!runInitialized)
255  throw cms::Exception("Generator|PartonShowerVeto")
256  << "Run not initialized in JetMatchingMGFastJet"
257  << std::endl;
258 
259  for (int i = 0; i < 3; i++)
260  {
261  typeIdx[i].clear();
262  }
263 
264  // Sort original process final state into light/heavy jets and 'other'.
265  // Criteria:
266  // 1 <= ID <= 5 and massless, or ID == 21 --> light jet (typeIdx[0])
267  // 4 <= ID <= 6 and massive --> heavy jet (typeIdx[1])
268  // All else --> other (typeIdx[2])
269  //
270  const lhef::HEPEUP& hepeup = *lhee->getHEPEUP();
271  int idx = 2;
272  for ( int i=0; i < hepeup.NUP; i++ )
273  {
274  if ( hepeup.ISTUP[i] < 0 ) continue;
275  if ( hepeup.MOTHUP[i].first != 1 && hepeup.MOTHUP[i].second !=2 ) continue; // this way we skip entries that come
276  // from resonance decays;
277  // we only take those that descent
278  // directly from "incoming partons"
279  idx = 2;
280  if ( hepeup.IDUP[i] == ID_GLUON || (fabs(hepeup.IDUP[i]) <= nQmatch) ) // light jet
281  // light jet
282  idx = 0;
283  else if ( fabs(hepeup.IDUP[i]) > nQmatch && fabs(hepeup.IDUP[i]) <= ID_TOP) // heavy jet
284  idx = 1;
285  // Store
286  typeIdx[idx].push_back(i);
287  }
288 
289  // NOTE: In principle, I should use exclusive, inclusive, or soup !!!
290  // should be like this:
291  if ( soup )
292  {
293  int NPartons = typeIdx[0].size();
294  fExcLocal = ( NPartons < nJetMax );
295  }
296  else
298 
299  return;
300 
301 }
int i
Definition: DBlmapReader.cc:9
const HEPEUP * getHEPEUP() const
Definition: LHEEvent.h:43
std::vector< std::pair< int, int > > MOTHUP
Definition: LesHouches.h:236
std::vector< int > ISTUP
Definition: LesHouches.h:230
std::vector< int > IDUP
Definition: LesHouches.h:225
std::vector< int > typeIdx[3]
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
virtual void gen::JetMatchingMGFastJet::beforeHadronisationExec ( )
inlineprotectedvirtual

Reimplemented from gen::JetMatching.

Definition at line 49 of file JetMatchingMGFastJet.h.

49 { return; }
double gen::JetMatchingMGFastJet::getJetEtaMax ( ) const
protectedvirtual

Implements gen::JetMatching.

Definition at line 201 of file JetMatchingMGFastJet.cc.

References gen::MEMAIN::etaclmax, and gen::memain_.

202 {
203  return memain_.etaclmax;
204 }
struct gen::MEMAIN memain_
template<typename T >
T gen::JetMatchingMGFastJet::getParameter ( const std::map< std::string, std::string > &  params,
const std::string &  var,
const T defValue = T() 
)
staticprivate

Definition at line 87 of file JetMatchingMGFastJet.cc.

Referenced by getParameter(), and updateOrDie().

90 {
91  std::map<std::string, std::string>::const_iterator pos =
92  params.find(var);
93  if (pos == params.end())
94  return defValue;
95  return parseParameter<T>(pos->second);
96 }
template<typename T >
T gen::JetMatchingMGFastJet::getParameter ( const std::string &  var,
const T defValue = T() 
) const
private

Definition at line 99 of file JetMatchingMGFastJet.cc.

References getParameter(), and mgParams.

101 {
102  return getParameter(mgParams, var, defValue);
103 }
static T getParameter(const std::map< std::string, std::string > &params, const std::string &var, const T &defValue=T())
std::map< std::string, std::string > mgParams
const std::vector<int>* gen::JetMatchingMGFastJet::getPartonList ( )
inlinevirtual

Reimplemented from gen::JetMatching.

Definition at line 41 of file JetMatchingMGFastJet.h.

References typeIdx.

41 { return typeIdx; }
std::vector< int > typeIdx[3]
void gen::JetMatchingMGFastJet::init ( const lhef::LHERunInfo runInfo)
protectedvirtual

Reimplemented from gen::JetMatching.

Definition at line 304 of file JetMatchingMGFastJet.cc.

References gather_cfg::cout, gen::MEMAIN::etaclmax, lhef::LHERunInfo::findHeader(), fIsInit, gen::UPPRIV::ickkw, initAfterBeams(), gen::MEMAIN::maxjets, gen::MEMAIN::mektsc, gen::memain_, mgParams, gen::MEMAIN::minjets, gen::MEMAIN::nqmatch, Parameters::parameters, gen::parseHeader(), gen::MEMAIN::qcut, runInitialized, gen::MEMAIN::showerkt, updateOrDie(), gen::uppriv_, and makeHLTPrescaleTable::values.

305 {
306  if (fIsInit) return;
307 
308  // read MGFastJet run card
309 
310  std::map<std::string, std::string> parameters;
311 
312  std::vector<std::string> header = runInfo->findHeader("MGRunCard");
313  if (header.empty())
314  throw cms::Exception("Generator|PartonShowerVeto")
315  << "In order to use MGFastJet jet matching, "
316  "the input file has to contain the corresponding "
317  "MadGraph headers." << std::endl;
318 
319  mgParams = parseHeader(header);
320 
321  // set variables in common block
322 
323  std::vector<Param> params;
324  std::vector<Param> values;
325  for(std::map<std::string, std::string>::const_iterator iter =
326  mgParams.begin(); iter != mgParams.end(); ++iter) {
327  params.push_back(" " + iter->first);
328  values.push_back(iter->second);
329 
330  }
331 
332  // set MG matching parameters
333 
334  uppriv_.ickkw = getParameter<int>("ickkw", 0);
335  memain_.mektsc = getParameter<int>("ktscheme", 0);
336 
337  header = runInfo->findHeader("MGParamCMS");
338 
339  std::map<std::string, std::string> mgInfoCMS = parseHeader(header);
340 
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;
344 
345  }
346 
347  updateOrDie(mgInfoCMS, memain_.etaclmax, "etaclmax");
348  updateOrDie(mgInfoCMS, memain_.qcut, "qcut");
349  updateOrDie(mgInfoCMS, memain_.minjets, "minjets");
350  updateOrDie(mgInfoCMS, memain_.maxjets, "maxjets");
351  updateOrDie(mgInfoCMS, memain_.showerkt, "showerkt");
352  updateOrDie(mgInfoCMS, memain_.nqmatch, "nqmatch");
353 
354  runInitialized = true;
355 
356  initAfterBeams();
357 
358  fIsInit=true;
359 
360  return;
361 }
struct gen::MEMAIN memain_
dictionary parameters
Definition: Parameters.py:2
static void updateOrDie(const std::map< std::string, std::string > &params, T &param, const std::string &name)
static std::map< std::string, std::string > parseHeader(const std::vector< std::string > &header)
std::map< std::string, std::string > mgParams
tuple cout
Definition: gather_cfg.py:121
std::vector< std::string > findHeader(const std::string &tag) const
Definition: LHERunInfo.cc:523
struct gen::UPPRIV uppriv_
bool gen::JetMatchingMGFastJet::initAfterBeams ( )
protectedvirtual

Reimplemented from gen::JetMatching.

Definition at line 207 of file JetMatchingMGFastJet.cc.

References clFact, coneRadius, doMerge, gen::MEMAIN::etaclmax, etaJetMax, fClusJets, fExcLocal, fIsInit, fJetFinder, fPtSortedJets, gen::UPPRIV::ickkw, jetAlgoPower, gen::MEMAIN::maxjets, gen::memain_, gen::MEMAIN::minjets, nJetMax, nJetMin, gen::MEMAIN::nqmatch, nQmatch, funct::pow(), gen::MEMAIN::qcut, qCut, qCutSq, and gen::uppriv_.

Referenced by init().

208 {
209 
210  if ( fIsInit ) return true;
211 
212  //
213  doMerge = uppriv_.ickkw;
214  // doMerge = true;
215  qCut = memain_.qcut; //
217  clFact = 1.; // Default value
218  // NOTE: ME2pythia seems to default to 1.5 - need to check !!!
219  // in general, needs to read key ALPSFACT from LHE file - fix CMSSW code !!!
222 
224 
225  coneRadius = 1.0;
226  jetAlgoPower = 1; // this is the kT algorithm !!!
227 
228  // Matching procedure
229  //
230  qCutSq = pow(qCut,2);
231  // this should be something like memaev_.iexc
232  fExcLocal = true;
233 
234 
235  // If not merging, then done (?????)
236  //
237  // if (!doMerge) return true;
238 
239  // Initialise chosen jet algorithm.
240  //
241  fJetFinder = new fastjet::JetDefinition( fastjet::kt_algorithm, coneRadius );
242  fClusJets.clear();
243  fPtSortedJets.clear();
244 
245  fIsInit = true;
246 
247  return true;
248 
249 }
struct gen::MEMAIN memain_
fastjet::JetDefinition * fJetFinder
std::vector< fastjet::PseudoJet > fClusJets
std::vector< fastjet::PseudoJet > fPtSortedJets
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
struct gen::UPPRIV uppriv_
int gen::JetMatchingMGFastJet::match ( const lhef::LHEEvent partonLevel,
const std::vector< fastjet::PseudoJet > *  jetInput 
)
protectedvirtual

Implements gen::JetMatching.

Definition at line 364 of file JetMatchingMGFastJet.cc.

References clFact, counter, edm::false, fClusJets, fDJROutFlag, fDJROutput, fExcLocal, fJetFinder, fPtSortedJets, lhef::LHEEvent::getHEPEUP(), i, customizeTrackingMonitorSeedNumber::idx, LESS_JETS, bookConverter::max, MORE_JETS, NONE, convertSQLitetoXML_cfg::output, funct::pow(), lhef::HEPEUP::PUP, qCut, qCutSq, mathSSE::sqrt(), typeIdx, and UNMATCHED_PARTON.

366 {
367 
368  // Number of hard partons
369  //
370  int NPartons = typeIdx[0].size();
371 
372  fClusJets.clear();
373  fPtSortedJets.clear();
374 
375  int ClusSeqNJets = 0;
376 
377  fastjet::ClusterSequence ClusSequence( *jetInput, *fJetFinder );
378 
379  if ( fExcLocal )
380  {
381  fClusJets = ClusSequence.exclusive_jets( qCutSq );
382  }
383  else
384  {
385  fClusJets = ClusSequence.inclusive_jets( qCut );
386  }
387 
388  ClusSeqNJets = fClusJets.size();
389 
390 
391  if ( ClusSeqNJets < NPartons ) return LESS_JETS;
392 
393  double localQcutSq = qCutSq;
394 
395  if ( fExcLocal ) // exclusive
396  {
397  if( ClusSeqNJets > NPartons ) return MORE_JETS;
398  }
399  else // inclusive
400  {
401  fPtSortedJets = fastjet::sorted_by_pt( *jetInput );
402  localQcutSq = std::max( qCutSq, fPtSortedJets[0].pt2() );
403  fClusJets = ClusSequence.exclusive_jets( NPartons ); // override
404  ClusSeqNJets = NPartons;
405  }
406 
407  if( clFact != 0 ) localQcutSq *= pow(clFact,2);
408 
409  std::vector<fastjet::PseudoJet> MatchingInput;
410 
411  std::vector<bool> jetAssigned;
412  jetAssigned.assign( fClusJets.size(), false );
413 
414  int counter = 0;
415 
416  const lhef::HEPEUP& hepeup = *partonLevel->getHEPEUP();
417 
418  while ( counter < NPartons )
419  {
420 
421  MatchingInput.clear();
422 
423  for ( int i=0; i<ClusSeqNJets; i++ )
424  {
425  if ( jetAssigned[i] ) continue;
426  if ( i == NPartons ) break;
427  //
428  // this looks "awkward" but this way we do NOT pass in cluster_hist_index
429  // which would actually indicate position in "history" of the "master" ClusSeq
430  //
431  fastjet::PseudoJet exjet = fClusJets[i];
432  MatchingInput.push_back( fastjet::PseudoJet( exjet.px(), exjet.py(), exjet.pz(), exjet.e() ) );
433  MatchingInput.back().set_user_index(i);
434  }
435 
436  int idx = typeIdx[0][counter];
437  MatchingInput.push_back( fastjet::PseudoJet( hepeup.PUP[idx][0],
438  hepeup.PUP[idx][1],
439  hepeup.PUP[idx][2],
440  hepeup.PUP[idx][3]) );
441 
442  //
443  // in principle, one can use ClusterSequence::n_particles()
444  //
445  int NBefore = MatchingInput.size();
446 
447  // Create new clustering object - which includes the 1st clustering run !!!
448  // NOTE-1: it better be a new object for each try, or history will be "too crowded".
449  // NOTE-2: when created, the it ALWAYS makes at least 1 clustering step, so in most
450  // cases at least 1 new jet object will be added; although in some cases
451  // no jet is added, if the system doesn't cluster to anything (for example,
452  // input jet(s) and parton(s) are going opposite directions)
453  //
454  fastjet::ClusterSequence ClusSeq( MatchingInput, *fJetFinder );
455 
456  const std::vector<fastjet::PseudoJet>& output = ClusSeq.jets();
457  int NClusJets = output.size() - NBefore;
458 
459  //
460  // JVY - I think this is the right idea:
461  // at least 1 (one) new jet needs to be added
462  // however, need to double check details and refine, especially for inclusive mode
463  //
464  //
465  if ( NClusJets < 1 )
466  {
467  return UNMATCHED_PARTON;
468  }
469  //
470  // very unlikely case but let's do it just to be safe
471  //
472  if ( NClusJets >= NBefore )
473  {
474  return MORE_JETS;
475  }
476 
477  // Now browse history and see how close the clustering distance
478  //
479  // NOTE: Remember, there maybe more than one new jet in the list (for example,
480  // for process=2,3,...);
481  // in this case we take the ** first ** one that satisfies the distance/cut,
482  // which is ** typically ** also the best one
483  //
484  bool matched = false;
485  const std::vector<fastjet::ClusterSequence::history_element>& history = ClusSeq.history();
486 
487  // for ( unsigned int i=nBefore; i<history.size(); i++ )
488  for ( unsigned int i=NBefore; i<output.size(); i++ )
489  {
490 
491  int hidx = output[i].cluster_sequence_history_index();
492  double dNext = history[hidx].dij;
493  if ( dNext < localQcutSq )
494  {
495  //
496  // the way we form input, parent1 is always jet, and parent2 can be any,
497  // but if it's a good match/cluster, it should point at the parton
498  //
499  int parent1 = history[hidx].parent1;
500  int parent2 = history[hidx].parent2;
501  if ( parent1 < 0 || parent2 < 0 ) break; // bad cluster, no match
502  //
503  // pull up jet's "global" index
504  //
505  int pidx = MatchingInput[parent1].user_index();
506  jetAssigned[pidx] = true;
507  matched = true;
508  break;
509  }
510  }
511  if ( !matched )
512  {
513  return UNMATCHED_PARTON;
514  }
515 
516  counter++;
517  }
518 
519  // Now save some info for DJR analysis (if requested).
520  // This mimics what is done in ME2pythia.f
521  // Basically, NJets and matching scale for these 4 cases:
522  // 1->0, 2->1, 3->2, and 4->3
523  //
524  if ( fDJROutFlag > 0 )
525  {
526  std::vector<double> MergingScale;
527  MergingScale.clear();
528  for ( int nj=0; nj<4; nj++ )
529  {
530  double dmscale2 = ClusSequence.exclusive_dmerge( nj );
531  double dmscale = sqrt( dmscale2 );
532  MergingScale.push_back( dmscale );
533  }
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;
540  fDJROutput.close();
541  }
542 
543  return NONE;
544 
545 }
int i
Definition: DBlmapReader.cc:9
fastjet::JetDefinition * fJetFinder
const HEPEUP * getHEPEUP() const
Definition: LHEEvent.h:43
T sqrt(T t)
Definition: SSEVec.h:48
std::vector< FiveVector > PUP
Definition: LesHouches.h:248
std::vector< fastjet::PseudoJet > fClusJets
std::vector< fastjet::PseudoJet > fPtSortedJets
std::vector< int > typeIdx[3]
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
static std::atomic< unsigned int > counter
volatile std::atomic< bool > shutdown_flag false
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
template<typename T >
T gen::JetMatchingMGFastJet::parseParameter ( const std::string &  value)
staticprivate

Definition at line 56 of file JetMatchingMGFastJet.cc.

References query::result, and contentValuesCheck::ss.

57 {
58  std::istringstream ss(value);
59  T result;
60  ss >> result;
61  return result;
62 }
tuple result
Definition: query.py:137
long double T
template<>
std::string gen::JetMatchingMGFastJet::parseParameter ( const std::string &  value)

Definition at line 65 of file JetMatchingMGFastJet.cc.

References query::result, and AlCaHLTBitMon_QueryRunRegistry::string.

66 {
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);
72  return result;
73 }
tuple result
Definition: query.py:137
template<>
bool gen::JetMatchingMGFastJet::parseParameter ( const std::string &  value_)

Definition at line 76 of file JetMatchingMGFastJet.cc.

References AlCaHLTBitMon_QueryRunRegistry::string, create_public_lumi_plots::transform, and relativeConstraints::value.

77 {
78  std::string value(value_);
79  std::transform(value.begin(), value.end(),
80  value.begin(), (int(*)(int))std::toupper);
81  return value == "T" || value == "Y" || value=="True" ||
82  value == "1" || value == ".TRUE.";
83 }
template<typename T >
void gen::JetMatchingMGFastJet::updateOrDie ( const std::map< std::string, std::string > &  params,
T param,
const std::string &  name 
)
staticprivate

Definition at line 137 of file JetMatchingMGFastJet.cc.

References Exception, and getParameter().

Referenced by init().

140 {
141  if (param < 0){
142  param = getParameter(params, name, param);
143  }
144  if (param < 0)
145  throw cms::Exception("Generator|PartonShowerVeto")
146  << "The MGParamCMS header does not specify the jet "
147  "matching parameter \"" << name << "\", but it "
148  "is requested by the CMSSW configuration."
149  << std::endl;
150 }
static T getParameter(const std::map< std::string, std::string > &params, const std::string &var, const T &defValue=T())

Member Data Documentation

double gen::JetMatchingMGFastJet::clFact
private

Definition at line 80 of file JetMatchingMGFastJet.h.

Referenced by initAfterBeams(), and match().

double gen::JetMatchingMGFastJet::coneRadius
private

Definition at line 91 of file JetMatchingMGFastJet.h.

Referenced by initAfterBeams().

bool gen::JetMatchingMGFastJet::doMerge
private

Definition at line 84 of file JetMatchingMGFastJet.h.

Referenced by initAfterBeams().

double gen::JetMatchingMGFastJet::etaJetMax
private

Definition at line 91 of file JetMatchingMGFastJet.h.

Referenced by initAfterBeams().

bool gen::JetMatchingMGFastJet::exclusive
private

Definition at line 103 of file JetMatchingMGFastJet.h.

Referenced by beforeHadronisation(), and JetMatchingMGFastJet().

std::vector<fastjet::PseudoJet> gen::JetMatchingMGFastJet::fClusJets
private

Definition at line 109 of file JetMatchingMGFastJet.h.

Referenced by initAfterBeams(), and match().

int gen::JetMatchingMGFastJet::fDJROutFlag
private

Definition at line 114 of file JetMatchingMGFastJet.h.

Referenced by JetMatchingMGFastJet(), and match().

std::ofstream gen::JetMatchingMGFastJet::fDJROutput
private

Definition at line 113 of file JetMatchingMGFastJet.h.

Referenced by match().

bool gen::JetMatchingMGFastJet::fExcLocal
private

Definition at line 96 of file JetMatchingMGFastJet.h.

Referenced by beforeHadronisation(), initAfterBeams(), and match().

bool gen::JetMatchingMGFastJet::fIsInit
private

Definition at line 116 of file JetMatchingMGFastJet.h.

Referenced by init(), and initAfterBeams().

fastjet::JetDefinition* gen::JetMatchingMGFastJet::fJetFinder
private

Definition at line 108 of file JetMatchingMGFastJet.h.

Referenced by initAfterBeams(), match(), and ~JetMatchingMGFastJet().

std::vector<fastjet::PseudoJet> gen::JetMatchingMGFastJet::fPtSortedJets
private

Definition at line 109 of file JetMatchingMGFastJet.h.

Referenced by initAfterBeams(), and match().

int gen::JetMatchingMGFastJet::jetAlgoPower
private

Definition at line 90 of file JetMatchingMGFastJet.h.

Referenced by initAfterBeams().

std::map<std::string, std::string> gen::JetMatchingMGFastJet::mgParams
private

Definition at line 72 of file JetMatchingMGFastJet.h.

Referenced by getParameter(), and init().

int gen::JetMatchingMGFastJet::nJetMax
private

Definition at line 87 of file JetMatchingMGFastJet.h.

Referenced by beforeHadronisation(), and initAfterBeams().

int gen::JetMatchingMGFastJet::nJetMin
private

Definition at line 87 of file JetMatchingMGFastJet.h.

Referenced by initAfterBeams().

int gen::JetMatchingMGFastJet::nQmatch
private

Definition at line 81 of file JetMatchingMGFastJet.h.

Referenced by beforeHadronisation(), and initAfterBeams().

double gen::JetMatchingMGFastJet::qCut
private

Definition at line 79 of file JetMatchingMGFastJet.h.

Referenced by initAfterBeams(), and match().

double gen::JetMatchingMGFastJet::qCutSq
private

Definition at line 79 of file JetMatchingMGFastJet.h.

Referenced by initAfterBeams(), and match().

bool gen::JetMatchingMGFastJet::runInitialized
private

Definition at line 101 of file JetMatchingMGFastJet.h.

Referenced by beforeHadronisation(), and init().

bool gen::JetMatchingMGFastJet::soup
private

Definition at line 102 of file JetMatchingMGFastJet.h.

Referenced by beforeHadronisation(), and JetMatchingMGFastJet().

std::vector< int > gen::JetMatchingMGFastJet::typeIdx[3]
private

Definition at line 99 of file JetMatchingMGFastJet.h.

Referenced by beforeHadronisation(), getPartonList(), and match().