13 #include <TLorentzVector.h> 14 #include <Math/IFunction.h> 15 #include <Math/Integrator.h> 22 #include "Framework/Conventions/GBuild.h" 47 using namespace genie;
91 const Kinematics & kinematics = interaction -> Kine();
92 const InitialState & init_state = interaction -> InitState();
103 double mNucleon = ( mNi + interaction->
RecoilNucleon()->Mass() ) / 2.;
108 double inNucleonOnShellEnergy = std::sqrt(
std::pow(mNi, 2)
116 TLorentzVector inNucleonMomOnShell( target.
HitNucP4().Vect(),
117 inNucleonOnShellEnergy );
128 TLorentzVector neutrinoMom = *tempNeutrino;
130 TLorentzVector inNucleonMom = target.
HitNucP4();
131 TLorentzVector leptonMom = kinematics.
FSLeptonP4();
132 TLorentzVector outNucleonMom = kinematics.
HadSystP4();
139 double pNf = outNucleonMom.P();
140 if ( pNf < kF )
return 0.;
148 *neutrinoMom.E()*outNucleonMom.E()*leptonMom.E()) / 4.0;
152 double ml2 = TMath::Power(ml, 2);
153 double coulombFactor = 1.0;
154 double plLocal = leptonMom.P();
161 double Vc =
vcr(& target, r);
164 int sign = is_neutrino ? 1 : -1;
165 double El = leptonMom.E();
166 double pl = leptonMom.P();
167 double ElLocal = El - sign*Vc;
169 if ( ElLocal - ml <= 0. ) {
170 LOG(
"Nieves",
pDEBUG) <<
"Event should be rejected. Coulomb effects" 171 <<
" push kinematics below threshold. Returning xsec = 0.0";
179 double KEl = El - ml;
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";
189 plLocal = TMath::Sqrt( ElLocal*ElLocal - ml2 );
192 coulombFactor= (plLocal * ElLocal) / (pl * El);
205 double q0Tilde = outNucleonMom.E() - inNucleonMomOnShell.E();
218 TVector3 neutrinoMom3 = neutrinoMom.Vect();
219 TVector3 leptonMom3 = leptonMom.Vect();
221 TVector3 inNucleonMom3 = inNucleonMom.Vect();
222 TVector3 outNucleonMom3 = outNucleonMom.Vect();
232 TVector3 leptonMomCoulomb3 = (!
fCoulomb ) ? leptonMom3
233 : plLocal * leptonMom3 * (1. / leptonMom3.Mag());
234 TVector3 q3VecTilde = neutrinoMom3 - leptonMomCoulomb3;
237 TVector3 zvec(0.0, 0.0, 1.0);
238 TVector3 rot = ( q3VecTilde.Cross(zvec) ).Unit();
240 double angle = zvec.Angle( q3VecTilde );
244 if ( q3VecTilde.Perp() == 0. && q3VecTilde.Z() < 0. ) {
245 rot = TVector3(0., 1., 0.);
252 neutrinoMom3.Rotate(angle,rot);
253 neutrinoMom.SetVect(neutrinoMom3);
255 leptonMom3.Rotate(angle,rot);
256 leptonMom.SetVect(leptonMom3);
258 inNucleonMom3.Rotate(angle,rot);
259 inNucleonMom.SetVect(inNucleonMom3);
260 inNucleonMomOnShell.SetVect(inNucleonMom3);
262 outNucleonMom3.Rotate(angle,rot);
263 outNucleonMom.SetVect(outNucleonMom3);
268 TLorentzVector qP4 = neutrinoMom - leptonMom;
269 TLorentzVector qTildeP4(0., 0., q3VecTilde.Mag(), q0Tilde);
271 double Q2 = -1. * qP4.Mag2();
272 double Q2tilde = -1. * qTildeP4.Mag2();
278 double q2 = -Q2tilde;
282 LOG(
"Nieves",
pWARN) <<
"q2 >= 0, returning xsec = 0.0";
299 double LmunuAnumuResult =
LmunuAnumu(neutrinoMom, inNucleonMomOnShell,
300 leptonMom, qTildeP4, mNucleon, is_neutrino, target,
304 double xsec = Gfactor*coulombFactor*LmunuAnumuResult;
315 <<
", q2 = " << q2 <<
", xsec = " << xsec;
325 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 327 <<
"Jacobian for transformation to: " 353 LOG(
"Nieves",
pFATAL) <<
"Splines for the Nieves CCQE model must be" 354 <<
" generated using genie::NewQELXSec";
376 bool prcok = proc_info.
IsWeakCC() && ((isP&&isnub) || (isN&&isnu));
377 if(!prcok)
return false;
396 bool good_config = true ;
399 fCos8c2 = TMath::Power(TMath::Cos(thc), 2);
409 this->
SubAlg(
"FormFactorsAlg") );
415 this->
SubAlg(
"XSec-Integrator") );
431 RgKey nuclkey =
"IntegralNuclearModel";
461 if ( temp_mode ==
"VertexGenerator" ) {
464 else if ( temp_mode ==
"Nieves" ) {
468 LOG(
"Nieves",
pFATAL) <<
"Unrecognized setting \"" << temp_mode
469 <<
"\" requested for the RmaxMode parameter in the" 470 <<
" configuration for NievesQELCCPXSec";
495 if(
GetConfig().Exists(
"QvalueShifterAlg") ) {
498 good_config = false ;
499 LOG(
"NievesQELCCPXSec",
pERROR) <<
"The required QvalueShifterAlg does not exist. AlgID is : " <<
SubAlg(
"QvalueShifterAlg")->
Id() ;
503 if( ! good_config ) {
504 LOG(
"NievesQELCCPXSec",
pERROR) <<
"Configuration has failed.";
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 517 if ( tgtIsNucleus && !assumeFreeNucleon ) {
518 double dq = qTildeP4.Vect().Mag();
519 double dq2 = TMath::Power(dq,2);
520 double q2 = 1 * qTildeP4.Mag2();
522 double hbarc2 = TMath::Power(
fhbarc,2);
529 double rho = rhop + rhon;
532 double fPrime = (0.33*rho/rho0+0.45*(1-rho/rho0))*c0;
538 kF1 = TMath::Power(3*
kPi2*rhop, 1.0/3.0) *
fhbarc;
539 kF2 = TMath::Power(3*
kPi2*rhon, 1.0/3.0) *
fhbarc;
541 kF1 = TMath::Power(3*
kPi2*rhon, 1.0/3.0) *
fhbarc;
542 kF2 = TMath::Power(3*
kPi2*rhop, 1.0/3.0) *
fhbarc;
554 double kF = TMath::Power(1.5*
kPi2*rho, 1.0/3.0) *
fhbarc;
556 std::complex<double> imU(
relLindhardIm(qTildeP4.E(),dq,kF1,kF2,
557 M,is_neutrino,
t0,r00));
559 imaginaryU = imag(imU);
561 std::complex<double> relLin(0,0),udel(0,0);
565 relLin =
relLindhard(qTildeP4.E(),dq,kF,M,is_neutrino,imU);
568 std::complex<double> relLinTot(relLin + udel);
575 (2* TMath::Power((6.25-0.5929)/(6.25-q2),2)*dq2/(q2-0.5929) + 0.63);
582 CN = 1.0/TMath::Power(
abs(1.0-fPrime*relLin/hbarc2),2);
584 CT = 1.0/TMath::Power(
abs(1.0-relLinTot*Vt),2);
585 CL = 1.0/TMath::Power(
abs(1.0-relLinTot*Vl),2);
599 double kFn,
double kFp,
605 double M2 = TMath::Power(M,2);
608 EF1 = TMath::Sqrt(M2+TMath::Power(kFn,2));
609 EF2 = TMath::Sqrt(M2+TMath::Power(kFp,2));
611 EF1 = TMath::Sqrt(M2+TMath::Power(kFp,2));
612 EF2 = TMath::Sqrt(M2+TMath::Power(kFn,2));
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);
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));
656 double dqgev,
double kFgev,
double M,
658 std::complex<double> )
const 666 LOG(
"Nieves",
pWARN) <<
"relLindhard() failed";
672 std::complex<double> ImLinRel(
relLindhardIm(q0gev,dqgev,kFgev,kFgev,M,isNeutrino,t0,r00));
674 return(RealLinRel*TMath::Power(
fhbarc,2) + 2.0*ImLinRel);
679 double kf,
double m)
const 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);
687 double ef = TMath::Sqrt(m2+kf2);
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))))));
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))));
706 std::complex<double> RlinrelX(-L1/(16.0*
kPi2)+
707 uL2*(2.0*ef+q0)/(32.0*
kPi2*qm)-
710 return RlinrelX*16.0*
m2;
739 double dq,
double rho,
double kF)
const 741 double q_zero = q0/
fhbarc;
743 double k_fermi = kF/
fhbarc;
751 double fdel_f = 2.13;
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);
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);
765 double fdel_f2 = TMath::Power(fdel_f, 2);
771 double s = m2+q_zero2-q_mod2+
772 2.0*q_zero *TMath::Sqrt(m2+3.0/5.0*k_fermi2);
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;
781 double sp = m2+q_zero2-q_mod2-
782 2.0*q_zero *TMath::Sqrt(m2+3.0/5.0*k_fermi2);
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;
793 const std::complex<double> iNum(0,1.0);
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));
800 std::complex<double> pzeta(0.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);
806 pzeta = z + (1.0-z*
z) * log((z+1.0)/(z-1.0))/2.0;
809 std::complex<double> pzetap(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);
815 pzetap = zp+ (1.0-zp*zp) * log((zp+1.0)/(zp-1.0))/2.0;
819 return 2.0/3.0 * rho * md/(q_mod*k_fermi) * (pzeta +pzetap) * fdel_f2 *
835 double c = TMath::Power(A,0.35),
z = 0.54;
840 Rmax = TMath::Sqrt(20.0)*1.75;
850 LOG(
"Nieves",
pFATAL) <<
"Unrecognized setting for fCoulombRmaxMode encountered" 851 <<
" in NievesQELCCPXSec::vcr()";
857 LOG(
"Nieves",
pNOTICE) <<
"Radius greater than maximum radius for coulomb corrections." 858 <<
" Integrating to max radius.";
862 ROOT::Math::IBaseFunctionOneDim *
func =
new 868 double reltol = 1
E-4;
869 int nmaxeval = 100000;
870 ROOT::Math::Integrator ig(*func,ig_type,abstol,reltol,nmaxeval);
871 double result = ig.Integral(0,Rmax);
885 int copy[4] = {input[0],input[1],input[2],input[3]};
886 int permutations = 0;
889 for(
int i=0;i<4;i++){
890 for(
int j=i+1;j<4;j++){
892 if(input[i] == input[j])
899 for(
int k=j;
k>i;
k--){
908 if(permutations % 2 == 0){
919 const TLorentzVector inNucleonMomOnShell,
const TLorentzVector leptonMom,
920 const TLorentzVector qTildeP4,
double M,
bool is_neutrino,
925 int tgt_pdgc = target.
Pdg();
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()};
935 double q2 = qTildeP4.Mag2();
937 const double q[4] = {qTildeP4.E(),qTildeP4.Px(),qTildeP4.Py(),qTildeP4.Pz()};
939 double dq = TMath::Sqrt(TMath::Power(q[1],2)+
940 TMath::Power(q[2],2)+TMath::Power(q[3],2));
942 int sign = (is_neutrino) ? 1 : -1;
953 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 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);
968 double CN=1.,CT=1.,CL=1.,imU=0;
970 tgt_pdgc, A, Z, N, hitNucIsProton, CN, CT, CL, imU,
971 t0, r00, assumeFreeNucleon);
976 if ( !
fRPA || assumeFreeNucleon ) {
982 double tulin[4] = {0.,0.,0.,0.};
983 double rulin[4][4] = { {0.,0.,0.,0.},
994 double t3 = (0.5*q2 + q0*
t0)/dq;
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;
1008 double r03 = (0.5*q2*t0 + q0*r00)/dq;
1012 rulin[1][1] =
gamma;
1013 rulin[2][2] =
gamma;
1015 rulin[3][3] = gamma+delta*dq2;
1019 tulin[0] = inNucleonMomOnShell.E();
1020 tulin[1] = inNucleonMomOnShell.Px();
1021 tulin[2] = inNucleonMomOnShell.Py();
1022 tulin[3] = inNucleonMomOnShell.Pz();
1024 for(
int i=0; i<4; i++)
1025 for(
int j=0; j<4; j++)
1026 rulin[i][j] = tulin[i]*tulin[j];
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;
1035 std::complex<double> sum(0.0,0.0);
1037 double kPrimek = k[0]*kPrime[0]-k[1]*kPrime[1]-k[2]*kPrime[2]-k[3]*kPrime[3];
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);
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++){
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;
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];
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;
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)+
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;
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])+
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];
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);
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;
1104 }
else if(mu == 2 && nu == 2){
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;
1112 }
else if(mu ==1 && nu == 2){
1113 Amunu = sign*16.0*iNum*FA*(xiF2V+F1V)*(-dq*tulin[0]*CT + q[0]*tulin[3]);
1116 sum += Lmunu*Anumu+Lnumu*Amunu;
1127 double tmugev = leptonMom.E() - leptonMom.Mag();
1129 std::ofstream ffstream;
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" 1159 LOG(
"Nieves",
pWARN) <<
"Imaginary part of tensor contraction is nonzero " 1160 <<
"in QEL XSec, real(sum) = " << real(sum)
1161 <<
"imag(sum) = " << imag(sum);
1170 double Rcurr,
int A,
int Z):
1171 ROOT::Math::IBaseFunctionOneDim()
1192 return rhop*rin*rin/
fRcurr;
1198 ROOT::Math::IBaseFunctionOneDim *
1223 double rmaxfrac = 0.25;
1225 bool carbon =
false;
1228 fTensorsOutFile =
"gen.RPA";
1230 fTensorsOutFile =
"gen.noRPA";
1251 rmax = TMath::Max(rp,rn) + 9.25*TMath::Max(ap,an);
1253 rmax = TMath::Sqrt(20.0)*TMath::Max(rp,rn);
1254 double r = rmax * rmaxfrac;
1258 const InitialState & init_state = interaction -> InitState();
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);
1274 double pin = TMath::Sqrt(ein*ein);
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;
1285 TLorentzVector qTildeP4(0, 0, dq, q0);
1286 TLorentzVector inNucleonMomOnShell(0,0,0,0);
1290 TLorentzVector neutrinoMom(0,0,pout+dq,eout+q0);
1291 TLorentzVector leptonMom(0,0,pout,eout);
1295 double Vc = vcr(& target, r);
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";
1307 double plLocal = TMath::Sqrt(ElLocal*ElLocal-ml*ml);
1310 double coulombFactor= plLocal*ElLocal/leptonMom.Vect().Mag()/El;
1311 fCoulombFactor = coulombFactor;
1316 fFormFactors.Calculate(interaction);
1317 LmunuAnumu(neutrinoMom, inNucleonMomOnShell, leptonMom, qTildeP4,
1318 M, is_neutrino, target,
false);
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.
std::complex< double > ruLinRelX(double q0, double qm, double kf, double m) const
double DoEval(double rin) const
ROOT::Math::IntegrationOneDim::Type Integration1DimTypeFromString(string type)
bool IsWeakCC(void) const
bool fRPA
use RPA corrections
bool IsNeutrino(int pdgc)
void Configure(const Registry &config)
double J(double q0, double q3, double Enu, double ml)
THE MAIN GENIE PROJECT NAMESPACE
static constexpr double g
Cross Section Integrator Interface.
double Q2(const Interaction *const i)
int HitNucPdg(void) const
void SetQ2(double Q2, bool selected=false)
Kinematics * KinePtr(void) const
double Density(double r, int A, double ring=0.)
TParticlePDG * RecoilNucleon(void) const
recoil nucleon
QELEvGen_BindingMode_t fIntegralNucleonBindingMode
const FermiMomentumTable * fKFTable
int RecoilNucleonPdg(void) const
recoil nucleon pdg
bool IsQuasiElastic(void) const
double HitNucMass(void) const
bool IsNucleus(void) const
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
static FermiMomentumTablePool * Instance(void)
Generated/set kinematical variables for an event.
Examines whether the generated event should be Pauli blocked. Is a concerete implementation of the Ev...
const XSecIntegratorI * fXSecIntegrator
const TLorentzVector & HadSystP4(void) const
Pure abstract base class. Defines the NuclearModelI interface to be implemented by any physics model ...
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
std::complex< double > deltaLindhard(double q0gev, double dqgev, double rho, double kFgev) const
virtual NuclearModel_t ModelType(const Target &) const =0
const QELFormFactorsModelI * fFormFactorsModel
double fhbarc
hbar*c in GeV*fm
double fCos8c2
cos^2(cabibbo angle)
Summary information for an interaction.
const TLorentzVector & HitNucP4(void) const
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.
const TLorentzVector & FSLeptonP4(void) const
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...
NievesQELvcrIntegrand(double Rcurr, int A, int Z)
static constexpr double m2
const FermiMomentumTable * GetTable(string name)
double EnergyDeltaFunctionSolutionQEL(const Interaction &inter)
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...
bool IsAntiNeutrino(int pdgc)
double Integral(const Interaction *i) const
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
virtual void Configure(const Registry &config)
bool fDoPauliBlocking
Whether to apply Pauli blocking in XSec()
static const double kASmallNum
virtual double Shift(const Target &target) const
unsigned int NDim(void) const
double fXSecScale
external xsec scaling factor
double gamma(double KE, const simb::MCParticle *part)
Misc GENIE control constants.
static const double kLightSpeed
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
double vcr(const Target *target, double r) const
TString fTensorsOutFile
file to print tensors to
virtual const AlgId & Id(void) const
Get algorithm ID.
virtual ~NievesQELCCPXSec()
A registry. Provides the container for algorithm configuration parameters.
const UInt_t kIAssumeFreeNucleon
double HitNucPosition(void) const
QELEvGen_BindingMode_t StringToQELBindingMode(const std::string &mode_str)
virtual double Integrate(const XSecAlgorithmI *model, const Interaction *interaction) const =0
double Jacobian(const Interaction *const i, KinePhaseSpace_t f, KinePhaseSpace_t t)
static constexpr double fermi
const InitialState & InitState(void) const
const ProcessInfo & ProcInfo(void) const
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?
double FindClosestKF(int target_pdgc, int nucleon_pdgc) const
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
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 const double kPlankConstant
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils...
Root of GENIE utility namespaces.
TLorentzVector * GetProbeP4(RefFrame_t rf=kRfHitNucRest) const
const UInt_t kISkipProcessChk
if set, skip process validity checks
std::complex< double > relLindhardIm(double q0gev, double dqgev, double kFngev, double kFpgev, double M, bool isNeutrino, double &t0, double &r00) const
Initial State information.
static const double kPionMass2
const Algorithm * SubAlg(const RgKey ®istry_key) const