cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.3
A generic central exclusive processes event generator
Loading...
Searching...
No Matches
Process.h
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#ifndef CepGen_Process_Process_h
20#define CepGen_Process_Process_h
21
22#include "CepGen/Event/Event.h"
27
28namespace cepgen {
30 namespace proc {
34 class Process : public NamedModule<Process> {
35 public:
36 explicit Process(const ParametersList&);
37 Process(const Process&);
38 virtual ~Process() = default;
39 Process& operator=(const Process&);
40
42
43 virtual std::unique_ptr<Process> clone() const;
44 virtual double computeWeight() = 0;
45
46 void clear();
47 void clearEvent();
48 void initialise();
49
50 const Kinematics& kinematics() const { return kin_; }
51 Kinematics& kinematics() { return kin_; }
52 void setKinematics();
53
54 // debugging utilities
55 double weight(const std::vector<double>&);
56 void dumpPoint(std::ostream* = nullptr) const;
57 void dumpVariables(std::ostream* = nullptr) const;
58
60 inline size_t ndim() const { return mapped_variables_.size(); }
61
62 bool hasEvent() const { return (bool)event_; }
63 const Event& event() const;
64 Event& event();
65 Event* eventPtr();
66
67 const Momentum& pA() const;
68 double mA2() const { return mA2_; }
69 double mA() const { return std::sqrt(mA2()); }
70 const Momentum& pB() const;
71 double mB2() const { return mB2_; }
72 double mB() const { return std::sqrt(mB2()); }
73 const Momentum& pX() const;
74 double mX2() const { return mX2_; }
75 double mX() const { return std::sqrt(mX2()); }
76 const Momentum& pY() const;
77 double mY2() const { return mY2_; }
78 double mY() const { return std::sqrt(mY2()); }
79
80 double x1() const { return x1_; }
81 double t1() const { return t1_; }
82 const Momentum& q1() const;
83 double x2() const { return x2_; }
84 double t2() const { return t2_; }
85 const Momentum& q2() const;
86
87 Momentum& q1();
88 Momentum& q2();
89
90 double wCM() const { return wcm_; }
91
92 Momentum& pA();
93 Momentum& pB();
94 Momentum& pX();
95 Momentum& pY();
96
97 double& mX2() { return mX2_; }
98 double& mY2() { return mY2_; }
99 double& t1() { return t1_; }
100 double& x1() { return x1_; }
101 double& t2() { return t2_; }
102 double& x2() { return x2_; }
103
104 Momentum& pc(size_t);
105 const Momentum& pc(size_t) const;
106
107 double s() const { return s_; }
108 double sqrtS() const { return sqs_; }
109 double inverseSqrtS() const { return inv_sqs_; }
110
112
113 protected:
114 virtual void addEventContent() = 0;
115 virtual void prepareKinematics() {}
116 virtual void fillKinematics() = 0;
117
118 //--- Mandelstam variables
119 double shat() const;
120
121 double mp_;
122 double mp2_;
123
124 public:
126 enum class Mapping {
127 linear = 0,
129 square,
130 power_law
131 };
133 friend std::ostream& operator<<(std::ostream&, const Mapping&);
141 Process& defineVariable(double& out,
142 const Mapping& type,
143 const Limits& lim,
144 const std::string& name,
145 const std::string& description = "");
147 double variableValue(size_t i, double x) const;
148
149 protected:
153 double generateVariables() const;
154
156 void setEventContent(const std::unordered_map<Particle::Role, spdgids_t>&);
157
158 double alphaEM(double q) const;
159 double alphaS(double q) const;
160
161 std::unique_ptr<utils::RandomGenerator> rnd_gen_;
162
163 inline const std::vector<double>& lastCoordinates() const { return point_coord_; }
164
165 private:
166 double s_{-1.};
167 double sqs_{-1.};
168 double inv_sqs_{-1.};
169 double wcm_{-1.};
170 double mA2_{-1.};
171 double mB2_{-1.};
172 double mX2_{-1.};
173 double mY2_{-1.};
174 double t1_{-1.};
175 double t2_{-1.};
176 double x1_{0.};
177 double x2_{0.};
178 std::unique_ptr<Coupling> alphaem_;
179 std::unique_ptr<Coupling> alphas_;
181 struct MappingVariable {
182 std::string name;
183 std::string description;
184 Limits limits;
185 double& value;
186 Mapping type;
187 size_t index;
188 };
190 std::vector<MappingVariable> mapped_variables_;
191 std::vector<double> point_coord_;
193 double base_jacobian_{1.};
194 Kinematics kin_{ParametersList()};
195 std::unique_ptr<Event> event_;
197 };
199 typedef std::unique_ptr<Process> ProcessPtr;
200 } // namespace proc
201} // namespace cepgen
202
203#endif
Container for the information on the in- and outgoing particles' kinematics.
Definition Event.h:28
List of kinematic constraints to apply on the process phase space.
Definition Kinematics.h:27
Validity interval for a variable.
Definition Limits.h:28
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
A description object for parameters collection.
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:264
void dumpPoint(std::ostream *=nullptr) const
Dump the coordinate of the phase-space point being evaluated.
Definition Process.cpp:364
virtual void prepareKinematics()
Compute the incoming state kinematics.
Definition Process.h:115
double wCM() const
Two-parton centre of mass energy.
Definition Process.h:90
double mB2() const
Negative-z incoming beam particle's squared mass.
Definition Process.h:71
const Kinematics & kinematics() const
Constant reference to the process kinematics.
Definition Process.h:50
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:287
const Momentum & pX() const
Positive-z outgoing beam particle's 4-momentum.
Definition Process.cpp:99
double t1() const
Positive-z incoming parton's squared mass.
Definition Process.h:81
double x2() const
Negative-z incoming parton's fractional momentum.
Definition Process.h:83
double mp_
Proton mass, in GeV/c .
Definition Process.h:121
double mX() const
Positive-z outgoing beam particle's mass.
Definition Process.h:75
const std::vector< double > & lastCoordinates() const
Last coordinates fed.
Definition Process.h:163
const Momentum & q2() const
Negative-z incoming parton's 4-momentum.
Definition Process.cpp:111
double inverseSqrtS() const
Inverse two-beam centre of mass energy.
Definition Process.h:109
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.
friend std::ostream & operator<<(std::ostream &, const Mapping &)
Human-friendly printout of the mapping type.
Definition Process.cpp:420
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:344
double & t1()
Positive-z incoming parton's squared mass.
Definition Process.h:99
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
Kinematics & kinematics()
Reference to the process kinematics.
Definition Process.h:51
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
virtual ~Process()=default
double s() const
Two-beam squared centre of mass energy.
Definition Process.h:107
double x1() const
Positive-z incoming parton's fractional momentum.
Definition Process.h:80
double & t2()
Negative-z incoming parton's squared mass.
Definition Process.h:101
double shat() const
Definition Process.cpp:127
double mA() const
Positive-z incoming beam particle's mass.
Definition Process.h:69
double mA2() const
Positive-z incoming beam particle's squared mass.
Definition Process.h:68
double mY2() const
Negative-z outgoing beam particle's squared mass.
Definition Process.h:77
double & x1()
Positive-z incoming parton's fractional momentum.
Definition Process.h:100
double sqrtS() const
Two-beam centre of mass energy.
Definition Process.h:108
double alphaEM(double q) const
Compute the electromagnetic running coupling algorithm at a given scale.
Definition Process.cpp:350
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:374
double mY() const
Negative-z outgoing beam particle's mass.
Definition Process.h:78
double mB() const
Negative-z incoming beam particle's mass.
Definition Process.h:72
bool hasEvent() const
Does the process contain (and hold) an event?
Definition Process.h:62
double & mX2()
Positive-z outgoing beam particle's squared mass.
Definition Process.h:97
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:236
Event * eventPtr()
Event pointer read/write accessor.
Definition Process.cpp:281
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:407
void clear()
Reset process prior to the phase space and variables definition.
Definition Process.cpp:129
double mX2() const
Positive-z outgoing beam particle's squared mass.
Definition Process.h:74
const Event & event() const
Handled particles objects and their relationships.
Definition Process.cpp:269
double & mY2()
Negative-z outgoing beam particle's squared mass.
Definition Process.h:98
double t2() const
Negative-z incoming parton's squared mass.
Definition Process.h:84
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:357
double generateVariables() const
Generate and initialise all variables handled by this process.
Definition Process.cpp:195
const Momentum & q1() const
Positive-z incoming parton's 4-momentum.
Definition Process.cpp:107
double & x2()
Negative-z incoming parton's fractional momentum.
Definition Process.h:102
double mp2_
Squared proton mass, in GeV /c .
Definition Process.h:122
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:190
A random number generator.
Kinematics
Type of scattering.
Definition Modes.h:28
std::unique_ptr< Process > ProcessPtr
Helper typedef for a Process unique pointer.
Definition Process.h:199
Common namespace for this Monte Carlo generator.