aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorThomas White <taw@bitwiz.org.uk>2010-12-17 19:53:39 -0800
committerThomas White <taw@physics.org>2012-02-22 15:27:10 +0100
commit32dd839ed7d10ef4b9dbcc5f11920dd6347e8ee4 (patch)
tree6a08165b16ebc29e47e6f486e9d33b5623a8d963 /src
parenta764d0a9512a7bd5c0e17af70b3c55e5fe7f9f71 (diff)
Remove row-crossing stuff since it should not be needed
Diffstat (limited to 'src')
-rw-r--r--src/hrs-scaling.c30
1 files changed, 8 insertions, 22 deletions
diff --git a/src/hrs-scaling.c b/src/hrs-scaling.c
index 3f03b8d8..ea87dec7 100644
--- a/src/hrs-scaling.c
+++ b/src/hrs-scaling.c
@@ -122,23 +122,19 @@ static double iterate_scale(struct image *images, int n,
gsl_vector *shifts;
int l;
double max_shift;
- int crossed = 0;
int n_ref;
- M = gsl_matrix_calloc(n-1, n-1);
- v = gsl_vector_calloc(n-1);
+ M = gsl_matrix_calloc(n, n);
+ v = gsl_vector_calloc(n);
n_ref = num_items(obs);
for ( l=0; l<n; l++ ) { /* "Equation number": one equation per frame */
int m; /* Frame (scale factor) number */
int h;
- int lcomp;
double vc_tot = 0.0;
struct image *imagel = &images[l];
- if ( l == crossed ) continue;
-
/* Determine the "solution" vector component */
for ( h=0; h<n_ref; h++ ) {
@@ -161,19 +157,15 @@ static double iterate_scale(struct image *images, int n,
}
- lcomp = l;
- if ( l > crossed ) lcomp--;
- gsl_vector_set(v, lcomp, vc_tot);
+ gsl_vector_set(v, l, vc_tot);
/* Now fill in the matrix components */
for ( m=0; m<n; m++ ) {
double mc_tot = 0.0;
- int mcomp;
struct image *imagem = &images[m];
if ( m > l ) continue; /* Matrix is symmetric */
- if ( m == crossed ) continue;
for ( h=0; h<n_ref; h++ ) {
@@ -204,29 +196,23 @@ static double iterate_scale(struct image *images, int n,
}
- mcomp = m;
- if ( m > crossed ) mcomp--;
- gsl_matrix_set(M, lcomp, mcomp, mc_tot);
- gsl_matrix_set(M, mcomp, lcomp, mc_tot);
+ gsl_matrix_set(M, l, m, mc_tot);
+ gsl_matrix_set(M, m, l, mc_tot);
}
}
- show_matrix_eqn(M, v, n-1);
+ show_matrix_eqn(M, v, n);
- shifts = gsl_vector_alloc(n-1);
+ shifts = gsl_vector_alloc(n);
gsl_linalg_HH_solve(M, v, shifts);
max_shift = 0.0;
for ( l=0; l<n-1; l++ ) {
double shift = gsl_vector_get(shifts, l);
- int limg;
-
- limg = l;
- if ( limg >= crossed ) limg++;
- images[limg].osf += shift;
+ images[l].osf += shift;
if ( fabs(shift) > fabs(max_shift) ) {
max_shift = fabs(shift);