MP-Gadget  5.0.1.dev1-76bc7d4726-dirty
omega_nu_single.h
Go to the documentation of this file.
1 #ifndef OMEGA_NU_SINGLE_H
2 #define OMEGA_NU_SINGLE_H
6 #include <gsl/gsl_interp.h>
7 
16 #define TNUCMB (pow(4/11.,1/3.)*1.00328)
20 #define NUSPECIES 3
21 
25  double * loga;
26  double * rhonu;
27  gsl_interp * interp;
28  gsl_interp_accel * acc;
29  /*Neutrino mass for this structure*/
30  double mnu;
31 };
32 typedef struct _rho_nu_single _rho_nu_single;
33 
39 void rho_nu_init(_rho_nu_single * rho_nu_tab, double a0, const double mnu, const double kBtnu);
40 
47 double rho_nu(const _rho_nu_single * rho_nu_tab, const double a, const double kT);
48 
53 struct _hybrid_nu {
54  /*True if hybrid neutrinos are enabled*/
55  int enabled;
56  /* This is the fraction of neutrino mass not followed by the analytic integrator.
57  The analytic method is cutoff at q < qcrit (specified using vcrit, below) and use
58  particles for the slower neutrinos.*/
60  /* Time at which to turn on the particle neutrinos.
61  * Ultimately we want something better than this.*/
62  double nu_crit_time;
63  /*Critical velocity as a fraction of lightspeed*/
64  double vcrit;
65 };
66 typedef struct _hybrid_nu _hybrid_nu;
67 
77 void init_hybrid_nu(_hybrid_nu * const hybnu, const double mnu[], const double vcrit, const double light, const double nu_crit_time, const double kBtnu);
78 
86 double particle_nu_fraction(const _hybrid_nu * const hybnu, const double a, int i);
87 
89 double nufrac_low(const double qc);
90 
96 struct _omega_nu {
97  /*Pointers to the array of structures we use to store rho_nu*/
99  /* Which species have the same mass and can thus be counted together.*/
101  /* Prefactor to turn density into matter density omega*/
102  double rhocrit;
103  /*neutrino temperature times Boltzmann constant*/
104  double kBtnu;
105  /*CMB temperature*/
106  double tcmb0;
107  /* Pointer to structure for hybrid neutrinos. */
109 };
110 typedef struct _omega_nu _omega_nu;
111 
118 void init_omega_nu(_omega_nu * const omnu, const double MNu[], const double a0, const double HubbleParam, const double tcmb0);
119 
121 double get_omega_nu(const _omega_nu * const omnu, const double a);
122 
124 double get_omega_nu_nopart(const _omega_nu * const omnu, const double a);
125 
127 double get_omegag(const _omega_nu * const omnu, const double a);
128 
133 double omega_nu_single(const _omega_nu * const rho_nu_tab, const double a, const int i);
134 
135 #endif
void init_omega_nu(_omega_nu *const omnu, const double MNu[], const double a0, const double HubbleParam, const double tcmb0)
double particle_nu_fraction(const _hybrid_nu *const hybnu, const double a, int i)
void init_hybrid_nu(_hybrid_nu *const hybnu, const double mnu[], const double vcrit, const double light, const double nu_crit_time, const double kBtnu)
double rho_nu(const _rho_nu_single *rho_nu_tab, const double a, const double kT)
#define NUSPECIES
double omega_nu_single(const _omega_nu *const rho_nu_tab, const double a, const int i)
double get_omegag(const _omega_nu *const omnu, const double a)
void rho_nu_init(_rho_nu_single *rho_nu_tab, double a0, const double mnu, const double kBtnu)
double nufrac_low(const double qc)
double get_omega_nu(const _omega_nu *const omnu, const double a)
double get_omega_nu_nopart(const _omega_nu *const omnu, const double a)
double nu_crit_time
double nufrac_low[NUSPECIES]
double rhocrit
_hybrid_nu hybnu
int nu_degeneracies[NUSPECIES]
_rho_nu_single RhoNuTab[NUSPECIES]
gsl_interp * interp
gsl_interp_accel * acc