NievesQELCCPXSec.cxx
Go to the documentation of this file.
1 //____________________________________________________________________________
2 /*
3  Copyright (c) 2003-2020, The GENIE Collaboration
4  For the full text of the license visit http://copyright.genie-mc.org
5 
6  Joe Johnston, University of Pittsburgh
7  Steven Dytman, University of Pittsburgh
8 */
9 //____________________________________________________________________________
10 
11 #include <TMath.h>
12 #include <TVector3.h>
13 #include <TLorentzVector.h>
14 #include <Math/IFunction.h>
15 #include <Math/Integrator.h>
16 #include <complex>
17 
22 #include "Framework/Conventions/GBuild.h"
42 
43 #include <iostream> // Used for testing code
44 #include <fstream> // Used for testing code
46 
47 using namespace genie;
48 using namespace genie::constants;
49 using namespace genie::controls;
50 using namespace genie::utils;
51 
52 //____________________________________________________________________________
54 XSecAlgorithmI("genie::NievesQELCCPXSec")
55 {
56 
57 }
58 //____________________________________________________________________________
60 XSecAlgorithmI("genie::NievesQELCCPXSec", config)
61 {
62 
63 }
64 //____________________________________________________________________________
66 {
67 
68  }
69 //____________________________________________________________________________
71  KinePhaseSpace_t kps) const
72 {
73  /*// TESTING CODE:
74  // The first time this method is called, output tensor elements and other
75  // kinmeatics variables for various kinematics. This can the be compared
76  // to Nieves' fortran code for validation purposes
77  if(fCompareNievesTensors){
78  LOG("Nieves",pNOTICE) << "Printing tensor elements for specific "
79  << "kinematics for testing purposes";
80  CompareNievesTensors(interaction);
81  fCompareNievesTensors = false;
82  exit(0);
83  }
84  // END TESTING CODE*/
85 
86 
87  if ( !this->ValidProcess (interaction) ) return 0.;
88  if ( !this->ValidKinematics(interaction) ) return 0.;
89 
90  // Get kinematics and init-state parameters
91  const Kinematics & kinematics = interaction -> Kine();
92  const InitialState & init_state = interaction -> InitState();
93  const Target & target = init_state.Tgt();
94 
95  // HitNucMass() looks up the PDGLibrary (on-shell) value for the initial
96  // struck nucleon
97  double mNi = target.HitNucMass();
98 
99  // Hadronic matrix element for CC neutrino interactions should really use
100  // the "nucleon mass," i.e., the mean of the proton and neutrino masses.
101  // This expression would also work for NC and EM scattering (since the
102  // initial and final on-shell nucleon masses would be the same)
103  double mNucleon = ( mNi + interaction->RecoilNucleon()->Mass() ) / 2.;
104 
105  // Create a copy of the struck nucleon 4-momentum that is forced
106  // to be on-shell (this will be needed later for the tensor contraction,
107  // in which the nucleon is treated in this way)
108  double inNucleonOnShellEnergy = std::sqrt( std::pow(mNi, 2)
109  + std::pow(target.HitNucP4().P(), 2) );
110 
111  // The Nieves CCQE model follows the de Forest prescription: free nucleon
112  // (i.e., on-shell) form factors and spinors are used, but an effective
113  // value of the 4-momentum transfer "qTilde" is used when computing the
114  // contraction of the hadronic tensor. See comments in the
115  // FullDifferentialXSec() method of LwlynSmithQELCCPXSec for more details.
116  TLorentzVector inNucleonMomOnShell( target.HitNucP4().Vect(),
117  inNucleonOnShellEnergy );
118 
119  // Get the four kinematic vectors and caluclate GFactor
120  // Create copies of all kinematics, so they can be rotated
121  // and boosted to the nucleon rest frame (because the tensor
122  // constraction below only applies for the initial nucleon
123  // at rest and q in the z direction)
124 
125  // All 4-momenta should already be stored, with the hit nucleon off-shell
126  // as appropriate
127  TLorentzVector* tempNeutrino = init_state.GetProbeP4(kRfLab);
128  TLorentzVector neutrinoMom = *tempNeutrino;
129  delete tempNeutrino;
130  TLorentzVector inNucleonMom = target.HitNucP4();
131  TLorentzVector leptonMom = kinematics.FSLeptonP4();
132  TLorentzVector outNucleonMom = kinematics.HadSystP4();
133 
134  // Apply Pauli blocking if enabled
135  if ( fDoPauliBlocking && target.IsNucleus() && !interaction->TestBit(kIAssumeFreeNucleon) ) {
136  int final_nucleon_pdg = interaction->RecoilNucleonPdg();
137  double kF = fPauliBlocker->GetFermiMomentum(target, final_nucleon_pdg,
138  target.HitNucPosition());
139  double pNf = outNucleonMom.P();
140  if ( pNf < kF ) return 0.;
141  }
142 
143  // Use the lab kinematics to calculate the Gfactor, in order to make
144  // the XSec differential in initial nucleon momentum and energy
145  // Divide by 4.0 because Nieves' conventions for the leptonic and hadronic
146  // tensor contraction differ from LwlynSmith by a factor of 4
147  double Gfactor = kGF2*fCos8c2 / (8.0*kPi*kPi*inNucleonOnShellEnergy
148  *neutrinoMom.E()*outNucleonMom.E()*leptonMom.E()) / 4.0;
149 
150  // Calculate Coulomb corrections
151  double ml = interaction->FSPrimLepton()->Mass();
152  double ml2 = TMath::Power(ml, 2);
153  double coulombFactor = 1.0;
154  double plLocal = leptonMom.P();
155 
156  bool is_neutrino = pdg::IsNeutrino(init_state.ProbePdg());
157  double r = target.HitNucPosition();
158 
159  if ( fCoulomb ) {
160  // Coulomb potential
161  double Vc = vcr(& target, r);
162 
163  // Outgoing lepton energy and momentum including Coulomb potential
164  int sign = is_neutrino ? 1 : -1;
165  double El = leptonMom.E();
166  double pl = leptonMom.P();
167  double ElLocal = El - sign*Vc;
168 
169  if ( ElLocal - ml <= 0. ) {
170  LOG("Nieves", pDEBUG) << "Event should be rejected. Coulomb effects"
171  << " push kinematics below threshold. Returning xsec = 0.0";
172  return 0.0;
173  }
174 
175  // The Coulomb correction factor blows up as pl -> 0. To guard against
176  // unphysically huge corrections here, require that the lepton kinetic energy
177  // (at infinity) is larger than the magnitude of the Coulomb potential
178  // (should be around a few MeV)
179  double KEl = El - ml;
180  if ( KEl <= std::abs(Vc) ) {
181  LOG("Nieves", pDEBUG) << "Outgoing lepton has a very small kinetic energy."
182  << " Protecting against near-singularities in the Coulomb correction"
183  << " factor by returning xsec = 0.0";
184  return 0.0;
185  }
186 
187  // Local value of the lepton 3-momentum magnitude for the Coulomb
188  // correction
189  plLocal = TMath::Sqrt( ElLocal*ElLocal - ml2 );
190 
191  // Correction factor
192  coulombFactor= (plLocal * ElLocal) / (pl * El);
193 
194  }
195 
196  // When computing the contraction of the leptonic and hadronic tensors,
197  // we need to use an effective value of the 4-momentum transfer q.
198  // The energy transfer (q0) needs to be modified to account for the binding
199  // energy of the struck nucleon, while the 3-momentum transfer needs to
200  // be corrected for Coulomb effects.
201  //
202  // See the original Valencia model paper:
203  // https://journals.aps.org/prc/abstract/10.1103/PhysRevC.70.055503
204 
205  double q0Tilde = outNucleonMom.E() - inNucleonMomOnShell.E();
206 
207  // Shift the q0Tilde if required:
208  if( fQvalueShifter ) q0Tilde += q0Tilde * fQvalueShifter->Shift(*interaction) ;
209 
210  // If binding energy effects pull us into an unphysical region, return
211  // zero for the differential cross section
212  if ( q0Tilde <= 0. && target.IsNucleus() && !interaction->TestBit(kIAssumeFreeNucleon) ) return 0.;
213 
214  // Note that we're working in the lab frame (i.e., the rest frame
215  // of the target nucleus). We can therefore use Nieves' explicit
216  // form of the Amunu tensor if we rotate the 3-momenta so that
217  // qTilde is in the +z direction
218  TVector3 neutrinoMom3 = neutrinoMom.Vect();
219  TVector3 leptonMom3 = leptonMom.Vect();
220 
221  TVector3 inNucleonMom3 = inNucleonMom.Vect();
222  TVector3 outNucleonMom3 = outNucleonMom.Vect();
223 
224  // If Coulomb corrections are being used, adjust the lepton 3-momentum used
225  // to get q3VecTilde so that its magnitude matches the local
226  // Coulomb-corrected value calculated earlier. Note that, although the
227  // treatment of Coulomb corrections by Nieves et al. doesn't change the
228  // direction of the lepton 3-momentum, it *does* change the direction of the
229  // 3-momentum transfer, and so the correction should be applied *before*
230  // rotating coordinates into a frame where q3VecTilde lies along the positive
231  // z axis.
232  TVector3 leptonMomCoulomb3 = (! fCoulomb ) ? leptonMom3
233  : plLocal * leptonMom3 * (1. / leptonMom3.Mag());
234  TVector3 q3VecTilde = neutrinoMom3 - leptonMomCoulomb3;
235 
236  // Find the rotation angle needed to put q3VecTilde along z
237  TVector3 zvec(0.0, 0.0, 1.0);
238  TVector3 rot = ( q3VecTilde.Cross(zvec) ).Unit(); // Vector to rotate about
239  // Angle between the z direction and q
240  double angle = zvec.Angle( q3VecTilde );
241 
242  // Handle the edge case where q3VecTilde is along -z, so the
243  // cross product above vanishes
244  if ( q3VecTilde.Perp() == 0. && q3VecTilde.Z() < 0. ) {
245  rot = TVector3(0., 1., 0.);
246  angle = kPi;
247  }
248 
249  // Rotate if the rotation vector is not 0
250  if ( rot.Mag() >= kASmallNum ) {
251 
252  neutrinoMom3.Rotate(angle,rot);
253  neutrinoMom.SetVect(neutrinoMom3);
254 
255  leptonMom3.Rotate(angle,rot);
256  leptonMom.SetVect(leptonMom3);
257 
258  inNucleonMom3.Rotate(angle,rot);
259  inNucleonMom.SetVect(inNucleonMom3);
260  inNucleonMomOnShell.SetVect(inNucleonMom3);
261 
262  outNucleonMom3.Rotate(angle,rot);
263  outNucleonMom.SetVect(outNucleonMom3);
264 
265  }
266 
267  // Calculate q and qTilde
268  TLorentzVector qP4 = neutrinoMom - leptonMom;
269  TLorentzVector qTildeP4(0., 0., q3VecTilde.Mag(), q0Tilde);
270 
271  double Q2 = -1. * qP4.Mag2();
272  double Q2tilde = -1. * qTildeP4.Mag2();
273 
274  // Store Q2tilde in the interaction so that we get the correct
275  // values of the form factors (according to the de Forest prescription)
276  interaction->KinePtr()->SetQ2(Q2tilde);
277 
278  double q2 = -Q2tilde;
279 
280  // Check that q2 < 0 (accounting for rounding errors)
281  if ( q2 >= kASmallNum ) {
282  LOG("Nieves", pWARN) << "q2 >= 0, returning xsec = 0.0";
283  return 0.0;
284  }
285 
286  // Calculate form factors
287  fFormFactors.Calculate( interaction );
288 
289  // Now that the form factors have been calculated, store Q2
290  // in the event instead of Q2tilde
291  interaction->KinePtr()->SetQ2( Q2 );
292 
293  // Do the contraction of the leptonic and hadronic tensors. See the
294  // RPA-corrected expressions for the hadronic tensor elements in appendix A
295  // of Phys. Rev. C 70, 055503 (2004). Note that the on-shell 4-momentum of
296  // the initial struck nucleon should be used in the calculation, as well as
297  // the effective 4-momentum transfer q tilde (corrected for the nucleon
298  // binding energy and Coulomb effects)
299  double LmunuAnumuResult = LmunuAnumu(neutrinoMom, inNucleonMomOnShell,
300  leptonMom, qTildeP4, mNucleon, is_neutrino, target,
301  interaction->TestBit( kIAssumeFreeNucleon ));
302 
303  // Calculate xsec
304  double xsec = Gfactor*coulombFactor*LmunuAnumuResult;
305 
306  // Apply the factor that arises from elimination of the energy-conserving
307  // delta function
308  xsec *= genie::utils::EnergyDeltaFunctionSolutionQEL( *interaction );
309 
310  // Apply given scaling factor
311  xsec *= fXSecScale;
312 
313  LOG("Nieves",pDEBUG) << "TESTING: RPA=" << fRPA
314  << ", Coulomb=" << fCoulomb
315  << ", q2 = " << q2 << ", xsec = " << xsec;
316 
317  //----- The algorithm computes dxsec/dQ2 or kPSQELEvGen
318  // Check whether variable tranformation is needed
319  if ( kps != kPSQELEvGen ) {
320 
321  // Compute the appropriate Jacobian for transformation to the requested
322  // phase space
323  double J = utils::kinematics::Jacobian(interaction, kPSQELEvGen, kps);
324 
325 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
326  LOG("Nieves", pDEBUG)
327  << "Jacobian for transformation to: "
328  << KinePhaseSpace::AsString(kps) << ", J = " << J;
329 #endif
330  xsec *= J;
331  }
332 
333  // Number of scattering centers in the target
334  int nucpdgc = target.HitNucPdg();
335  int NNucl = (pdg::IsProton(nucpdgc)) ? target.Z() : target.N();
336 
337  xsec *= NNucl; // nuclear xsec
338 
339  return xsec;
340 }
341 //____________________________________________________________________________
342 double NievesQELCCPXSec::Integral(const Interaction * in) const
343 {
344  // If we're using the new spline generation method (which integrates
345  // over the kPSQELEvGen phase space used by QELEventGenerator) then
346  // let the cross section integrator do all of the work. It's smart
347  // enough to handle free nucleon vs. nuclear targets, different
348  // nuclear models (including the local Fermi gas model), etc.
349  if ( fXSecIntegrator->Id().Name() == "genie::NewQELXSec" ) {
350  return fXSecIntegrator->Integrate(this, in);
351  }
352  else {
353  LOG("Nieves", pFATAL) << "Splines for the Nieves CCQE model must be"
354  << " generated using genie::NewQELXSec";
355  std::exit(1);
356  }
357 }
358 //____________________________________________________________________________
360 {
361  if(interaction->TestBit(kISkipProcessChk)) return true;
362 
363  const InitialState & init_state = interaction->InitState();
364  const ProcessInfo & proc_info = interaction->ProcInfo();
365 
366  if(!proc_info.IsQuasiElastic()) return false;
367 
368  int nuc = init_state.Tgt().HitNucPdg();
369  int nu = init_state.ProbePdg();
370 
371  bool isP = pdg::IsProton(nuc);
372  bool isN = pdg::IsNeutron(nuc);
373  bool isnu = pdg::IsNeutrino(nu);
374  bool isnub = pdg::IsAntiNeutrino(nu);
375 
376  bool prcok = proc_info.IsWeakCC() && ((isP&&isnub) || (isN&&isnu));
377  if(!prcok) return false;
378 
379  return true;
380 }
381 //____________________________________________________________________________
383 {
384  Algorithm::Configure(config);
385  this->LoadConfig();
386 }
387 //____________________________________________________________________________
389 {
390  Algorithm::Configure(config);
391  this->LoadConfig();
392 }
393 //____________________________________________________________________________
395 {
396  bool good_config = true ;
397  double thc;
398  GetParam( "CabibboAngle", thc ) ;
399  fCos8c2 = TMath::Power(TMath::Cos(thc), 2);
400 
401  // Cross section scaling factor
402  GetParam( "QEL-CC-XSecScale", fXSecScale ) ;
403 
404  // hbarc for unit conversion, GeV*fm
406 
407  // load QEL form factors model
408  fFormFactorsModel = dynamic_cast<const QELFormFactorsModelI *> (
409  this->SubAlg("FormFactorsAlg") );
410  assert(fFormFactorsModel);
411  fFormFactors.SetModel( fFormFactorsModel ); // <-- attach algorithm
412 
413  // load XSec Integrator
414  fXSecIntegrator = dynamic_cast<const XSecIntegratorI*>(
415  this->SubAlg("XSec-Integrator") );
416  assert(fXSecIntegrator);
417 
418  // Load settings for RPA and Coulomb effects
419 
420  // RPA corrections will not affect a free nucleon
421  GetParamDef("RPA", fRPA, true ) ;
422 
423  // Coulomb Correction- adds a correction factor, and alters outgoing lepton
424  // 3-momentum magnitude (but not direction)
425  // Correction only becomes sizeable near threshold and/or for heavy nuclei
426  GetParamDef( "Coulomb", fCoulomb, true ) ;
427 
428  LOG("Nieves", pNOTICE) << "RPA=" << fRPA << ", useCoulomb=" << fCoulomb;
429 
430  // Get nuclear model for use in Integral()
431  RgKey nuclkey = "IntegralNuclearModel";
432  fNuclModel = dynamic_cast<const NuclearModelI *> (this->SubAlg(nuclkey));
433  assert(fNuclModel);
434 
435  // Check if the model is a local Fermi gas
437 
438  if(!fLFG){
439  // get the Fermi momentum table for relativistic Fermi gas
440  GetParam( "FermiMomentumTable", fKFTableName ) ;
441 
442  fKFTable = 0;
444  fKFTable = kftp->GetTable( fKFTableName );
445  assert( fKFTable );
446  }
447 
448  // TESTING CODE
449  GetParamDef( "PrintDebugData", fCompareNievesTensors, false ) ;
450  // END TESTING CODE
451 
452  // Nuclear radius parameter (R = R0*A^(1/3)) to use when computing
453  // the maximum radius to use to integrate the Coulomb potential
454  GetParam("NUCL-R0", fR0) ; // fm
455 
456  std::string temp_mode;
457  GetParamDef( "RmaxMode", temp_mode, std::string("VertexGenerator") ) ;
458 
459  // Translate the string setting the Rmax mode to the appropriate
460  // enum value, or complain if one couldn't be found
461  if ( temp_mode == "VertexGenerator" ) {
463  }
464  else if ( temp_mode == "Nieves" ) {
466  }
467  else {
468  LOG("Nieves", pFATAL) << "Unrecognized setting \"" << temp_mode
469  << "\" requested for the RmaxMode parameter in the"
470  << " configuration for NievesQELCCPXSec";
471  gAbortingInErr = true;
472  std::exit(1);
473  }
474 
475  // Method to use to calculate the binding energy of the initial hit nucleon when
476  // generating splines
477  std::string temp_binding_mode;
478  GetParamDef( "IntegralNucleonBindingMode", temp_binding_mode, std::string("UseNuclearModel") );
480 
481  // Cutoff energy for integrating over nucleon momentum distribution (above this
482  // lab-frame probe energy, the effects of Fermi motion and binding energy
483  // are taken to be negligible for computing the total cross section)
484  GetParamDef("IntegralNuclearInfluenceCutoffEnergy", fEnergyCutOff, 2.5 ) ;
485 
486  // Get PauliBlocker for possible use in XSec()
487  fPauliBlocker = dynamic_cast<const PauliBlocker*>( this->SubAlg("PauliBlockerAlg") );
488  assert( fPauliBlocker );
489 
490  // Decide whether or not it should be used in XSec()
491  GetParamDef( "DoPauliBlocking", fDoPauliBlocking, true );
492 
493  // Read optional QvalueShifter:
494  fQvalueShifter = nullptr;
495  if( GetConfig().Exists("QvalueShifterAlg") ) {
496  fQvalueShifter = dynamic_cast<const QvalueShifter *> ( this->SubAlg("QvalueShifterAlg") );
497  if( !fQvalueShifter ) {
498  good_config = false ;
499  LOG("NievesQELCCPXSec", pERROR) << "The required QvalueShifterAlg does not exist. AlgID is : " << SubAlg("QvalueShifterAlg")->Id() ;
500  }
501  }
502 
503  if( ! good_config ) {
504  LOG("NievesQELCCPXSec", pERROR) << "Configuration has failed.";
505  exit(78) ;
506  }
507 
508  // Scaling factor for the Coulomb potential
509  GetParamDef( "CoulombScale", fCoulombScale, 1.0 );
510 }
511 //___________________________________________________________________________
512 void NievesQELCCPXSec::CNCTCLimUcalc(TLorentzVector qTildeP4,
513  double M, double r, bool is_neutrino, bool tgtIsNucleus, int tgt_pdgc,
514  int A, int Z, int N, bool hitNucIsProton, double & CN, double & CT, double & CL,
515  double & imaginaryU, double & t0, double & r00, bool assumeFreeNucleon) const
516 {
517  if ( tgtIsNucleus && !assumeFreeNucleon ) {
518  double dq = qTildeP4.Vect().Mag();
519  double dq2 = TMath::Power(dq,2);
520  double q2 = 1 * qTildeP4.Mag2();
521  //Terms for polarization coefficients CN,CT, and CL
522  double hbarc2 = TMath::Power(fhbarc,2);
523  double c0 = 0.380/fhbarc;//Constant for CN in natural units
524 
525  //Density gives the nuclear density, normalized to 1
526  //Input radius r must be in fm
527  double rhop = nuclear::Density(r,A)*Z;
528  double rhon = nuclear::Density(r,A)*N;
529  double rho = rhop + rhon;
530  double rho0 = A*nuclear::Density(0,A);
531 
532  double fPrime = (0.33*rho/rho0+0.45*(1-rho/rho0))*c0;
533 
534  // Get Fermi momenta
535  double kF1, kF2;
536  if(fLFG){
537  if(hitNucIsProton){
538  kF1 = TMath::Power(3*kPi2*rhop, 1.0/3.0) *fhbarc;
539  kF2 = TMath::Power(3*kPi2*rhon, 1.0/3.0) *fhbarc;
540  }else{
541  kF1 = TMath::Power(3*kPi2*rhon, 1.0/3.0) *fhbarc;
542  kF2 = TMath::Power(3*kPi2*rhop, 1.0/3.0) *fhbarc;
543  }
544  }else{
545  if(hitNucIsProton){
546  kF1 = fKFTable->FindClosestKF(tgt_pdgc, kPdgProton);
547  kF2 = fKFTable->FindClosestKF(tgt_pdgc, kPdgNeutron);
548  }else{
549  kF1 = fKFTable->FindClosestKF(tgt_pdgc, kPdgNeutron);
550  kF2 = fKFTable->FindClosestKF(tgt_pdgc, kPdgProton);
551  }
552  }
553 
554  double kF = TMath::Power(1.5*kPi2*rho, 1.0/3.0) *fhbarc;
555 
556  std::complex<double> imU(relLindhardIm(qTildeP4.E(),dq,kF1,kF2,
557  M,is_neutrino,t0,r00));
558 
559  imaginaryU = imag(imU);
560 
561  std::complex<double> relLin(0,0),udel(0,0);
562 
563  // By comparison with Nieves' fortran code
564  if(imaginaryU < 0.){
565  relLin = relLindhard(qTildeP4.E(),dq,kF,M,is_neutrino,imU);
566  udel = deltaLindhard(qTildeP4.E(),dq,rho,kF);
567  }
568  std::complex<double> relLinTot(relLin + udel);
569 
570  /* CRho = 2
571  DeltaRho = 2500 MeV, (2.5 GeV)^2 = 6.25 GeV^2
572  mRho = 770 MeV, (0.770 GeV)^2 = 0.5929 GeV^2
573  g' = 0.63 */
574  double Vt = 0.08*4*kPi/kPionMass2 *
575  (2* TMath::Power((6.25-0.5929)/(6.25-q2),2)*dq2/(q2-0.5929) + 0.63);
576  /* f^2/4/Pi = 0.08
577  DeltaSubPi = 1200 MeV, (1.2 GeV)^2 = 1.44 GeV^2
578  g' = 0.63 */
579  double Vl = 0.08*4*kPi/kPionMass2 *
580  (TMath::Power((1.44-kPionMass2)/(1.44-q2),2)*dq2/(q2-kPionMass2)+0.63);
581 
582  CN = 1.0/TMath::Power(abs(1.0-fPrime*relLin/hbarc2),2);
583 
584  CT = 1.0/TMath::Power(abs(1.0-relLinTot*Vt),2);
585  CL = 1.0/TMath::Power(abs(1.0-relLinTot*Vl),2);
586  }
587  else {
588  //Polarization Coefficients: all equal to 1.0 for free nucleon
589  CN = 1.0;
590  CT = 1.0;
591  CL = 1.0;
592  imaginaryU = 0.0;
593  }
594 }
595 //____________________________________________________________________________
596 // Gives the imaginary part of the relativistic lindhard function in GeV^2
597 // and sets the values of t0 and r00
598 std::complex<double> NievesQELCCPXSec::relLindhardIm(double q0, double dq,
599  double kFn, double kFp,
600  double M,
601  bool isNeutrino,
602  double & t0,
603  double & r00) const
604 {
605  double M2 = TMath::Power(M,2);
606  double EF1,EF2;
607  if(isNeutrino){
608  EF1 = TMath::Sqrt(M2+TMath::Power(kFn,2)); //EFn
609  EF2 = TMath::Sqrt(M2+TMath::Power(kFp,2)); //EFp
610  }else{
611  EF1 = TMath::Sqrt(M2+TMath::Power(kFp,2)); //EFp
612  EF2 = TMath::Sqrt(M2+TMath::Power(kFn,2)); //EFn
613  }
614 
615  double q2 = TMath::Power(q0,2) - TMath::Power(dq,2);
616  double a = (-q0+dq*TMath::Sqrt(1-4.0*M2/q2))/2.0;
617  double epsRP = TMath::Max(TMath::Max(M,EF2-q0),a);
618 
619  // Other theta functions for q are handled by nuclear suppression
620  // That is, q0>0 and -q2>0 are always handled, and q0>EF2-EF1 is
621  // handled if pauli blocking is on, because otherwise the final
622  // nucleon would be below the fermi sea
623  //if(fNievesSuppression && !interaction->TestBit(kIAssumeFreeNucleon )
624  //&& !EF1-epsRP<0){
625  //LOG("Nieves", pINFO) << "Average value of E(p) above Fermi sea";
626  //return 0;
627  //}else{
628  t0 = 0.5*(EF1+epsRP);
629  r00 = (TMath::Power(EF1,2)+TMath::Power(epsRP,2)+EF1*epsRP)/3.0;
630  std::complex<double> result(0.0,-M2/2.0/kPi/dq*(EF1-epsRP));
631  return result;
632  //}
633 }
634 //____________________________________________________________________________
635 //Following obtained from fortran code by J Nieves, which contained the following comment:
636 /*
637  NUCLEON relativistic Lindhard Function
638  Same normalization as ULIN
639  Real part
640  taken from Eur.Phys.J.A25:299-318,2005 (Barbaro et al)
641  Eq. 61
642 
643  Im. part: Juan.
644 
645  INPUT: Real*8
646  q0:Energy [fm]
647  qm: modulus 3mom [fm]
648  kf: Fermi mom [fm]
649 
650  OUTPUT: Complex*16 [fm]
651 
652  USES: ruLinRelX, relLindhardIm
653  */
654 //Takes inputs in GeV (with imU in GeV^2), and gives output in GeV^2
655 std::complex<double> NievesQELCCPXSec::relLindhard(double q0gev,
656  double dqgev, double kFgev, double M,
657  bool isNeutrino,
658  std::complex<double> /* relLindIm */) const
659 {
660  double q0 = q0gev/fhbarc;
661  double qm = dqgev/fhbarc;
662  double kf = kFgev/fhbarc;
663  double m = M/fhbarc;
664 
665  if(q0>qm){
666  LOG("Nieves", pWARN) << "relLindhard() failed";
667  return 0.0;
668  }
669 
670  std::complex<double> RealLinRel(ruLinRelX(q0,qm,kf,m)+ruLinRelX(-q0,qm,kf,m));
671  double t0,r00;
672  std::complex<double> ImLinRel(relLindhardIm(q0gev,dqgev,kFgev,kFgev,M,isNeutrino,t0,r00));
673  //Units of GeV^2
674  return(RealLinRel*TMath::Power(fhbarc,2) + 2.0*ImLinRel);
675 }
676 //____________________________________________________________________________
677 //Inputs assumed to be in natural units
678 std::complex<double> NievesQELCCPXSec::ruLinRelX(double q0, double qm,
679  double kf, double m) const
680 {
681  double q02 = TMath::Power(q0, 2);
682  double qm2 = TMath::Power(qm, 2);
683  double kf2 = TMath::Power(kf, 2);
684  double m2 = TMath::Power(m, 2);
685  double m4 = TMath::Power(m, 4);
686 
687  double ef = TMath::Sqrt(m2+kf2);
688  double q2 = q02-qm2;
689  double q4 = TMath::Power(q2,2);
690  double ds = TMath::Sqrt(1.0-4.0*m2/q2);
691  double L1 = log((kf+ef)/m);
692  std::complex<double> uL2(
693  TMath::Log(TMath::Abs(
694  (ef + q0 - TMath::Sqrt(m2+TMath::Power(kf-qm,2)))/
695  (ef + q0 - TMath::Sqrt(m2 + TMath::Power(kf + qm,2))))) +
696  TMath::Log(TMath::Abs(
697  (ef + q0 + TMath::Sqrt(m2 + TMath::Power(kf - qm,2)))/
698  (ef + q0 + TMath::Sqrt(m2 + TMath::Power(kf + qm,2))))));
699 
700  std::complex<double> uL3(
701  TMath::Log(TMath::Abs((TMath::Power(2*kf + q0*ds,2)-qm2)/
702  (TMath::Power(2*kf - q0*ds,2)-qm2))) +
703  TMath::Log(TMath::Abs((TMath::Power(kf-ef*ds,2) - (4*m4*qm2)/q4)/
704  (TMath::Power(kf+ef*ds,2) - (4*m4*qm2)/q4))));
705 
706  std::complex<double> RlinrelX(-L1/(16.0*kPi2)+
707  uL2*(2.0*ef+q0)/(32.0*kPi2*qm)-
708  uL3*ds/(64.0*kPi2));
709 
710  return RlinrelX*16.0*m2;
711 }
712 //____________________________________________________________________________
713 //Following obtained from fortran code by J Nieves, which contained the following comment:
714 /*
715  complex Lindhard function for symmetric nuclear matter:
716  from Appendix of
717  E.Oset et al Phys. Rept. 188:79, 1990
718  formula A.4
719 
720  input variables:
721  q_zero [fm^-1] : Energy
722  q_mod [fm^-1] : Momentum
723  rho [fm^3] : Nuclear density
724  k_fermi[fm^-1] : Fermi momentum
725 
726  All variables are real*8
727 
728  output variable:
729  delta_lind [fm^-2]
730 
731  ATTENTION!!!
732  Only works properly for real q_zero,
733  if q_zero has an imaginary part calculates the L. function
734  assuming Gamma= 0.
735  Therefore this subroutine provides two different functions
736  depending on whether q_zero is real or not!!!!!!!!!!!
737 */
738 std::complex<double> NievesQELCCPXSec::deltaLindhard(double q0,
739  double dq, double rho, double kF) const
740 {
741  double q_zero = q0/fhbarc;
742  double q_mod = dq/fhbarc;
743  double k_fermi = kF/fhbarc;
744  //Divide by hbarc in order to use natural units (rho is already in the correct units)
745 
746  //m = 939/197.3, md = 1232/197.3, mpi = 139/197.3
747  double m = 4.7592;
748  double md = 6.2433;
749  double mpi = 0.7045;
750 
751  double fdel_f = 2.13;
752  double wr = md-m;
753  double gamma = 0;
754  double gammap = 0;
755 
756  double q_zero2 = TMath::Power(q_zero, 2);
757  double q_mod2 = TMath::Power(q_mod, 2);
758  double k_fermi2 = TMath::Power(k_fermi, 2);
759 
760  double m2 = TMath::Power(m, 2);
761  double m4 = TMath::Power(m, 4);
762  double mpi2 = TMath::Power(mpi, 2);
763  double mpi4 = TMath::Power(mpi, 4);
764 
765  double fdel_f2 = TMath::Power(fdel_f, 2);
766 
767  //For the current code q_zero is always real
768  //If q_zero can have an imaginary part then only the real part is used
769  //until z and zp are calculated
770 
771  double s = m2+q_zero2-q_mod2+
772  2.0*q_zero *TMath::Sqrt(m2+3.0/5.0*k_fermi2);
773 
774  if(s>TMath::Power(m+mpi,2)){
775  double srot = TMath::Sqrt(s);
776  double qcm = TMath::Sqrt(TMath::Power(s,2)+mpi4+m4-2.0*(s*mpi2+s*m2+
777  mpi2*m2)) /(2.0*srot);
778  gamma = 1.0/3.0 * 1.0/(4.0*kPi) * fdel_f2*
779  TMath::Power(qcm,3)/srot*(m+TMath::Sqrt(m2+TMath::Power(qcm,2)))/mpi2;
780  }
781  double sp = m2+q_zero2-q_mod2-
782  2.0*q_zero *TMath::Sqrt(m2+3.0/5.0*k_fermi2);
783 
784 
785  if(sp > TMath::Power(m+mpi,2)){
786  double srotp = TMath::Sqrt(sp);
787  double qcmp = TMath::Sqrt(TMath::Power(sp,2)+mpi4+m4-2.0*(sp*mpi2+sp*m2+
788  mpi2*m2))/(2.0*srotp);
789  gammap = 1.0/3.0 * 1.0/(4.0*kPi) * fdel_f2*
790  TMath::Power(qcmp,3)/srotp*(m+TMath::Sqrt(m2+TMath::Power(qcmp,2)))/mpi2;
791  }
792  //}//End if statement
793  const std::complex<double> iNum(0,1.0);
794 
795  std::complex<double> z(md/(q_mod*k_fermi)*(q_zero-q_mod2/(2.0*md)
796  -wr +iNum*gamma/2.0));
797  std::complex<double> zp(md/(q_mod*k_fermi)*(-q_zero-q_mod2/(2.0*md)
798  -wr +iNum*gammap/2.0));
799 
800  std::complex<double> pzeta(0.0);
801  if(abs(z) > 50.0){
802  pzeta = 2.0/(3.0*z)+2.0/(15.0*z*z*z);
803  }else if(abs(z) < TMath::Power(10.0,-2)){
804  pzeta = 2.0*z-2.0/3.0*z*z*z-iNum*kPi/2.0*(1.0-z*z);
805  }else{
806  pzeta = z + (1.0-z*z) * log((z+1.0)/(z-1.0))/2.0;
807  }
808 
809  std::complex<double> pzetap(0);
810  if(abs(zp) > 50.0){
811  pzetap = 2.0/(3.0*zp)+2.0/(15.0*zp*zp*zp);
812  }else if(abs(zp) < TMath::Power(10.0,-2)){
813  pzetap = 2.0*zp-2.0/3.0*zp*zp*zp-iNum*kPi/2.0*(1.0-zp*zp);
814  }else{
815  pzetap = zp+ (1.0-zp*zp) * log((zp+1.0)/(zp-1.0))/2.0;
816  }
817 
818  //Multiply by hbarc^2 to give answer in units of GeV^2
819  return 2.0/3.0 * rho * md/(q_mod*k_fermi) * (pzeta +pzetap) * fdel_f2 *
820  TMath::Power(fhbarc,2);
821 }
822 
823 //____________________________________________________________________________
824 // Gives coulomb potential in units of GeV
825 double NievesQELCCPXSec::vcr(const Target * target, double Rcurr) const{
826  if(target->IsNucleus()){
827  int A = target->A();
828  int Z = target->Z();
829  double Rmax = 0.;
830 
831  if ( fCoulombRmaxMode == kMatchNieves ) {
832  // Rmax calculated using formula from Nieves' fortran code and default
833  // charge and neutron matter density parameters from NuclearUtils.cxx
834  if (A > 20) {
835  double c = TMath::Power(A,0.35), z = 0.54;
836  Rmax = c + 9.25*z;
837  }
838  else {
839  // c = 1.75 for A <= 20
840  Rmax = TMath::Sqrt(20.0)*1.75;
841  }
842  }
844  // TODO: This solution is fragile. If the formula used by VertexGenerator
845  // changes, then this one will need to change too. Switch to using
846  // a common function to get Rmax for both.
847  Rmax = 3. * fR0 * std::pow(A, 1./3.);
848  }
849  else {
850  LOG("Nieves", pFATAL) << "Unrecognized setting for fCoulombRmaxMode encountered"
851  << " in NievesQELCCPXSec::vcr()";
852  gAbortingInErr = true;
853  std::exit(1);
854  }
855 
856  if(Rcurr >= Rmax){
857  LOG("Nieves",pNOTICE) << "Radius greater than maximum radius for coulomb corrections."
858  << " Integrating to max radius.";
859  Rcurr = Rmax;
860  }
861 
862  ROOT::Math::IBaseFunctionOneDim * func = new
866 
867  double abstol = 1; // We mostly care about relative tolerance;
868  double reltol = 1E-4;
869  int nmaxeval = 100000;
870  ROOT::Math::Integrator ig(*func,ig_type,abstol,reltol,nmaxeval);
871  double result = ig.Integral(0,Rmax);
872  delete func;
873 
874  // Multiply by Z to normalize densities to number of protons
875  // Multiply by hbarc to put result in GeV instead of fm
876  // Multiply by an extra configurable scaling factor that defaults to unity
877  return -kAem*4*kPi*result*fhbarc*fCoulombScale;
878  }else{
879  // If target is not a nucleus the potential will be 0
880  return 0.0;
881  }
882 }
883 //____________________________________________________________________________
885  int copy[4] = {input[0],input[1],input[2],input[3]};
886  int permutations = 0;
887  int temp;
888 
889  for(int i=0;i<4;i++){
890  for(int j=i+1;j<4;j++){
891  //If any two elements are equal return 0
892  if(input[i] == input[j])
893  return 0;
894  //If a larger number is before a smaller one, use permutations
895  //(exchanges of two adjacent elements) to move the smaller element
896  //so it is before the larger element, eg 2341->2314->2134->1234
897  if(copy[i]>copy[j]){
898  temp = copy[j];
899  for(int k=j;k>i;k--){
900  copy[k] = copy[k-1];
901  permutations++;
902  }
903  copy[i] = temp;
904  }
905  }
906  }
907 
908  if(permutations % 2 == 0){
909  return 1;
910  }else{
911  return -1;
912  }
913 }
914 //____________________________________________________________________________
915 // Calculates the constraction of the leptonic and hadronic tensors. The
916 // expressions used here are valid in a frame in which the
917 // initial nucleus is at rest, and qTilde must be in the z direction.
918 double NievesQELCCPXSec::LmunuAnumu(const TLorentzVector neutrinoMom,
919 const TLorentzVector inNucleonMomOnShell, const TLorentzVector leptonMom,
920 const TLorentzVector qTildeP4, double M, bool is_neutrino,
921 const Target& target, bool assumeFreeNucleon) const
922 {
923  double r = target.HitNucPosition();
924  bool tgtIsNucleus = target.IsNucleus();
925  int tgt_pdgc = target.Pdg();
926  int A = target.A();
927  int Z = target.Z();
928  int N = target.N();
929  bool hitNucIsProton = pdg::IsProton( target.HitNucPdg() );
930 
931  const double k[4] = {neutrinoMom.E(),neutrinoMom.Px(),neutrinoMom.Py(),neutrinoMom.Pz()};
932  const double kPrime[4] = {leptonMom.E(),leptonMom.Px(),
933  leptonMom.Py(),leptonMom.Pz()};
934 
935  double q2 = qTildeP4.Mag2();
936 
937  const double q[4] = {qTildeP4.E(),qTildeP4.Px(),qTildeP4.Py(),qTildeP4.Pz()};
938  double q0 = q[0];
939  double dq = TMath::Sqrt(TMath::Power(q[1],2)+
940  TMath::Power(q[2],2)+TMath::Power(q[3],2));
941 
942  int sign = (is_neutrino) ? 1 : -1;
943 
944  // Get the QEL form factors (were calculated before this method was called)
945  double F1V = 0.5*fFormFactors.F1V();
946  double xiF2V = 0.5*fFormFactors.xiF2V();
947  double FA = -fFormFactors.FA();
948  // According to Nieves' paper, Fp = 2.0*M*FA/(kPionMass2-q2), but Llewelyn-
949  // Smith uses Fp = 2.0*M^2*FA/(kPionMass2-q2), so I divide by M
950  // This gives units of GeV^-1
951  double Fp = -1.0/M*fFormFactors.Fp();
952 
953 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
954  LOG("Nieves", pDEBUG) << "\n" << fFormFactors;
955 #endif
956 
957  // Calculate auxiliary parameters
958  double M2 = TMath::Power(M, 2);
959  double FA2 = TMath::Power(FA, 2);
960  double Fp2 = TMath::Power(Fp, 2);
961  double F1V2 = TMath::Power(F1V, 2);
962  double xiF2V2 = TMath::Power(xiF2V, 2);
963  double q02 = TMath::Power(q[0], 2);
964  double dq2 = TMath::Power(dq, 2);
965  double q4 = TMath::Power(q2, 2);
966 
967  double t0,r00;
968  double CN=1.,CT=1.,CL=1.,imU=0;
969  CNCTCLimUcalc(qTildeP4, M, r, is_neutrino, tgtIsNucleus,
970  tgt_pdgc, A, Z, N, hitNucIsProton, CN, CT, CL, imU,
971  t0, r00, assumeFreeNucleon);
972 
973  if ( imU > kASmallNum )
974  return 0.;
975 
976  if ( !fRPA || assumeFreeNucleon ) {
977  CN = 1.0;
978  CT = 1.0;
979  CL = 1.0;
980  }
981 
982  double tulin[4] = {0.,0.,0.,0.};
983  double rulin[4][4] = { {0.,0.,0.,0.},
984  {0.,0.,0.,0.},
985  {0.,0.,0.,0.},
986  {0.,0.,0.,0.} };
987 
988  // TESTING CODE:
990  // Use average values for initial momentum to calculate A, as given
991  // in Appendix B of Nieves' paper. T gives average values of components
992  // of p, and R gives the average value of two components multiplied
993  // together
994  double t3 = (0.5*q2 + q0*t0)/dq; // Average pz
995 
996  // Vector of p
997 
998  tulin[0] = t0;
999  tulin[3] = t3;
1000 
1001  // R is a 4x4 matrix, with R[mu][nu] is the average
1002  // value of p[mu]*p[nu]
1003  double aR = r00-M2;
1004  double bR = (q4+4.0*r00*q02+4.0*q2*q0*t0)/(4.0*dq2);
1005  double gamma = (aR-bR)/2.0;
1006  double delta = (-aR+3.0*bR)/2.0/dq2;
1007 
1008  double r03 = (0.5*q2*t0 + q0*r00)/dq; // Average E(p)*pz
1009 
1010  rulin[0][0] = r00;
1011  rulin[0][3] = r03;
1012  rulin[1][1] = gamma;
1013  rulin[2][2] = gamma;
1014  rulin[3][0] = r03;
1015  rulin[3][3] = gamma+delta*dq2; // END TESTING CODE
1016  }
1017  else {
1018  // For normal code execulation, tulin is the initial nucleon momentum
1019  tulin[0] = inNucleonMomOnShell.E();
1020  tulin[1] = inNucleonMomOnShell.Px();
1021  tulin[2] = inNucleonMomOnShell.Py();
1022  tulin[3] = inNucleonMomOnShell.Pz();
1023 
1024  for(int i=0; i<4; i++)
1025  for(int j=0; j<4; j++)
1026  rulin[i][j] = tulin[i]*tulin[j];
1027  }
1028 
1029  //Additional constants and variables
1030  const int g[4][4] = {{1,0,0,0},{0,-1,0,0},{0,0,-1,0},{0,0,0,-1}};
1031  const std::complex<double> iNum(0,1);
1032  int leviCivitaIndexArray[4];
1033  double imaginaryPart = 0;
1034 
1035  std::complex<double> sum(0.0,0.0);
1036 
1037  double kPrimek = k[0]*kPrime[0]-k[1]*kPrime[1]-k[2]*kPrime[2]-k[3]*kPrime[3];
1038 
1039  std::complex<double> Lmunu(0.0,0.0),Lnumu(0.0,0.0);
1040  std::complex<double> Amunu(0.0,0.0),Anumu(0.0,0.0);
1041 
1042  // Calculate LmunuAnumu by iterating over mu and nu
1043  // In each iteration, add LmunuAnumu to sum if mu=nu, and add
1044  // LmunuAnumu + LnumuAmunu if mu != nu, since we start nu at mu
1045  double axx=0.,azz=0.,a0z=0.,a00=0.,axy=0.;
1046  for(int mu=0;mu<4;mu++){
1047  for(int nu=mu;nu<4;nu++){
1048  imaginaryPart = 0;
1049  if(mu == nu){
1050  //if mu==nu then levi-civita = 0, so imaginary part = 0
1051  Lmunu = g[mu][mu]*kPrime[mu]*g[nu][nu]*k[nu]+g[nu][nu]*kPrime[nu]*g[mu][mu]*k[mu]-g[mu][nu]*kPrimek;
1052  }else{
1053  //if mu!=nu, then g[mu][nu] = 0
1054  //This same leviCivitaIndex array can be used in the else portion when
1055  //calculating Anumu
1056  leviCivitaIndexArray[0] = mu;
1057  leviCivitaIndexArray[1] = nu;
1058  for(int a=0;a<4;a++){
1059  for(int b=0;b<4;b++){
1060  leviCivitaIndexArray[2] = a;
1061  leviCivitaIndexArray[3] = b;
1062  imaginaryPart += - leviCivita(leviCivitaIndexArray)*kPrime[a]*k[b];
1063  }
1064  }
1065  //real(Lmunu) is symmetric, and imag(Lmunu) is antisymmetric
1066  //std::complex<double> num(g[mu][mu]*kPrime[mu]*g[nu][nu]*k[nu]+g[nu][nu]*kPrime[nu]*g[mu][mu]*k[mu],imaginaryPart);
1067  Lmunu = g[mu][mu]*kPrime[mu]*g[nu][nu]*k[nu]+g[nu][nu]*kPrime[nu]*g[mu][mu]*k[mu] + iNum*imaginaryPart;
1068  Lnumu = g[nu][nu]*kPrime[nu]*g[mu][mu]*k[mu]+g[mu][mu]*kPrime[mu]*g[nu][nu]*k[nu ]- iNum*imaginaryPart;
1069  } // End Lmunu calculation
1070 
1071  if(mu ==0 && nu == 0){
1072  Amunu = 16.0*F1V2*(2.0*rulin[0][0]*CN+2.0*q[0]*tulin[0]+q2/2.0)+
1073  2.0*q2*xiF2V2*
1074  (4.0-4.0*rulin[0][0]/M2-4.0*q[0]*tulin[0]/M2-q02*(4.0/q2+1.0/M2)) +
1075  4.0*FA2*(2.0*rulin[0][0]+2.0*q[0]*tulin[0]+(q2/2.0-2.0*M2))-
1076  (2.0*CL*Fp2*q2+8.0*FA*Fp*CL*M)*q02-16.0*F1V*xiF2V*(-q2+q02)*CN;
1077  a00 = real(Amunu); // TESTING CODE
1078  sum += Lmunu*Amunu;
1079  }else if(mu == 0 && nu == 3){
1080  Amunu = 16.0*F1V2*((2.0*rulin[0][3]+tulin[0]*dq)*CN+tulin[3]*q[0])+
1081  2.0*q2*xiF2V2*
1082  (-4.0*rulin[0][3]/M2-2.0*(dq*tulin[0]+q[0]*tulin[3])/M2-dq*q[0]*(4.0/q2+1.0/M2))+
1083  4.0*FA2*((2.0*rulin[0][3]+dq*tulin[0])*CL+q[0]*tulin[3])-
1084  (2.0*CL*Fp2*q2+8.0*FA*Fp*CL*M)*dq*q[0]-
1085  16.0*F1V*xiF2V*dq*q[0];
1086  a0z= real(Amunu); // TESTING CODE
1087  Anumu = Amunu;
1088  sum += Lmunu*Anumu + Lnumu*Amunu;
1089  }else if(mu == 3 && nu == 3){
1090  Amunu = 16.0*F1V2*(2.0*rulin[3][3]+2.0*dq*tulin[3]-q2/2.0)+
1091  2.0*q2*xiF2V2*(-4.0-4.0*rulin[3][3]/M2-4.0*dq*tulin[3]/M2-dq2*(4.0/q2+1.0/M2))+
1092  4.0*FA2*(2.0*rulin[3][3]+2.0*dq*tulin[3]-(q2/2.0-2.0*CL*M2))-
1093  (2.0*CL*Fp2*q2+8.0*FA*Fp*CL*M)*dq2-
1094  16.0*F1V*xiF2V*(q2+dq2);
1095  azz = real(Amunu); // TESTING CODE
1096  sum += Lmunu*Amunu;
1097  }else if(mu ==1 && nu == 1){
1098  Amunu = 16.0*F1V2*(2.0*rulin[1][1]-q2/2.0)+
1099  2.0*q2*xiF2V2*(-4.0*CT-4.0*rulin[1][1]/M2) +
1100  4.0*FA2*(2.0*rulin[1][1]-(q2/2.0-2.0*CT*M2))-
1101  16.0*F1V*xiF2V*CT*q2;
1102  axx = real(Amunu); // TESTING CODE
1103  sum += Lmunu*Amunu;
1104  }else if(mu == 2 && nu == 2){
1105  // Ayy not explicitly listed in paper. This is included so rotating the
1106  // coordinates of k and k' about the z-axis does not change the xsec.
1107  Amunu = 16.0*F1V2*(2.0*rulin[2][2]-q2/2.0)+
1108  2.0*q2*xiF2V2*(-4.0*CT-4.0*rulin[2][2]/M2) +
1109  4.0*FA2*(2.0*rulin[2][2]-(q2/2.0-2.0*CT*M2))-
1110  16.0*F1V*xiF2V*CT*q2;
1111  sum += Lmunu*Amunu;
1112  }else if(mu ==1 && nu == 2){
1113  Amunu = sign*16.0*iNum*FA*(xiF2V+F1V)*(-dq*tulin[0]*CT + q[0]*tulin[3]);
1114  Anumu = -Amunu; // Im(A) is antisymmetric
1115  axy = imag(Amunu); // TESTING CODE
1116  sum += Lmunu*Anumu+Lnumu*Amunu;
1117  }
1118  // All other terms will be 0 because the initial nucleus is at rest and
1119  // qTilde is in the z direction
1120 
1121  } // End loop over nu
1122  } // End loop over mu
1123 
1124  // TESTING CODE
1126  // get tmu
1127  double tmugev = leptonMom.E() - leptonMom.Mag();
1128  // Print Q2, form factors, and tensor elts
1129  std::ofstream ffstream;
1130  ffstream.open(fTensorsOutFile, std::ios_base::app);
1131  if(q0 > 0){
1132  ffstream << -q2 << "\t" << q[0] << "\t" << dq
1133  << "\t" << axx << "\t" << azz << "\t" << a0z
1134  << "\t" << a00 << "\t" << axy << "\t"
1135  << CT << "\t" << CL << "\t" << CN << "\t"
1136  << tmugev << "\t" << imU << "\t"
1137  << F1V << "\t" << xiF2V << "\t"
1138  << FA << "\t" << Fp << "\t"
1139  << tulin[0] << "\t"<< tulin[1] << "\t"
1140  << tulin[2] << "\t"<< tulin[3] << "\t"
1141  << rulin[0][0]<< "\t"<< rulin[0][1]<< "\t"
1142  << rulin[0][2]<< "\t"<< rulin[0][3]<< "\t"
1143  << rulin[1][0]<< "\t"<< rulin[1][1]<< "\t"
1144  << rulin[1][2]<< "\t"<< rulin[1][3]<< "\t"
1145  << rulin[2][0]<< "\t"<< rulin[2][1]<< "\t"
1146  << rulin[2][2]<< "\t"<< rulin[2][3]<< "\t"
1147  << rulin[3][0]<< "\t"<< rulin[3][1]<< "\t"
1148  << rulin[3][2]<< "\t"<< rulin[3][3]<< "\t"
1149  << fVc << "\t" << fCoulombFactor << "\t";
1150  ffstream << "\n";
1151  }
1152  ffstream.close();
1153  }
1154  // END TESTING CODE
1155 
1156  // Since the real parts of A and L are both symmetric and the imaginary
1157  // parts are antisymmetric, the contraction should be real
1158  if ( imag(sum) > kASmallNum )
1159  LOG("Nieves",pWARN) << "Imaginary part of tensor contraction is nonzero "
1160  << "in QEL XSec, real(sum) = " << real(sum)
1161  << "imag(sum) = " << imag(sum);
1162 
1163  return real(sum);
1164 }
1165 
1166 //___________________________________________________________________________
1167 // Auxiliary scalar function for internal integration
1168 //____________________________________________________________________________
1170  double Rcurr, int A, int Z):
1171 ROOT::Math::IBaseFunctionOneDim()
1172 {
1173  fRcurr = Rcurr;
1174  fA = A;
1175  fZ = Z;
1176 }
1177 //____________________________________________________________________________
1179 {
1180 
1181 }
1182 //____________________________________________________________________________
1184 {
1185  return 1;
1186 }
1187 //____________________________________________________________________________
1189 {
1190  double rhop = fZ*nuclear::Density(rin,fA);
1191  if(rin<fRcurr){
1192  return rhop*rin*rin/fRcurr;
1193  }else{
1194  return rhop*rin;
1195  }
1196 }
1197 //____________________________________________________________________________
1198 ROOT::Math::IBaseFunctionOneDim *
1200 {
1202 }
1203 //____________________________________________________________________________
1204 
1205 //____________________________________________________________________________
1206 //
1207 // NOTE: THE REMAINING IS TESTING CODE
1208 //
1209 // This method prints the tensor elements (as well as various inputs) for
1210 // different kinematics. The tensor elements can then be compared to the
1211 // output of Nieves' fortran code.
1212 //
1213 // The results of this code will only agree exactlly with Nieves' fortran
1214 // if Dipole form factors are set (in UserPhysicsOptions).
1215 //
1217  const {
1218  Interaction * interaction = new Interaction(*in); // copy in
1219 
1220  // Set input values here
1221  double ein = 0.2;
1222  double ctl = 0.5;
1223  double rmaxfrac = 0.25;
1224 
1225  bool carbon = false; // true -> C12, false -> Pb208
1226 
1227  if(fRPA)
1228  fTensorsOutFile = "gen.RPA";
1229  else
1230  fTensorsOutFile = "gen.noRPA";
1231 
1232  // Calculate radius
1233  bool klave;
1234  double rp,ap,rn,an;
1235  if(carbon){
1236  klave = true;
1237  rp = 1.692;
1238  ap = 1.082;
1239  rn = 1.692;
1240  an = 1.082;
1241  }else{
1242  // Pb208
1243  klave = false;
1244  rp = 6.624;
1245  ap = 0.549;
1246  rn = 6.890;
1247  an = 0.549;
1248  }
1249  double rmax;
1250  if(!klave)
1251  rmax = TMath::Max(rp,rn) + 9.25*TMath::Max(ap,an);
1252  else
1253  rmax = TMath::Sqrt(20.0)*TMath::Max(rp,rn);
1254  double r = rmax * rmaxfrac;
1255 
1256  // Relevant objects and parameters
1257  //const Kinematics & kinematics = interaction -> Kine();
1258  const InitialState & init_state = interaction -> InitState();
1259  const Target & target = init_state.Tgt();
1260 
1261  // Parameters required for LmunuAnumu
1262  double M = target.HitNucMass();
1263  double ml = interaction->FSPrimLepton()->Mass();
1264  bool is_neutrino = pdg::IsNeutrino(init_state.ProbePdg());
1265 
1266  // Iterate over lepton energy (which then affects q, which is passed to
1267  // LmunuAnumu using in and out NucleonMom
1268  double delta = (ein-0.025)/100.0;
1269  for(int it=0;it<100;it++){
1270  double tmu = it*delta;
1271  double eout = ml + tmu;
1272  double pout = TMath::Sqrt(eout*eout-ml*ml);
1273 
1274  double pin = TMath::Sqrt(ein*ein); // Assume massless neutrinos
1275 
1276  double q0 = ein-eout;
1277  double dq = TMath::Sqrt(pin*pin+pout*pout-2.0*ctl*pin*pout);
1278  double q2 = q0*q0-dq*dq;
1279  interaction->KinePtr()->SetQ2(-q2);
1280 
1281  // When this method is called, inNucleonMomOnShell is unused.
1282  // I can thus provide the calculated values using a null vector for
1283  // inNucleonMomOnShell. I also need to put qTildeP4 in the z direction, as
1284  // Nieves does in his paper.
1285  TLorentzVector qTildeP4(0, 0, dq, q0);
1286  TLorentzVector inNucleonMomOnShell(0,0,0,0);
1287 
1288  // neutrinoMom and leptonMom only directly affect the leptonic tensor, which
1289  // we are not calculating now. Use them to transfer q.
1290  TLorentzVector neutrinoMom(0,0,pout+dq,eout+q0);
1291  TLorentzVector leptonMom(0,0,pout,eout);
1292 
1293  if(fCoulomb){ // Use same steps as in XSec()
1294  // Coulomb potential
1295  double Vc = vcr(& target, r);
1296  fVc = Vc;
1297 
1298  // Outgoing lepton energy and momentum including coulomb potential
1299  int sign = is_neutrino ? 1 : -1;
1300  double El = leptonMom.E();
1301  double ElLocal = El - sign*Vc;
1302  if(ElLocal - ml <= 0.0){
1303  LOG("Nieves",pINFO) << "Event should be rejected. Coulomb effects "
1304  << "push kinematics below threshold";
1305  return;
1306  }
1307  double plLocal = TMath::Sqrt(ElLocal*ElLocal-ml*ml);
1308 
1309  // Correction factor
1310  double coulombFactor= plLocal*ElLocal/leptonMom.Vect().Mag()/El;
1311  fCoulombFactor = coulombFactor; // Store and print
1312  }
1313 
1314  // TODO: apply Coulomb correction to 3-momentum transfer dq
1315 
1316  fFormFactors.Calculate(interaction);
1317  LmunuAnumu(neutrinoMom, inNucleonMomOnShell, leptonMom, qTildeP4,
1318  M, is_neutrino, target, false);
1319  }
1320  return;
1321 } // END TESTING CODE
1322 //____________________________________________________________________________
code to link reconstructed objects back to the MC truth information
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.
Basic constants.
std::complex< double > ruLinRelX(double q0, double qm, double kf, double m) const
ROOT::Math::IntegrationOneDim::Type Integration1DimTypeFromString(string type)
Definition: GSLUtils.cxx:23
bool IsWeakCC(void) const
bool fRPA
use RPA corrections
bool IsNeutrino(int pdgc)
Definition: PDGUtils.cxx:107
void Configure(const Registry &config)
double J(double q0, double q3, double Enu, double ml)
Definition: MECUtils.cxx:147
THE MAIN GENIE PROJECT NAMESPACE
Definition: AlgCmp.h:25
static constexpr double g
Definition: Units.h:144
static QCString result
#define pERROR
Definition: Messenger.h:59
Cross Section Integrator Interface.
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1064
int Type
Definition: 018_def.c:12
int HitNucPdg(void) const
Definition: Target.cxx:304
void SetQ2(double Q2, bool selected=false)
Definition: Kinematics.cxx:255
Kinematics * KinePtr(void) const
Definition: Interaction.h:76
double Density(double r, int A, double ring=0.)
TParticlePDG * RecoilNucleon(void) const
recoil nucleon
QELEvGen_BindingMode_t fIntegralNucleonBindingMode
std::string string
Definition: nybbler.cc:12
const FermiMomentumTable * fKFTable
int RecoilNucleonPdg(void) const
recoil nucleon pdg
int A(void) const
Definition: Target.h:70
void SetModel(const QELFormFactorsModelI *model)
Attach an algorithm.
bool IsQuasiElastic(void) const
Definition: ProcessInfo.cxx:69
double HitNucMass(void) const
Definition: Target.cxx:233
#define pFATAL
Definition: Messenger.h:56
bool IsNucleus(void) const
Definition: Target.cxx:272
constexpr T pow(T x)
Definition: pow.h:72
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
int Pdg(void) const
Definition: Target.h:71
static FermiMomentumTablePool * Instance(void)
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
Examines whether the generated event should be Pauli blocked. Is a concerete implementation of the Ev...
Definition: PauliBlocker.h:36
const XSecIntegratorI * fXSecIntegrator
const TLorentzVector & HadSystP4(void) const
Definition: Kinematics.h:66
Pure abstract base class. Defines the NuclearModelI interface to be implemented by any physics model ...
Definition: NuclearModelI.h:46
QELFormFactors fFormFactors
double fCoulombScale
Scaling factor for the Coulomb potential.
bool fCoulomb
use Coulomb corrections
enum genie::EKinePhaseSpace KinePhaseSpace_t
Nieves_Coulomb_Rmax_t fCoulombRmaxMode
ROOT::Math::IBaseFunctionOneDim * Clone(void) const
virtual const Registry & GetConfig(void) const
Definition: Algorithm.cxx:246
std::complex< double > deltaLindhard(double q0gev, double dqgev, double rho, double kFgev) const
virtual NuclearModel_t ModelType(const Target &) const =0
const QELFormFactorsModelI * fFormFactorsModel
static const double kAem
Definition: Constants.h:56
double fhbarc
hbar*c in GeV*fm
double fCos8c2
cos^2(cabibbo angle)
bool IsNeutron(int pdgc)
Definition: PDGUtils.cxx:338
Summary information for an interaction.
Definition: Interaction.h:56
T abs(T value)
const TLorentzVector & HitNucP4(void) const
Definition: Target.h:91
double GetFermiMomentum(const Target &tgt, int pdg_Nf, double radius=0.0) const
Get the Fermi momentum needed to check Pauli blocking.
virtual bool ValidKinematics(const Interaction *i) const
Is the input kinematical point a physically allowed one?
const NuclearModelI * fNuclModel
Nuclear Model for integration.
bool IsProton(int pdgc)
Definition: PDGUtils.cxx:333
const TLorentzVector & FSLeptonP4(void) const
Definition: Kinematics.h:65
Singleton class to load & serve tables of Fermi momentum constants.
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
Pure abstract base class. Defines the QELFormFactorsModelI interface to be implemented by any algorit...
static int input(void)
Definition: code.cpp:15695
string Name(void) const
Definition: AlgId.h:44
static constexpr double m2
Definition: Units.h:72
const FermiMomentumTable * GetTable(string name)
double EnergyDeltaFunctionSolutionQEL(const Interaction &inter)
Definition: QELUtils.cxx:50
static Config * config
Definition: config.cpp:1054
static string AsString(KinePhaseSpace_t kps)
bool fCompareNievesTensors
print tensors
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
const double a
bool IsAntiNeutrino(int pdgc)
Definition: PDGUtils.cxx:115
double Integral(const Interaction *i) const
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
int ProbePdg(void) const
Definition: InitialState.h:64
bool fDoPauliBlocking
Whether to apply Pauli blocking in XSec()
int Z(void) const
Definition: Target.h:68
static const double kASmallNum
Definition: Controls.h:40
virtual double Shift(const Target &target) const
#define pINFO
Definition: Messenger.h:62
double fXSecScale
external xsec scaling factor
double gamma(double KE, const simb::MCParticle *part)
Misc GENIE control constants.
double xiF2V(void) const
Get the computed form factor xi*F2V.
#define pWARN
Definition: Messenger.h:60
static const double kLightSpeed
Definition: Constants.h:31
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
void Calculate(const Interaction *interaction)
Compute the form factors for the input interaction using the attached model.
double vcr(const Target *target, double r) const
int sign(double val)
Definition: UtilFunc.cxx:104
TString fTensorsOutFile
file to print tensors to
virtual const AlgId & Id(void) const
Get algorithm ID.
Definition: Algorithm.h:97
int N(void) const
Definition: Target.h:69
string RgKey
A registry. Provides the container for algorithm configuration parameters.
Definition: Registry.h:65
const UInt_t kIAssumeFreeNucleon
Definition: Interaction.h:49
double HitNucPosition(void) const
Definition: Target.h:89
QELEvGen_BindingMode_t StringToQELBindingMode(const std::string &mode_str)
Definition: QELUtils.cxx:194
virtual double Integrate(const XSecAlgorithmI *model, const Interaction *interaction) const =0
def func()
Definition: docstring.py:7
T copy(T const &v)
#define A
Definition: memgrp.cpp:38
static bool * b
Definition: config.cpp:1043
double Jacobian(const Interaction *const i, KinePhaseSpace_t f, KinePhaseSpace_t t)
Definition: KineUtils.cxx:130
static constexpr double fermi
Definition: Units.h:55
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
int leviCivita(int input[]) const
void CompareNievesTensors(const Interaction *i) const
const int kPdgProton
Definition: PDGCodes.h:81
double F1V(void) const
Get the computed form factor F1V.
bool ValidProcess(const Interaction *i) const
Can this cross section algorithm handle the input process?
double FindClosestKF(int target_pdgc, int nucleon_pdgc) const
#define pNOTICE
Definition: Messenger.h:61
bool GetParamDef(const RgKey &name, T &p, const T &def) const
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
const Target & Tgt(void) const
Definition: InitialState.h:66
std::complex< double > relLindhard(double q0gev, double dqgev, double kFgev, double M, bool isNeutrino, std::complex< double > relLindIm) const
static const double kGF2
Definition: Constants.h:59
double Fp(void) const
Get the computed form factor Fp.
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
bool gAbortingInErr
Definition: Messenger.cxx:34
const genie::PauliBlocker * fPauliBlocker
The PauliBlocker instance to use to apply that correction.
static const double kPlankConstant
Definition: Constants.h:32
const int kPdgNeutron
Definition: PDGCodes.h:83
static const double kPi
Definition: Constants.h:37
double FA(void) const
Get the computed form factor FA.
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils...
double Density(double r)
Definition: PREM.cxx:18
static QCString * s
Definition: config.cpp:1042
static const double kPi2
Definition: Constants.h:38
Root of GENIE utility namespaces.
TLorentzVector * GetProbeP4(RefFrame_t rf=kRfHitNucRest) const
const UInt_t kISkipProcessChk
if set, skip process validity checks
Definition: Interaction.h:47
std::complex< double > relLindhardIm(double q0gev, double dqgev, double kFngev, double kFpgev, double M, bool isNeutrino, double &t0, double &r00) const
TFile * ef
Definition: doAna.cpp:25
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63
static const double kPionMass2
Definition: Constants.h:86
const Algorithm * SubAlg(const RgKey &registry_key) const
Definition: Algorithm.cxx:345