cepgen is hosted by Hepforge, IPPP Durham

Processes development

Contents

Processes development#

Two possible tracks are given to the CepGen process developer: the implementation of the full matrix element, starting from the beams kinematics definition down to the event weight and event content, giving a full flexibility to the whole process definition, or using a parton-from-beam-factorised simplified implementation.

The earlier relies on a sub-classing of the cepgen::Proc::Process object, already providing some useful kinematics variables of interest:

class Process : public NamedModule<Process>#

Class template to define any process to compute using this MC integrator/events generator.

Author

Laurent Forthomme laurent.forthomme@cern.ch

Date

Jan 2014

Subclassed by LPAIR, FactorisedProcess

Public Types

enum class Mapping#

Type of mapping to apply on the variable.

Values:

enumerator linear#

a linear \({\rm d}x\) mapping

enumerator exponential#

an exponential \(\frac{\dot{x}}{x} = \dot{\log x}\) mapping

enumerator square#

a square \({\rm d}x^2=2x\cdot\dot{x}\) mapping

enumerator power_law#

a power-law mapping inherited from LPAIR

Public Functions

Process(const Process&)#

Copy constructor for a user process.

Process &operator=(const Process&)#

Assignment operator.

virtual std::unique_ptr<Process> clone() const#

Copy all process attributes into a new object.

virtual double computeWeight() = 0#

Compute the phase space point weight.

void clear()#

Reset process prior to the phase space and variables definition.

void clearEvent()#

Restore the event object to its initial state.

void initialise()#

Initialise the process once the kinematics has been set.

inline const Kinematics &kinematics() const#

Constant reference to the process kinematics.

inline Kinematics &kinematics()#

Reference to the process kinematics.

double weight(const std::vector<double>&)#

Compute the weight for a phase-space point.

void dumpPoint(std::ostream* = nullptr) const#

Dump the coordinate of the phase-space point being evaluated.

void dumpVariables(std::ostream* = nullptr) const#

List all variables handled by this generic process.

inline size_t ndim() const#

Number of dimensions on which the integration is performed.

inline bool hasEvent() const#

Does the process contain (and hold) an event?

const Event &event() const#

Handled particles objects and their relationships.

Event &event()#

Event object read/write accessor.

Event *eventPtr()#

Event pointer read/write accessor.

const Momentum &pA() const#

Positive-z incoming beam particle’s 4-momentum.

inline double mA2() const#

Positive-z incoming beam particle’s squared mass.

inline double mA() const#

Positive-z incoming beam particle’s mass.

const Momentum &pB() const#

Negative-z incoming beam particle’s 4-momentum.

inline double mB2() const#

Negative-z incoming beam particle’s squared mass.

inline double mB() const#

Negative-z incoming beam particle’s mass.

const Momentum &pX() const#

Positive-z outgoing beam particle’s 4-momentum.

inline double mX2() const#

Positive-z outgoing beam particle’s squared mass.

inline double mX() const#

Positive-z outgoing beam particle’s mass.

const Momentum &pY() const#

Negative-z outgoing beam particle’s 4-momentum.

inline double mY2() const#

Negative-z outgoing beam particle’s squared mass.

inline double mY() const#

Negative-z outgoing beam particle’s mass.

inline double x1() const#

Positive-z incoming parton’s fractional momentum.

inline double t1() const#

Positive-z incoming parton’s squared mass.

const Momentum &q1() const#

Positive-z incoming parton’s 4-momentum.

inline double x2() const#

Negative-z incoming parton’s fractional momentum.

inline double t2() const#

Negative-z incoming parton’s squared mass.

const Momentum &q2() const#

Negative-z incoming parton’s 4-momentum.

Momentum &q1()#

Positive-z incoming parton’s 4-momentum.

Momentum &q2()#

Negative-z incoming parton’s 4-momentum.

inline double wCM() const#

Two-parton centre of mass energy.

Momentum &pA()#

Positive-z incoming beam particle’s 4-momentum.

Momentum &pB()#

Negative-z incoming beam particle’s 4-momentum.

Momentum &pX()#

Positive-z outgoing beam particle’s 4-momentum.

Momentum &pY()#

Negative-z outgoing beam particle’s 4-momentum.

inline double &mX2()#

Positive-z outgoing beam particle’s squared mass.

inline double &mY2()#

Negative-z outgoing beam particle’s squared mass.

inline double &t1()#

Positive-z incoming parton’s squared mass.

