00001
00002
00003
00004
00005
00006
00007
00008
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056 #include <regress.h>
00057 #include <misc.h>
00058
00061 int
00062 regress(double *coef, double *x, double *y, double *cte, double *yreg, double *yerr,
00063 double *chisq, double *rsq, double *vif, double *autocor, int nterm, int npts) {
00081 gsl_vector_view row;
00082
00083 gsl_matrix *xx;
00084 gsl_matrix *cov;
00085
00086 gsl_multifit_linear_workspace *work;
00087
00088 gsl_vector *yy;
00089 gsl_vector *ccoef;
00090 gsl_vector *res;
00091
00092 int istat;
00093 int term;
00094 int vterm;
00095 int cterm;
00096 int pts;
00097
00098 size_t stride = 1;
00099
00100 double vchisq;
00101 double *ytmp;
00102 double ystd;
00103
00104
00105 xx = gsl_matrix_alloc(npts, nterm+1);
00106 yy = gsl_vector_alloc(npts);
00107
00108 ccoef = gsl_vector_alloc(nterm+1);
00109 cov = gsl_matrix_alloc(nterm+1, nterm+1);
00110
00111
00112 for (term=0; term<nterm; term++)
00113 for (pts=0; pts<npts; pts++)
00114 (void) gsl_matrix_set(xx, pts, term+1, x[pts+term*npts]);
00115
00116
00117 for (pts=0; pts<npts; pts++)
00118 (void) gsl_matrix_set(xx, pts, 0, 1.0);
00119
00120
00121 for (pts=0; pts<npts; pts++)
00122 (void) gsl_vector_set(yy, pts, y[pts]);
00123
00124
00125 work = gsl_multifit_linear_alloc(npts, nterm+1);
00126
00127
00128 istat = gsl_multifit_linear(xx, yy, ccoef, cov, chisq, work);
00129 if (istat != GSL_SUCCESS) {
00130 (void) fprintf(stderr, "%s: Line %d: Error %d in multifitting algorithm!\n", __FILE__, __LINE__, istat);
00131 (void) gsl_multifit_linear_free(work);
00132 (void) gsl_matrix_free(xx);
00133 (void) gsl_matrix_free(cov);
00134 (void) gsl_vector_free(yy);
00135 (void) gsl_vector_free(ccoef);
00136 return istat;
00137 }
00138
00139 (void) gsl_multifit_linear_free(work);
00140
00141
00142 #if DEBUG >= 7
00143 #define C(i) (gsl_vector_get(ccoef,(i)))
00144 #define COV(i,j) (gsl_matrix_get(cov,(i),(j)))
00145 #define X(i,j) (gsl_matrix_get(xx,(i),(j)))
00146
00147 {
00148 printf ("# best fit: Y = %g + %g X + %g X^2\n",
00149 C(0), C(1), C(2));
00150
00151 printf ("# x matrix:\n");
00152
00153 for (term=0; term<nterm+1; term++)
00154 for (pts=0; pts<npts; pts++)
00155 printf("term=%d pts=%d x=%lf\n", term, pts, X(pts,term));
00156
00157 printf ("# covariance matrix:\n");
00158 printf ("[ %+.5e, %+.5e, %+.5e \n",
00159 COV(0,0), COV(0,1), COV(0,2));
00160 printf (" %+.5e, %+.5e, %+.5e \n",
00161 COV(1,0), COV(1,1), COV(1,2));
00162 printf (" %+.5e, %+.5e, %+.5e ]\n",
00163 COV(2,0), COV(2,1), COV(2,2));
00164 printf ("# chisq = %g\n", *chisq);
00165 }
00166 #endif
00167
00168
00169 for (term=0; term<nterm; term++)
00170 coef[term] = gsl_vector_get(ccoef, term+1);
00171
00172 *cte = gsl_vector_get(ccoef, 0);
00173
00174
00175 *rsq = 1.0 - ( (*chisq) / gsl_stats_tss(y, stride, (size_t) pts) );
00176
00177
00178 #if DEBUG >= 7
00179 (void) fprintf(stdout, "%s: Vector reconstruction\n", __FILE__);
00180 #endif
00181 for (pts=0; pts<npts; pts++) {
00182 row = gsl_matrix_row(xx, pts);
00183 istat = gsl_multifit_linear_est(&row.vector, ccoef, cov, &(yreg[pts]), &ystd);
00184 if (istat != GSL_SUCCESS) {
00185 (void) fprintf(stderr, "%s: Line %d: Error %d in multifitting linear values estimation!\n", __FILE__, __LINE__, istat);
00186 (void) gsl_matrix_free(xx);
00187 (void) gsl_matrix_free(cov);
00188 (void) gsl_vector_free(yy);
00189 (void) gsl_vector_free(ccoef);
00190 return istat;
00191 }
00192 }
00193
00194 #if DEBUG >= 7
00195 (void) fprintf(stdout, "%s: Residuals calculation\n", __FILE__);
00196 #endif
00197 res = gsl_vector_alloc(npts);
00198 istat = gsl_multifit_linear_residuals(xx, yy, ccoef, res);
00199 if (istat != GSL_SUCCESS) {
00200 (void) fprintf(stderr, "%s: Line %d: Error %d in multifitting linear values residuals estimation!\n", __FILE__, __LINE__, istat);
00201 (void) gsl_matrix_free(xx);
00202 (void) gsl_matrix_free(cov);
00203 (void) gsl_vector_free(yy);
00204 (void) gsl_vector_free(res);
00205 (void) gsl_vector_free(ccoef);
00206 return istat;
00207 }
00208
00209 for (pts=0; pts<npts; pts++)
00210 yerr[pts] = gsl_vector_get(res, pts);
00211
00212
00213 *autocor = gsl_stats_lag1_autocorrelation(yerr, stride, npts);
00214
00215
00216 (void) gsl_matrix_free(xx);
00217 (void) gsl_matrix_free(cov);
00218 (void) gsl_vector_free(yy);
00219 (void) gsl_vector_free(res);
00220 (void) gsl_vector_free(ccoef);
00221
00222
00225
00226 #if DEBUG >= 7
00227 (void) fprintf(stdout, "%s: VIF calculation\n", __FILE__);
00228 #endif
00229 xx = gsl_matrix_alloc(npts, nterm);
00230 yy = gsl_vector_alloc(npts);
00231
00232 ccoef = gsl_vector_alloc(nterm);
00233 cov = gsl_matrix_alloc(nterm, nterm);
00234
00235 ytmp = (double *) malloc(npts * sizeof(double));
00236 if (ytmp == NULL) alloc_error(__FILE__, __LINE__);
00237
00238
00239 for (vterm=0; vterm<nterm; vterm++) {
00240
00241 (void) gsl_matrix_set_zero(xx);
00242 (void) gsl_matrix_set_zero(cov);
00243
00244 (void) gsl_vector_set_zero(yy);
00245 (void) gsl_vector_set_zero(ccoef);
00246
00247
00248 cterm = 0;
00249 for (term=0; term<nterm; term++) {
00250 if (term != vterm) {
00251 for (pts=0; pts<npts; pts++)
00252 (void) gsl_matrix_set(xx, pts, cterm+1, x[pts+term*npts]);
00253 cterm++;
00254 }
00255 }
00256
00257
00258 for (pts=0; pts<npts; pts++)
00259 (void) gsl_matrix_set(xx, pts, 0, 1.0);
00260
00261
00262 for (pts=0; pts<npts; pts++) {
00263 ytmp[pts] = x[pts+vterm*npts];
00264 (void) gsl_vector_set(yy, pts, ytmp[pts]);
00265 }
00266
00267
00268 work = gsl_multifit_linear_alloc(npts, nterm);
00269
00270
00271 istat = gsl_multifit_linear(xx, yy, ccoef, cov, &vchisq, work);
00272 if (istat != GSL_SUCCESS) {
00273 (void) fprintf(stderr, "%s: Line %d: Error %d in multifitting algorithm!\n", __FILE__, __LINE__, istat);
00274 (void) gsl_multifit_linear_free(work);
00275 (void) gsl_matrix_free(xx);
00276 (void) gsl_matrix_free(cov);
00277 (void) gsl_vector_free(yy);
00278 (void) gsl_vector_free(ccoef);
00279 (void) free(ytmp);
00280 return istat;
00281 }
00282
00283
00284 (void) gsl_multifit_linear_free(work);
00285
00286
00287
00288 vif[vterm] = 1.0 / (1.0 - (1.0 - ( vchisq / gsl_stats_tss(ytmp, stride, (size_t) pts) ) ));
00289 }
00290
00291
00292 (void) gsl_matrix_free(xx);
00293 (void) gsl_matrix_free(cov);
00294 (void) gsl_vector_free(yy);
00295 (void) gsl_vector_free(ccoef);
00296
00297 (void) free(ytmp);
00298
00299
00300 return 0;
00301 }