1 |
#include <stdio.h> |
2 |
#include <stdlib.h> |
3 |
#include <math.h> |
4 |
// #include "rdutil.h" |
5 |
|
6 |
int read_fit_v(FILE *fpt, int *npts, int **n, double **l, double **f, |
7 |
double **ef, double **ux, double **eux, double **uy, double **euy) { |
8 |
/* |
9 |
Function - read_fit |
10 |
Takes an open file pointer, and reads n, l, frequencies, and velocities |
11 |
*/ |
12 |
int i, nlines; |
13 |
char buffer[8192]; |
14 |
if( ferror(fpt) ) { |
15 |
//perror(fpt); |
16 |
return -1; |
17 |
} |
18 |
if( feof(fpt)) rewind(fpt); |
19 |
|
20 |
nlines = 0; |
21 |
while(!feof(fpt)) { |
22 |
fgets(buffer, 8192, fpt); |
23 |
if(buffer[0] != '#' && !feof(fpt)) nlines++; |
24 |
} |
25 |
if(nlines == 0) return 1; |
26 |
nlines++; |
27 |
*n = (int*) malloc(nlines*sizeof(int)); |
28 |
*l = (double *)calloc(nlines, sizeof(double)); |
29 |
*f = (double*) malloc(nlines*sizeof(double)); |
30 |
*ef = (double*) malloc(nlines*sizeof(double)); |
31 |
*ux = (double*) malloc(nlines*sizeof(double)); |
32 |
*eux = (double*) malloc(nlines*sizeof(double)); |
33 |
*uy = (double*) malloc(nlines*sizeof(double)); |
34 |
*euy = (double*) malloc(nlines*sizeof(double)); |
35 |
nlines--; |
36 |
rewind(fpt); |
37 |
|
38 |
for(i=0; i<nlines; i++) { |
39 |
fgets(buffer, 8192, fpt); |
40 |
while(buffer[0] == '#') fgets(buffer, 8192, fpt); |
41 |
sscanf(buffer, "%i %lf %*f %lf %lf %lf %lf %lf %lf", &(*n)[i], &(*l)[i], &(*f)[i], |
42 |
&(*ef)[i], &(*ux)[i], &(*eux)[i], &(*uy)[i], &(*euy)[i]); |
43 |
} |
44 |
*npts = nlines; |
45 |
//fprintf (stderr, "read_fit_v(): returning npts = %d\n", nlines); |
46 |
|
47 |
return 0; |
48 |
} |
49 |
|
50 |
int nearst(double xb, double x[], int ntab) { |
51 |
int low, igh, mid; |
52 |
|
53 |
low=0; igh=ntab-1; |
54 |
if((xb < x[low]) != (xb < x[igh]) ) { |
55 |
|
56 |
/* If the point is within the range of table, then locate it by bisection */ |
57 |
|
58 |
while(igh-low > 1) { |
59 |
mid=(low+igh)/2; |
60 |
if((xb < x[mid]) == (xb < x[low])) low=mid; |
61 |
else igh=mid; |
62 |
} |
63 |
} |
64 |
|
65 |
if(fabs(xb-x[low]) < fabs(xb-x[igh])) return low; |
66 |
else return igh; |
67 |
} |
68 |
|
69 |
int read_fit(FILE *fpt, int *npts, int **n, double **l, double **f, double **ef) { |
70 |
/* |
71 |
Function - read_fit |
72 |
Takes an open file pointer, and reads n, l, frequencies, and velocities |
73 |
*/ |
74 |
int i, nlines; |
75 |
char buffer[8192]; |
76 |
if( ferror(fpt) ) { |
77 |
//perror(fpt); |
78 |
return -1; |
79 |
} |
80 |
if( feof(fpt)) rewind(fpt); |
81 |
|
82 |
nlines = 0; |
83 |
while(!feof(fpt)) { |
84 |
fgets(buffer, 8192, fpt); |
85 |
if(buffer[0] != '#' && !feof(fpt)) nlines++; |
86 |
} |
87 |
if(nlines == 0) return 1; |
88 |
|
89 |
*n = (int*) malloc(nlines*sizeof(int)); |
90 |
*l = (double*) malloc(nlines*sizeof(double)); |
91 |
*f = (double*) malloc(nlines*sizeof(double)); |
92 |
*ef = (double*) malloc(nlines*sizeof(double)); |
93 |
|
94 |
rewind(fpt); |
95 |
|
96 |
for(i=0; i<nlines; i++) { |
97 |
fgets(buffer, 8192, fpt); |
98 |
while(buffer[0] == '#') fgets(buffer, 8192, fpt); |
99 |
sscanf(buffer, "%i %lf %*f %lf %lf", &(*n)[i], &(*l)[i], &(*f)[i], |
100 |
&(*ef)[i]); |
101 |
} |
102 |
*npts = nlines; |
103 |
|
104 |
return 0; |
105 |
} |
106 |
|
107 |
int divdif(double xb, double x[], double f[], int *nuse, int ntab, |
108 |
double fb[], double aeps, double *dfb, double *ddfb) { |
109 |
int i,j,k,next,in,ip,nit,ier, nmax=10; |
110 |
double err,px,dpx,ddpx,xn[11],xd[11]; |
111 |
|
112 |
/* Find the nearest point */ |
113 |
|
114 |
next=nearst(xb,x,ntab); |
115 |
fb[1]=f[next]; |
116 |
xd[1]=f[next]; |
117 |
xn[1]=x[next]; |
118 |
ier=0; |
119 |
px=1.0; |
120 |
|
121 |
/* Initialisation for the derivatives */ |
122 |
|
123 |
*dfb=0.0; *ddfb=0.0; |
124 |
dpx=0.0; ddpx=0.0; |
125 |
|
126 |
/* Points between IN and IP are used for interpolation */ |
127 |
|
128 |
ip=next; in=next; |
129 |
|
130 |
/* Maximum number of points to be used for interpolation */ |
131 |
nit=*nuse; if(nmax<nit) nit=nmax; if(ntab<nit) nit=ntab; |
132 |
if(*nuse>nmax || *nuse>ntab) ier=22; |
133 |
if(*nuse<1) { |
134 |
ier=21; |
135 |
nit=6; if(nmax<nit) nit=nmax; if(ntab<nit) nit=ntab; |
136 |
} |
137 |
*nuse=1; |
138 |
|
139 |
/* Calculate successive interpolation polynomial */ |
140 |
for(j=2; j<=nit; ++j) { |
141 |
/* Choose the next nearest point to XB */ |
142 |
if(in<=0 ) { |
143 |
ip=ip+1; next=ip; |
144 |
} |
145 |
else if(ip >= ntab-1) { |
146 |
in=in-1; next=in; |
147 |
} |
148 |
else if(fabs(xb-x[ip+1]) < fabs(xb-x[in-1]) ) { |
149 |
ip=ip+1; next=ip; |
150 |
} |
151 |
else { |
152 |
in=in-1; next=in; |
153 |
} |
154 |
|
155 |
/* Calculating the divided differences */ |
156 |
xd[j]=f[next]; |
157 |
xn[j]=x[next]; |
158 |
for(k=j-1; k>=1; --k) xd[k]=(xd[k+1]-xd[k])/(xn[j]-xn[k]); |
159 |
|
160 |
/* Calculating the derivatives */ |
161 |
ddpx=ddpx*(xb-xn[j-1])+2.*dpx; |
162 |
dpx=dpx*(xb-xn[j-1])+px; |
163 |
*dfb = *dfb+dpx*xd[1]; |
164 |
*ddfb = *ddfb+ddpx*xd[1]; |
165 |
|
166 |
px=px*(xb-xn[j-1]); |
167 |
err=xd[1]*px; |
168 |
fb[j]=fb[j-1]+err; |
169 |
*nuse=j; |
170 |
|
171 |
if(fabs(err) < aeps) return ier; |
172 |
} |
173 |
return 23; |
174 |
} |
175 |
|
176 |
int interp_vel(int *n, double *l, double *f, double *ef, double *ux, |
177 |
double *eux, double *uy, double *euy, int npts) { |
178 |
/* |
179 |
Function - interp |
180 |
Takes power spectra fit parameters, and interpolates them to integer l. |
181 |
Data are returned in the input arrays - **data is overwritten!** |
182 |
*/ |
183 |
int i, j, nuse, ierr, offset=0; |
184 |
int n_num[13]; |
185 |
double fb[10], ll; |
186 |
double dfb, ddfb, aeps=1e-7; |
187 |
double *inp_l, *inp_f, *inp_ef, *inp_ux, *inp_uy, *inp_eux, *inp_euy; |
188 |
FILE * outBug; |
189 |
|
190 |
// count number of frequencies for each n |
191 |
// note that it is assumed that all frequencies for each |
192 |
// n are grouped together in the arrays |
193 |
for(i=0; i<13; i++) n_num[i]=0; |
194 |
for(i=0; i < npts; i++) { |
195 |
n_num[n[i]] ++; |
196 |
} |
197 |
|
198 |
inp_l=(double*) malloc(npts*sizeof(double)); |
199 |
inp_f=(double*) malloc(npts*sizeof(double)); |
200 |
inp_ef=(double*) malloc(npts*sizeof(double)); |
201 |
inp_ux=(double*) malloc(npts*sizeof(double)); |
202 |
inp_uy=(double*) malloc(npts*sizeof(double)); |
203 |
inp_eux=(double*) malloc(npts*sizeof(double)); |
204 |
inp_euy=(double*) malloc(npts*sizeof(double)); |
205 |
// loop over n |
206 |
// outBug=fopen("tesbug","w"); |
207 |
for(i=0; i < 10; i++) { |
208 |
for(j=0;j<npts;j++) { |
209 |
inp_l[j]=0.0; |
210 |
inp_f[j]=0.0; |
211 |
inp_ef[j]=0.0; |
212 |
inp_ux[j]=0.0; |
213 |
inp_uy[j]=0.0; |
214 |
inp_eux[j]=0.0; |
215 |
inp_euy[j]=0.0; |
216 |
} |
217 |
for(j=0;j<n_num[i];j++) { |
218 |
inp_l[j]=l[j+offset]; |
219 |
inp_f[j]=f[j+offset]; |
220 |
inp_ef[j]=ef[j+offset]; |
221 |
inp_ux[j]=ux[j+offset]; |
222 |
inp_uy[j]=uy[j+offset]; |
223 |
inp_eux[j]=eux[j+offset]; |
224 |
inp_euy[j]=euy[j+offset]; |
225 |
} |
226 |
for(j=0;j<n_num[i];j++) { |
227 |
ll=(int)(inp_l[j]+0.5); |
228 |
nuse=3; |
229 |
ierr=divdif(ll, inp_l, inp_f, &nuse, |
230 |
n_num[i], fb, aeps, &dfb, &ddfb); |
231 |
f[j+offset]=fb[nuse]; |
232 |
nuse=3; |
233 |
ierr=divdif(ll, inp_l, inp_ef, &nuse, |
234 |
n_num[i], fb, aeps, &dfb, &ddfb); |
235 |
ef[j+offset]=fb[nuse]; |
236 |
nuse=3; |
237 |
ierr=divdif(ll, inp_l, inp_ux, &nuse, |
238 |
n_num[i], fb, aeps, &dfb, &ddfb); |
239 |
ux[j+offset]=fb[nuse]; |
240 |
nuse=3; |
241 |
ierr=divdif(ll, inp_l, inp_uy, &nuse, |
242 |
n_num[i], fb, aeps, &dfb, &ddfb); |
243 |
uy[j+offset]=fb[nuse]; |
244 |
nuse=3; |
245 |
ierr=divdif(ll, inp_l, inp_eux, &nuse, |
246 |
n_num[i], fb, aeps, &dfb, &ddfb); |
247 |
eux[j+offset]=fb[nuse]; |
248 |
nuse=3; |
249 |
ierr=divdif(ll, inp_l, inp_euy, &nuse, |
250 |
n_num[i], fb, aeps, &dfb, &ddfb); |
251 |
euy[j+offset]=fb[nuse]; |
252 |
l[j+offset]=ll; |
253 |
} |
254 |
offset+=n_num[i]; |
255 |
} |
256 |
|
257 |
return 0; |
258 |
} |
259 |
|
260 |
int interp_freq(int *n, double *l, double *f, double *ef, int npts, |
261 |
int *ni, int *li, double *fi, double *efi, int *npts_out) { |
262 |
int i, j, nn, ll, nmin, nmax, this_npts, nuse; |
263 |
double *this_l, *this_f, *this_ef, all, fb[5]; |
264 |
double reps, dfb, ddfb, ff, error; |
265 |
|
266 |
this_l = (double*) malloc(npts*sizeof(double)); |
267 |
this_f = (double*) malloc(npts*sizeof(double)); |
268 |
this_ef = (double*) malloc(npts*sizeof(double)); |
269 |
|
270 |
*npts_out = 0; |
271 |
nmin = 100; |
272 |
nmax = 0; |
273 |
for(i=0; i<npts; i++) { |
274 |
if(n[i] < nmin) nmin = n[i]; |
275 |
if(n[i] > nmax) nmax = n[i]; |
276 |
} |
277 |
for(nn=nmin; nn<=nmax; nn++) { |
278 |
j = 0; |
279 |
for(i=0; i<npts; i++) { |
280 |
if(n[i]==nn) { |
281 |
this_l[j] = l[i]; |
282 |
this_f[j] = f[i]; |
283 |
this_ef[j] = ef[i]; |
284 |
j++; |
285 |
} |
286 |
} |
287 |
this_npts = j; |
288 |
j = *npts_out; |
289 |
//printf("%i %i\n",this_npts, *npts_out); |
290 |
for(i=0; i<this_npts; i++) { |
291 |
ll = (int) this_l[i] + 0.5; |
292 |
all = (double) ll; // cast back to double for function argument |
293 |
nuse = 3; |
294 |
reps = 1.0e-7; |
295 |
//printf("got here\n"); |
296 |
divdif(all, this_l, this_f, &nuse, this_npts, fb, reps, &dfb, &ddfb); |
297 |
ff = fb[nuse]; |
298 |
nuse = 3; |
299 |
divdif(all, this_l, this_ef, &nuse, this_npts, fb, reps, &dfb, &ddfb); |
300 |
error = fb[nuse]; |
301 |
//printf("%i %i %f %f\n",nn,ll,ff,error); |
302 |
if(error>0) { |
303 |
ni[j] = nn; |
304 |
li[j] = ll; |
305 |
fi[j] = ff; |
306 |
efi[j] = error; |
307 |
j++; |
308 |
} |
309 |
} |
310 |
//printf("%i\n",j); |
311 |
*npts_out = j; |
312 |
} |
313 |
|
314 |
return 0; |
315 |
} |
316 |
|
317 |
int freqdif(int *n1, int *n2, int *l1, double *l2, double *f1, double *f2, |
318 |
double *ef1, double *ef2, int npts1, int npts2, int *n, int *l, double *f, |
319 |
double *df, double *edf, int *npts) { |
320 |
int i, j, nn, ll, *index, nmin, nmax, lmin, lmax, this_npts, ind; |
321 |
int nuse; |
322 |
double *this_l, *this_f, *this_ef, all, fb[5], dfb, ddfb, reps; |
323 |
double ff, error; |
324 |
|
325 |
*npts = 0; |
326 |
|
327 |
this_l = (double*) malloc(npts2*sizeof(double)); |
328 |
this_f = (double*) malloc(npts2*sizeof(double)); |
329 |
this_ef = (double*) malloc(npts2*sizeof(double)); |
330 |
|
331 |
nmin = n1[0]; |
332 |
nmax = n1[0]; |
333 |
// lmin = (int) l2[0]+0.5; |
334 |
// lmax = (int) l2[0]+0.5; |
335 |
for(i=1; i<npts1; i++) { |
336 |
if(n1[i] < nmin) nmin = n1[i]; |
337 |
if(n1[i] > nmax) nmax = n1[i]; |
338 |
// if(l2[i] < lmin) lmin = (int) l2[i]+0.5; |
339 |
// if(l2[i] > lmax) lmax = (int) l2[i]+0.5; |
340 |
} |
341 |
|
342 |
// index = (int*) malloc((nmax-nmin)*(lmax-lmin)*sizeof(int)); |
343 |
// for(i=0; i<(nmax-nmin)*(lmax-lmin); i++) index[i] = -1; |
344 |
// for(i=0; i<npts2; i++) { |
345 |
// ll = (int) l2[i]+0.5; |
346 |
// index[n2[i]*(lmax-lmin)+ll-lmin] = i; |
347 |
// } |
348 |
|
349 |
for(nn=nmin; nn<=nmax; nn++) { |
350 |
j = 0; |
351 |
lmin = 10000; |
352 |
lmax = 0; |
353 |
for(i=0; i<npts2; i++) { |
354 |
if(n2[i] == nn) { |
355 |
this_l[j] = l2[i]; |
356 |
this_f[j] = f2[i]; |
357 |
this_ef[j] = ef2[i]; |
358 |
if(lmin > l2[i]) lmin = (int) l2[i]+0.5; |
359 |
if(lmax < l2[i]) lmax = (int) l2[i]+0.5; |
360 |
j++; |
361 |
} |
362 |
} |
363 |
this_npts = j; |
364 |
//printf("this_npts = %i\n", this_npts); |
365 |
j = *npts; |
366 |
for(i=0; i<npts1; i++) { |
367 |
// ind = index[n1[i]*(lmax-lmin)+l1[i]-lmin]; |
368 |
// if(ind>=0) { |
369 |
if(n1[i] == nn && l1[i] >= lmin && l1[i] <= lmax) { |
370 |
//ll = (int) this_l[i] + 0.5; |
371 |
ll = l1[i]; |
372 |
all = (double) ll; // cast back to double for function argument |
373 |
nuse = 3; |
374 |
reps = 1.0e-7; |
375 |
divdif(all, this_l, this_f, &nuse, this_npts, fb, reps, &dfb, &ddfb); |
376 |
ff = fb[nuse]; |
377 |
nuse = 3; |
378 |
divdif(all, this_l, this_ef, &nuse, this_npts, fb, reps, &dfb, &ddfb); |
379 |
error = fb[nuse]; |
380 |
if(error>0.0) { |
381 |
n[j] = nn; |
382 |
l[j] = ll; |
383 |
f[j] = f1[i]; |
384 |
df[j] = f1[i] - ff; |
385 |
edf[j] = sqrt(ef1[i]*ef1[i]+error*error); |
386 |
//printf("%i %i %i %f %f %f\n",j, n[j], l[j], f[j], df[j], edf[j]); |
387 |
j++; |
388 |
} |
389 |
} |
390 |
} |
391 |
*npts = j; |
392 |
//printf("%i\n",*npts); |
393 |
} |
394 |
|
395 |
return 0; |
396 |
} |
397 |
|
398 |
int autoweed_vel(int* n, double* l, double *ux, double *uy, |
399 |
int *mask, int npts) { |
400 |
/* |
401 |
Function - autoweed_vel |
402 |
Function to automatically weed velocities. An int array |
403 |
'mask' is returned of length npts - mask[i] is False for |
404 |
rejected modes. Ux and Uy are weeded together. |
405 |
*/ |
406 |
const int maxn=8; |
407 |
const double tol=5.0; |
408 |
|
409 |
int i, j, offset; |
410 |
int n_num[maxn]; |
411 |
double num; |
412 |
double sumx, sumxx, sumy, sumyy, meanx, meany, stdx, stdy; |
413 |
double range, uxmin, uxmax, uymin, uymax; |
414 |
double lmin,lmax; |
415 |
|
416 |
for(i=0; i<maxn; i++) { |
417 |
n_num[i]=0; |
418 |
} |
419 |
for(i=0; i<npts; i++) { |
420 |
if(n[i] < maxn) n_num[n[i]]++; |
421 |
} |
422 |
|
423 |
offset=0; |
424 |
for(i=0; i<maxn; i++) { |
425 |
num=0.0; |
426 |
sumx=0.0; |
427 |
sumy=0.0; |
428 |
sumxx=0.0; |
429 |
sumyy=0.0; |
430 |
// get l limits |
431 |
//fprintf (stderr, "autoweed_vel(): i = %d: lmin = lmax = l[%d] = %f\n", |
432 |
//i, offset, l[offset]); |
433 |
lmin=lmax=l[offset]; |
434 |
for(j=offset; j<n_num[i]+offset; j++) { |
435 |
if(l[j] < lmin) lmin=l[j]; |
436 |
if(l[j] > lmax) lmax= l[j]; |
437 |
} |
438 |
range=(lmax-lmin)/5.0; |
439 |
//printf("%f %f %f",lmin,lmax,range); |
440 |
lmin+=2.0*range; |
441 |
lmax-=2.0*range; |
442 |
for(j=offset; j<n_num[i]+offset; j++) { |
443 |
if(l[j] <= lmax && l[j] >= lmin) { |
444 |
sumx+=ux[j]; |
445 |
sumxx+=ux[j]*ux[j]; |
446 |
sumy+=uy[j]; |
447 |
sumyy+=uy[j]*uy[j]; |
448 |
num++; |
449 |
} |
450 |
} |
451 |
meanx=sumx/num; |
452 |
meany=sumy/num; |
453 |
stdx=num*sumxx-sumx*sumx; |
454 |
stdy=num*sumyy-sumy*sumy; |
455 |
stdx/=(num*(num-1)); |
456 |
stdy/=(num*(num-1)); |
457 |
stdx=sqrt(stdx); |
458 |
stdy=sqrt(stdy); |
459 |
//printf("%f %f %f %f\n",stdx,stdy,lmin,lmax); |
460 |
uxmin=meanx-5.0*stdx; |
461 |
uxmax=meanx+5.0*stdx; |
462 |
uymin=meany-5.0*stdy; |
463 |
uymax=meany+5.0*stdy; |
464 |
//printf("%f %f %f %f\n",uxmin,uxmax,uymin,uymax); |
465 |
|
466 |
for(j=offset; j<n_num[i]+offset; j++) { |
467 |
if(ux[j] <= uxmax && ux[j] >= uxmin && |
468 |
uy[j] <= uymax && uy[j] >= uymin) mask[j]=1; |
469 |
else mask[j]=0; |
470 |
} |
471 |
offset+=n_num[i]; |
472 |
} |
473 |
return 0; |
474 |
} |
475 |
|
476 |
void line(int m, double *x, double *y) { |
477 |
y[0] = x[0]; |
478 |
y[1] = 1.0; |
479 |
return; |
480 |
} |
481 |
|
482 |
int svd(int n, int m, double *a, double *v, double sigma[], int la, int lv) |
483 |
|
484 |
{ |
485 |
int i,j,k,l,itr,ier, itmax=30; |
486 |
double f, g, h, rmax, s, r1, r2, c, x, y, z, aeps, eps=1.e-16; |
487 |
// double *e; |
488 |
double e[2]; |
489 |
|
490 |
if(n>m || n<=0 || m<=0 || n>la || n>lv) return 105; |
491 |
ier=0; |
492 |
|
493 |
/* Reduction to Bidiagonal form using Householder transformations */ |
494 |
g=0.0; rmax=0.0; |
495 |
// e=(double *) calloc(n, sizeof(double)); |
496 |
|
497 |
for(i=0; i<n; ++i) { |
498 |
/* Off-diagonal elements of bidiagonal form */ |
499 |
e[i]=g; |
500 |
s=0.0; |
501 |
for(j=i; j<m; ++j) s=s+a[i+j*la]*a[i+j*la]; |
502 |
if(s <= 0.0) { |
503 |
/* transformation not required */ |
504 |
g=0.0; |
505 |
} |
506 |
else { |
507 |
f= a[i+i*la]; |
508 |
g=sqrt(s); |
509 |
if(f>=0.0) g=-g; |
510 |
h=f*g-s; |
511 |
a[i+i*la] = f-g; |
512 |
|
513 |
for(j=i+1; j<n; ++j) { |
514 |
s=0.0; |
515 |
for(k=i; k<m; ++k) s=s+a[i+k*la]*a[j+k*la]; |
516 |
f=s/h; |
517 |
for(k=i; k<m; ++k) a[j+k*la]= a[j+k*la]+f*a[i+k*la]; |
518 |
} |
519 |
} |
520 |
|
521 |
/* Diagonal elements of bidiagonal form */ |
522 |
sigma[i]=g; |
523 |
s=0.0; |
524 |
for(j=i+1; j<n; ++j) s=s+a[j+i*la]*a[j+i*la]; |
525 |
|
526 |
if(s<= 0.0) g=0.0; |
527 |
else { |
528 |
f= a[i*la+(i+1)]; |
529 |
g=sqrt(s); |
530 |
if(f>= 0.0) g=-g; |
531 |
h=f*g-s; |
532 |
a[i*la+(i+1)]=f-g; |
533 |
for(j=i+1; j<n; ++j) e[j]=a[j+i*la]/h; |
534 |
|
535 |
for(j=i+1; j<m; ++j) { |
536 |
s=0.0; |
537 |
for(k=i+1; k<n; ++k) s=s+a[k+j*la]*a[k+i*la]; |
538 |
for(k=i+1; k<n; ++k) a[k+j*la] = a[k+j*la]+s*e[k]; |
539 |
} |
540 |
} |
541 |
r1=fabs(sigma[i])+fabs(e[i]); |
542 |
if(r1 > rmax) rmax=r1; |
543 |
} |
544 |
|
545 |
/* Accumulation of right hand transformation in array V */ |
546 |
for(i=n-1; i>=0; --i) { |
547 |
if(g != 0.0) { |
548 |
h=g*a[i*la+(i+1)]; |
549 |
for(j=i+1; j<n; ++j) v[i+j*lv]=a[j+i*la]/h; |
550 |
|
551 |
for(j=i+1; j<n; ++j) { |
552 |
s=0.0; |
553 |
for(k=i+1; k<n; ++k) s=s+a[k+i*la]*v[j+k*lv]; |
554 |
for(k=i+1; k<n; ++k) v[j+k*lv]=v[j+k*lv]+s*v[i+k*lv]; |
555 |
} |
556 |
} |
557 |
|
558 |
for(j=i+1; j<n; ++j) { |
559 |
v[j+i*lv]=0.0; v[i+j*lv]=0.0; |
560 |
} |
561 |
v[i+i*lv]=1; |
562 |
g= e[i]; |
563 |
} |
564 |
|
565 |
/* Accumulation of left hand transformation overwritten on matrix A */ |
566 |
for(i=n-1; i>=0; --i) { |
567 |
g=sigma[i]; |
568 |
for(j=i+1; j<n; ++j) a[j+i*la]=0.0; |
569 |
if(g != 0.0) { |
570 |
h=g*a[i+i*la]; |
571 |
|
572 |
for(j=i+1; j<n; ++j) { |
573 |
s=0.0; |
574 |
for(k=i+1; k<m; ++k) s=s+a[i+k*la]*a[j+k*la]; |
575 |
f=s/h; |
576 |
for(k=i; k<m; ++k) a[j+k*la]=a[j+k*la]+f*a[i+k*la]; |
577 |
} |
578 |
|
579 |
for(j=i; j<m; ++j) a[i+j*la]=a[i+j*la]/g; |
580 |
} |
581 |
else { |
582 |
for(j=i; j<m; ++j) a[i+j*la]=0.0; |
583 |
} |
584 |
a[i+i*la] = a[i+i*la]+1; |
585 |
} |
586 |
|
587 |
/* Diagonalisation of the bidiagonal form */ |
588 |
aeps=eps*rmax; |
589 |
/* Loop over the singular values */ |
590 |
for(k=n-1; k>=0; --k) { |
591 |
/* The QR transformation */ |
592 |
for(itr=1; itr<=itmax; ++itr) { |
593 |
|
594 |
/* Test for splitting */ |
595 |
for(l=k; l>=0; --l) { |
596 |
if(fabs(e[l]) < aeps) goto split; |
597 |
if(fabs(sigma[l-1]) < aeps) break; |
598 |
} |
599 |
|
600 |
/* cancellation of E[L] if L>1 */ |
601 |
c=0.0; s=1.0; |
602 |
for(i=l; i<=k; ++i) { |
603 |
f=s*e[i]; |
604 |
e[i] = c*e[i]; |
605 |
if(fabs(f) < aeps) goto split; |
606 |
g=sigma[i]; |
607 |
sigma[i]=sqrt(f*f+g*g); |
608 |
c=g/sigma[i]; |
609 |
s=-f/sigma[i]; |
610 |
|
611 |
for(j=0; j<m; ++j) { |
612 |
r1= a[j*la+(l-1)]; |
613 |
r2= a[i+j*la]; |
614 |
a[j*la+(l-1)]=r1*c+r2*s; |
615 |
a[i+j*la]=c*r2-s*r1; |
616 |
} |
617 |
} |
618 |
|
619 |
split: z=sigma[k]; |
620 |
if(l == k) { |
621 |
/* QR iteration has converged */ |
622 |
if(z < 0.0) { |
623 |
sigma[k] = -z; |
624 |
for(j=0; j<n; ++j) v[k+j*lv]=-v[k+j*lv]; |
625 |
} |
626 |
break; |
627 |
} |
628 |
|
629 |
if(itr==itmax) {ier=12; break;} |
630 |
|
631 |
/* calculating shift from bottom 2x2 minor */ |
632 |
x=sigma[l]; |
633 |
y=sigma[k-1]; |
634 |
g=e[k-1]; |
635 |
h=e[k]; |
636 |
f=((y-z)*(y+z)+(g-h)*(g+h))/(2.*h*y); |
637 |
g=sqrt(1.+f*f); |
638 |
if(f < 0.0) g=-g; |
639 |
f=((x-z)*(x+z)+h*(y/(f+g)-h))/x; |
640 |
|
641 |
/* next QR transformation */ |
642 |
c=1.0; s=1.0; |
643 |
/* Given's rotation */ |
644 |
for(i=l+1; i<=k; ++i) { |
645 |
g=e[i]; |
646 |
y=sigma[i]; |
647 |
h=s*g; |
648 |
g=c*g; |
649 |
e[i-1]=sqrt(f*f+h*h); |
650 |
c=f/e[i-1]; |
651 |
s=h/e[i-1]; |
652 |
f=c*x+s*g; |
653 |
g=c*g-s*x; |
654 |
h=s*y; |
655 |
y=c*y; |
656 |
|
657 |
for(j=0; j<n; ++j) { |
658 |
x=v[j*lv+(i-1)]; |
659 |
z=v[i+j*lv]; |
660 |
v[j*lv+(i-1)]=c*x+s*z; |
661 |
v[i+j*lv]=c*z-s*x; |
662 |
} |
663 |
|
664 |
sigma[i-1]=sqrt(f*f+h*h); |
665 |
if(sigma[i-1] != 0.0) { |
666 |
c=f/sigma[i-1]; |
667 |
s=h/sigma[i-1]; |
668 |
} |
669 |
f=c*g+s*y; |
670 |
x=c*y-s*g; |
671 |
for(j=0; j<m; ++j) { |
672 |
y= a[j*la+(i-1)]; |
673 |
z= a[i+j*la]; |
674 |
a[j*la+(i-1)] = c*y+s*z; |
675 |
a[i+j*la] = c*z-s*y; |
676 |
} |
677 |
} |
678 |
|
679 |
e[l]=0.0; |
680 |
e[k]=f; |
681 |
sigma[k]=x; |
682 |
} |
683 |
} |
684 |
// free(e); |
685 |
return ier; |
686 |
} |
687 |
|
688 |
int svdevl(int n, int m, double *u, double *v, double sigma[], int lu, |
689 |
int lv, double b[], double reps) |
690 |
|
691 |
{ |
692 |
int i,j; |
693 |
double smax, aeps, s; |
694 |
double *wk; |
695 |
|
696 |
/* Finding the largest singular value */ |
697 |
smax=0.0; |
698 |
for(i=0; i<n; ++i) |
699 |
if(sigma[i] > smax) smax=sigma[i]; |
700 |
|
701 |
aeps=smax*reps; |
702 |
wk=(double *)calloc(n, sizeof(double)); |
703 |
for(i=0; i<n; ++i) { |
704 |
s=0.0; |
705 |
/* Only SIGMA[I] > AEPS contribute to the solution */ |
706 |
if(sigma[i] > aeps) { |
707 |
for(j=0; j<m; ++j) s=s+b[j]*u[i+j*lu]; |
708 |
s=s/sigma[i]; |
709 |
} |
710 |
wk[i]=s; |
711 |
} |
712 |
|
713 |
for(i=0; i<n; ++i) { |
714 |
s=0.0; |
715 |
for(j=0; j<n; ++j) s=s+v[j+i*lv]*wk[j]; |
716 |
b[i]=s; |
717 |
} |
718 |
free(wk); |
719 |
return 0; |
720 |
} |
721 |
|
722 |
int llsq(int n, int m, int k, double *x, int ix, double f[], double ef[], |
723 |
double a[], double *u, double *v, int iu, int iv, double sigma[], |
724 |
double y[], void (* phi) (int , double * , double * ), double reps, |
725 |
double *chisq) |
726 |
|
727 |
{ |
728 |
int i,j,ier; |
729 |
double s1; |
730 |
// double *wk; |
731 |
double wk[2]; |
732 |
|
733 |
if(m>n || m<=0 || n<=0 || k>ix) return 606; |
734 |
|
735 |
// wk=(double *) calloc(m, sizeof(double)); |
736 |
/* Setting up the design matrix and the RHS */ |
737 |
for(i=0; i<n; ++i) { |
738 |
if(ef[i]<=0.0) {/*free(wk); */return 607;} |
739 |
a[i]=f[i]/ef[i]; |
740 |
phi(m,&x[i*ix],wk); |
741 |
for(j=0; j<m; ++j) u[j+i*iu]=wk[j]/ef[i]; |
742 |
} |
743 |
|
744 |
ier=svd(m,n,u,v,sigma,iu,iv); |
745 |
if(ier>100) {/*free(wk);*/ return ier;} |
746 |
|
747 |
/* Calculate the least squares solution */ |
748 |
ier=svdevl(m,n,u,v,sigma,iu,iv,a,reps); |
749 |
|
750 |
/* Computing the \chi^2 from fitted coefficients */ |
751 |
*chisq=0.0; |
752 |
for(i=0; i<n; ++i) { |
753 |
phi(m,&x[i*ix],wk); |
754 |
s1=0.0; |
755 |
for(j=0; j<m; ++j) s1=s1+a[j]*wk[j]; |
756 |
y[i]=s1; |
757 |
s1=(f[i]-y[i])/ef[i]; |
758 |
*chisq=(*chisq)+s1*s1; |
759 |
} |
760 |
|
761 |
// free(wk); |
762 |
return 0; |
763 |
} |
764 |
|
765 |
int autoweed(int *l, int *n, double *f, double *df, double *edf, int *msk, |
766 |
int npts) { |
767 |
int avglen = 4; |
768 |
double tol = 4.0; |
769 |
double delnu = 50.0; |
770 |
int i, j, nn, this_npts, llim, ulim, offset; |
771 |
int data_range[2], *this_l; |
772 |
double *this_f, *this_df, *this_edf, *slopes, *a, *u, *y; |
773 |
double v[4], sigma[2], chisq; |
774 |
double mean, std; |
775 |
|
776 |
this_l = (int*) malloc(npts*sizeof(int)); |
777 |
this_f = (double*) malloc(npts*sizeof(double)); |
778 |
this_df = (double*) malloc(npts*sizeof(double)); |
779 |
this_edf = (double*) malloc(npts*sizeof(double)); |
780 |
slopes = (double*) malloc(npts*sizeof(double)); |
781 |
a = (double*) malloc(avglen*sizeof(double)); |
782 |
u = (double*) malloc(2*avglen*sizeof(double)); |
783 |
y = (double*) malloc(avglen*sizeof(double)); |
784 |
offset = 0; |
785 |
for(i=0; i<npts; i++) msk[i] = 0; |
786 |
for(nn=0; nn<7; nn++) { |
787 |
j = 0; |
788 |
mean = 0.0; |
789 |
std = 0.0; |
790 |
for(i=0; i<npts; i++) { |
791 |
if(n[i] == nn) { |
792 |
this_l[j] = l[i]; |
793 |
this_f[j] = f[i]; |
794 |
this_df[j] = df[i]; |
795 |
this_edf[j] = edf[i]; |
796 |
j++; |
797 |
} |
798 |
} |
799 |
this_npts = j; |
800 |
if(this_npts > avglen + 5) { |
801 |
data_range[0] = (int) 2*this_npts/10.; |
802 |
data_range[1] = (int) 8*this_npts/10.; |
803 |
for(i=0; i<this_npts-avglen; i++) { |
804 |
llsq(avglen, 2, 1, &this_f[i], 1, &this_df[i], &this_edf[i], a, u, |
805 |
v, 2, 2, sigma, y, *line, 1.0e-6, &chisq); |
806 |
slopes[i] = a[0]; |
807 |
//printf("%i %f\n",nn,a[0]); |
808 |
mean += a[0]; |
809 |
} |
810 |
mean /= (double) this_npts; |
811 |
for(i=0; i<this_npts-avglen; i++) std += (slopes[i]-mean)*(slopes[i]-mean); |
812 |
std /= (double) this_npts; |
813 |
std = sqrt(std); |
814 |
// printf("%i %i %i\n",data_range[0], data_range[1], this_npts/2); |
815 |
llim = data_range[0]; |
816 |
ulim = data_range[1]; |
817 |
while(llim > 0) { |
818 |
if( fabs(slopes[llim]-mean)/std > tol) break; |
819 |
llim--; |
820 |
} |
821 |
while(ulim<this_npts-avglen-1) { |
822 |
if( fabs(slopes[ulim]-mean)/std > tol) break; |
823 |
ulim++; |
824 |
} |
825 |
i = this_npts/2; |
826 |
while(i>0) { |
827 |
if(this_f[i]-this_f[i-1] > delnu) break; |
828 |
i--; |
829 |
} |
830 |
if(i > llim) llim = i; |
831 |
i = this_npts/2; |
832 |
while(i<this_npts) { |
833 |
if(this_f[i]-this_f[i-1] > delnu) break; |
834 |
i++; |
835 |
} |
836 |
if(i < ulim) ulim = i; |
837 |
printf("%i %i %i %i\n",nn,this_npts,llim,ulim); |
838 |
for(i=offset+llim; i<offset+ulim; i++) msk[i] = 1;//i<this_npts+offset-ulim; i++) msk[i] = 1; |
839 |
offset += this_npts; |
840 |
} |
841 |
} |
842 |
|
843 |
return 0; |
844 |
} |
845 |
|