inline double &x1()#

Positive-z incoming parton’s fractional momentum.

inline double &t2()#

Negative-z incoming parton’s squared mass.

inline double &x2()#

Negative-z incoming parton’s fractional momentum.

Momentum &pc(size_t)#

Central particle’s 4-momentum.

const Momentum &pc(size_t) const#

Central particle’s 4-momentum.

inline double s() const#

Two-beam squared centre of mass energy.

inline double sqrtS() const#

Two-beam centre of mass energy.

inline double inverseSqrtS() const#

Inverse two-beam centre of mass energy.

utils::RandomGenerator &randomGenerator() const#

Accessor for this process’ random number generator.

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 calculated.

Note

To be run once per generation (before any point computation)

Parameters:
  • out[out] Reference to the variable to be mapped

  • type[in] Type of mapping to apply

  • lim[in] Integration limits

  • name[in] Computer-readable variable name

  • description[in] Human-readable description of the variable

double variableValue(size_t i, double x) const#

Retrieve the physical value for one variable.

inline const std::string &name() const#

Module unique indexing name.

inline bool operator==(const SteeredObject &oth) const#

Equality operator.

inline bool operator!=(const SteeredObject &oth) const#

Inequality operator.

inline virtual const ParametersList &parameters() const override#

Module user-defined parameters.

inline virtual void setParameters(const ParametersList &params) override#

Set module parameters.

inline void setDescribedParameters(const ParametersList &params_orig)#

Set (documented) module parameters.

Protected Functions

virtual void addEventContent() = 0#

Set the incoming and outgoing state to be expected in the process.

inline virtual void prepareKinematics()#

Compute the incoming state kinematics.

virtual void fillKinematics() = 0#

Fill the Event object with the particles’ kinematics.

double shat() const#

\(\hat s=(p_1+p_2)^2=(p_3+...)^2\)

double generateVariables() const#

Generate and initialise all variables handled by this process.

Note

To be run at each point computation (therefore, to be optimised!)

Returns:

Phase space point-dependent component of the Jacobian weight of the point in the phase space for integration

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 accordingly)

double alphaEM(double q) const#

Compute the electromagnetic running coupling algorithm at a given scale.

double alphaS(double q) const#

Compute the strong coupling algorithm at a given scale.

inline const std::vector<double> &lastCoordinates() const#

Last coordinates fed.

template<typename T>
inline T steer(const std::string &key) const#

Retrieve a parameters as previously steered.

template<typename T, typename U>
inline U steerAs(const std::string &key) const#

Retrieve a recasted parameters as previously steered.

inline std::string steerName() const#

Retrieve module name from parameters.

std::string steerPath(const std::string &key) const#

Retrieve a path from common search paths.

Protected Attributes

double mp_#

Proton mass, in GeV/c \(^2\).

double mp2_#

Squared proton mass, in GeV \(^2\)/c \(^4\).

std::unique_ptr<utils::RandomGenerator> rnd_gen_#

Process-local random number generator engine.

const std::string name_#

Module unique indexing name.

mutable ParametersList params_#

Module parameters.

Friends

friend std::ostream &operator<<(std::ostream&, const Mapping&)#

Human-friendly printout of the mapping type.

As noted above, to avoid code redundancy, the following helper accessors can be used to retrieve and set the various event content kinematics at the level of the process:

However, in most cases, such an access to the full beam dynamics and kinematics content is not necessarily required. Therefore, the so-called factorised processes definition can be useful for most of user cases.

Factorised processes#

New in version 0.9.

We recommend the implementation of all new processes in C++, as described below. However, we conveniently provide an interface to Fortran implementations for historical reasons. A general description of this latter can be seen down this page.

C++ interface#

The cepgen::proc::FactorisedProcess helper derivated-class of the earlier is introduced to allow the parton emission part to be transparent to the process developper.

class FactorisedProcess : public Process#

Subclassed by MadGraphProcessBuilder, PPtoFF, PPtoWW, FortranFactorisedProcess

class FactorisedProcess : public Process

Generic parton emission-factorised process.

Author

Laurent Forthomme laurent.forthomme@cern.ch

Date

Jul 2023

Note

0 to 2 dimensions may be used for the scattered diffractive system(s)’ invariant mass definition.

Subclassed by MadGraphProcessBuilder, PPtoFF, PPtoWW, FortranFactorisedProcess

Public Types

enum class Mapping#

Type of mapping to apply on the variable.

Values:

