eic-smear  1.0.3
A collection of ROOT classes for Monte Carlo events and a fast-smearing code simulating detector effects for the Electron-Ion Collider task force
erhic::EventPythia Class Reference

Pythia 6 DIS event. More...

#include <EventPythia.h>

Inheritance diagram for erhic::EventPythia:
erhic::EventMC erhic::EventDis erhic::VirtualEvent erhic::EventBeagle

Public Member Functions

 EventPythia (const std::string &str="")
 Constructor. More...
 
virtual ~EventPythia ()
 Destructor. More...
 
virtual bool Parse (const std::string &)
 Parses the event information from a text string. More...
 
virtual void SetNucleon (int n)
 Sets the nucleon species. More...
 
virtual void SetTargetParton (int n)
 Sets the target parton species. More...
 
virtual void SetBeamParton (int n)
 Sets the beam parton species. More...
 
virtual void SetGenEvent (int n)
 Sets the number of trials required to generate this event. More...
 
virtual void SetTargetPartonX (double xB)
 Sets the x of the target parton. More...
 
virtual void SetBeamPartonX (double xB)
 Sets the x of the beam parton. More...
 
virtual void SetBeamPartonTheta (double radians)
 Sets the polar angle of the beam parton in the cm frame. More...
 
virtual void SetLeptonPhi (double radians)
 Azimuthal angle of the scattered lepton in the cm frame. More...
 
virtual void SetF1 (double f1)
 Used for radiative corrections. More...
 
virtual void SetF2 (double f2)
 Used for radiative corrections. More...
 
virtual void SetSigmaRad (double sr)
 Used for radiative corrections. More...
 
virtual void SetHardS (double s)
 Sets the Mandelstamm s of the hard interaction. More...
 
virtual void SetHardT (double t)
 Sets the Mandelstamm t of the hard interaction. More...
 
virtual void SetHardU (double u)
 Sets the Mandelstamm u of the hard interaction. More...
 
virtual void SetHardQ2 (double Q2)
 Sets the Q2 of the hard interaction. More...
 
virtual void SetHardPt2 (double pt2)
 Sets the squared pT of the hard interaction. More...
 
virtual void SetSigRadCor (double s)
 Used for radiative corrections. More...
 
virtual void SetEBrems (double e)
 Sets the energy per radiative photon in the nuclear rest frame. More...
 
virtual void SetPhotonFlux (double f)
 Flux factor, see VINT(319) in PYTHIA 6. More...
 
virtual void SetTrueY (double inelasticity)
 Sets the true (not reconstructed) value for inelasticity. More...
 
virtual void SetTrueQ2 (double Q2)
 Sets the true (not reconstructed) value for Q2. More...
 
virtual void SetTrueX (double x)
 Sets the true (not reconstructed) value for x. More...
 
virtual void SetTrueW2 (double W2)
 Sets the true (not reconstructed) value for W2. More...
 
virtual void SetTrueNu (double Nu)
 Sets the true (not reconstructed) value for nu. More...
 
virtual void SetR (double r)
 Used for radiative corrections. More...
 
virtual int GetGenEvent () const
 Returns the number of trials required to generate this event. More...
 
virtual double GetTargetPartonX () const
 Returns the x of the target parton. More...
 
virtual double GetBeamPartonX () const
 Returns the x of the beam parton. More...
 
virtual double GetBeamPartonTheta () const
 Returns the polar angle of the beam parton in the cm frame, in radians in the range [0, pi]. More...
 
virtual double GetLeptonPhi () const
 Returns the azimuthal angle of the scattered lepton. More...
 
virtual double GetF1 () const
 Used for radiative corrections. More...
 
virtual double GetF2 () const
 Used for radiative corrections. More...
 
virtual double GetSigmaRad () const
 Used for radiative corrections. More...
 
virtual double GetHardS () const
 Returns the Mandelstamm s of the hard interaction. More...
 
virtual double GetHardT () const
 Returns the Mandelstamm t of the hard interaction. More...
 
virtual double GetHardU () const
 Returns the Mandelstamm u of the hard interaction. More...
 
virtual double GetHardQ2 () const
 Returns the Q2 of the hard interaction. More...
 
virtual double GetHardPt2 () const
 Returns the squared pT of the hard interaction. More...
 
virtual double GetSigRadCor () const
 Used for radiative corrections. More...
 
virtual double GetEBrems () const
 Returnss the energy per radiative photon in the nuclear rest frame. More...
 
