gibuu is hosted by Hepforge, IPPP Durham
GiBUU

TABLE OF CONTENTS


/eventGenerator_eN_lowEnergy [ Modules ]

[ Top ] [ Modules ]

NAME

module eventGenerator_eN_lowEnergy

PURPOSE

This module includes initilization routines for low energetic electron induced events (Resonance region and quasi-elastic regime, i.e. 0.7<W<2 GeV)

INPUTS

no namelist. All parameters are given as function arguments.


eventGenerator_eN_lowEnergy/iC_XXX [ Global module-variables ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Global module-variables ]

PURPOSE

Constants defined for ease of coding: possible Channels

SOURCE

  integer, parameter  :: iC_QE     = 1
  integer, parameter  :: iC_Res    = 2
  integer, parameter  :: iC_1Pi    = 3
  integer, parameter  :: iC_2Pi    = 4
  integer, parameter  :: iC_DIS    = 5
  integer, parameter  :: iC_VMDrho = 6
  integer, parameter  :: iC_2p2hQE = 7
  integer, parameter  :: iC_2p2hDelta = 8

  integer, parameter  :: nC        = 8

eventGenerator_eN_lowEnergy/eventGen_eN_lowEnergy [ Subroutines ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Subroutines ]

NAME

subroutine eventGen_eN_lowEnergy(eN,doC,whichRes,DoPauli,OutPart,channel,flagOK,XS,XS_Arr)

PURPOSE

Initializes one e^- N -> X events for 0.7<W<2 GeV . The weight each event is given by (total Xsection of event=sum over all possible channels). We make a Monte Carlo decision to choose a specific channel, which is returned as "OutPart".

The particles are produced at the place of the nucleon target.

INPUTS

Scattering particles:

  • type(electronNucleon_event) :: eN -- The incoming electron and nucleon
  • logical, dimension(:) :: doC -- switches, see below
  • logical :: DoPauli -- if .true., every event is checked for pauli blocking. if it is blocked, the corresponding cross section will be set to 0 and the Monte Carlo decision will neglect this event class. Thus the total cross section is reduced. Maybe all channel can be closed.
  • type(particle),dimension(:,:) :: realParticles -- real particle vector, only needed for pauli blocking

Switches (cf. iC_'module Electron_origin' and XXX):

  • QE -- Switch on/off quasi-elastic events
  • Res -- Switch on/off resonance production events
  • 1Pi -- Switch on/off single pion production events
  • 2Pi -- Switch on/off 2 pion background events
  • DIS -- Switch on/off DIS events
  • VMDrho -- Switch on/off special treatment of gamma N -> rho0 N
  • 2p2hQE -- Switch on/off gamma N N --> N' N'
  • 2p2hDelta -- Switch on/off gamma N N --> Delta N'

Special Switch for resonance production:

  • logical, dimension(2:nres+1), intent(in) :: whichRes

With this switch special resonances can be selected, if res_flag=.true., and whichRes is defined by the lines...

   whichRes=.false.
   whichRes(Delta)=.true.

... then only the Delta is included as a possible resonance channel.

OUTPUT

  • integer :: channel -- value according to naming scheme defined in module "Electron_origin". For all 1-Body final states the value of "channel" is given by the ID of the produced particle.
  • type(particle), dimension(:) :: OutPart -- FinalState particles
  • logical :: flagOK -- .true. if OutPart was created
  • real, OPTIONAL :: XS -- cross section in mub
  • real,dimension(...),OPTIONAL :: XS_Arr -- cross sections according the internal weights (in mub, maybe also negative)


eventGenerator_eN_lowEnergy/checkEvent [ Functions ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Functions ]

NAME

logical function checkEvent(eN,f,channel)

PURPOSE

Check Charge and momentum conservation (.true. -> all is ok!)

INPUTS

OUTPUT

  • function value

NOTES


eventGenerator_eN_lowEnergy/init_2p2hQE [ Subroutines ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Subroutines ]

NAME

subroutine init_2p2hQE(eN,OutPart,XS)

PURPOSE

Generate a e^- N1 N2 -> N1' N2' event.

NOTES

just a wrapper for el_2p2h_DoQE


eventGenerator_eN_lowEnergy/init_2p2hDelta [ Subroutines ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Subroutines ]

NAME

subroutine init_2p2hDelta(eN,OutPart,XS)

PURPOSE

Generate a e^- N1 N2 -> N Delta event.

NOTES

just a wrapper for el_2p2h_DoDelta


eventGenerator_eN_lowEnergy/init_VMDrho [ Subroutines ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Subroutines ]

NAME

subroutine init_VMDrho(eN,OutPart,XS)

PURPOSE

Generate a e^- N -> rho0 event.

NOTES

This can be easily generalized to all other vector mesons.


eventGenerator_eN_lowEnergy/init_DIS [ Subroutines ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Subroutines ]

NAME

subroutine init_DIS(eN,OutPart,XS)

PURPOSE

Generate a e^- N -> DIS event.

NOTES

This is mainly just a wrapper around DoColl_gammaN_Py.

The returned cross section (and also the weights of the particles) is 'flux * sigma^*',


eventGenerator_eN_lowEnergy/init_2Pi [ Subroutines ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Subroutines ]

NAME

subroutine init_2Pi(eN,OutPart,XS)

PURPOSE

Generate a e^- N -> e^- N pion pion event.

NOTES

