Git commit db6b7da941de2069914d1f76e2a9b5e177a30975 by Stefan Gerlach. Committed on 31/05/2017 at 22:22. Pushed by sgerlach into branch 'master'. 1. finished distribution fit models with F distribution 2. clean up code M +2 -2 src/backend/gsl/functions.h M +89 -76 src/backend/nsl/nsl_fit.c M +1 -0 src/backend/nsl/nsl_fit.h M +2 -2 src/backend/nsl/nsl_sf_stats.c M +17 -2 src/backend/worksheet/plots/cartesian/XYFitCurve.cpp M +10 -10 src/kdefrontend/dockwidgets/XYFitCurveDock.cpp D +- -- src/pics/gsl_distributions/F.jpg A +- -- src/pics/gsl_distributions/fdist.jpg https://commits.kde.org/labplot/db6b7da941de2069914d1f76e2a9b5e177a30975 diff --git a/src/backend/gsl/functions.h b/src/backend/gsl/functions.h index da2f9b6c..4caeae1c 100644 --- a/src/backend/gsl/functions.h +++ b/src/backend/gsl/functions.h @@ -44,8 +44,8 @@ double my_random() { return random(); } double my_drand() { return random()/(double)RAND_MAX; } /* math.h */ #ifndef _WIN32 -double my_jn(double n, double x) { return jn((int)n, x); } -double my_yn(double n,double x) { return yn((int)n, x); } +double my_jn(double n, double x) { return jn(round(n), x); } +double my_yn(double n,double x) { return yn(round(n), x); } double my_sgn(double x) { return copysign(1.0, x); } #else double my_sgn(double x) { diff --git a/src/backend/nsl/nsl_fit.c b/src/backend/nsl/nsl_fit.c index 2888f6fc..f5049093 100644 --- a/src/backend/nsl/nsl_fit.c +++ b/src/backend/nsl/nsl_fit.c @@ -183,49 +183,49 @@ double nsl_fit_model_fourier_param_deriv(int param, i= nt degree, double x, double } = /* peak */ -double nsl_fit_model_gaussian_param_deriv(int param, double x, double s, d= ouble mu, double a, double weight) { +double nsl_fit_model_gaussian_param_deriv(int param, double x, double s, d= ouble mu, double A, double weight) { double s2 =3D s*s, norm =3D weight/sqrt(2.*M_PI)/s, efactor =3D exp(-(x-m= u)*(x-mu)/(2.*s2)); = if (param =3D=3D 0) - return a * norm/(s*s2) * ((x-mu)*(x-mu) - s2) * efactor; + return A * norm/(s*s2) * ((x-mu)*(x-mu) - s2) * efactor; if (param =3D=3D 1) - return a * norm/s2 * (x-mu) * efactor; + return A * norm/s2 * (x-mu) * efactor; if (param =3D=3D 2) return norm * efactor; = return 0; } -double nsl_fit_model_lorentz_param_deriv(int param, double x, double s, do= uble t, double a, double weight) { +double nsl_fit_model_lorentz_param_deriv(int param, double x, double s, do= uble t, double A, double weight) { double norm =3D weight/M_PI, denom =3D s*s+(x-t)*(x-t); = if (param =3D=3D 0) - return a * norm * ((x-t)*(x-t) - s*s)/(denom*denom); + return A * norm * ((x-t)*(x-t) - s*s)/(denom*denom); if (param =3D=3D 1) - return a * norm * 2.*s*(x-t)/(denom*denom); + return A * norm * 2.*s*(x-t)/(denom*denom); if (param =3D=3D 2) return norm * s/denom; = return 0; } -double nsl_fit_model_sech_param_deriv(int param, double x, double s, doubl= e mu, double a, double weight) { +double nsl_fit_model_sech_param_deriv(int param, double x, double s, doubl= e mu, double A, double weight) { double y =3D (x-mu)/s, norm =3D weight/M_PI/s; = if (param =3D=3D 0) - return a/s * norm * (y*tanh(y)-1.)/cosh(y); + return A/s * norm * (y*tanh(y)-1.)/cosh(y); if (param =3D=3D 1) - return a/s * norm * tanh(y)/cosh(y); + return A/s * norm * tanh(y)/cosh(y); if (param =3D=3D 2) return norm/cosh(y); = return 0; } -double nsl_fit_model_logistic_param_deriv(int param, double x, double s, d= ouble mu, double a, double weight) { +double nsl_fit_model_logistic_param_deriv(int param, double x, double s, d= ouble mu, double A, double weight) { double y =3D (x-mu)/2./s, norm =3D weight/4./s; = if (param =3D=3D 0) - return a/s * norm * (2.*y*tanh(y)-1.)/cosh(y); + return A/s * norm * (2.*y*tanh(y)-1.)/cosh(y); if (param =3D=3D 1) - return a/s * norm * tanh(y)/cosh(y)/cosh(y); + return A/s * norm * tanh(y)/cosh(y)/cosh(y); if (param =3D=3D 2) return norm/cosh(y)/cosh(y); = @@ -233,67 +233,67 @@ double nsl_fit_model_logistic_param_deriv(int param, = double x, double s, double } = /* growth */ -double nsl_fit_model_atan_param_deriv(int param, double x, double s, doubl= e mu, double a, double weight) { +double nsl_fit_model_atan_param_deriv(int param, double x, double s, doubl= e mu, double A, double weight) { double norm =3D weight, y =3D (x-mu)/s; if (param =3D=3D 0) - return -a/s * norm * y/(1.+y*y); + return -A/s * norm * y/(1.+y*y); if (param =3D=3D 1) - return -a/s * norm * 1./(1+y*y); + return -A/s * norm * 1./(1+y*y); if (param =3D=3D 2) return norm * atan(y); = return 0; } -double nsl_fit_model_tanh_param_deriv(int param, double x, double s, doubl= e mu, double a, double weight) { +double nsl_fit_model_tanh_param_deriv(int param, double x, double s, doubl= e mu, double A, double weight) { double norm =3D weight, y =3D (x-mu)/s; if (param =3D=3D 0) - return -a/s * norm * y/cosh(y)/cosh(y); + return -A/s * norm * y/cosh(y)/cosh(y); if (param =3D=3D 1) - return -a/s * norm * 1./cosh(y)/cosh(y); + return -A/s * norm * 1./cosh(y)/cosh(y); if (param =3D=3D 2) return norm * tanh(y); = return 0; } -double nsl_fit_model_algebraic_sigmoid_param_deriv(int param, double x, do= uble s, double mu, double a, double weight) { +double nsl_fit_model_algebraic_sigmoid_param_deriv(int param, double x, do= uble s, double mu, double A, double weight) { double norm =3D weight, y =3D (x-mu)/s, y2 =3D y*y; if (param =3D=3D 0) - return -a/s * norm * y/pow(1.+y2, 1.5); + return -A/s * norm * y/pow(1.+y2, 1.5); if (param =3D=3D 1) - return -a/s * norm * 1./pow(1.+y2, 1.5); + return -A/s * norm * 1./pow(1.+y2, 1.5); if (param =3D=3D 2) return norm * y/sqrt(1.+y2); = return 0; } -double nsl_fit_model_sigmoid_param_deriv(int param, double x, double k, do= uble mu, double a, double weight) { +double nsl_fit_model_sigmoid_param_deriv(int param, double x, double k, do= uble mu, double A, double weight) { double norm =3D weight, y =3D k*(x-mu); if (param =3D=3D 0) - return a/k * norm * y*exp(-y)/gsl_pow_2(1. + exp(-y)); + return A/k * norm * y*exp(-y)/gsl_pow_2(1. + exp(-y)); if (param =3D=3D 1) - return -a*k * norm * exp(-y)/gsl_pow_2(1. + exp(-y)); + return -A*k * norm * exp(-y)/gsl_pow_2(1. + exp(-y)); if (param =3D=3D 2) return norm/(1. + exp(-y)); = return 0; } -double nsl_fit_model_erf_param_deriv(int param, double x, double s, double= mu, double a, double weight) { +double nsl_fit_model_erf_param_deriv(int param, double x, double s, double= mu, double A, double weight) { double norm =3D weight, y =3D (x-mu)/(sqrt(2.)*s); if (param =3D=3D 0) - return -a/sqrt(M_PI)/s * norm * y*exp(-y*y); + return -A/sqrt(M_PI)/s * norm * y*exp(-y*y); if (param =3D=3D 1) - return -a/sqrt(2.*M_PI)/s * norm * exp(-y*y); + return -A/sqrt(2.*M_PI)/s * norm * exp(-y*y); if (param =3D=3D 2) return norm/2. * gsl_sf_erf(y); = return 0; } -double nsl_fit_model_hill_param_deriv(int param, double x, double s, doubl= e n, double a, double weight) { +double nsl_fit_model_hill_param_deriv(int param, double x, double s, doubl= e n, double A, double weight) { double norm =3D weight, y =3D x/s; if (param =3D=3D 0) - return -a*n/s * norm * pow(y, n)/gsl_pow_2(1.+pow(y, n)); + return -A*n/s * norm * pow(y, n)/gsl_pow_2(1.+pow(y, n)); if (param =3D=3D 1) - return a * norm * log(y)*pow(y, n)/gsl_pow_2(1.+pow(y, n)); + return A * norm * log(y)*pow(y, n)/gsl_pow_2(1.+pow(y, n)); if (param =3D=3D 2) return norm * pow(y, n)/(1.+pow(y, n)); = @@ -309,12 +309,12 @@ double nsl_fit_model_gompertz_param_deriv(int param, = double x, double a, double = return 0; } -double nsl_fit_model_gudermann_param_deriv(int param, double x, double s, = double mu, double a, double weight) { +double nsl_fit_model_gudermann_param_deriv(int param, double x, double s, = double mu, double A, double weight) { double norm =3D weight, y =3D (x-mu)/s; if (param =3D=3D 0) - return -a/s * norm * y/cosh(y); + return -A/s * norm * y/cosh(y); if (param =3D=3D 1) - return -a/s * norm * 1./cosh(y); + return -A/s * norm * 1./cosh(y); if (param =3D=3D 2) return -asin(tanh(y)); = @@ -338,27 +338,27 @@ double nsl_fit_model_gaussian_tail_param_deriv(int pa= ram, double x, double s, do = return 0; } -double nsl_fit_model_exponential_param_deriv(int param, double x, double l= , double mu, double a, double weight) { +double nsl_fit_model_exponential_param_deriv(int param, double x, double l= , double mu, double A, double weight) { if (x < mu) return 0; double y =3D l*(x-mu), efactor =3D exp(-y); = if (param =3D=3D 0) - return weight * a * (1. - y) * efactor; + return weight * A * (1. - y) * efactor; if (param =3D=3D 1) - return weight * a * gsl_pow_2(l) * efactor; + return weight * A * gsl_pow_2(l) * efactor; if (param =3D=3D 2) return weight * l * efactor; = return 0; } -double nsl_fit_model_laplace_param_deriv(int param, double x, double s, do= uble mu, double a, double weight) { +double nsl_fit_model_laplace_param_deriv(int param, double x, double s, do= uble mu, double A, double weight) { double norm =3D weight/(2.*s), y =3D fabs((x-mu)/s), efactor =3D exp(-y); = if (param =3D=3D 0) - return a/s*norm * (y-1.) * efactor; + return A/s*norm * (y-1.) * efactor; if (param =3D=3D 1) - return a/(s*s)*norm * (x-mu)/y * efactor; + return A/(s*s)*norm * (x-mu)/y * efactor; if (param =3D=3D 2) return norm * efactor; = @@ -388,35 +388,35 @@ double nsl_fit_model_maxwell_param_deriv(int param, d= ouble x, double a, double c = return 0; } -double nsl_fit_model_poisson_param_deriv(int param, double x, double l, do= uble a, double weight) { +double nsl_fit_model_poisson_param_deriv(int param, double x, double l, do= uble A, double weight) { double norm =3D weight*pow(l, x)/gsl_sf_gamma(x+1.); = if (param =3D=3D 0) - return a/l * norm *(x-l)*exp(-l); + return A/l * norm *(x-l)*exp(-l); if (param =3D=3D 1) return norm * exp(-l); = return 0; } -double nsl_fit_model_lognormal_param_deriv(int param, double x, double s, = double mu, double a, double weight) { +double nsl_fit_model_lognormal_param_deriv(int param, double x, double s, = double mu, double A, double weight) { double norm =3D weight/sqrt(2.*M_PI)/(x*s), y =3D log(x)-mu, efactor =3D = exp(-(y/s)*(y/s)/2.); = if (param =3D=3D 0) - return a * norm * (y*y - s*s) * efactor; + return A * norm * (y*y - s*s) * efactor; if (param =3D=3D 1) - return a * norm * y/(s*s) * efactor; + return A * norm * y/(s*s) * efactor; if (param =3D=3D 2) return norm * efactor; = return 0; } -double nsl_fit_model_gamma_param_deriv(int param, double x, double t, doub= le k, double a, double weight) { +double nsl_fit_model_gamma_param_deriv(int param, double x, double t, doub= le k, double A, double weight) { double factor =3D weight*pow(x, k-1.)/pow(t, k)/gsl_sf_gamma(k), efactor = =3D exp(-x/t); = if (param =3D=3D 0) - return a * factor/t * (x/t-k) * efactor; + return A * factor/t * (x/t-k) * efactor; if (param =3D=3D 1) - return a * factor * (log(x/t) - gsl_sf_psi(k)) * efactor; + return A * factor * (log(x/t) - gsl_sf_psi(k)) * efactor; if (param =3D=3D 2) return factor * efactor; = @@ -435,35 +435,35 @@ double nsl_fit_model_flat_param_deriv(int param, doub= le x, double a, double b, d = return 0; } -double nsl_fit_model_rayleigh_param_deriv(int param, double x, double s, d= ouble a, double weight) { +double nsl_fit_model_rayleigh_param_deriv(int param, double x, double s, d= ouble A, double weight) { double y=3Dx/s, norm =3D weight*y/s, efactor =3D exp(-y*y/2.); = if (param =3D=3D 0) - return a*y/(s*s) * (y*y-2.)*efactor; + return A*y/(s*s) * (y*y-2.)*efactor; if (param =3D=3D 1) return norm * efactor; = return 0; } -double nsl_fit_model_rayleigh_tail_param_deriv(int param, double x, double= s, double mu, double a, double weight) { +double nsl_fit_model_rayleigh_tail_param_deriv(int param, double x, double= s, double mu, double A, double weight) { double norm =3D weight*x/(s*s), y =3D (mu*mu - x*x)/2./(s*s); = if (param =3D=3D 0) - return -2. * a * norm/s * (1. + y) * exp(y); + return -2. * A * norm/s * (1. + y) * exp(y); if (param =3D=3D 1) - return a * mu * norm/(s*s) * exp(y); + return A * mu * norm/(s*s) * exp(y); if (param =3D=3D 2) return norm * exp(y); = return 0; = } -double nsl_fit_model_levy_param_deriv(int param, double x, double g, doubl= e mu, double a, double weight) { +double nsl_fit_model_levy_param_deriv(int param, double x, double g, doubl= e mu, double A, double weight) { double y=3Dx-mu, norm =3D weight*sqrt(g/(2.*M_PI))/pow(y, 1.5), efactor = =3D exp(-g/2./y); = if (param =3D=3D 0) - return a/2.*norm/g/y * (y - g) * efactor; + return A/2.*norm/g/y * (y - g) * efactor; if (param =3D=3D 1) - return a/2.*norm/y/y * (3.*y - g) * efactor; + return A/2.*norm/y/y * (3.*y - g) * efactor; if (param =3D=3D 2) return norm * efactor; = @@ -475,25 +475,38 @@ double nsl_fit_model_landau_param_deriv(int param, do= uble x, double weight) { = return 0; } -double nsl_fit_model_chi_square_param_deriv(int param, double x, double n,= double a, double weight) { +double nsl_fit_model_chi_square_param_deriv(int param, double x, double n,= double A, double weight) { double y=3Dn/2., norm =3D weight*pow(x, y-1.)/pow(2., y)/gsl_sf_gamma(y),= efactor =3D exp(-x/2.); = if (param =3D=3D 0) - return a/2. * norm * (log(x/2.) - gsl_sf_psi(y)) * efactor; + return A/2. * norm * (log(x/2.) - gsl_sf_psi(y)) * efactor; if (param =3D=3D 1) return norm * efactor; = return 0; } -double nsl_fit_model_students_t_param_deriv(int param, double x, double n,= double a, double weight) { +double nsl_fit_model_students_t_param_deriv(int param, double x, double n,= double A, double weight) { if (param =3D=3D 0) - return weight * a * gsl_sf_gamma((n+1.)/2.)/2./pow(n, 1.5)/sqrt(M_PI)/gs= l_sf_gamma(n/2.) * pow(1.+x*x/n, - (n+3.)/2.) + return weight * A * gsl_sf_gamma((n+1.)/2.)/2./pow(n, 1.5)/sqrt(M_PI)/gs= l_sf_gamma(n/2.) * pow(1.+x*x/n, - (n+3.)/2.) * (x*x - 1. - (n+x*x)*log(1.+x*x/n) + (n+x*x)*(gsl_sf_psi((n+1.)/2.) -= gsl_sf_psi(n/2.)) ) ; if (param =3D=3D 1) return weight * gsl_ran_tdist_pdf(x, n); = return 0; } +double nsl_fit_model_fdist_param_deriv(int param, double x, double n1, dou= ble n2, double A, double weight) { + double norm =3D weight * gsl_sf_gamma((n1+n2)/2.)/gsl_sf_gamma(n1/2.)/gsl= _sf_gamma(n2/2.) * pow(n1, n1/2.) * pow(n2, n2/2.) * pow(x, n1/2.-1.); + double y =3D n2+n1*x; + + if (param =3D=3D 0) + return A/2. * norm * pow(y, -(n1+n2+2.)/2.) * (n2*(1.-x) + y*(log(n1) + = log(x) - log(y) + gsl_sf_psi((n1+n2)/2.) - gsl_sf_psi(n1/2.))); + if (param =3D=3D 1) + return A/2. * norm * pow(y, -(n1+n2+2.)/2.) * (n1*(x-1.) + y*(log(n2) - = log(y) + gsl_sf_psi((n1+n2)/2.) - gsl_sf_psi(n2/2.))); + if (param =3D=3D 2) + return weight * gsl_ran_fdist_pdf(x, n1, n2); + + return 0; +} double nsl_fit_model_beta_param_deriv(int param, double x, double a, doubl= e b, double A, double weight) { double norm =3D weight * A * gsl_sf_gamma(a+b)/gsl_sf_gamma(a)/gsl_sf_gam= ma(b) * pow(x, a-1.) * pow(1.-x, b-1.); = @@ -520,43 +533,43 @@ double nsl_fit_model_pareto_param_deriv(int param, do= uble x, double a, double b, = return 0; } -double nsl_fit_model_weibull_param_deriv(int param, double x, double k, do= uble l, double mu, double a, double weight) { +double nsl_fit_model_weibull_param_deriv(int param, double x, double k, do= uble l, double mu, double A, double weight) { double y =3D (x-mu)/l, z =3D pow(y, k), efactor =3D exp(-z); = if (param =3D=3D 0) - return weight*a/l * z/y*(k*log(y)*(1.-z) + 1.) * efactor; + return weight * A/l * z/y*(k*log(y)*(1.-z) + 1.) * efactor; if (param =3D=3D 1) - return weight*a*k*k/l/l * z/y*(z-1.) * efactor; + return weight * A*k*k/l/l * z/y*(z-1.) * efactor; if (param =3D=3D 2) - return weight*a*k/l/l * z/y/y*(k*z + 1. - k) * efactor; + return weight * A*k/l/l * z/y/y*(k*z + 1. - k) * efactor; if (param =3D=3D 3) - return weight*k/l * z/y * efactor; + return weight * k/l * z/y * efactor; = return 0; } -double nsl_fit_model_frechet_param_deriv(int param, double x, double g, do= uble mu, double s, double a, double weight) { +double nsl_fit_model_frechet_param_deriv(int param, double x, double g, do= uble mu, double s, double A, double weight) { double y =3D (x-mu)/s, efactor =3D exp(-pow(y, -g)); = if (param =3D=3D 0) - return weight*a/s * pow(y, -2.*g-1.) * (g*log(y)*(1.-pow(y, g))+pow(y, g= )) * efactor; + return weight * A/s * pow(y, -2.*g-1.) * (g*log(y)*(1.-pow(y, g))+pow(y,= g)) * efactor; if (param =3D=3D 1) - return a*weight * g/(s*s)*pow(y, -g-2.) * (g+1.-g*pow(y, -g)) * efactor; + return A * weight * g/(s*s)*pow(y, -g-2.) * (g+1.-g*pow(y, -g)) * efacto= r; if (param =3D=3D 2) - return a*weight * gsl_pow_2(g/s)*pow(y, -2.*g-1.) * (pow(y, g)-1.) * efa= ctor; + return A * weight * gsl_pow_2(g/s)*pow(y, -2.*g-1.) * (pow(y, g)-1.) * e= factor; if (param =3D=3D 3) - return g*weight/s * pow(y, -g-1.) * efactor; + return g * weight/s * pow(y, -g-1.) * efactor; = return 0; } -double nsl_fit_model_gumbel1_param_deriv(int param, double x, double s, do= uble b, double mu, double a, double weight) { +double nsl_fit_model_gumbel1_param_deriv(int param, double x, double s, do= uble b, double mu, double A, double weight) { double norm =3D weight/s, y =3D (x-mu)/s, efactor =3D exp(-y - b*exp(-y)); = if (param =3D=3D 0) - return a/s * norm * (y - 1. - b*exp(-y)) * efactor; + return A/s * norm * (y - 1. - b*exp(-y)) * efactor; if (param =3D=3D 1) - return -a * norm * exp(-y) * efactor; + return -A * norm * exp(-y) * efactor; if (param =3D=3D 2) - return a/s * norm * (1. - b*exp(-y)) * efactor; + return A/s * norm * (1. - b*exp(-y)) * efactor; if (param =3D=3D 3) return norm * efactor; = @@ -641,13 +654,13 @@ double nsl_fit_model_logarithmic_param_deriv(int para= m, double k, double p, doub = return 0; } -double nsl_fit_model_sech_dist_param_deriv(int param, double x, double s, = double mu, double a, double weight) { +double nsl_fit_model_sech_dist_param_deriv(int param, double x, double s, = double mu, double A, double weight) { double norm =3D weight/2./s, y =3D M_PI/2.*(x-mu)/s; = if (param =3D=3D 0) - return -a/s * norm * (y*tanh(y)+1.)/cosh(y); + return -A/s * norm * (y*tanh(y)+1.)/cosh(y); if (param =3D=3D 1) - return a*M_PI/2./s * norm * tanh(y)/cosh(y); + return A*M_PI/2./s * norm * tanh(y)/cosh(y); if (param =3D=3D 2) return norm * 1./cosh(y); = diff --git a/src/backend/nsl/nsl_fit.h b/src/backend/nsl/nsl_fit.h index 78e52a64..6be9a5b7 100644 --- a/src/backend/nsl/nsl_fit.h +++ b/src/backend/nsl/nsl_fit.h @@ -109,6 +109,7 @@ double nsl_fit_model_rayleigh_tail_param_deriv(int para= m, double x, double s, do double nsl_fit_model_landau_param_deriv(int param, double x, double weight= ); double nsl_fit_model_chi_square_param_deriv(int param, double x, double n,= double a, double weight); double nsl_fit_model_students_t_param_deriv(int param, double x, double n,= double a, double weight); +double nsl_fit_model_fdist_param_deriv(int param, double x, double n1, dou= ble n2, double a, double weight); double nsl_fit_model_beta_param_deriv(int param, double x, double a, doubl= e b, double A, double weight); double nsl_fit_model_pareto_param_deriv(int param, double x, double a, dou= ble b, double A, double weight); double nsl_fit_model_weibull_param_deriv(int param, double x, double k, do= uble l, double mu, double a, double weight); diff --git a/src/backend/nsl/nsl_sf_stats.c b/src/backend/nsl/nsl_sf_stats.c index 4c8eb6d1..de41ea70 100644 --- a/src/backend/nsl/nsl_sf_stats.c +++ b/src/backend/nsl/nsl_sf_stats.c @@ -37,7 +37,7 @@ const char* nsl_sf_stats_distribution_name[] =3D {i18n("G= aussian (Normal)"), i18n( i18n("Hypergeometric"), i18n("Logarithmic"), i18n("Maxwell-Boltzmann"), i= 18n("Hyperbolic secant (sech)"), i18n("Levy"), i18n("Frechet (inverse Weibu= ll)")}; const char* nsl_sf_stats_distribution_pic_name[] =3D { "gaussian", "gaussian_tail", "exponential", "laplace", "exponential_power= ", "cauchy_lorentz", "rayleigh", "rayleigh_tail", "landau", - "levy_alpha_stable", "levy_skew_alpha_stable","gamma", "flat", "lognormal= ", "chi_squared", "F", "students_t", "beta", "logistic", + "levy_alpha_stable", "levy_skew_alpha_stable","gamma", "flat", "lognormal= ", "chi_squared", "fdist", "students_t", "beta", "logistic", "pareto", "weibull", "gumbel1", "gumbel2", "poisson", "bernoulli", "binom= ial", "negative_binomial", "pascal", "geometric", "hypergeometric", "logarithmic", "maxwell_boltzmann", "sech", "levy", "fr= echet"}; const char* nsl_sf_stats_distribution_equation[] =3D { @@ -46,7 +46,7 @@ const char* nsl_sf_stats_distribution_equation[] =3D { "a * x/s^2 * exp(-(x/s)^2/2)", "a*x/s^2 * exp((mu^2-x^2)/2/s^2)", "a*land= au(x)", "Levy alpha-stable", "Levy-skew", "a/gamma(k)/t^k * x^(k-1)*exp(-x/t)", "A/(b-a)*theta(b-x)*theta(x-a)", "a/sqrt(2*pi)/x/s * exp(-( (log(x)-mu)/s= )^2/2)", "a * x^(n/2.-1.)/2^(n/2.)/gamma(n/2.) * exp(-x/2.)", - "F", "a*gamma((n+1)/2)/sqrt(pi*n)/gamma(n/2) * (1+x^2/n)^(-(n+1)/2)", "ga= mma(a+b)/gamma(a)/gamma(b) * x^(a-1) * (1-x)^(b-1)", + "a * fdist(x, n1, n2)", "a*gamma((n+1)/2)/sqrt(pi*n)/gamma(n/2) * (1+x^2/= n)^(-(n+1)/2)", "gamma(a+b)/gamma(a)/gamma(b) * x^(a-1) * (1-x)^(b-1)", "a/4/s * sech((x-mu)/2/s)**2", "a/b * (b/x)^(a+1) * theta(x-b)", "a * k/l= * ((x-mu)/l)^(k-1) * exp(-((x-mu)/l)^k)", "a/s * exp(-(x-mu)/s - b*exp(-(x-mu)/s))", "A*a*b * (x-mu)^(-a-1) * exp(-= b*(x-mu)^(-a))", "a * l^x/gamma(x+1) * exp(-l)", "Bernoulli", "a * binomial(x, p, n)", "a * negative_binomial(x, p, n)", diff --git a/src/backend/worksheet/plots/cartesian/XYFitCurve.cpp b/src/bac= kend/worksheet/plots/cartesian/XYFitCurve.cpp index 96559535..ba43d45f 100644 --- a/src/backend/worksheet/plots/cartesian/XYFitCurve.cpp +++ b/src/backend/worksheet/plots/cartesian/XYFitCurve.cpp @@ -697,6 +697,22 @@ int func_df(const gsl_vector* paramValues, void* param= s, gsl_matrix* J) { } break; } + case nsl_sf_stats_fdist: { + double n1 =3D nsl_fit_map_bound(gsl_vector_get(paramValues, 0), min[0],= max[0]); + double n2 =3D nsl_fit_map_bound(gsl_vector_get(paramValues, 1), min[1],= max[1]); + double a =3D nsl_fit_map_bound(gsl_vector_get(paramValues, 2), min[2], = max[2]); + for (size_t i =3D 0; i < n; i++) { + x =3D xVector[i]; + + for (int j =3D 0; j < 3; j++) { + if (fixed[j]) + gsl_matrix_set(J, i, j, 0.); + else + gsl_matrix_set(J, i, j, nsl_fit_model_fdist_param_deriv(j, x, n1, n2= , a, weight[i])); + } + } + break; + } case nsl_sf_stats_beta: case nsl_sf_stats_pareto: { double a =3D nsl_fit_map_bound(gsl_vector_get(paramValues, 0), min[0], = max[0]); @@ -905,10 +921,9 @@ int func_df(const gsl_vector* paramValues, void* param= s, gsl_matrix* J) { } break; } - // TODO: not implemented yet: + // unused distributions case nsl_sf_stats_levy_alpha_stable: case nsl_sf_stats_levy_skew_alpha_stable: - case nsl_sf_stats_fdist: case nsl_sf_stats_bernoulli: break; } diff --git a/src/kdefrontend/dockwidgets/XYFitCurveDock.cpp b/src/kdefronte= nd/dockwidgets/XYFitCurveDock.cpp index aeed2280..ae4ebbe7 100644 --- a/src/kdefrontend/dockwidgets/XYFitCurveDock.cpp +++ b/src/kdefrontend/dockwidgets/XYFitCurveDock.cpp @@ -410,13 +410,12 @@ void XYFitCurveDock::categoryChanged(int index) { for(int i =3D 0; i < NSL_SF_STATS_DISTRIBUTION_COUNT; i++) uiGeneralTab.cbModel->addItem(nsl_sf_stats_distribution_name[i]); = - // non-used items are disabled here + // not-used items are disabled here const QStandardItemModel* model =3D qobject_cast(uiGeneralTab.cbModel->model()); = for(int i =3D 1; i < NSL_SF_STATS_DISTRIBUTION_COUNT; i++) { - //TODO: not implemented yet: - if (i =3D=3D nsl_sf_stats_levy_alpha_stable || i =3D=3D nsl_sf_stats_le= vy_skew_alpha_stable || - i =3D=3D nsl_sf_stats_fdist || i =3D=3D nsl_sf_stats_bernoulli) { + // unused distributions + if (i =3D=3D nsl_sf_stats_levy_alpha_stable || i =3D=3D nsl_sf_stats_le= vy_skew_alpha_stable || i =3D=3D nsl_sf_stats_bernoulli) { QStandardItem* item =3D model->item(i); item->setFlags(item->flags() & ~(Qt::ItemIsSelectable|Qt::ItemIsEnabl= ed)); } @@ -767,7 +766,6 @@ void XYFitCurveDock::updateModelEquation() { break; case nsl_fit_model_distribution: switch ((nsl_sf_stats_distribution)m_fitData.modelType) { - // TODO: add missing GSL distributions (see nsl_sf_stats.c) case nsl_sf_stats_gaussian: case nsl_sf_stats_laplace: case nsl_sf_stats_rayleigh_tail: @@ -802,8 +800,9 @@ void XYFitCurveDock::updateModelEquation() { m_fitData.paramNames << "a"; m_fitData.paramNamesUtf8 << "A"; break; - case nsl_sf_stats_levy_alpha_stable: // TODO - case nsl_sf_stats_levy_skew_alpha_stable: // TODO + case nsl_sf_stats_levy_alpha_stable: // unused distributions + case nsl_sf_stats_levy_skew_alpha_stable: + case nsl_sf_stats_bernoulli: break; case nsl_sf_stats_gamma: m_fitData.paramNames << "t" << "k" << "a"; @@ -816,7 +815,10 @@ void XYFitCurveDock::updateModelEquation() { m_fitData.paramNames << "n" << "a"; m_fitData.paramNamesUtf8 << "n" << "A"; break; - case nsl_sf_stats_fdist: // TODO + case nsl_sf_stats_fdist: + m_fitData.paramNames << "n1" << "n2" << "a"; + m_fitData.paramNamesUtf8 << QString::fromUtf8("\u03bd") + QString::from= Utf8("\u2081") + << QString::fromUtf8("\u03bd") + QString::fromUtf8("\u2082") << "A"; break; case nsl_sf_stats_tdist: m_fitData.paramNames << "n" << "a"; @@ -842,8 +844,6 @@ void XYFitCurveDock::updateModelEquation() { m_fitData.paramNames << "l" << "a"; m_fitData.paramNamesUtf8 << QString::fromUtf8("\u03bb") << "A"; break; - case nsl_sf_stats_bernoulli: // TODO - break; case nsl_sf_stats_binomial: case nsl_sf_stats_negative_binomial: case nsl_sf_stats_pascal: diff --git a/src/pics/gsl_distributions/F.jpg b/src/pics/gsl_distributions/= F.jpg deleted file mode 100644 index c2174bd6..00000000 Binary files a/src/pics/gsl_distributions/F.jpg and /dev/null differ diff --git a/src/pics/gsl_distributions/fdist.jpg b/src/pics/gsl_distributi= ons/fdist.jpg new file mode 100644 index 00000000..cf1267a9 Binary files /dev/null and b/src/pics/gsl_distributions/fdist.jpg differ