fisher
Functions
-
void tape_waveform_gsl_subroutine(gsl_subroutine *params_packed)
-
void tape_time_gsl_subroutine(gsl_subroutine *params_packed)
-
void tape_phase_gsl_subroutine(gsl_subroutine *params_packed)
-
void prep_gsl_subroutine(gsl_subroutine *params_packed)
-
void fisher_numerical(double *frequency, int length, string generation_method, string detector, string reference_detector, double **output, int dimension, gen_params_base<double> *parameters, int order, int *amp_tapes = NULL, int *phase_tapes = NULL, double *noise = NULL, Quadrature *quadMethod = NULL)
Calculates the fisher matrix for the given arguments.
Utilizes numerical derivatives — non-skyaveraged supports up to 4th order finite difference (sky averaged supports second order only)
NOT SAFE FOR LISA YET — see hessian of phase for time derivative
- Parameters:
length – if 0, standard frequency range for the detector is used
output – double [dimension][dimension]
order – Order of the numerical derivative (2 or 4)*
amp_tapes – if speed is required, precomputed tapes can be used - assumed the user knows what they’re doing, no checks done here to make sure that the number of tapes matches the requirement by the generation_method — if using numerical derivatives or speed isn’t that important, just set to NULL
phase_tapes – if speed is required, precomputed tapes can be used - assumed the user knows what they’re doing, no checks done here to make sure that the number of tapes matches the requirement by the generation_method
quadMethod – Quadrature method
-
void calculate_fisher_elements(double *frequency, int length, int dimension, std::complex<double> **response_deriv, double **output, double *psd, std::string integration_method, double *weights, bool log10_f)
-
void calculate_fisher_elements(double **output, std::complex<double> **response_deriv, double *psd, int dimension, Quadrature *quadMethod)
-
void calculate_fisher_elements_batch(double *frequency, int length, int base_dimension, int full_dimension, std::complex<double> **response_deriv, double **output, double *psd, std::string integration_method, double *weights, bool log10_f)
Subroutine to calculate fisher elements for a subset of the fisher.
Skips elements that have dimensions (i,j) for i!=j && i>base_dim && j>base_dim
Sets non-computed elements to zero
-
void calculate_derivatives(std::complex<double> **response_deriv, double *frequencies, int length, int dimension, string detector, string reference_detector, string gen_method, gen_params_base<double> *parameters, int order)
-
void fisher_autodiff(double *frequency, int length, string generation_method, string detector, string reference_detector, double **output, int dimension, gen_params *parameters, std::string integration_method = "SIMPSONS", double *weights = NULL, bool log10_f = false, double *noise = NULL, int *amp_tapes = NULL, int *phase_tapes = NULL)
-
void fisher_autodiff_gq_internal(double *frequency, int length, string generation_method, string detector, double **output, int dimension, gen_params *parameters, double *weights, bool log_freq, int *amp_tapes = NULL, int *phase_tapes = NULL, double *noise = NULL)
-
void fisher_autodiff_interp(double *frequency, int length, string generation_method, string detector, string reference_detector, double **output, int dimension, gen_params *parameters, int downsampling_factor, int *amp_tapes = NULL, int *phase_tapes = NULL, double *noise = NULL)
-
void fisher_autodiff_batch_mod(double *frequency, int length, string generation_method, string detector, string reference_detector, double **output, int base_dimension, int full_dimension, gen_params *parameters, std::string integration_method = "SIMPSONS", double *weights = NULL, bool log10_f = false, double *noise = NULL, int *amp_tapes = NULL, int *phase_tapes = NULL)
-
void PhenomP_fisher(double *frequency, int length, gen_params *parameters, std::complex<double> **waveform_derivative, int *amp_tapes, int *phase_tapes, std::string detector)
-
void construct_waveform_derivative(double *frequency, int length, int dimension, std::complex<double> **waveform_deriv, source_parameters<double> *input_params, int *waveform_tapes)
-
void calculate_derivatives_autodiff(double *frequency, int length, int dimension, std::string generation_method, gen_params *parameters, std::complex<double> **waveform_deriv, int *waveform_tapes, std::string detector, bool autodiff_time_deriv, std::string reference_detector)
Calculates the derivatives of the detector response using automatic differentiation.
Possibly slower than the numerical derivative, but not susceptible to truncation error from finite difference
Higher dimensional fishers actually could be faster
NOTE: dimension parameter ALWAYS refers to the dimension of the fisher (ie the length of the source parameter vector), even though the derivatives are computed wrt dimension +1 or dimension + 2 — the +1(+2) are for the frequency deriv(time deriv)
autodiff_time_deriv indicates if the time derivatives for space detectors should be done, either with the full autodiff hessian (slow) or with the numerical approximation (autodiff then numerical)
-
void time_phase_corrected_derivative_autodiff_full_hess(double **dt, int length, double *frequencies, gen_params_base<double> *params, std::string generation_method, int dimension, bool correct_time)
Computes the derivative of the phase w.r.t. source parameters AS DEFINED BY FISHER FILE — hessian of the phase.
If specific derivatives need to taken, take this routine as a template and write it yourself.
The dt array has shape [dimension+1][length] (dimension + 1 for the frequency derivative, so dimension should only include the source parameters)
Computes the derivative of the phase w.r.t. source parameters AS DEFINED BY FISHER FILE — hessian of the phase
If specific derivatives need to taken, take this routine as a template and write it yourself.
The dt array has shape [dimension+1][length] (dimension + 1 for the frequency derivative, so dimension should only include the source parameters)
-
void time_phase_corrected_derivative_autodiff(double **dt, int length, double *frequencies, gen_params_base<double> *params, std::string generation_method, int dimension, bool correct_time)
Computes the derivative of the phase w.r.t. source parameters AS DEFINED BY FISHER FILE — hessian of the phase.
If specific derivatives need to taken, take this routine as a template and write it yourself.
The dt array has shape [dimension+1][length] (dimension + 1 for the frequency derivative, so dimension should only include the source parameters)
-
template<class T>
void time_phase_corrected_derivative_numerical(T **dt, int length, T *frequencies, gen_params_base<T> *params, std::string generation_method, int dimension, bool correct_time) Computes the derivative of the phase w.r.t. source parameters AS DEFINED BY FISHER FILE — hessian of the phase — numerical.
IN PROGRESS — DO NOT USE
If specific derivatives need to taken, take this routine as a template and write it yourself.
The dt array has shape [dimension+1][length] (dimension + 1 for the frequency derivative, so dimension should only include the (full) source parameters)
-
void time_phase_corrected_derivative_autodiff_numerical(double **dt, int length, double *frequencies, gen_params_base<double> *params, std::string generation_method, int dimension, bool correct_time)
Computes the derivative of the phase w.r.t. source parameters AS DEFINED BY FISHER FILE — hessian of the phase.
If specific derivatives need to taken, take this routine as a template and write it yourself.
The dt array has shape [dimension+1][length] (dimension + 1 for the frequency derivative, so dimension should only include the source parameters)
This takes the autodiff derivative wrt source parameters, and a numerical derivative for the frequency
-
std::string local_generation_method(std::string generation_method)
Utility for mapping generation method string to one accepted by the waveform_generation routines.
Certain combinations of parameters are labeled by generation method strings not under the waveform_generation routines, so a transformation is needed
-
void prep_fisher_calculation(double *parameters, bool*, double*, double*, int, gen_params_base<double> *input_params, std::string generation_method, int dim)
-
void detect_adjust_parameters(double *freq_boundaries, double *grad_freqs, int *boundary_num, gen_params_base<double> *input_params, std::string generation_method, std::string detector, int dim)
Adjust parameters for detector specific configurations (namely, LISA introduces extra transitions that needs to be accounted for)
This is kept separate to improve the modularity of the code. Waveform specific parameters are taken care of in prep_fisher_calculation and detector specific parameters are taken care of here.
-
void unpack_parameters(double *parameters, gen_params_base<double> *input_params, std::string generation_method, int dimension, bool *log_factors)
Unpacks the input gen_params object into a double array for use with the fisher routines.
-
template<class T>
void repack_parameters(T *avec_parameters, gen_params_base<T> *a_params, std::string generation_method, int dim, gen_params_base<double> *original_params) Repack the parameters from an adouble vector to a gen_params_base<adouble> object and freqeuncy.
This is one of the places where the generation-method/dimension/sky_average specific modifications should go
-
template<class T>
void repack_non_parameter_options(gen_params_base<T> *waveform_params, gen_params_base<double> *input_params, std::string gen_method) Utilitiy to transfer non-parameter options from one gen_params structure to another.
If ppE waveform ALLOCATES MEMORY — MUST be deallocated
-
template<class T>
void deallocate_non_param_options(gen_params_base<T> *waveform_params, gen_params_base<double> *input_params, std::string gen_method)
-
double calculate_integrand_autodiff_gsl_subroutine(double frequency, void *params_in)
Calculates the derivatives of the detector response using automatic differentiation — one frequency for gsl_integration.
Possibly slower than the numerical derivative, but not susceptible to truncation error from finite difference
Higher dimensional fishers actually could be faster
NOTE: dimension parameter ALWAYS refers to the dimension of the fisher (ie the length of the source parameter vector), even though the derivatives are computed wrt dimension +1 or dimension + 2 — the +1(+2) are for the frequency deriv(time deriv)
-
void fisher_autodiff_gsl_integration(double *frequency_bounds, string generation_method, string sensitivity_curve, string detector, string reference_detector, double **output, double **error, int dimension, gen_params *parameters, double abserr, double relerr)
Routine that implements GSL numerical integration to calculate the Fishers.
This can be faster than brute force calculations in fisher_autodiff, but that depends
Trade offs:
Every element is calculated independently, so no information is retained between elements. In the brute force calculation, there is information reused.
However, time can be saved by spending less time on trivial elements (identically 0 elements, etc) and better spent on complicated elements
Does not have a direct interpretation in terms of integration time, as the scheme is adaptative. Sampling frequency and integration time are `as good as they need to be’ to calculate the fisher
Implements (GSL_INTEG_GAUSS15)
- Parameters:
frequency_bounds – Bounds of integration in fourier space
generation_method – Method of waveform generation
sensitivity_curve – Sensitivity curve to be used for the PSD — MUST BE ANALYTIC
detector – Detector to use for the response function
reference_detector – Detector to use for the response function
output – [out] Output Fisher — must be preallocated — shape [dimension][dimension]
error – [out] Estimated error, as specified by GSL’s integration — must be preallocated — shape [dimension][dimension]
dimension – Dimension of the Fisher
parameters – Generation parameters specifying source parameters and waveform options
abserr – Target absolute error (0 if this should be ignored — ONE type of error must be specified)
relerr – Target relative error (0 if this should be ignored — ONE type of error must be specified)
-
void fisher_autodiff_gsl_integration(double *frequency_bounds, string generation_method, string sensitivity_curve, string detector, string reference_detector, double **output, double **error, int dimension, gen_params *parameters, double abserr, double relerr, std::string error_log, bool logerr)
Routine that implements GSL numerical integration to calculate the Fishers.
This can be faster than brute force calculations in fisher_autodiff, but that depends
Trade offs:
Every element is calculated independently, so no information is retained between elements. In the brute force calculation, there is information reused.
However, time can be saved by spending less time on trivial elements (identically 0 elements, etc) and better spent on complicated elements
Does not have a direct interpretation in terms of integration time, as the scheme is adaptative. Sampling frequency and integration time are `as good as they need to be’ to calculate the fisher
Implements (GSL_INTEG_GAUSS15)
Now includes option to log error instead of ending program for certain types of errors
- Parameters:
frequency_bounds – Bounds of integration in fourier space
generation_method – Method of waveform generation
sensitivity_curve – Sensitivity curve to be used for the PSD — MUST BE ANALYTIC
detector – Detector to use for the response function
reference_detector – Detector to use for the response function
output – [out] Output Fisher — must be preallocated — shape [dimension][dimension]
error – [out] Estimated error, as specified by GSL’s integration — must be preallocated — shape [dimension][dimension]
dimension – Dimension of the Fisher
parameters – Generation parameters specifying source parameters and waveform options
abserr – Target absolute error (0 if this should be ignored — ONE type of error must be specified)
relerr – Target relative error (0 if this should be ignored — ONE type of error must be specified)
error_log – File to write non-critical error codes to (roundoff error)
logerr – Whether or not to end program with certain error codes, or to log them and continue
-
void fisher_autodiff_gsl_integration_batch_mod(double *frequency_bounds, string generation_method, string sensitivity_curve, string detector, string reference_detector, double **output, double **error, int base_dimension, int full_dimension, gen_params *parameters, double abserr, double relerr)
Routine that implements GSL numerical integration to calculate the Fishers — batch modifications version.
Calculates Fisher for multiple modifications at a time, neglecting covariance between modifications (set to 0 in Fisher)
Modifications MUST BE evaluated at 0 for this routine to calculate correct results
This can be faster than brute force calculations in fisher_autodiff, but that depends
Trade offs:
Every element is calculated independently, so no information is retained between elements. In the brute force calculation, there is information reused.
However, time can be saved by spending less time on trivial elements (identically 0 elements, etc) and better spent on complicated elements
Does not have a direct interpretation in terms of integration time, as the scheme is adaptative. Sampling frequency and integration time are `as good as they need to be’ to calculate the fisher
Implements (GSL_INTEG_GAUSS15)
- Parameters:
frequency_bounds – Bounds of integration in fourier space
generation_method – Method of waveform generation
sensitivity_curve – Sensitivity curve to be used for the PSD — MUST BE ANALYTIC
detector – Detector to use for the response function
reference_detector – Detector to use for the response function
output – [out] Output Fisher — must be preallocated — shape [full_dimension][full_dimension]
error – [out] Estimated error, as specified by GSL’s integration — must be preallocated — shape [full_dimension][full_dimension]
base_dimension – Dimension of base model (ie GR dimension)
full_dimension – Full dimension (GR dimension + Nmod)
parameters – Generation parameters specifying source parameters and waveform options
abserr – Target absolute error (0 if this should be ignored — ONE type of error must be specified)
relerr – Target relative error (0 if this should be ignored — ONE type of error must be specified)
-
void fisher_autodiff_gsl_integration_batch_mod(double *frequency_bounds, string generation_method, string sensitivity_curve, string detector, string reference_detector, double **output, double **error, int base_dimension, int full_dimension, gen_params *parameters, double abserr, double relerr, std::string error_log, bool logerr)
Routine that implements GSL numerical integration to calculate the Fishers — batch modifications version.
Calculates Fisher for multiple modifications at a time, neglecting covariance between modifications (set to 0 in Fisher)
Modifications MUST BE evaluated at 0 for this routine to calculate correct results
This can be faster than brute force calculations in fisher_autodiff, but that depends
Trade offs:
Every element is calculated independently, so no information is retained between elements. In the brute force calculation, there is information reused.
However, time can be saved by spending less time on trivial elements (identically 0 elements, etc) and better spent on complicated elements
Does not have a direct interpretation in terms of integration time, as the scheme is adaptative. Sampling frequency and integration time are `as good as they need to be’ to calculate the fisher
Implements (GSL_INTEG_GAUSS15)
Now includes option to log error instead of ending program for certain types of errors
- Parameters:
frequency_bounds – Bounds of integration in fourier space
generation_method – Method of waveform generation
sensitivity_curve – Sensitivity curve to be used for the PSD — MUST BE ANALYTIC
detector – Detector to use for the response function
reference_detector – Detector to use for the response function
output – [out] Output Fisher — must be preallocated — shape [full_dimension][full_dimension]
error – [out] Estimated error, as specified by GSL’s integration — must be preallocated — shape [full_dimension][full_dimension]
base_dimension – Dimension of base model (ie GR dimension)
full_dimension – Full dimension (GR dimension + Nmod)
parameters – Generation parameters specifying source parameters and waveform options
abserr – Target absolute error (0 if this should be ignored — ONE type of error must be specified)
relerr – Target relative error (0 if this should be ignored — ONE type of error must be specified)
error_log – File to write non-critical error codes to (roundoff error)
logerr – Whether or not to end program with certain error codes, or to log them and continue
-
void ppE_theory_fisher_transformation(std::string original_method, std::string new_method, int dimension, gen_params_base<double> *param, double **old_fisher, double **new_fisher)
Transform a generic ppE Fisher matrix to a theory specific Fisher matrix.
See ppE_utilities.cpp for a list of supported theories.
The beta value in the input parameter structure must be the value of the theory-specific parameter, ie \alpha_dCS^2 and not \beta_2PN
-
void ppE_theory_transformation_jac(std::string original_method, std::string new_method, int dimension, gen_params_base<double> *param, double **jac)
Calculates the jacobian matrix going from ppE to a specific theory.
d \theta_i / d \theta_j where \theta_i is the ppE parameterrs, and \theta_j are the theory specific parameters
-
void ppE_theory_fisher_transformation_calculate_derivatives(std::string original_method, std::string new_method, int dimension, int base_dim, gen_params_base<double> *param, double **derivatives)
-
void ppE_theory_covariance_transformation(std::string original_method, std::string new_method, int dimension, gen_params_base<double> *param, double **old_cov, double **new_cov)
Transform a generic ppE covariance matrix to a theory specific covariance matrix.
See ppE_utilities.cpp for a list of supported theories.
The covariance matrix is defined as (\Gamma)^-1 where, \Gamma is the Fisher matrix
-
struct gsl_subroutine
- #include <fisher.h>
Public Members
-
string detector
-
string reference_detector
-
string generation_method
-
string sensitivity_curve
-
gen_params *gen_params_in
-
int dim
-
int id1
-
int id2
-
int *waveform_tapes
-
int *time_tapes
-
int *phase_tapes
-
double *freq_boundaries
-
double *grad_freqs
-
int boundary_num
-
bool *log_factors
-
string detector
-
namespace std