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
EventGmcTrans.h
Go to the documentation of this file.
1 
10 #ifndef INCLUDE_EICSMEAR_ERHIC_EVENTGMCTRANS_H_
11 #define INCLUDE_EICSMEAR_ERHIC_EVENTGMCTRANS_H_
12 
13 #include <string>
14 
15 #include <TLorentzVector.h>
16 
17 #include "eicsmear/erhic/EventMC.h"
18 
19 namespace erhic {
23 class EventGmcTrans : public EventMC {
24  public:
30  explicit EventGmcTrans(const std::string& = "");
31 
35  virtual ~EventGmcTrans() { }
36 
49  virtual bool Parse(const std::string&);
50 
54  Double_t GetPhiSpin() const { return mPhiSpin; }
55 
59  Double_t GetPhiHadron() const { return mPhiHadron; }
60 
64  Double_t GetHadronZ() const { return mZ; }
65 
66  Double_t GetHadronPt() const { return mHadronPt; }
67 
68  Double_t GetF1() const { return mF1; }
69 
70  Double_t GetG1() const { return mG1; }
71 
72  Double_t GetH1() const { return mH1; }
73 
74  Double_t GetD1() const { return mD1; }
75 
76  Double_t GetF1TPerp() const { return mF1TPerp; }
77 
78  Double_t GetF1TPerp1() const { return mF1TPerp1; }
79 
80  Double_t GetH1Perp() const { return mH1Perp; }
81 
82  Double_t GetH1Perp1() const { return mH1Perp1; }
83 
84  Double_t GetH1Perp12() const { return mH1Perp12; }
85 
86  Double_t GetSivers() const { return mAutSivAllQ; }
87 
88  Double_t GetSiversWeight() const { return mAutWtSivAllQ; }
89 
90  Double_t GetSiversStruckQuark() const { return mAutSiv; }
91 
92  Double_t GetSiversStruckQuarkWeight() const { return mAutWtSiv; }
93 
94  Double_t GetSiversPiDifference() const { return mAutSivPiDiff; }
95 
96  Double_t GetSiversPiDifferenceWeight() const { return mAutWtSivPiDiff; }
97 
98  Double_t GetCollins() const { return mAutColAllQ; }
99 
100  Double_t GetCollinsWeight() const { return mAutWtColAllQ; }
101 
102  Double_t GetCollinsStruckQuark() const { return mAutCol; }
103 
104  Double_t GetCollinsStructQuarkWeight() const { return mAutWtCol; }
105 
106  Double_t GetCollinsTwist3() const { return mAutTw3Col; }
107 
108  Double_t GetCollinsTwist3Weight() const { return mAutWtTw3Col; }
109 
110  Double_t GetCrossSectionUnpolarised() const { return mXUnpolarised; }
111 
112  Double_t GetCrossSectionSivers() const { return mXSivers; }
113 
114  Double_t GetCrossSectionCollins() const { return mXCollins; }
115 
131  TLorentzVector GetHadronPolarisation() const { return TLorentzVector(); }
132 
140  virtual double GetHermesPhiS() const { return 0.; }
141 
142  protected:
143  Int_t mStruckQuark;
144  Double32_t mQSquared;
145  Double32_t mBjorkenX;
147  Double32_t mInelasticity; // or mY?
149  Double32_t mWSquared;
150  Double32_t mNu;
151  Double32_t mS;
152  Double32_t mZ;
153  Double32_t mHadronPt;
154  Double32_t mLeptonTheta;
155  Double32_t mLeptonPhi;
156  Double32_t mPhiSpin;
157  Double32_t mPhiHadron;
159  Double32_t mF1;
160  Double32_t mG1;
161  Double32_t mH1;
162  Double32_t mD1;
163  Double32_t mF1TPerp;
164  Double32_t mF1TPerp1;
165  Double32_t mF1TPerp12;
166  Double32_t mH1Perp;
167  Double32_t mH1Perp1;
168  Double32_t mH1Perp12;
169  Double32_t mAutSiv;
170  Double32_t mAutWtSiv;
171  Double32_t mAutSivAllQ;
172  Double32_t mAutWtSivAllQ;
173  Double32_t mAutSivPiDiff;
174  Double32_t mAutWtSivPiDiff;
175  Double32_t mAutCol;
176  Double32_t mAutWtCol;
177  Double32_t mAutTw3Col;
178  Double32_t mAutWtTw3Col;
179  Double32_t mAutColAllQ;
180  Double32_t mAutWtColAllQ;
181  Double32_t mXUnpolarised;
182  Double32_t mXSivers;
183  Double32_t mXCollins;
184 
185  ClassDef(erhic::EventGmcTrans, 1)
186 };
187 
188 } // namespace erhic
189 
190 #endif // INCLUDE_EICSMEAR_ERHIC_EVENTGMCTRANS_H_
Double32_t mZ
z of the produced hadron
Double32_t mNu
Energy of the exchanged boson in the hadron rest frame.
virtual double GetHermesPhiS() const
Double32_t mXCollins
Cross sections.
virtual bool Parse(const std::string &)
Double32_t mLeptonTheta
Polar angle of the scattered lepton.
Double_t GetHadronZ() const
Definition: EventGmcTrans.h:64
Int_t mStruckQuark
Flavour of struck quark.
Double32_t mLeptonPhi
Azimuthal angle of the scattered lepton.
Double32_t mHadronPt
pT of the produced hadron
Double32_t mWSquared
Invariant mass of the hadronic final state.
Double32_t mS
Square of the centre of mass energy.
virtual ~EventGmcTrans()
Definition: EventGmcTrans.h:35
EventGmcTrans(const std::string &="")
Double_t GetPhiHadron() const
Definition: EventGmcTrans.h:59
TLorentzVector GetHadronPolarisation() const
Double_t GetPhiSpin() const
Definition: EventGmcTrans.h:54