/home/docs/checkouts/readthedocs.org/user_builds/ratpac/checkouts/latest/src/gen/include/RAT/ESCrossSec.hh Source File

Ratpac-two: /home/docs/checkouts/readthedocs.org/user_builds/ratpac/checkouts/latest/src/gen/include/RAT/ESCrossSec.hh Source File
Ratpac-two
ESCrossSec.hh
1 
45 #ifndef __RAT_ESCrossSec__
46 #define __RAT_ESCrossSec__
47 
48 // G4 headers
49 #include <globals.hh>
50 // RAT headers
51 #include <RAT/GLG4StringUtil.hh>
52 #include <RAT/LinearInterp.hh>
53 // ROOT headers
54 #include <TGraph.h>
55 
56 namespace RAT {
57 
58 // Forward declarations within the namespace
59 class ESCrossSecMessenger;
60 
61 class ESCrossSec {
62  public:
63  enum NuEType { nue, nuebar, numu, numubar };
64 
65  ESCrossSec(const char *flavor = "nue");
66 
67  ~ESCrossSec();
68 
75  double Sigma(const double Enu) const;
76 
86  double dSigmadT(const double Enu, const double Te) const;
87 
97  double IntegraldSigmadT(const double Enu, const double T1, const double T2) const;
98 
108  double dSigmadCosTh(const double Enu, const double CosTh) const;
109 
110  //-----------------------------------------------------------------------
123  void SetRadiativeCorrection(const int ii);
124 
130  inline int GetRadiativeCorrection() const { return fRadiativeCorrection; };
131 
138  void SetSinThetaW(const double &sintw);
139 
146  double GetSinThetaW() const { return fsinthetaW2; };
147 
153  void SetReaction(const std::string &reaction);
154 
160  const std::string &GetReaction() const { return fReactionStr; };
161 
169  TGraph *DrawdSigmadT(const double Enu) const;
170 
179  double CrossSecNorm() const { return 1e-42; };
180 
181  protected:
185  void LoadTablesDB();
186 
187  void CalcG();
188 
189  private:
190  // Sets the defaults for the calculation
191  void Defaults();
192 
193  // add by y.t. 16-JAN-2003
194  double fL(const double x) const;
195 
196  NuEType fReaction;
197  std::string fReactionStr;
198 
199  double fEmin, fEmax;
200 
201  // Some constants
202  static const double fGf;
203  static const double fhbarc;
204  static const double fhbarc2;
205  static const double falpha;
206 
212  static double fsinthetaW2;
213 
214  double fMe;
215  double fSigmaOverMe;
216 
217  double fgL;
218  double fgR;
219 
220  //-----------------------------------------------------------------------
221  int fRadiativeCorrection; // flag to use radiative correction or not
222 
225  static const double fTableTeMin;
226  static const double fTableTeMax;
227 
228  // Total cross section table
229  bool fExistTableTot;
230  int fNDataTot;
231  double fEnuStepTot;
232  LinearInterp<double> fTableTot_e;
233 
234  // Differential cross section table
235  bool fExistTableDif;
236  int fNDataDif;
237  double fEnuStepDif;
238  LinearInterp<double> fTableDif_e;
239  std::vector<double> fTableDif;
240 
241  ESCrossSecMessenger *fMessenger;
242 };
243 
244 //
245 // Some inline definitions to make the code a bit faster
246 //
247 
248 inline void ESCrossSec::CalcG() {
249  // calculate the gL and gR for
250  // the reaction type
251  if (fReaction == nue) {
252  fgL = 0.5 + fsinthetaW2;
253  fgR = fsinthetaW2;
254  } else if (fReaction == nuebar) {
255  fgL = fsinthetaW2;
256  fgR = 0.5 + fsinthetaW2;
257  } else if (fReaction == numu) {
258  fgL = -0.5 + fsinthetaW2;
259  fgR = fsinthetaW2;
260  } else if (fReaction == numubar) {
261  fgL = fsinthetaW2;
262  fgR = -0.5 + fsinthetaW2;
263  } else {
264  throw std::invalid_argument("Unknown reaction " + fReactionStr);
265  }
266 }
267 } // namespace RAT
268 
269 #endif
Messenger class to control cross section options.
Definition: ESCrossSecMessenger.hh:50
Calculates neutrino-electron elastic scattering. (based on original QSNO code by F....
Definition: ESCrossSec.hh:61
int GetRadiativeCorrection() const
Getter of the cross section calculation type being used.
Definition: ESCrossSec.hh:130
double IntegraldSigmadT(const double Enu, const double T1, const double T2) const
Definition: ESCrossSec.cc:426
void SetSinThetaW(const double &sintw)
Setter for the Weak mixing angle.
Definition: ESCrossSec.cc:568
double dSigmadT(const double Enu, const double Te) const
Calculate the differential cross section for the neutrino energy Enu.
Definition: ESCrossSec.cc:204
void LoadTablesDB()
Internal function to load the data from the DB.
Definition: ESCrossSec.cc:488
void SetRadiativeCorrection(const int ii)
Setter for which calculation to use.
Definition: ESCrossSec.cc:476
TGraph * DrawdSigmadT(const double Enu) const
Return a TGraph with the differential cross section for an incoming neutrino with energy Enu.
Definition: ESCrossSec.cc:575
double CrossSecNorm() const
Definition: ESCrossSec.hh:179
double GetSinThetaW() const
Getter for the Weak mixing angle.
Definition: ESCrossSec.hh:146
const std::string & GetReaction() const
Getter for the neutrino type to be used.
Definition: ESCrossSec.hh:160
void SetReaction(const std::string &reaction)
Setter for the neutrino type to be used.
Definition: ESCrossSec.cc:534
double Sigma(const double Enu) const
Calculate the total cross section for the neutrino energy Enu.
Definition: ESCrossSec.cc:121
double dSigmadCosTh(const double Enu, const double CosTh) const
Calculate the differential cross section as a function of the recoil angle.
Definition: ESCrossSec.cc:442
Definition: CCCrossSecMessenger.hh:29