enumerator linear#

a linear \({\rm d}x\) mapping

enumerator exponential#

an exponential \(\frac{\dot{x}}{x} = \dot{\log x}\) mapping

enumerator square#

a square \({\rm d}x^2=2x\cdot\dot{x}\) mapping

enumerator power_law#

a power-law mapping inherited from LPAIR

Public Functions

explicit FactorisedProcess(const ParametersList &params, const spdgids_t &output)#

Class constructor.

Parameters:
  • params[in] Parameters list

  • output[in] Produced final state particles

virtual double computeWeight() override#

Compute the phase space point weight.

virtual void fillKinematics() final override#

Fill the Event object with the particles’ kinematics.

virtual std::unique_ptr<Process> clone() const#

Copy all process attributes into a new object.

void clear()#

Reset process prior to the phase space and variables definition.

void clearEvent()#

Restore the event object to its initial state.

void initialise()#

Initialise the process once the kinematics has been set.

inline const Kinematics &kinematics() const#

Constant reference to the process kinematics.

inline Kinematics &kinematics()#

Reference to the process kinematics.

double weight(const std::vector<double>&)#

Compute the weight for a phase-space point.

void dumpPoint(std::ostream* = nullptr) const#

Dump the coordinate of the phase-space point being evaluated.

void dumpVariables(std::ostream* = nullptr) const#

List all variables handled by this generic process.

inline size_t ndim() const#

Number of dimensions on which the integration is performed.

inline bool hasEvent() const#

Does the process contain (and hold) an event?

const Event &event() const#

Handled particles objects and their relationships.

Event &event()#

Event object read/write accessor.

Event *eventPtr()#

Event pointer read/write accessor.

const Momentum &pA() const#

Positive-z incoming beam particle’s 4-momentum.

Momentum &pA()#

Positive-z incoming beam particle’s 4-momentum.

inline double mA2() const#

Positive-z incoming beam particle’s squared mass.

inline double mA() const#

Positive-z incoming beam particle’s mass.

const Momentum &pB() const#

Negative-z incoming beam particle’s 4-momentum.

Momentum &pB()#

Negative-z incoming beam particle’s 4-momentum.

inline double mB2() const#

Negative-z incoming beam particle’s squared mass.

inline double mB() const#

Negative-z incoming beam particle’s mass.

const Momentum &pX() const#

Positive-z outgoing beam particle’s 4-momentum.

Momentum &pX()#

Positive-z outgoing beam particle’s 4-momentum.

inline double mX2() const#

Positive-z outgoing beam particle’s squared mass.

inline double &mX2()#

Positive-z outgoing beam particle’s squared mass.

inline double mX() const#

Positive-z outgoing beam particle’s mass.

const Momentum &pY() const#

Negative-z outgoing beam particle’s 4-momentum.

Momentum &pY()#

Negative-z outgoing beam particle’s 4-momentum.

inline double mY2() const#

Negative-z outgoing beam particle’s squared mass.

inline double &mY2()#

Negative-z outgoing beam particle’s squared mass.

inline double mY() const#

Negative-z outgoing beam particle’s mass.

inline double x1() const#

Positive-z incoming parton’s fractional momentum.

inline double &x1()#

Positive-z incoming parton’s fractional momentum.

inline double t1() const#

Positive-z incoming parton’s squared mass.

inline double &t1()#

Positive-z incoming parton’s squared mass.

const Momentum &q1() const#

Positive-z incoming parton’s 4-momentum.

Momentum &q1()#

Positive-z incoming parton’s 4-momentum.

inline double x2() const#

Negative-z incoming parton’s fractional momentum.

inline double &x2()#

Negative-z incoming parton’s fractional momentum.

inline double t2() const#

Negative-z incoming parton’s squared mass.

inline double &t2()#

Negative-z incoming parton’s squared mass.

const Momentum &q2() const#

Negative-z incoming parton’s 4-momentum.

Momentum &q2()#

Negative-z incoming parton’s 4-momentum.

inline double wCM() const#

Two-parton centre of mass energy.

Momentum &pc(size_t)#

Central particle’s 4-momentum.

const Momentum &pc(size_t) const#

Central particle’s 4-momentum.

inline double s() const#

Two-beam squared centre of mass energy.

inline double sqrtS() const#

Two-beam centre of mass energy.

inline double inverseSqrtS() const#

Inverse two-beam centre of mass energy.

utils::RandomGenerator &randomGenerator() const#

