/* * hrs-scaling.c * * Intensity scaling using generalised HRS target function * * (c) 2006-2011 Thomas White * * Part of CrystFEL - crystallography with a FEL * */ #ifdef HAVE_CONFIG_H #include #endif #include #include #include #include #include #include #include #include "image.h" #include "peaks.h" #include "symmetry.h" #include "geometry.h" #include "cell.h" #include "utils.h" /* Maximum number of iterations of NLSq scaling per macrocycle. */ #define MAX_CYCLES (30) static void s_uhavha(signed int hat, signed int kat, signed int lat, struct image *images, int n, const char *sym, int a, double *uha, double *vha) { int k; double uha_val = 0.0; double vha_val = 0.0; for ( k=0; kcpeaks; if ( k != a ) continue; for ( hi=0; hin_cpeaks; hi++ ) { double ic, sigi; if ( !spots[hi].scalable ) continue; if ( spots[hi].h != hat ) continue; if ( spots[hi].k != kat ) continue; if ( spots[hi].l != lat ) continue; ic = spots[hi].intensity / spots[hi].p; sigi = sqrt(fabs(ic)); uha_val += 1.0 / pow(sigi, 2.0); vha_val += ic / pow(sigi, 2.0); } } if ( uha != NULL ) *uha = uha_val; if ( vha != NULL ) *vha = vha_val; } static double s_uh(struct image *images, int n, signed int h, signed int k, signed int l, const char *sym) { int a; double val = 0.0; for ( a=0; ah, it->k, it->l, sym); vh = s_vh(images, n, it->h, it->k, it->l, sym); set_intensity(uh_arr, it->h, it->k, it->l, uh); set_intensity(vh_arr, it->h, it->k, it->l, vh); } uha_arr = malloc(n*sizeof(double)); vha_arr = malloc(n*sizeof(double)); for ( h=0; hh; const signed int k = it->k; const signed int l = it->l; /* For this reflection, calculate all the possible * values of uha and vha */ for ( a=0; aosf * uha * Ih; vc = Ih * rha; /* Determine the matrix component */ for ( b=0; b a ) continue; uhb = uha_arr[b]; vhb = vha_arr[b]; rhb = vhb - image_b->osf * uhb * Ih; mc = (rha*vhb + vha*rhb - vha*vhb) / uh; if ( isnan(mc) ) mc = 0.0; /* 0 / 0 = 0 */ if ( a == b ) { mc += pow(Ih, 2.0) * uha; } tval = gsl_matrix_get(M, a, b); gsl_matrix_set(M, a, b, tval+mc); gsl_matrix_set(M, b, a, tval+mc); } vval = gsl_vector_get(v, a); gsl_vector_set(v, a, vval+vc); } } free(uh_arr); free(vh_arr); free(uha_arr); free(vha_arr); /* Fox and Holmes method */ gsl_eigen_symmv_workspace *work; gsl_vector *e_val; gsl_matrix *e_vec; int val; /* Diagonalise */ work = gsl_eigen_symmv_alloc(n); e_val = gsl_vector_alloc(n); e_vec = gsl_matrix_alloc(n, n); val = gsl_eigen_symmv(M, e_val, e_vec, work); gsl_eigen_symmv_free(work); val = gsl_eigen_symmv_sort(e_val, e_vec, GSL_EIGEN_SORT_ABS_DESC); /* Rotate the "solution vector" */ gsl_vector *rprime; rprime = gsl_vector_alloc(n); val = gsl_blas_dgemv(CblasTrans, 1.0, e_vec, v, 0.0, rprime); /* Solve (now easy) */ gsl_vector *sprime; sprime = gsl_vector_alloc(n); for ( frame=0; frame fabs(max_shift) ) { max_shift = fabs(shift); } } gsl_vector_free(shifts); gsl_vector_free(sprime); gsl_matrix_free(e_vec); gsl_vector_free(e_val); gsl_matrix_free(M); gsl_vector_free(v); return max_shift; } static double *lsq_intensities(struct image *images, int n, ReflItemList *obs, const char *sym) { double *I_full; int i; I_full = new_list_intensity(); for ( i=0; ih, it->k, it->l, &h, &k, &l, sym); /* For each frame */ for ( m=0; m 0.01) && (i < MAX_CYCLES) ); I_full = lsq_intensities(images, n, obs, sym); return I_full; }