QELUtils.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  Steven Gardiner <gardiner \at fnal.gov>
7  Fermi National Accelerator Laboratory
8 */
9 //____________________________________________________________________________
10 
11 // standard library includes
12 #include <cmath>
13 #include <limits>
14 #include <sstream>
15 
16 // ROOT includes
17 #include "TDatabasePDG.h"
18 #include "TLorentzVector.h"
19 
20 // GENIE includes
30 
31 
32 // Helper functions local to this source file
33 namespace {
34 
35  TVector3 COMframe2Lab(const genie::InitialState& initialState)
36  {
37  TLorentzVector* k4 = initialState.GetProbeP4( genie::kRfLab );
38  TLorentzVector* p4 = initialState.TgtPtr()->HitNucP4Ptr();
39  TLorentzVector totMom = *k4 + *p4;
40 
41  TVector3 beta = totMom.BoostVector();
42 
43  delete k4;
44 
45  return beta;
46  }
47 
48 }
49 
51  const genie::Interaction& inter)
52 {
53  // Get the final-state lepton and struck nucleon 4-momenta in the lab frame
54  const TLorentzVector& lepton = inter.Kine().FSLeptonP4();
55  const TLorentzVector& outNucleon = inter.Kine().HadSystP4();
56 
57  // Get beta for a Lorentz boost from the lab frame to the COM frame and
58  // vice-versa
59  TLorentzVector* probe = inter.InitStatePtr()->GetProbeP4( kRfLab );
60  const TLorentzVector& hit_nucleon = inter.InitStatePtr()->TgtPtr()
61  ->HitNucP4();
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;
65 
66  // Square of the Lorentz gamma factor for the boosts
67  double gamma2 = std::pow(total_p4.Gamma(), 2);
68 
69  // Get the final-state lepton 4-momentum in the COM frame
70  TLorentzVector leptonCOM = TLorentzVector( lepton );
71  leptonCOM.Boost( beta_lab_to_COM );
72 
73  // Angle between COM outgoing lepton and lab-frame velocity of COM frame
74  double theta0 = leptonCOM.Angle( beta_COM_to_lab );
75  double cos_theta0_squared = std::pow(std::cos( theta0 ), 2);
76  //double sin_theta0_squared = 1. - cos_theta0_squared; // sin^2(x) + cos^2(x) = 1
77 
78  // Difference in lab frame velocities between lepton and outgoing nucleon
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();
82 
83  // Compute the factor in the kPSQELEvGen phase space that arises from
84  // eliminating the energy-conserving delta function
85  double factor = std::sqrt( 1. + (1. - cos_theta0_squared)*(gamma2 - 1.) )
86  * std::pow(leptonCOM.P(), 2) / vel_diff;
87 
88  delete probe;
89 
90  return factor;
91 }
92 
94  const genie::NuclearModelI* nucl_model, const genie::XSecAlgorithmI* xsec_model,
95  double cos_theta_0, double phi_0, double& Eb,
96  genie::QELEvGen_BindingMode_t hitNucleonBindingMode, double min_angle_EM,
97  bool bind_nucleon)
98 {
99  // If requested, set the initial hit nucleon 4-momentum to be off-shell
100  // according to the binding mode specified in the function call
101  if ( bind_nucleon ) {
102  genie::utils::BindHitNucleon(*interaction, *nucl_model, Eb,
103  hitNucleonBindingMode);
104  }
105 
106  // Mass of the outgoing lepton
107  double lepMass = interaction->FSPrimLepton()->Mass();
108 
109  // Look up the (on-shell) mass of the final nucleon
110  TDatabasePDG *tb = TDatabasePDG::Instance();
111  double mNf = tb->GetParticle( interaction->RecoilNucleonPdg() )->Mass();
112 
113  // Mandelstam s for the probe/hit nucleon system
114  double s = std::pow( interaction->InitState().CMEnergy(), 2 );
115 
116  // Return a differential cross section of zero if we're below threshold (and
117  // therefore need to sample a new event)
118  if ( std::sqrt(s) < lepMass + mNf ) return 0.;
119 
120  double outLeptonEnergy = ( s - mNf*mNf + lepMass*lepMass ) / (2 * std::sqrt(s));
121 
122  if (outLeptonEnergy*outLeptonEnergy - lepMass*lepMass < 0.) return 0.;
123  double outMomentum = TMath::Sqrt(outLeptonEnergy*outLeptonEnergy - lepMass*lepMass);
124 
125  // Compute the boost vector for moving from the COM frame to the
126  // lab frame, i.e., the velocity of the COM frame as measured
127  // in the lab frame.
128  TVector3 beta = COMframe2Lab( interaction->InitState() );
129 
130  // FullDifferentialXSec depends on theta_0 and phi_0, the lepton COM
131  // frame angles with respect to the direction of the COM frame velocity
132  // as measured in the lab frame. To generate the correct dependence
133  // here, first set the lepton COM frame angles with respect to +z
134  // (via TVector3::SetTheta() and TVector3::SetPhi()).
135  TVector3 lepton3Mom(0., 0., outMomentum);
136  lepton3Mom.SetTheta( TMath::ACos(cos_theta_0) );
137  lepton3Mom.SetPhi( phi_0 );
138 
139  // Then rotate the lepton 3-momentum so that the old +z direction now
140  // points along the COM frame velocity (beta)
141  TVector3 zvec(0., 0., 1.);
142  TVector3 rot = ( zvec.Cross(beta) ).Unit();
143  double angle = beta.Angle( zvec );
144 
145  // Handle the edge case where beta is along -z, so the
146  // cross product above vanishes
147  if ( beta.Perp() == 0. && beta.Z() < 0. ) {
148  rot = TVector3(0., 1., 0.);
149  angle = genie::constants::kPi;
150  }
151 
152  // Rotate if the rotation vector is not 0
153  if ( rot.Mag() >= genie::controls::kASmallNum ) {
154  lepton3Mom.Rotate(angle, rot);
155  }
156 
157  // Construct the lepton 4-momentum in the COM frame
158  TLorentzVector lepton(lepton3Mom, outLeptonEnergy);
159 
160  // The final state nucleon will have an equal and opposite 3-momentum
161  // in the COM frame and will be on the mass shell
162  TLorentzVector outNucleon(-1*lepton.Px(),-1*lepton.Py(),-1*lepton.Pz(),
163  TMath::Sqrt(outMomentum*outMomentum + mNf*mNf));
164 
165  // Boost the 4-momenta for both particles into the lab frame
166  lepton.Boost(beta);
167  outNucleon.Boost(beta);
168 
169  // Check if event is at a low angle - if so return 0 and stop wasting time
170  if (180 * lepton.Theta() / genie::constants::kPi < min_angle_EM && interaction->ProcInfo().IsEM()) {
171  return 0;
172  }
173 
174  TLorentzVector * nuP4 = interaction->InitState().GetProbeP4( genie::kRfLab );
175  TLorentzVector qP4 = *nuP4 - lepton;
176  delete nuP4;
177  double Q2 = -1 * qP4.Mag2();
178 
179  interaction->KinePtr()->SetFSLeptonP4( lepton );
180  interaction->KinePtr()->SetHadSystP4( outNucleon );
181  interaction->KinePtr()->SetQ2( Q2 );
182 
183  // Check the Q2 range. If we're outside of it, don't bother
184  // with the rest of the calculation.
185  Range1D_t Q2lim = interaction->PhaseSpace().Q2Lim();
186  if (Q2 < Q2lim.min || Q2 > Q2lim.max) return 0.;
187 
188  // Compute the QE cross section for the current kinematics
189  double xsec = xsec_model->XSec(interaction, genie::kPSQELEvGen);
190 
191  return xsec;
192 }
193 
195  const std::string& binding_mode)
196 {
197  // Translate the string setting the binding mode to the appropriate
198  // enum value, or complain if one couldn't be found
199  if ( binding_mode == "UseGroundStateRemnant" ) {
200  return kUseGroundStateRemnant;
201  }
202  else if ( binding_mode == "UseNuclearModel" ) {
203  return kUseNuclearModel;
204  }
205  else if ( binding_mode == "OnShell" ) {
206  return kOnShell;
207  }
208  else if ( binding_mode == "ValenciaStyleQValue" ) {
209  return kValenciaStyleQValue;
210  }
211  else {
212  LOG("QELEvent", pFATAL) << "Unrecognized setting \"" << binding_mode
213  << "\" requested in genie::utils::StringToQELBindingMode()";
214  gAbortingInErr = true;
215  std::exit(1);
216  }
217 }
218 
220 
221  // q0 > 0 only needs to be enforced (indirectly via a calculation of
222  // CosTheta0Max) for bound nucleons. The Q2 limits should take care of valid
223  // kinematics for free nucleons.
224  if ( !interaction.InitState().Tgt().IsNucleus()
225  || interaction.TestBit(kIAssumeFreeNucleon) ) return 1.;
226 
227  double probe_E_lab = interaction.InitState().ProbeE( genie::kRfLab );
228 
229  TVector3 beta = COMframe2Lab( interaction.InitState() );
230  double gamma = 1. / std::sqrt(1. - beta.Mag2());
231 
232  double sqrt_s = interaction.InitState().CMEnergy();
233  double mNf = interaction.RecoilNucleon()->Mass();
234  double ml = interaction.FSPrimLepton()->Mass();
235  double lepton_E_COM = (sqrt_s*sqrt_s + ml*ml - mNf*mNf) / (2.*sqrt_s);
236 
237  // If there isn't enough available energy to create an on-shell
238  // final lepton, then don't bother with the rest of the calculation
239  // NOTE: C++11 would allow use to use lowest() here instead
240  if ( lepton_E_COM <= ml ) return -std::numeric_limits<double>::max();
241 
242  double lepton_p_COM = std::sqrt( std::max(lepton_E_COM*lepton_E_COM
243  - ml*ml, 0.) );
244 
245  // Possibly off-shell initial struck nucleon total energy
246  // (BindHitNucleon() should have been called previously if needed)
247  const TLorentzVector& p4Ni = interaction.InitState().Tgt().HitNucP4();
248  double ENi = p4Ni.E();
249  // On-shell mass of initial struck nucleon
250  double mNi = interaction.InitState().Tgt().HitNucMass();
251  // On-shell initial struck nucleon energy
252  double ENi_on_shell = std::sqrt( mNi*mNi + p4Ni.Vect().Mag2() );
253  // Energy needed to put initial nucleon on the mass shell
254  double epsilon_B = ENi_on_shell - ENi;
255 
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;
259 }
260 
262  const genie::NuclearModelI& nucl_model, double& Eb,
263  genie::QELEvGen_BindingMode_t hitNucleonBindingMode)
264 {
265  genie::Target* tgt = interaction.InitState().TgtPtr();
266  TLorentzVector* p4Ni = tgt->HitNucP4Ptr();
267 
268  // Initial nucleon 3-momentum (lab frame)
269  TVector3 p3Ni = nucl_model.Momentum3();
270 
271  // Look up the (on-shell) mass of the initial nucleon
272  TDatabasePDG* tb = TDatabasePDG::Instance();
273  double mNi = tb->GetParticle( tgt->HitNucPdg() )->Mass();
274 
275  // Set the (possibly off-shell) initial nucleon energy based on
276  // the selected binding energy mode. Always put the initial nucleon
277  // on shell if it is not part of a composite nucleus
278  double ENi = 0.;
279  if ( tgt->IsNucleus() && hitNucleonBindingMode != genie::kOnShell ) {
280 
281  // For a nuclear target with a bound initial struck nucleon, take binding
282  // energy effects and Pauli blocking into account when computing QE
283  // differential cross sections
284  interaction.ResetBit( kIAssumeFreeNucleon );
285 
286  // Initial nucleus mass
287  double Mi = tgt->Mass();
288 
289  // Final nucleus mass
290  double Mf = 0.;
291 
292  // If we use the removal energy reported by the nuclear
293  // model, then it implies a certain value for the final
294  // nucleus mass
295  if ( hitNucleonBindingMode == genie::kUseNuclearModel ) {
296  Eb = nucl_model.RemovalEnergy();
297  // This equation is the definition that we assume
298  // here for the "removal energy" (Eb) returned by the
299  // nuclear model. It matches GENIE's convention for
300  // the Bodek/Ritchie Fermi gas model.
301  Mf = Mi + Eb - mNi;
302  }
303  // We can also assume that the final nucleus is in its
304  // ground state. In this case, we can just look up its
305  // mass from the standard table. This implies a particular
306  // binding energy for the hit nucleon.
307  else if ( hitNucleonBindingMode == genie::kUseGroundStateRemnant ) {
308  // Determine the mass and proton numbers for the remnant nucleus
309  int Af = tgt->A() - 1;
310  int Zf = tgt->Z();
311  if ( genie::pdg::IsProton( tgt->HitNucPdg()) ) --Zf;
312  Mf = genie::PDGLibrary::Instance()->Find( genie::pdg::IonPdgCode(Af, Zf) )->Mass();
313 
314  // Deduce the binding energy from the final nucleus mass
315  Eb = Mf - Mi + mNi;
316  }
317  // A third option is to assign the binding energy and final nuclear
318  // mass in a way that is equivalent to the original Valencia model
319  // treatment
320  else if ( hitNucleonBindingMode == genie::kValenciaStyleQValue ) {
321  // Compute the Q-value needed for the ground-state-to-ground-state
322  // transition without nucleon removal (see Sec. III B of
323  // https://arxiv.org/abs/nucl-th/0408005). Note that this will be
324  // zero for NC and EM scattering since the proton and nucleon numbers
325  // are unchanged. Also get Q_LFG, the difference in Fermi energies
326  // between the final and initial struck nucleon species. This will
327  // likewise be zero for NC and EM interactions.
328  const genie::ProcessInfo& info = interaction.ProcInfo();
329  double Qvalue = 0.;
330  double Q_LFG = 0.;
331  if ( info.IsWeakCC() ) {
332  // Without nucleon removal, the final nucleon number remains the same
333  int Af = tgt->A();
334  // For CC interactions, the proton number will change by one. Choose
335  // the right change by checking whether we're working with a neutrino
336  // or an antineutrino.
337  int Zf = tgt->Z();
338  int probe_pdg = interaction.InitState().ProbePdg();
339  if ( genie::pdg::IsAntiNeutrino(probe_pdg) ) {
340  --Zf;
341  }
342  else if ( genie::pdg::IsNeutrino(probe_pdg) ) {
343  ++Zf;
344  }
345  else {
346  LOG( "QELEvent", pFATAL ) << "Unhandled probe PDG code " << probe_pdg
347  << " encountered for a CC interaction in"
348  << " genie::utils::BindHitNucleon()";
349  gAbortingInErr = true;
350  std::exit( 1 );
351  }
352 
353  // We have what we need to get the Q-value. Get the final nuclear
354  // mass (without nucleon removal)
355  double mf_keep_nucleon = genie::PDGLibrary::Instance()
356  ->Find( genie::pdg::IonPdgCode(Af, Zf) )->Mass();
357 
358  Qvalue = mf_keep_nucleon - Mi;
359 
360  // Get the Fermi energies for the initial and final nucleons. Include
361  // the radial dependence if using the LFG.
362  double hit_nucleon_radius = tgt->HitNucPosition();
363 
364  // Average of the proton and neutron masses. It may actually be better
365  // to use the exact on-shell masses here. However, the original paper
366  // uses the same value for protons and neutrons when computing the
367  // difference in Fermi energies. For consistency with the Valencia
368  // model publication, I'll do the same.
369  const double mN = genie::constants::kNucleonMass;
370 
371  double kF_Ni = nucl_model.LocalFermiMomentum( *tgt,
372  tgt->HitNucPdg(), hit_nucleon_radius );
373  double EFermi_Ni = std::sqrt( std::max(0., mN*mN + kF_Ni*kF_Ni) );
374 
375  double kF_Nf = nucl_model.LocalFermiMomentum( *tgt,
376  interaction.RecoilNucleonPdg(), hit_nucleon_radius );
377  // (On-shell) final nucleon mass
378  //double mNf = interaction.RecoilNucleon()->Mass();
379  double EFermi_Nf = std::sqrt( std::max(0., mN*mN + kF_Nf*kF_Nf) );
380 
381  // The difference in Fermi energies is Q_LFG
382  Q_LFG = EFermi_Nf - EFermi_Ni;
383  }
384  // Fail here for interactions that aren't one of EM/NC/CC. This is
385  // intended to help avoid silently doing the wrong thing in the future.
386  else if ( !info.IsEM() && !info.IsWeakNC() ) {
387  LOG( "QELEvent", pFATAL ) << "Unhandled process type \"" << info
388  << "\" encountered in genie::utils::BindHitNucleon()";
389  gAbortingInErr = true;
390  std::exit( 1 );
391  }
392 
393  // On-shell total energy of the initial struck nucleon
394  double ENi_OnShell = std::sqrt( std::max(0., mNi*mNi + p3Ni.Mag2()) );
395 
396  // Total energy of the remnant nucleus (with the hit nucleon removed)
397  double Ef = Mi - ENi_OnShell + Qvalue - Q_LFG;
398 
399  // Mass and kinetic energy of the remnant nucleus (with the hit nucleon
400  // removed)
401  Mf = std::sqrt( std::max(0., Ef*Ef - p3Ni.Mag2()) );
402 
403  // Deduce the binding energy from the final nucleus mass
404  Eb = Mf - Mi + mNi;
405 
406  LOG( "QELEvent", pDEBUG ) << "Qvalue = " << Qvalue
407  << ", Q_LFG = " << Q_LFG << " at radius = " << tgt->HitNucPosition();
408  }
409 
410  // The (lab-frame) off-shell initial nucleon energy is the difference
411  // between the lab frame total energies of the initial and remnant nuclei
412  ENi = Mi - std::sqrt( Mf*Mf + p3Ni.Mag2() );
413  }
414  else {
415  // Keep the struck nucleon on shell either because
416  // hitNucleonBindingMode == kOnShell or because
417  // the target is a single nucleon
418  ENi = std::sqrt( p3Ni.Mag2() + std::pow(mNi, 2) );
419  Eb = 0.;
420 
421  // If we're dealing with a nuclear target but using the on-shell
422  // binding mode, set the "assume free nucleon" flag. This turns off
423  // Pauli blocking and the requirement that q0 > 0 in the QE cross section
424  // models (an on-shell nucleon *is* a free nucleon)
425  if ( tgt->IsNucleus() ) interaction.SetBit( kIAssumeFreeNucleon );
426  }
427 
428  LOG( "QELEvent", pDEBUG ) << "Eb = " << Eb << ", pNi = " << p3Ni.Mag()
429  << ", ENi = " << ENi;
430 
431  // Update the initial nucleon lab-frame 4-momentum in the interaction with
432  // its current components
433  p4Ni->SetVect( p3Ni );
434  p4Ni->SetE( ENi );
435 
436 }
Cross Section Calculation Interface.
const KPhaseSpace & PhaseSpace(void) const
Definition: Interaction.h:73
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)
Definition: PDGUtils.cxx:107
double RemovalEnergy(void) const
Definition: NuclearModelI.h:65
static const double kNucleonMass
Definition: Constants.h:77
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1064
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
TParticlePDG * RecoilNucleon(void) const
recoil nucleon
std::string string
Definition: nybbler.cc:12
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
int A(void) const
Definition: Target.h:70
A simple [min,max] interval for doubles.
Definition: Range1.h:42
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
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
double Mass(Resonance_t res)
resonance mass (GeV)
void BindHitNucleon(Interaction &interaction, const NuclearModelI &nucl_model, double &Eb, DMELEvGen_BindingMode_t hitNucleonBindingMode)
Definition: DMELUtils.cxx:259
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)
Definition: QELUtils.cxx:93
Range1D_t Q2Lim(void) const
Q2 limits.
double Mass(void) const
Definition: Target.cxx:224
enum genie::EQELEvGenBindingMode QELEvGen_BindingMode_t
const TVector3 & Momentum3(void) const
Definition: NuclearModelI.h:75
double Qvalue(int targetpdg, int nupdg)
Definition: MECUtils.cxx:164
Summary information for an interaction.
Definition: Interaction.h:56
const TLorentzVector & HitNucP4(void) const
Definition: Target.h:91
bool IsProton(int pdgc)
Definition: PDGUtils.cxx:333
bool IsWeakNC(void) const
const TLorentzVector & FSLeptonP4(void) const
Definition: Kinematics.h:65
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void SetFSLeptonP4(const TLorentzVector &p4)
Definition: Kinematics.cxx:297
double EnergyDeltaFunctionSolutionQEL(const Interaction &inter)
Definition: QELUtils.cxx:50
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
bool IsAntiNeutrino(int pdgc)
Definition: PDGUtils.cxx:115
const Kinematics & Kine(void) const
Definition: Interaction.h:71
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
int ProbePdg(void) const
Definition: InitialState.h:64
int Z(void) const
Definition: Target.h:68
static const double kASmallNum
Definition: Controls.h:40
static int max(int a, int b)
double gamma(double KE, const simb::MCParticle *part)
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
bool IsEM(void) const
double max
Definition: Range1.h:53
TLorentzVector * HitNucP4Ptr(void) const
Definition: Target.cxx:247
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:57
double CosTheta0Max(const genie::Interaction &interaction)
Definition: DMELUtils.cxx:217
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
double CMEnergy() const
centre-of-mass energy (sqrt s)
Target * TgtPtr(void) const
Definition: InitialState.h:67
int IonPdgCode(int A, int Z)
Definition: PDGUtils.cxx:68
void SetHadSystP4(const TLorentzVector &p4)
Definition: Kinematics.cxx:307
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:75
const Target & Tgt(void) const
Definition: InitialState.h:66
bool gAbortingInErr
Definition: Messenger.cxx:34
double ProbeE(RefFrame_t rf) const
static const double kPi
Definition: Constants.h:37
static QCString * s
Definition: config.cpp:1042
TLorentzVector * GetProbeP4(RefFrame_t rf=kRfHitNucRest) const
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63