Accessor for this process’ random number generator.

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 calculated.

Note

To be run once per generation (before any point computation)

Parameters:
  • out[out] Reference to the variable to be mapped

  • type[in] Type of mapping to apply

  • lim[in] Integration limits

  • name[in] Computer-readable variable name

  • description[in] Human-readable description of the variable

double variableValue(size_t i, double x) const#

Retrieve the physical value for one variable.

inline const std::string &name() const#

Module unique indexing name.

inline bool operator==(const SteeredObject &oth) const#

Equality operator.

inline bool operator!=(const SteeredObject &oth) const#

Inequality operator.

inline virtual const ParametersList &parameters() const override#

Module user-defined parameters.

inline virtual void setParameters(const ParametersList &params) override#

Set module parameters.

inline void setDescribedParameters(const ParametersList &params_orig)#

Set (documented) module parameters.

Protected Functions

virtual void addEventContent() override#

Set the incoming and outgoing state to be expected in the process.

virtual void prepareKinematics() final override#

Compute the incoming state kinematics.

virtual void prepareFactorisedPhaseSpace() = 0#

Prepare central part of the Jacobian after kinematics is set.

virtual double computeFactorisedMatrixElement() = 0#

Factorised matrix element (event weight)

double that() const#

\(\hat t=\frac{1}{2}\left[(p_1-p_3)^2+(p_2-p_4)^2\right]\)

double uhat() const#

\(\hat u=\frac{1}{2}\left[(p_1-p_4)^2+(p_2-p_3)^2\right]\)

double shat() const#

\(\hat s=(p_1+p_2)^2=(p_3+...)^2\)

double generateVariables() const#

Generate and initialise all variables handled by this process.

Note

To be run at each point computation (therefore, to be optimised!)

Returns:

Phase space point-dependent component of the Jacobian weight of the point in the phase space for integration

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 accordingly)

double alphaEM(double q) const#

Compute the electromagnetic running coupling algorithm at a given scale.

double alphaS(double q) const#

Compute the strong coupling algorithm at a given scale.

inline const std::vector<double> &lastCoordinates() const#

Last coordinates fed.

template<typename T>
inline T steer(const std::string &key) const#

Retrieve a parameters as previously steered.

template<typename T, typename U>
inline U steerAs(const std::string &key) const#

Retrieve a recasted parameters as previously steered.

inline std::string steerName() const#

Retrieve module name from parameters.

std::string steerPath(const std::string &key) const#

Retrieve a path from common search paths.

Protected Attributes

const std::unique_ptr<PhaseSpaceGenerator> psgen_#

Kinematic variables generator for the phase space coverage.

double mp_#

Proton mass, in GeV/c \(^2\).

double mp2_#

Squared proton mass, in GeV \(^2\)/c \(^4\).

std::unique_ptr<utils::RandomGenerator> rnd_gen_#

Process-local random number generator engine.

const std::string name_#

Module unique indexing name.

mutable ParametersList params_#

Module parameters.

As this object is pure virtual, the two following members have to be defined for any factorised process:

The two-parton (or more generally the full central kinematics) can be generated automatically by CepGen according to the two scenarios enumerated above: collinear, or transverse momentum-dependent parton emission.

An experimental interfacing to other phase space mappers (such as Rambo) was recently made technically possible, and is currently under development.

The basic object of interest for this mapping is:

class PhaseSpaceGenerator : public NamedModule<PhaseSpaceGenerator>#

Class template to define any phase space helper process.

Author

Laurent Forthomme laurent.forthomme@cern.ch

Date

Feb 2024

Subclassed by PhaseSpaceGenerator2to4< T >

Public Functions

inline virtual void setCentralCuts(const cuts::Central&)#

Set cuts on central particles.

virtual void initialise(proc::FactorisedProcess*) = 0#

Set all process parameters.

virtual bool generate() = 0#

Generate a kinematics combination, and return a success flag.

virtual double weight() const = 0#

Return the event weight for a kinematics combination.

virtual pdgids_t partons() const = 0#

List of incoming partons in kinematics.

virtual void setCentral(const std::vector<int>&) = 0#

Override the central particles list.

virtual std::vector<int> central() const = 0#

List of outgoing central particles in kinematics.

inline const std::string &name() const#

Module unique indexing name.

inline bool operator==(const SteeredObject &oth) const#

Equality operator.

inline bool operator!=(const SteeredObject &oth) const#

Inequality operator.