virtual double GetPhotonFlux () const
 Returns the flux factor, see VINT(319) in PYTHIA 6. More...
 
virtual double GetTrueY () const
 Sets the true (not reconstructed) value for inelasticity. More...
 
virtual double GetTrueQ2 () const
 Sets the true (not reconstructed) value for Q2. More...
 
virtual double GetTrueX () const
 Sets the true (not reconstructed) value for x. More...
 
virtual double GetTrueW2 () const
 Sets the true (not reconstructed) value for W2. More...
 
virtual double GetTrueNu () const
 Sets the true (not reconstructed) value for ν. More...
 
virtual double GetR () const
 Used for radiative corrections. More...
 
virtual const ParticleMCScatteredLepton () const
 Returns the scattered lepton. More...
 
- Public Member Functions inherited from erhic::EventMC
 EventMC ()
 Constructor. More...
 
virtual ~EventMC ()
 Destructor. More...
 
virtual bool RequiresEaParticleFields ()
 
virtual ULong64_t GetN () const
 Returns a unique identifier for this event. More...
 
virtual Int_t GetProcess () const
 Returns a code describing the production process of this event. More...
 
virtual UInt_t GetNTracks () const
 Returns the number of tracks in the event. More...
 
virtual const ParticleMCGetTrack (UInt_t) const
 Returns the nth track. More...
 
virtual ParticleMCGetTrack (UInt_t)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts. More...
 
virtual const ParticleMCBeamLepton () const
 Returns a pointer to the incident lepton beam particle. More...
 
virtual const ParticleMCBeamHadron () const
 Returns a pointer to the incident hadron beam particle. More...
 
virtual const ParticleMCExchangeBoson () const
 Returns a pointer to the exchanged boson. More...
 
virtual void AddLast (ParticleMC *track)
 Add a copy of a track argument to the end of the track list. More...
 
virtual void Reset ()
 Resets event properties to defaults. More...
 
virtual void Clear (Option_t *="")
 Clears event contents. More...
 
virtual void SetProcess (int code)
 Sets the code describing the production process of this event. More...
 
virtual void SetN (int n)
 Sets the unique identifier for this event. More...
 
virtual void SetNTracks (int n)
 Sets the track count for this event. More...
 
virtual void SetELeptonInNuclearFrame (double energy)
 Set incident lepton energy in the nuclear rest frame. More...
 
virtual void SetEScatteredInNuclearFrame (double energy)
 Set scattered lepton energy in the nuclear rest frame. More...
 
void FinalState (ParticlePtrList &particles) const
 Stores pointers to all final state particles in the list. More...
 
void HadronicFinalState (ParticlePtrList &) const
 Yields all particles that belong to the hadronic final state. More...
 
TLorentzVector FinalStateMomentum () const
 Returns the total momentum of the final state in GeV/c. More...
 
TLorentzVector HadronicFinalStateMomentum () const
 Returns the total momentum of the hadronic final state in GeV/c. More...
 
Double_t FinalStateCharge () const
 Returns the total charge of the final state in units of e. More...
 
std::vector< const VirtualParticle * > GetTracks () const
 Returns pointers to all tracks in the event. More...
 
- Public Member Functions inherited from erhic::EventDis
virtual ~EventDis ()
 Destructor. More...
 
 EventDis ()
 Default constructor. More...
 
 EventDis (const EventDis &)
 Constructor copying another event's kinematics. More...
 
EventDisoperator= (const EventDis &)
 Assign another event's kinematics to this EventDis. More...
 
virtual Double_t GetX () const
 Returns Bjorken-x of the event. More...
 
virtual Double_t GetQ2 () const
 Returns the four-momentum transfer (exchange boson mass) Q2. More...
 
virtual Double_t GetY () const
 Returns the event inelasticity. More...
 
virtual Double_t GetYPlus () const
 Returns Y+ = y2 / (1 + (1-y)2) More...
 
virtual Double_t GetW2 () const
 Returns the invariant mass of the hadronic final state. More...
 
virtual Double_t GetNu () const
 Returns the exchange boson energy in the beam hadron rest frame. More...
 
virtual double GetXDoubleAngle () const
 Returns Bjorken x computed via the double-angle method. More...
 
virtual double GetQ2DoubleAngle () const
 Returns Q-squared computed via the double-angle method. More...
 
