/* * 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" #include "reflist.h" /* Maximum number of iterations of NLSq scaling per macrocycle. */ #define MAX_CYCLES (10) char *find_common_reflections(struct image *images, int n) { int i; char *cref; cref = calloc(n*n, sizeof(char)); for ( i=0; ireflections; Reflection *refl; for ( refl = find_refl(reflections, hat, kat, lat); refl != NULL; refl = next_found_refl(refl) ) { double ic, sigi; if ( !get_scalable(refl) ) continue; ic = get_intensity(refl) / get_partiality(refl); sigi = ic / 10.0; /* FIXME */ uha_val += 1.0 / pow(sigi, 2.0); vha_val += ic / pow(sigi, 2.0); } *uha = uha_val; *vha = vha_val; } static void s_uhvh(struct image *images, int n, signed int h, signed int k, signed int l, const char *sym, double *uhp, double *vhp) { int a; double uh = 0.0; double vh = 0.0; for ( a=0; asize; if ( v->size != M->size1 ) return NULL; if ( v->size != M->size2 ) return NULL; /* 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); if ( val ) { ERROR("Couldn't diagonalise matrix.\n"); return NULL; } 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; framesize; if ( v->size != M->size1 ) return NULL; if ( v->size != M->size2 ) return NULL; shifts = gsl_vector_alloc(n); if ( shifts == NULL ) return NULL; for ( frame=0; frameh, it->k, it->l, sym, &uh, &vh); 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 ( refidx=0; refidxh; 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; /* Off-diagonals zero if reference available */ if ( (reference != NULL) && (a!=b) ) continue; /* Zero if no common reflections */ if ( cref[a+n*b] != 0 ) 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 ( reference != NULL ) mc = 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); //show_matrix_eqn(M, v, n); if ( reference == NULL ) { shifts = solve_by_eigenvalue_filtration(v, M); } else { shifts = solve_diagonal(v, M); } gsl_matrix_free(M); gsl_vector_free(v); if ( shifts == NULL ) return 0.0; /* Apply shifts */ max_shift = 0.0; for ( frame=0; frame %5.2f\n", // frame, shift, images[frame].osf); if ( fabs(shift) > fabs(max_shift) ) { max_shift = fabs(shift); } } gsl_vector_free(shifts); return max_shift; } static RefList *lsq_intensities(struct image *images, int n, ReflItemList *obs, const char *sym) { RefList *full; int i; full = reflist_new(); for ( i=0; ih, it->k, it->l); refl != NULL; refl = next_found_refl(refl) ) { double p; if ( !get_scalable(refl) ) continue; p = get_partiality(refl); num += get_intensity(refl) * p * G; den += pow(p, 2.0) * pow(G, 2.0); redundancy++; } } if ( !isnan(num/den) ) { new = add_refl(full, it->h, it->k, it->l); set_int(new, num/den); set_redundancy(new, redundancy); } else { ERROR("Couldn't calculate LSQ full intensity for" "%3i %3i %3i\n", it->h, it->k, it->l); /* Doom is probably impending */ } } return full; } static void normalise_osfs(struct image *images, int n) { int m; double tot = 0.0; double norm; for ( m=0; m 0.01) && (i < MAX_CYCLES) ); full = lsq_intensities(images, n, obs, sym); return full; }