CMS 3D CMS Logo

JetMatchingAlpgen.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <cstring>
3 #include <vector>
4 #include <memory>
5 #include <string>
6 
7 #include <HepMC/GenEvent.h>
8 
11 
15 
16 namespace gen {
17 
18  extern "C" {
19  // Put here all the functions for interfacing Fortran code.
20 
21  // extern void alinit_();
22  extern void alsetp_();
23  // extern void alevnt_();
24  extern void alveto_(int* ipveto);
25 
26  extern void dbpart_();
27  extern void pyupre_();
28 
29  void alshcd_(char csho[3]);
30  void alshen_();
31  } // extern "C"
32 
33  // Constructor
35  : JetMatching(params), applyMatching(params.getParameter<bool>("applyMatching")), runInitialized(false) {
36  ahopts_.etclus = params.getParameter<double>("etMin");
37  ahopts_.rclus = params.getParameter<double>("drMin");
38  ahopts_.iexc = params.getParameter<bool>("exclusive");
39  }
40 
41  // Destructor
43 
44  std::set<std::string> JetMatchingAlpgen::capabilities() const {
45  std::set<std::string> result;
46  result.insert("psFinalState");
47  result.insert("hepevt");
48  result.insert("pythia6");
49  return result;
50  }
51 
52  // Implements the Alpgen MLM method - use alpg_match.F
53 
55  // Read Alpgen run card stored in the LHERunInfo object.
56  std::vector<std::string> headerLines = runInfo->findHeader("AlpgenUnwParFile");
57  if (headerLines.empty())
58  throw cms::Exception("Generator|PartonShowerVeto") << "In order to use Alpgen jet matching, "
59  "the input file has to contain the corresponding "
60  "Alpgen headers."
61  << std::endl;
62 
63  // Parse the header using its bultin function.
64  header.parse(headerLines.begin(), headerLines.end());
65 
66  // I don't want to print this right now.
67  // std::cout << "Alpgen header" << std::endl;
68  // std::cout << "========================" << std::endl;
69  // std::cout << "\tihrd = " << header.ihrd << std::endl;
70  // std::cout << "\tmc = " << header.masses[AlpgenHeader::mc]
71  // << ", mb = " << header.masses[AlpgenHeader::mb]
72  // << ", mt = " << header.masses[AlpgenHeader::mt]
73  // << ", mw = " << header.masses[AlpgenHeader::mw]
74  // << ", mz = " << header.masses[AlpgenHeader::mz]
75  // << ", mh = " << header.masses[AlpgenHeader::mh]
76  // << std::endl;
77  // for(std::map<AlpgenHeader::Parameter, double>::const_iterator iter =
78  // header.params.begin(); iter != header.params.end(); ++iter)
79  // std::cout << "\t" << AlpgenHeader::parameterName(iter->first)
80  // << " = " << iter->second << std::endl;
81  // std::cout << "\txsec = " << header.xsec
82  // << " +-" << header.xsecErr << std::endl;
83  // std::cout << "========================" << std::endl;
84 
85  // Here we pass a few header variables to common block and
86  // call Alpgen init routine to do the rest.
87  // The variables passed first are the ones directly that
88  // need to be set up "manually": IHRD and the masses.
89  // (ebeam is set just to mimic the original code)
90  // Afterwards, we pass the full spectrum of Alpgen
91  // parameters directly into the AHPARS structure, to be
92  // treated by AHSPAR which is called inside alsetp_().
93 
97 
98  for (std::map<AlpgenHeader::Parameter, double>::const_iterator iter = header.params.begin();
99  iter != header.params.end();
100  ++iter) {
101  if (iter->first <= 0 || iter->first >= (int)AHPARS::nparam - 1)
102  continue;
103  ahpars_.parval[(int)iter->first - 1] = iter->second;
104  }
105 
106  // Run the rest of the setup.
107  alsetp_();
108 
109  // When we reach this point, the run is fully initialized.
110  runInitialized = true;
111  }
112 
114  // We can't continue if the run has not been initialized.
115  if (!runInitialized)
116  throw cms::Exception("Generator|PartonShowerVeto") << "Run not initialized in JetMatchingAlpgen" << std::endl;
117 
118  // We are called just after LHEInterface has filled in
119  // the Fortran common block (and Pythia6 called UPEVNT).
120 
121  // Possibly not interesting for us.
122  // (except perhaps for debugging?)
123  // pyupre_();
124  // dbpart_();
125  eventInitialized = true;
126  }
127 
128  /*
129 int JetMatchingAlpgen::match(const HepMC::GenEvent *partonLevel,
130  const HepMC::GenEvent *finalState,
131  bool showeredFinalState)
132 */
133  int JetMatchingAlpgen::match(const lhef::LHEEvent* partonLevel, const std::vector<fastjet::PseudoJet>* jetInput) {
134  /*
135  if (!showeredFinalState)
136  throw cms::Exception("Generator|PartonShowerVeto")
137  << "Alpgen matching expected parton shower "
138  "final state." << std::endl;
139 */
140 
141  if (!runInitialized)
142  throw cms::Exception("Generator|PartonShowerVeto") << "Run not initialized in JetMatchingAlpgen" << std::endl;
143 
144  if (!eventInitialized)
145  throw cms::Exception("Generator|PartonShowerVeto") << "Event not initialized in JetMatchingAlpgen" << std::endl;
146 
147  // If matching not required (e.g., icckw = 0), don't run the
148  // FORTRAN veto code.
149  if (!applyMatching)
150  return 0;
151 
152  // Call the Fortran veto code.
153  int veto = 0;
154  alveto_(&veto);
155 
156  eventInitialized = false;
157 
158  // If event was vetoed, the variable veto will contain the number 1.
159  // In this case, we must return 1 - that will be used as the return value from UPVETO.
160  // If event was accepted, the variable veto will contain the number 0.
161  // In this case, we must return 0 - that will be used as the return value from UPVETO.
162  return veto ? 1 : 0;
163  }
164 
165  void alshcd_(char csho[3]) {
166  std::memcpy(csho, "PYT", 3); // or "HER"
167  }
168 
169  void alshen_() {}
170 
171 } // end namespace gen
AHOPTS::etclus
double etclus
Definition: AlpgenCommonBlocks.h:8
AHPARS::nparam
static const unsigned int nparam
Definition: AlpgenCommonBlocks.h:31
electrons_cff.bool
bool
Definition: electrons_cff.py:372
MessageLogger.h
funct::false
false
Definition: Factorize.h:34
filterCSVwithJSON.copy
copy
Definition: filterCSVwithJSON.py:36
ahppara_
struct AHPPARA ahppara_
CalibrationSummaryClient_cfi.params
params
Definition: CalibrationSummaryClient_cfi.py:14
ahpars_
struct AHPARS ahpars_
AHPPARA::ihrd
int ihrd
Definition: AlpgenCommonBlocks.h:23
gen::JetMatchingAlpgen::capabilities
std::set< std::string > capabilities() const override
Definition: JetMatchingAlpgen.cc:44
gen::JetMatching
Definition: JetMatching.h:27
AlpgenHeader::ebeam
Definition: AlpgenHeader.h:13
AlpgenHeader::params
std::map< Parameter, double > params
Definition: AlpgenHeader.h:65
gen::JetMatchingAlpgen::init
void init(const lhef::LHERunInfo *runInfo) override
Definition: JetMatchingAlpgen.cc:54
AlpgenHeader::ihrd
unsigned int ihrd
Definition: AlpgenHeader.h:66
gen::JetMatchingAlpgen::runInitialized
bool runInitialized
Definition: JetMatchingAlpgen.h:29
gen::dbpart_
void dbpart_()
gen::JetMatchingAlpgen::applyMatching
bool applyMatching
Definition: JetMatchingAlpgen.h:28
gen::JetMatchingAlpgen::header
AlpgenHeader header
Definition: JetMatchingAlpgen.h:32
JetMatchingAlpgen.h
gen
Definition: PythiaDecays.h:13
gen::alshen_
void alshen_()
Definition: JetMatchingAlpgen.cc:169
AHOPTS::iexc
int iexc
Definition: AlpgenCommonBlocks.h:11
AHPPARA::masses
double masses[6]
Definition: AlpgenCommonBlocks.h:21
lhef::LHERunInfo
Definition: LHERunInfo.h:25
edm::ParameterSet
Definition: ParameterSet.h:36
AHPPARA::ebeam
double ebeam
Definition: AlpgenCommonBlocks.h:22
gen::JetMatchingAlpgen::~JetMatchingAlpgen
~JetMatchingAlpgen() override
Definition: JetMatchingAlpgen.cc:42
generator_cfi.applyMatching
applyMatching
Definition: generator_cfi.py:23
LHERunInfo.h
gen::alveto_
void alveto_(int *ipveto)
lhef::LHEEvent
Definition: LHEEvent.h:23
gen::JetMatchingAlpgen::beforeHadronisation
void beforeHadronisation(const lhef::LHEEvent *event) override
Definition: JetMatchingAlpgen.cc:113
createfilelist.int
int
Definition: createfilelist.py:10
gen::alsetp_
void alsetp_()
gen::alshcd_
void alshcd_(char csho[3])
Definition: JetMatchingAlpgen.cc:165
Exception
Definition: hltDiff.cc:246
AHPARS::parval
double parval[nparam]
Definition: AlpgenCommonBlocks.h:33
LHEEvent.h
gen::pyupre_
void pyupre_()
mps_fire.result
result
Definition: mps_fire.py:303
cms::Exception
Definition: Exception.h:70
ParameterSet.h
AHOPTS::rclus
double rclus
Definition: AlpgenCommonBlocks.h:9
AlpgenHeader::masses
double masses[MASS_MAX]
Definition: AlpgenHeader.h:71
gen::JetMatchingAlpgen::eventInitialized
bool eventInitialized
Definition: JetMatchingAlpgen.h:30
event
Definition: event.py:1
lhef::LHERunInfo::findHeader
std::vector< std::string > findHeader(const std::string &tag) const
Definition: LHERunInfo.cc:436
PbPb_ZMuSkimMuonDPG_cff.veto
veto
Definition: PbPb_ZMuSkimMuonDPG_cff.py:61
ahopts_
struct AHOPTS ahopts_
gen::JetMatchingAlpgen::match
int match(const lhef::LHEEvent *partonLevel, const std::vector< fastjet::PseudoJet > *jetInput) override
Definition: JetMatchingAlpgen.cc:133
AlpgenHeader::parse
bool parse(const std::vector< std::string >::const_iterator &begin, const std::vector< std::string >::const_iterator &end)
Definition: AlpgenHeader.cc:62
AlpgenHeader::MASS_MAX
Definition: AlpgenHeader.h:58
gen::JetMatchingAlpgen::JetMatchingAlpgen
JetMatchingAlpgen(const edm::ParameterSet &params)
Definition: JetMatchingAlpgen.cc:34