40 #include <TLorentzVector.h> 48 #include "Framework/Conventions/GBuild.h" 72 using std::ostringstream;
73 using namespace genie;
80 int pdgc,
const TLorentzVector & x4,
const TLorentzVector & p4,
81 double A,
double Z,
double nRpi,
double nRnuc,
const bool useOset,
const bool altOset,
const bool xsecNNCorr,
string INukeMode)
95 bool is_kaon = pdgc ==
kPdgKP;
98 if(!is_pion && !is_nucleon && !is_kaon && !is_gamma)
return 0.;
109 double ring = (momentum>0) ? 1.240/momentum : 0;
111 if(A<=20) { ring /= 2.; }
118 if(INukeMode==
"hN2018")
120 if (is_pion ) { ring *= nRpi; }
121 else if (is_nucleon ) { ring *= nRnuc; }
122 else if (is_gamma || is_kaon || useOset) { ring = 0.;}
126 if (is_pion || is_kaon ) { ring *= nRpi; }
127 else if (is_nucleon ) { ring *= nRnuc; }
128 else if (is_gamma ) { ring = 0.; }
132 double rnow = x4.Vect().Mag();
138 double ke = (p4.Energy() - p4.M()) /
units::MeV;
145 double ppcnt = (double) Z/ (
double)
A;
148 if (is_pion and (INukeMode ==
"hN2018") and useOset and ke < 350.0)
151 { sigtot = fHadroData2018 -> XSecPipp_Tot() -> Evaluate(ke)*ppcnt;
152 sigtot+= fHadroData2018 -> XSecPipn_Tot() -> Evaluate(ke)*(1-ppcnt);}
154 { sigtot = fHadroData2018 -> XSecPi0p_Tot() -> Evaluate(ke)*ppcnt;
155 sigtot+= fHadroData2018 -> XSecPi0n_Tot() -> Evaluate(ke)*(1-ppcnt);}
157 { sigtot = fHadroData2018 -> XSecPipn_Tot() -> Evaluate(ke)*ppcnt;
158 sigtot+= fHadroData2018 -> XSecPipp_Tot() -> Evaluate(ke)*(1-ppcnt);}
161 sigtot = fHadroData2018 -> XSecPp_Tot() -> Evaluate(ke)*ppcnt;
166 double R0 = 1.25 * TMath::Power(A,1./3.) + 2.0 * 0.65;
167 double Mp = pLib->
Find(2212)->Mass();
168 double M = pLib->
Find(pdgc)->Mass();
171 if (Z*hc/137./x4.Vect().Mag() >
E)
174 double Bc = z*Z*hc/137./R0;
176 double f = TMath::ACos(TMath::Power(x,0.5)) - TMath::Power(x*(1-x),0.5);
177 double B = 0.63*z*Z*TMath::Power((M/Mp)/E,0.5);
178 double Pc = TMath::Exp(-B*f);
181 sigtot+= fHadroData2018 -> XSecPn_Tot() -> Evaluate(ke)*(1-ppcnt);
183 double E0 = TMath::Power(A,0.2)*12.;
184 if (INukeMode==
"hN2018"){
if(ke<E0){sigtot=0.0;}}
189 sigtot = fHadroData2018 -> XSecPn_Tot() -> Evaluate(ke)*ppcnt;
190 sigtot+= fHadroData2018 -> XSecNn_Tot() -> Evaluate(ke)*(1-ppcnt);
191 double E0 = TMath::Power(A,0.2)*12.;
192 if (INukeMode==
"hN2018"){
if(ke<E0){sigtot=0.0;}}
196 { sigtot = fHadroData2018 -> XSecKpN_Tot() -> Evaluate(ke);
200 { sigtot = fHadroData2018 -> XSecGamp_fs() -> Evaluate(ke)*ppcnt;
201 sigtot+= fHadroData2018 -> XSecGamn_fs() -> Evaluate(ke)*(1-ppcnt);}
213 if (xsecNNCorr and is_nucleon)
218 if(sigtot<1
E-6){sigtot=1
E-6;}
221 double lamda = 1. / (rho * sigtot);
224 if( ! TMath::Finite(lamda) ) {
237 int pdgc,
const TLorentzVector & x4,
const TLorentzVector & p4,
double A)
252 if(!is_deltapp)
return 0.;
255 double rnow = x4.Vect().Mag();
260 double ke = (p4.Energy() - p4.M()) /
units::MeV;
261 ke = TMath::Max(1500., ke);
262 ke = TMath::Min( 0., ke);
267 else if (ke<1000) sig=40;
274 double lamda = 1. / (rho * sig);
277 if( ! TMath::Finite(lamda) ) {
285 int pdgc,
const TLorentzVector & x4,
const TLorentzVector & p4,
double A,
double ,
310 double NR = fsi_model.
GetNR();
313 double R0 = fsi_model.
GetR0();
329 double R = NR * R0 * TMath::Power(A, 1./3.);
331 TVector3 dr3 = p4.Vect().Unit();
332 TLorentzVector dr4(dr3,0);
338 <<
"Calculating survival probability for hadron with PDG code = " << pdgc
339 <<
" and momentum = " << p4.P() <<
" GeV";
341 <<
"mfp scale = " << mfp_scale_factor
342 <<
", nRpi = " << nRpi <<
", nRnuc = " << nRnuc <<
", NR = " << NR
343 <<
", R0 = " << R0 <<
" fm";
345 TLorentzVector x4_curr(x4);
348 double rnow = x4_curr.Vect().Mag();
349 if (rnow > (R+step))
break;
351 x4_curr += (step*dr4);
352 rnow = x4_curr.Vect().Mag();
355 remnA, remnZ, nRpi, nRnuc, useOset, altOset, xsecNNCorr, inuke_mode );
356 double mfp_twk = mfp * mfp_scale_factor;
358 double dprob = (mfp_twk>0) ? TMath::Exp(-step/mfp_twk) : 0.;
368 LOG(
"INukeUtils",
pDEBUG) <<
"Psurv = " << prob;
374 const TLorentzVector & x4,
const TLorentzVector & p4,
375 double A,
double NR,
double R0)
381 double R = NR * R0 * TMath::Power(A, 1./3.);
384 TVector3 dr3 = p4.Vect().Unit();
385 TLorentzVector dr4(dr3,0);
387 TLorentzVector x4_curr(x4);
391 double rnow = x4_curr.Vect().Mag();
392 x4_curr += (step*dr4);
394 rnow = x4_curr.Vect().Mag();
401 int pdgc,
const TLorentzVector & x4,
const TLorentzVector & p4,
402 double A,
double Z,
double NR,
double R0)
411 double R = NR * R0 * TMath::Power(A, 1./3.);
414 TVector3 dr3 = p4.Vect().Unit();
415 TLorentzVector dr4(dr3,0);
417 TLorentzVector x4_curr(x4);
422 double rnow = x4_curr.Vect().Mag();
423 x4_curr += (step*dr4);
425 rnow = x4_curr.Vect().Mag();
428 d_mfp += (step/lambda);
446 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 448 <<
"Stepping particle [" << p->
Name() <<
"] by dr = " << step <<
" fm";
452 TVector3 dr = p->
P4()->Vect().Unit();
455 TLorentzVector dx4(dr,dt);
456 TLorentzVector x4new = *(p->
X4()) + dx4;
458 if(nuclear_radius > 0.) {
463 double r = x4new.Vect().Mag();
464 double rmax = nuclear_radius+epsilon;
467 <<
"Particle was stepped too far away (r = " << r <<
" fm)";
469 <<
"Placing it " << epsilon
470 <<
" fm outside the nucleus (r' = " << rmax <<
" fm)";
471 double scale = rmax/
r;
476 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 493 int &RemnA,
int &RemnZ, TLorentzVector &RemnP4,
498 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 500 <<
"PreEquilibrium() is invoked for a : " << p->
Name()
501 <<
" whose kinetic energy is : " << p->
KinE();
508 bool allow_dup =
true;
511 double ppcnt = (double) RemnZ / (
double) RemnA;
520 if(rnd->
RndFsi().Rndm()<ppcnt)
529 ppcnt = (double) RemnZ / (
double) RemnA;
560 if(success)
LOG(
"INukeUtils2018",
pINFO) <<
"Successful phase space decay for pre-equilibrium nucleus FSI event";
564 exception.
SetReason(
"Phase space generation of pre-equilibrium nucleus final state failed, details above");
569 while(p_loc<ev->GetEntries())
584 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 586 <<
"Particle at: " << p_loc;
593 int f_loc = p_loc + 1;
613 for(
unsigned int j=0;j<descendants->
size();j++)
615 loc = (*descendants)[j];
644 int &RemnA,
int &RemnZ, TLorentzVector &RemnP4,
649 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 651 <<
"Equilibrium() is invoked for a : " << p->
Name()
652 <<
" whose kinetic energy is : " << p->
KinE();
659 bool allow_dup =
true;
663 double ppcnt = (double) RemnZ / (
double) RemnA;
673 if(rnd->
RndFsi().Rndm()<ppcnt)
682 ppcnt = (double) RemnZ / (
double) RemnA;
685 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__ 687 <<
"Remnant nucleus (A,Z) = (" << RemnA <<
", " << RemnZ <<
")";
717 if (success)
LOG(
"INukeUtils",
pINFO) <<
"successful equilibrium interaction";
721 exception.
SetReason(
"Phase space generation of compound nucleus final state failed, details above");
730 GHepRecord* ev,
int pcode,
int tcode,
int scode,
int s2code,
double C3CM,
744 double E3L, P3L, E4L, P4L;
745 TVector3 tP1L, tPtot, tbeta, tbetadir, tTrans, tVect;
746 TVector3 tP1zCM, tP2zCM, tP3L, tP4L;
758 M1 = pLib->
Find(pcode)->Mass();
760 M3 = pLib->
Find(scode)->Mass();
761 M4 = pLib->
Find(s2code)->Mass();
764 TLorentzVector t4P1L = *p->
P4();
765 TLorentzVector t4P2L = *t->
P4();
768 double bindE = 0.025;
771 LOG(
"TwoBodyCollision",
pNOTICE) <<
"M1 = " << t4P1L.M() <<
" , M2 = " << t4P2L.M();
772 LOG(
"TwoBodyCollision",
pNOTICE) <<
"E1 = " << t4P1L.E() <<
" , E2 = " << t4P2L.E();
774 if ( (p->
Energy()-p->
Mass()) < bindE ){bindE = 0.;}
777 if((pcode==2112||pcode==2212)&&(t4P1L.E()-M1)<.1) bindE = 0.0;
778 if((pcode==211||pcode==-211||pcode==111)&&(t4P1L.E()-M1)<.08) bindE = 0.0;
779 if((pcode==321)&&(t4P1L.E()-M1)<.1) bindE = 0.0;
782 TLorentzVector t4P3L, t4P4L;
785 P3L = t4P3L.Vect().Mag();
786 P4L = t4P4L.Vect().Mag();
791 <<
"TwoBodyKinematics fails: C3CM, P3 = " << C3CM <<
" " 792 << P3L <<
" " << E3L <<
"\n" <<
" P4 = " 793 << P4L <<
" " << E4L ;
798 P3L = t4P3L.Vect().Mag();
799 P4L = t4P4L.Vect().Mag();
803 <<
"C3CM, P3 = " << C3CM <<
" " 804 << P3L <<
" " << E3L <<
"\n" <<
" P4 = " 805 << P4L <<
" " << E4L ;
808 if(!(TMath::Finite(P3L)) || P3L<.001)
811 <<
"Particle 3 momentum small or non-finite: " << P3L
812 <<
"\n" <<
"--> Assigning .001 as new momentum";
814 E3L = TMath::Sqrt(P3L*P3L + M3*M3);
816 if(!(TMath::Finite(P4L)) || P4L<.001)
819 <<
"Particle 4 momentum small or non-finite: " << P4L
820 <<
"\n" <<
"--> Assigning .001 as new momentum";
822 E4L = TMath::Sqrt(P4L*P4L + M4*M4);
841 <<
"t4P2L= " << t4P2L.E() <<
" " << t4P2L.Z()
842 <<
" RemnP4= " << RemnP4.E() <<
" " << RemnP4.Z() ;
868 LOG(
"INukeUtils",
pINFO) <<
"Successful 2 body collision";
874 double M3,
double M4, TLorentzVector t4P1L, TLorentzVector t4P2L,
875 TLorentzVector &t4P3L, TLorentzVector &t4P4L,
double C3CM, TLorentzVector &RemnP4,
double bindE)
895 double E1L, E2L, P1L, P2L, E3L, P3L;
899 double E1CM, E2CM, E3CM, P3CM;
902 double theta1, theta2, theta5, P1zL, P2zL, P1tL, P2tL;
903 TVector3 tbeta, tbetadir, tTrans, tVect;
904 TVector3 tP1zCM, tP2zCM, vP3L;
905 TLorentzVector t4P1buf, t4P2buf, t4Ptot;
914 if (C3CM < -1. || C3CM > 1.)
return false;
917 S3CM = TMath::Sqrt(1.0 - C3CM*C3CM);
928 t4Ptot = t4P1buf + t4P2buf;
930 LOG(
"INukeUtils",
pINFO) <<
"M1 "<<t4P1L.M()<<
", M2 "<<t4P2L.M();
931 LOG(
"INukeUtils",
pINFO) <<
"E1L "<<E1L<<
", E1CM "<<E1CM;
932 LOG(
"INukeUtils",
pINFO) <<
"bindE = " << bindE;
942 LOG(
"INukeUtils",
pNOTICE) <<
"TwoBodyKinematics Failed: Forbidden by binding energy";
943 LOG(
"INukeUtils",
pNOTICE) <<
"E1L, E2L, M3, M4 : "<< E1L <<
", "<< E2L <<
", "<< M3 <<
", "<< M4;
944 t4P3L.SetPxPyPzE(0,0,0,0);
945 t4P4L.SetPxPyPzE(0,0,0,0);
951 tbeta = t4Ptot.Vect() * (1.0 / (E1L + E2L));
952 tbetadir = tbeta.Unit();
954 gm = 1.0 / TMath::Sqrt(1.0 - beta*beta);
957 theta1 = t4P1buf.Angle(tbeta);
958 theta2 = t4P2buf.Angle(tbeta);
959 P1zL = P1L*TMath::Cos(theta1);
960 P2zL = P2L*TMath::Cos(theta2);
961 P1tL = TMath::Sqrt(P1L*P1L - P1zL*P1zL);
962 P2tL = -TMath::Sqrt(P2L*P2L - P2zL*P2zL);
964 if(TMath::Abs((tVect - tbetadir).Mag())<.01) tVect.SetXYZ(0,1,0);
965 theta5 = tVect.Angle(tbetadir);
966 tTrans = (tVect - TMath::Cos(theta5)*tbetadir).Unit();
969 E1CM = gm*E1L - gm*beta*P1zL;
970 tP1zCM = gm*P1zL*tbetadir - gm*tbeta*E1L;
971 E2CM = gm*E2L - gm*beta*P2zL;
972 tP2zCM = gm*P2zL*tbetadir - gm*tbeta*E2L;
976 LOG(
"INukeUtils",
pINFO) <<
"M1 "<<t4P1L.M()<<
", M2 "<<t4P2L.M();
977 LOG(
"INukeUtils",
pINFO) <<
"E1L "<<E1L<<
", E1CM "<<E1CM;
978 LOG(
"INukeUtils",
pINFO) <<
"P1zL "<<P1zL<<
", P1zCM "<<tP1zCM.Mag()<<
", P1tL "<<P1tL;
979 LOG(
"INukeUtils",
pINFO) <<
"E2L "<<E2L<<
", E2CM "<<E2CM;
980 LOG(
"INukeUtils",
pINFO) <<
"P2zL "<<P2zL<<
", P2zCM "<<tP2zCM.Mag()<<
", P2tL "<<P2tL;
981 LOG(
"INukeUtils",
pINFO) <<
"C3CM "<<C3CM;
984 E3CM = (Et*Et + M3*M3 - M4*M4) / (2.0 * Et);
987 if(E3CM*E3CM - M3*M3<0 || E3CM<0 || Et<0)
989 if (Et<0)
LOG(
"INukeUtils",
pNOTICE) <<
"TwoBodyKinematics Failed: Total energy is negative";
990 if (E3CM<M3)
LOG(
"INukeUtils",
pNOTICE) <<
"TwoBodyKinematics Failed: Scattered Particle 3 CM energy is too small";
991 if (E3CM*E3CM - M3*M3<0)
LOG(
"INukeUtils",
pNOTICE) <<
"TwoBodyKinematics Failed: Scattered Particle 3 CM momentum is nonreal";
992 t4P3L.SetPxPyPzE(0,0,0,0);
993 t4P4L.SetPxPyPzE(0,0,0,0);
996 P3CM = TMath::Sqrt(E3CM*E3CM - M3*M3);
999 P3zL = gm*beta*E3CM + gm*P3CM*C3CM;
1002 P3L = TMath::Sqrt(P3zL*P3zL + P3tL*P3tL);
1003 E3L = TMath::Sqrt(P3L*P3L + M3*M3);
1007 double E4CM = Et-E3CM;
1008 double P4zL = gm*beta*E4CM - gm*P3CM*C3CM;
1009 double P4tL = -1.*P3tL;
1010 double P4L = TMath::Sqrt(P4zL*P4zL + P4tL*P4tL);
1011 double E4L = TMath::Sqrt(P4L*P4L + M4*M4);
1013 LOG(
"INukeUtils",
pINFO) <<
"M3 "<< M3 <<
", M4 "<< M4;
1014 LOG(
"INukeUtils",
pINFO) <<
"E3L "<<E3L<<
", E3CM "<<E3CM;
1015 LOG(
"INukeUtils",
pINFO) <<
"P3zL "<<P3zL<<
", P3tL "<<P3tL;
1016 LOG(
"INukeUtils",
pINFO) <<
"C3L "<<P3zL/P3L;
1017 LOG(
"INukeUtils",
pINFO) <<
"Check:";
1018 LOG(
"INukeUtils",
pINFO) <<
"E4L "<<E4L<<
", E4CM "<<E4CM;
1019 LOG(
"INukeUtils",
pINFO) <<
"P4zL "<<P4zL<<
", P4tL "<<P4tL;
1020 LOG(
"INukeUtils",
pINFO) <<
"P4L "<<P4L;
1021 LOG(
"INukeUtils",
pINFO) <<
"C4L "<<P4zL/P4L;
1023 double echeck = E1L + E2L - (E3L + E4L);
1024 double pzcheck = P1zL+ P2zL - (P3zL + P4zL);
1025 double ptcheck = P1tL+ P2tL - (P3tL + P4tL);
1027 LOG(
"INukeUtils",
pINFO) <<
"Check 4-momentum conservation - Energy "<<echeck<<
", z momentum "<<pzcheck <<
", transverse momentum " << ptcheck ;
1032 if(!(TMath::Finite(P3L)) || P3L<.001)
1035 <<
"Particle 3 momentum small or non-finite: " << P3L
1036 <<
"\n" <<
"--> Assigning .001 as new momentum";
1040 E3L = TMath::Sqrt(P3L*P3L + M3*M3);
1046 vP3L = P3zL*tbetadir + P3tL*tTrans;
1047 vP3L.Rotate(PHI3,tbetadir);
1049 t4P3L.SetVect(vP3L);
1052 t4P4L = t4P1buf + t4P2buf - t4P3L;
1053 t4P4L-= TLorentzVector(0,0,0,bindE);
1060 if(t4P4L.Mag2()<0 || t4P4L.E()<0)
1062 LOG(
"INukeUtils",
pNOTICE)<<
"TwoBodyKinematics Failed: Target mass or energy is negative";
1063 t4P3L.SetPxPyPzE(0,0,0,0);
1064 t4P4L.SetPxPyPzE(0,0,0,0);
1068 if (bindE!=0) RemnP4 += TLorentzVector(0,0,0,bindE);
1074 bool DoFermi,
double FermiFac,
double FermiMomentum,
const NuclearModelI* Nuclmodel)
1086 double M1, M2, M3, M4, M5;
1087 double P1L, P2L, P3L, P4L, P5L;
1088 double E1L, E2L, E3L, E4L, E5L;
1089 double E1CM, E2CM, P3tL;
1090 double PizL, PitL, PiL, EiL;
1091 double EiCM, P4CM2, E4CM2, E5CM2, P3CM, E3CM;
1092 double beta, gm, beta2, gm2;
1093 double P3zL, P4zL, P4tL, P5zL, P5tL;
1094 double Et, M, theta1, theta2;
1096 double theta3, theta4, phi3, phi4, theta5;
1097 TVector3 tP2L, tP1L, tPtot, tbeta, tbetadir, tTrans, tP4L, tP5L;
1098 TVector3 tP1zCM, tP2zCM, tP3L, tPiL, tbeta2, tbetadir2, tVect, tTrans2;
1107 M2 = pLib->
Find(tcode)->Mass();
1118 target.SetHitNucPdg(tcode);
1120 tP2L = FermiFac * Nuclmodel->
Momentum3();
1122 E2L = TMath::Sqrt(tP2L.Mag2() + M2*M2);
1126 tP2L.SetXYZ(0.0, 0.0, 0.0);
1134 P1L = TMath::Sqrt(E1L*E1L - M1*M1);
1135 tP1L = p->
P4()->Vect();
1136 tPtot = tP1L + tP2L;
1138 tbeta = tPtot * (1.0 / (E1L + E2L));
1139 tbetadir = tbeta.Unit();
1141 gm = 1.0 / TMath::Sqrt(1.0 - beta*beta);
1143 theta1 = tP1L.Angle(tbeta);
1144 theta2 = tP2L.Angle(tbeta);
1145 P1zL = P1L*TMath::Cos(theta1);
1146 P2zL = P2L*TMath::Cos(theta2);
1147 tVect.SetXYZ(1,0,0);
1148 if(TMath::Abs((tVect - tbetadir).Mag())<.01) tVect.SetXYZ(0,1,0);
1149 theta5 = tVect.Angle(tbetadir);
1150 tTrans = (tVect - TMath::Cos(theta5)*tbetadir).Unit();
1152 E1CM = gm*E1L - gm*beta*P1zL;
1153 tP1zCM = gm*P1zL*tbetadir - gm*tbeta*E1L;
1154 E2CM = gm*E2L - gm*beta*P2zL;
1155 tP2zCM = gm*P2zL*tbetadir - gm*tbeta*E2L;
1157 M = (rnd->
RndFsi().Rndm()*(Et - M3 - M4 - M5)) + (M4 + M5);
1158 E3CM = (Et*Et + M3*M3 - M*M)/(2*Et);
1160 if(E3CM*E3CM - M3*M3<0)
1163 <<
"PionProduction P3 has non-real momentum - retry kinematics";
1164 LOG(
"INukeUtils",
pNOTICE) <<
"Energy, masses of 3 fs particales:" 1165 << E3CM <<
" " << M3 <<
" " <<
" " << M4 <<
" " << M5;
1167 exception.
SetReason(
"PionProduction particle 3 has non-real momentum");
1171 P3CM = TMath::Sqrt(E3CM*E3CM - M3*M3);
1178 P3zL = gm*beta*E3CM + gm*P3CM*TMath::Cos(theta3);
1179 P3tL = P3CM*TMath::Sin(theta3);
1180 PizL = gm*beta*EiCM - gm*P3CM*TMath::Cos(theta3);
1181 PitL = -P3CM*TMath::Sin(theta3);
1183 P3L = TMath::Sqrt(P3zL*P3zL + P3tL*P3tL);
1184 PiL = TMath::Sqrt(PizL*PizL + PitL*PitL);
1185 E3L = TMath::Sqrt(P3L*P3L + M3*M3);
1186 EiL = TMath::Sqrt(PiL*PiL + M*M);
1189 if(!(TMath::Finite(P3L)) || P3L < .001)
1192 <<
"Particle 3 " << M3 <<
" momentum small or non-finite: " << P3L
1193 <<
"\n" <<
"--> Assigning .001 as new momentum";
1197 E3L = TMath::Sqrt(P3L*P3L + M3*M3);
1200 tP3L = P3zL*tbetadir + P3tL*tTrans;
1201 tPiL = PizL*tbetadir + PitL*tTrans;
1202 tP3L.Rotate(phi3,tbetadir);
1203 tPiL.Rotate(phi3,tbetadir);
1206 tbeta2 = tPiL * (1.0 / EiL);
1207 tbetadir2 = tbeta2.Unit();
1208 beta2 = tbeta2.Mag();
1209 gm2 = 1.0 / TMath::Sqrt(1.0 - beta2*beta2);
1211 E4CM2 = (M*M + M4*M4 - M5*M5) / (2*M);
1213 P4CM2 = TMath::Sqrt(E4CM2*E4CM2 - M4*M4);
1215 tVect.SetXYZ(1,0,0);
1216 if(TMath::Abs((tVect - tbetadir2).Mag())<.01) tVect.SetXYZ(0,1,0);
1217 theta5 = tVect.Angle(tbetadir2);
1218 tTrans2 = (tVect - TMath::Cos(theta5)*tbetadir2).Unit();
1220 P4zL = gm2*beta2*E4CM2 + gm2*P4CM2*TMath::Cos(theta4);
1221 P4tL = P4CM2*TMath::Sin(theta4);
1222 P5zL = gm2*beta2*E5CM2 - gm2*P4CM2*TMath::Cos(theta4);
1225 P4L = TMath::Sqrt(P4zL*P4zL + P4tL*P4tL);
1226 P5L = TMath::Sqrt(P5zL*P5zL + P5tL*P5tL);
1227 E4L = TMath::Sqrt(P4L*P4L + M4*M4);
1228 E5L = TMath::Sqrt(P5L*P5L + M5*M5);
1231 if(!(TMath::Finite(P4L)) || P4L < .001)
1234 <<
"Particle 4 " << M4 <<
" momentum small or non-finite: " << P4L
1235 <<
"\n" <<
"--> Assigning .001 as new momentum";
1239 E4L = TMath::Sqrt(P4L*P4L + M4*M4);
1241 if(!(TMath::Finite(P5L)) || P5L < .001)
1244 <<
"Particle 5 " << M5 <<
" momentum small or non-finite: " << P5L
1245 <<
"\n" <<
"--> Assigning .001 as new momentum";
1249 E5L = TMath::Sqrt(P5L*P5L + M5*M5);
1252 tP4L = P4zL*tbetadir2 + P4tL*tTrans2;
1253 tP5L = P5zL*tbetadir2 + P5tL*tTrans2;
1254 tP4L.Rotate(phi4,tbetadir2);
1255 tP5L.Rotate(phi4,tbetadir2);
1261 <<
"PionProduction fails because of Pauli blocking - retry kinematics";
1263 exception.
SetReason(
"PionProduction final state not determined");
1277 TLorentzVector(tP3L,E3L);
1278 TLorentzVector(tP4L,E4L);
1279 TLorentzVector(tP5L,E5L);
1285 LOG (
"INukeUtils",
pDEBUG) <<
"in Pi Prod, mode = " << mode;
1303 TLorentzVector &RemnP4,
bool DoFermi,
double FermiFac,
double FermiMomentum,
const NuclearModelI* Nuclmodel)
1324 int pcode = p->
Pdg();
1326 int p1code = p->
Pdg();
1328 int p3code = 0, p4code = 0, p5code = 0;
1344 double kine = 1000*p->
KinE();
1352 double xsec2pi0n = TMath::Max(0.,TMath::Exp(8.82 - (573.2/(kine-107.3))));
1355 double xsecpippimn = TMath::Max(0.,TMath::Exp(11.06 - (985.9/(kine-88.2))));
1358 double xsecpimpi0p = TMath::Max(0.,TMath::Exp(9.58 - (1229.4/(kine-60.5))));
1359 double totpimp = xsec2pi0n + xsecpippimn + xsecpimpi0p;
1365 double xsec2pipn = TMath::Max(0.,TMath::Exp(5.64 - (222.6/(kine-150.))));
1368 double xsecpippi0p = TMath::Max(0.,TMath::Exp(7.95 - (852.6/(kine-77.8))));
1369 double totpipp = xsec2pipn + xsecpippi0p;
1371 if (totpimp<=0 && totpipp<=0) {
1372 LOG(
"INukeUtils",
pNOTICE) <<
"InelasticHN called below threshold energy";
1378 double xsecp, xsecn;
1380 case kPdgPi0: xsecp = 0.5 * (totpimp + totpipp); xsecn = xsecp;
break;
1381 case kPdgPiP: xsecp = totpipp; xsecn = totpimp;
break;
1382 case kPdgPiM: xsecp = totpimp; xsecn = totpipp;
break;
1384 LOG(
"INukeUtils",
pWARN) <<
"InelasticHN cannot handle probe: " 1387 exception.
SetReason(
"PionProduction final state not determined");
1396 xsecn *= RemnA-RemnZ;
1400 double rand = rnd->
RndFsi().Rndm() * (xsecp + xsecn);
1402 { rand /= RemnZ; ptarg =
true;}
1404 { rand -= xsecp; rand /= RemnA-RemnZ; ptarg =
false;}
1409 if (((ptarg==
true)&&(p1code==
kPdgPiP))
1410 || ((ptarg==
false)&&(p1code==
kPdgPiM)))
1412 if (rand < xsec2pipn)
1424 else if (((ptarg==
false)&&(p1code==
kPdgPiP))
1425 || ((ptarg==
true)&&(p1code==
kPdgPiM)))
1427 if (rand < xsec2pi0n)
1433 else if (rand < (xsec2pi0n + xsecpippimn))
1448 rand = rnd->
RndFsi().Rndm();
1449 if (rand < 191./270.)
1455 else if (rand < 7./135.)
1472 exception.
SetReason(
"PionProduction final state not determined");
1480 double tote = p->
Energy();
1481 double pMass = pLib->
Find(2212)->Mass();
1482 double nMass = pLib->
Find(2112)->Mass();
1483 double etapp2ppPi0 =
1485 double etapp2pnPip =
1487 pMass+nMass,pLib->
Find(211)->Mass());
1488 double etapn2nnPip =
1490 double etapn2ppPim =
1493 if ((etapp2ppPi0<=0.)&&(etapp2pnPip<=0.)&&(etapn2nnPip<=0.)&&(etapn2ppPim<=0.)) {
1494 LOG(
"INukeUtils",
pNOTICE) <<
"PionProduction() called below threshold energy";
1496 exception.
SetReason(
"PionProduction final state not possible - below threshold");
1502 double xsecppPi0=0,xsecpnPiP=0,xsecnnPiP=0,xsecppPiM=0;
1504 xsecppPi0 = 4511*(1-.91*TMath::Exp(-TMath::Power((etapp2ppPi0-.705),2)));
1505 xsecppPi0 *= (1-TMath::Exp(-TMath::Power((.556*etapp2ppPi0),3.5)));
1506 xsecppPi0 *= (1-TMath::Exp(-56.897/(etapp2ppPi0*etapp2ppPi0)));
1507 xsecppPi0 = TMath::Max(0.,xsecppPi0);}
1510 xsecpnPiP = 18840*(1-.808*TMath::Exp(-TMath::Power((etapp2pnPip-.371),2)));
1511 xsecpnPiP *= (1-TMath::Exp(-TMath::Power((.568*etapp2pnPip),3.2)));
1512 xsecpnPiP *= (1-TMath::Exp(-39.818/(etapp2pnPip*etapp2pnPip)));
1513 xsecpnPiP = TMath::Max(0.,xsecpnPiP);}
1516 xsecnnPiP = 7670*(1-.479*TMath::Exp(-TMath::Power((etapn2nnPip-.947),2)));
1517 xsecnnPiP *= (1-TMath::Exp(-TMath::Power((.35*etapn2nnPip),3.2)));
1518 xsecnnPiP *= (1-TMath::Exp(-71.53/(etapn2nnPip*etapn2nnPip)));
1519 xsecnnPiP = TMath::Max(0.,xsecnnPiP);}
1522 xsecppPiM = 7670*(1-.479*TMath::Exp(-TMath::Power((etapn2ppPim-.947),2)));
1523 xsecppPiM *= (1-TMath::Exp(-TMath::Power((.35*etapn2ppPim),3.2)));
1524 xsecppPiM *= (1-TMath::Exp(-71.53/(etapn2ppPim*etapn2ppPim)));
1525 xsecppPiM = TMath::Max(0.,xsecppPiM);}
1528 double sigma10 = TMath::Max(0.,xsecpnPiP - xsecppPi0);
1529 double sigma01 = TMath::Max(0.,xsecppPiM + xsecnnPiP - xsecppPi0);
1531 double xsecpnPi0 = .5*(sigma10 + sigma01);
1532 xsecpnPi0 = TMath::Max(xsecpnPi0,0.);
1534 LOG(
"INukeUtils",
pDEBUG) <<
'\n' <<
"Cross section values: "<<
'\n' 1535 << xsecppPi0 <<
" PP pi0" <<
'\n' 1536 << xsecpnPiP <<
" PN pi+" <<
'\n' 1537 << xsecnnPiP <<
" NN pi+" <<
'\n' 1538 << xsecpnPi0 <<
" PN pi0";
1540 double xsecp=0,xsecn=0;
1542 case kPdgProton: xsecp=xsecppPi0+xsecpnPiP; xsecn=xsecppPiM+xsecnnPiP+xsecpnPi0;
break;
1543 case kPdgNeutron: xsecp=xsecppPiM+xsecnnPiP+xsecpnPi0; xsecn=xsecppPi0+xsecpnPiP;
break;
1545 LOG(
"INukeUtils",
pWARN) <<
"InelasticHN cannot handle probe: " 1554 xsecn *= RemnA-RemnZ;
1558 double rand = rnd->
RndFsi().Rndm() * (xsecp + xsecn);
1560 { rand /= RemnZ; ptarg =
true;}
1562 { rand -= xsecp; rand /= RemnA-RemnZ; ptarg =
false;}
1596 <<
"Unable to handle probe (=" << p1code <<
") in InelasticHN()";
1603 LOG(
"INukeUtils",
pNOTICE) <<
"PionProduction() failed : no nucleons to produce pions";
1611 LOG(
"INukeUtils",
pNOTICE) <<
"PionProduction() failed : too few protons in nucleus";
1613 exception.
SetReason(
"PionProduction fails - too few protons available");
1639 LOG(
"INukeUtils",
pDEBUG) <<
"Remnant (A,Z) = (" <<RemnA<<
','<<RemnZ<<
')';
1641 RemnP4 -= *s1->
P4() + *s2->
P4() + *s3->
P4() - *p->
P4();
1646 exception.
SetReason(
"PionProduction final state not determined");
1653 double Mtwopart,
double Mpi)
1663 double Scm = Minc*Minc + Mtarg*Mtarg + 2*Mtarg*nrg;
1665 eta= TMath::Power(Scm,2) + TMath::Power(Mtwopart,4) + TMath::Power(Mpi,4);
1666 eta-= 2*TMath::Power(Mtwopart*Mpi,2);
1667 eta-= 2*Scm*TMath::Power(Mtwopart,2);
1668 eta-= 2*Scm*TMath::Power(Mpi,2);
1669 eta = TMath::Power(eta/Scm,1./2.);
1672 eta=TMath::Max(eta,0.);
1680 TLorentzVector &RemnP4,
double NucRmvE,
EINukeMode mode)
1686 LOG(
"INukeUtils",
pINFO) <<
"*** Performing a Phase Space Decay";
1687 assert(pdgv.
size() > 1);
1689 LOG(
"INukeUtils",
pINFO) <<
"probe mass: M = " << p->
Mass();
1693 ostringstream state_sstream;
1694 state_sstream <<
"( ";
1697 double * mass =
new double[pdgv.
size()];
1698 double mass_sum = 0;
1699 for(pdg_iter = pdgv.begin(); pdg_iter != pdgv.end(); ++pdg_iter) {
1700 int pdgc = *pdg_iter;
1705 state_sstream << nm <<
" ";
1707 state_sstream <<
")";
1709 TLorentzVector * pd = p->
GetP4();
1716 double availE = pd->Energy() + mass_sum;
1717 if(is_nuc||is_kaon) availE -= p->
Mass();
1721 <<
"size, mass_sum, availE, pd mass, energy = " << pdgv.
size() <<
" " 1722 << mass_sum <<
" " << availE <<
" " << p->
Mass() <<
" " << p->
Energy() ;
1725 double dE = mass_sum;
1726 if(is_nuc||is_kaon) dE -= p->
Mass();
1727 TLorentzVector premnsub(0,0,0,dE);
1731 <<
"Final state = " << state_sstream.str() <<
" has N = " << pdgv.
size()
1732 <<
" particles / total mass = " << mass_sum;
1737 TGenPhaseSpace GenPhaseSpace;
1738 bool permitted = GenPhaseSpace.SetDecay(*pd, pdgv.
size(), mass);
1741 <<
" *** Phase space decay is not permitted \n" 1742 <<
" Total particle mass = " << mass_sum <<
"\n" 1759 for(
int k=0;
k<200;
k++) {
1760 double w = GenPhaseSpace.Generate();
1761 wmax = TMath::Max(wmax,w);
1766 <<
"Max phase space gen. weight @ current hadronic interaction: " << wmax;
1773 bool accept_decay=
false;
1774 unsigned int itry=0;
1776 while(!accept_decay)
1783 <<
"Couldn't generate an unweighted phase space decay after " 1784 << itry <<
" attempts";
1790 double w = GenPhaseSpace.Generate();
1791 double gw = wmax * rnd->
RndFsi().Rndm();
1795 <<
"Decay weight = " << w <<
" > max decay weight = " << wmax;
1798 LOG(
"INukeUtils",
pNOTICE) <<
"Decay weight = " << w <<
" / R = " << gw;
1799 accept_decay = (gw<=
w);
1807 LOG(
"INukeUtils",
pNOTICE) <<
"mother index = " << mom;
1811 TLorentzVector * v4 = p->
GetX4();
1813 double checkpx = p->
Px();
1814 double checkpy = p->
Py();
1815 double checkpz = p->
Pz();
1816 double checkE = p->
E();
1820 for(pdg_iter = pdgv.begin(); pdg_iter != pdgv.end(); ++pdg_iter) {
1823 int pdgc = *pdg_iter;
1827 TLorentzVector * p4fin = GenPhaseSpace.GetDecay(i++);
1834 double En = p4fin->Energy();
1841 double dE_leftover = TMath::Min(NucRmvE, KE);
1844 double pmag_old = p4fin->P();
1845 double pmag_new = TMath::Sqrt(TMath::Max(0.,En*En-M*M));
1846 double scale = pmag_new / pmag_old;
1847 double pxn = scale * p4fin->Px();
1848 double pyn = scale * p4fin->Py();
1849 double pzn = scale * p4fin->Pz();
1851 TLorentzVector p4n(pxn,pyn,pzn,En);
1862 if (p4n.Vect().Mag()>=0.001)
1864 GHepParticle new_particle(pdgc, ist_pi, mom,-1,-1,-1, p4n, *v4);
1872 LOG(
"INukeUtils",
pINFO)<<
"Momentum too small; assigning 0.001 as new momentum";
1874 double phi = 2*
kPi*rnd->
RndFsi().Rndm();
1875 double omega = 2*rnd->
RndFsi().Rndm();
1878 double E4n = TMath::Sqrt(0.001*0.001+M*M);
1879 p4n.SetPxPyPzE(0.001,0,0,E4n);
1880 p4n.Rotate(TMath::ACos(1-omega),TVector3(0,0,1));
1881 p4n.Rotate(phi,TVector3(1,0,0));
1883 RemnP4 -= (p4n - TLorentzVector(0,0,0,M));
1885 GHepParticle new_particle(pdgc, ist, mom,-1,-1,-1, p4n, *v4);
1891 GHepParticle new_particle(pdgc, ist, mom,-1,-1,-1, p4n, *v4);
1897 double dpx = (1-scale)*p4fin->Px();
1898 double dpy = (1-scale)*p4fin->Py();
1899 double dpz = (1-scale)*p4fin->Pz();
1900 TLorentzVector premnadd(dpx,dpy,dpz,dE_leftover);
1904 LOG(
"INukeUtils",
pNOTICE) <<
"check conservation: Px = " << checkpx <<
" Py = " << checkpy
1905 <<
" Pz = " << checkpz <<
" E = " << checkE;
1916 const double &pionKineticEnergy,
1917 const double &density,
1919 const double &protonFraction,
1920 const bool &isTableChosen
1926 if (iNukeOset == NULL)
1931 static const std::string dataDir = (gSystem->Getenv(
"GINUKEHADRONDATA")) ?
1932 string(gSystem->Getenv(
"GINUKEHADRONDATA")) :
1933 string(gSystem->Getenv(
"GENIE")) +
1934 string(
"/data/evgen/intranuke/");
1936 static const std::string dataFile = dataDir +
"tot_xsec/" 1937 "intranuke-xsection-pi+N-Oset.dat";
1946 iNukeOset->
setupOset (density, pionKineticEnergy, pionPDG, protonFraction);
double MeanFreePath_Delta(int pdgc, const TLorentzVector &x4, const TLorentzVector &p4, double A)
Mean free path (Delta++ test)
const int kPdgP33m1232_DeltaPP
An exception thrown by SimulateHadronState for kinematics problems. TwoBodyCollision/Kinematics used ...
void SetFirstMother(int m)
virtual GHepParticle * Particle(int position) const
void PreEquilibrium(GHepRecord *ev, GHepParticle *p, int &RemnA, int &RemnZ, TLorentzVector &RemnP4, bool DoFermi, double FermiFac, const NuclearModelI *Nuclmodel, double NucRmvE, EINukeMode mode=kIMdHN)
virtual int ParticlePosition(int pdg, GHepStatus_t i, int start=0) const
TLorentzVector * GetX4(void) const
TRandom3 & RndFsi(void) const
rnd number generator used by intranuclear cascade monte carlos
double beta(double KE, const simb::MCParticle *part)
static INukeNucleonCorr * getInstance()
get single instance of INukeNucleonCorr; create if necessary
THE MAIN GENIE PROJECT NAMESPACE
double E(void) const
Get energy.
static RandomGen * Instance()
Access instance.
const TLorentzVector * P4(void) const
bool PhaseSpaceDecay(GHepRecord *ev, GHepParticle *p, const PDGCodeList &pdgv, TLorentzVector &RemnP4, double NucRmvE, EINukeMode mode=kIMdHA)
general phase space decay method
double Density(double r, int A, double ring=0.)
virtual string GetINukeMode() const
Table-based implementation of Oset model.
void StepParticle(GHepParticle *p, double step, double nuclear_radius=-1.)
Step particle.
bool ThreeBodyKinematics(GHepRecord *ev, GHepParticle *p, int tcode, GHepParticle *s1, GHepParticle *s2, GHepParticle *s3, bool DoFermi=false, double FermiFac=0, double FermiMomentum=0, const NuclearModelI *Nuclmodel=(const NuclearModelI *) 0)
bool PionProduction(GHepRecord *ev, GHepParticle *p, GHepParticle *s1, GHepParticle *s2, GHepParticle *s3, int &RemnA, int &RemnZ, TLorentzVector &RemnP4, bool DoFermi, double FermiFac, double FermiMomentum, const NuclearModelI *Nuclmodel)
string P4AsString(const TLorentzVector *p)
static const unsigned int kMaxUnweightDecayIterations
Pure abstract base class. Defines the NuclearModelI interface to be implemented by any physics model ...
double Mass(Resonance_t res)
resonance mass (GeV)
void SetMomentum(const TLorentzVector &p4)
double GetHadStep() const
double Mass(void) const
Mass that corresponds to the PDG code.
static constexpr double MeV
A singleton holding random number generator classes. All random number generation in GENIE should tak...
double Pz(void) const
Get Pz.
const TVector3 & Momentum3(void) const
double Energy(void) const
Get energy.
double GetDelRPion() const
double Px(void) const
Get Px.
bool CompareMomentum(const GHepParticle *p) const
int LastMother(void) const
double Dist2Exit(const TLorentzVector &x4, const TLorentzVector &p4, double A, double NR=3, double R0=1.4)
Distance to exit.
int FirstMother(void) const
string Name(void) const
Name that corresponds to the PDG code.
const int kPdgCompNuclCluster
void SetPosition(const TLorentzVector &v4)
double getTotalCrossSection() const
return total = (qel+cex+abs) cross section
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
static constexpr double mb
void SetReason(string reason)
TLorentzVector * GetP4(void) const
bool ComparePdgCodes(const GHepParticle *p) const
static constexpr double fm2
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
bool TwoBodyKinematics(double M3, double M4, TLorentzVector tP1L, TLorentzVector tP2L, TLorentzVector &tP3L, TLorentzVector &tP4L, double C3CM, TLorentzVector &RemnP4, double bindE=0)
bool TwoBodyKinematics(double M3, double M4, TLorentzVector tP1L, TLorentzVector tP2L, TLorentzVector &tP3L, TLorentzVector &tP4L, double C3CM, TLorentzVector &RemnP4, double bindE=0)
bool TwoBodyCollision(GHepRecord *ev, int pcode, int tcode, int scode, int s2code, double C3CM, GHepParticle *p, GHepParticle *t, int &RemnA, int &RemnZ, TLorentzVector &RemnP4, EINukeMode mode=kIMdHA)
Intranuke utility functions.
virtual GHepParticle * TargetNucleus(void) const
virtual vector< int > * GetStableDescendants(int position) const
double CalculateEta(double Minc, double ke, double Mtarg, double Mtwopart, double Mpi)
double MeanFreePath(int pdgc, const TLorentzVector &x4, const TLorentzVector &p4, double A, double Z, double nRpi=0.5, double nRnuc=1.0, const bool useOset=false, const bool altOset=false, const bool xsecNNCorr=false, string INukeMode="XX2018")
Mean free path (pions, nucleons)
void Equilibrium(GHepRecord *ev, GHepParticle *p, int &RemnA, int &RemnZ, TLorentzVector &RemnP4, bool DoFermi, double FermiFac, const NuclearModelI *Nuclmodel, double NucRmvE, EINukeMode mode=kIMdHN)
Misc GENIE control constants.
virtual void setupOset(const double &density, const double &pionTk, const int &pionPDG, const double &protonFraction)=0
use to set up Oset class (assign pion Tk, nuclear density etc)
static double fMinKinEnergy
bool CompareStatusCodes(const GHepParticle *p) const
void SetRemovalEnergy(double Erm)
double KinE(bool mass_from_pdg=false) const
Get kinetic energy.
static PDGLibrary * Instance(void)
static INukeHadroData2018 * Instance(void)
void SetLastMother(int m)
double ProbSurvival(int pdgc, const TLorentzVector &x4, const TLorentzVector &p4, double A, double Z, double mfp_scale_factor, const Intranuke2018 &fsi_model)
Hadron survival probability.
Singleton class to load & serve a TDatabasePDG.
const TLorentzVector * X4(void) const
void SetStatus(GHepStatus_t s)
bool IsNeutronOrProton(int pdgc)
virtual bool GenerateNucleon(const Target &) const =0
virtual void AddParticle(const GHepParticle &p)
static double fMaxKinEnergyHN
TParticlePDG * Find(int pdgc, bool must_exist=true)
string X4AsString(const TLorentzVector *x)
def momentum(x1, x2, x3, scale=1.)
double sigmaTotalOset(const double &pionKineticEnergy, const double &density, const int &pionPDG, const double &protonFraction, const bool &isTableChosen=true)
GENIE's GHEP MC event record.
double GetDelRNucleon() const
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils...
STDHEP-like event record entry that can fit a particle or a nucleus.
double Dist2ExitMFP(int pdgc, const TLorentzVector &x4, const TLorentzVector &p4, double A, double Z, double NR=3, double R0=1.4)
Distance to exit.
bool GetXsecNNCorr() const
string Vec3AsString(const TVector3 *vec)
Root of GENIE utility namespaces.
cet::coded_exception< error, detail::translate > exception
void push_back(int pdg_code)
enum genie::EGHepStatus GHepStatus_t
double Py(void) const
Get Py.