gibuu is hosted by Hepforge, IPPP Durham
GiBUU

TABLE OF CONTENTS


/formFactor_ResProd [ Modules ]

[ Top ] [ Modules ]

NAME

module formFactor_ResProd

PURPOSE

Provides the form factors for electro-weak resonance excitation.

One might choose whether the form factors are calculated from MAID's helicity amplitudes directly or whether the fit of Lalakulich (PRD 74, 014009 (2006)) is used.


formFactor_ResProd/debugFlag [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  logical, parameter :: debugFlag=.false.

PURPOSE

Switch for debug information


formFactor_ResProd/FF_ResProd [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  integer, save :: FF_ResProd=0

PURPOSE

select how the form factors are calculated:

  • 0: MAID's helicity amplitudes (Luis' helicity expressions - CM frame)
  • 1: fit of Lalakulich (PRD 74, 014009 (2006))
  • 2: MAID's helicity amplitudes (Lalakulich's helicity expressions - LAB frame)


formFactor_ResProd/DeltaAxFF [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  integer, save :: DeltaAxFF=1

PURPOSE

choose between different axial form factors for the Delta:

  • 1: Adler
  • 2: Paschos
  • 3: dipol


formFactor_ResProd/MA [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  real, save :: MA=1.05   ! 0.95 tuned to ANL  ! 1.3 tuned to BNL

PURPOSE

Delta resonance axial mass parameter. Wilkinson et al have shown in Phys.Rev.D 90 (2014) that the ANL pion production data are more reliable than the BNL ones.


formFactor_ResProd/C5A0corr [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  real, save :: C5A0corr = 0.85

PURPOSE

fit parameter for C_5^A (Adler), adjusts strength of C5A0 for the Delta


formFactor_ResProd/aDelta [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  real, save :: aDelta= 0

PURPOSE

fit parameter for C_5^A (Adler)


formFactor_ResProd/bDelta [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  real, save :: bDelta= 0

PURPOSE

fit parameter for C_5^A (Adler)


formFactor_ResProd/cDelta [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  real, save :: cDelta=3.

PURPOSE

fit parameter for C_5^A (Paschos)


formFactor_ResProd/DeltaCouplrelErr [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  real, save :: DeltaCouplrelErr=0.

PURPOSE

error in percent for C_5^A(0) for the Delta


formFactor_ResProd/HNV_axialFF [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  logical, save :: HNV_axialFF=.false.

PURPOSE

choose axial form factors for the Delta:

  • .true. is Hernandez-Nieves-Valverde fit with C5A=0.867,MA=0.985 (PRD 76)
  • .false. is as it was used by Lalakulich et al in PRD 74


formFactor_ResProd/nenner_C5A_Lalakulich [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  real, save :: nenner_C5A_Lalakulich=3.0

PURPOSE

Factor which appears in the Lalakulich parameterization of the axial C_5^A form factor:

  • 3.0 was fitted to BNL and used in Lalakulich PRD71 and PRD 74
  • 0.5 is given by a fit of ANL


formFactor_ResProd/refit_barnu_axialFF [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  logical, save :: refit_barnu_axialFF=.false.

PURPOSE

if .true., axial form factors are refitted to explain the low value of antineutrino cross section (exp. data Bolognese PLB81,393 (1979))


formFactor_ResProd/W_cutOff_switch [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  logical, save :: W_cutOff_switch=.false.

PURPOSE

Switch to include a W-dependent cut-off function for the vector form factor of the Delta


formFactor_ResProd/W_cutOff_switchAll [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  logical, save :: W_cutOff_switchAll=.false.

PURPOSE

Switch to include a W-dependent cut-off function for the vector and the axial form factor of all resonances

NOTES

For the Delta, the cut off W_cutOff_lambda is used. For all other resonances the pole mass is used, i.e. lambda = mass


formFactor_ResProd/W_cutOff_lambda [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  real, save :: W_cutOff_lambda=1.071

PURPOSE

Value for lambda in the W-dependent cut-off function.

NOTES

cf. https://gibuu.hepforge.org/trac/wiki/Xsections/Compton and J.Weil, PhD thesis, fig. 8


formFactor_ResProd/vector_FF_switch [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  logical, save :: vector_FF_switch=.true.

PURPOSE

Switch to turn off the vector form factors


formFactor_ResProd/axial_FF_switch [ Global module-variables ]

[ Top ] [ formFactor_ResProd ] [ Global module-variables ]

SOURCE

  logical, save :: axial_FF_switch=.true.

PURPOSE

Switch to turn off the axial form factors


formFactor_ResProd/input_FF_ResProd [ Namelists ]

[ Top ] [ formFactor_ResProd ] [ Namelists ]

NAME

NAMELIST /input_FF_ResProd/

PURPOSE

This Namelist for module "formFactor_ResProd" includes:


formFactor_ResProd/getFormfactor_Res [ Functions ]

[ Top ] [ formFactor_ResProd ] [ Functions ]

NAME

function getFormfactor_Res(Qs, bare_mass, resID, targetCharge, process, flagOK)

PURPOSE

This function serves as an interface to the subroutines which calculate the form factors.

INPUTS

  • real :: Qs -- momentum transfer
  • real :: bare_mass -- bare mass of the resonance (no potentials included!)
  • integer :: resID -- ID of resonance (see IDtable)
  • integer :: targetCharge -- charge of target (nucleon)
  • integer :: process -- EM, CC or NC (or anti if negative)
  • logical :: flagOK -- control flag whether form factors have been set successfully

OUTPUT

real, dimension(1:8) :: getFormfactor_Res -- first 4 entries contain the vector form factors, last 4 entries contain the axial form factors


formFactor_ResProd/vecFF_MAID_CM [ Subroutines ]

[ Top ] [ formFactor_ResProd ] [ Subroutines ]

NAME

subroutine vecFF_MAID_CM(FF,Qs,resID,targetCharge,process,flagOK)

PURPOSE

returns the vector form factors calculated using the MAID helicity amplitudes, helicity amplitudes from Luis are used (CM frame!!)

INPUTS

  • real :: Qs -- momentum transfer
  • integer :: resID -- ID of resonance (see IDtable)
  • integer :: targetCharge -- charge of target (nucleon)
  • integer :: process -- EM, CC or NC (or anti if negative)
  • logical :: flagOK -- control flag whether form factors have been set successfully

OUTPUT

  • real, dimension(1:4):: FF -- contains the vector FF


formFactor_ResProd/vecFF_MAIDLala [ Subroutines ]

[ Top ] [ formFactor_ResProd ] [ Subroutines ]

NAME

subroutine vecFF_MAIDLala(FF,Qs,resID,targetCharge,process,flagOK)

PURPOSE

returns the vector form factors calculated using the MAID helicity amplitudes, helicity amplitudes from Lalakulich are used

INPUTS

  • real :: Qs -- momentum transfer
  • integer :: resID -- ID of resonance (see IDtable)
  • integer :: targetCharge -- charge of target (nucleon)
  • integer :: process -- EM, CC or NC (or anti if negative)
  • logical :: flagOK -- controll flag whether form factors have been set successfully

OUTPUT

  • real, dimension(1:4):: FF -- contains the vector FF


formFactor_ResProd/vecFF_Lala [ Subroutines ]

[ Top ] [ formFactor_ResProd ] [ Subroutines ]

NAME

subroutine vecFF_Lala(FF,Qs,resID,targetCharge,process,flagOK)

PURPOSE

returns the vector form factors according to a fit by Lalakulich et al. PRD 74, 014009 (2006).

INPUTS

  • real :: Qs -- momentum transfer
  • integer :: resID -- ID of resonance (see IDtable)
  • integer :: targetCharge -- charge of target (nucleon)
  • integer :: process -- EM, CC or NC (or anti if negative)
  • logical :: flagOK -- controll flag whether form factors have been set successfully

OUTPUT

  • real, dimension(1:4):: FF -- contains the vector FF


formFactor_ResProd/axialFF [ Subroutines ]

[ Top ] [ formFactor_ResProd ] [ Subroutines ]

NAME

subroutine axialFF(FF,Qs,resID,targetCharge,process,flagOK)

PURPOSE

returns the axial form factors, see Tina's notes and Mathematica notebook 'axialff_coupling.nb' for the calculation of the axial couplings

INPUTS

  • real :: Qs -- momentum transfer
  • integer :: resID -- ID of resonance (see IDtable)
  • integer :: targetCharge -- charge of target (nucleon)
  • integer :: process -- EM, CC or NC (or anti if negative)
  • logical :: flagOK -- controll flag whether form factors have been set successfully

OUTPUT

  • real, dimension(1:4):: FF -- contains the axial FF


formFactor_ResProd/axialFF_Lala [ Subroutines ]

[ Top ] [ formFactor_ResProd ] [ Subroutines ]

NAME

subroutine axialFF_Lala(FF,Qs,resID,targetCharge,process,flagOK)

PURPOSE

returns the axial form factors according to a fit by Lalakulich et al. PRD 74, 014009 (2006).

INPUTS

  • real :: Qs -- momentum transfer
  • integer :: resID -- ID of resonance (see IDtable)
  • integer :: targetCharge -- charge of target (nucleon)
  • integer :: process -- EM, CC or NC (or anti if negative)
  • logical :: flagOK -- controll flag whether form factors have been set successfully

OUTPUT

  • real, dimension(1:4):: FF -- contains the axial FF


formFactor_ResProd/ff_resProd_setCutoff [ Subroutines ]

[ Top ] [ formFactor_ResProd ] [ Subroutines ]

NAME

subroutine ff_resProd_setCutoff(x)

INPUTS

real, intent(in) :: x -- new value for lambda in the W-dependence of the form factor

PURPOSE

Sets a new cut-off for the W-dependent form factor, in particular it sets W_cutoff_lambda=x.