cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
Process.cpp
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 2013-2024 Laurent Forthomme
4 *
5 * This program is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19#include <iomanip>
20
26#include "CepGen/Physics/PDG.h"
28#include "CepGen/Utils/Math.h"
30#include "CepGen/Utils/String.h"
31
32namespace cepgen {
33 namespace proc {
34 auto compute_value = [](double in, const Process::Mapping& type) -> double {
35 switch (type) {
37 case Process::Mapping::power_law: //FIXME
38 default:
39 return in;
41 return in * in;
43 return std::exp(in);
44 }
45 };
46
48 : NamedModule(params),
49 mp_(PDG::get().mass(PDG::proton)),
50 mp2_(mp_ * mp_),
51 rnd_gen_(RandomGeneratorFactory::get().build(steer<ParametersList>("randomGenerator"))) {
52 if (const auto& kin = steer<ParametersList>("kinematics"); !kin.empty())
53 kin_.setParameters(kin);
54 if (steer<bool>("hasEvent"))
55 event_.reset(new Event);
56 }
57
58 Process::Process(const Process& proc) : Process(proc.parameters()) { *this = proc; }
59
61 s_ = proc.s_;
62 sqs_ = proc.sqs_;
63 inv_sqs_ = proc.inv_sqs_;
64 mA2_ = proc.mA2_;
65 mB2_ = proc.mB2_;
66 mX2_ = proc.mX2_;
67 mY2_ = proc.mY2_;
68 point_coord_ = proc.point_coord_;
69 base_jacobian_ = proc.base_jacobian_;
70 if (proc.event_)
71 event_.reset(new Event(*proc.event_));
72 CG_DEBUG("Process").log([&](auto& log) {
73 log << "Process " << name_ << " cloned with "
74 << utils::s("integration variable", mapped_variables_.size(), true) << ":";
75 for (const auto& var : mapped_variables_)
76 log << "\n\t" << var.index << ") " << var.description << " (type: " << var.type << ", limits: " << var.limits
77 << ").";
78 if (event_)
79 log << "\n\t" << *event_;
80 });
81 kin_ = proc.kin_;
82 return *this;
83 }
84
85 std::unique_ptr<Process> Process::clone() const {
86 throw CG_FATAL("Process:clone") << "Process \"" << name_ << "\" has no cloning method implementation!";
87 }
88
90
92
94
96
98
100
102
104
106
108
110
112
114 if (event()[Particle::CentralSystem].size() <= i)
115 throw CG_FATAL("Process:pc") << "Trying to retrieve central particle #" << i << " while only "
116 << event()[Particle::CentralSystem].size() << " is/are registered.";
117 return event()[Particle::CentralSystem].at(i).get().momentum();
118 }
119
120 const Momentum& Process::pc(size_t i) const {
121 if (event()(Particle::CentralSystem).size() <= i)
122 throw CG_FATAL("Process:pc") << "Trying to retrieve central particle #" << i << " while only "
123 << event()(Particle::CentralSystem).size() << " is/are registered.";
124 return event()(Particle::CentralSystem).at(i).momentum();
125 }
126
127 double Process::shat() const { return (q1() + q2()).mass2(); }
128
131 //--- initialise the "constant" (wrt x) part of the Jacobian
132 base_jacobian_ = 1.;
133 mapped_variables_.clear();
134 CG_DEBUG("Process:clear") << "Process event content, and integration variables cleared.";
135 }
136
137 void Process::dumpVariables(std::ostream* os) const {
138 std::ostringstream ss;
139 ss << "List of variables handled by this process:";
140 for (const auto& var : mapped_variables_)
141 ss << "\n\t(" << var.index << ") " << var.type << " mapping (" << var.description << ")"
142 << " in range " << var.limits;
143 if (os)
144 (*os) << ss.str();
145 else
146 CG_LOG << ss.str();
147 }
148
150 double& out, const Mapping& type, const Limits& lim, const std::string& name, const std::string& descr) {
151 if (lim.min() == lim.max()) {
152 if (lim.hasMin()) {
153 out = compute_value(lim.min(), type);
154 CG_DEBUG("Process:defineVariable")
155 << "Quantity " << descr << " is set to be constant with a value " << out << ".";
156 return *this;
157 } else
158 throw CG_FATAL("Process:defineVariable")
159 << "The limits for '" << descr << "' (" << lim << ") could not be retrieved from the user configuration.";
160 }
161
162 double jacob_weight = 1.; // initialise the local weight for this variable
163 switch (type) {
164 case Mapping::linear:
166 jacob_weight = lim.range();
167 break;
168 case Mapping::square:
169 jacob_weight = 2. * lim.range();
170 break;
172 jacob_weight = log(lim.max() / lim.min());
173 break;
174 }
175 const auto var_desc =
176 (!descr.empty() ? descr : (!name.empty() ? name : utils::format("var%z", mapped_variables_.size())));
177 mapped_variables_.emplace_back(MappingVariable{name, var_desc, lim, out, type, mapped_variables_.size()});
178 point_coord_.emplace_back(0.);
179 base_jacobian_ *= jacob_weight;
180 CG_DEBUG("Process:defineVariable") << "\n\t" << descr << " has been mapped to variable "
181 << mapped_variables_.size() << ".\n\t"
182 << "Allowed range for integration: " << lim << ".\n\t"
183 << "Variable integration mode: " << type << ".\n\t"
184 << "Weight in the Jacobian: " << jacob_weight << ".";
185 return *this;
186 }
187
188 double Process::variableValue(size_t i, double x) const {
189 const auto& var = mapped_variables_.at(i);
190 return compute_value(var.limits.x(x), var.type);
191 }
192
194 if (mapped_variables_.size() == 0)
195 throw CG_FATAL("Process:vars") << "No variables are mapped for this process!";
196 if (base_jacobian_ == 0.)
197 throw CG_FATAL("Process:vars") << "Point-independent component of the Jacobian for this "
198 << "process is null.\n\t"
199 << "Please check the validity of the phase space!";
200
201 double jacobian = 1.;
202 for (const auto& var : mapped_variables_) {
203 if (!var.limits.valid())
204 continue;
205 if (var.index >= point_coord_.size())
206 throw CG_FATAL("Process:x") << "Failed to retrieve coordinate " << var.index << " from "
207 << "a dimension-" << ndim() << " process!";
208 const auto& xv = point_coord_.at(var.index); // between 0 and 1
209 var.value = compute_value(var.limits.x(xv), var.type);
210 switch (var.type) {
211 case Mapping::linear: {
212 // jacobian *= 1
213 } break;
215 jacobian *= var.value;
216 } break;
217 case Mapping::square: {
218 jacobian *= var.limits.x(xv);
219 } break;
220 case Mapping::power_law: {
221 const double y = var.limits.max() / var.limits.min();
222 var.value = var.limits.min() * std::pow(y, xv);
223 jacobian *= var.value;
224 } break;
225 }
226 CG_DEBUG_LOOP("Process:vars") << "\n\tvariable " << var.index << std::left << std::setw(60)
227 << (!var.description.empty() ? " (" + var.description + ")" : "") << " in range "
228 << std::setw(20) << var.limits << " has value " << std::setw(20) << var.value
229 << " (x=" << point_coord_.at(var.index) << std::right << ")";
230 }
231 return jacobian;
232 }
233
234 double Process::weight(const std::vector<double>& x) {
235 point_coord_ = x;
236
237 //--- generate and initialise all variables and generate auxiliary
238 // (x-dependent) part of the Jacobian for this phase space point.
239 const auto aux_jacobian = generateVariables();
240
241 CG_DEBUG_LOOP("Process:weight").log([&](auto& log) {
242 log << "Jacobian: " << base_jacobian_ << " * " << aux_jacobian << " = " << (base_jacobian_ * aux_jacobian)
243 << ".\n\t";
244 dumpPoint(&log.stream());
245 });
246
247 if (!utils::positive(aux_jacobian))
248 return 0.;
249
250 //--- compute the integrand
251 const auto me_integrand = computeWeight();
252 CG_DEBUG_LOOP("Process:weight") << "Integrand = " << me_integrand << "\n\t"
253 << "Proc.-specific integrand * Jacobian (excl. global Jacobian) = "
254 << (me_integrand * aux_jacobian) << ".";
255 if (!utils::positive(me_integrand))
256 return 0.;
257
258 //--- combine every component into a single weight for this point
259 return (base_jacobian_ * aux_jacobian) * me_integrand * constants::GEVM2_TO_PB;
260 }
261
263 if (event_)
264 event_->restore();
265 }
266
267 const Event& Process::event() const {
268 if (!event_)
269 throw CG_FATAL("Process:event") << "Process does not have an event object!";
270 return *event_;
271 }
272
274 if (!event_)
275 throw CG_FATAL("Process:event") << "Process does not have an event object!";
276 return *event_;
277 }
278
280 if (!event_)
281 throw CG_FATAL("Process:event") << "Process does not have an event object!";
282 return event_.get();
283 }
284
286 CG_DEBUG("Process:initialise") << "Preparing to set the kinematics parameters. Input parameters: "
287 << ParametersDescription(kin_.parameters()) << ".";
288
289 clear(); // also resets the "first run" flag
290
291 // build the coupling objects
292 if (const auto& alpha_em = steer<ParametersList>("alphaEM"); !alpha_em.empty())
293 alphaem_ = AlphaEMFactory::get().build(alpha_em);
294 if (const auto& alpha_s = steer<ParametersList>("alphaS"); !alpha_s.empty())
295 alphas_ = AlphaSFactory::get().build(alpha_s);
296
297 const auto& p1 = kin_.incomingBeams().positive().momentum();
298 const auto& p2 = kin_.incomingBeams().negative().momentum();
299 //--- define incoming system
300 if (event_) {
301 auto& ib1 = event_->oneWithRole(Particle::IncomingBeam1);
302 ib1.setIntegerPdgId(kin_.incomingBeams().positive().integerPdgId());
303 ib1.setMomentum(p1);
304 auto& ib2 = event_->oneWithRole(Particle::IncomingBeam2);
305 ib2.setIntegerPdgId(kin_.incomingBeams().negative().integerPdgId());
306 ib2.setMomentum(p2);
307 auto& ob1 = event_->oneWithRole(Particle::OutgoingBeam1);
308 ob1.setIntegerPdgId(kin_.incomingBeams().positive().integerPdgId());
311 auto& ob2 = event_->oneWithRole(Particle::OutgoingBeam2);
312 ob2.setIntegerPdgId(kin_.incomingBeams().negative().integerPdgId());
315 for (auto& cp : (*event_)[Particle::CentralSystem])
316 cp.get().setIntegerPdgId(cp.get().integerPdgId());
317 }
318 s_ = kin_.incomingBeams().s();
319 sqs_ = std::sqrt(s_);
320 inv_sqs_ = 1. / sqs_;
321
322 mA2_ = mX2_ = p1.mass2();
323 mB2_ = mY2_ = p2.mass2();
324 wcm_ = 0.5 * (1. + std::sqrt(1. - 4. * std::sqrt(mA2_ * mB2_) / s_));
325
327
328 if (event_) {
329 CG_DEBUG("Process:initialise").log([this, &p1, &p2](auto& log) {
330 log << "Kinematics successfully set!\n"
331 << " sqrt(s) = " << sqs_ * 1.e-3 << " TeV,\n"
332 << " p1=" << p1 << ",\tmass=" << p1.mass() << " GeV\n"
333 << " p2=" << p2 << ",\tmass=" << p2.mass() << " GeV.\n";
334 dumpVariables(&log.stream());
335 });
336 clearEvent();
337 }
338 }
339
341 if (!rnd_gen_)
342 throw CG_FATAL("Process:randomGenerator") << "Process-local random generator was not yet initialised.";
343 return *rnd_gen_;
344 }
345
346 double Process::alphaEM(double q) const {
347 if (!alphaem_)
348 throw CG_FATAL("Process:alphaEM")
349 << "Trying to compute the electromagnetic running coupling while it is not initialised.";
350 return (*alphaem_)(q);
351 }
352
353 double Process::alphaS(double q) const {
354 if (!alphas_)
355 throw CG_FATAL("Process:alphaS")
356 << "Trying to compute the strong running coupling while it is not initialised.";
357 return (*alphas_)(q);
358 }
359
360 void Process::dumpPoint(std::ostream* os) const {
361 std::ostringstream oss;
362 oss << "Number of integration parameters: " << mapped_variables_.size() << ", point: {"
363 << utils::merge(point_coord_, ", ") << "}.";
364 if (!os)
365 CG_INFO("Process") << oss.str();
366 else
367 (*os) << oss.str();
368 }
369
370 void Process::setEventContent(const std::unordered_map<Particle::Role, spdgids_t>& part_ids) {
371 if (!event_)
372 return;
373 if (part_ids.count(Particle::Role::CentralSystem) == 0)
374 throw CG_FATAL("Process") << "The central system was not specified for this process.";
375
376 *event_ = Event::minimal(part_ids.at(Particle::Role::CentralSystem).size());
377 for (const auto& role_vs_parts : part_ids) {
378 auto evt_parts = (*event_)[role_vs_parts.first];
379 if (evt_parts.size() != role_vs_parts.second.size())
380 throw CG_FATAL("Process") << "Invalid number of '" << role_vs_parts.first << "' given. "
381 << "Expecting " << evt_parts.size() << ", got " << role_vs_parts.second.size()
382 << ".";
383 for (size_t i = 0; i < evt_parts.size(); ++i) {
384 auto& evt_part = evt_parts.at(i).get();
385 const auto user_evt_part_pdgid = role_vs_parts.second.at(i);
386 evt_part.setIntegerPdgId(user_evt_part_pdgid);
387 evt_part.momentum().setMass(PDG::get().mass(user_evt_part_pdgid));
388 }
389 }
390 event_->freeze(); // freeze the event as it is
391 }
392
395 if (event().hasRole(Particle::Intermediate)) {
396 Momentum interm_mom;
397 for (size_t i = 0; i < event()[Particle::CentralSystem].size(); ++i)
398 interm_mom += pc(i);
400 }
401 }
402
404 auto desc = ParametersDescription();
405 desc.add<ParametersDescription>("alphaEM", AlphaEMFactory::get().describeParameters("fixed"))
406 .setDescription("electromagnetic coupling evolution algorithm");
407 desc.add<ParametersDescription>("alphaS", AlphaSFactory::get().describeParameters("pegasus"))
408 .setDescription("strong coupling evolution algorithm");
409 desc.add<bool>("hasEvent", true).setDescription("does the process carry an event definition");
410 desc.add<ParametersDescription>("randomGenerator", RandomGeneratorFactory::get().describeParameters("stl"))
411 .setDescription("random number generator engine");
412 desc.add("kinematics", Kinematics::description());
413 return desc;
414 }
415
416 std::ostream& operator<<(std::ostream& os, const Process::Mapping& type) {
417 switch (type) {
419 return os << "linear";
421 return os << "exponential";
423 return os << "squared";
425 return os << "power law";
426 }
427 return os;
428 }
429 } // namespace proc
430} // namespace cepgen
#define CG_FATAL(mod)
Definition Exception.h:61
#define CG_DEBUG_LOOP(mod)
Definition Message.h:224
#define CG_LOG
Definition Message.h:212
#define CG_DEBUG(mod)
Definition Message.h:220
#define CG_INFO(mod)
Definition Message.h:216
spdgid_t integerPdgId() const
Beam particle PDG id Set the beam particle PDG id.
Definition Beam.h:47
const Momentum & momentum() const
Beam particle 4-momentum Set the beam particle 4-momentum.
Definition Beam.h:54
bool elastic() const
Does the beam remain on-shell after parton emission? Specify if the beam remains on-shell after parto...
Definition Beam.h:40
Container for the information on the in- and outgoing particles' kinematics.
Definition Event.h:28
Particle & oneWithRole(Particle::Role role)
First Particle object with a given role in the event.
Definition Event.cpp:190
size_t size() const
Number of particles in the event.
Definition Event.cpp:345
static Event minimal(size_t num_out_particles=1)
Build a trivial event with the minimal information.
Definition Event.cpp:52
const Beam & positive() const
Reference to the positive-z beam information.
double s() const
Incoming beams squared centre of mass energy (in GeV^2)
const Beam & negative() const
Reference to the negative-z beam information.
void setParameters(const ParametersList &) override
Set module parameters.
IncomingBeams & incomingBeams()
Beam/primary particle kinematics.
Definition Kinematics.h:36
static ParametersDescription description()
const ParametersList & parameters() const override
List containing all parameters handled.
Validity interval for a variable.
Definition Limits.h:28
double range() const
Full variable range allowed.
Definition Limits.cpp:65
bool hasMin() const
Have a lower limit?
Definition Limits.cpp:73
double min() const
Lower limit to apply on the variable.
Definition Limits.h:52
double max() const
Upper limit to apply on the variable.
Definition Limits.h:54
Container for a particle's 4-momentum, along with useful methods to ease the development of any matri...
Definition Momentum.h:33
Base runtime module object.
Definition NamedModule.h:28
const std::string & name() const
Module unique indexing name.
Definition NamedModule.h:42
const std::string name_
Module unique indexing name.
Definition NamedModule.h:45
A singleton holding all physics constants associated to particles.
Definition PDG.h:28
double mass(spdgid_t) const
Particle mass (in GeV)
Definition PDG.cpp:90
static PDG & get()
Retrieve a unique instance of this particles info collection.
Definition PDG.cpp:41
A description object for parameters collection.
Momentum & momentum()
Retrieve the momentum object associated with this particle Retrieve the momentum object associated wi...
Definition Particle.h:123
Particle & setMomentum(const Momentum &, bool offshell=false)
Associate a momentum object to this particle.
Definition Particle.cpp:77
@ FinalState
Stable, final state particle.
@ Unfragmented
Particle to be hadronised externally.
@ IncomingBeam2
incoming beam particle
Definition Particle.h:53
@ Parton2
beam incoming parton
Definition Particle.h:59
@ OutgoingBeam1
outgoing beam state/particle
Definition Particle.h:54
@ IncomingBeam1
incoming beam particle
Definition Particle.h:52
@ OutgoingBeam2
outgoing beam state/particle
Definition Particle.h:55
@ Parton1
beam incoming parton
Definition Particle.h:58
@ CentralSystem
Central particles system.
Definition Particle.h:56
@ Intermediate
Intermediate two-parton system.
Definition Particle.h:57
Class template to define any process to compute using this MC integrator/events generator.
Definition Process.h:34
void clearEvent()
Restore the event object to its initial state.
Definition Process.cpp:262
void dumpPoint(std::ostream *=nullptr) const
Dump the coordinate of the phase-space point being evaluated.
Definition Process.cpp:360
virtual void prepareKinematics()
Compute the incoming state kinematics.
Definition Process.h:115
const Momentum & pA() const
Positive-z incoming beam particle's 4-momentum.
Definition Process.cpp:91
void initialise()
Initialise the process once the kinematics has been set.
Definition Process.cpp:285
const Momentum & pX() const
Positive-z outgoing beam particle's 4-momentum.
Definition Process.cpp:99
const Momentum & q2() const
Negative-z incoming parton's 4-momentum.
Definition Process.cpp:111
virtual void addEventContent()=0
Set the incoming and outgoing state to be expected in the process.
virtual void fillKinematics()=0
Fill the Event object with the particles' kinematics.
Process & defineVariable(double &out, const Mapping &type, const Limits &lim, const std::string &name, const std::string &description="")
Register a variable to be handled and populated whenever a new phase space point weight is to be calc...
Definition Process.cpp:149
utils::RandomGenerator & randomGenerator() const
Accessor for this process' random number generator.
Definition Process.cpp:340
virtual std::unique_ptr< Process > clone() const
Copy all process attributes into a new object.
Definition Process.cpp:85
void dumpVariables(std::ostream *=nullptr) const
List all variables handled by this generic process.
Definition Process.cpp:137
Mapping
Type of mapping to apply on the variable.
Definition Process.h:126
@ exponential
an exponential mapping
@ power_law
a power-law mapping inherited from LPAIR
double shat() const
Definition Process.cpp:127
Process(const ParametersList &)
Definition Process.cpp:47
double alphaEM(double q) const
Compute the electromagnetic running coupling algorithm at a given scale.
Definition Process.cpp:346
const Momentum & pB() const
Negative-z incoming beam particle's 4-momentum.
Definition Process.cpp:95
void setEventContent(const std::unordered_map< Particle::Role, spdgids_t > &)
Set the incoming and outgoing states to be defined in this process (and prepare the Event object acco...
Definition Process.cpp:370
Momentum & pc(size_t)
Central particle's 4-momentum.
Definition Process.cpp:113
double weight(const std::vector< double > &)
Compute the weight for a phase-space point.
Definition Process.cpp:234
Event * eventPtr()
Event pointer read/write accessor.
Definition Process.cpp:279
const Momentum & pY() const
Negative-z outgoing beam particle's 4-momentum.
Definition Process.cpp:103
std::unique_ptr< utils::RandomGenerator > rnd_gen_
Process-local random number generator engine.
Definition Process.h:161
virtual double computeWeight()=0
Compute the phase space point weight.
static ParametersDescription description()
Definition Process.cpp:403
void clear()
Reset process prior to the phase space and variables definition.
Definition Process.cpp:129
const Event & event() const
Handled particles objects and their relationships.
Definition Process.cpp:267
size_t ndim() const
Number of dimensions on which the integration is performed.
Definition Process.h:60
double alphaS(double q) const
Compute the strong coupling algorithm at a given scale.
Definition Process.cpp:353
double generateVariables() const
Generate and initialise all variables handled by this process.
Definition Process.cpp:193
const Momentum & q1() const
Positive-z incoming parton's 4-momentum.
Definition Process.cpp:107
Process & operator=(const Process &)
Assignment operator.
Definition Process.cpp:60
double variableValue(size_t i, double x) const
Retrieve the physical value for one variable.
Definition Process.cpp:188
A random number generator.
constexpr double GEVM2_TO_PB
Conversion factor between GeV and barn i.e. in GeV .
Definition Constants.h:37
auto compute_value
Definition Process.cpp:34
std::ostream & operator<<(std::ostream &os, const Process::Mapping &type)
Definition Process.cpp:416
std::string format(const std::string &fmt, Args... args)
Format a string using a printf style format descriptor.
Definition String.h:61
std::string s(const std::string &word, float num, bool show_number)
Add a trailing "s" when needed.
Definition String.cpp:228
bool positive(const T &val)
Check if a number is positive and finite.
Definition Math.cpp:26
std::string merge(const std::vector< T > &vec, const std::string &delim)
Merge a collection of a printable type in a single string.
Definition String.cpp:248
Common namespace for this Monte Carlo generator.