Jpp test-rotations-old-533-g2bdbdb559
the software that should make you happy
Loading...
Searching...
No Matches
JPHYSICS::JDeltaRays Struct Reference

Auxiliary data structure for delta-rays. More...

#include <JDeltaRays.hh>

Classes

struct  FIT
 Auxiliary data structure to encapsulate energy loss methods based on fit. More...
 
struct  JFormFactor
 Auxiliary class for 2nd order polynomial form factor. More...
 

Static Public Member Functions

static constexpr double getZoverA (const JSeaWater::atom_type &parameters)
 Get ratio charge to mass number for given atomic parameters.
 
static constexpr double getZoverA ()
 Get average ratio charge to mass number for sea water.
 
static double getTmin ()
 Get minimum delta-ray kinetic energy.
 
static double getTmax (const double E, const double M)
 Get maximum delta-ray kinetic energy for given lepton energy and mass.
 
template<class JFormFactor_t >
static double getCount (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 number of delta-rays per unit track length for an ionising particle with given energy and given mass.
 
template<class JFormFactor_t >
static double getEnergyLoss (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 loss due to delta-rays per unit track length
for an ionising particle with given energy and given mass and for a given form factor.
 
static double getCount (const double E, const double M, const double Tmin, const double Tmax, const double Z, const double A)
 Get number of delta-rays per unit track length for an ionising particle with given energy and given mass.
 
static double getEnergyLoss (const double E, const double M, const double Tmin, const double Tmax, const double Z, const double A)
 Get equivalent EM-shower energy loss due to delta-rays per unit track length
for an ionising particle with given energy and given mass.
 
static double getCount (const double E, const double M, const double Tmin, const double Tmax)
 Get number of delta-rays per unit track length for an ionising particle with given energy and given mass in sea water.
 
static double getEnergyLoss (const double E, const double M, const double Tmin, const double Tmax)
 Get equivalent EM-shower energy loss due to delta-rays per unit track length
for an ionising particle with given energy and given mass in sea water.
 
static double getEnergylossFromElectron (const double E, const JEnergyRange T_GeV=JEnergyRange(TMIN_GEV, TMAX_GEV))
 Equivalent EM-shower energy loss due to delta-rays per unit electron track length in sea water.
 
static double getEnergyLossFromMuon (const double E, const JEnergyRange T_GeV=JEnergyRange(TMIN_GEV, TMAX_GEV))
 Equivalent EM-shower energy loss due to delta-rays per unit muon track length in sea water.
 
static double getEnergyLossFromTau (const double E, const JEnergyRange T_GeV=JEnergyRange(TMIN_GEV, TMAX_GEV))
 Equivalent EM-shower energy loss due to delta-rays per unit tau track length in sea water.
 
static double getProbability (const double x)
 Emission profile of photons from delta-rays.
 

Static Public Attributes

static constexpr double TMIN_GEV = 0.000915499
 Minimum allowed delta-ray kinetic energy [GeV].
 
static constexpr double TMAX_GEV = 1.0e10
 Maximum allowed delta-ray kinetic energy [GeV].
 
static constexpr double K = 0.307075
 internal constant [MeV mol^-1 cm^2]
 

Detailed Description

Auxiliary data structure for delta-rays.

Definition at line 21 of file JDeltaRays.hh.

Member Function Documentation

◆ getZoverA() [1/2]

static constexpr double JPHYSICS::JDeltaRays::getZoverA ( const JSeaWater::atom_type & parameters)
inlinestaticconstexpr

Get ratio charge to mass number for given atomic parameters.

Parameters
parametersatomic parameters
Returns
ratio charge to mass number

Definition at line 34 of file JDeltaRays.hh.

35 {
36 return parameters() * parameters.Z / parameters.A;
37 }

◆ getZoverA() [2/2]

static constexpr double JPHYSICS::JDeltaRays::getZoverA ( )
inlinestaticconstexpr

Get average ratio charge to mass number for sea water.

Returns
ratio charge to mass number

Definition at line 45 of file JDeltaRays.hh.

46 {
47 return (getZoverA(JSeaWater::H) +
51 }
static constexpr double getZoverA()
Get average ratio charge to mass number for sea water.
Definition JDeltaRays.hh:45
static constexpr atom_type Cl
Definition JSeaWater.hh:32
static constexpr atom_type H
Definition JSeaWater.hh:29
static constexpr atom_type O
Definition JSeaWater.hh:30
static constexpr atom_type Na
Definition JSeaWater.hh:31

◆ getTmin()

static double JPHYSICS::JDeltaRays::getTmin ( )
inlinestatic

Get minimum delta-ray kinetic energy.

Returns
minimum delta-ray kinetic energy [GeV]

Definition at line 59 of file JDeltaRays.hh.

60 {
61 const double Emin = MASS_ELECTRON / getSinThetaC(); // delta-ray Cherenkov threshold [GeV]
62
63 return getKineticEnergy(Emin, MASS_ELECTRON);
64 }
static const double MASS_ELECTRON
electron mass [GeV]
double getSinThetaC()
Get average sine of Cherenkov angle of water corresponding to group velocity.
double getKineticEnergy(const double E, const double m)
Get kinetic energy of particle with given energy and mass.

◆ getTmax()

static double JPHYSICS::JDeltaRays::getTmax ( const double E,
const double M )
inlinestatic

Get maximum delta-ray kinetic energy for given lepton energy and mass.


This formula is taken from reference https://pdg.lbl.gov/2020/reviews/rpp2020-rev-passage-particles-matter.pdf

N.B.: This function should not be used for electrons.
For electrons use 0.5 * getKineticEnergy(E, MASS_ELECTRON) instead.

Parameters
Eparticle energy [GeV]
Mparticle mass [GeV]
Returns
maximum delta-ray kinetic energy [GeV]

Definition at line 79 of file JDeltaRays.hh.

81 {
82 const double ratio = MASS_ELECTRON / M;
83
84 const double gamma = E / M;
85 const double beta = getBeta(gamma);
86
87 return ( (2.0*MASS_ELECTRON*beta*beta*gamma*gamma) /
88 (1.0 + 2.0*gamma*ratio + ratio*ratio) );
89 }
double getBeta(const double gamma)
Get relative velocity given Lorentz boost.

◆ getCount() [1/3]

template<class JFormFactor_t >
static double JPHYSICS::JDeltaRays::getCount ( 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 )
inlinestatic

Get number of delta-rays per unit track length for an ionising particle with given energy and given mass.

The template parameter corresponds to a class which contains an operator()(const double)
to compute the form factor corresponding to a given delta-ray kinetic energy.

Parameters
Eparticle energy [GeV]
Mparticle mass [GeV]
Tminminimum delta-ray kinetic energy [GeV]
Tmaxmaximum delta-ray kinetic energy [GeV]
Zatomic number [unit]
Aatomic mass [g/mol]
Fform factor functor
Nnumber of points for numeric integration
Returns
delta-ray count [g^-1 cm^2]

Definition at line 149 of file JDeltaRays.hh.

157 {
158 if (Tmax > Tmin) {
159
160 const double gamma = E / M;
161 const double beta = getBeta(gamma);
162
163 const double W = 0.5 * K * (Z/A) * (1.0/(beta*beta)) * 1.0e-3; // [GeV g^-1 cm^2]
164
165 const double xmin = 1.0/Tmax;
166 const double xmax = 1.0/Tmin;
167 const double dx = (xmax - xmin) / ((double) N);
168
169 double count = 0.0;
170
171 for (double x = xmin; x <= xmax; x += dx) {
172
173 const double T = 1.0/x;
174 const double y = W * F(T) * dx; // [g^-1 cm^2]
175
176 count += y;
177 }
178
179 return count; // [g^-1 cm^2]
180
181 } else {
182
183 return 0.0;
184 }
185 }
static constexpr double K
internal constant [MeV mol^-1 cm^2]
Definition JDeltaRays.hh:25

◆ getEnergyLoss() [1/3]

template<class JFormFactor_t >
static double JPHYSICS::JDeltaRays::getEnergyLoss ( 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 )
inlinestatic

Get equivalent EM-shower energy loss due to delta-rays per unit track length
for an ionising particle with given energy and given mass and for a given form factor.

The template parameter corresponds to a class which contains an operator()(const double)
to compute the form factor corresponding to a given delta-ray kinetic energy.

Parameters
Eparticle energy [GeV]
Mparticle mass [GeV]
Tminminimum delta-ray kinetic energy [GeV]
Tmaxmaximum delta-ray kinetic energy [GeV]
Zatomic number [unit]
Aatomic mass [g/mol]
Fform factor functor
Nnumber of points for numeric integration
Returns
equivalent energy loss [GeV g^-1 cm^2]

Definition at line 206 of file JDeltaRays.hh.

214 {
215 if (Tmax > Tmin) {
216
217 const double gamma = E / M;
218 const double beta = getBeta(gamma);
219
220 const double W = 0.5 * K * (Z/A) * (1.0/(beta*beta)); // [MeV g^-1 cm^2]
221
222 const double xmin = log(Tmin);
223 const double xmax = log(Tmax);
224 const double dx = (xmax - xmin) / ((double) N);
225
226 double weight = 0.0;
227
228 for (double x = xmin; x <= xmax; x += dx) {
229
230 const double T = exp(x);
231 const double y = W * F(T) * dx; // [MeV g^-1 cm^2]
232
233 weight += y;
234 }
235
236 return weight * 1.0e-3; // [GeV g^-1 cm^2]
237
238 } else {
239
240 return 0.0;
241 }
242 }

◆ getCount() [2/3]

static double JPHYSICS::JDeltaRays::getCount ( const double E,
const double M,
const double Tmin,
const double Tmax,
const double Z,
const double A )
inlinestatic

Get number of delta-rays per unit track length for an ionising particle with given energy and given mass.

N.B: For this function a form factor $F(T) = \frac{1}{2E^2}T^2 - \frac{\beta^2}{T_{\rm max}} + 1$ is assumed, where:

  • $T$ corresponds to the delta-ray kinetic energy,
  • $T_{\rm max}$ to the maximum delta-ray kinetic energy,
  • $E$ to the ionising particle's energy and
  • $\beta$ to the corresponding particle velocity, relative to the speed of light.
Parameters
Eparticle energy [GeV]
Mparticle mass [GeV]
Tminminimum delta-ray kinetic energy [GeV]
Tmaxmaximum delta-ray kinetic energy [GeV]
Zatomic number [unit]
Aatomic mass [g/mol]
Returns
delta-ray count [g^-1 cm^2]

Definition at line 262 of file JDeltaRays.hh.

268 {
269 if (Tmax > Tmin) {
270
271 const double gamma = E / M;
272 const double beta = getBeta(gamma);
273
274 const double a = 0.5/(E*E);
275 const double b = beta*beta/Tmax;
276 const double c = 1.0;
277
278 const double W = 0.5 * K * (Z/A) * (1.0/(beta*beta)); // [MeV g^-1 cm^2]
279
280 const double dT = Tmax - Tmin;
281 const double rT = Tmax / Tmin;
282 const double mT = Tmax * Tmin;
283
284 return W * (a*dT - b*log(rT) + c*dT/mT) * 1.0e-3; // [g^-1 cm^2]
285
286 } else {
287
288 return 0.0;
289 }
290 }

◆ getEnergyLoss() [2/3]

static double JPHYSICS::JDeltaRays::getEnergyLoss ( const double E,
const double M,
const double Tmin,
const double Tmax,
const double Z,
const double A )
inlinestatic

Get equivalent EM-shower energy loss due to delta-rays per unit track length
for an ionising particle with given energy and given mass.



N.B: For this function a form factor $F(T) = \frac{1}{2E^2}T^2 - \frac{\beta^2}{T_{\rm max}} + 1$ is assumed, where:

  • $T$ corresponds to the delta-ray kinetic energy,
  • $T_{\rm max}$ to the maximum delta-ray kinetic energy,
  • $E$ to the ionising particle's energy and
  • $\beta$ to the corresponding particle velocity, relative to the speed of light.
Parameters
Eparticle energy [GeV]
Mparticle mass [GeV]
Tminminimum delta-ray kinetic energy [GeV]
Tmaxmaximum delta-ray kinetic energy [GeV]
Zatomic number [unit]
Aatomic mass [g/mol]
Returns
equivalent energy loss [GeV g^-1 cm^2]

Definition at line 311 of file JDeltaRays.hh.

317 {
318 if (Tmax > Tmin) {
319
320 const double gamma = E / M;
321 const double beta = getBeta(gamma);
322
323 const double a = 0.25/(E*E);
324 const double b = beta*beta/Tmax;
325 const double c = 1.0;
326
327 const double W = 0.5 * K * (Z/A) * (1.0/(beta*beta)); // [MeV g^-1 cm^2]
328
329 const double sT = Tmax + Tmin;
330 const double dT = Tmax - Tmin;
331 const double rT = Tmax / Tmin;
332
333 return W * (a*sT*dT - b*dT + c*log(rT)) * 1.0e-3; // [GeV g^-1 cm^2]
334
335 } else {
336
337 return 0.0;
338 }
339 }

◆ getCount() [3/3]

static double JPHYSICS::JDeltaRays::getCount ( const double E,
const double M,
const double Tmin,
const double Tmax )
inlinestatic

Get number of delta-rays per unit track length for an ionising particle with given energy and given mass in sea water.

Parameters
Eparticle energy [GeV]
Mparticle mass [GeV]
Tminminimum delta-ray kinetic energy [GeV]
Tmaxmaximum delta-ray kinetic energy [GeV]
Returns
delta-ray count [g^-1 cm^2]

Definition at line 351 of file JDeltaRays.hh.

355 {
356 return getCount(E, M, Tmin, Tmax, 1.0, 1.0) * getZoverA();
357 }
static double getCount(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 number of delta-rays per unit track length for an ionising particle with given energy and given m...

◆ getEnergyLoss() [3/3]

static double JPHYSICS::JDeltaRays::getEnergyLoss ( const double E,
const double M,
const double Tmin,
const double Tmax )
inlinestatic

Get equivalent EM-shower energy loss due to delta-rays per unit track length
for an ionising particle with given energy and given mass in sea water.

Parameters
Eparticle energy [GeV]
Mparticle mass [GeV]
Tminminimum delta-ray kinetic energy [GeV]
Tmaxmaximum delta-ray kinetic energy [GeV]
Returns
equivalent energy loss [GeV g^-1 cm^2]

Definition at line 370 of file JDeltaRays.hh.

374 {
375 return getEnergyLoss(E, M, Tmin, Tmax, 1.0, 1.0) * getZoverA();
376 }
static double getEnergyLoss(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 loss due to delta-rays per unit track length for an ionising particle...

◆ getEnergylossFromElectron()

static double JPHYSICS::JDeltaRays::getEnergylossFromElectron ( const double E,
const JEnergyRange T_GeV = JEnergyRange(TMIN_GEVTMAX_GEV) )
inlinestatic

Equivalent EM-shower energy loss due to delta-rays per unit electron track length in sea water.

Parameters
Eelectron energy [GeV]
T_GeVallowed kinetic energy range of delta-rays [GeV]
Returns
equivalent energy loss [GeV/m]

Definition at line 386 of file JDeltaRays.hh.

388 {
389 const double Tmin = T_GeV.constrain(getTmin());
390 const double Tmax = T_GeV.constrain(0.5 * getKineticEnergy(E, MASS_ELECTRON));
391
392 return getEnergyLoss(E, MASS_ELECTRON, Tmin, Tmax) * DENSITY_SEA_WATER * 1.0e2;
393 }
static const double DENSITY_SEA_WATER
Fixed environment values.
static double getTmin()
Get minimum delta-ray kinetic energy.
Definition JDeltaRays.hh:59

◆ getEnergyLossFromMuon()

static double JPHYSICS::JDeltaRays::getEnergyLossFromMuon ( const double E,
const JEnergyRange T_GeV = JEnergyRange(TMIN_GEVTMAX_GEV) )
inlinestatic

Equivalent EM-shower energy loss due to delta-rays per unit muon track length in sea water.

Parameters
Emuon energy [GeV]
T_GeVallowed kinetic energy range of delta-rays [GeV]
Returns
equivalent energy loss [GeV/m]

Definition at line 403 of file JDeltaRays.hh.

405 {
406 const double Tmin = T_GeV.constrain(getTmin());
407 const double Tmax = T_GeV.constrain(getTmax(E, MASS_MUON));
408
409 return getEnergyLoss(E, MASS_MUON, Tmin, Tmax) * DENSITY_SEA_WATER * 1.0e2;
410 }
static double getTmax(const double E, const double M)
Get maximum delta-ray kinetic energy for given lepton energy and mass.
Definition JDeltaRays.hh:79

◆ getEnergyLossFromTau()

static double JPHYSICS::JDeltaRays::getEnergyLossFromTau ( const double E,
const JEnergyRange T_GeV = JEnergyRange(TMIN_GEVTMAX_GEV) )
inlinestatic

Equivalent EM-shower energy loss due to delta-rays per unit tau track length in sea water.

Parameters
Etau energy [GeV]
T_GeVallowed kinetic energy range of delta-rays [GeV]
Returns
equivalent energy loss [GeV/m]

Definition at line 420 of file JDeltaRays.hh.

422 {
423 const double Tmin = T_GeV.constrain(getTmin());
424 const double Tmax = T_GeV.constrain(getTmax(E, MASS_TAU));
425
426 return getEnergyLoss(E, MASS_TAU, Tmin, Tmax) * DENSITY_SEA_WATER * 1.0e2;
427 }
static const double MASS_TAU
tau mass [GeV]

◆ getProbability()

static double JPHYSICS::JDeltaRays::getProbability ( const double x)
inlinestatic

Emission profile of photons from delta-rays.

Profile is taken from reference ANTARES-SOFT-2002-015, J. Brunner (fig. 3).

Parameters
xcosine emission angle
Returns
probability

Definition at line 500 of file JDeltaRays.hh.

501 {
502 //return 1.0 / (4.0 * PI);
503 return 0.188 * exp(-1.25 * pow(fabs(x - 0.90), 1.30));
504 }
T pow(const T &x, const double y)
Power .
Definition JMath.hh:97

Member Data Documentation

◆ TMIN_GEV

double JPHYSICS::JDeltaRays::TMIN_GEV = 0.000915499
staticconstexpr

Minimum allowed delta-ray kinetic energy [GeV].

Definition at line 23 of file JDeltaRays.hh.

◆ TMAX_GEV

double JPHYSICS::JDeltaRays::TMAX_GEV = 1.0e10
staticconstexpr

Maximum allowed delta-ray kinetic energy [GeV].

Definition at line 24 of file JDeltaRays.hh.

◆ K

double JPHYSICS::JDeltaRays::K = 0.307075
staticconstexpr

internal constant [MeV mol^-1 cm^2]

Definition at line 25 of file JDeltaRays.hh.


The documentation for this struct was generated from the following file: