17 #include "TDatabasePDG.h" 18 #include "TLorentzVector.h" 39 TLorentzVector totMom = *k4 + *p4;
41 TVector3
beta = totMom.BoostVector();
62 TLorentzVector total_p4 = (*probe) + hit_nucleon;
63 TVector3 beta_COM_to_lab = total_p4.BoostVector();
64 TVector3 beta_lab_to_COM = -beta_COM_to_lab;
67 double gamma2 =
std::pow(total_p4.Gamma(), 2);
70 TLorentzVector leptonCOM = TLorentzVector( lepton );
71 leptonCOM.Boost( beta_lab_to_COM );
74 double theta0 = leptonCOM.Angle( beta_COM_to_lab );
75 double cos_theta0_squared =
std::pow(std::cos( theta0 ), 2);
79 TVector3 lepton_vel = ( lepton.Vect() * (1. / lepton.E()) );
80 TVector3 outNucleon_vel = ( outNucleon.Vect() * (1. / outNucleon.E()) );
81 double vel_diff = ( lepton_vel - outNucleon_vel ).Mag();
85 double factor = std::sqrt( 1. + (1. - cos_theta0_squared)*(gamma2 - 1.) )
86 *
std::pow(leptonCOM.P(), 2) / vel_diff;
95 double cos_theta_0,
double phi_0,
double& Eb,
101 if ( bind_nucleon ) {
103 hitNucleonBindingMode);
110 TDatabasePDG *tb = TDatabasePDG::Instance();
118 if ( std::sqrt(s) < lepMass + mNf )
return 0.;
120 double outLeptonEnergy = ( s - mNf*mNf + lepMass*lepMass ) / (2 * std::sqrt(s));
122 if (outLeptonEnergy*outLeptonEnergy - lepMass*lepMass < 0.)
return 0.;
123 double outMomentum = TMath::Sqrt(outLeptonEnergy*outLeptonEnergy - lepMass*lepMass);
135 TVector3 lepton3Mom(0., 0., outMomentum);
136 lepton3Mom.SetTheta( TMath::ACos(cos_theta_0) );
137 lepton3Mom.SetPhi( phi_0 );
141 TVector3 zvec(0., 0., 1.);
142 TVector3 rot = ( zvec.Cross(beta) ).Unit();
143 double angle = beta.Angle( zvec );
147 if ( beta.Perp() == 0. && beta.Z() < 0. ) {
148 rot = TVector3(0., 1., 0.);
154 lepton3Mom.Rotate(angle, rot);
158 TLorentzVector lepton(lepton3Mom, outLeptonEnergy);
162 TLorentzVector outNucleon(-1*lepton.Px(),-1*lepton.Py(),-1*lepton.Pz(),
163 TMath::Sqrt(outMomentum*outMomentum + mNf*mNf));
167 outNucleon.Boost(beta);
175 TLorentzVector qP4 = *nuP4 - lepton;
177 double Q2 = -1 * qP4.Mag2();
186 if (Q2 < Q2lim.min || Q2 > Q2lim.
max)
return 0.;
199 if ( binding_mode ==
"UseGroundStateRemnant" ) {
202 else if ( binding_mode ==
"UseNuclearModel" ) {
205 else if ( binding_mode ==
"OnShell" ) {
208 else if ( binding_mode ==
"ValenciaStyleQValue" ) {
212 LOG(
"QELEvent",
pFATAL) <<
"Unrecognized setting \"" << binding_mode
213 <<
"\" requested in genie::utils::StringToQELBindingMode()";
230 double gamma = 1. / std::sqrt(1. - beta.Mag2());
235 double lepton_E_COM = (sqrt_s*sqrt_s + ml*ml - mNf*mNf) / (2.*sqrt_s);
242 double lepton_p_COM = std::sqrt(
std::max(lepton_E_COM*lepton_E_COM
248 double ENi = p4Ni.E();
252 double ENi_on_shell = std::sqrt( mNi*mNi + p4Ni.Vect().Mag2() );
254 double epsilon_B = ENi_on_shell - ENi;
256 double cos_theta0_max = ( probe_E_lab - gamma*lepton_E_COM - epsilon_B )
257 / ( gamma * lepton_p_COM * beta.Mag() );
258 return cos_theta0_max;
272 TDatabasePDG* tb = TDatabasePDG::Instance();
287 double Mi = tgt->
Mass();
309 int Af = tgt->
A() - 1;
346 LOG(
"QELEvent",
pFATAL ) <<
"Unhandled probe PDG code " << probe_pdg
347 <<
" encountered for a CC interaction in" 348 <<
" genie::utils::BindHitNucleon()";
358 Qvalue = mf_keep_nucleon - Mi;
373 double EFermi_Ni = std::sqrt(
std::max(0., mN*mN + kF_Ni*kF_Ni) );
379 double EFermi_Nf = std::sqrt(
std::max(0., mN*mN + kF_Nf*kF_Nf) );
382 Q_LFG = EFermi_Nf - EFermi_Ni;
387 LOG(
"QELEvent",
pFATAL ) <<
"Unhandled process type \"" << info
388 <<
"\" encountered in genie::utils::BindHitNucleon()";
394 double ENi_OnShell = std::sqrt(
std::max(0., mNi*mNi + p3Ni.Mag2()) );
397 double Ef = Mi - ENi_OnShell + Qvalue - Q_LFG;
401 Mf = std::sqrt(
std::max(0., Ef*Ef - p3Ni.Mag2()) );
406 LOG(
"QELEvent",
pDEBUG ) <<
"Qvalue = " << Qvalue
407 <<
", Q_LFG = " << Q_LFG <<
" at radius = " << tgt->
HitNucPosition();
412 ENi = Mi - std::sqrt( Mf*Mf + p3Ni.Mag2() );
418 ENi = std::sqrt( p3Ni.Mag2() +
std::pow(mNi, 2) );
428 LOG(
"QELEvent",
pDEBUG ) <<
"Eb = " << Eb <<
", pNi = " << p3Ni.Mag()
429 <<
", ENi = " << ENi;
433 p4Ni->SetVect( p3Ni );
Cross Section Calculation Interface.
const KPhaseSpace & PhaseSpace(void) const
virtual double LocalFermiMomentum(const Target &, int nucleon_pdg, double radius) const
bool IsWeakCC(void) const
double beta(double KE, const simb::MCParticle *part)
bool IsNeutrino(int pdgc)
double RemovalEnergy(void) const
static const double kNucleonMass
double Q2(const Interaction *const i)
int HitNucPdg(void) const
void SetQ2(double Q2, bool selected=false)
Kinematics * KinePtr(void) const
TParticlePDG * RecoilNucleon(void) const
recoil nucleon
virtual double XSec(const Interaction *i, KinePhaseSpace_t k=kPSfE) const =0
Compute the cross section for the input interaction.
int RecoilNucleonPdg(void) const
recoil nucleon pdg
A simple [min,max] interval for doubles.
double HitNucMass(void) const
bool IsNucleus(void) const
const TLorentzVector & HadSystP4(void) const
Pure abstract base class. Defines the NuclearModelI interface to be implemented by any physics model ...
double Mass(Resonance_t res)
resonance mass (GeV)
void BindHitNucleon(Interaction &interaction, const NuclearModelI &nucl_model, double &Eb, DMELEvGen_BindingMode_t hitNucleonBindingMode)
double ComputeFullQELPXSec(Interaction *interaction, const NuclearModelI *nucl_model, const XSecAlgorithmI *xsec_model, double cos_theta_0, double phi_0, double &Eb, QELEvGen_BindingMode_t hitNucleonBindingMode, double min_angle_EM=0., bool bind_nucleon=true)
Range1D_t Q2Lim(void) const
Q2 limits.
enum genie::EQELEvGenBindingMode QELEvGen_BindingMode_t
const TVector3 & Momentum3(void) const
double Qvalue(int targetpdg, int nupdg)
Summary information for an interaction.
const TLorentzVector & HitNucP4(void) const
bool IsWeakNC(void) const
const TLorentzVector & FSLeptonP4(void) const
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
void SetFSLeptonP4(const TLorentzVector &p4)
double EnergyDeltaFunctionSolutionQEL(const Interaction &inter)
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
bool IsAntiNeutrino(int pdgc)
const Kinematics & Kine(void) const
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
static const double kASmallNum
static int max(int a, int b)
double gamma(double KE, const simb::MCParticle *part)
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
TLorentzVector * HitNucP4Ptr(void) const
static PDGLibrary * Instance(void)
double CosTheta0Max(const genie::Interaction &interaction)
const UInt_t kIAssumeFreeNucleon
double HitNucPosition(void) const
QELEvGen_BindingMode_t StringToQELBindingMode(const std::string &mode_str)
double CMEnergy() const
centre-of-mass energy (sqrt s)
Target * TgtPtr(void) const
int IonPdgCode(int A, int Z)
void SetHadSystP4(const TLorentzVector &p4)
InitialState * InitStatePtr(void) const
const InitialState & InitState(void) const
const ProcessInfo & ProcInfo(void) const
TParticlePDG * Find(int pdgc, bool must_exist=true)
const Target & Tgt(void) const
double ProbeE(RefFrame_t rf) const
TLorentzVector * GetProbeP4(RefFrame_t rf=kRfHitNucRest) const
Initial State information.