#include <scatteractionphoton.h>
ScatterActionPhoton is a special action which takes two incoming particles and performs a perturbative electromagnetic scattering.
The final state particles are not further propagated, only written to the output.
Definition at line 27 of file scatteractionphoton.h.
Public Member Functions | |
ScatterActionPhoton (const ParticleList &in, const double time, const int n_frac_photons, const double hadronic_cross_section_input) | |
Construct a ScatterActionPhoton object. More... | |
void | perform_photons (const OutputsList &outputs) |
Create the photon final state and write to output. More... | |
void | generate_final_state () override |
Generate the final-state for the photon scatter process. More... | |
double | get_total_weight () const override |
Return the weight of the last created photon. More... | |
double | hadronic_cross_section () const |
Return the total cross section of the underlying hadronic scattering. More... | |
double | sample_out_hadron_mass (const ParticleTypePtr out_type) |
Sample the mass of the outgoing hadron. More... | |
void | add_dummy_hadronic_process (double reaction_cross_section) |
Adds one hadronic process with a given cross-section. More... | |
void | add_single_process () |
Add the photonic process. More... | |
![]() | |
ScatterAction (const ParticleData &in_part1, const ParticleData &in_part2, double time, bool isotropic=false, double string_formation_time=1.0) | |
Construct a ScatterAction object. More... | |
void | add_collision (CollisionBranchPtr p) |
Add a new collision channel. More... | |
void | add_collisions (CollisionBranchList pv) |
Add several new collision channels at once. More... | |
double | transverse_distance_sqr () const |
Calculate the transverse distance of the two incoming particles in their local rest frame. More... | |
double | get_partial_weight () const override |
Get the partial cross section of the chosen channel. More... | |
void | sample_angles (std::pair< double, double > masses, double kinetic_energy_cm) override |
Sample final-state angles in a 2->2 collision (possibly anisotropic). More... | |
void | add_all_scatterings (double elastic_parameter, bool two_to_one, ReactionsBitSet included_2to2, double low_snn_cut, bool strings_switch, bool use_AQM, bool strings_with_probability, NNbarTreatment nnbar_treatment) |
Add all possible scattering subprocesses for this action object. More... | |
const CollisionBranchList & | collision_channels () |
Get list of possible collision channels. More... | |
void | set_string_interface (StringProcess *str_proc) |
Set the StringProcess object to be used. More... | |
virtual double | cross_section () const |
Get the total cross section of the scattering particles. More... | |
![]() | |
Action (const ParticleList &in_part, double time) | |
Construct an action object with incoming particles and relative time. More... | |
Action (const ParticleData &in_part, const ParticleData &out_part, double time, ProcessType type) | |
Construct an action object with the incoming particles, relative time, and the already known outgoing particles and type of the process. More... | |
Action (const ParticleList &in_part, const ParticleList &out_part, double absolute_execution_time, ProcessType type) | |
Construct an action object with the incoming particles, absolute time, and the already known outgoing particles and type of the process. More... | |
Action (const Action &)=delete | |
Copying is disabled. Use pointers or create a new Action. More... | |
virtual | ~Action () |
Virtual Destructor. More... | |
bool | operator< (const Action &rhs) const |
Determine whether one action takes place before another in time. More... | |
virtual ProcessType | get_type () const |
Get the process type. More... | |
template<typename Branch > | |
void | add_process (ProcessBranchPtr< Branch > &p, ProcessBranchList< Branch > &subprocesses, double &total_weight) |
Add a new subprocess. More... | |
template<typename Branch > | |
void | add_processes (ProcessBranchList< Branch > pv, ProcessBranchList< Branch > &subprocesses, double &total_weight) |
Add several new subprocesses at once. More... | |
virtual void | perform (Particles *particles, uint32_t id_process) |
Actually perform the action, e.g. More... | |
bool | is_valid (const Particles &particles) const |
Check whether the action still applies. More... | |
bool | is_pauli_blocked (const Particles &particles, const PauliBlocker &p_bl) const |
Check if the action is Pauli-blocked. More... | |
const ParticleList & | incoming_particles () const |
Get the list of particles that go into the action. More... | |
void | update_incoming (const Particles &particles) |
Update the incoming particles that are stored in this action to the state they have in the global particle list. More... | |
const ParticleList & | outgoing_particles () const |
Get the list of particles that resulted from the action. More... | |
double | time_of_execution () const |
Get the time at which the action is supposed to be performed. More... | |
void | check_conservation (const uint32_t id_process) const |
Check various conservation laws. More... | |
double | sqrt_s () const |
Determine the total energy in the center-of-mass frame [GeV]. More... | |
FourVector | total_momentum_of_outgoing_particles () const |
Calculate the total kinetic momentum of the outgoing particles. More... | |
FourVector | get_interaction_point () const |
Get the interaction point. More... | |
std::pair< FourVector, FourVector > | get_potential_at_interaction_point () const |
Get the skyrme and asymmetry potential at the interaction point. More... | |
Static Public Member Functions | |
static ReactionType | photon_reaction_type (const ParticleList &in) |
Determine photon process from incoming particles. More... | |
static bool | is_photon_reaction (const ParticleList &in) |
Check if particles can undergo an implemented photon process. More... | |
static ParticleTypePtr | outgoing_hadron_type (const ParticleList &in) |
Return ParticleTypePtr of hadron in the out channel, given the incoming particles. More... | |
static ParticleTypePtr | outgoing_hadron_type (const ReactionType reaction) |
Return ParticleTypePtr of hadron in the out channel, given the ReactionType. More... | |
static bool | is_kinematically_possible (const double s_sqrt, const ParticleList &in) |
Check if CM-energy is sufficient to produce hadron in final state. More... | |
Private Types | |
enum | MediatorType { MediatorType::SUM, MediatorType::PION, MediatorType::OMEGA } |
Compile-time switch for setting the handling of processes which can happen via different mediating particles. More... | |
Private Member Functions | |
double | diff_cross_section (const double t, const double m_rho, MediatorType mediator=default_mediator_) const |
Calculate the differential cross section of photon process. More... | |
double | rho_mass () const |
Find the mass of the participating rho-particle. More... | |
CollisionBranchList | photon_cross_sections (MediatorType mediator=default_mediator_) |
Computes the total cross section of the photon process. More... | |
std::pair< double, double > | diff_cross_section_single (const double t, const double m_rho) |
For processes which can happen via (pi, a1, rho) and omega exchange, return the differential cross section for the (pi, a1, rho) process in the first argument, for the omega process in the second. More... | |
std::pair< double, double > | form_factor_single (const double E_photon) |
For processes which can happen via (pi, a1, rho) and omega exchange, return the form factor for the (pi, a1, rho) process in the first argument, for the omega process in the second. More... | |
double | form_factor_pion (const double E_photon) const |
Compute the form factor for a process with a pion as the lightest exchange particle. More... | |
double | form_factor_omega (const double E_photon) const |
Compute the form factor for a process with a omega as the lightest exchange particle. More... | |
double | diff_cross_section_w_ff (const double t, const double m_rho, const double E_photon) |
Compute the differential cross section with form factors included. More... | |
Private Attributes | |
CollisionBranchList | collision_processes_photons_ |
Holds the photon branch. More... | |
const ReactionType | reac_ |
Photonic process as determined from incoming particles. More... | |
const int | number_of_fractional_photons_ |
Number of photons created for each hadronic scattering, needed for correct weighting. More... | |
const ParticleTypePtr | hadron_out_t_ |
ParticleTypePtr to the type of the outgoing hadron. More... | |
const double | hadron_out_mass_ |
Mass of outgoing hadron. More... | |
double | weight_ = 0.0 |
Weight of the produced photon. More... | |
double | cross_section_photons_ = 0.0 |
Total cross section of photonic process. More... | |
const double | hadronic_cross_section_ |
Total hadronic cross section. More... | |
Static Private Attributes | |
static constexpr MediatorType | default_mediator_ = MediatorType::SUM |
Value used for default exchange particle. See MediatorType. More... | |
Additional Inherited Members | |
![]() | |
double | mandelstam_s () const |
Determine the Mandelstam s variable,. More... | |
double | cm_momentum () const |
Get the momentum of the center of mass of the incoming particles in the calculation frame. More... | |
double | cm_momentum_squared () const |
Get the squared momentum of the center of mass of the incoming particles in the calculation frame. More... | |
ThreeVector | beta_cm () const |
Get the velocity of the center of mass of the scattering particles in the calculation frame. More... | |
double | gamma_cm () const |
Get the gamma factor corresponding to a boost to the center of mass frame of the colliding particles. More... | |
void | elastic_scattering () |
Perform an elastic two-body scattering, i.e. just exchange momentum. More... | |
void | inelastic_scattering () |
Perform an inelastic two-body scattering, i.e. new particles are formed. More... | |
void | string_excitation () |
Todo(ryu): document better - it is not really UrQMD-based, isn't it? Perform the UrQMD-based string excitation and decay. More... | |
void | format_debug_output (std::ostream &out) const override |
Writes information about this scatter action to the out stream. More... | |
![]() | |
FourVector | total_momentum () const |
Sum of 4-momenta of incoming particles. More... | |
template<typename Branch > | |
const Branch * | choose_channel (const ProcessBranchList< Branch > &subprocesses, double total_weight) |
Decide for a particular final-state channel via Monte-Carlo and return it as a ProcessBranch. More... | |
virtual std::pair< double, double > | sample_masses (double kinetic_energy_cm) const |
Sample final-state masses in general X->2 processes (thus also fixing the absolute c.o.m. More... | |
void | sample_2body_phasespace () |
Sample the full 2-body phase-space (masses, momenta, angles) in the center-of-mass frame for the final state particles. More... | |
virtual void | sample_3body_phasespace () |
Sample the full 3-body phase-space (masses, momenta, angles) in the center-of-mass frame for the final state particles. More... | |
![]() | |
CollisionBranchList | collision_channels_ |
List of possible collisions. More... | |
double | total_cross_section_ |
Total hadronic cross section. More... | |
double | partial_cross_section_ |
Partial cross-section to the chosen outgoing channel. More... | |
bool | isotropic_ = false |
Do this collision isotropically? More... | |
double | string_formation_time_ = 1.0 |
Time fragments take to be fully formed in hard string excitation. More... | |
![]() | |
ParticleList | incoming_particles_ |
List with data of incoming particles. More... | |
ParticleList | outgoing_particles_ |
Initially this stores only the PDG codes of final-state particles. More... | |
const double | time_of_execution_ |
Time at which the action is supposed to be performed (absolute time in the lab frame in fm/c). More... | |
ProcessType | process_type_ |
type of process More... | |
|
strong |
Enum for encoding the photon process.
It is uniquely determined by the incoming particles. The naming scheme is : Incoming_1__Incoming_2__Outgoing_hadron. The photon is omitted in the naming.
Enumerator | |
---|---|
no_reaction | |
pi_z_pi_p_rho_p | |
pi_z_pi_m_rho_m | |
pi_p_rho_z_pi_p | |
pi_m_rho_z_pi_m | |
pi_m_rho_p_pi_z | |
pi_p_rho_m_pi_z | |
pi_z_rho_p_pi_p | |
pi_z_rho_m_pi_m | |
pi_p_pi_m_rho_z | |
pi_z_rho_z_pi_z |
Definition at line 112 of file scatteractionphoton.h.
|
strongprivate |
Compile-time switch for setting the handling of processes which can happen via different mediating particles.
Relevant only for the processes pi0 + rho => pi + y and pi + rho => pi0 + gamma, which both can happen via exchange of (rho, a1, pi) or omega. If MediatorType::SUM is set, the cross section for both processes is added. If MediatorType::PION/ OMEGA is set, only the respective processes are computed.
Enumerator | |
---|---|
SUM | |
PION | |
OMEGA |
Definition at line 215 of file scatteractionphoton.h.
smash::ScatterActionPhoton::ScatterActionPhoton | ( | const ParticleList & | in, |
const double | time, | ||
const int | n_frac_photons, | ||
const double | hadronic_cross_section_input | ||
) |
Construct a ScatterActionPhoton object.
[in] | in | ParticleList of incoming particles. |
[in] | time | Time relative to underlying hadronic action. |
[in] | n_frac_photons | Number of photons to produce for each hadronic scattering. |
[in] | hadronic_cross_section_input | Cross-section of underlying hadronic cross-section. |
Definition at line 29 of file scatteractionphoton.cc.
void smash::ScatterActionPhoton::perform_photons | ( | const OutputsList & | outputs | ) |
Create the photon final state and write to output.
[in] | outputs | List of all outputs. Does not have to be a specific photon output, the function will take care of this. |
Definition at line 92 of file scatteractionphoton.cc.
|
overridevirtual |
Generate the final-state for the photon scatter process.
Generates only one photon / hadron pair
Reimplemented from smash::ScatterAction.
Definition at line 193 of file scatteractionphoton.cc.
|
inlineoverridevirtual |
Return the weight of the last created photon.
Reimplemented from smash::ScatterAction.
Definition at line 64 of file scatteractionphoton.h.
|
inline |
Return the total cross section of the underlying hadronic scattering.
Definition at line 71 of file scatteractionphoton.h.
double smash::ScatterActionPhoton::sample_out_hadron_mass | ( | const ParticleTypePtr | out_type | ) |
Sample the mass of the outgoing hadron.
Returns the pole mass if particle is stable.
[in] | out_type | TypePtr of the outgoing hadron. |
Definition at line 291 of file scatteractionphoton.cc.
void smash::ScatterActionPhoton::add_dummy_hadronic_process | ( | double | reaction_cross_section | ) |
Adds one hadronic process with a given cross-section.
The intended use is to add the hadronic cross-section from the already performed hadronic action without recomputing it.
[in] | reaction_cross_section | Total cross-section of underlying hadronic process [mb] |
Definition at line 283 of file scatteractionphoton.cc.
|
inline |
Add the photonic process.
Also compute the total cross section as a side effect.
Definition at line 100 of file scatteractionphoton.h.
|
static |
Determine photon process from incoming particles.
If incoming particles are not part of any implemented photonic process, return no_reaction.
[in] | in | ParticleList of incoming particles. |
Definition at line 39 of file scatteractionphoton.cc.
|
inlinestatic |
Check if particles can undergo an implemented photon process.
This function does not check the involved kinematics.
[in] | in | ParticleList of incoming particles. |
Definition at line 145 of file scatteractionphoton.h.
|
static |
Return ParticleTypePtr of hadron in the out channel, given the incoming particles.
This function is overloaded since we need the hadron type in different places.
[in] | in | ParticleList of incoming particles. |
Definition at line 150 of file scatteractionphoton.cc.
|
static |
Return ParticleTypePtr of hadron in the out channel, given the ReactionType.
This function is overloaded since we need the hadron type in different places.
[in] | reaction | ReactionType, determined from incoming particles. |
Definition at line 104 of file scatteractionphoton.cc.
|
static |
Check if CM-energy is sufficient to produce hadron in final state.
[in] | s_sqrt | CM-energy [GeV] |
[in] | in | ParticleList of incoming hadrons |
Definition at line 156 of file scatteractionphoton.cc.
|
private |
Calculate the differential cross section of photon process.
Formfactors are not included
[in] | t | Mandelstam-t [GeV^2]. |
[in] | m_rho | Mass of the incoming or outgoing rho-particle [GeV] |
[in] | mediator | Switch for determing which mediating particle to use |
Definition at line 417 of file scatteractionphoton.cc.
|
private |
Find the mass of the participating rho-particle.
In case of a rho in the incoming channel it is the mass of the incoming rho, in case of an rho in the outgoing channel it is the mass sampled in the constructor. When an rho acts in addition as a mediator, its mass is the same as the incoming / outgoing rho. This function returns the alrady sampled mass or the mass of the incoming rho, depending on the process.
Definition at line 307 of file scatteractionphoton.cc.
|
private |
Computes the total cross section of the photon process.
[in] | mediator | Switch for determing which mediating particle to use. |
Definition at line 333 of file scatteractionphoton.cc.
|
private |
For processes which can happen via (pi, a1, rho) and omega exchange, return the differential cross section for the (pi, a1, rho) process in the first argument, for the omega process in the second.
If only one process exists, both values are the same.
[in] | t | Mandelstam-t [GeV^2] |
[in] | m_rho | Mass of the incoming or outgoing rho-particle [GeV] |
Definition at line 574 of file scatteractionphoton.cc.
|
private |
For processes which can happen via (pi, a1, rho) and omega exchange, return the form factor for the (pi, a1, rho) process in the first argument, for the omega process in the second.
If only one process exists, both values are the same.
[in] | E_photon | Energy of the photon [GeV] |
Definition at line 568 of file scatteractionphoton.cc.
|
private |
Compute the form factor for a process with a pion as the lightest exchange particle.
See wiki for details how form factors are handled.
[in] | E_photon | Energy of photon [GeV] |
Definition at line 543 of file scatteractionphoton.cc.
|
private |
Compute the form factor for a process with a omega as the lightest exchange particle.
See wiki for details how form factors are handled.
[in] | E_photon | Energy of photon [GeV] |
Definition at line 555 of file scatteractionphoton.cc.
|
private |
Compute the differential cross section with form factors included.
Takes care of correct handling of reactions with multiple processes by reading the default_mediator_ member variable.
[in] | t | Mandelstam-t [GeV^2] |
[in] | m_rho | Mass of the incoming or outgoing rho-particle [GeV] |
[in] | E_photon | of outgoing photon [GeV] |
The form factor is assumed to be a hadronic dipole form factor which takes the shape: FF = (2*Lambda^2/(2*Lambda^2 - t))^2 with Lambda = 1.0 GeV. t depends on the lightest possible exchange particle in the different channels. This could either be a pion or an omega meson. For the computation the parametrizations given in (Turbide:2006) are used.
Definition at line 478 of file scatteractionphoton.cc.
|
private |
Holds the photon branch.
As of now, this will always hold only one branch.
Definition at line 188 of file scatteractionphoton.h.
|
private |
Photonic process as determined from incoming particles.
Definition at line 191 of file scatteractionphoton.h.
|
private |
Number of photons created for each hadronic scattering, needed for correct weighting.
Note that in generate_final_state() only one photon + hadron is created.
Definition at line 198 of file scatteractionphoton.h.
|
private |
ParticleTypePtr to the type of the outgoing hadron.
Definition at line 201 of file scatteractionphoton.h.
|
private |
Mass of outgoing hadron.
Definition at line 204 of file scatteractionphoton.h.
|
staticprivate |
Value used for default exchange particle. See MediatorType.
Definition at line 217 of file scatteractionphoton.h.
|
private |
Weight of the produced photon.
Definition at line 220 of file scatteractionphoton.h.
|
private |
Total cross section of photonic process.
Definition at line 223 of file scatteractionphoton.h.
|
private |
Total hadronic cross section.
Definition at line 226 of file scatteractionphoton.h.