icemc
|
#include <Taumodel.hh>
Public Member Functions | |
Taumodel () | |
double | GetTauWeight (Primaries *primary1, Settings *settings1, IceModel *antarctica1, Interaction *interaction1, double pnu, int nu_nubar, double &ptauf, int &crust_entered) |
GetTauWeight is the function that will calculate the probability that a tau neutrino will interact along its path through the earth,and the tau will survive the rest of the journey and decay in the ice. This probability is calculated for final energies from 10^15.5 to the energy of the neutrino. More... | |
ClassDef (Taumodel, 1) | |
Public Attributes | |
double | ptauf |
Final energy of the tau. More... | |
double | weight_tau_prob |
Weight for tau neutrino to interact, create a tau, tau survives and decays in the ice. More... | |
int | inu |
double | weight_nu_prob |
Private Member Functions | |
void | GetDensityVectors (IceModel *antarctica1, Interaction *interaction1, Vector nchord, double step, double Distance, int &totalnusteps, int &crust_entered) |
Get Density Vectors sets two density vectors. One has the density at each step along the path, the other has an average density from the starting point to the current step. More... | |
void | GetEnergyVector (double Etau_final, double step, int totalnusteps, int &totalsteps, double &totaltaudistance, double pnu) |
Get Energy Vector sets the energy of tau particle at every step along the path. It starts from the final energy and works back towards the nuetrino interaction point. More... | |
void | GetTauSurvVector (double step, int totalsteps) |
Get Tau Surv Vector gets a vector that is filled with the probability the tau will survive from neutrino interaction point to the ice. More... | |
Private Attributes | |
double | B0 |
double | B1 |
double | E0 |
double | mT |
for B, the tau elecromagnetic energy loss parameter. More... | |
double | cT |
double | A |
Used in Connolly Calc 2011.(d_dzPsurvNu()) More... | |
double | Mn |
vector< double > | mydensityvector |
vector< double > | myavgdensityvector |
vector< double > | myenergyvector |
vector< double > | myPsurvvector |
vector< double > | etaufarray |
vector< double > | PDFarray |
Taumodel::Taumodel | ( | ) |
For Total Tau Survival probability equation n.b. not in SI units. //from Tau neutrino propagaiton and tau energy loss 2005 Dutta, Huang, & Reno. Equation 16 & used in Equation 30.
| m^2/kg |
| m^2/kg | }parameterization using a logarithmic dependence on energy for B,
| eV | the tau elecromagnetic energy loss parameter.
| eV |Mass of Tau
| m |Tau Decay length (86.93 microMeters) | |
| g |nucleon/ proton mass in grams,also equal to 0.938 GeV.
| none |constant that sets the total probability to unity
these last two constants from Connolly Calc 2011, used in d_dzPsurvNu().
Taumodel::ClassDef | ( | Taumodel | , |
1 | |||
) |
|
private |
Get Density Vectors sets two density vectors. One has the density at each step along the path, the other has an average density from the starting point to the current step.
<Filled with the total weight(a running sum) up to the current final energy.>
Get Density Vectors sets two density vectors. One has the density at each step along the path, the other has an average density from the starting point to the current step.
filled with density at that step
avgdensity
avgdensity up to that step for neutrino
total number of steps neutrino will take through earth
|
private |
Get Energy Vector sets the energy of tau particle at every step along the path. It starts from the final energy and works back towards the nuetrino interaction point.
Get Energy Vector sets the energy of tau particle at every step along the path. It starts from the final energy and works back towards the nuetrino interaction point.
calculate the initial energy needed at the step so the tau will end at the correct final energy
start at nuexit, work backwards
get the density at this step
E_i=E_last +dE/dz*dz
as long as the energy is less than the neutrino energy
distance tau can travel
number of steps tau can take
fill energy vector
Initial energy cannot go above pnu
|
private |
Get Tau Surv Vector gets a vector that is filled with the probability the tau will survive from neutrino interaction point to the ice.
Get Tau Surv Vector gets a vector that is filled with the probability the tau will survive from neutrino interaction point to the ice.
tau surv vector
energy vector starts at the endpoint and goes to earth_in;
calculate chance for tau to survive
is the tau runs out of energy, it wont make it to the ice to decay
fill the vector
double Taumodel::GetTauWeight | ( | Primaries * | primary1, |
Settings * | settings1, | ||
IceModel * | antarctica1, | ||
Interaction * | interaction1, | ||
double | pnu, | ||
int | nu_nubar, | ||
double & | ptauf, | ||
int & | crust_entered | ||
) |
GetTauWeight is the function that will calculate the probability that a tau neutrino will interact along its path through the earth,and the tau will survive the rest of the journey and decay in the ice. This probability is calculated for final energies from 10^15.5 to the energy of the neutrino.
GetTauWeight is the function that will calculate the probability that a tau neutrino will interact along its path through the earth,and the tau will survive the rest of the journey and decay in the ice. This probability is calculated for final energies from 10^15.5 to the energy of the neutrino.
vector from earth_in to "interaction point". Sets the path direction
normalized chord3
Bring in useful variables from other classes
Find the chord, its length and its unit vector.
posnu-earth_in;
normalized chord3
total tau survival probability.
how big is the step size
set up stuff to be used later.
where the particle enters the ice.
Get the density vectors this function needs
set energy vector
set tau surv vector
Set the starting position for tau. First possible interaction point for neutrino
step number to get the correct initial energy at that point;
vector pointing to the step we are currently on.
average density neutrino has seen to that step
density of the current step
Energy vector is filled backwards (final to initial), pull out correct energy
chance tau survives from creation to ice
inelasticity
if neutrino interacts in ice, the tau will already be in ice.
prob neutrino survives to this point
how the initial tau energy is related to final
probability that tau survives to ice
chance to decay in ice
chance to decay in ice if already in ice
total prob for neutrino to make tau and for tau to decay in ice at this energy
dP/dE ->dP/d(log10(E))
dP/dlog10(E) * dlog10(E)
We have the weight. Now use a PDF to find the final energy of the tau.///
|
private |
Used in Connolly Calc 2011.(d_dzPsurvNu())
constant that sets the total probability to unity>
|
private |
For Tau Weight & Survival probability equations. n.b. not in SI units. from Tau neutrino propagaiton and tau energy loss 2005 Dutta, Huang, & Reno. Equation 16 & used in Equation 30.
|
private |
|
private |
Tau Decay length in cm>
|
private |
parameterization using a logarithmic dependence on energy>
|
private |
int Taumodel::inu |
|
private |
|
private |
for B, the tau elecromagnetic energy loss parameter.
mass of the Tau in Gev>
|
private |
Filled with avg density from earth-in to interaction point.>
|
private |
< nucleon/ proton mass in grams,also equal to 0.938 GeV.> Filled with density at each step>
|
private |
vector to hold initial energies.>
|
private |
vector to hold the chance tau would survive from that step.>
|
private |
<filled with all the possible tau final energies.>
double Taumodel::ptauf |
Final energy of the tau.
double Taumodel::weight_nu_prob |
double Taumodel::weight_tau_prob |
Weight for tau neutrino to interact, create a tau, tau survives and decays in the ice.