gibuu is hosted by Hepforge, IPPP Durham
GiBUU

TABLE OF CONTENTS


/FF_QE_nucleonScattering [ Modules ]

[ Top ] [ Modules ]

NAME

module FF_QE_nucleonScattering

PURPOSE

Provides the electro-weak form factors for the process lepton nucleon -> lepton' nucleon' and the Sachs form factors.


FF_QE_nucleonScattering/parametrization [ Modules ]

[ Top ] [ FF_QE_nucleonScattering ] [ Modules ]

SOURCE

  integer,save :: parametrization=3

PURPOSE

  • 0 = dipole approximation
  • 1 = BBA03 parametrization
  • 2 = BBBA05 parametrization
  • 3 = BBBA07 parametrization


FF_QE_nucleonScattering/useNonStandardMA [ Modules ]

[ Top ] [ FF_QE_nucleonScattering ] [ Modules ]

SOURCE

  logical,save::useNonStandardMA=.false.

PURPOSE

if one wants to use a specific axial mass, set this to true and choose value for MA_in


FF_QE_nucleonScattering/MA_in [ Modules ]

[ Top ] [ FF_QE_nucleonScattering ] [ Modules ]

SOURCE

  real,save :: MA_in=1.0

PURPOSE

axial mass (only if useNonStandardMA=.true.)


FF_QE_nucleonScattering/MV2 [ Modules ]

[ Top ] [ FF_QE_nucleonScattering ] [ Modules ]

SOURCE

  real,save :: MV2=0.71

PURPOSE

vector mass squared in the dipole parametrization of the vector form factors


FF_QE_nucleonScattering/deltas [ Modules ]

[ Top ] [ FF_QE_nucleonScattering ] [ Modules ]

SOURCE

  real,save :: deltas=-0.15

PURPOSE

strange contribution to the axial ff.


FF_QE_nucleonScattering/axialMonopole [ Modules ]

[ Top ] [ FF_QE_nucleonScattering ] [ Modules ]

SOURCE

  logical,save :: axialMonopole=.false.

PURPOSE

use axial ff. of Gari, Kaulfuss PLB 138 (1984)


FF_QE_nucleonScattering/ff_QE [ Namelists ]

[ Top ] [ FF_QE_nucleonScattering ] [ Namelists ]

NAME

NAMELIST ff_QE

PURPOSE

Includes parameters for neutrino matrix elements:


FF_QE_nucleonScattering/formfactors_QE [ Subroutines ]

[ Top ] [ FF_QE_nucleonScattering ] [ Subroutines ]

NAME

subroutine formfactors_QE(QSquared,processID,initialState_charge,F1,F2,FA,FP)

PURPOSE

Provides the electro-weak form factors for the process lepton nucleon -> lepton' nucleon'.

INPUTS

  • real, intent(in) :: QSQuared ! =Q^2 : virtuality of gauge boson in units of GeV**2 =(-1)*Mandelstam-t
  • integer, intent(in) :: initialState_charge ! Specifies the charge of the incoming nucleon (either 1 or 0)
  • integer, intent(in) :: processID

"processID" specifies the reaction type: EM, CC, NC

OUTPUT

  • real, intent(out) :: F1,F2 ! Vector Form factors
  • real, intent(out),optional :: FA,FP ! Axial Form factors (not necessary as input if processID=3)


FF_QE_nucleonScattering/DipoleFF [ Subroutines ]

[ Top ] [ FF_QE_nucleonScattering ] [ Subroutines ]

NAME

subroutine DipoleFF(QSquared, GEp, GMp, GEn, GMn)

PURPOSE

  • Provides the Sachs form factors in dipole approximation.

INPUTS

  • real, intent(in) :: QSquared ! = Q^2 : virtuality of gauge boson in units of GeV**2

OUTPUT

  • real, intent(out) :: GEp, GMp, GEn, GMn ! Sachs Form factors


FF_QE_nucleonScattering/BBA2003 [ Subroutines ]

[ Top ] [ FF_QE_nucleonScattering ] [ Subroutines ]

NAME

subroutine BBA2003(tau, GEp, GMp, GEn, GMn)

PURPOSE

  • Provides the BBA2003 formfactors according to hep-ex/0308005,

valid up to Q**2= 6 GeV**2

INPUTS

  • real, intent(in) :: QSquared ! = Q^2 : virtuality of gauge boson in units of GeV**2

OUTPUT

  • real, intent(out) :: GEp, GMp, GEn, GMn ! Sachs Form factors


FF_QE_nucleonScattering/BBBA2005 [ Subroutines ]

[ Top ] [ FF_QE_nucleonScattering ] [ Subroutines ]

NAME

subroutine BBBA2005(tau, GEp, GMp, GEn, GMn)

PURPOSE

  • Provides the BBBA2005 formfactors according to hep-ex/0602017 v3 (2006)

INPUTS

  • real, intent(in) :: tau

OUTPUT

  • real, intent(out) :: GEp, GMp, GEn, GMn ! Sachs Form factors


FF_QE_nucleonScattering/BBBA2007 [ Subroutines ]

[ Top ] [ FF_QE_nucleonScattering ] [ Subroutines ]

NAME

subroutine BBBA2007(tau, GEp, GMp, GEn, GMn)

PURPOSE

  • Provides the BBBA2007 formfactors according to arXiv:0708.1827 [hep-ex] and arXiv:0708.1946 [hep-ex]
  • c++ file download from www.pas.rochester.edu/~bodek/FF/

INPUTS

  • real, intent(in) :: tau

OUTPUT

  • real, intent(out) :: GEp, GMp, GEn, GMn ! Sachs Form factors