/* * hrs-scaling.c * * Intensity scaling using generalised HRS target function * * (c) 2006-2010 Thomas White * * Part of CrystFEL - crystallography with a FEL * */ #ifdef HAVE_CONFIG_H #include #endif #include #include #include #include #include #include "image.h" #include "peaks.h" #include "symmetry.h" #include "geometry.h" #include "cell.h" /* Maximum number of iterations of NLSq scaling per macrocycle. */ #define MAX_CYCLES (30) static void show_matrix_eqn(gsl_matrix *M, gsl_vector *v, int r) { int i, j; for ( i=0; icpeaks; int found = 0; for ( hi=0; hin_cpeaks; hi++ ) { double ic; signed int ha, ka, la; if ( !spots[hi].valid ) continue; if ( spots[hi].p < 0.1 ) continue; get_asymm(spots[hi].h, spots[hi].k, spots[hi].l, &ha, &ka, &la, sym); if ( ha != hat ) continue; if ( ka != kat ) continue; if ( la != lat ) continue; ic = spots[hi].intensity / spots[hi].p; *sigma_GI += ic * image->osf; found = 1; } if ( found ) { *sigma_Gsq += pow(image->osf, 2.0); } } } static double find_occurrances(struct image *image, const char *sym, signed int h, signed int k, signed int l) { double Ihl = 0.0; int find; struct cpeak *spots = image->cpeaks; for ( find=0; findn_cpeaks; find++ ) { signed int ha, ka, la; if ( !spots[find].valid ) continue; if ( spots[find].p < 0.1 ) continue; get_asymm(spots[find].h, spots[find].k, spots[find].l, &ha, &ka, &la, sym); if ( ha != h ) continue; if ( ka != k ) continue; if ( la != l ) continue; Ihl += spots[find].intensity / spots[find].p; } return Ihl; } static double iterate_scale(struct image *images, int n, ReflItemList *obs, const char *sym) { gsl_matrix *M; gsl_vector *v; gsl_vector *shifts; int l; double max_shift; int n_ref; M = gsl_matrix_calloc(n, n); v = gsl_vector_calloc(n); n_ref = num_items(obs); for ( l=0; lh, it->k, it->l, &sigma_GI, &sigma_Gsq); /* Add up symmetric equivalents within the pattern */ Ihl = find_occurrances(imagel, sym, it->h, it->k, it->l); vc = Ihl * sigma_GI / sigma_Gsq; vc -= imagel->osf * pow(sigma_GI, 2.0) / sigma_Gsq; vc_tot += vc; } gsl_vector_set(v, l, vc_tot); /* Now fill in the matrix components */ for ( m=0; m l ) continue; /* Matrix is symmetric */ for ( h=0; hh, it->k, it->l, &sigma_GI, &sigma_Gsq); if ( l == m ) { mc += pow(sigma_GI, 2.0) / pow(sigma_Gsq, 2.0); } Ihl = find_occurrances(imagel, sym, it->h, it->k, it->l); Ihm = find_occurrances(imagem, sym, it->h, it->k, it->l); mc += Ihl * Ihm / sigma_Gsq; mc += (sigma_GI / pow(sigma_Gsq, 2.0) ) * ( imagel->osf*Ihm + imagem->osf * Ihl); mc_tot += mc; } gsl_matrix_set(M, l, m, mc_tot); gsl_matrix_set(M, m, l, mc_tot); } } show_matrix_eqn(M, v, n); shifts = gsl_vector_alloc(n); gsl_linalg_HH_solve(M, v, shifts); max_shift = 0.0; for ( l=0; l fabs(max_shift) ) { max_shift = fabs(shift); } } gsl_vector_free(shifts); 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) ); for ( m=0; m