CLASS MANUAL
|
#include "input.h"
Functions | |
int | input_init_from_arguments (int argc, char **argv, struct precision *ppr, struct background *pba, struct thermo *pth, struct perturbs *ppt, struct transfers *ptr, struct primordial *ppm, struct spectra *psp, struct nonlinear *pnl, struct lensing *ple, struct output *pop, ErrorMsg errmsg) |
int | input_init (struct file_content *pfc, struct precision *ppr, struct background *pba, struct thermo *pth, struct perturbs *ppt, struct transfers *ptr, struct primordial *ppm, struct spectra *psp, struct nonlinear *pnl, struct lensing *ple, struct output *pop, ErrorMsg errmsg) |
int | input_read_parameters (struct file_content *pfc, struct precision *ppr, struct background *pba, struct thermo *pth, struct perturbs *ppt, struct transfers *ptr, struct primordial *ppm, struct spectra *psp, struct nonlinear *pnl, struct lensing *ple, struct output *pop, ErrorMsg errmsg) |
int | input_default_params (struct background *pba, struct thermo *pth, struct perturbs *ppt, struct transfers *ptr, struct primordial *ppm, struct spectra *psp, struct nonlinear *pnl, struct lensing *ple, struct output *pop) |
int | input_default_precision (struct precision *ppr) |
int | get_machine_precision (double *smallest_allowed_variation) |
int | class_fzero_ridder (int(*func)(double x, void *param, double *y, ErrorMsg error_message), double x1, double x2, double xtol, void *param, double *Fx1, double *Fx2, double *xzero, int *fevals, ErrorMsg error_message) |
int | input_try_unknown_parameters (double *unknown_parameter, int unknown_parameters_size, void *voidpfzw, double *output, ErrorMsg errmsg) |
int | input_get_guess (double *xguess, double *dxdy, struct fzerofun_workspace *pfzw, ErrorMsg errmsg) |
int | input_find_root (double *xzero, int *fevals, struct fzerofun_workspace *pfzw, ErrorMsg errmsg) |
Documented input module.
Julien Lesgourgues, 27.08.2010
int input_init_from_arguments | ( | int | argc, |
char ** | argv, | ||
struct precision * | ppr, | ||
struct background * | pba, | ||
struct thermo * | pth, | ||
struct perturbs * | ppt, | ||
struct transfers * | ptr, | ||
struct primordial * | ppm, | ||
struct spectra * | psp, | ||
struct nonlinear * | pnl, | ||
struct lensing * | ple, | ||
struct output * | pop, | ||
ErrorMsg | errmsg | ||
) |
Use this routine to extract initial parameters from files 'xxx.ini' and/or 'xxx.pre'. They can be the arguments of the main() routine.
If class is embedded into another code, you will probably prefer to call directly input_init() in order to pass input parameters through a 'file_content' structure.
Summary:
int input_init | ( | struct file_content * | pfc, |
struct precision * | ppr, | ||
struct background * | pba, | ||
struct thermo * | pth, | ||
struct perturbs * | ppt, | ||
struct transfers * | ptr, | ||
struct primordial * | ppm, | ||
struct spectra * | psp, | ||
struct nonlinear * | pnl, | ||
struct lensing * | ple, | ||
struct output * | pop, | ||
ErrorMsg | errmsg | ||
) |
Initialize each parameter, first to its default values, and then from what can be interpreted from the values passed in the input 'file_content' structure. If its size is null, all parameters keep their default values.
These two arrays must contain the strings of names to be searched for and the corresponding new parameter
int input_read_parameters | ( | struct file_content * | pfc, |
struct precision * | ppr, | ||
struct background * | pba, | ||
struct thermo * | pth, | ||
struct perturbs * | ppt, | ||
struct transfers * | ptr, | ||
struct primordial * | ppm, | ||
struct spectra * | psp, | ||
struct nonlinear * | pnl, | ||
struct lensing * | ple, | ||
struct output * | pop, | ||
ErrorMsg | errmsg | ||
) |
Summary:
Knowing the gauge from the very beginning is useful (even if this could be a run not requiring perturbations at all: even in that case, knowing the gauge is important e.g. for fixing the sampling in momentum space for non-cold dark matter)
(a) background parameters
(b) assign values to thermodynamics cosmological parameters
(c) define which perturbations and sources should be computed, and down to which scale
(d) define the primordial spectrum
(e) parameters for final spectra
(f) parameter related to the non-linear spectra computation
(g) amount of information sent to standard output (none if all set to zero)
(h) all precision parameters
(i) Write values in file
int input_default_params | ( | struct background * | pba, |
struct thermo * | pth, | ||
struct perturbs * | ppt, | ||
struct transfers * | ptr, | ||
struct primordial * | ppm, | ||
struct spectra * | psp, | ||
struct nonlinear * | pnl, | ||
struct lensing * | ple, | ||
struct output * | pop | ||
) |
All default parameter values (for input parameters)
pba | Input: pointer to background structure |
pth | Input: pointer to thermodynamics structure |
ppt | Input: pointer to perturbation structure |
ptr | Input: pointer to transfer structure |
ppm | Input: pointer to primordial structure |
psp | Input: pointer to spectra structure |
pnl | Input: pointer to nonlinear structure |
ple | Input: pointer to lensing structure |
pop | Input: pointer to output structure |
Define all default parameter values (for input parameters) for each structure:
int input_default_precision | ( | struct precision * | ppr | ) |
Initialize the precision parameter structure.
All precision parameters used in the other modules are listed here and assigned here a default value.
ppr | Input/Output: a precision_params structure pointer |
Initialize presicion parameters for different structures:
int get_machine_precision | ( | double * | smallest_allowed_variation | ) |
Automatically computes the machine precision.
smallest_allowed_variation | a pointer to the smallest allowed variation |
Returns the smallest allowed variation (minimum epsilon * TOLVAR)
int class_fzero_ridder | ( | int(*)(double x, void *param, double *y, ErrorMsg error_message) | func, |
double | x1, | ||
double | x2, | ||
double | xtol, | ||
void * | param, | ||
double * | Fx1, | ||
double * | Fx2, | ||
double * | xzero, | ||
int * | fevals, | ||
ErrorMsg | error_message | ||
) |
Using Ridders' method, return the root of a function func known to lie between x1 and x2. The root, returned as zriddr, will be found to an approximate accuracy xtol.
int input_try_unknown_parameters | ( | double * | unknown_parameter, |
int | unknown_parameters_size, | ||
void * | voidpfzw, | ||
double * | output, | ||
ErrorMsg | errmsg | ||
) |
Summary:
int input_get_guess | ( | double * | xguess, |
double * | dxdy, | ||
struct fzerofun_workspace * | pfzw, | ||
ErrorMsg | errmsg | ||
) |
Summary:
xguess[index_guess] = 1.77835*pow(ba.Omega0_scf,-2./7.); dxdy[index_guess] = -0.5081*pow(ba.Omega0_scf,-9./7.)
;int input_find_root | ( | double * | xzero, |
int * | fevals, | ||
struct fzerofun_workspace * | pfzw, | ||
ErrorMsg | errmsg | ||
) |
Summary: