Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2024-11-16 09:02:05

0001 /**
0002  \file
0003  Declaration of class erhic::EventRapgap.
0004  
0005  \author    Thomas Burton
0006  \date      2011-08-31
0007  \copyright 2011 Brookhaven National Lab
0008  */
0009 
0010 #ifndef INCLUDE_EICSMEAR_ERHIC_EVENTPYTHIA_H_
0011 #define INCLUDE_EICSMEAR_ERHIC_EVENTPYTHIA_H_
0012 
0013 #include <string>
0014 
0015 #include <Rtypes.h>
0016 
0017 #include "eicsmear/erhic/EventMC.h"
0018 
0019 namespace erhic {
0020 
0021 /**
0022  \brief Pythia 6 DIS event
0023  
0024  Describes an event from the EIC PYTHIA6 implementation.
0025  
0026  \todo Change sHat/t_hat/u_hat naming to be consistent
0027  */
0028 class EventPythia : public EventMC {
0029  public:
0030   /**
0031    Constructor.
0032    @param [in] str A text string setting event-wise quantities. The string
0033    format should be (no newlines):
0034    "I ievent genevent subprocess nucleon targetparton xtargparton
0035    beamparton xbeamparton thetabeamprtn truey trueQ2 truex trueW2 trueNu
0036    leptonphi s_hat t_hat u_hat pt2_hat Q2_hat F2 F1 R sigma_rad SigRadCor
0037    EBrems photonflux nrTracks"
0038    */
0039   explicit EventPythia(const std::string& str = "");
0040 
0041   /**
0042    Destructor
0043    */
0044   virtual ~EventPythia();
0045 
0046   /**
0047    Parses the event information from a text string.
0048    See the constructor for the string format.
0049    Returns true in the event of a successful read operation,
0050    false in case of an error.
0051    */
0052   virtual bool Parse(const std::string&);
0053 
0054   /**
0055    Sets the nucleon species.
0056    @param n [in] PDG code of the hadron beam, see MSTI(12)
0057    */
0058   virtual void SetNucleon(int n);
0059 
0060   /**
0061    Sets the target parton species.
0062    @param n [in] PDG code of the struck parton in the hadron beam,
0063    see MSTI(16)
0064    */
0065   virtual void SetTargetParton(int n);
0066 
0067   /**
0068    Sets the beam parton species.
0069    @param n [in] PDG code of the parton interacting with the hadron beam in
0070    the case of resolved photon processes and soft VMD, see MSTI(15)
0071    */
0072   virtual void SetBeamParton(int n);
0073 
0074   /**
0075    Sets the number of trials required to generate this event
0076    @param n [in] The number of trials
0077    */
0078   virtual void SetGenEvent(int n);
0079 
0080   /**
0081    Sets the x of the target parton.
0082    */
0083   virtual void SetTargetPartonX(double xB);
0084 
0085   /**
0086    Sets the x of the beam parton.
0087    */
0088   virtual void SetBeamPartonX(double xB);
0089 
0090   /**
0091    Sets the polar angle of the beam parton in the cm frame
0092    \todo enforce range [0, pi] when setting
0093    */
0094   virtual void SetBeamPartonTheta(double radians);
0095 
0096   /**
0097    Azimuthal angle of the scattered lepton in the cm frame
0098    \todo enforce range [0, 2pi) when setting
0099    */
0100   virtual void SetLeptonPhi(double radians);
0101 
0102   /**
0103    Used for radiative corrections.
0104    */
0105   virtual void SetF1(double f1);
0106 
0107   /**
0108    Used for radiative corrections.
0109    */
0110   virtual void SetF2(double f2);
0111 
0112   /**
0113    Used for radiative corrections.
0114    */
0115   virtual void SetSigmaRad(double sr);
0116 
0117   /**
0118    Sets the Mandelstamm s of the hard interaction.
0119    */
0120   virtual void SetHardS(double s);
0121 
0122   /**
0123    Sets the Mandelstamm t of the hard interaction.
0124    */
0125   virtual void SetHardT(double t);
0126 
0127   /**
0128    Sets the Mandelstamm u of the hard interaction.
0129    */
0130   virtual void SetHardU(double u);
0131 
0132   /**
0133    Sets the Q<sup>2</sup> of the hard interaction.
0134    */
0135   virtual void SetHardQ2(double Q2);
0136 
0137   /**
0138    Sets the squared p<sub>T</sub> of the hard interaction.
0139    */
0140   virtual void SetHardPt2(double pt2);
0141 
0142   /**
0143    Used for radiative corrections.
0144    */
0145   virtual void SetSigRadCor(double s);
0146 
0147   /**
0148    Sets the energy per radiative photon in the nuclear rest frame.
0149    */
0150   virtual void SetEBrems(double e);
0151 
0152   /**
0153    Flux factor, see VINT(319) in PYTHIA 6.
0154    */
0155   virtual void SetPhotonFlux(double f);
0156 
0157   /**
0158    Sets the true (not reconstructed) value for inelasticity.
0159    */
0160   virtual void SetTrueY(double inelasticity);
0161 
0162   /**
0163    Sets the true (not reconstructed) value for Q<sup>2</sup>.
0164    */
0165   virtual void SetTrueQ2(double Q2);
0166 
0167   /**
0168    Sets the true (not reconstructed) value for x.
0169    */
0170   virtual void SetTrueX(double x);
0171 
0172   /**
0173    Sets the true (not reconstructed) value for W<sup>2</sup>.
0174    */
0175   virtual void SetTrueW2(double W2);
0176 
0177   /**
0178    Sets the true (not reconstructed) value for nu.
0179    */
0180   virtual void SetTrueNu(double Nu);
0181 
0182   /**
0183    Used for radiative corrections.
0184    */
0185   virtual void SetR(double r);
0186 
0187   /**
0188    Returns the number of trials required to generate this event.
0189    */
0190   virtual int GetGenEvent() const;
0191 
0192   /**
0193    Returns the x of the target parton.
0194    */
0195   virtual double GetTargetPartonX() const;
0196 
0197   /**
0198    Returns the x of the beam parton.
0199    */
0200   virtual double GetBeamPartonX() const;
0201 
0202   /**
0203    Returns the polar angle of the beam parton in the cm frame,
0204    in radians in the range [0, pi]
0205    */
0206   virtual double GetBeamPartonTheta() const;
0207 
0208   /**
0209    Returns the azimuthal angle of the scattered lepton.
0210    
0211    Angle is given in the centre-of-mass frame, in radians in the range [0, 2pi).
0212    */
0213   virtual double GetLeptonPhi() const;
0214 
0215   /**
0216    Used for radiative corrections.
0217    */
0218   virtual double GetF1() const;
0219 
0220   /**
0221    Used for radiative corrections.
0222    */
0223   virtual double GetF2() const;
0224 
0225   /**
0226    Used for radiative corrections.
0227    */
0228   virtual double GetSigmaRad() const;
0229 
0230   /**
0231    Returns the Mandelstamm s of the hard interaction.
0232    */
0233   virtual double GetHardS() const;
0234 
0235   /**
0236    Returns the Mandelstamm t of the hard interaction.
0237    */
0238   virtual double GetHardT() const;
0239 
0240   /**
0241    Returns the Mandelstamm u of the hard interaction.
0242    */
0243   virtual double GetHardU() const;
0244 
0245   /**
0246    Returns the Q<sup>2</sup> of the hard interaction.
0247    */
0248   virtual double GetHardQ2() const;
0249 
0250   /**
0251    Returns the squared p<sub>T</sub> of the hard interaction.
0252    */
0253   virtual double GetHardPt2() const;
0254 
0255   /**
0256    Used for radiative corrections.
0257    */
0258   virtual double GetSigRadCor() const;
0259 
0260   /**
0261    Returnss the energy per radiative photon in the nuclear rest frame.
0262    */
0263   virtual double GetEBrems() const;
0264 
0265   /**
0266    Returns the flux factor, see VINT(319) in PYTHIA 6.
0267    */
0268   virtual double GetPhotonFlux() const;
0269 
0270   /**
0271    Sets the true (not reconstructed) value for inelasticity.
0272    */
0273   virtual double GetTrueY() const;
0274 
0275   /**
0276    Sets the true (not reconstructed) value for Q<sup>2</sup>.
0277    */
0278   virtual double GetTrueQ2() const;
0279 
0280   /**
0281    Sets the true (not reconstructed) value for x.
0282    */
0283   virtual double GetTrueX() const;
0284 
0285   /**
0286    Sets the true (not reconstructed) value for W<sup>2</sup>.
0287    */
0288   virtual double GetTrueW2() const;
0289 
0290   /**
0291    Sets the true (not reconstructed) value for &nu;.
0292    */
0293   virtual double GetTrueNu() const;
0294 
0295   /**
0296    Used for radiative corrections.
0297    */
0298   virtual double GetR() const;
0299 
0300   /**
0301    Returns the scattered lepton.
0302    This is the first final state particle with the same species
0303    as the beam lepton and parent index equal to three
0304    (counting index from 1).
0305    */
0306   virtual const ParticleMC* ScatteredLepton() const;
0307 
0308   // Let them all be public; this access method dances for POD does not make sense;
0309   //protected:
0310   // Inline comments after field names will appear in ROOT
0311   // when EventPythia::Dump() is called.
0312   Int_t       nucleon;          ///< PDG code of the hadron beam,
0313                                 ///< see MSTI(12)
0314   Int_t       tgtparton;        ///< PDG code of the struck parton
0315                                 ///< in the hadron beam, see MSTI(16)
0316   Int_t       beamparton;       ///< Parton interacting with the hadron
0317                                 ///< beam in the case of resolved photon
0318                                 ///< processes and soft VMD, see MSTI(15)
0319   Int_t       genevent;         ///< Trials required for this event
0320   Double32_t  xtgtparton;       ///< Momentum fraction taken by the
0321                                 ///< target parton, see PARI(34)
0322   Double32_t  xbeamparton;      ///< Momentum fraction taken by the
0323                                 ///< beam parton, see PARI(33)
0324   Double32_t  thetabeamparton;  ///< Polar angle of the beam parton in
0325                                 ///< the cm frame, between 0 and pi
0326                                 ///< radians, see PARI(53)
0327   Double32_t  leptonphi;        ///< Azimuthal angle of the scattered
0328                                 ///< lepton in the cm frame
0329   Double32_t  F1;               ///< Value used for radiative corrections
0330   Double32_t  sigma_rad;        ///< Value used for radiative corrections
0331   Double32_t  t_hat;            ///< Mandelstam t of the hard subprocess,
0332                                 ///< see PARI(15)
0333   Double32_t  u_hat;            ///< Mandelstam u of the hard subprocess,
0334                                 ///< see PARI(16)
0335   Double32_t  Q2_hat;           ///< Q<sup>2</sup> of the hard subprocess,
0336                                 ///< see PARI(22)
0337   Double32_t  SigRadCor;        ///< Value used for radiative corrections
0338   Double32_t  EBrems;           ///< Energy per radiative photon in the
0339                                 ///< nuclear rest frame.
0340   Double32_t  photonflux;       ///< Flux factor, see VINT(319)
0341   Double32_t  trueY;            ///< Generated y of the event,
0342                                 ///< see VINT(309)
0343   Double32_t  trueQ2;           ///< Generated Q<sup>2</sup> of the event,
0344                                 ///< see VINT(307)
0345   Double32_t  trueX;            ///< Generated x of the event
0346   Double32_t  trueW2;           ///< Generated W<sup>2</sup> of the event
0347   Double32_t  trueNu;           ///< Generated nu of the event
0348   Double32_t  F2;               ///< Value used for radiative corrections
0349   Double32_t  R;                ///< Value used for radiative corrections
0350   Double32_t  pt2_hat;          ///< Squared p<sub>T</sub> of the hard
0351                                 ///< subprocess, see PARI(18)
0352   Double32_t  sHat;             ///< Mandelstam s of the hard subprocess,
0353                                 ///< see PARI(14)
0354   ClassDef(erhic::EventPythia, 2)
0355 };
0356 
0357 inline void EventPythia::SetNucleon(int n) {
0358   nucleon = n;
0359 }
0360 
0361 inline void EventPythia::SetTargetParton(int n) {
0362   tgtparton = n;
0363 }
0364 
0365 inline void EventPythia::SetBeamParton(int n) {
0366   beamparton = n;
0367 }
0368 
0369 inline void EventPythia::SetGenEvent(int n) {
0370   genevent = n;
0371 }
0372 
0373 inline void EventPythia::SetTargetPartonX(double xB) {
0374   xtgtparton = xB;
0375 }
0376 
0377 inline void EventPythia::SetBeamPartonX(double xB) {
0378   xbeamparton = xB;
0379 }
0380 
0381 inline void EventPythia::SetBeamPartonTheta(double radians) {
0382   thetabeamparton = radians;
0383 }
0384 
0385 inline void EventPythia::SetLeptonPhi(double radians) {
0386   leptonphi = radians;
0387 }
0388 
0389 inline void EventPythia::SetF1(double f1) {
0390   F1 = f1;
0391 }
0392 
0393 inline void EventPythia::SetF2(double f2) {
0394   F2 = f2;
0395 }
0396 
0397 inline void EventPythia::SetSigmaRad(double sr) {
0398   sigma_rad = sr;
0399 }
0400 
0401 inline void EventPythia::SetHardS(double s) {
0402   sHat = s;
0403 }
0404 
0405 inline void EventPythia::SetHardT(double t) {
0406   t_hat = t;
0407 }
0408 
0409 inline void EventPythia::SetHardU(double u) {
0410   u_hat = u;
0411 }
0412 
0413 inline void EventPythia::SetHardQ2(double Q2) {
0414   Q2_hat = Q2;
0415 }
0416 
0417 inline void EventPythia::SetHardPt2(double pt2) {
0418   pt2_hat = pt2;
0419 }
0420 
0421 inline void EventPythia::SetSigRadCor(double s) {
0422   SigRadCor = s;
0423 }
0424 
0425 inline void EventPythia::SetEBrems(double e) {
0426   EBrems = e;
0427 }
0428 
0429 inline void EventPythia::SetPhotonFlux(double f) {
0430   photonflux = f;
0431 }
0432 
0433 inline void EventPythia::SetTrueY(double inelasticity) {
0434   trueY = inelasticity;
0435 }
0436 
0437 inline void EventPythia::SetTrueQ2(double Q2) {
0438   trueQ2 = Q2;
0439 }
0440 
0441 inline void EventPythia::SetTrueX(double xB) {
0442   trueX = xB;
0443 }
0444 
0445 inline void EventPythia::SetTrueW2(double W2) {
0446   trueW2 = W2;
0447 }
0448 
0449 inline void EventPythia::SetTrueNu(double Nu) {
0450   trueNu = Nu;
0451 }
0452 
0453 inline void EventPythia::SetR(double r) {
0454   R = r;
0455 }
0456 
0457 inline int EventPythia::GetGenEvent() const {
0458   return genevent;
0459 }
0460 
0461 inline double EventPythia::GetTargetPartonX() const {
0462   return xtgtparton;
0463 }
0464 
0465 inline double EventPythia::GetBeamPartonX() const {
0466   return xbeamparton;
0467 }
0468 
0469 inline double EventPythia::GetBeamPartonTheta() const {
0470   return thetabeamparton;
0471 }
0472 
0473 inline double EventPythia::GetLeptonPhi() const {
0474   return leptonphi;
0475 }
0476 
0477 inline double EventPythia::GetF1() const {
0478   return F1;
0479 }
0480 
0481 inline double EventPythia::GetF2() const {
0482   return F2;
0483 }
0484 
0485 inline double EventPythia::GetSigmaRad() const {
0486   return sigma_rad;
0487 }
0488 
0489 inline double EventPythia::GetHardS() const {
0490   return sHat;
0491 }
0492 
0493 inline double EventPythia::GetHardT() const {
0494   return t_hat;
0495 }
0496 
0497 inline double EventPythia::GetHardU() const {
0498   return u_hat;
0499 }
0500 
0501 inline double EventPythia::GetHardQ2() const {
0502   return Q2_hat;
0503 }
0504 
0505 inline double EventPythia::GetHardPt2() const {
0506   return pt2_hat;
0507 }
0508 
0509 inline double EventPythia::GetSigRadCor() const {
0510   return SigRadCor;
0511 }
0512 
0513 inline double EventPythia::GetEBrems() const {
0514   return EBrems;
0515 }
0516 
0517 inline double EventPythia::GetPhotonFlux() const {
0518   return photonflux;
0519 }
0520 
0521 inline double EventPythia::GetTrueY() const {
0522   return trueY;
0523 }
0524 
0525 inline double EventPythia::GetTrueQ2() const {
0526   return trueQ2;
0527 }
0528 
0529 inline double EventPythia::GetTrueX() const {
0530   return trueX;
0531 }
0532 
0533 inline double EventPythia::GetTrueW2() const {
0534   return trueW2;
0535 }
0536 
0537 inline double EventPythia::GetTrueNu() const {
0538   return trueNu;
0539 }
0540 
0541 inline double EventPythia::GetR() const {
0542   return R;
0543 }
0544 
0545 
0546 class EventBeagle : public EventPythia {
0547  public:
0548   explicit EventBeagle(const std::string& str = "");
0549 
0550   ~EventBeagle();
0551 
0552   bool RequiresEaParticleFields() { return true; };
0553 
0554   bool Parse(const std::string&);
0555 
0556   ///=================additional variables for BeAGLE==================
0557     //put in the public region to be easily accessed
0558   Int_t lepton;  ///< lepton beam ID
0559   Int_t Atarg;   ///< mass number of target beam
0560   Int_t Ztarg;   ///< charge number of target beam
0561   Double32_t pzlep;     ///< lepton beam momentum
0562   Double32_t pztarg;        ///< target beam momentum
0563   Double32_t pznucl;        ///< target nucleon momentum
0564     Double32_t crang;           ///< crossing angle (mr). crang=1000*atan(px/pz), one beam px=py=0, the other py=0
0565     Double32_t crori;           ///< crossing angle orientation, +-1 lepton beam along +-z, +-2 hadron beam along +-z, 0 meas no crossing angle
0566   Double32_t b;     ///< impact parameter value
0567   Double32_t Phib;      ///< phi of impact parameter vector
0568     Double32_t Thickness; ///< T(b) in nucleons/fm^2
0569     Double32_t ThickScl; ///< T(b)/rho0 in fm
0570     Int_t Ncollt; ///< Number of collisions in target
0571     Int_t Ncolli; ///< Number of inelastic collisions in target
0572     Int_t Nwound; ///< Number of wounded nucleon including those in INC
0573     Int_t Nwdch; ///< Number of wounded proton including those in INC
0574     Int_t Nnevap; ///< Number of neutrons from evaporation
0575     Int_t Npevap; ///< Number of protons from evaporation
0576     Int_t Aremn; ///< A of the nuclear remnant after evaporation and breakup
0577   Int_t NINC; ///< Number of stable hadrons from intranuclear cascade
0578     Int_t NINCch; ///< Number of charged stable hadrons from intranuclear cascade
0579     Double32_t d1st; ///< density-weighted distance from first collision to the edge of the nucleus (amount of material traversed / rho0)
0580     Double32_t davg; ///< Average density-weighted distance from all inelastic collisions to the edge of the nucleus
0581     Double32_t pxf; ///< Sum fermi momentum of all inelastic participant in target rest frame z along gamma*
0582     Double32_t pyf; ///< Sum fermi momentum of all inelastic participant in target rest frame z along gamma*
0583     Double32_t pzf; ///< Sum fermi momentum of all inelastic participant in target rest frame z along gamma*
0584     Double32_t Eexc; ///< Excitation energy in the nuclear remnant before evaporation and breakup 
0585     Double32_t RAevt; ///< Nuclear PDF ratio for the up sea for the given event kinematics 
0586     Double32_t User1; ///< User variables to prevent/delay future format changes 
0587     Double32_t User2; ///< User variables to prevent/delay future format changes 
0588     Double32_t User3; ///< User variables to prevent/delay future format changes 
0589 
0590   ClassDef(erhic::EventBeagle, 1)
0591 };
0592 
0593 }  // namespace erhic
0594 
0595 #endif  // INCLUDE_EICSMEAR_ERHIC_EVENTPYTHIA_H_