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.h
Go to the documentation of this file.
1 
10 #ifndef INCLUDE_EICSMEAR_ERHIC_EVENTPYTHIA_H_
11 #define INCLUDE_EICSMEAR_ERHIC_EVENTPYTHIA_H_
12 
13 #include <string>
14 
15 #include <Rtypes.h>
16 
17 #include "eicsmear/erhic/EventMC.h"
18 
19 namespace erhic {
20 
28 class EventPythia : public EventMC {
29  public:
39  explicit EventPythia(const std::string& str = "");
40 
44  virtual ~EventPythia();
45 
52  virtual bool Parse(const std::string&);
53 
58  virtual void SetNucleon(int n);
59 
65  virtual void SetTargetParton(int n);
66 
72  virtual void SetBeamParton(int n);
73 
78  virtual void SetGenEvent(int n);
79 
83  virtual void SetTargetPartonX(double xB);
84 
88  virtual void SetBeamPartonX(double xB);
89 
94  virtual void SetBeamPartonTheta(double radians);
95 
100  virtual void SetLeptonPhi(double radians);
101 
105  virtual void SetF1(double f1);
106 
110  virtual void SetF2(double f2);
111 
115  virtual void SetSigmaRad(double sr);
116 
120  virtual void SetHardS(double s);
121 
125  virtual void SetHardT(double t);
126 
130  virtual void SetHardU(double u);
131 
135  virtual void SetHardQ2(double Q2);
136 
140  virtual void SetHardPt2(double pt2);
141 
145  virtual void SetSigRadCor(double s);
146 
150  virtual void SetEBrems(double e);
151 
155  virtual void SetPhotonFlux(double f);
156 
160  virtual void SetTrueY(double inelasticity);
161 
165  virtual void SetTrueQ2(double Q2);
166 
170  virtual void SetTrueX(double x);
171 
175  virtual void SetTrueW2(double W2);
176 
180  virtual void SetTrueNu(double Nu);
181 
185  virtual void SetR(double r);
186 
190  virtual int GetGenEvent() const;
191 
195  virtual double GetTargetPartonX() const;
196 
200  virtual double GetBeamPartonX() const;
201 
206  virtual double GetBeamPartonTheta() const;
207 
213  virtual double GetLeptonPhi() const;
214 
218  virtual double GetF1() const;
219 
223  virtual double GetF2() const;
224 
228  virtual double GetSigmaRad() const;
229 
233  virtual double GetHardS() const;
234 
238  virtual double GetHardT() const;
239 
243  virtual double GetHardU() const;
244 
248  virtual double GetHardQ2() const;
249 
253  virtual double GetHardPt2() const;
254 
258  virtual double GetSigRadCor() const;
259 
263  virtual double GetEBrems() const;
264 
268  virtual double GetPhotonFlux() const;
269 
273  virtual double GetTrueY() const;
274 
278  virtual double GetTrueQ2() const;
279 
283  virtual double GetTrueX() const;
284 
288  virtual double GetTrueW2() const;
289 
293  virtual double GetTrueNu() const;
294 
298  virtual double GetR() const;
299 
306  virtual const ParticleMC* ScatteredLepton() const;
307 
308  // Let them all be public; this access method dances for POD does not make sense;
309  //protected:
310  // Inline comments after field names will appear in ROOT
311  // when EventPythia::Dump() is called.
312  Int_t nucleon;
313  Int_t tgtparton;
315  Int_t beamparton;
317  Int_t genevent;
320  Double32_t xtgtparton;
321  Double32_t xbeamparton;
323  Double32_t thetabeamparton;
325  Double32_t leptonphi;
328  Double32_t F1;
330  Double32_t sigma_rad;
331  Double32_t t_hat;
332  Double32_t u_hat;
334  Double32_t Q2_hat;
336  Double32_t SigRadCor;
338  Double32_t EBrems;
339  Double32_t photonflux;
341  Double32_t trueY;
342  Double32_t trueQ2;
344  Double32_t trueX;
346  Double32_t trueW2;
347  Double32_t trueNu;
348  Double32_t F2;
349  Double32_t R;
350  Double32_t pt2_hat;
351  Double32_t sHat;
353  ClassDef(erhic::EventPythia, 2)
355 };
356 
357 inline void EventPythia::SetNucleon(int n) {
358  nucleon = n;
359 }
360 
361 inline void EventPythia::SetTargetParton(int n) {
362  tgtparton = n;
363 }
364 
365 inline void EventPythia::SetBeamParton(int n) {
366  beamparton = n;
367 }
368 
369 inline void EventPythia::SetGenEvent(int n) {
370  genevent = n;
371 }
372 
373 inline void EventPythia::SetTargetPartonX(double xB) {
374  xtgtparton = xB;
375 }
376 
377 inline void EventPythia::SetBeamPartonX(double xB) {
378  xbeamparton = xB;
379 }
380 
381 inline void EventPythia::SetBeamPartonTheta(double radians) {
382  thetabeamparton = radians;
383 }
384 
385 inline void EventPythia::SetLeptonPhi(double radians) {
386  leptonphi = radians;
387 }
388 
389 inline void EventPythia::SetF1(double f1) {
390  F1 = f1;
391 }
392 
393 inline void EventPythia::SetF2(double f2) {
394  F2 = f2;
395 }
396 
397 inline void EventPythia::SetSigmaRad(double sr) {
398  sigma_rad = sr;
399 }
400 
401 inline void EventPythia::SetHardS(double s) {
402  sHat = s;
403 }
404 
405 inline void EventPythia::SetHardT(double t) {
406  t_hat = t;
407 }
408 
409 inline void EventPythia::SetHardU(double u) {
410  u_hat = u;
411 }
412 
413 inline void EventPythia::SetHardQ2(double Q2) {
414  Q2_hat = Q2;
415 }
416 
417 inline void EventPythia::SetHardPt2(double pt2) {
418  pt2_hat = pt2;
419 }
420 
421 inline void EventPythia::SetSigRadCor(double s) {
422  SigRadCor = s;
423 }
424 
425 inline void EventPythia::SetEBrems(double e) {
426  EBrems = e;
427 }
428 
429 inline void EventPythia::SetPhotonFlux(double f) {
430  photonflux = f;
431 }
432 
433 inline void EventPythia::SetTrueY(double inelasticity) {
434  trueY = inelasticity;
435 }
436 
437 inline void EventPythia::SetTrueQ2(double Q2) {
438  trueQ2 = Q2;
439 }
440 
441 inline void EventPythia::SetTrueX(double xB) {
442  trueX = xB;
443 }
444 
445 inline void EventPythia::SetTrueW2(double W2) {
446  trueW2 = W2;
447 }
448 
449 inline void EventPythia::SetTrueNu(double Nu) {
450  trueNu = Nu;
451 }
452 
453 inline void EventPythia::SetR(double r) {
454  R = r;
455 }
456 
457 inline int EventPythia::GetGenEvent() const {
458  return genevent;
459 }
460 
461 inline double EventPythia::GetTargetPartonX() const {
462  return xtgtparton;
463 }
464 
465 inline double EventPythia::GetBeamPartonX() const {
466  return xbeamparton;
467 }
468 
469 inline double EventPythia::GetBeamPartonTheta() const {
470  return thetabeamparton;
471 }
472 
473 inline double EventPythia::GetLeptonPhi() const {
474  return leptonphi;
475 }
476 
477 inline double EventPythia::GetF1() const {
478  return F1;
479 }
480 
481 inline double EventPythia::GetF2() const {
482  return F2;
483 }
484 
485 inline double EventPythia::GetSigmaRad() const {
486  return sigma_rad;
487 }
488 
489 inline double EventPythia::GetHardS() const {
490  return sHat;
491 }
492 
493 inline double EventPythia::GetHardT() const {
494  return t_hat;
495 }
496 
497 inline double EventPythia::GetHardU() const {
498  return u_hat;
499 }
500 
501 inline double EventPythia::GetHardQ2() const {
502  return Q2_hat;
503 }
504 
505 inline double EventPythia::GetHardPt2() const {
506  return pt2_hat;
507 }
508 
509 inline double EventPythia::GetSigRadCor() const {
510  return SigRadCor;
511 }
512 
513 inline double EventPythia::GetEBrems() const {
514  return EBrems;
515 }
516 
517 inline double EventPythia::GetPhotonFlux() const {
518  return photonflux;
519 }
520 
521 inline double EventPythia::GetTrueY() const {
522  return trueY;
523 }
524 
525 inline double EventPythia::GetTrueQ2() const {
526  return trueQ2;
527 }
528 
529 inline double EventPythia::GetTrueX() const {
530  return trueX;
531 }
532 
533 inline double EventPythia::GetTrueW2() const {
534  return trueW2;
535 }
536 
537 inline double EventPythia::GetTrueNu() const {
538  return trueNu;
539 }
540 
541 inline double EventPythia::GetR() const {
542  return R;
543 }
544 
545 
546 class EventBeagle : public EventPythia {
547  public:
548  explicit EventBeagle(const std::string& str = "");
549 
550  ~EventBeagle();
551 
552  bool RequiresEaParticleFields() { return true; };
553 
554  bool Parse(const std::string&);
555 
557  //put in the public region to be easily accessed
558  Int_t lepton;
559  Int_t Atarg;
560  Int_t Ztarg;
561  Double32_t pzlep;
562  Double32_t pztarg;
563  Double32_t pznucl;
564  Double32_t crang;
565  Double32_t crori;
566  Double32_t b;
567  Double32_t Phib;
568  Double32_t Thickness;
569  Double32_t ThickScl;
570  Int_t Ncollt;
571  Int_t Ncolli;
572  Int_t Nwound;
573  Int_t Nwdch;
574  Int_t Nnevap;
575  Int_t Npevap;
576  Int_t Aremn;
577  Int_t NINC;
578  Int_t NINCch;
579  Double32_t d1st;
580  Double32_t davg;
581  Double32_t pxf;
582  Double32_t pyf;
583  Double32_t pzf;
584  Double32_t Eexc;
585  Double32_t RAevt;
586  Double32_t User1;
587  Double32_t User2;
588  Double32_t User3;
589 
590  ClassDef(erhic::EventBeagle, 1)
591 };
592 
593 } // namespace erhic
594 
595 #endif // INCLUDE_EICSMEAR_ERHIC_EVENTPYTHIA_H_
erhic::EventPythia::t_hat
Double32_t t_hat
Mandelstam t of the hard subprocess, see PARI(15)
Definition: erhic/EventPythia.h:331
erhic::EventPythia::xtgtparton
Double32_t xtgtparton
Momentum fraction taken by the target parton, see PARI(34)
Definition: erhic/EventPythia.h:320
erhic::EventBeagle::Ncollt
Int_t Ncollt
Number of collisions in target.
Definition: erhic/EventPythia.h:570
erhic::EventPythia::SetBeamPartonX
virtual void SetBeamPartonX(double xB)
Sets the x of the beam parton.
Definition: erhic/EventPythia.h:377
erhic::EventBeagle::RAevt
Double32_t RAevt
Nuclear PDF ratio for the up sea for the given event kinematics.
Definition: erhic/EventPythia.h:585
erhic::EventPythia::sigma_rad
Double32_t sigma_rad
Value used for radiative corrections.
Definition: erhic/EventPythia.h:330
erhic::EventBeagle::lepton
Int_t lepton
=================additional variables for BeAGLE==================
Definition: erhic/EventPythia.h:558
erhic::EventPythia::SetBeamParton
virtual void SetBeamParton(int n)
Sets the beam parton species.
Definition: erhic/EventPythia.h:365
erhic::EventPythia::EventPythia
EventPythia(const std::string &str="")
Constructor.
Definition: erhic/EventPythia.cxx:20
erhic
Definition: EventDis.cxx:14
erhic::EventPythia::SetTrueX
virtual void SetTrueX(double x)
Sets the true (not reconstructed) value for x.
Definition: erhic/EventPythia.h:441
erhic::EventBeagle::Npevap
Int_t Npevap
Number of protons from evaporation.
Definition: erhic/EventPythia.h:575
erhic::EventPythia::SetSigRadCor
virtual void SetSigRadCor(double s)
Used for radiative corrections.
Definition: erhic/EventPythia.h:421
erhic::EventPythia::sHat
Double32_t sHat
Mandelstam s of the hard subprocess, see PARI(14)
Definition: erhic/EventPythia.h:352
erhic::EventBeagle::pyf
Double32_t pyf
Sum fermi momentum of all inelastic participant in target rest frame z along gamma*.
Definition: erhic/EventPythia.h:582
erhic::EventPythia::nucleon
Int_t nucleon
PDG code of the hadron beam, see MSTI(12)
Definition: erhic/EventPythia.h:312
erhic::EventPythia::u_hat
Double32_t u_hat
Mandelstam u of the hard subprocess, see PARI(16)
Definition: erhic/EventPythia.h:333
erhic::EventBeagle::Ztarg
Int_t Ztarg
charge number of target beam
Definition: erhic/EventPythia.h:560
erhic::EventPythia::GetBeamPartonX
virtual double GetBeamPartonX() const
Returns the x of the beam parton.
Definition: erhic/EventPythia.h:465
erhic::EventPythia::SetEBrems
virtual void SetEBrems(double e)
Sets the energy per radiative photon in the nuclear rest frame.
Definition: erhic/EventPythia.h:425
erhic::EventPythia::Q2_hat
Double32_t Q2_hat
Q2 of the hard subprocess, see PARI(22)
Definition: erhic/EventPythia.h:335
erhic::EventPythia::GetEBrems
virtual double GetEBrems() const
Returnss the energy per radiative photon in the nuclear rest frame.
Definition: erhic/EventPythia.h:513
erhic::EventPythia::GetLeptonPhi
virtual double GetLeptonPhi() const
Returns the azimuthal angle of the scattered lepton.
Definition: erhic/EventPythia.h:473
erhic::EventPythia::GetR
virtual double GetR() const
Used for radiative corrections.
Definition: erhic/EventPythia.h:541
erhic::EventBeagle::crang
Double32_t crang
crossing angle (mr). crang=1000*atan(px/pz), one beam px=py=0, the other py=0
Definition: erhic/EventPythia.h:564
erhic::EventPythia::GetHardQ2
virtual double GetHardQ2() const
Returns the Q2 of the hard interaction.
Definition: erhic/EventPythia.h:501
erhic::EventBeagle::NINC
Int_t NINC
Number of stable hadrons from intranuclear cascade.
Definition: erhic/EventPythia.h:577
erhic::EventPythia::GetSigRadCor
virtual double GetSigRadCor() const
Used for radiative corrections.
Definition: erhic/EventPythia.h:509
erhic::EventPythia::SetTrueNu
virtual void SetTrueNu(double Nu)
Sets the true (not reconstructed) value for nu.
Definition: erhic/EventPythia.h:449
erhic::EventPythia::~EventPythia
virtual ~EventPythia()
Destructor.
Definition: erhic/EventPythia.cxx:48
EventMC.h
erhic::EventPythia::GetTrueQ2
virtual double GetTrueQ2() const
Sets the true (not reconstructed) value for Q2.
Definition: erhic/EventPythia.h:525
erhic::EventBeagle::Atarg
Int_t Atarg
mass number of target beam
Definition: erhic/EventPythia.h:559
erhic::EventBeagle::d1st
Double32_t d1st
density-weighted distance from first collision to the edge of the nucleus (amount of material travers...
Definition: erhic/EventPythia.h:579
erhic::EventPythia::trueQ2
Double32_t trueQ2
Generated Q2 of the event, see VINT(307)
Definition: erhic/EventPythia.h:343
erhic::EventBeagle::User1
Double32_t User1
User variables to prevent/delay future format changes.
Definition: erhic/EventPythia.h:586
erhic::EventBeagle::~EventBeagle
~EventBeagle()
Definition: erhic/EventPythia.cxx:135
erhic::EventBeagle::NINCch
Int_t NINCch
Number of charged stable hadrons from intranuclear cascade.
Definition: erhic/EventPythia.h:578
erhic::EventPythia::R
Double32_t R
Value used for radiative corrections.
Definition: erhic/EventPythia.h:349
erhic::EventPythia::thetabeamparton
Double32_t thetabeamparton
Polar angle of the beam parton in the cm frame, between 0 and pi radians, see PARI(53)
Definition: erhic/EventPythia.h:324
erhic::EventPythia::ScatteredLepton
virtual const ParticleMC * ScatteredLepton() const
Returns the scattered lepton.
Definition: erhic/EventPythia.cxx:71
erhic::EventPythia::GetHardPt2
virtual double GetHardPt2() const
Returns the squared pT of the hard interaction.
Definition: erhic/EventPythia.h:505
erhic::EventPythia::trueX
Double32_t trueX
Generated x of the event.
Definition: erhic/EventPythia.h:345
erhic::EventPythia::SetGenEvent
virtual void SetGenEvent(int n)
Sets the number of trials required to generate this event.
Definition: erhic/EventPythia.h:369
erhic::EventPythia::SetTrueQ2
virtual void SetTrueQ2(double Q2)
Sets the true (not reconstructed) value for Q2.
Definition: erhic/EventPythia.h:437
erhic::EventBeagle::RequiresEaParticleFields
bool RequiresEaParticleFields()
Definition: erhic/EventPythia.h:552
erhic::EventPythia
Pythia 6 DIS event.
Definition: erhic/EventPythia.h:28
erhic::EventPythia::SigRadCor
Double32_t SigRadCor
Value used for radiative corrections.
Definition: erhic/EventPythia.h:337
erhic::EventBeagle::Nwdch
Int_t Nwdch
Number of wounded proton including those in INC.
Definition: erhic/EventPythia.h:573
erhic::EventBeagle::User2
Double32_t User2
User variables to prevent/delay future format changes.
Definition: erhic/EventPythia.h:587
erhic::EventBeagle::b
Double32_t b
impact parameter value
Definition: erhic/EventPythia.h:566
erhic::EventBeagle::pxf
Double32_t pxf
Sum fermi momentum of all inelastic participant in target rest frame z along gamma*.
Definition: erhic/EventPythia.h:581
erhic::EventPythia::genevent
Int_t genevent
Trials required for this event.
Definition: erhic/EventPythia.h:319
erhic::EventBeagle::Nnevap
Int_t Nnevap
Number of neutrons from evaporation.
Definition: erhic/EventPythia.h:574
erhic::EventPythia::SetHardQ2
virtual void SetHardQ2(double Q2)
Sets the Q2 of the hard interaction.
Definition: erhic/EventPythia.h:413
erhic::EventPythia::GetGenEvent
virtual int GetGenEvent() const
Returns the number of trials required to generate this event.
Definition: erhic/EventPythia.h:457
erhic::EventPythia::SetSigmaRad
virtual void SetSigmaRad(double sr)
Used for radiative corrections.
Definition: erhic/EventPythia.h:397
erhic::EventPythia::SetF2
virtual void SetF2(double f2)
Used for radiative corrections.
Definition: erhic/EventPythia.h:393
erhic::EventPythia::F2
Double32_t F2
Value used for radiative corrections.
Definition: erhic/EventPythia.h:348
erhic::EventBeagle::crori
Double32_t crori
crossing angle orientation, +-1 lepton beam along +-z, +-2 hadron beam along +-z, 0 meas no crossing ...
Definition: erhic/EventPythia.h:565
erhic::EventBeagle::davg
Double32_t davg
Average density-weighted distance from all inelastic collisions to the edge of the nucleus.
Definition: erhic/EventPythia.h:580
erhic::EventPythia::trueNu
Double32_t trueNu
Generated nu of the event.
Definition: erhic/EventPythia.h:347
erhic::EventPythia::GetHardT
virtual double GetHardT() const
Returns the Mandelstamm t of the hard interaction.
Definition: erhic/EventPythia.h:493
erhic::EventBeagle::User3
Double32_t User3
User variables to prevent/delay future format changes.
Definition: erhic/EventPythia.h:588
erhic::EventPythia::GetTrueY
virtual double GetTrueY() const
Sets the true (not reconstructed) value for inelasticity.
Definition: erhic/EventPythia.h:521
erhic::EventPythia::xbeamparton
Double32_t xbeamparton
Momentum fraction taken by the beam parton, see PARI(33)
Definition: erhic/EventPythia.h:322
erhic::EventBeagle::pznucl
Double32_t pznucl
target nucleon momentum
Definition: erhic/EventPythia.h:563
erhic::EventBeagle::pztarg
Double32_t pztarg
target beam momentum
Definition: erhic/EventPythia.h:562
erhic::EventPythia::leptonphi
Double32_t leptonphi
Azimuthal angle of the scattered lepton in the cm frame.
Definition: erhic/EventPythia.h:327
erhic::EventDis::x
Double32_t x
Bjorken scaling variable.
Definition: EventDis.h:181
erhic::EventPythia::EBrems
Double32_t EBrems
Energy per radiative photon in the nuclear rest frame.
Definition: erhic/EventPythia.h:338
erhic::EventPythia::tgtparton
Int_t tgtparton
PDG code of the struck parton in the hadron beam, see MSTI(16)
Definition: erhic/EventPythia.h:314
erhic::EventPythia::GetTrueNu
virtual double GetTrueNu() const
Sets the true (not reconstructed) value for ν.
Definition: erhic/EventPythia.h:537
erhic::EventPythia::SetF1
virtual void SetF1(double f1)
Used for radiative corrections.
Definition: erhic/EventPythia.h:389
erhic::ParticleMC
Definition: erhic/ParticleMC.h:403
erhic::EventPythia::SetHardT
virtual void SetHardT(double t)
Sets the Mandelstamm t of the hard interaction.
Definition: erhic/EventPythia.h:405
erhic::EventPythia::SetR
virtual void SetR(double r)
Used for radiative corrections.
Definition: erhic/EventPythia.h:453
erhic::EventPythia::GetHardS
virtual double GetHardS() const
Returns the Mandelstamm s of the hard interaction.
Definition: erhic/EventPythia.h:489
erhic::EventBeagle::Phib
Double32_t Phib
phi of impact parameter vector
Definition: erhic/EventPythia.h:567
erhic::EventBeagle::Ncolli
Int_t Ncolli
Number of inelastic collisions in target.
Definition: erhic/EventPythia.h:571
erhic::EventPythia::GetF1
virtual double GetF1() const
Used for radiative corrections.
Definition: erhic/EventPythia.h:477
erhic::EventBeagle::ThickScl
Double32_t ThickScl
T(b)/rho0 in fm.
Definition: erhic/EventPythia.h:569
erhic::EventPythia::SetLeptonPhi
virtual void SetLeptonPhi(double radians)
Azimuthal angle of the scattered lepton in the cm frame.
Definition: erhic/EventPythia.h:385
erhic::EventBeagle::Nwound
Int_t Nwound
Number of wounded nucleon including those in INC.
Definition: erhic/EventPythia.h:572
erhic::EventPythia::GetBeamPartonTheta
virtual double GetBeamPartonTheta() const
Returns the polar angle of the beam parton in the cm frame, in radians in the range [0,...
Definition: erhic/EventPythia.h:469
erhic::EventPythia::photonflux
Double32_t photonflux
Flux factor, see VINT(319)
Definition: erhic/EventPythia.h:340
erhic::EventPythia::GetPhotonFlux
virtual double GetPhotonFlux() const
Returns the flux factor, see VINT(319) in PYTHIA 6.
Definition: erhic/EventPythia.h:517
erhic::EventPythia::SetHardPt2
virtual void SetHardPt2(double pt2)
Sets the squared pT of the hard interaction.
Definition: erhic/EventPythia.h:417
erhic::EventPythia::Parse
virtual bool Parse(const std::string &)
Parses the event information from a text string.
Definition: erhic/EventPythia.cxx:50
erhic::EventPythia::SetHardU
virtual void SetHardU(double u)
Sets the Mandelstamm u of the hard interaction.
Definition: erhic/EventPythia.h:409
erhic::EventPythia::GetTargetPartonX
virtual double GetTargetPartonX() const
Returns the x of the target parton.
Definition: erhic/EventPythia.h:461
erhic::EventPythia::trueY
Double32_t trueY
Generated y of the event, see VINT(309)
Definition: erhic/EventPythia.h:341
erhic::EventPythia::SetTargetPartonX
virtual void SetTargetPartonX(double xB)
Sets the x of the target parton.
Definition: erhic/EventPythia.h:373
erhic::EventBeagle::Aremn
Int_t Aremn
A of the nuclear remnant after evaporation and breakup.
Definition: erhic/EventPythia.h:576
erhic::EventBeagle::Eexc
Double32_t Eexc
Excitation energy in the nuclear remnant before evaporation and breakup.
Definition: erhic/EventPythia.h:584
erhic::EventBeagle::EventBeagle
EventBeagle(const std::string &str="")
Definition: erhic/EventPythia.cxx:98
erhic::EventPythia::GetF2
virtual double GetF2() const
Used for radiative corrections.
Definition: erhic/EventPythia.h:481
erhic::EventPythia::SetBeamPartonTheta
virtual void SetBeamPartonTheta(double radians)
Sets the polar angle of the beam parton in the cm frame.
Definition: erhic/EventPythia.h:381
erhic::EventPythia::F1
Double32_t F1
Value used for radiative corrections.
Definition: erhic/EventPythia.h:329
erhic::EventPythia::SetHardS
virtual void SetHardS(double s)
Sets the Mandelstamm s of the hard interaction.
Definition: erhic/EventPythia.h:401
erhic::EventPythia::GetTrueX
virtual double GetTrueX() const
Sets the true (not reconstructed) value for x.
Definition: erhic/EventPythia.h:529
erhic::EventBeagle::pzf
Double32_t pzf
Sum fermi momentum of all inelastic participant in target rest frame z along gamma*.
Definition: erhic/EventPythia.h:583
erhic::EventPythia::SetTargetParton
virtual void SetTargetParton(int n)
Sets the target parton species.
Definition: erhic/EventPythia.h:361
erhic::EventBeagle::Thickness
Double32_t Thickness
T(b) in nucleons/fm^2.
Definition: erhic/EventPythia.h:568
erhic::EventPythia::SetNucleon
virtual void SetNucleon(int n)
Sets the nucleon species.
Definition: erhic/EventPythia.h:357
erhic::EventBeagle::Parse
bool Parse(const std::string &)
Parses the event information from a text string.
Definition: erhic/EventPythia.cxx:137
erhic::EventPythia::SetTrueY
virtual void SetTrueY(double inelasticity)
Sets the true (not reconstructed) value for inelasticity.
Definition: erhic/EventPythia.h:433
erhic::EventPythia::SetPhotonFlux
virtual void SetPhotonFlux(double f)
Flux factor, see VINT(319) in PYTHIA 6.
Definition: erhic/EventPythia.h:429
erhic::EventPythia::GetHardU
virtual double GetHardU() const
Returns the Mandelstamm u of the hard interaction.
Definition: erhic/EventPythia.h:497
erhic::EventPythia::trueW2
Double32_t trueW2
Generated W2 of the event.
Definition: erhic/EventPythia.h:346
erhic::EventPythia::pt2_hat
Double32_t pt2_hat
Squared pT of the hard subprocess, see PARI(18)
Definition: erhic/EventPythia.h:350
erhic::EventPythia::beamparton
Int_t beamparton
Parton interacting with the hadron beam in the case of resolved photon processes and soft VMD,...
Definition: erhic/EventPythia.h:316
erhic::EventPythia::GetSigmaRad
virtual double GetSigmaRad() const
Used for radiative corrections.
Definition: erhic/EventPythia.h:485
erhic::EventBeagle::pzlep
Double32_t pzlep
lepton beam momentum
Definition: erhic/EventPythia.h:561
erhic::EventMC
Abstract base class for DIS Monte Carlo events.
Definition: erhic/EventMC.h:30
erhic::EventPythia::SetTrueW2
virtual void SetTrueW2(double W2)
Sets the true (not reconstructed) value for W2.
Definition: erhic/EventPythia.h:445
erhic::EventPythia::GetTrueW2
virtual double GetTrueW2() const
Sets the true (not reconstructed) value for W2.
Definition: erhic/EventPythia.h:533
erhic::EventBeagle
Definition: erhic/EventPythia.h:546