ESCrossSec Class Reference
|
Ratpac-two
|
Calculates neutrino-electron elastic scattering. (based on original QSNO code by F. Duncan, M. Chen and Y. Takeuchi). More...
#include <ESCrossSec.hh>
Public Types | |
| enum | NuEType { nue , nuebar , numu , numubar } |
Public Member Functions | |
| ESCrossSec (const char *flavor="nue") | |
| double | Sigma (const double Enu) const |
| Calculate the total cross section for the neutrino energy Enu. More... | |
| double | dSigmadT (const double Enu, const double Te) const |
| Calculate the differential cross section for the neutrino energy Enu. More... | |
| double | IntegraldSigmadT (const double Enu, const double T1, const double T2) const |
| double | dSigmadCosTh (const double Enu, const double CosTh) const |
| Calculate the differential cross section as a function of the recoil angle. More... | |
| void | SetRadiativeCorrection (const int ii) |
| Setter for which calculation to use. More... | |
| int | GetRadiativeCorrection () const |
| Getter of the cross section calculation type being used. More... | |
| void | SetSinThetaW (const double &sintw) |
| Setter for the Weak mixing angle. More... | |
| double | GetSinThetaW () const |
| Getter for the Weak mixing angle. More... | |
| void | SetReaction (const std::string &reaction) |
| Setter for the neutrino type to be used. More... | |
| const std::string & | GetReaction () const |
| Getter for the neutrino type to be used. More... | |
| TGraph * | DrawdSigmadT (const double Enu) const |
| Return a TGraph with the differential cross section for an incoming neutrino with energy Enu. More... | |
| double | CrossSecNorm () const |
Protected Member Functions | |
| void | LoadTablesDB () |
| Internal function to load the data from the DB. | |
| void | CalcG () |
Detailed Description
Calculates neutrino-electron elastic scattering. (based on original QSNO code by F. Duncan, M. Chen and Y. Takeuchi).
Code contributed by the SNO+ collaboration
REVISION HISTORY:
23-NOV-2010 - Nuno Barros - Imported and adapted into SNO+ RAT
- Shamelessly taken from SNO QPhysics class PNuE
22-JUN-2012 - Nuno Barros - Cleaned up code and solved a small problem with unit conversion.
07-JUL-2012 - Nuno Barros - Removed all Geant4 streamers replacing them by RAT log objects.
- Changed geant4 data types by C++ data types (suggestion by CIC).
- Revised constness of members.
Calculates the neutrino-electron ES cross section. The calculation is performed either for the ES_e (W+Z) channel, or the ES_mu,tau (Z) channel. Some remarks concerning the calculations.
- Parameters
-
E,Enu : Neutrino energy (MeV). T,Te : Recoil electron energy (MeV).
- Returns
- Cross section in units of $10^{-42}cm^{2}$
Available strategies for the ES cross section calculation (set by messenger):
- 1 : Original routine from QSNO::PNuE (Bahcall).
- 2 : Improved routine from QSNO::PNuE (without rad. corrections).
- 3 : Improved routine from QSNO::PNuE (with rad. corrections - analytical).
- 4 (default) : Improved routine from QSNO::PNuE (with rad. corrections - table).
Member Function Documentation
◆ CrossSecNorm()
|
inline |
Returns the global normalization of the cross section calculation. For precision reasons, the cross-section is performed on a different scale, and therefore any result returned by the calculation is missing this scale, which has to be applied separately.
- Returns
- cross section scaling factor (1e-42)
◆ DrawdSigmadT()
| TGraph * RAT::ESCrossSec::DrawdSigmadT | ( | const double | Enu | ) | const |
Return a TGraph with the differential cross section for an incoming neutrino with energy Enu.
- Parameters
-
Enu Incoming neutrino energy (MeV)
- Returns
- TGraph with the shape of
in units of
.
◆ dSigmadCosTh()
| double RAT::ESCrossSec::dSigmadCosTh | ( | const double | Enu, |
| const double | CosTh | ||
| ) | const |
Calculate the differential cross section as a function of the recoil angle.
- Parameters
-
Enu Incoming neutrino energy. CosTh Cosine of laboratory recoil angle of the electron.
- Returns
in units of
.
◆ dSigmadT()
| double RAT::ESCrossSec::dSigmadT | ( | const double | Enu, |
| const double | Te | ||
| ) | const |
Calculate the differential cross section for the neutrino energy Enu.
- Parameters
-
Enu Incoming neutrino energy (MeV). Te Recoil electron energy (MeV).
- Returns
- Differencial cross section
in units of
.
- Todo:
- Optimize the interpolation.
List of relevant variables at this point:
- k1 : Index of Enu part of the table.
- k2 : Index of the Te part of the table.
- e1,e1 : Closest points in Enu in the table.
- t1,t2 : Closest points in Te in the table.
- sig_e1_t1,sig_e2_t2,sig_e2_t1,sig_e1_t2 :

- Te1,Te2 : Limits on the recoil energy in the table.
◆ GetRadiativeCorrection()
|
inline |
Getter of the cross section calculation type being used.
- Returns
- index of calculation type.
◆ GetReaction()
|
inline |
Getter for the neutrino type to be used.
- Returns
- reaction neutrino type (nue,numu,nuebar,numubar)
◆ GetSinThetaW()
|
inline |
Getter for the Weak mixing angle.
- Returns
- Weak mixing angle :
◆ IntegraldSigmadT()
| double RAT::ESCrossSec::IntegraldSigmadT | ( | const double | Enu, |
| const double | T1, | ||
| const double | T2 | ||
| ) | const |
Integrate the differential cross section between recoil energies T1 and T2.
- Parameters
-
Enu Incoming neutrino energy (MeV). T1 Lower limit of recoil energy interval (MeV). T2 Upper limit of recoil energy interval (MeV).
- Returns
in units of
.
◆ SetRadiativeCorrection()
| void RAT::ESCrossSec::SetRadiativeCorrection | ( | const int | ii | ) |
Setter for which calculation to use.
The available calculations are:
- 1 : No radiative corrections (Bahcall first calculation).
- 2 : No radiative corrections from table.
- 3 : Full analytical calculation with radiative corrections.
- 4 : Full calculations with radiative corrections from table (default).
- Parameters
-
ii Option for calculation to use.
◆ SetReaction()
| void RAT::ESCrossSec::SetReaction | ( | const std::string & | reaction | ) |
Setter for the neutrino type to be used.
- Parameters
-
reaction neutrino type (nue,numu,nuebar,numubar)
◆ SetSinThetaW()
| void RAT::ESCrossSec::SetSinThetaW | ( | const double & | sintw | ) |
Setter for the Weak mixing angle.
- Parameters
-
sintw Weak mixing angle :
◆ Sigma()
| double RAT::ESCrossSec::Sigma | ( | const double | Enu | ) | const |
Calculate the total cross section for the neutrino energy Enu.
Calculates total cross section for a given neutrino energy.
- Parameters
-
Enu
- Returns
- total cross section in units of
.
The documentation for this class was generated from the following files:
- /home/docs/checkouts/readthedocs.org/user_builds/ratpac/checkouts/latest/src/gen/include/RAT/ESCrossSec.hh
- /home/docs/checkouts/readthedocs.org/user_builds/ratpac/checkouts/latest/src/gen/src/ESCrossSec.cc
Generated by