inline virtual const ParametersList &parameters() const override#

Module user-defined parameters.

inline virtual void setParameters(const ParametersList &params) override#

Set module parameters.

inline void setDescribedParameters(const ParametersList &params_orig)#

Set (documented) module parameters.

Public Static Functions

static inline ParametersDescription description()#

Describe all steering parameters for this module.

Protected Functions

template<typename T>
inline T steer(const std::string &key) const#

Retrieve a parameters as previously steered.

template<typename T, typename U>
inline U steerAs(const std::string &key) const#

Retrieve a recasted parameters as previously steered.

inline std::string steerName() const#

Retrieve module name from parameters.

std::string steerPath(const std::string &key) const#

Retrieve a path from common search paths.

Protected Attributes

const std::string name_#

Module unique indexing name.

mutable ParametersList params_#

Module parameters.

Which currently contains one instantiation, for 2-to-4 processes (such as pptoff and pptoww):

template<typename T>
class PhaseSpaceGenerator2to4 : public PhaseSpaceGenerator#

A 2-to-4 (or 2-to-2 central) phase space generator.

Protected Functions

template<typename T>
inline T steer(const std::string &key) const#

Retrieve a parameters as previously steered.

template<typename T, typename U>
inline U steerAs(const std::string &key) const#

Retrieve a recasted parameters as previously steered.

inline std::string steerName() const#

Retrieve module name from parameters.

std::string steerPath(const std::string &key) const#

Retrieve a path from common search paths.

Protected Attributes

const std::string name_#

Module unique indexing name.

mutable ParametersList params_#

Module parameters.

Again, the same accessors as for cepgen::proc::Process can be used for the retrieval of the event kinematics used in computing the event weight.

Fortran interface#

For development or testing purposes (and increased flexibility in general), Fortran definitions of physics processes can be fed and used by CepGen for cross section computation and event generation. In this page a summary and hands-on example of such a Fortran implementation and linkage is described.

Before other things, you will have to provide the definition of your process and its topology for CepGen to handle it properly. This requires you to fill in a set of predefined common blocks shared between your process definition and the core CepGen instance.

All these new processes are then linked to CepGen through the construction and association of each matrix element definition and interface to the cepgen::proc::FortranFactorisedProcess object.

class FortranFactorisedProcess : public FactorisedProcess#

Compute the matrix element for a generic factorised process defined in a Fortran weighting function.

Public Types

enum class Mapping#

Type of mapping to apply on the variable.

Values:

enumerator linear#

a linear \({\rm d}x\) mapping

enumerator exponential#

an exponential \(\frac{\dot{x}}{x} = \dot{\log x}\) mapping

enumerator square#

a square \({\rm d}x^2=2x\cdot\dot{x}\) mapping

enumerator power_law#

a power-law mapping inherited from LPAIR

Public Functions

explicit FortranFactorisedProcess(const ParametersList&, const std::function<double(void)> &func)#

Construct a Fortran-CepGen interface object using a double precision argument-less F77 function.

Parameters:

func[in] a double precision argument-less Fortran function returning the event weight

inline virtual ProcessPtr clone() const override#

Copy all process attributes into a new object.

virtual double computeWeight() override#

Compute the phase space point weight.

virtual void fillKinematics() final override#

Fill the Event object with the particles’ kinematics.

void clear()#

Reset process prior to the phase space and variables definition.

void clearEvent()#

Restore the event object to its initial state.

void initialise()#

Initialise the process once the kinematics has been set.

inline const Kinematics &kinematics() const#

Constant reference to the process kinematics.

inline Kinematics &kinematics()#

Reference to the process kinematics.

double weight(const std::vector<double>&)#

Compute the weight for a phase-space point.

void dumpPoint(std::ostream* = nullptr) const#

Dump the coordinate of the phase-space point being evaluated.

void dumpVariables(std::ostream* = nullptr) const#

List all variables handled by this generic process.

inline size_t ndim() const#

Number of dimensions on which the integration is performed.

inline bool hasEvent() const#

Does the process contain (and hold) an event?

const Event &event() const#

Handled particles objects and their relationships.

Event &event()#

Event object read/write accessor.

Event *eventPtr()#

Event pointer read/write accessor.

const Momentum &pA() const#

Positive-z incoming beam particle’s 4-momentum.

Momentum &pA()#

Positive-z incoming beam particle’s 4-momentum.

inline double mA2() const#

Positive-z incoming beam particle’s squared mass.

