CosmicFish
Reference documentation for version 1.0
Looking into future Cosmology
|
This module contains the code that computes the Fisher matrix for Cls. More...
Functions/Subroutines | |
subroutine, public | dimension_cl_covariance (P, FP, num) |
This subroutine computes the number of entrance of the cls covariance matrix The matrix of the Cls is made up by T, E, CMBWL, W1,..., Wn, B. More... | |
subroutine, public | cl_covariance_out (P, FP, cl_dim, l_min, l_max, cl_out, err) |
This subroutine that gives the covariance matrix of the Cls for the desired model. Calls internally CAMB_get results. More... | |
subroutine, public | save_cl_covariance_to_file (cl_cov_in, cl_dim, l_min, l_max, filename) |
This subroutine saves to file the input Cls covariance matrix. More... | |
subroutine, public | cls_derivative (P, FP, param_num, cl_initial, cl_derivative, cl_error, initial_step, out_dim, l_min, l_max, num_param, error_code) |
This subroutine computes the derivative of the Cls covariance matrix with a fixed stepsize finite difference formula or by the Richardson’s deferred approach to the limit. More... | |
subroutine, public | add_noise_cls_to_fiducial (P, FP, cl_fiducial, out_dim, l_min, l_max, error_code) |
This subroutine computes and adds the noise cls to the fiducial cls. On output the fiducial cls will be modified. More... | |
subroutine, public | fisher_cls (P, FP, num_param, Fisher_Matrix, outroot) |
This subroutine computes the Fisher matrix for Cls. More... | |
This module contains the code that computes the Fisher matrix for Cls.
subroutine, public fisher_calculator_cls::add_noise_cls_to_fiducial | ( | type(cambparams), intent(in) | P, |
type(cosmicfish_params), intent(in) | FP, | ||
real(dl), dimension(out_dim, out_dim, l_max-l_min+1), intent(inout) | cl_fiducial, | ||
integer, intent(in) | out_dim, | ||
integer, intent(in) | l_min, | ||
integer, intent(in) | l_max, | ||
integer, intent(out) | error_code | ||
) |
This subroutine computes and adds the noise cls to the fiducial cls. On output the fiducial cls will be modified.
[in] | p | Input CAMBparams |
[in] | fp | Input Cosmicfish params |
[in] | out_dim | Input dimension of the Cls matrix |
[in] | l_min | Input minimum value of l |
[in] | l_max | Input maximum value of l |
[in,out] | cl_fiducial | Fiducial cls covariance. On output this contains the fiducial plus the noise cls. |
[out] | error_code | Output error code. When run the code will not stop on error but report it 0 = everything was fine computation exited correctly 1 = error on input parameters 2 = error on computation 3 = error on quality of the results |
Definition at line 681 of file 008_Fisher_calculator_Cls.f90.
subroutine, public fisher_calculator_cls::cl_covariance_out | ( | type(cambparams) | P, |
type(cosmicfish_params) | FP, | ||
integer, intent(in) | cl_dim, | ||
integer, intent(in) | l_min, | ||
integer, intent(in) | l_max, | ||
real(dl), dimension(cl_dim, cl_dim, l_max-l_min+1) | cl_out, | ||
integer, intent(out) | err | ||
) |
This subroutine that gives the covariance matrix of the Cls for the desired model. Calls internally CAMB_get results.
p | Input CAMBparams | |
fp | Input Cosmicfish params | |
[in] | cl_dim | Input dimension of the Cl matrix |
[in] | l_min | Input min l computed |
[in] | l_max | Input max l computed |
cl_out | Output matrix containing the cls | |
[out] | err | Output error code: 0 = all fine 1 = error in input 2 = error in computation 3 = error in quality 4 = routine specific |
Definition at line 87 of file 008_Fisher_calculator_Cls.f90.
subroutine, public fisher_calculator_cls::cls_derivative | ( | type(cambparams), intent(in) | P, |
type(cosmicfish_params), intent(in) | FP, | ||
integer, intent(in) | param_num, | ||
real(dl), dimension(out_dim, out_dim, l_max-l_min+1), intent(in) | cl_initial, | ||
real(dl), dimension(out_dim, out_dim, l_max-l_min+1), intent(out) | cl_derivative, | ||
real(dl), dimension(out_dim, out_dim, l_max-l_min+1), intent(out) | cl_error, | ||
real(dl), intent(in) | initial_step, | ||
integer, intent(in) | out_dim, | ||
integer, intent(in) | l_min, | ||
integer, intent(in) | l_max, | ||
integer, intent(in) | num_param, | ||
integer, intent(out) | error_code | ||
) |
This subroutine computes the derivative of the Cls covariance matrix with a fixed stepsize finite difference formula or by the Richardson’s deferred approach to the limit.
[in] | p | Input CAMBparams |
[in] | fp | Input Cosmicfish params |
[in] | param_num | Input number of the parameter wrt the derivative is computed |
[in] | out_dim | Input dimension of the Cls matrix |
[in] | l_min | Input minimum value of l |
[in] | l_max | Input maximum value of l |
[in] | num_param | Input number of parameters involved in the calculation |
[in] | cl_initial | Cls covariance matrix computed for h=0 |
[out] | cl_derivative | Output returns the derivative of the Cls |
[out] | cl_error | Output returns the error on the Cls derivative at a given l |
[in] | initial_step | Input initial value of h. Does not need to be small for the adaptive algorithm. |
[out] | error_code | Output error code. When run the code will not stop on error but report it 0 = everything was fine computation exited correctly 1 = error on input parameters 2 = error on computation 3 = error on quality of the results |
Definition at line 427 of file 008_Fisher_calculator_Cls.f90.
subroutine, public fisher_calculator_cls::dimension_cl_covariance | ( | type(cambparams), intent(in) | P, |
type(cosmicfish_params), intent(in) | FP, | ||
integer, intent(out) | num | ||
) |
This subroutine computes the number of entrance of the cls covariance matrix The matrix of the Cls is made up by T, E, CMBWL, W1,..., Wn, B.
[in] | p | Input CAMBparams |
[in] | fp | Input Cosmicfish params |
[out] | num | Output size of the Cls covariance matrix |
Definition at line 56 of file 008_Fisher_calculator_Cls.f90.
subroutine, public fisher_calculator_cls::fisher_cls | ( | type(cambparams), intent(in) | P, |
type(cosmicfish_params), intent(in) | FP, | ||
integer, intent(in) | num_param, | ||
real(dl), dimension(num_param,num_param), intent(out) | Fisher_Matrix, | ||
character(len=*), intent(in), optional | outroot | ||
) |
This subroutine computes the Fisher matrix for Cls.
[in] | p | Input CAMBparams |
[in] | fp | Input Cosmicfish params |
[in] | num_param | Input dimension of the Fisher matrix |
[out] | fisher_matrix | Output Fisher matrix |
[in] | outroot | Optional input: filename that will be used to dump to file optional output if feedback is greater than 1. |
Definition at line 946 of file 008_Fisher_calculator_Cls.f90.
subroutine, public fisher_calculator_cls::save_cl_covariance_to_file | ( | real(dl), dimension(cl_dim, cl_dim, l_max-l_min+1), intent(in) | cl_cov_in, |
integer, intent(in) | cl_dim, | ||
integer, intent(in) | l_min, | ||
integer, intent(in) | l_max, | ||
character(len=*), intent(in) | filename | ||
) |
This subroutine saves to file the input Cls covariance matrix.
[in] | cl_dim | Number of Cls in a row of the matrix |
[in] | l_min | Minimum l of the Cls |
[in] | l_max | Maximum l of the Cls |
[in] | cl_cov_in | Input Cls covariance matrix |
[in] | filename | Name of the file where the Cls covariance will be saved |
Definition at line 386 of file 008_Fisher_calculator_Cls.f90.