mirror of
https://github.com/cosmo-sims/monofonIC.git
synced 2024-09-19 17:03:45 +02:00
updated default cosmological model to Planck2018, added neutrinos in class module (not in backscaling model)
This commit is contained in:
parent
d7cec75495
commit
c437e6f720
5 changed files with 167 additions and 151 deletions
23
example.conf
23
example.conf
|
@ -30,21 +30,26 @@ ParticleLoad = sc # particle load, can be 'sc' (1x), 'bcc' (2x) or 'fcc
|
|||
## transfer = ... specifies the Einstein-Boltzmann plugin module
|
||||
|
||||
# main cosmological parameters
|
||||
Omega_m = 0.302
|
||||
Omega_b = 0.045
|
||||
Omega_L = 0.698
|
||||
H0 = 70.3
|
||||
nspec = 0.961
|
||||
sigma_8 = 0.811
|
||||
Omega_m = 0.3158
|
||||
Omega_b = 0.0494
|
||||
Omega_L = 0.6842
|
||||
H0 = 67.321
|
||||
n_s = 0.9661
|
||||
sigma_8 = 0.8102
|
||||
|
||||
ZeroRadiation = false # For Back-scaling only: set to true if your simulation code
|
||||
# cannot deal with Omega_r!=0 in its background FLRW model
|
||||
|
||||
## Additional cosmological parameters (set by default to the given values)
|
||||
# w0 = -1.0 # not supported yet!
|
||||
# wa = 0.0 # not supported yet!
|
||||
# Tcmb = 2.7255
|
||||
# Neff = 3.046
|
||||
# k_p = 0.05
|
||||
# N_ur = 2.046
|
||||
# m_nu1 = 0.06
|
||||
# m_nu2 = 0.0
|
||||
# m_nu3 = 0.0
|
||||
# w_0 = -1.0 # not supported yet!
|
||||
# w_a = 0.0 # not supported yet!
|
||||
|
||||
|
||||
## Use below for anisotropic large scale tidal field ICs up to 2LPT
|
||||
## see Stuecker+2020 (https://arxiv.org/abs/2003.06427)
|
||||
|
|
4
external/class.cmake
vendored
4
external/class.cmake
vendored
|
@ -2,7 +2,7 @@ cmake_minimum_required(VERSION 3.11)
|
|||
include(FetchContent)
|
||||
FetchContent_Declare(
|
||||
class
|
||||
GIT_REPOSITORY https://github.com/michaelbuehlmann/class_public.git
|
||||
GIT_REPOSITORY https://github.com/ohahn/class_public.git
|
||||
GIT_TAG master
|
||||
)
|
||||
|
||||
|
@ -10,4 +10,4 @@ FetchContent_GetProperties(class)
|
|||
if(NOT class_POPULATED)
|
||||
FetchContent_Populate(class)
|
||||
add_subdirectory(${class_SOURCE_DIR} ${class_BINARY_DIR})
|
||||
endif()
|
||||
endif()
|
||||
|
|
|
@ -152,7 +152,7 @@ public:
|
|||
*/
|
||||
explicit calculator(config_file &cf)
|
||||
: cosmo_param_(cf), astart_( 1.0/(1.0+cf.get_value<double>("setup","zstart")) ),
|
||||
atarget_( 1.0/(1.0+cf.get_value_safe<double>("cosmology","ztarget",1./astart_-1.)))
|
||||
atarget_( 1.0/(1.0+cf.get_value_safe<double>("cosmology","ztarget",0.0)) )
|
||||
{
|
||||
// pre-compute growth factors and store for interpolation
|
||||
std::vector<double> tab_a, tab_D, tab_f;
|
||||
|
@ -170,9 +170,9 @@ public:
|
|||
// set up transfer functions and compute normalisation
|
||||
transfer_function_ = std::move(select_TransferFunction_plugin(cf, cosmo_param_));
|
||||
transfer_function_->intialise();
|
||||
if( !transfer_function_->tf_isnormalised_ )
|
||||
if( !transfer_function_->tf_isnormalised_ ){
|
||||
cosmo_param_.set("pnorm", this->compute_pnorm_from_sigma8() );
|
||||
else{
|
||||
}else{
|
||||
cosmo_param_.set("pnorm", 1.0/Dplus_target_/Dplus_target_);
|
||||
auto sigma8 = this->compute_sigma8();
|
||||
music::ilog << "Measured sigma_8 for given PS normalisation is " << sigma8 << std::endl;
|
||||
|
@ -323,13 +323,11 @@ public:
|
|||
the transfer function and the window function of 8 Mpc/h at wave number k */
|
||||
static double dSigma8(double k, void *pParams)
|
||||
{
|
||||
if (k <= 0.0)
|
||||
return 0.0f;
|
||||
|
||||
cosmology::calculator *pcc = reinterpret_cast<cosmology::calculator *>(pParams);
|
||||
|
||||
double x = k * 8.0;
|
||||
double w = 3.0 * (sin(x) - x * cos(x)) / (x * x * x);
|
||||
const double x = k * 8.0;
|
||||
const double w = (x < 0.01)? 1.0-0.1*x*x : 3.0 * (std::sin(x) - x * std::cos(x)) / (x * x * x);
|
||||
|
||||
static double nspect = (double)pcc->cosmo_param_["n_s"];
|
||||
double tf = pcc->transfer_function_->compute(k, total);
|
||||
|
||||
|
@ -342,13 +340,11 @@ public:
|
|||
the transfer function and the window function of 8 Mpc/h at wave number k */
|
||||
static double dSigma8_0(double k, void *pParams)
|
||||
{
|
||||
if (k <= 0.0)
|
||||
return 0.0f;
|
||||
|
||||
cosmology::calculator *pcc = reinterpret_cast<cosmology::calculator *>(pParams);
|
||||
|
||||
double x = k * 8.0;
|
||||
double w = 3.0 * (std::sin(x) - x * std::cos(x)) / (x * x * x);
|
||||
const double x = k * 8.0;
|
||||
const double w = (x < 0.01)? 1.0-0.1*x*x : 3.0 * (std::sin(x) - x * std::cos(x)) / (x * x * x);
|
||||
|
||||
static double nspect = static_cast<double>(pcc->cosmo_param_["n_s"]);
|
||||
double tf = pcc->transfer_function_->compute(k, total0);
|
||||
|
||||
|
|
|
@ -1,17 +1,17 @@
|
|||
// This file is part of monofonIC (MUSIC2)
|
||||
// A software package to generate ICs for cosmological simulations
|
||||
// Copyright (C) 2020 by Oliver Hahn
|
||||
//
|
||||
//
|
||||
// monofonIC is free software: you can redistribute it and/or modify
|
||||
// it under the terms of the GNU General Public License as published by
|
||||
// the Free Software Foundation, either version 3 of the License, or
|
||||
// (at your option) any later version.
|
||||
//
|
||||
//
|
||||
// monofonIC is distributed in the hope that it will be useful,
|
||||
// but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
// GNU General Public License for more details.
|
||||
//
|
||||
//
|
||||
// You should have received a copy of the GNU General Public License
|
||||
// along with this program. If not, see <http://www.gnu.org/licenses/>.
|
||||
#pragma once
|
||||
|
@ -24,10 +24,10 @@
|
|||
|
||||
namespace cosmology
|
||||
{
|
||||
//! singleton structure for cosmological parameters
|
||||
class parameters
|
||||
{
|
||||
/*double
|
||||
//! singleton structure for cosmological parameters
|
||||
class parameters
|
||||
{
|
||||
/*double
|
||||
Omega_m, //!< baryon+dark matter density
|
||||
Omega_b, //!< baryon matter density
|
||||
Omega_DE, //!< dark energy density (cosmological constant or parameterised)
|
||||
|
@ -51,128 +51,135 @@ class parameters
|
|||
sqrtpnorm, //!< sqrt of power spectrum normalisation factor
|
||||
vfact; //!< velocity<->displacement conversion factor in Zel'dovich approx.
|
||||
*/
|
||||
private:
|
||||
std::map<std::string,double> pmap_;
|
||||
private:
|
||||
std::map<std::string, double> pmap_;
|
||||
|
||||
public:
|
||||
|
||||
double get( const std::string& key ) const
|
||||
{
|
||||
auto it = pmap_.find(key);
|
||||
if( it == pmap_.end() ){
|
||||
auto errmsg = std::string("Cosmological parameter \'")+key+std::string("\' does not exist in internal list.");
|
||||
music::elog << errmsg << std::endl;
|
||||
throw std::runtime_error(errmsg.c_str());
|
||||
}
|
||||
return it->second;
|
||||
}
|
||||
|
||||
void set( const std::string& key, const double value )
|
||||
{
|
||||
auto it = pmap_.find(key);
|
||||
if( it != pmap_.end() ){
|
||||
pmap_[key] = value;
|
||||
}else{
|
||||
auto errmsg = std::string("Cosmological parameter \'")+key+std::string("\' does not exist in internal list. Needs to be defaulted before it can be set!");
|
||||
music::elog << errmsg << std::endl;
|
||||
throw std::runtime_error(errmsg.c_str());
|
||||
}
|
||||
}
|
||||
|
||||
inline double operator[]( const std::string& key ) const { return this->get( key ); }
|
||||
|
||||
parameters() = delete;
|
||||
|
||||
parameters( const parameters& ) = default;
|
||||
|
||||
explicit parameters( config_file& cf )
|
||||
{
|
||||
pmap_["H0"] = cf.get_value<double>("cosmology", "H0");
|
||||
pmap_["h"] = cf.get_value<double>("cosmology", "H0") / 100.0;
|
||||
pmap_["n_s"] = cf.get_value<double>("cosmology", "nspec");
|
||||
pmap_["Omega_b"] = cf.get_value<double>("cosmology", "Omega_b");
|
||||
pmap_["Omega_m"] = cf.get_value<double>("cosmology", "Omega_m");
|
||||
pmap_["Omega_c"] = this->get("Omega_m") - this->get("Omega_b");
|
||||
pmap_["Omega_DE"] = cf.get_value<double>("cosmology", "Omega_L");
|
||||
pmap_["w_0"] = cf.get_value_safe<double>("cosmology", "w0", -1.0);
|
||||
pmap_["w_a"] = cf.get_value_safe<double>("cosmology", "wa", 0.0);
|
||||
pmap_["Tcmb"] = cf.get_value_safe<double>("cosmology", "Tcmb", 2.7255);
|
||||
pmap_["YHe"] = cf.get_value_safe<double>("cosmology", "YHe", 0.248);
|
||||
pmap_["Neff"] = cf.get_value_safe<double>("cosmology", "Neff", 3.046);
|
||||
pmap_["sigma_8"] = cf.get_value_safe<double>("cosmology", "sigma_8",-1.0);
|
||||
pmap_["A_s"] = cf.get_value_safe<double>("cosmology", "A_s",-1.0);
|
||||
pmap_["Omega_gamma"] = 4 * phys_const::sigma_SI / std::pow(phys_const::c_SI, 3) * std::pow(this->get("Tcmb"), 4.0) / phys_const::rhocrit_h2_SI / (this->get("h") * this->get("h"));
|
||||
pmap_["Omega_nu"] = this->get("Neff") * this->get("Omega_gamma") * 7. / 8. * std::pow(4. / 11., 4. / 3.);
|
||||
pmap_["Omega_r"] = this->get("Omega_gamma") + this->get("Omega_nu");
|
||||
// H0 = cf.get_value<double>("cosmology", "H0");
|
||||
// h = H0 / 100.0;
|
||||
|
||||
// nspect = cf.get_value<double>("cosmology", "nspec");
|
||||
|
||||
// Omega_b = cf.get_value<double>("cosmology", "Omega_b");
|
||||
|
||||
// Omega_m = cf.get_value<double>("cosmology", "Omega_m");
|
||||
|
||||
// Omega_DE = cf.get_value<double>("cosmology", "Omega_L");
|
||||
|
||||
// w_0 = cf.get_value_safe<double>("cosmology", "w0", -1.0);
|
||||
|
||||
// w_a = cf.get_value_safe<double>("cosmology", "wa", 0.0);
|
||||
|
||||
// Tcmb = cf.get_value_safe<double>("cosmology", "Tcmb", 2.7255);
|
||||
|
||||
// YHe = cf.get_value_safe<double>("cosmology", "YHe", 0.248);
|
||||
|
||||
// Neff = cf.get_value_safe<double>("cosmology", "Neff", 3.046);
|
||||
|
||||
// sigma8 = cf.get_value_safe<double>("cosmology", "sigma_8",-1.0);
|
||||
|
||||
// calculate energy density in ultrarelativistic species from Tcmb and Neff
|
||||
// double Omega_gamma = 4 * phys_const::sigma_SI / std::pow(phys_const::c_SI, 3) * std::pow(Tcmb, 4.0) / phys_const::rhocrit_h2_SI / (h * h);
|
||||
// double Omega_nu = Neff * Omega_gamma * 7. / 8. * std::pow(4. / 11., 4. / 3.);
|
||||
// Omega_r = Omega_gamma + Omega_nu;
|
||||
|
||||
if (cf.get_value_safe<bool>("cosmology", "ZeroRadiation", false))
|
||||
public:
|
||||
double get(const std::string &key) const
|
||||
{
|
||||
pmap_["Omega_r"] = 0.0;
|
||||
auto it = pmap_.find(key);
|
||||
if (it == pmap_.end())
|
||||
{
|
||||
auto errmsg = std::string("Cosmological parameter \'") + key + std::string("\' does not exist in internal list.");
|
||||
music::elog << errmsg << std::endl;
|
||||
throw std::runtime_error(errmsg.c_str());
|
||||
}
|
||||
return it->second;
|
||||
}
|
||||
|
||||
pmap_["f_b"] = this->get("Omega_b") / this->get("Omega_m");
|
||||
pmap_["f_c"] = 1.0 - this->get("f_b");
|
||||
void set(const std::string &key, const double value)
|
||||
{
|
||||
auto it = pmap_.find(key);
|
||||
if (it != pmap_.end())
|
||||
{
|
||||
pmap_[key] = value;
|
||||
}
|
||||
else
|
||||
{
|
||||
auto errmsg = std::string("Cosmological parameter \'") + key + std::string("\' does not exist in internal list. Needs to be defaulted before it can be set!");
|
||||
music::elog << errmsg << std::endl;
|
||||
throw std::runtime_error(errmsg.c_str());
|
||||
}
|
||||
}
|
||||
|
||||
inline double operator[](const std::string &key) const { return this->get(key); }
|
||||
|
||||
parameters() = delete;
|
||||
|
||||
parameters(const parameters &) = default;
|
||||
|
||||
explicit parameters(config_file &cf)
|
||||
{
|
||||
// CMB
|
||||
pmap_["Tcmb"] = cf.get_value_safe<double>("cosmology", "Tcmb", 2.7255);
|
||||
pmap_["YHe"] = cf.get_value_safe<double>("cosmology", "YHe", 0.2454006);
|
||||
|
||||
// H0
|
||||
pmap_["H0"] = cf.get_value<double>("cosmology", "H0");
|
||||
pmap_["h"] = cf.get_value<double>("cosmology", "H0") / 100.0;
|
||||
const double h = pmap_["h"];
|
||||
|
||||
// primordial and normalisation
|
||||
if(!cf.contains_key("cosmology/n_s"))
|
||||
pmap_["n_s"] = cf.get_value<double>("cosmology", "nspec");
|
||||
else
|
||||
pmap_["n_s"] = cf.get_value<double>("cosmology", "n_s");
|
||||
pmap_["A_s"] = cf.get_value_safe<double>("cosmology", "A_s", -1.0);
|
||||
pmap_["k_p"] = cf.get_value_safe<double>("cosmology", "k_p", 0.05);
|
||||
pmap_["sigma_8"] = cf.get_value_safe<double>("cosmology", "sigma_8", -1.0);
|
||||
|
||||
// baryon and non-relativistic matter content
|
||||
pmap_["Omega_b"] = cf.get_value<double>("cosmology", "Omega_b");
|
||||
pmap_["Omega_m"] = cf.get_value<double>("cosmology", "Omega_m");
|
||||
|
||||
// massive neutrino species
|
||||
pmap_["m_nu1"] = cf.get_value_safe<double>("cosmology", "m_nu1", 0.06);
|
||||
pmap_["m_nu2"] = cf.get_value_safe<double>("cosmology", "m_nu2", 0.0);
|
||||
pmap_["m_nu3"] = cf.get_value_safe<double>("cosmology", "m_nu3", 0.0);
|
||||
pmap_["N_nu_massive"] = int(this->get("m_nu1") > 1e-9) + int(this->get("m_nu2") > 1e-9) + int(this->get("m_nu3") > 1e-9);
|
||||
const double sum_m_nu = this->get("m_nu1") + this->get("m_nu2") + this->get("m_nu3");
|
||||
|
||||
// number ultrarelativistic neutrinos
|
||||
pmap_["N_ur"] = cf.get_value_safe<double>("cosmology", "N_ur", 3.046 - this->get("N_nu_massive"));
|
||||
pmap_["Omega_nu_massive"] = sum_m_nu / (93.14 * h * h); // Omega_nu_m = \sum_i m_i / (93.14 eV h^2)
|
||||
|
||||
// compute amount of cold dark matter as the rest
|
||||
pmap_["Omega_c"] = this->get("Omega_m") - this->get("Omega_b") - this->get("Omega_nu_massive");
|
||||
|
||||
// calculate energy density in ultrarelativistic species from Tcmb and Neff
|
||||
// photons
|
||||
pmap_["Omega_gamma"] = 4 * phys_const::sigma_SI / std::pow(phys_const::c_SI, 3) * std::pow(this->get("Tcmb"), 4.0)
|
||||
/ phys_const::rhocrit_h2_SI / (this->get("h") * this->get("h"));
|
||||
// massless neutrinos
|
||||
pmap_["Omega_nu_massless"] = this->get("N_ur") * this->get("Omega_gamma") * 7. / 8. * std::pow(4. / 11., 4. / 3.);
|
||||
// total relativistic
|
||||
pmap_["Omega_r"] = this->get("Omega_gamma") + this->get("Omega_nu_massless");
|
||||
|
||||
// dark energy
|
||||
pmap_["Omega_DE"] = cf.get_value<double>("cosmology", "Omega_L");
|
||||
pmap_["w_0"] = cf.get_value_safe<double>("cosmology", "w_0", -1.0);
|
||||
pmap_["w_a"] = cf.get_value_safe<double>("cosmology", "w_a", 0.0);
|
||||
|
||||
if (cf.get_value_safe<bool>("cosmology", "ZeroRadiation", false))
|
||||
{
|
||||
pmap_["Omega_r"] = 0.0;
|
||||
}
|
||||
|
||||
pmap_["f_b"] = this->get("Omega_b") / this->get("Omega_m");
|
||||
pmap_["f_c"] = 1.0 - this->get("f_b"); // this means we add massive neutrinos to CDM here
|
||||
|
||||
#if 1
|
||||
// assume zero curvature, take difference from dark energy
|
||||
pmap_["Omega_DE"] += 1.0 - this->get("Omega_m") - this->get("Omega_DE") - this->get("Omega_r");
|
||||
// Omega_DE += 1.0 - Omega_m - Omega_DE - Omega_r;
|
||||
pmap_["Omega_k"] = 0.0;
|
||||
// assume zero curvature, take difference from dark energy
|
||||
pmap_["Omega_DE"] += 1.0 - this->get("Omega_m") - this->get("Omega_DE") - this->get("Omega_r");
|
||||
// Omega_DE += 1.0 - Omega_m - Omega_DE - Omega_r;
|
||||
pmap_["Omega_k"] = 0.0;
|
||||
#else
|
||||
// allow for curvature
|
||||
Omega_k = 1.0 - Omega_m - Omega_DE - Omega_r;
|
||||
// allow for curvature
|
||||
Omega_k = 1.0 - Omega_m - Omega_DE - Omega_r;
|
||||
#endif
|
||||
|
||||
pmap_["dplus"] = 0.0;
|
||||
pmap_["pnorm"] = 0.0;
|
||||
pmap_["sqrtpnorm"] = 0.0;
|
||||
pmap_["vfact"] = 0.0;
|
||||
pmap_["dplus"] = 0.0;
|
||||
pmap_["pnorm"] = 0.0;
|
||||
pmap_["sqrtpnorm"] = 0.0;
|
||||
pmap_["vfact"] = 0.0;
|
||||
|
||||
music::ilog << "-------------------------------------------------------------------------------" << std::endl;
|
||||
music::ilog << "Cosmological parameters are: " << std::endl;
|
||||
music::ilog << " H0 = " << std::setw(16) << this->get("H0") << "sigma_8 = " << std::setw(16) << this->get("sigma_8") << std::endl;
|
||||
music::ilog << " Omega_c = " << std::setw(16) << this->get("Omega_c") << "Omega_b = " << std::setw(16) << this->get("Omega_b") << std::endl;
|
||||
if (!cf.get_value_safe<bool>("cosmology", "ZeroRadiation", false)){
|
||||
music::ilog << " Omega_g = " << std::setw(16) << this->get("Omega_gamma") << "Omega_nu = " << std::setw(16) << this->get("Omega_nu") << std::endl;
|
||||
}else{
|
||||
music::ilog << " Omega_r = " << std::setw(16) << this->get("Omega_r") << std::endl;
|
||||
music::ilog << "-------------------------------------------------------------------------------" << std::endl;
|
||||
music::ilog << "Cosmological parameters are: " << std::endl;
|
||||
music::ilog << " h = " << std::setw(16) << this->get("h");
|
||||
if( this->get("A_s") > 0.0 )
|
||||
music::ilog << "A_s = " << std::setw(16) << this->get("A_s");
|
||||
else
|
||||
music::ilog << "sigma_8 = " << std::setw(16) << this->get("sigma_8");
|
||||
music::ilog << "nspec = " << std::setw(16) << this->get("n_s") << std::endl;
|
||||
music::ilog << " Omega_c = " << std::setw(16) << this->get("Omega_c") << "Omega_b = " << std::setw(16) << this->get("Omega_b") << "Omega_m = " << std::setw(16) << this->get("Omega_m") << std::endl;
|
||||
music::ilog << " Omega_r = " << std::setw(16) << this->get("Omega_r") << "Omega_nu = " << std::setw(16) << this->get("Omega_nu_massive") << "∑m_nu = " << std::setw(11) << sum_m_nu << "eV" << std::endl;
|
||||
music::ilog << " Omega_DE = " << std::setw(16) << this->get("Omega_DE") << "w_0 = " << std::setw(16) << this->get("w_0") << "w_a = " << std::setw(16) << this->get("w_a") << std::endl;
|
||||
|
||||
if (this->get("Omega_r") > 0.0)
|
||||
{
|
||||
music::wlog << " Radiation enabled, using Omega_r=" << this->get("Omega_r") << " internally for backscaling." << std::endl;
|
||||
music::wlog << " Make sure your sim code supports this, otherwise set [cosmology] / ZeroRadiation=true." << std::endl;
|
||||
}
|
||||
}
|
||||
music::ilog << " Omega_DE = " << std::setw(16) << this->get("Omega_DE") << "nspect = " << std::setw(16) << this->get("n_s") << std::endl;
|
||||
music::ilog << " w0 = " << std::setw(16) << this->get("w_0") << "w_a = " << std::setw(16) << this->get("w_a") << std::endl;
|
||||
|
||||
if( this->get("Omega_r") > 0.0 )
|
||||
{
|
||||
music::wlog << "Radiation enabled, using Omega_r=" << this->get("Omega_r") << " internally."<< std::endl;
|
||||
music::wlog << "Make sure your sim code supports this..." << std::endl;
|
||||
}
|
||||
}
|
||||
|
||||
};
|
||||
};
|
||||
} // namespace cosmology
|
|
@ -59,7 +59,7 @@ private:
|
|||
{
|
||||
//--- general parameters ------------------------------------------
|
||||
add_class_parameter("z_max_pk", std::max(std::max(zstart_, ztarget_),199.0)); // use 1.2 as safety
|
||||
add_class_parameter("P_k_max_h/Mpc", kmax_);
|
||||
add_class_parameter("P_k_max_h/Mpc", std::max(2.0,kmax_));
|
||||
add_class_parameter("output", "dTk,vTk");
|
||||
add_class_parameter("extra metric transfer functions","yes");
|
||||
// add_class_parameter("lensing", "no");
|
||||
|
@ -84,18 +84,26 @@ private:
|
|||
// add_class_parameter("cs2_fld", 1);
|
||||
|
||||
//--- massive neutrinos -------------------------------------------
|
||||
#if 1
|
||||
#if 0
|
||||
//default off
|
||||
// add_class_parameter("Omega_ur",0.0);
|
||||
add_class_parameter("N_eff", cosmo_params_.get("Neff"));
|
||||
add_class_parameter("N_ur", cosmo_params_.get("N_ur"));
|
||||
add_class_parameter("N_ncdm", 0);
|
||||
|
||||
#else
|
||||
|
||||
add_class_parameter("N_ur", cosmo_params_.get("N_ur"));
|
||||
add_class_parameter("N_ncdm", cosmo_params_.get("N_nu_massive"));
|
||||
std::stringstream sstr;
|
||||
if( cosmo_params_.get("m_nu1") > 1e-9 ) sstr << cosmo_params_.get("m_nu1");
|
||||
if( cosmo_params_.get("m_nu2") > 1e-9 ) sstr << ", " << cosmo_params_.get("m_nu2");
|
||||
if( cosmo_params_.get("m_nu3") > 1e-9 ) sstr << ", " << cosmo_params_.get("m_nu3");
|
||||
add_class_parameter("m_ncdm", sstr.str().c_str());
|
||||
|
||||
// change above to enable
|
||||
add_class_parameter("N_ur", 0);
|
||||
add_class_parameter("N_ncdm", 1);
|
||||
add_class_parameter("m_ncdm", "0.4");
|
||||
add_class_parameter("T_ncdm", 0.71611);
|
||||
//add_class_parameter("omega_ncdm", 0.0006451439);
|
||||
//add_class_parameter("m_ncdm", "0.4");
|
||||
//add_class_parameter("T_ncdm", 0.71611);
|
||||
#endif
|
||||
|
||||
//--- cosmological parameters, primordial -------------------------
|
||||
|
|
Loading…
Reference in a new issue