inline double mA() const#

Positive-z incoming beam particle’s mass.

const Momentum &pB() const#

Negative-z incoming beam particle’s 4-momentum.

Momentum &pB()#

Negative-z incoming beam particle’s 4-momentum.

inline double mB2() const#

Negative-z incoming beam particle’s squared mass.

inline double mB() const#

Negative-z incoming beam particle’s mass.

const Momentum &pX() const#

Positive-z outgoing beam particle’s 4-momentum.

Momentum &pX()#

Positive-z outgoing beam particle’s 4-momentum.

inline double mX2() const#

Positive-z outgoing beam particle’s squared mass.

inline double &mX2()#

Positive-z outgoing beam particle’s squared mass.

inline double mX() const#

Positive-z outgoing beam particle’s mass.

const Momentum &pY() const#

Negative-z outgoing beam particle’s 4-momentum.

Momentum &pY()#

Negative-z outgoing beam particle’s 4-momentum.

inline double mY2() const#

Negative-z outgoing beam particle’s squared mass.

inline double &mY2()#

Negative-z outgoing beam particle’s squared mass.

inline double mY() const#

Negative-z outgoing beam particle’s mass.

inline double x1() const#

Positive-z incoming parton’s fractional momentum.

inline double &x1()#

Positive-z incoming parton’s fractional momentum.

inline double t1() const#

Positive-z incoming parton’s squared mass.

inline double &t1()#

Positive-z incoming parton’s squared mass.

const Momentum &q1() const#

Positive-z incoming parton’s 4-momentum.

Momentum &q1()#

Positive-z incoming parton’s 4-momentum.

inline double x2() const#

Negative-z incoming parton’s fractional momentum.

inline double &x2()#

Negative-z incoming parton’s fractional momentum.

inline double t2() const#

Negative-z incoming parton’s squared mass.

inline double &t2()#

Negative-z incoming parton’s squared mass.

const Momentum &q2() const#

Negative-z incoming parton’s 4-momentum.

Momentum &q2()#

Negative-z incoming parton’s 4-momentum.

inline double wCM() const#

Two-parton centre of mass energy.

Momentum &pc(size_t)#

Central particle’s 4-momentum.

const Momentum &pc(size_t) const#

Central particle’s 4-momentum.

inline double s() const#

Two-beam squared centre of mass energy.

inline double sqrtS() const#

Two-beam centre of mass energy.

inline double inverseSqrtS() const#

Inverse two-beam centre of mass energy.

utils::RandomGenerator &randomGenerator() const#

Accessor for this process’ random number generator.

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 calculated.

Note

To be run once per generation (before any point computation)

Parameters:
  • out[out] Reference to the variable to be mapped

  • type[in] Type of mapping to apply

  • lim[in] Integration limits

  • name[in] Computer-readable variable name

  • description[in] Human-readable description of the variable

double variableValue(size_t i, double x) const#

Retrieve the physical value for one variable.

inline const std::string &name() const#

Module unique indexing name.

inline bool operator==(const SteeredObject &oth) const#

Equality operator.

inline bool operator!=(const SteeredObject &oth) const#

Inequality operator.

inline virtual const ParametersList &parameters() const override#

Module user-defined parameters.

inline virtual void setParameters(const ParametersList &params) override#

Set module parameters.

inline void setDescribedParameters(const ParametersList &params_orig)#

Set (documented) module parameters.

Public Static Attributes

static ParametersList kProcParameters#

List of parameters to steer the process.

In a later paragraph, this linking recipe will be described.

Output event kinematics#

This common block lets CepGen access the whole event structure and kinematics information. Beside the outgoing beam particles 4-momenta px and py, it contains the PDG id and 4-momentum of all central system particles in pc.

For this latter, a maximum multiplicity of 10 particles is handled by default in CepGen.

common/evtkin/nout,ipdg,idum2,pc,px,py
integer nout,idum2,ipdg(10)
double precision pc(4,10),px(4),py(4)

Should this maximal value not be sufficient for your implementation purposes, please get in touch with us.

Process definition#

This block specifies all useful run-level parameters, such as:

  • icontri, the kinematics mode considered (1 = elastic, 2-3 = single-dissociative, 4 = double-dissociative),

  • iflux1, iflux2, the type of \(\kt\)-dependent flux used to parameterise the incoming partons’ kinematics,

  • sfmod, the structure functions modelling used [1],

  • inp1, inp2, the longitudinal beam momenta (in GeV/c).