This is done by taking the (vacuum) 2pi-background in photoproduction and scaling it with Q2 as the total XS. (many better scalings possible!)

returns XS = dsigma/(dE' dOmega) in mb/GeV

Restricted to Q^2 < 5 GeV^2 and W < 2 GeV, otherwise set to zero.


eventGenerator_eN_lowEnergy/init_2Pi_getBG [ Subroutines ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Subroutines ]

NAME

subroutine init_2Pi_getBG(nucleon_free, Wfree, sig2pi)

PURPOSE

Calculate the background contribution for an e^- N -> e^- N pion pion event at Q^2=0.

Since te resonance cross section is only calculated up to W=2 GeV, it is also only possible to get some 2pi BG up to this value.

NOTES

This is a helper routine for init_2Pi, but also used in the neutrino case


eventGenerator_eN_lowEnergy/init_1Pi [ Subroutines ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Subroutines ]

NAME

subroutine init_1Pi(eN,OutPart,XS,theta_k,phi_k,modeBckGrnd)

PURPOSE

Generate one e^- N -> e^- N pion event.

The weight of the each event is given by (total Xsection of event=sum over all possible pion charges). We make a Monte Carlo decision to determine the pion charge.

The particles are produced at the place of the nucleon target.

If one or more of the input angles theta_k, phi_k are negative, then we make a Monte-Carlo decision for those angles which are negative. In this procedure, we distribute:

  • phi_k flat in [0,2*pi]
  • cos(theta_k) flat in [-1,1]

Note that those angles are defined in the CM-frame of the outgoing pion and nucleon.

If a Monte-Carlo decision on the angle is performed, then the perweight of each event includes the following integral measure:

  • phi and theta decision: measure=int dphi_k dtheta_k=4*pi --> perweight=4*pi*dsigma_dOmega_pion(phi_k,theta_k)
  • phi integration : measure=int dphi_k dtheta_k=2*pi --> perweight=2*pi*dsigma_dOmega_pion(phi_k,theta_k
  • theta integration : measure=int dtheta_k =2 --> perweight=2 *dsigma_dOmega_pion(phi_k,theta_k

INPUTS

  • type(electronNucleon_event) :: eN -- The underlying electron and nucleon event
  • real :: theta_k,phi_k -- Outgoing pion angles (if negative then Monte Carlo Integration). The angles are defined in the CM-Frame of the outgoing pion and nucleon
  • logical :: modeBckGrnd -- true: mode==background, false: mode=normal

OUPTUT * type(particle), dimension(1:2) :: OutPart * logical :: flagOK


eventGenerator_eN_lowEnergy/init_QE [ Subroutines ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Subroutines ]

NAME

subroutine init_QE(eN,OutPart,XS)

PURPOSE

Initializes one e^- N -> e^- N' events. The weight of the event is given by the Xsection for QE scattering.

The particles are produced at the place of the nucleon target.

INPUTS

OUPTUT * real :: XS - The resulting cross section (=0 for failure) * type(particle) :: OutPart


eventGenerator_eN_lowEnergy/init_Res [ Subroutines ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Subroutines ]

NAME

subroutine init_Res(eN,OutPart,XS)

PURPOSE

Initializes an e^- N -> Resonance event. The (per)weight of the event is given by (total Xsection of event=sum over all resonances). We make a Monte Carlo decision to determine the resonance type.

The particles are produced at the place of the nucleon target.

INPUTS

  • type(electronNucleon_event) :: eN -- The underlying electron and nucleon event
  • logical, dimension(2:nres+1) :: useRes -- Switch on/off each resonance

OUPTUT * real :: XS

   -- The resulting cross section (=0 for failure)

* type(particle) :: OutPart


eventGenerator_eN_lowEnergy/generateEvent_1Pi [ Subroutines ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Subroutines ]

NAME

subroutine generateEvent_1Pi(OutPart,initNuc,kf,pf,xSection,pionCharge)

PURPOSE

Given the evaluated kinematics and cross sections, final state particles are initialized for pion nucleon production.

INPUTS

  • integer :: pionCharge -- Charge of outgoing pion
  • type(particle) :: initNuc -- initial Nucleon
  • real, dimension(0:3) :: kf,pf -- pion momentum and nucleon momentum
  • real :: xSection -- Xsection for producing this event, including e.g. d(Omega)

OUPTUT * type(particle),dimension(1:2) :: OutPart -- final state particles


eventGenerator_eN_lowEnergy/generateEvent_1Body [ Subroutines ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Subroutines ]

NAME

subroutine generateEvent_1Body(OutPart,initNuc,pf,xSection,ID,mass)

PURPOSE

Given the evaluated kinematics and cross sections, final state particles are initialized for 1-body final states.

INPUTS

  • integer :: ID -- ID of produced particle
  • real :: mass -- mass of produced particle
  • type(particle) :: initNuc -- initial Nucleon
  • real, dimension(0:3) :: pf -- final state
  • real :: xSection -- Xsection for producing this event, including e.g. d(Omega)

OUPTUT * type(particle), intent(out) :: OutPart -- finalstate particle


eventGenerator_eN_lowEnergy/setOutPartDefaults [ Subroutines ]

[ Top ] [ eventGenerator_eN_lowEnergy ] [ Subroutines ]

NAME

subroutine setOutPartDefaults(OutPart, initNuc)

PURPOSE

set some fields of the outgoing particles to default values