gibuu is hosted by Hepforge, IPPP Durham
GiBUU

TABLE OF CONTENTS


/Dilepton_Analysis [ Modules ]

[ Top ] [ Modules ]

NAME

module Dilepton_Analysis

PURPOSE

This analysis module produces various dilepton spectra. Currently only di-electrons are supported.


Dilepton_Analysis/Enable [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  logical, save :: Enable = .false.

PURPOSE

If .true. the dilepton analysis will be performed, otherwise not.


Dilepton_Analysis/Extra [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  logical, save :: Extra = .false.

PURPOSE

If .true. an extended analysis will be performed, writing out many extra histograms (beyond the basic ones: mass, pT and rapidity).


Dilepton_Analysis/DeltaDalitz [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  integer, save :: DeltaDalitz = 2

PURPOSE

Choose between different parametrizations of the Delta Dalitz decay width (Delta -> N e+e-):


Dilepton_Analysis/DeltaDalitzFF [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  integer, save :: DeltaDalitzFF = 1

PURPOSE

Choose a parametrization of the electromagnetic N-Delta transition form factor for the Delta Dalitz decay (only used for DeltaDalitz = 2):


Dilepton_Analysis/omegaDalitzFF [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  integer, save :: omegaDalitzFF = 1

PURPOSE

Choose between different parametrizations of the omega Dalitz decay (omega -> pi^0 e+e-) form factor:

  • 0 = constant
  • 1 = Effenberger/Bratkovskaya (default)
  • 2 = standard VMD
  • 3 = Terschluesen/Leupold


Dilepton_Analysis/b_pi [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  real, save :: b_pi = 5.5

PURPOSE

This constant represents the b parameter in the form factor of the pi0 Dalitz decay (in GeV^-2), cf. Effenberger Diss. eq. (2.141). Originally taken from L.G. Landsberg, Phys. Rep. 128, 301 (1985).


Dilepton_Analysis/lambda_eta [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  real, save :: lambda_eta = 0.716

PURPOSE

This constant represents the Lambda parameter in the form factor of the eta Dalitz decay in GeV. Values:

  • L.G. Landsberg, Phys. Rep. 128, 301 (1985): Lambda = 720 MeV
  • HADES pp@2.2, B. Spruck, Diss. (2008): Lambda = 676 MeV
  • NA60, Arnaldi et al, PLB 677 (2009): Lambda = 716 MeV (default)
  • CB/TAPS, Berghäuser et al, PLB 701 (2011): Lambda = 722 MeV


Dilepton_Analysis/etaPrimeDalitzFF [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  integer, save :: etaPrimeDalitzFF = 0

PURPOSE

Choose between different parametrizations of the eta' Dalitz decay (eta' -> e+e- gamma) form factor:

  • 0 = constant (default)
  • 1 = eta FF (cf. lambda_eta)
  • 2 = generic VMD
  • 3 = Genesis / Lepton-G
  • 4 = standard VMD (Terschluesen)


Dilepton_Analysis/angDist [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  integer, save :: angDist = 1

PURPOSE

This switch determines the angular distribution of the pseudoscalar Dalitz decays P -> e+ e- gamma (with P=pi0,eta,etaPrime):

  • 0 = isotropic decay
  • 1 = anisotropic decay according to 1 + B*cos**2(theta) with B=1
  • 2 = the Dalitz decays of pi0 and eta will be done via Pythia.


Dilepton_Analysis/brems [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  integer, save :: brems = 1

PURPOSE

This switch determines how the bremsstrahlung contribution is obtained:

  • 0 = none
  • 1 = soft-photon approximation (SPA)
  • 2 = according to the one-boson-exchange (OBE) model by R. Shyam (for NN bremstrahlung only, no em. form factors)
  • 3 = as 2, but with pion em. form factor (for pn)


Dilepton_Analysis/nEvent [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  integer, save :: nEvent = 10

PURPOSE

Number of events to generate for each dilepton decay (to enhance statistics).


Dilepton_Analysis/nEvent_BH [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  integer, save :: nEvent_BH = 1000

PURPOSE

Number of events for Bethe-Heitler simulation. BH typically needs a lot more statistics than the other dilepton channels. Therefore nEvent_BH should be much bigger than nEvent.


Dilepton_Analysis/kp_cut [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  logical, save :: kp_cut = .false.

PURPOSE

Perform a cut on (k*p) in the dilepton analysis, where k is the photon 4-momentum, and p is the electron or positron 4-momentum. This is useful for suppressing the BH contribution. Cf. "kp_min".


Dilepton_Analysis/kp_min [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  real, save :: kp_min = 0.01

PURPOSE

If kp_cut=.true. a cut on (k*p) is performed. kp_min determines the position of this cut. Only events with (k*p)>kp_min are taken into account.


Dilepton_Analysis/filter [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  integer, save :: filter = 0

PURPOSE

If filter is nonzero, a filtering algorithm will be applied to the dilepton pairs, otherwise they will be written to the histograms unfiltered. For details on the filtering parameters see routine 'CS'. Choices:

  • 0 = no filter
  • 1 = DLS
  • 2 = HADES (simple cuts on polar angle, absolute momentum and opening angle)
  • 3 = HADES (full acceptance filter, using pair acceptance)
  • 4 = HADES (full acceptance filter, using single-particle acceptance)
  • 5 = g7/CLAS @ JLab
  • 6 = KEK E325 (cuts on rapidity, transverse momentum and opening angle)
  • 7 = JPARC E16

NOTES

For filtering modes 3 and 4, the file containing the acceptance matrices must be specified (cf. hadesFilterFile).


Dilepton_Analysis/hadesFilterFile [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  character(len=1000), save :: hadesFilterFile = ""

PURPOSE

This character string determines the location of the file containing the HADES acceptance matrices (filename with absolute or relative path). It has to be set for filtering modes 3 and 4.


Dilepton_Analysis/p_lep_min [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  real, save :: p_lep_min = 0.

PURPOSE

This switch sets a lower bound on the lepton momentum. Only leptons with momenta larger than this threshold will pass the filter. This switch is only used for filtering mode 5 (JLab).


Dilepton_Analysis/beta_gamma_cut [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  real, save :: beta_gamma_cut = 1.25

PURPOSE

This is an upper bound on the beta*gamma value of the lepton pair. Since beta*gamma = p/m, it cuts on slow pairs.


Dilepton_Analysis/binsz [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  real, save :: binsz = 0.01

PURPOSE

This determines the bin size of the dilepton mass spectrum in GeV. Default is 10 MeV.


Dilepton_Analysis/WriteEvents [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  integer, save :: WriteEvents = 0

PURPOSE

This switch decides whether we write out the simulated events. Possible values:

  • 0: Don't write events (default).
  • 1: We write out only the lepton pair information (including charge, four-momentum, perturbative weight, source channel and filter result). All of this will be written to a file called 'Dilepton_Events.dat'.
  • 2: As 1, but only writing exclusive events (NN->NNe+e-).
  • 3: We write out all produced particles in the event (including the lepton pair) to a file called 'Dilepton_FullEvents.dat'.
  • 4: As 3, but only writing exclusive events (NN->NNe+e-).
  • 5: As 4, but only writing out R->Ne+e- events (with R=N*,Delta*).


Dilepton_Analysis/massBinning [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  real, dimension(1:4), save :: massBinning = (/ 0.150, 0.550, 9.999, 9.999 /)

PURPOSE

We produce several histograms (e.g. p,pT,mT,y,theta_cm) not only for the full mass range, but also for (up to 5) different mass bins. The borders of these bins are given by this array.


Dilepton_Analysis/particle_source [ Global module-variables ]

[ Top ] [ Dilepton_Analysis ] [ Global module-variables ]

SOURCE

  logical, save :: particle_source = .true.

PURPOSE

This switch determines whether the mass spectrum will contain separate channels for different sources of particles, such as decays (R->rho N) or collisions (pi pi -> rho, K K -> phi). Currently this is only done for the rho and phi mesons. Note: If using this switch, the "sum" channel in the mass histogram should not be used, since the rho and phi contributions will enter twice.


Dilepton_Analysis/Dilep_Init [ Subroutines ]

[ Top ] [ Dilepton_Analysis ] [ Subroutines ]

NAME

subroutine Dilep_Init (Ekin_max, N, cms, beta)

PURPOSE

This callback routine passes some eventtype-specific information to the Dilepton_Analysis module.

INPUTS

  • real, intent(in) :: Ekin_max - maximum energy to be expected during the simulation
  • integer, intent(in), optional :: N - number of projectile particles
  • logical, intent(in), optional :: cms - is the collision performed in the center-of-mass frame?
  • real, dimension(1:3), intent(in), optional :: beta - boost vector to the laboratory frame


Dilepton_Analysis/readInput [ Subroutines ]

[ Top ] [ Dilepton_Analysis ] [ Subroutines ]

NAME

subroutine readInput

PURPOSE

Reads input in jobcard out of namelist "DileptonAnalysis" and initializes all the histograms.


Dilepton_Analysis/DileptonAnalysis [ Namelists ]

[ Top ] [ Dilepton_Analysis ] [ Namelists ]

NAME

NAMELIST /DileptonAnalysis/

PURPOSE

Includes the switches:


Dilepton_Analysis/Dilep_Decays [ Subroutines ]

[ Top ] [ Dilepton_Analysis ] [ Subroutines ]

NAME

subroutine Dilep_Decays (t, part, last)

PURPOSE

Subroutine for calculation of dilepton spectra (all cross sections are calculated in microbarn/GeV).

INPUTS

  • real, intent(in) :: t - time [fm]
  • type(particle) :: part(:,:) - particle vector (real or perturbative)
  • integer,intent(in) :: last - are we in the last timestep yet? 0=during time evolution, before last timestep, 1=last timestep, before/during forced decays, 2=last timestep, after forced decays

OUTPUT

  • the dilepton spectra are stored in the histograms msigma,psigma,esigma etc
  • at then end of the run, these histograms are written to files 'Dilepton*.dat' etc (e.g. 'DileptonMass.dat'), cf. Dilep_Write_CS


Dilepton_Analysis/dGamma_dM_omegaDalitz [ Functions ]

[ Top ] [ Dilepton_Analysis ] [ Functions ]

NAME

real function dGamma_dM_omegaDalitz(m_D,M)

PURPOSE

This function calculates the mass differential decay width dGamma/dM of omega -> pi0 e+e-, neglecting the electron mass. See:

INPUTS

  • real,intent(in) :: m_D ! mass of the omega
  • real,intent(in) :: M ! invariant mass of the gamma*


Dilepton_Analysis/dGamma_dM_etaPrimeDalitz [ Functions ]

[ Top ] [ Dilepton_Analysis ] [ Functions ]

NAME

real function dGamma_dM_etaPrimeDalitz (m_etap, M)

PURPOSE

This function calculates the mass differential decay width dGamma/dM of eta' -> e+e- gamma, neglecting the electron mass. See:

  • Effenberger, PhD thesis, eq. (2.140).

INPUTS

  • real,intent(in) :: m_etap ! mass of the eta prime
  • real,intent(in) :: M ! invariant mass of the gamma*


Dilepton_Analysis/dGamma_dM_DeltaDalitz [ Functions ]

[ Top ] [ Dilepton_Analysis ] [ Functions ]

NAME

real function dGamma_dM_DeltaDalitz (W, q, charge)

PURPOSE

This function calculates the mass differential decay width dGamma/dM of Delta -> N e+e-. The charge of the Delta must be +1 or 0. See e.g. Krivoruchenko/Faessler, Phys. Rev. D65 (2001), 017502, equation (3) and (4).

INPUTS

  • real,intent(in) :: W ! mass of the Delta
  • real,intent(in) :: q ! invariant mass of the gamma*
  • integer,intent(in) :: charge ! charge of the Delta


Dilepton_Analysis/DeltaWidth_gammaN [ Functions ]

[ Top ] [ Dilepton_Analysis ] [ Functions ]

NAME

real function DeltaWidth_gammaN (W, q, charge)

PURPOSE

This function calculates the decay width of Delta -> N gamma*. The charge of the Delta must be +1 or 0. See e.g. Krivoruchenko/Faessler, Phys. Rev. D65 (2001), 017502, equation (2).

INPUTS

  • real,intent(in) :: W ! mass of the Delta
  • real,intent(in) :: q ! invariant mass of the gamma*
  • integer,intent(in) :: charge ! charge of the Delta


Dilepton_Analysis/Gamma0_DeltaDalitz_Wolf [ Functions ]

[ Top ] [ Dilepton_Analysis ] [ Functions ]

NAME

real function Gamma0_DeltaDalitz_Wolf (W, q)

PURPOSE

This function calculates the total decay rate of a Delta resonance going into a nucleon and a gamma*. See Wolf et al., Nucl. Phys. A517 (1990) 615-638, http://inspirehep.net/record/306273 , equation (4.9). Effenberger gives a slightly modified formula in his thesis, equation (2.139). REMARKS The coupling constant 'g' in Wolf's paper is wrong. Effenberger has the correct value of g = 5.44. With a real photon and an on-shell Delta one should get the photonic decay width of Gamma_0 (0) = 0.72 MeV.

INPUTS

  • real, intent(in) :: W ! mass of the Delta
  • real, intent(in) :: q ! invariant mass of the gamma*

OUTPUT

Gamma_0 in GeV.


Dilepton_Analysis/Gamma0_DeltaDalitz_Ernst [ Functions ]

[ Top ] [ Dilepton_Analysis ] [ Functions ]

NAME

real function Gamma0_DeltaDalitz_Ernst (W, q)

PURPOSE

This function calculates the total decay rate of a Delta resonance going into a nucleon and a gamma*, assuming a constant form factor. See Ernst et al., Phys. Rev. C 58 (1998) 447-456, http://inspirehep.net/record/452782 , equations (3),(12),(13).

INPUTS

  • real,intent(in) :: W ! mass of the Delta
  • real,intent(in) :: q ! invariant mass of the gamma* (dilepton)

OUTPUT

Gamma_0 in GeV.


Dilepton_Analysis/Gamma0_DeltaDalitz_Krivo [ Functions ]

[ Top ] [ Dilepton_Analysis ] [ Functions ]

NAME

real function Gamma0_DeltaDalitz_Krivo (W, q)

PURPOSE

This function calculates the total decay rate of a Delta resonance going into a nucleon and a gamma*. See Krivoruchenko/Faessler, Phys. Rev. D65 (2001) 017502, http://inspirehep.net/record/555421 , equation (2). REMARKS With a real photon and an on-shell Delta one should get the PDG value of Gamma_0 (0) = 0.66 MeV.

INPUTS

  • real, intent(in) :: W ! mass of the Delta
  • real, intent(in) :: q ! invariant mass of the gamma* (dilepton)

OUTPUT

Gamma_0 in GeV.


Dilepton_Analysis/Delta_FF_MAID [ Functions ]

[ Top ] [ Dilepton_Analysis ] [ Functions ]

NAME

real function Delta_FF_MAID (W,q2)

PURPOSE

This function uses the electromagnetic N-Delta transition form factors from MAID2005, as used in GiBUU electroproduction but with off-shell W and continued to negative Q^2. References:

  • D. Drechsel, O. Hanstein, S.S. Kamalov, L. Tiator; Nucl. Phys. A 645 (1999) 145, eq. (21) - (23).
  • M. Warns, H. Schroeder, W. Pfeil, H. Rollnik; Z. Phys. C 45 (1990) 627, eq. (2.16) - (2.17)

INPUTS

  • real,intent(in) :: W ! mass of the Delta
  • real,intent(in) :: q2 ! q^2 = m_ee^2 = - Q^2

OUTPUT

Squared form factor.


Dilepton_Analysis/Delta_FF_Iachello [ Functions ]

[ Top ] [ Dilepton_Analysis ] [ Functions ]

NAME

real function Delta_FF_Iachello (W, q2)

PURPOSE

This function calculates the electromagnetic N-Delta transition form factor according to the two-component quark model by Wan/Iachello. References:

INPUTS

  • real, intent(in) :: W ! off-shell mass of the Delta [GeV]
  • real, intent(in) :: q2 ! q^2 = m_ee^2 [GeV^2]

OUTPUT

Squared form factor.


Dilepton_Analysis/dGamma_dM_N1520Dalitz [ Functions ]

[ Top ] [ Dilepton_Analysis ] [ Functions ]

NAME

real function dGamma_dM_N1520Dalitz (W, q)

PURPOSE

This function calculates the mass differential decay width dGamma/dM of N*(1520) -> N e+e-. See Krivoruchenko/Martemyanov/Faessler/Fuchs, Ann. Phys. 296 (2002), 299-346, equation (III.22).

INPUTS

  • real,intent(in) :: W ! mass of the N*
  • real,intent(in) :: q ! invariant mass of the gamma*


Dilepton_Analysis/CountParticles [ Subroutines ]

[ Top ] [ Dilepton_Analysis ] [ Subroutines ]

NAME

subroutine CountParticles(nt,pertPart)

PURPOSE

This subroutine simply counts the number of particles at each timestep. Currently only rho^0, omega, phi, pi^0, eta and Delta^(0,+) are counted. All these particles are counted separately, and in addition the particles which are not in formation are counted (for each species). The numbers are per ensemble, and averaged over all the runs.

INPUTS

OUTPUT

The data is stored in two multi-channel histograms: "partnum" and "partnum_noform".


Dilepton_Analysis/writePartNum [ Subroutines ]

[ Top ] [ Dilepton_Analysis ] [ Subroutines ]

NAME

subroutine writePartNum()

PURPOSE

Writes the data collected by "CountParticles" to disk, creating four output files:


Dilepton_Analysis/PartNum.dat [ Output files ]

[ Top ] [ Dilepton_Analysis ] [ Output files ]

NAME

file PartNum.dat

PURPOSE

Total particle numbers per timestep, including those which are still in formation.


Dilepton_Analysis/PartNumInt.dat [ Output files ]

[ Top ] [ Dilepton_Analysis ] [ Output files ]

NAME

file PartNumInt.dat

PURPOSE

Time-integrated particle numbers.


Dilepton_Analysis/PartNum_NoForm.dat [ Output files ]

[ Top ] [ Dilepton_Analysis ] [ Output files ]

NAME

file PartNum_NoForm.dat

PURPOSE

Numbers of particles which are not in formation.


Dilepton_Analysis/VMmass.dat [ Output files ]

[ Top ] [ Dilepton_Analysis ] [ Output files ]

NAME

file VMmass.dat

PURPOSE

Mass spectrum of all neutral vector mesons (rho0, omega and phi) produced in the simulation (independent of decay mode, integrated over time).


Dilepton_Analysis/WriteFullEvent [ Subroutines ]

[ Top ] [ Dilepton_Analysis ] [ Subroutines ]

NAME

subroutine WriteFullEvent (parts, dil, pnr, pw, iso)

PURPOSE

Print out a full dilepton event, including all hadrons and the lepton pair. The event will be written to the file "Dilepton_FullEvents.dat" in a format similar to the LesHouches format.


Dilepton_Analysis/Dilepton_FullEvents.dat [ Output files ]

[ Top ] [ Dilepton_Analysis ] [ Output files ]

NAME

file Dilepton_FullEvents.dat

PURPOSE

This file contains dilepton events (including all hadrons). It will only be written if enabled by the WriteEvents switch. Attention: This file can get very large for long runs! The format is XML-like and similar to the LesHouches format: Each event is contained in an <event> .... </event> block, with a header in the <event> line that contains the perturbative weight, source type (channel) and filter result (1=accept, 0=reject). Each line in the event describes one particle and has the following entries:

  • Column #1: Particle ID.
  • Column #2: Charge.
  • Column #3-6: 4-momentum (E,p_x,p_y,p_z) in GeV.
  • Column #7: Parent particle ID (if the particle came from a decay, 0 otherwise).


Dilepton_Analysis/Dilepton_Events.dat [ Output files ]

[ Top ] [ Dilepton_Analysis ] [ Output files ]

NAME

file Dilepton_Events.dat

PURPOSE

This file contains dilepton pair events (leptons only, no hadrons). It will only be written if enabled by the WriteEvents switch. Attention: This file can get very large for long runs! Entries:

  • Column #1: Projectile energy in GeV.
  • Column #2: Lepton charge (+ or -).
  • Column #3-5: Lepton 3-momentum (p_x,p_y,p_z) in GeV.
  • Column #6: Perturbative weight.
  • Column #7: Source type (channel).
  • Column #8: Filter result (1=accept, 0=reject).


Dilepton_Analysis/CS [ Subroutines ]

[ Top ] [ Dilepton_Analysis ] [ Subroutines ]

NAME

subroutine CS (pout_, iso, pw, dens, gen, part, enr, pnr)

PURPOSE

Calculates differential cross sections (in microbarn) and puts them into various histograms.

INPUTS

  • real, dimension(0:3,1:2), intent(in) :: pout --- 4-momenta of the lepton pair (1=positron, 2=electron)
  • integer, intent(in) :: iso --- source type (channel)
  • real, intent(in) :: pw --- perturbative weight
  • real, intent(in) :: dens --- density at decay vertex
  • integer, intent(in) :: gen --- generation
  • type(particle), intent(in), optional :: part(:,:) --- particle vector
  • integer, intent(in), optional :: enr --- ensemble nr
  • integer, intent(in), optional :: pnr --- particle nr

OUTPUT

  • cross section data is stored in histograms msigma, psigma, esigma, ...


Dilepton_Analysis/applyFilter [ Functions ]

[ Top ] [ Dilepton_Analysis ] [ Functions ]

NAME

integer function applyFilter (p1,p2)

PURPOSE

Performs acceptance filtering and detector resolution smearing based on the value of the 'filter' switch.

INPUTS

  • real,dimension(0:3)::p1,p2 --- 4-momenta of positron and electron

OUTPUT

  • returns "1" if the events is accepted, "0" if rejected
  • p1 and p2 will be modified if resolution smearing is applied


Dilepton_Analysis/DileptonMass.dat [ Output files ]

[ Top ] [ Dilepton_Analysis ] [ Output files ]

NAME

file DileptonMass.dat

PURPOSE

Dilepton mass spectrum.


Dilepton_Analysis/DileptonPlab.dat [ Output files ]

[ Top ] [ Dilepton_Analysis ] [ Output files ]

NAME

file DileptonPlab.dat

PURPOSE

Dilepton momentum spectrum.


Dilepton_Analysis/DileptonPt.dat [ Output files ]

[ Top ] [ Dilepton_Analysis ] [ Output files ]

NAME

file DileptonPt.dat

PURPOSE

Dilepton transverse momentum spectrum.


Dilepton_Analysis/DileptonMt.dat [ Output files ]

[ Top ] [ Dilepton_Analysis ] [ Output files ]

NAME

file DileptonMt.dat

PURPOSE

Dilepton transverse mass spectrum.


Dilepton_Analysis/DileptonY.dat [ Output files ]

[ Top ] [ Dilepton_Analysis ] [ Output files ]

NAME

file DileptonY.dat

PURPOSE

Dilepton rapidity spectrum.


Dilepton_Analysis/DileptonEkin.dat [ Output files ]

[ Top ] [ Dilepton_Analysis ] [ Output files ]

NAME

file DileptonEkin.dat

PURPOSE

Dilepton kinetic energy spectrum.


Dilepton_Analysis/Dilep_write_CS [ Subroutines ]

[ Top ] [ Dilepton_Analysis ] [ Subroutines ]

NAME

subroutine Dilep_write_CS

PURPOSE

Writes out the cross section spectra (called at the end each run).

INPUTS

OUTPUT

The spectra are written to files called 'DileptonMass.dat', 'DileptonPlab.dat', etc


Dilepton_Analysis/DilepSim2 [ Functions ]

[ Top ] [ Dilepton_Analysis ] [ Functions ]

NAME

function DilepSim2 (ptot, massdil) result(plep)

PURPOSE

Simulates a dilepton decay with 2-particle final state, using isotropic angular distribution.

INPUTS

  • real, intent(in) :: ptot(0:3) - four momentum of decaying particle
  • real, intent(in) :: massdil - invariant mass of the dilepton pair

OUTPUT

  • real, dimension(0:3,1:2) :: plep - 4-vectors of the generated lepton pair


Dilepton_Analysis/DilepSim3 [ Functions ]

[ Top ] [ Dilepton_Analysis ] [ Functions ]

NAME

function DilepSim3 (ptot, massdil, mass1, ang) result (plep)

PURPOSE

Simulates a dilepton decay with 3-particle final state (Dalitz decay).

INPUTS

  • real, intent(in) :: ptot(0:3) - four momentum of decaying particle
  • real, intent(in) :: massdil - invariant mass of the dilepton pair
  • real, intent(in) :: mass1 - mass of 2nd particle
  • integer, intent(in) :: ang - angular distribution of gamma* decay: 0=isotropic, 1,2=anisotropic

OUTPUT

  • real, dimension(0:3,1:2) :: plep - 4-vectors of the generated lepton pair


Dilepton_Analysis/PseudoscalarDalitzPythia [ Functions ]

[ Top ] [ Dilepton_Analysis ] [ Functions ]

NAME

function PseudoscalarDalitzPythia (ID, ptot) result (pdil)

PURPOSE

Performs the Dalitz decay of a pseudoscalar meson (pi0 or eta) via Pythia. Cf. Pythia 6.4 manual, chapter 13.3.1 ("Strong and electromagnetic decays"), eq. (275) and (276).

INPUTS

  • integer, intent(in) :: ID --- ID of decaying particle (pi0 or eta)
  • real, dimension(0:3), intent(in) --- four momentum of decaying particle

OUTPUT

  • real,dimension(0:3,1:2) :: pdil --- four vectors of the generated lepton pair (1=positron,2=electron)


Dilepton_Analysis/Dilep_Brems [ Subroutines ]

[ Top ] [ Dilepton_Analysis ] [ Subroutines ]

NAME

subroutine Dilep_Brems (part, srtfree, XS_elast, XS_tot)

PURPOSE

Calculate Dilepton Production via Nucleon-Nucleon (or pion-Nucleon) Bremsstrahlung in Soft-Photon-Approximation (SPA). See: Wolf et al., Nucl. Phys. A517 (1990) 615-638, chapter 4.2

INPUTS

  • type(particle), intent(in) :: part(1:2) --- colliding particles
  • real, intent(in) :: srtfree --- CoM energy
  • real, intent(in) :: XS_elast --- elastic cross section in mb (NN or piN)
  • real, intent(in) :: XS_tot --- total cross section in mb (NN or piN)

TODO

  • Pauli blocking ?


Dilepton_Analysis/Dilep_BH [ Subroutines ]

[ Top ] [ Dilepton_Analysis ] [ Subroutines ]

NAME

subroutine Dilep_BH(k_in,nuc)

PURPOSE

Calculate Dilepton Production via the Bethe-Heitler Process. See Diploma Thesis J.Weil.

INPUTS