!-*-f90-*- ! ! Interfaces: Statistics ! function gsl_stats_mean (data, stride, n) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double) :: gsl_stats_mean end function gsl_stats_mean function gsl_stats_variance (data, stride, n) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double) :: gsl_stats_variance end function gsl_stats_variance function gsl_stats_variance_m (data, stride, n, mean) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double), value :: mean real(c_double) :: gsl_stats_variance_m end function gsl_stats_variance_m function gsl_stats_sd (data, stride, n) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double) :: gsl_stats_sd end function gsl_stats_sd function gsl_stats_sd_m (data, stride, n, mean) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double), value :: mean real(c_double) :: gsl_stats_sd_m end function gsl_stats_sd_m function gsl_stats_variance_with_fixed_mean (data, stride, n, mean) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double), value :: mean real(c_double) :: gsl_stats_variance_with_fixed_mean end function gsl_stats_variance_with_fixed_mean function gsl_stats_sd_with_fixed_mean (data, stride, n, mean) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double), value :: mean real(c_double) :: gsl_stats_sd_with_fixed_mean end function gsl_stats_sd_with_fixed_mean function gsl_stats_absdev (data, stride, n) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double) :: gsl_stats_absdev end function gsl_stats_absdev function gsl_stats_absdev_m (data, stride, n, mean) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double), value :: mean real(c_double) :: gsl_stats_absdev_m end function gsl_stats_absdev_m function gsl_stats_skew (data, stride, n) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double) :: gsl_stats_skew end function gsl_stats_skew function gsl_stats_skew_m_sd (data, stride, n, mean, sd) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double), value :: mean, sd real(c_double) :: gsl_stats_skew_m_sd end function gsl_stats_skew_m_sd function gsl_stats_kurtosis (data, stride, n) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double) :: gsl_stats_kurtosis end function gsl_stats_kurtosis function gsl_stats_kurtosis_m_sd (data, stride, n, mean, sd) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double), value :: mean, sd real(c_double) :: gsl_stats_kurtosis_m_sd end function gsl_stats_kurtosis_m_sd function gsl_stats_lag1_autocorrelation (data, stride, n) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double) :: gsl_stats_lag1_autocorrelation end function gsl_stats_lag1_autocorrelation function gsl_stats_lag1_autocorrelation_m (data, stride, n, mean) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double), value :: mean real(c_double) :: gsl_stats_lag1_autocorrelation_m end function gsl_stats_lag1_autocorrelation_m function gsl_stats_covariance(data1, stride1, data2, stride2, n) bind(c) import type(c_ptr), value :: data1, data2 integer(c_size_t), value :: stride1, stride2, n real(c_double) :: gsl_stats_covariance end function gsl_stats_covariance function gsl_stats_covariance_m(data1, stride1, data2, stride2, n, & mean1, mean2) bind(c) import type(c_ptr), value :: data1, data2 integer(c_size_t), value :: stride1, stride2, n real(c_double), value :: mean1, mean2 real(c_double) :: gsl_stats_covariance_m end function gsl_stats_covariance_m function gsl_stats_correlation(data1, stride1, data2, stride2, n) bind(c) import type(c_ptr), value :: data1, data2 integer(c_size_t), value :: stride1, stride2, n real(c_double) :: gsl_stats_correlation end function gsl_stats_correlation function gsl_stats_spearman(data1, stride1, data2, stride2, n, work) bind(c) import type(c_ptr), value :: data1, data2 type(c_ptr), value :: work integer(c_size_t), value :: stride1, stride2, n real(c_double) :: gsl_stats_spearman end function gsl_stats_spearman function gsl_stats_wmean(w, wstride, data, stride, n) bind(c) import type(c_ptr), value :: w, data integer(c_size_t), value :: wstride, stride, n real(c_double) :: gsl_stats_wmean end function gsl_stats_wmean function gsl_stats_wvariance(w, wstride, data, stride, n) bind(c) import type(c_ptr), value :: w, data integer(c_size_t), value :: wstride, stride, n real(c_double) :: gsl_stats_wvariance end function gsl_stats_wvariance function gsl_stats_wvariance_m(w, wstride, data, stride, n, mean) bind(c) import type(c_ptr), value :: w, data integer(c_size_t), value :: wstride, stride, n real(c_double), value :: mean real(c_double) :: gsl_stats_wvariance_m end function gsl_stats_wvariance_m function gsl_stats_wsd(w, wstride, data, stride, n) bind(c) import type(c_ptr), value :: w, data integer(c_size_t), value :: wstride, stride, n real(c_double) :: gsl_stats_wsd end function gsl_stats_wsd function gsl_stats_wsd_m(w, wstride, data, stride, n, mean) bind(c) import type(c_ptr), value :: w, data integer(c_size_t), value :: wstride, stride, n real(c_double), value :: mean real(c_double) :: gsl_stats_wsd_m end function gsl_stats_wsd_m function gsl_stats_wvariance_with_fixed_mean(w, wstride, data, stride, n, mean) bind(c) import type(c_ptr), value :: w, data integer(c_size_t), value :: wstride, stride, n real(c_double), value :: mean real(c_double) :: gsl_stats_wvariance_with_fixed_mean end function gsl_stats_wvariance_with_fixed_mean function gsl_stats_wsd_with_fixed_mean(w, wstride, data, stride, n, mean) bind(c) import type(c_ptr), value :: w, data integer(c_size_t), value :: wstride, stride, n real(c_double), value :: mean real(c_double) :: gsl_stats_wsd_with_fixed_mean end function gsl_stats_wsd_with_fixed_mean function gsl_stats_wabsdev(w, wstride, data, stride, n) bind(c) import type(c_ptr), value :: w, data integer(c_size_t), value :: wstride, stride, n real(c_double) :: gsl_stats_wabsdev end function gsl_stats_wabsdev function gsl_stats_wabsdev_m(w, wstride, data, stride, n, mean) bind(c) import type(c_ptr), value :: w, data integer(c_size_t), value :: wstride, stride, n real(c_double), value :: mean real(c_double) :: gsl_stats_wabsdev_m end function gsl_stats_wabsdev_m function gsl_stats_wskew(w, wstride, data, stride, n) bind(c) import type(c_ptr), value :: w, data integer(c_size_t), value :: wstride, stride, n real(c_double) :: gsl_stats_wskew end function gsl_stats_wskew function gsl_stats_wskew_m_sd(w, wstride, data, stride, n, mean, sd) bind(c) import type(c_ptr), value :: w, data integer(c_size_t), value :: wstride, stride, n real(c_double), value :: mean, sd real(c_double) :: gsl_stats_wskew_m_sd end function gsl_stats_wskew_m_sd function gsl_stats_wkurtosis(w, wstride, data, stride, n) bind(c) import type(c_ptr), value :: w, data integer(c_size_t), value :: wstride, stride, n real(c_double) :: gsl_stats_wkurtosis end function gsl_stats_wkurtosis function gsl_stats_wkurtosis_m_sd(w, wstride, data, stride, n, mean, sd) bind(c) import type(c_ptr), value :: w, data integer(c_size_t), value :: wstride, stride, n real(c_double), value :: mean, sd real(c_double) :: gsl_stats_wkurtosis_m_sd end function gsl_stats_wkurtosis_m_sd function gsl_stats_max(data, stride, n) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double) :: gsl_stats_max end function gsl_stats_max function gsl_stats_min(data, stride, n) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double) :: gsl_stats_min end function gsl_stats_min subroutine gsl_stats_minmax(min, max, data, stride, n) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double) :: min, max end subroutine gsl_stats_minmax function gsl_stats_max_index(data, stride, n) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n integer(c_size_t) :: gsl_stats_max_index end function gsl_stats_max_index function gsl_stats_min_index(data, stride, n) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n integer(c_size_t) :: gsl_stats_min_index end function gsl_stats_min_index subroutine gsl_stats_minmax_index(min_index, max_index, data, stride, n) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n integer(c_size_t) :: min_index, max_index end subroutine gsl_stats_minmax_index function gsl_stats_median_from_sorted_data(data, stride, n) bind(c) import type(c_ptr), value :: data integer(c_size_t), value :: stride, n real(c_double) :: gsl_stats_median_from_sorted_data end function gsl_stats_median_from_sorted_data function fgsl_stats_median(data, stride, n) bind(c, name='gsl_stats_median') import :: fgsl_size_t, fgsl_double real(fgsl_double), dimension(*), intent(in) :: data integer(fgsl_size_t), value :: stride, n real(fgsl_double) :: fgsl_stats_median end function fgsl_stats_median function gsl_stats_quantile_from_sorted_data(data, stride, n, f) bind(c) import type(c_ptr), value :: data real(c_double), value :: f integer(c_size_t), value :: stride, n real(c_double) :: gsl_stats_quantile_from_sorted_data end function gsl_stats_quantile_from_sorted_data subroutine fgsl_stats_select(data, stride, n, k) bind(c, name='gsl_stats_select') import :: fgsl_size_t, fgsl_double real(fgsl_double), dimension(*), intent(inout) :: data integer(fgsl_size_t), value :: stride, n, k end subroutine fgsl_stats_select function fgsl_stats_trmean_from_sorted_data(alpha, data, stride, n) & bind(c, name='gsl_stats_trmean_from_sorted_data') import :: fgsl_size_t, fgsl_double real(fgsl_double), value :: alpha real(fgsl_double), dimension(*), intent(in) :: data integer(fgsl_size_t), value :: stride, n real(fgsl_double) :: fgsl_stats_trmean_from_sorted_data end function fgsl_stats_trmean_from_sorted_data function fgsl_stats_gastwirth_from_sorted_data(data, stride, n) & bind(c, name='gsl_stats_gastwirth_from_sorted_data') import :: fgsl_size_t, fgsl_double real(fgsl_double), dimension(*), intent(in) :: data integer(fgsl_size_t), value :: stride, n real(fgsl_double) :: fgsl_stats_gastwirth_from_sorted_data end function fgsl_stats_gastwirth_from_sorted_data function fgsl_stats_mad0(data, stride, n, work) bind(c, name='gsl_stats_mad0') import :: fgsl_size_t, fgsl_double real(fgsl_double), dimension(*), intent(in) :: data integer(fgsl_size_t), value :: stride, n real(fgsl_double), dimension(*), intent(inout) :: work real(fgsl_double) :: fgsl_stats_mad0 end function fgsl_stats_mad0 function fgsl_stats_mad(data, stride, n, work) bind(c, name='gsl_stats_mad') import :: fgsl_size_t, fgsl_double real(fgsl_double), dimension(*), intent(in) :: data integer(fgsl_size_t), value :: stride, n real(fgsl_double), dimension(*), intent(inout) :: work real(fgsl_double) :: fgsl_stats_mad end function fgsl_stats_mad function fgsl_stats_sn0_from_sorted_data(data, stride, n, work) & bind(c, name='gsl_stats_Sn0_from_sorted_data') import :: fgsl_size_t, fgsl_double real(fgsl_double), dimension(*), intent(in) :: data integer(fgsl_size_t), value :: stride, n real(fgsl_double), dimension(*), intent(inout) :: work real(fgsl_double) :: fgsl_stats_sn0_from_sorted_data end function fgsl_stats_sn0_from_sorted_data function fgsl_stats_sn_from_sorted_data(data, stride, n, work) & bind(c, name='gsl_stats_Sn_from_sorted_data') import :: fgsl_size_t, fgsl_double real(fgsl_double), dimension(*), intent(in) :: data integer(fgsl_size_t), value :: stride, n real(fgsl_double), dimension(*), intent(inout) :: work real(fgsl_double) :: fgsl_stats_sn_from_sorted_data end function fgsl_stats_sn_from_sorted_data function fgsl_stats_qn0_from_sorted_data(data, stride, n, work, work_int) & bind(c, name='gsl_stats_Qn0_from_sorted_data') import :: fgsl_size_t, fgsl_double, fgsl_int real(fgsl_double), dimension(*), intent(in) :: data integer(fgsl_size_t), value :: stride, n real(fgsl_double), dimension(*), intent(inout) :: work integer(fgsl_int), dimension(*), intent(inout) :: work_int real(fgsl_double) :: fgsl_stats_qn0_from_sorted_data end function fgsl_stats_qn0_from_sorted_data function fgsl_stats_qn_from_sorted_data(data, stride, n, work, work_int) & bind(c, name='gsl_stats_Qn_from_sorted_data') import :: fgsl_size_t, fgsl_double, fgsl_int real(fgsl_double), dimension(*), intent(in) :: data integer(fgsl_size_t), value :: stride, n real(fgsl_double), dimension(*), intent(inout) :: work integer(fgsl_int), dimension(*), intent(inout) :: work_int real(fgsl_double) :: fgsl_stats_qn_from_sorted_data end function fgsl_stats_qn_from_sorted_data