Additionally, for heavy ions initial states, the following parameters are used:

  • a_nuc1, a_nuc2, the mass numbers of positive- and negative-z incoming beam particles,

  • z_nuc1, z_nuc2, the atomic numbers of positive- and negative-z incoming beam particles,

For increased modularity, a set of process parameters directly parsed from the input cards (see e.g. the processParameters block) may also be introduced. At the Fortran process definition level, all parameters may be accessed through a set of interfacing functions:

int cepgen_param_int_(char *pname, int &def)#

Retrieve an integer process parameter from runtime parameters collection.

Parameters:
  • pname[in] Parameter name string

  • def[in] Default parameter value if not found in collection

double cepgen_param_real_(char *pname, double &def)#

Retrieve a double precision floating point process parameter from runtime parameters collection.

Parameters:
  • pname[in] Parameter name string

  • def[in] Default parameter value if not found in collection

Note

These latter may be translated into the following Fortran signature:

external CepGen_param_int,CepGen_param_real
integer CepGen_param_int
double precision CepGen_param_real

Those two can be called from Fortran using e.g.

integer ival
real*8 rval
! [...]
ival = CepGen_param_int('int_parameter', 0)      ! first argument is parameter name
rval = CepGen_param_real('real_parameter', 0.d0) ! second argument is default value

\(k_{\rm T}\)-factorisation kinematics#

This common block holds the following \(k_{\rm T}\)-factorisation processes kinematics variables (generated and filled by the core CepGen instance):

  • q1t and q2t, the 2-norm of transverse incoming parton virtualities,

  • phiq1t and phiq2t, the equivalent azimutal angle \(\phi\) of these virtualities,

  • y1 and y2 the outgoing system particles’ rapidities,

  • ptdiff and phiptdiff the 2-norm and azimutal angle of the outgoing system’s transverse momentum difference,

  • am_x and am_y the dissociative outgoing beams’ invariant masses (or \(m_p\) for elastic parton emission from proton).

   common/ktkin/q1t,q2t,phiq1t,phiq2t,y1,y2,ptdiff,phiptdiff,
&     am_x,am_y
   double precision q1t,q2t,phiq1t,phiq2t,y1,y2,ptdiff,phiptdiff,
&     am_x,am_y

Phase space cuts#

A limited set of pre-registered phase space cuts (one boolean switch, a lower and an upper value to apply) may be found in this interfacing library, such as:

  • single particles \(\pt\), energy, pseudo-rapidity \(\eta\),

  • central system invariant mass, \(\pt\),

  • central system correlations: \(\Delta y\).

   common/constants/am_p,units,pi
   common/kincuts/ipt,iene,ieta,iinvm,iptsum,idely,
&     pt_min,pt_max,ene_min,ene_max,eta_min,eta_max,
&     invm_min,invm_max,ptsum_min,ptsum_max,
&     dely_min,dely_max
   logical ipt,iene,ieta,iinvm,iptsum,idely
   double precision pt_min,pt_max,ene_min,ene_max,eta_min,eta_max,
&     invm_min,invm_max,ptsum_min,ptsum_max,
&     dely_min,dely_max

Should this collection not be sufficient for your purposes, please contact us.

Physics constants#

This block introduces all useful and standard physical constants in double precision to help the process definition: \(m_p\), GeV\({}^2\)/barn, or \(pi\).

common/constants/am_p,units,pi
double precision am_p,units,pi

Matrix element definition#

The core definition of your process is to be implemented within a new .f file you will store in the External/Processes/ directory of your local CepGen install.

Currently, the only following Fortran process is defined: nucl_to_ff

For illustration in this documentation page, we will be using a dummy process name, such as dummy_f77_process. For the sake of simplicity in your processes definition bookkeeping we suggest you to match the filename to the matrix element weighting function name. In our example, this corresponds to External/Processes/dummy_f77_process.f

We expect your weighting function to follow the C-equivalent signature:

double dummy_f77_process_(void);

This might get translated, for instance, into the following minimal working example (here using the free form Fortran coding convention):

 1      function dummy_f77_process()
 2      implicit none
 3      double precision dummy_f77_process
 4      !--------------------------------------------------------------------------
 5      ! CepGen overhead
 6      !--------------------------------------------------------------------------
 7      include 'KTBlocks.inc' ! mandatory, include the kinematics common blocks
 8      call CepGen_print      ! optional, display some run parameters information
 9      !--------------------------------------------------------------------------