virtual double GetYDoubleAngle () const
 Returns inelasticity computed via the double-angle method. More...
 
virtual double GetW2DoubleAngle () const
 Returns W-squared computed via the double-angle method. More...
 
virtual double GetXJacquetBlondel () const
 Returns Bjorken x computed via the Jacquet-Blondel method. More...
 
virtual double GetQ2JacquetBlondel () const
 Returns Q-squared computed via the Jacquet-Blondel method. More...
 
virtual double GetYJacquetBlondel () const
 Returns inelasticity computed via the Jacquet-Blondel method. More...
 
virtual double GetW2JacquetBlondel () const
 Returns W-squared computed via the Jacquet-Blondel method. More...
 
virtual void SetLeptonKinematics (const DisKinematics &)
 Set the kinematics computed from the scattered lepton. More...
 
virtual void SetJacquetBlondelKinematics (const DisKinematics &)
 Set the kinematics computed from the Jacquet-Blondel method. More...
 
virtual void SetDoubleAngleKinematics (const DisKinematics &)
 Set the kinematics computed from the double-angle method. More...
 
virtual void CopyKinematics (const EventDis &)
 Copy the kinematics from another event to this event. More...
 
- Public Member Functions inherited from erhic::VirtualEvent
virtual ~VirtualEvent ()
 Destructor. More...
 

Public Attributes

Int_t nucleon
 PDG code of the hadron beam, see MSTI(12) More...
 
Int_t tgtparton
 PDG code of the struck parton in the hadron beam, see MSTI(16) More...
 
Int_t beamparton
 Parton interacting with the hadron beam in the case of resolved photon processes and soft VMD, see MSTI(15) More...
 
Int_t genevent
 Trials required for this event. More...
 
Double32_t xtgtparton
 Momentum fraction taken by the target parton, see PARI(34) More...
 
Double32_t xbeamparton
 Momentum fraction taken by the beam parton, see PARI(33) More...
 
Double32_t thetabeamparton
 Polar angle of the beam parton in the cm frame, between 0 and pi radians, see PARI(53) More...
 
Double32_t leptonphi
 Azimuthal angle of the scattered lepton in the cm frame. More...
 
Double32_t F1
 Value used for radiative corrections. More...
 
Double32_t sigma_rad
 Value used for radiative corrections. More...
 
Double32_t t_hat
 Mandelstam t of the hard subprocess, see PARI(15) More...
 
Double32_t u_hat
 Mandelstam u of the hard subprocess, see PARI(16) More...
 
Double32_t Q2_hat
 Q2 of the hard subprocess, see PARI(22) More...
 
Double32_t SigRadCor
 Value used for radiative corrections. More...
 
Double32_t EBrems
 Energy per radiative photon in the nuclear rest frame. More...
 
Double32_t photonflux
 Flux factor, see VINT(319) More...
 
Double32_t trueY
 Generated y of the event, see VINT(309) More...
 
Double32_t trueQ2
 Generated Q2 of the event, see VINT(307) More...
 
Double32_t trueX
 Generated x of the event. More...
 
Double32_t trueW2
 Generated W2 of the event. More...
 
Double32_t trueNu
 Generated nu of the event. More...
 
Double32_t F2
 Value used for radiative corrections. More...
 
Double32_t R
 Value used for radiative corrections. More...
 
Double32_t pt2_hat
 Squared pT of the hard subprocess, see PARI(18) More...
 
Double32_t sHat
 Mandelstam s of the hard subprocess, see PARI(14) More...
 
- Public Attributes inherited from erhic::EventDis
Double32_t x
 Bjorken scaling variable. More...
 
Double32_t QSquared
 Q2 calculated from scattered electron. More...
 
Double32_t y
 Inelasticity. More...
 
Double32_t WSquared
 Invariant mass of the hadronic system. More...
 
Double32_t nu
 Energy transfer from the electron. More...
 
Double32_t yJB
 y calculated via the Jacquet-Blondel method More...
 
Double32_t QSquaredJB
 Q2 calculated via the Jacquet-Blondel method. More...
 
Double32_t xJB
 x calculated via the Jacquet-Blondel method More...
 
Double32_t WSquaredJB
 W2 calculated via the Jacquet-Blondel method. More...
 
Double32_t yDA
 y calculated via the double-angle method More...
 
Double32_t QSquaredDA
 Q2 calculated via the double-angle method. More...
 
Double32_t xDA
 x calculated via the double-angle method More...
 
