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.

INPUTS

Via FF_ResProd (namelist "input_FF_ResProd" in the Jobcard) 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

With FF_ResProd (namelist "input_FF_ResProd" in the Jobcard) one can choose 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=0.95   ! 0.95 tuned to ANL  ! 1.3 tuned to BNL

PURPOSE

delta resonance axial mass parameter.


formFactor_ResProd/aDelta [ Global module-variables ]

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

SOURCE

  real, save :: aDelta=-0.25

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

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

With .true. or .false. HNV_axialFF (namelist "input_FF_ResProd" in the Jobcard) one can choose which axial form factors to use for Delta-resonance:

  • .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 wich appear 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 fit of ANL gave 0.5


formFactor_ResProd/refit_barnu_axialFF [ Global module-variables ]

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

SOURCE

  logical, save :: refit_barnu_axialFF=.false.

PURPOSE

With .true. refit_barnu_axialFF (namelist "input_FF_ResProd" in the Jobcard) means that the axial form factors are refitted to explain the low value of antineutrino cross section ( exper 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:

  • false = excluded
  • true = included


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:

  • false = excluded
  • true = included

NOTES

we assume the same dependence as for the Delta vector form factor


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.


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:

  • false = off
  • true = on


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:

  • false = off
  • true = on


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, FF_set)

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 :: FF_set -- controll 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/vectorformfactor_MAID_CM [ Subroutines ]

[ Top ] [ formFactor_ResProd ] [ Subroutines ]

NAME

subroutine vectorformfactor_MAID_CM(vecformfactor,Qs,resID,targetCharge,process,FF_set)

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 :: FF_set -- controll flag whether form factors have been set successfully

OUTPUT

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


formFactor_ResProd/vectorformfactor_MAIDLala [ Subroutines ]

[ Top ] [ formFactor_ResProd ] [ Subroutines ]

NAME

subroutine vectorformfactor_MAIDLala(vecformfactor,Qs,resID,targetCharge,process,FF_set)

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 :: FF_set -- controll flag whether form factors have been set successfully

OUTPUT

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


formFactor_ResProd/vectorformfactor_Lalakulich [ Subroutines ]

[ Top ] [ formFactor_ResProd ] [ Subroutines ]

NAME

subroutine vectorformfactor_Lalakulich(vecformfactor,Qs,resID,targetCharge,process,FF_set)

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 :: FF_set -- controll flag whether form factors have been set successfully

OUTPUT

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


formFactor_ResProd/axialformfactor [ Subroutines ]

[ Top ] [ formFactor_ResProd ] [ Subroutines ]

NAME

subroutine axialformfactor(axialformfactor,Qs,resID,targetCharge,process,FF_set)

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 :: FF_set -- controll flag whether form factors have been set successfully

OUTPUT

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


formFactor_ResProd/axialformfactor_Lalakulich [ Subroutines ]

[ Top ] [ formFactor_ResProd ] [ Subroutines ]

NAME

subroutine axialformfactor_Lalakulich(axialformfactor,Qs,resID,targetCharge,process,FF_set)

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 :: FF_set -- controll flag whether form factors have been set successfully

OUTPUT


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.