22#ifndef _NIEVES_QELCC_CROSS_SECTION_H_
23#define _NIEVES_QELCC_CROSS_SECTION_H_
29#include <Math/IFunction.h>
45class QELFormFactorsModelI;
121 void CNCTCLimUcalc(TLorentzVector qTildeP4,
double M,
double r,
122 bool is_neutrino,
bool tgtIsNucleus,
int tgt_pdgc,
int A,
int Z,
int N,
123 bool hitNucIsProton,
double & CN,
double & CT,
double & CL,
double & imU,
124 double & t0,
double & r00,
bool assumeFreeNucleon)
const;
127 std::complex<double>
relLindhardIm(
double q0gev,
double dqgev,
128 double kFngev,
double kFpgev,
129 double M,
bool isNeutrino,
130 double & t0,
double & r00)
const;
131 std::complex<double>
relLindhard(
double q0gev,
double dqgev,
132 double kFgev,
double M,
134 std::complex<double> relLindIm)
const;
135 std::complex<double>
ruLinRelX(
double q0,
double qm,
136 double kf,
double m)
const;
137 std::complex<double>
deltaLindhard(
double q0gev,
double dqgev,
138 double rho,
double kFgev)
const;
141 double vcr(
const Target * target,
double r)
const;
148 double LmunuAnumu(
const TLorentzVector neutrinoMom,
149 const TLorentzVector inNucleonMom,
const TLorentzVector leptonMom,
150 const TLorentzVector outNucleonMom,
double M,
bool is_neutrino,
151 const Target& target,
bool assumeFreeNucleon)
const;
189 unsigned int NDim (
void)
const;
190 double DoEval (
double rin)
const;
191 ROOT::Math::IBaseFunctionOneDim *
Clone (
void)
const;
A table of Fermi momentum constants.
Summary information for an interaction.
bool fCompareNievesTensors
print tensors
bool fCoulomb
use Coulomb corrections
Nieves_Coulomb_Rmax_t fCoulombRmaxMode
QELEvGen_BindingMode_t fIntegralNucleonBindingMode
int leviCivita(int input[]) const
std::complex< double > relLindhard(double q0gev, double dqgev, double kFgev, double M, bool isNeutrino, std::complex< double > relLindIm) const
void Configure(const Registry &config)
double LmunuAnumu(const TLorentzVector neutrinoMom, const TLorentzVector inNucleonMom, const TLorentzVector leptonMom, const TLorentzVector outNucleonMom, double M, bool is_neutrino, const Target &target, bool assumeFreeNucleon) const
double fCoulombScale
Scaling factor for the Coulomb potential.
void CompareNievesTensors(const Interaction *i) const
const NuclearModelI * fNuclModel
Nuclear Model for integration.
virtual ~NievesQELCCPXSec()
void CNCTCLimUcalc(TLorentzVector qTildeP4, double M, double r, bool is_neutrino, bool tgtIsNucleus, int tgt_pdgc, int A, int Z, int N, bool hitNucIsProton, double &CN, double &CT, double &CL, double &imU, double &t0, double &r00, bool assumeFreeNucleon) const
std::complex< double > ruLinRelX(double q0, double qm, double kf, double m) const
const XSecIntegratorI * fXSecIntegrator
double fXSecNCScale
external xsec scaling factor for NC
double vcr(const Target *target, double r) const
double fhbarc
hbar*c in GeV*fm
bool fDoPauliBlocking
Whether to apply Pauli blocking in XSec()
std::complex< double > deltaLindhard(double q0gev, double dqgev, double rho, double kFgev) const
std::complex< double > relLindhardIm(double q0gev, double dqgev, double kFngev, double kFpgev, double M, bool isNeutrino, double &t0, double &r00) const
bool fRPA
use RPA corrections
bool ValidProcess(const Interaction *i) const
Can this cross section algorithm handle the input process?
const genie::PauliBlocker * fPauliBlocker
The PauliBlocker instance to use to apply that correction.
const QELFormFactorsModelI * fFormFactorsModel
double fXSecEMScale
external xsec scaling factor for EM
const FermiMomentumTable * fKFTable
TString fTensorsOutFile
file to print tensors to
double XSec(const Interaction *i, KinePhaseSpace_t k) const
Compute the cross section for the input interaction.
double fXSecCCScale
external xsec scaling factor for CC
const QvalueShifter * fQvalueShifter
Optional algorithm to retrieve the qvalue shift for a given target.
QELFormFactors fFormFactors
double Integral(const Interaction *i) const
double fCos8c2
cos^2(cabibbo angle)
Pure abstract base class. Defines the NuclearModelI interface to be implemented by any physics model ...
Examines whether the generated event should be Pauli blocked. Is a concerete implementation of the Ev...
A registry. Provides the container for algorithm configuration parameters.
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Cross Section Integrator Interface.
unsigned int NDim(void) const
NievesQELvcrIntegrand(double Rcurr, int A, int Z)
double DoEval(double rin) const
ROOT::Math::IBaseFunctionOneDim * Clone(void) const
THE MAIN GENIE PROJECT NAMESPACE
enum genie::EQELRmax Nieves_Coulomb_Rmax_t
enum genie::EQELEvGenBindingMode QELEvGen_BindingMode_t
enum genie::EKinePhaseSpace KinePhaseSpace_t
@ kMatchVertexGeneratorRmax