Double32_t WSquaredDA
 W2 calculated via the double-angle method. More...
 

Additional Inherited Members

- Public Types inherited from erhic::VirtualEvent
typedef std::vector< const erhic::VirtualParticle * > ParticlePtrList
 typedef for a track pointer collection. More...
 
- Protected Attributes inherited from erhic::EventMC
Int_t number
 Event number. More...
 
Int_t process
 PYTHIA code for the physics process producing the event. More...
 
Int_t nTracks
 Number of Particles in the event (intermediate + final) More...
 
Double32_t ELeptonInNucl
 Incident lepton energy in the nuclear rest frame. More...
 
Double32_t ELeptonOutNucl
 Scattered lepton energy in the nuclear rest frame. More...
 
TClonesArray particles
 Particle list. More...
 

Detailed Description

Pythia 6 DIS event.

Describes an event from the EIC PYTHIA6 implementation.

Todo:
Change sHat/t_hat/u_hat naming to be consistent

Definition at line 28 of file erhic/EventPythia.h.

Constructor & Destructor Documentation

◆ EventPythia()

erhic::EventPythia::EventPythia ( const std::string &  str = "")
explicit

Constructor.

Parameters
[in]strA text string setting event-wise quantities. The string format should be (no newlines): "I ievent genevent subprocess nucleon targetparton xtargparton beamparton xbeamparton thetabeamprtn truey trueQ2 truex trueW2 trueNu leptonphi s_hat t_hat u_hat pt2_hat Q2_hat F2 F1 R sigma_rad SigRadCor EBrems photonflux nrTracks"

Definition at line 20 of file erhic/EventPythia.cxx.

◆ ~EventPythia()

erhic::EventPythia::~EventPythia ( )
virtual

Destructor.

Definition at line 48 of file erhic/EventPythia.cxx.

Member Function Documentation

◆ GetBeamPartonTheta()

double erhic::EventPythia::GetBeamPartonTheta ( ) const
inlinevirtual

Returns the polar angle of the beam parton in the cm frame, in radians in the range [0, pi].

Definition at line 469 of file erhic/EventPythia.h.

◆ GetBeamPartonX()

double erhic::EventPythia::GetBeamPartonX ( ) const
inlinevirtual

Returns the x of the beam parton.

Definition at line 465 of file erhic/EventPythia.h.

◆ GetEBrems()

double erhic::EventPythia::GetEBrems ( ) const
inlinevirtual

Returnss the energy per radiative photon in the nuclear rest frame.

Definition at line 513 of file erhic/EventPythia.h.

◆ GetF1()

double erhic::EventPythia::GetF1 ( ) const
inlinevirtual

Used for radiative corrections.

Definition at line 477 of file erhic/EventPythia.h.

◆ GetF2()

double erhic::EventPythia::GetF2 ( ) const
inlinevirtual

Used for radiative corrections.

Definition at line 481 of file erhic/EventPythia.h.

◆ GetGenEvent()

int erhic::EventPythia::GetGenEvent ( ) const
inlinevirtual

Returns the number of trials required to generate this event.

Definition at line 457 of file erhic/EventPythia.h.

◆ GetHardPt2()

double erhic::EventPythia::GetHardPt2 ( ) const
inlinevirtual

Returns the squared pT of the hard interaction.

Definition at line 505 of file erhic/EventPythia.h.

◆ GetHardQ2()

double erhic::EventPythia::GetHardQ2 ( ) const
inlinevirtual

Returns the Q2 of the hard interaction.

Definition at line 501 of file erhic/EventPythia.h.

◆ GetHardS()

double erhic::EventPythia::GetHardS ( ) const
inlinevirtual

Returns the Mandelstamm s of the hard interaction.

Definition at line 489 of file erhic/EventPythia.h.

◆ GetHardT()

double erhic::EventPythia::GetHardT ( ) const
inlinevirtual

Returns the Mandelstamm t of the hard interaction.

Definition at line 493 of file erhic/EventPythia.h.

◆ GetHardU()

double erhic::EventPythia::GetHardU ( ) const
inlinevirtual

Returns the Mandelstamm u of the hard interaction.

Definition at line 497 of file erhic/EventPythia.h.

◆ GetLeptonPhi()

double erhic::EventPythia::GetLeptonPhi ( ) const
inlinevirtual

Returns the azimuthal angle of the scattered lepton.

Angle is given in the centre-of-mass frame, in radians in the range [0, 2pi).