10      ! end of overhead, beginning of process definition
11      !--------------------------------------------------------------------------
12      dummy_f77_process = 1.D0 ! placeholder, your actual definition is to be
13                               ! implemented here
14      !--------------------------------------------------------------------------
15      ! end of process definition
16      !--------------------------------------------------------------------------
17      return
18      end

With the kinematics common blocks defined through the include statement described above.

Helper tools#

To ease the work of the process developer, we provide utilitaries for the evaluation of unintegrated parton fluxes and other physics quantities through the CepGen core C++ library methods. The external methods exposed to the Fortran process through the include statement above are:

double cepgen_kt_flux_(int &fmode, double &x, double &kt2, int &sfmode, double &min, double &mout)#

Compute a \(k_{\rm T}\)-dependent flux for single nucleons.

Parameters:
  • fmode[in] Flux mode

  • x[in] Fractional momentum loss

  • kt2[in] The \(k_{\rm T}\) transverse momentum norm

  • sfmode[in] Structure functions set for dissociative emission

  • min[in] Incoming particle mass

  • mout[in] Diffractive state mass for dissociative emission

double cepgen_kt_flux_hi_(int &fmode, double &x, double &kt2, int &a, int &z)#

Compute a \(k_{\rm T}\)-dependent flux for heavy ions.

Parameters:
  • fmode[in] Flux mode

  • x[in] Fractional momentum loss

  • kt2[in] The \(k_{\rm T}\) transverse momentum norm

  • a[in] Mass number for the heavy ion

  • z[in] Atomic number for the heavy ion

double cepgen_particle_charge_(int &pdg_id)#

Charge of a particle, in e.

double cepgen_particle_mass_(int &pdg_id)#

Mass of a particle, in GeV/c^2.

These can be translated in the following Fortran subroutines/functions definitions:

external CepGen_kT_flux,CepGen_kT_flux_HI
double precision CepGen_kT_flux,CepGen_kT_flux_HI
external CepGen_particle_charge,CepGen_particle_mass
double precision CepGen_particle_charge,CepGen_particle_mass

Overall linking#

To interface your process to CepGen, edit the CepGenProcesses/ProcessesWrapper.cpp file to link your function implementation to the core processes module.

The following macros are used to declare the function name and link it to CepGen.

Warning

doxygendefine: Cannot find define “DECLARE_FORTRAN_FUNCTION” in doxygen xml output for project “CepGen” from directory: /Documentation/build/xml

REGISTER_FORTRAN_PROCESS(name, descr, f77_func)#

Add the Fortran process definition to the list of handled processes.

Note

No single, nor double trailing underscore (_) is required for this name registration.

For the latter, the signature is the following:

REGISTER_FORTRAN_PROCESS(my_dummy_process, dummy_f77_process,
                         "this process is only for illustration purposes")

Or, following this nomenclature:

  • the first argument corresponds to the process name (i.e. may be used with the my_dummy_process string in steering cards),

  • the second argument is the function name as defined above,

  • the third one is a one-line, human-readable process description that will be shown at the run start.

For this version, the following process is already registered:

List of registered external Fortran processes, as imported from CepGenProcesses/ProcessesWrapper.cpp.#
/*
 *  CepGen: a central exclusive processes event generator
 *  Copyright (C) 2013-2024   Laurent Forthomme
 *
 *  This program is free software: you can redistribute it and/or modify
 *  it under the terms of the GNU General Public License as published by
 *  the Free Software Foundation, either version 3 of the License, or
 *  any later version.
 *
 *  This program is distributed in the hope that it will be useful,
 *  but WITHOUT ANY WARRANTY; without even the implied warranty of
 *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 *  GNU General Public License for more details.
 *
 *  You should have received a copy of the GNU General Public License
 *  along with this program.  If not, see <http://www.gnu.org/licenses/>.
 */

#include "CepGen/Event/Event.h"
#include "CepGen/Modules/ProcessFactory.h"
#include "CepGen/Process/FortranFactorisedProcess.h"

//=============================================================================
// START THE MAPPING name -> Fortran matrix element evaluation function
// usage:
//  REGISTER_FORTRAN_PROCESS(name, "description", function_name)
// with the Fortran function name written in lowercase (no trailing '_')
//=============================================================================

REGISTER_FORTRAN_PROCESS(pptoff_f77, "(p/A)(p/A) ↝ (g/γ)γ → f⁺f¯", nucl_to_ff);