GENIEGenerator
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
NievesQELCCPXSec.h
Go to the documentation of this file.
1 //____________________________________________________________________________
2 /*!
3 
4 \class genie::NievesQELCCPXSec
5 
6 \brief Computes neutrino-nucleon(nucleus) QELCC differential cross section
7  with RPA corrections
8  Is a concrete implementation of the XSecAlgorithmI interface. \n
9 
10 \ref Physical Review C 70, 055503 (2004)
11 
12 \author Joe Johnston, University of Pittsburgh
13  Steven Dytman, University of Pittsburgh
14 
15 \created April 2016
16 
17 \cpright Copyright (c) 2003-2024, The GENIE Collaboration
18  For the full text of the license visit http://copyright.genie-mc.org
19 */
20 //____________________________________________________________________________
21 
22 #ifndef _NIEVES_QELCC_CROSS_SECTION_H_
23 #define _NIEVES_QELCC_CROSS_SECTION_H_
24 
28 #include <complex>
29 #include <Math/IFunction.h>
34 
35 namespace genie {
36 
37 typedef enum EQELRmax {
38  // Use the same maximum radius as VertexGenerator (3*R0*A^(1/3))
40 
41  // Use the method for calculting Rmax from Nieves' Fortran code
44 
45 class QELFormFactorsModelI;
46 class XSecIntegratorI;
47 
49 
50 public:
52  NievesQELCCPXSec(string config);
53  virtual ~NievesQELCCPXSec();
54 
55  // XSecAlgorithmI interface implementation
56  double XSec (const Interaction * i, KinePhaseSpace_t k) const;
57  double Integral (const Interaction * i) const;
58  bool ValidProcess (const Interaction * i) const;
59 
60  // Override the Algorithm::Configure methods to load configuration
61  // data to private data members
62  void Configure (const Registry & config);
63  void Configure (string param_set);
64 
65 private:
66  void LoadConfig (void);
67 
71  double fCos8c2; ///< cos^2(cabibbo angle)
72 
73  double fXSecCCScale; ///< external xsec scaling factor for CC
74  double fXSecNCScale; ///< external xsec scaling factor for NC
75  double fXSecEMScale; ///< external xsec scaling factor for EM
76  const QvalueShifter * fQvalueShifter ; ///< Optional algorithm to retrieve the qvalue shift for a given target
77 
78  double fhbarc; ///< hbar*c in GeV*fm
79 
80  // mutable for testing purposes only!
81  mutable bool fRPA; ///< use RPA corrections
82  bool fCoulomb; ///< use Coulomb corrections
83 
84  const NuclearModelI* fNuclModel; ///< Nuclear Model for integration
85  // Detect whether the nuclear model is local Fermi gas, and store
86  // the relativistic Fermi momentum table if not
87  bool fLFG;
89  string fKFTableName;
90 
91  /// Enum specifying the method to use when calculating the binding energy of
92  /// the initial hit nucleon during spline generation
94 
95  /// Cutoff lab-frame probe energy above which the effects of Fermi motion and
96  /// binding energy are ignored when computing the total cross section
97  double fEnergyCutOff;
98 
99  /// Whether to apply Pauli blocking in XSec()
101  /// The PauliBlocker instance to use to apply that correction
103 
104  /// Nuclear radius parameter r = R0*A^(1/3) used to compute the
105  /// maximum radius for integration of the Coulomb potential
106  /// when matching the VertexGenerator method
107  double fR0;
108 
109  /// Scaling factor for the Coulomb potential
111 
112  /// Enum variable describing which method of computing Rmax should be used
113  /// for integrating the Coulomb potential
115 
116  //Functions needed to calculate XSec:
117 
118  // Calculates values of CN, CT, CL, and imU, and stores them in the provided
119  // variables. If target is not a nucleus, then CN, CN, and CL are all 1.0.
120  // r must be in units of fm.
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;
125 
126  //Equations to calculate the relativistic Lindhard function for Amunu
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,
133  bool isNeutrino,
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;
139 
140  // Potential for coulomb correction
141  double vcr(const Target * target, double r) const;
142 
143  //input must be length 4. Returns 1 if input is an even permutation of 0123,
144  //-1 if input is an odd permutation of 0123, and 0 if any two elements
145  //are equal
146  int leviCivita(int input[]) const;
147 
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;
152 
153  // NOTE: THE FOLLOWING CODE IS FOR TESTING PURPOSES ONLY
154  // Used to print tensor elements and various inputs for comparison to Nieves'
155  // fortran code
156  mutable bool fCompareNievesTensors; ///< print tensors
157  mutable TString fTensorsOutFile; ///< file to print tensors to
158  mutable double fVc,fCoulombFactor;
159  void CompareNievesTensors(const Interaction* i) const;
160  // END TESTING CODE
161 };
162 } // genie namespace
163 
164 //____________________________________________________________________________
165 /*!
166 \class genie::utils::gsl::wrap::NievesQELIntegrand
167 
168 \brief Auxiliary scalar function for integration over the nuclear density
169  when calculaing the Coulomb correction in the Nieves QEL xsec model
170 
171 \author Joe Johnston, University of Pittsburgh
172  Steven Dytman, University of Pittsburgh
173 
174 \created June 03, 2016
175 */
176 //____________________________________________________________________________
177 
178 namespace genie {
179  namespace utils {
180  namespace gsl {
181  namespace wrap {
182 
183  class NievesQELvcrIntegrand : public ROOT::Math::IBaseFunctionOneDim
184  {
185  public:
186  NievesQELvcrIntegrand(double Rcurr, int A, int Z);
188  // ROOT::Math::IBaseFunctionOneDim interface
189  unsigned int NDim (void) const;
190  double DoEval (double rin) const;
191  ROOT::Math::IBaseFunctionOneDim * Clone (void) const;
192  private:
193  double fRcurr;
194  double fA;
195  double fZ;
196  };
197 
198  } // wrap namespace
199  } // gsl namespace
200  } // utils namespace
201 } // genie namespace
202 
203 
204 #endif
Cross Section Calculation Interface.
double XSec(const Interaction *i, KinePhaseSpace_t k) const
Compute the cross section for the input interaction.
const QvalueShifter * fQvalueShifter
Optional algorithm to retrieve the qvalue shift for a given target.
std::complex< double > ruLinRelX(double q0, double qm, double kf, double m) const
bool fRPA
use RPA corrections
double fXSecCCScale
external xsec scaling factor for CC
A class holding Quasi Elastic (QEL) Form Factors.
void Configure(const Registry &config)
Cross Section Integrator Interface.
enum genie::EQELRmax Nieves_Coulomb_Rmax_t
QELEvGen_BindingMode_t fIntegralNucleonBindingMode
const FermiMomentumTable * fKFTable
double fXSecEMScale
external xsec scaling factor for EM
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
Examines whether the generated event should be Pauli blocked. Is a concerete implementation of the Ev...
Definition: PauliBlocker.h:36
const XSecIntegratorI * fXSecIntegrator
Pure abstract base class. Defines the NuclearModelI interface to be implemented by any physics model ...
Definition: NuclearModelI.h:46
A table of Fermi momentum constants.
QELFormFactors fFormFactors
double fCoulombScale
Scaling factor for the Coulomb potential.
bool fCoulomb
use Coulomb corrections
enum genie::EKinePhaseSpace KinePhaseSpace_t
enum genie::EQELEvGenBindingMode QELEvGen_BindingMode_t
Nieves_Coulomb_Rmax_t fCoulombRmaxMode
ROOT::Math::IBaseFunctionOneDim * Clone(void) const
std::complex< double > deltaLindhard(double q0gev, double dqgev, double rho, double kFgev) const
const QELFormFactorsModelI * fFormFactorsModel
double fhbarc
hbar*c in GeV*fm
double fCos8c2
cos^2(cabibbo angle)
Summary information for an interaction.
Definition: Interaction.h:56
const NuclearModelI * fNuclModel
Nuclear Model for integration.
Pure abstract base class. Defines the QELFormFactorsModelI interface to be implemented by any algorit...
static constexpr double A
Definition: Units.h:74
bool fCompareNievesTensors
print tensors
double Integral(const Interaction *i) const
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
bool fDoPauliBlocking
Whether to apply Pauli blocking in XSec()
double fXSecNCScale
external xsec scaling factor for NC
Computes neutrino-nucleon(nucleus) QELCC differential cross section with RPA corrections Is a concret...
double vcr(const Target *target, double r) const
TString fTensorsOutFile
file to print tensors to
A registry. Provides the container for algorithm configuration parameters.
Definition: Registry.h:65
int leviCivita(int input[]) const
void CompareNievesTensors(const Interaction *i) const
bool ValidProcess(const Interaction *i) const
Can this cross section algorithm handle the input process?
std::complex< double > relLindhard(double q0gev, double dqgev, double kFgev, double M, bool isNeutrino, std::complex< double > relLindIm) const
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
const genie::PauliBlocker * fPauliBlocker
The PauliBlocker instance to use to apply that correction.
static constexpr double m
Definition: Units.h:71
std::complex< double > relLindhardIm(double q0gev, double dqgev, double kFngev, double kFpgev, double M, bool isNeutrino, double &t0, double &r00) const