Definition at line 473 of file erhic/EventPythia.h.

◆ GetPhotonFlux()

double erhic::EventPythia::GetPhotonFlux ( ) const
inlinevirtual

Returns the flux factor, see VINT(319) in PYTHIA 6.

Definition at line 517 of file erhic/EventPythia.h.

◆ GetR()

double erhic::EventPythia::GetR ( ) const
inlinevirtual

Used for radiative corrections.

Definition at line 541 of file erhic/EventPythia.h.

◆ GetSigmaRad()

double erhic::EventPythia::GetSigmaRad ( ) const
inlinevirtual

Used for radiative corrections.

Definition at line 485 of file erhic/EventPythia.h.

◆ GetSigRadCor()

double erhic::EventPythia::GetSigRadCor ( ) const
inlinevirtual

Used for radiative corrections.

Definition at line 509 of file erhic/EventPythia.h.

◆ GetTargetPartonX()

double erhic::EventPythia::GetTargetPartonX ( ) const
inlinevirtual

Returns the x of the target parton.

Definition at line 461 of file erhic/EventPythia.h.

◆ GetTrueNu()

double erhic::EventPythia::GetTrueNu ( ) const
inlinevirtual

Sets the true (not reconstructed) value for ν.

Definition at line 537 of file erhic/EventPythia.h.

◆ GetTrueQ2()

double erhic::EventPythia::GetTrueQ2 ( ) const
inlinevirtual

Sets the true (not reconstructed) value for Q2.

Definition at line 525 of file erhic/EventPythia.h.

◆ GetTrueW2()

double erhic::EventPythia::GetTrueW2 ( ) const
inlinevirtual

Sets the true (not reconstructed) value for W2.

Definition at line 533 of file erhic/EventPythia.h.

◆ GetTrueX()

double erhic::EventPythia::GetTrueX ( ) const
inlinevirtual

Sets the true (not reconstructed) value for x.

Definition at line 529 of file erhic/EventPythia.h.

◆ GetTrueY()

double erhic::EventPythia::GetTrueY ( ) const
inlinevirtual

Sets the true (not reconstructed) value for inelasticity.

Definition at line 521 of file erhic/EventPythia.h.

◆ Parse()

bool erhic::EventPythia::Parse ( const std::string &  line)
virtual

Parses the event information from a text string.

See the constructor for the string format. Returns true in the event of a successful read operation, false in case of an error.

Implements erhic::EventMC.

Reimplemented in erhic::EventBeagle.

Definition at line 50 of file erhic/EventPythia.cxx.

◆ ScatteredLepton()

const ParticleMC * erhic::EventPythia::ScatteredLepton ( ) const
virtual

Returns the scattered lepton.

This is the first final state particle with the same species as the beam lepton and parent index equal to three (counting index from 1).

Reimplemented from erhic::EventMC.

Definition at line 71 of file erhic/EventPythia.cxx.

◆ SetBeamParton()

void erhic::EventPythia::SetBeamParton ( int  n)
inlinevirtual

Sets the beam parton species.

Parameters
n[in] PDG code of the parton interacting with the hadron beam in the case of resolved photon processes and soft VMD, see MSTI(15)

Definition at line 365 of file erhic/EventPythia.h.

◆ SetBeamPartonTheta()

void erhic::EventPythia::SetBeamPartonTheta ( double  radians)
inlinevirtual

Sets the polar angle of the beam parton in the cm frame.

Todo:
enforce range [0, pi] when setting

Definition at line 381 of file erhic/EventPythia.h.

◆ SetBeamPartonX()

void erhic::EventPythia::SetBeamPartonX ( double  xB)
inlinevirtual

Sets the x of the beam parton.

Definition at line 377 of file erhic/EventPythia.h.

◆ SetEBrems()

void erhic::EventPythia::SetEBrems ( double  e)
inlinevirtual

Sets the energy per radiative photon in the nuclear rest frame.

Definition at line 425 of file erhic/EventPythia.h.

◆ SetF1()

void erhic::EventPythia::SetF1 ( double  f1)
inlinevirtual

Used for radiative corrections.

Definition at line 389 of file erhic/EventPythia.h.

◆ SetF2()

void erhic::EventPythia::SetF2 ( double  f2)
inlinevirtual

Used for radiative corrections.

Definition at line 393 of file erhic/EventPythia.h.

◆ SetGenEvent()

