1 #ifndef __JPHYSICS__JPDFTOOLKIT__ 
    2 #define __JPHYSICS__JPDFTOOLKIT__ 
   62     const double x = 
n*lambda;
 
   98     const double gamma = E / M;
 
   99     const double beta  = sqrt((1.0 + 1.0/gamma) * (1.0 - 1.0/gamma));    
 
  102              (1.0 + 2.0*gamma*ratio + ratio*ratio) );
 
  133       const double K = 0.307075;                            
 
  135       const double gamma = E / M;
 
  136       const double beta  = sqrt((1.0 + 1.0/gamma) * (1 - 1.0/gamma));
 
  138       const double a = 0.25/(E*E);
 
  139       const double b = beta*beta/Tmax;
 
  140       const double c = 1.0;
 
  142       const double W = 0.5 * 
K * (Z/A) * (1.0/(beta*beta)); 
 
  144       const double sT = Tmax + Tmin;
 
  145       const double dT = Tmax - Tmin;
 
  146       const double rT = Tmax / Tmin;
 
  148       const double weight = W * (
a*sT*dT - b*dT + c*log(rT));
 
  175   template<
class JFormFactor_t>
 
  182                              const JFormFactor_t& F,
 
  183                              const int            N = 1000000)
 
  189       const double K = 0.307075;                            
 
  191       const double gamma = E / M;
 
  192       const double beta  = sqrt((1.0 + 1.0/gamma) * (1 - 1.0/gamma));
 
  194       const double W = 0.5 * 
K * (Z/A) * (1.0/(beta*beta)); 
 
  196       const double xmin = log(Tmin);
 
  197       const double xmax = log(Tmax);
 
  198       const double dx   = (
xmax - 
xmin) / ((
double) N);
 
  204         const double T = exp(
x);
 
  205         const double y = W * F(T) * dx;                     
 
  236     static const double Z = 10.0; 
 
  237     static const double A = 18.0; 
 
  263     static const double Z = 10.0; 
 
  264     static const double A = 18.0; 
 
  288     static const double a =  3.195e-01;
 
  289     static const double b =  3.373e-01;
 
  290     static const double c = -2.731e-02;
 
  291     static const double d =  1.610e-03;
 
  292     static const double Emin = 0.13078; 
 
  296       const double x = log10(E);                      
 
  297       const double y = 
a + 
x*(b + 
x*(c + 
x*(d)));     
 
  325     static const double Z = 10.0; 
 
  326     static const double A = 18.0; 
 
  350     static const double a = -2.178e-01;
 
  351     static const double b =  4.995e-01;
 
  352     static const double c = -3.906e-02;
 
  353     static const double d =  1.615e-03;
 
  354     static const double Emin = 2.19500; 
 
  358       const double x = log10(E);                      
 
  359       const double y = 
a + 
x*(b + 
x*(c + 
x*(d)));     
 
  381     return 0.188 * exp(-1.25 * 
pow(fabs(
x - 0.90), 1.30));
 
  395     static const double d = 0.36;                    
 
  396     static const double U = 
PI*
PI*
PI*
PI*
PI*2.0/3.0;
 
  397     static const double V = d*d*d*d*d*d;
 
  399     const double W     = (
n*
n - 1.0) / (
n*
n + 2.0);
 
  400     const double sigma = 1.0e-14 * U*V*W*W / (lambda*lambda*lambda*lambda);   
 
  416     static const double amu = 18.01528; 
 
Auxiliary class to define a range between two values.
 
static const uint32_t K[64]
 
double getKineticEnergy(const Trk &trk)
Get track kinetic energy.
 
T pow(const T &x, const double y)
Power .
 
static const double PI
Mathematical constants.
 
Auxiliary methods for light properties of deep-sea water.
 
double getDeltaRaysFromMuon(const double E, const JRange< double > T_GeV=JRange< double >(DELTARAY_TMIN, DELTARAY_TMAX))
Equivalent EM-shower energy due to delta-rays per unit muon track length.
 
static const double ALPHA_ELECTRO_MAGNETIC
Electro-Magnetic coupling constant.
 
static const double DELTARAY_TMIN
Minimum allowed delta-ray kinetic energy [GeV].
 
static const double DELTARAY_TMAX
Maximum allowed delta-ray kinetic energy [GeV].
 
double getDeltaRayTmin()
Get minimum delta-ray kinetic energy.
 
double getDeltaRays(const double E, const double M, const double Tmin, const double Tmax, const double Z, const double A, const JFormFactor_t &F, const int N=1000000)
Get equivalent EM-shower energy due to delta-rays per unit track length for an ionising particle with...
 
double getDeltaRaysFromMuonFit(const double E)
Equivalent EM-shower energy due to delta-rays per unit muon track length.
 
static const double DENSITY_SEA_WATER
Fixed environment values.
 
static const double MASS_MUON
muon mass [GeV]
 
const double getRayleighScatteringLength(const double n, const double lambda)
Rayleigh scattering length.
 
double getDeltaRaysFromTauFit(const double E)
Equivalent EM-shower energy due to delta-rays per unit tau track length.
 
static const double MASS_ELECTRON
electron mass [GeV]
 
double getDeltaRayProbability(const double x)
Emission profile of photons from delta-rays.
 
const double getRayleighCrossSection(const double n, const double lambda)
Rayleigh cross section.
 
double cherenkov(const double lambda, const double n)
Number of Cherenkov photons per unit track length and per unit wavelength.
 
double getMinimalWavelength()
Get minimal wavelength for PDF evaluations.
 
double getDeltaRayTmax(const double E, const double M)
Get maximum delta-ray kinetic energy for given lepton energy and mass.
 
static const double MASS_TAU
tau mass [GeV]
 
static const double AVOGADRO
Avogadro's number.
 
double getSinThetaC()
Get average sine of Cherenkov angle of water corresponding to group velocity.
 
double getMaximalWavelength()
Get maximal wavelength for PDF evaluations.
 
double getDeltaRaysFromTau(const double E, const JRange< double > T_GeV=JRange< double >(DELTARAY_TMIN, DELTARAY_TMAX))
Equivalent EM-shower energy due to delta-rays per unit tau track length.
 
double getDeltaRaysFromElectron(const double E, const JRange< double > T_GeV=JRange< double >(DELTARAY_TMIN, DELTARAY_TMAX))
Equivalent EM-shower energy due to delta-rays per unit electron track length.
 
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
 
Auxiliary data structure to list files in directory.