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
|
Go to the documentation of this file.
10 #ifndef INCLUDE_EICSMEAR_ERHIC_EVENTPYTHIA_H_
11 #define INCLUDE_EICSMEAR_ERHIC_EVENTPYTHIA_H_
52 virtual bool Parse(
const std::string&);
105 virtual void SetF1(
double f1);
110 virtual void SetF2(
double f2);
160 virtual void SetTrueY(
double inelasticity);
185 virtual void SetR(
double r);
218 virtual double GetF1()
const;
223 virtual double GetF2()
const;
298 virtual double GetR()
const;
434 trueY = inelasticity;
554 bool Parse(
const std::string&);
595 #endif // INCLUDE_EICSMEAR_ERHIC_EVENTPYTHIA_H_
Double32_t t_hat
Mandelstam t of the hard subprocess, see PARI(15)
Double32_t xtgtparton
Momentum fraction taken by the target parton, see PARI(34)
Int_t Ncollt
Number of collisions in target.
virtual void SetBeamPartonX(double xB)
Sets the x of the beam parton.
Double32_t RAevt
Nuclear PDF ratio for the up sea for the given event kinematics.
Double32_t sigma_rad
Value used for radiative corrections.
Int_t lepton
=================additional variables for BeAGLE==================
virtual void SetBeamParton(int n)
Sets the beam parton species.
EventPythia(const std::string &str="")
Constructor.
virtual void SetTrueX(double x)
Sets the true (not reconstructed) value for x.
Int_t Npevap
Number of protons from evaporation.
virtual void SetSigRadCor(double s)
Used for radiative corrections.
Double32_t sHat
Mandelstam s of the hard subprocess, see PARI(14)
Double32_t pyf
Sum fermi momentum of all inelastic participant in target rest frame z along gamma*.
Int_t nucleon
PDG code of the hadron beam, see MSTI(12)
Double32_t u_hat
Mandelstam u of the hard subprocess, see PARI(16)
Int_t Ztarg
charge number of target beam
virtual double GetBeamPartonX() const
Returns the x of the beam parton.
virtual void SetEBrems(double e)
Sets the energy per radiative photon in the nuclear rest frame.
Double32_t Q2_hat
Q2 of the hard subprocess, see PARI(22)
virtual double GetEBrems() const
Returnss the energy per radiative photon in the nuclear rest frame.
virtual double GetLeptonPhi() const
Returns the azimuthal angle of the scattered lepton.
virtual double GetR() const
Used for radiative corrections.
Double32_t crang
crossing angle (mr). crang=1000*atan(px/pz), one beam px=py=0, the other py=0
virtual double GetHardQ2() const
Returns the Q2 of the hard interaction.
Int_t NINC
Number of stable hadrons from intranuclear cascade.
virtual double GetSigRadCor() const
Used for radiative corrections.
virtual void SetTrueNu(double Nu)
Sets the true (not reconstructed) value for nu.
virtual ~EventPythia()
Destructor.
virtual double GetTrueQ2() const
Sets the true (not reconstructed) value for Q2.
Int_t Atarg
mass number of target beam
Double32_t d1st
density-weighted distance from first collision to the edge of the nucleus (amount of material travers...
Double32_t trueQ2
Generated Q2 of the event, see VINT(307)
Double32_t User1
User variables to prevent/delay future format changes.
Int_t NINCch
Number of charged stable hadrons from intranuclear cascade.
Double32_t R
Value used for radiative corrections.
Double32_t thetabeamparton
Polar angle of the beam parton in the cm frame, between 0 and pi radians, see PARI(53)
virtual const ParticleMC * ScatteredLepton() const
Returns the scattered lepton.
virtual double GetHardPt2() const
Returns the squared pT of the hard interaction.
Double32_t trueX
Generated x of the event.
virtual void SetGenEvent(int n)
Sets the number of trials required to generate this event.
virtual void SetTrueQ2(double Q2)
Sets the true (not reconstructed) value for Q2.
bool RequiresEaParticleFields()
Double32_t SigRadCor
Value used for radiative corrections.
Int_t Nwdch
Number of wounded proton including those in INC.
Double32_t User2
User variables to prevent/delay future format changes.
Double32_t b
impact parameter value
Double32_t pxf
Sum fermi momentum of all inelastic participant in target rest frame z along gamma*.
Int_t genevent
Trials required for this event.
Int_t Nnevap
Number of neutrons from evaporation.
virtual void SetHardQ2(double Q2)
Sets the Q2 of the hard interaction.
virtual int GetGenEvent() const
Returns the number of trials required to generate this event.
virtual void SetSigmaRad(double sr)
Used for radiative corrections.
virtual void SetF2(double f2)
Used for radiative corrections.
Double32_t F2
Value used for radiative corrections.
Double32_t crori
crossing angle orientation, +-1 lepton beam along +-z, +-2 hadron beam along +-z, 0 meas no crossing ...
Double32_t davg
Average density-weighted distance from all inelastic collisions to the edge of the nucleus.
Double32_t trueNu
Generated nu of the event.
virtual double GetHardT() const
Returns the Mandelstamm t of the hard interaction.
Double32_t User3
User variables to prevent/delay future format changes.
virtual double GetTrueY() const
Sets the true (not reconstructed) value for inelasticity.
Double32_t xbeamparton
Momentum fraction taken by the beam parton, see PARI(33)
Double32_t pznucl
target nucleon momentum
Double32_t pztarg
target beam momentum
Double32_t leptonphi
Azimuthal angle of the scattered lepton in the cm frame.
Double32_t x
Bjorken scaling variable.
Double32_t EBrems
Energy per radiative photon in the nuclear rest frame.
Int_t tgtparton
PDG code of the struck parton in the hadron beam, see MSTI(16)
virtual double GetTrueNu() const
Sets the true (not reconstructed) value for ν.
virtual void SetF1(double f1)
Used for radiative corrections.
virtual void SetHardT(double t)
Sets the Mandelstamm t of the hard interaction.
virtual void SetR(double r)
Used for radiative corrections.
virtual double GetHardS() const
Returns the Mandelstamm s of the hard interaction.
Double32_t Phib
phi of impact parameter vector
Int_t Ncolli
Number of inelastic collisions in target.
virtual double GetF1() const
Used for radiative corrections.
Double32_t ThickScl
T(b)/rho0 in fm.
virtual void SetLeptonPhi(double radians)
Azimuthal angle of the scattered lepton in the cm frame.
Int_t Nwound
Number of wounded nucleon including those in INC.
virtual double GetBeamPartonTheta() const
Returns the polar angle of the beam parton in the cm frame, in radians in the range [0,...
Double32_t photonflux
Flux factor, see VINT(319)
virtual double GetPhotonFlux() const
Returns the flux factor, see VINT(319) in PYTHIA 6.
virtual void SetHardPt2(double pt2)
Sets the squared pT of the hard interaction.
virtual bool Parse(const std::string &)
Parses the event information from a text string.
virtual void SetHardU(double u)
Sets the Mandelstamm u of the hard interaction.
virtual double GetTargetPartonX() const
Returns the x of the target parton.
Double32_t trueY
Generated y of the event, see VINT(309)
virtual void SetTargetPartonX(double xB)
Sets the x of the target parton.
Int_t Aremn
A of the nuclear remnant after evaporation and breakup.
Double32_t Eexc
Excitation energy in the nuclear remnant before evaporation and breakup.
EventBeagle(const std::string &str="")
virtual double GetF2() const
Used for radiative corrections.
virtual void SetBeamPartonTheta(double radians)
Sets the polar angle of the beam parton in the cm frame.
Double32_t F1
Value used for radiative corrections.
virtual void SetHardS(double s)
Sets the Mandelstamm s of the hard interaction.
virtual double GetTrueX() const
Sets the true (not reconstructed) value for x.
Double32_t pzf
Sum fermi momentum of all inelastic participant in target rest frame z along gamma*.
virtual void SetTargetParton(int n)
Sets the target parton species.
Double32_t Thickness
T(b) in nucleons/fm^2.
virtual void SetNucleon(int n)
Sets the nucleon species.
bool Parse(const std::string &)
Parses the event information from a text string.
virtual void SetTrueY(double inelasticity)
Sets the true (not reconstructed) value for inelasticity.
virtual void SetPhotonFlux(double f)
Flux factor, see VINT(319) in PYTHIA 6.
virtual double GetHardU() const
Returns the Mandelstamm u of the hard interaction.
Double32_t trueW2
Generated W2 of the event.
Double32_t pt2_hat
Squared pT of the hard subprocess, see PARI(18)
Int_t beamparton
Parton interacting with the hadron beam in the case of resolved photon processes and soft VMD,...
virtual double GetSigmaRad() const
Used for radiative corrections.
Double32_t pzlep
lepton beam momentum
Abstract base class for DIS Monte Carlo events.
virtual void SetTrueW2(double W2)
Sets the true (not reconstructed) value for W2.
virtual double GetTrueW2() const
Sets the true (not reconstructed) value for W2.