void erhic::EventPythia::SetGenEvent ( int  n)
inlinevirtual

Sets the number of trials required to generate this event.

Parameters
n[in] The number of trials

Definition at line 369 of file erhic/EventPythia.h.

◆ SetHardPt2()

void erhic::EventPythia::SetHardPt2 ( double  pt2)
inlinevirtual

Sets the squared pT of the hard interaction.

Definition at line 417 of file erhic/EventPythia.h.

◆ SetHardQ2()

void erhic::EventPythia::SetHardQ2 ( double  Q2)
inlinevirtual

Sets the Q2 of the hard interaction.

Definition at line 413 of file erhic/EventPythia.h.

◆ SetHardS()

void erhic::EventPythia::SetHardS ( double  s)
inlinevirtual

Sets the Mandelstamm s of the hard interaction.

Definition at line 401 of file erhic/EventPythia.h.

◆ SetHardT()

void erhic::EventPythia::SetHardT ( double  t)
inlinevirtual

Sets the Mandelstamm t of the hard interaction.

Definition at line 405 of file erhic/EventPythia.h.

◆ SetHardU()

void erhic::EventPythia::SetHardU ( double  u)
inlinevirtual

Sets the Mandelstamm u of the hard interaction.

Definition at line 409 of file erhic/EventPythia.h.

◆ SetLeptonPhi()

void erhic::EventPythia::SetLeptonPhi ( double  radians)
inlinevirtual

Azimuthal angle of the scattered lepton in the cm frame.

Todo:
enforce range [0, 2pi) when setting

Definition at line 385 of file erhic/EventPythia.h.

◆ SetNucleon()

void erhic::EventPythia::SetNucleon ( int  n)
inlinevirtual

Sets the nucleon species.

Parameters
n[in] PDG code of the hadron beam, see MSTI(12)

Definition at line 357 of file erhic/EventPythia.h.

◆ SetPhotonFlux()

void erhic::EventPythia::SetPhotonFlux ( double  f)
inlinevirtual

Flux factor, see VINT(319) in PYTHIA 6.

Definition at line 429 of file erhic/EventPythia.h.

◆ SetR()

void erhic::EventPythia::SetR ( double  r)
inlinevirtual

Used for radiative corrections.

Definition at line 453 of file erhic/EventPythia.h.

◆ SetSigmaRad()

void erhic::EventPythia::SetSigmaRad ( double  sr)
inlinevirtual

Used for radiative corrections.

Definition at line 397 of file erhic/EventPythia.h.

◆ SetSigRadCor()

void erhic::EventPythia::SetSigRadCor ( double  s)
inlinevirtual

Used for radiative corrections.

Definition at line 421 of file erhic/EventPythia.h.

◆ SetTargetParton()

void erhic::EventPythia::SetTargetParton ( int  n)
inlinevirtual

Sets the target parton species.

Parameters
n[in] PDG code of the struck parton in the hadron beam, see MSTI(16)

Definition at line 361 of file erhic/EventPythia.h.

◆ SetTargetPartonX()

void erhic::EventPythia::SetTargetPartonX ( double  xB)
inlinevirtual

Sets the x of the target parton.

Definition at line 373 of file erhic/EventPythia.h.

◆ SetTrueNu()

void erhic::EventPythia::SetTrueNu ( double  Nu)
inlinevirtual

Sets the true (not reconstructed) value for nu.

Definition at line 449 of file erhic/EventPythia.h.

◆ SetTrueQ2()

void erhic::EventPythia::SetTrueQ2 ( double  Q2)
inlinevirtual

Sets the true (not reconstructed) value for Q2.

Definition at line 437 of file erhic/EventPythia.h.

◆ SetTrueW2()

void erhic::EventPythia::SetTrueW2 ( double  W2)
inlinevirtual

Sets the true (not reconstructed) value for W2.

Definition at line 445 of file erhic/EventPythia.h.

◆ SetTrueX()

void erhic::EventPythia::SetTrueX ( double  x)
inlinevirtual

Sets the true (not reconstructed) value for x.

Definition at line 441 of file erhic/EventPythia.h.

◆ SetTrueY()

void erhic::EventPythia::SetTrueY ( double  inelasticity)
inlinevirtual

Sets the true (not reconstructed) value for inelasticity.

Definition at line 433 of file erhic/EventPythia.h.

Member Data Documentation

◆ beamparton

Int_t erhic::EventPythia::beamparton

