![]() |
CosmoBolognaLib
Free Software C++/Python libraries for cosmological calculations
|
The namespace of the CAMB wrappers More...
Functions | |
| void * | GetCAMBparams () |
| wrapper of the function GetCAMBparams in CAMB/CAMBinterface.f90, creating an instance of the CAMBparams type contained in External/CAMB/fortran/model.f90 More... | |
| void * | GetCAMBdata () |
| void | SetCAMBparams (void *params, const double ombh2, const double omch2, const double omnuh2, const double massless_nu, const int massive_nu, const int neutrino_hierarchy, const double omk, const double H0, const int dark_energy_model, const double w, const double wa, const double tau, const double cs2_lam, const double T_cmb, const double helium_fraction) |
| wrapper of the subroutine SetCAMBparams in CAMB/CAMBinterface.f90, setting CAMBparams parameters More... | |
| void | SetCAMBPk (void *params, const double redshift, const double ns, const double As, const double pivot_scalar, const bool accurate_massive_nu, const double kmax, const bool nonlinear) |
| wrapper of the subroutine SetCAMBPk in CAMB/CAMBinterface.f90, setting the power spectrum parameters More... | |
| void | GetCAMBresults (void *params, void *data) |
| wrapper of the subroutine GetCAMBresults in CAMB/CAMBinterface.f90, calling the subroutine CAMB_GetResults contained in External/CAMB/fortran/camb.f90 More... | |
| void | GetCAMBPk (void *data, double *Pk, const double mink, const double dlnkh, const int npoints) |
| wrapper of the subroutine GetCAMBPk in CAMB/CAMBinterface.f90, calling the subroutine Transfer_GetMatterPowerD contained in External/CAMB/fortran/results.f90 More... | |
| void * | ReleaseCAMBparams (void *params) |
| wrapper of the subroutine ReleaseCAMBparams in CAMB/CAMBinterface.f90, deallocating CAMBparams instance More... | |
| void * | ReleaseCAMBdata (void *data) |
| wrapper of the subroutine ReleaseCAMBdata in CAMB/CAMBinterface.f90, deallocating CAMBdata instance More... | |
| std::vector< double > | Pk_CAMB (const bool nonlinear, const double redshift, const double kmin, const double kmax, const int npoints, const double ombh2, const double omch2, const double omnuh2, const double massless_nu, const int massive_nu, const double omk, const double H0, const double ns, const double As, const double pivot_scalar=0.05, const double w=-1., const double wa=0., const double tau=2.1e-9, const bool accurate_massive_nu=false, const int neutrino_hierarchy=3, const int dark_energy_model=1, const double cs2_lam=1., const double T_cmb=2.7255, const double helium_fraction=0.24) |
| Get the matter power spectrum from CAMB. More... | |
The namespace of the CAMB wrappers
The camb namespace contains all the wrapper functions of the CAMB routines
| void* cbl::wrapper::camb::GetCAMBdata | ( | ) |
wrapper of the function GetCAMBdata in CAMB/CAMBinterface.f90, creating an instance of the CAMBdata type contained in External/CAMB/fortran/results.f90
| void* cbl::wrapper::camb::GetCAMBparams | ( | ) |
wrapper of the function GetCAMBparams in CAMB/CAMBinterface.f90, creating an instance of the CAMBparams type contained in External/CAMB/fortran/model.f90
| void cbl::wrapper::camb::GetCAMBPk | ( | void * | data, |
| double * | Pk, | ||
| const double | mink, | ||
| const double | dlnkh, | ||
| const int | npoints | ||
| ) |
wrapper of the subroutine GetCAMBPk in CAMB/CAMBinterface.f90, calling the subroutine Transfer_GetMatterPowerD contained in External/CAMB/fortran/results.f90
| data | pointer to CAMBdata instance |
| Pk | output power spectrum |
| mink | minimum \(k\) |
| dlnkh | log space between each \(k\) |
| npoints | number of points where \(P(k)\) is computed |
| void cbl::wrapper::camb::GetCAMBresults | ( | void * | params, |
| void * | data | ||
| ) |
wrapper of the subroutine GetCAMBresults in CAMB/CAMBinterface.f90, calling the subroutine CAMB_GetResults contained in External/CAMB/fortran/camb.f90
| params | pointer to CAMBparams instance |
| data | pointer to CAMBdata instance |
| std::vector< double > cbl::wrapper::camb::Pk_CAMB | ( | const bool | nonlinear, |
| const double | redshift, | ||
| const double | kmin, | ||
| const double | kmax, | ||
| const int | npoints, | ||
| const double | ombh2, | ||
| const double | omch2, | ||
| const double | omnuh2, | ||
| const double | massless_nu, | ||
| const int | massive_nu, | ||
| const double | omk, | ||
| const double | H0, | ||
| const double | ns, | ||
| const double | As, | ||
| const double | pivot_scalar = 0.05, |
||
| const double | w = -1., |
||
| const double | wa = 0., |
||
| const double | tau = 2.1e-9, |
||
| const bool | accurate_massive_nu = false, |
||
| const int | neutrino_hierarchy = 3, |
||
| const int | dark_energy_model = 1, |
||
| const double | cs2_lam = 1., |
||
| const double | T_cmb = 2.7255, |
||
| const double | helium_fraction = 0.24 |
||
| ) |
Get the matter power spectrum from CAMB.
| nonlinear | if true, \(P(k)\) is nonlinear |
| redshift | the redshift |
| kmin | minimum \(k\) |
| kmax | maximum \(k\) |
| npoints | number of points where \(P(k)\) is computed |
| ombh2 | \(\Omega_{\rm b}h^2\) |
| omch2 | \(\Omega_{\rm c}h^2\) |
| omnuh2 | \(\Omega_{\nu}h^2\) |
| massless_nu | contribution of massless neutrinos to \(N_{\rm eff}\) |
| massive_nu | number of massive neutrinos |
| omk | \(\Omega_{\rm k}\) |
| H0 | \(H_0\) |
| ns | \(n_s\) |
| As | \(A_s\) |
| pivot_scalar | the pivot scalar |
| w | \(w_0\) |
| wa | \(w_{\rm a}\) |
| tau | \(\tau\) |
| accurate_massive_nu | whether you care about accuracy of the neutrino transfers themselves |
| neutrino_hierarchy | 1 \(\rightarrow\) normal hierarchy; 2 \(\rightarrow\) inverted hierarchy; 3 \(\rightarrow\) degenerate hierarchy |
| dark_energy_model | 0 \(\rightarrow\) fluid; 1 \(\rightarrow\) PPF; 2 \(\rightarrow\) axion effective fluid; 3 \(\rightarrow\) early quintessence |
| cs2_lam | constant comoving sound speed of the dark energy (1=quintessence), for fluid model |
| T_cmb | CMB temperature |
| helium_fraction | helium fraction |
| void* cbl::wrapper::camb::ReleaseCAMBdata | ( | void * | data | ) |
wrapper of the subroutine ReleaseCAMBdata in CAMB/CAMBinterface.f90, deallocating CAMBdata instance
| data | pointer to CAMBdata instance |
| void* cbl::wrapper::camb::ReleaseCAMBparams | ( | void * | params | ) |
wrapper of the subroutine ReleaseCAMBparams in CAMB/CAMBinterface.f90, deallocating CAMBparams instance
| params | pointer to CAMBparams instance |
| void cbl::wrapper::camb::SetCAMBparams | ( | void * | params, |
| const double | ombh2, | ||
| const double | omch2, | ||
| const double | omnuh2, | ||
| const double | massless_nu, | ||
| const int | massive_nu, | ||
| const int | neutrino_hierarchy, | ||
| const double | omk, | ||
| const double | H0, | ||
| const int | dark_energy_model, | ||
| const double | w, | ||
| const double | wa, | ||
| const double | tau, | ||
| const double | cs2_lam, | ||
| const double | T_cmb, | ||
| const double | helium_fraction | ||
| ) |
wrapper of the subroutine SetCAMBparams in CAMB/CAMBinterface.f90, setting CAMBparams parameters
| params | pointer to CAMBparams instance |
| ombh2 | \(\Omega_{\rm b}h^2\) |
| omch2 | \(\Omega_{\rm c}h^2\) |
| omnuh2 | \(\Omega_{\nu}h^2\) |
| massless_nu | contribution of massless neutrinos to \(N_{\rm eff}\) |
| massive_nu | number of massive neutrinos |
| neutrino_hierarchy | 1 \(\rightarrow\) normal hierarchy; 2 \(\rightarrow\) inverted hierarchy; 3 \(\rightarrow\) degenerate hierarchy |
| omk | \(\Omega_{\rm k}\) |
| H0 | \(H_0\) |
| dark_energy_model | 0 \(\rightarrow\) fluid; 1 \(\rightarrow\) PPF; 2 \(\rightarrow\) axion effective fluid; 3 \(\rightarrow\) early quintessence |
| w | \(w_0\) |
| wa | \(w_{\rm a}\) |
| tau | \(\tau\) |
| cs2_lam | constant comoving sound speed of the dark energy (1=quintessence), for fluid model |
| T_cmb | CMB temperature |
| helium_fraction | helium fraction |
| void cbl::wrapper::camb::SetCAMBPk | ( | void * | params, |
| const double | redshift, | ||
| const double | ns, | ||
| const double | As, | ||
| const double | pivot_scalar, | ||
| const bool | accurate_massive_nu, | ||
| const double | kmax, | ||
| const bool | nonlinear | ||
| ) |
wrapper of the subroutine SetCAMBPk in CAMB/CAMBinterface.f90, setting the power spectrum parameters
| params | pointer to CAMBparams instance |
| redshift | the redshift |
| ns | \(n_s\) |
| As | \(A_s\) |
| pivot_scalar | the pivot scalar |
| accurate_massive_nu | whether you care about accuracy of the neutrino transfers themselves |
| kmax | maximum \(k\) |
| nonlinear | if true, \(P(k)\) is nonlinear |