Parton interacting with the hadron beam in the case of resolved photon processes and soft VMD, see MSTI(15)

Definition at line 316 of file erhic/EventPythia.h.

◆ EBrems

Double32_t erhic::EventPythia::EBrems

Energy per radiative photon in the nuclear rest frame.

Definition at line 338 of file erhic/EventPythia.h.

◆ F1

Double32_t erhic::EventPythia::F1

Value used for radiative corrections.

Definition at line 329 of file erhic/EventPythia.h.

◆ F2

Double32_t erhic::EventPythia::F2

Value used for radiative corrections.

Definition at line 348 of file erhic/EventPythia.h.

◆ genevent

Int_t erhic::EventPythia::genevent

Trials required for this event.

Definition at line 319 of file erhic/EventPythia.h.

◆ leptonphi

Double32_t erhic::EventPythia::leptonphi

Azimuthal angle of the scattered lepton in the cm frame.

Definition at line 327 of file erhic/EventPythia.h.

◆ nucleon

Int_t erhic::EventPythia::nucleon

PDG code of the hadron beam, see MSTI(12)

Definition at line 312 of file erhic/EventPythia.h.

◆ photonflux

Double32_t erhic::EventPythia::photonflux

Flux factor, see VINT(319)

Definition at line 340 of file erhic/EventPythia.h.

◆ pt2_hat

Double32_t erhic::EventPythia::pt2_hat

Squared pT of the hard subprocess, see PARI(18)

Definition at line 350 of file erhic/EventPythia.h.

◆ Q2_hat

Double32_t erhic::EventPythia::Q2_hat

Q2 of the hard subprocess, see PARI(22)

Definition at line 335 of file erhic/EventPythia.h.

◆ R

Double32_t erhic::EventPythia::R

Value used for radiative corrections.

Definition at line 349 of file erhic/EventPythia.h.

◆ sHat

Double32_t erhic::EventPythia::sHat

Mandelstam s of the hard subprocess, see PARI(14)

Definition at line 352 of file erhic/EventPythia.h.

◆ sigma_rad

Double32_t erhic::EventPythia::sigma_rad

Value used for radiative corrections.

Definition at line 330 of file erhic/EventPythia.h.

◆ SigRadCor

Double32_t erhic::EventPythia::SigRadCor

Value used for radiative corrections.

Definition at line 337 of file erhic/EventPythia.h.

◆ t_hat

Double32_t erhic::EventPythia::t_hat

Mandelstam t of the hard subprocess, see PARI(15)

Definition at line 331 of file erhic/EventPythia.h.

◆ tgtparton

Int_t erhic::EventPythia::tgtparton

PDG code of the struck parton in the hadron beam, see MSTI(16)

Definition at line 314 of file erhic/EventPythia.h.

◆ thetabeamparton

Double32_t erhic::EventPythia::thetabeamparton

Polar angle of the beam parton in the cm frame, between 0 and pi radians, see PARI(53)

Definition at line 324 of file erhic/EventPythia.h.

◆ trueNu

Double32_t erhic::EventPythia::trueNu

Generated nu of the event.

Definition at line 347 of file erhic/EventPythia.h.

◆ trueQ2

Double32_t erhic::EventPythia::trueQ2

Generated Q2 of the event, see VINT(307)

Definition at line 343 of file erhic/EventPythia.h.

◆ trueW2

Double32_t erhic::EventPythia::trueW2

Generated W2 of the event.

Definition at line 346 of file erhic/EventPythia.h.

◆ trueX

Double32_t erhic::EventPythia::trueX

Generated x of the event.

Definition at line 345 of file erhic/EventPythia.h.

◆ trueY

Double32_t erhic::EventPythia::trueY

Generated y of the event, see VINT(309)

Definition at line 341 of file erhic/EventPythia.h.

◆ u_hat

Double32_t erhic::EventPythia::u_hat

Mandelstam u of the hard subprocess, see PARI(16)

Definition at line 333 of file erhic/EventPythia.h.

◆ xbeamparton

Double32_t erhic::EventPythia::xbeamparton

Momentum fraction taken by the beam parton, see PARI(33)

Definition at line 322 of file erhic/EventPythia.h.

◆ xtgtparton

Double32_t erhic::EventPythia::xtgtparton

Momentum fraction taken by the target parton, see PARI(34)

Definition at line 320 of file erhic/EventPythia.h.


The documentation for this class was generated from the following files: