/* * predict-refine.c * * Prediction refinement * * Copyright © 2012-2020 Deutsches Elektronen-Synchrotron DESY, * a research centre of the Helmholtz Association. * * Authors: * 2010-2016 Thomas White * 2016 Valerio Mariani * * This file is part of CrystFEL. * * CrystFEL is free software: you can redistribute it and/or modify * it under the terms of the GNU General Public License as published by * the Free Software Foundation, either version 3 of the License, or * (at your option) any later version. * * CrystFEL is distributed in the hope that it will be useful, * but WITHOUT ANY WARRANTY; without even the implied warranty of * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * GNU General Public License for more details. * * You should have received a copy of the GNU General Public License * along with CrystFEL. If not, see . * */ #ifdef HAVE_CONFIG_H #include #endif #include #include #include #include #include "image.h" #include "geometry.h" #include "cell-utils.h" /** \file predict-refine.h */ /* Maximum number of iterations of NLSq to do for each image per macrocycle. */ #define MAX_CYCLES (10) /* Weighting of excitation error term (m^-1) compared to position term (m) */ #define EXC_WEIGHT (4e-20) /* Parameters to refine */ static const enum gparam rv[] = { GPARAM_ASX, GPARAM_ASY, GPARAM_ASZ, GPARAM_BSX, GPARAM_BSY, GPARAM_BSZ, GPARAM_CSX, GPARAM_CSY, GPARAM_CSZ, GPARAM_DETX, GPARAM_DETY, }; static const int num_params = 11; struct reflpeak { Reflection *refl; struct imagefeature *peak; double Ih; /* normalised */ struct panel *panel; /* panel the reflection appears on * (we assume this never changes) */ }; static void twod_mapping(double fs, double ss, double *px, double *py, struct panel *p) { double xs, ys; xs = fs*p->fsx + ss*p->ssx; ys = fs*p->fsy + ss*p->ssy; *px = (xs + p->cnx) / p->res; *py = (ys + p->cny) / p->res; } static double r_dev(struct reflpeak *rp) { /* Excitation error term */ return get_exerr(rp->refl); } static double x_dev(struct reflpeak *rp, struct detector *det) { /* Peak position term */ double xpk, ypk, xh, yh; double fsh, ssh; twod_mapping(rp->peak->fs, rp->peak->ss, &xpk, &ypk, rp->panel); get_detector_pos(rp->refl, &fsh, &ssh); twod_mapping(fsh, ssh, &xh, &yh, rp->panel); return xh-xpk; } static double y_dev(struct reflpeak *rp, struct detector *det) { /* Peak position term */ double xpk, ypk, xh, yh; double fsh, ssh; twod_mapping(rp->peak->fs, rp->peak->ss, &xpk, &ypk, rp->panel); get_detector_pos(rp->refl, &fsh, &ssh); twod_mapping(fsh, ssh, &xh, &yh, rp->panel); return yh-ypk; } static void UNUSED write_pairs(const char *filename, struct reflpeak *rps, int n, struct detector *det) { int i; FILE *fh; fh = fopen(filename, "w"); if ( fh == NULL ) { ERROR("Failed to open '%s'\n", filename); return; } for ( i=0; ifs; ss = rps[i].peak->ss; p = rps[i].panel; get_indices(rps[i].refl, &h, &k, &l); write_fs = fs + p->orig_min_fs; write_ss = ss + p->orig_min_ss; fprintf(fh, "%7.2f %7.2f dev r,x,y: %9f %9f %9f %9f\n", write_fs, write_ss, r_dev(&rps[i])/1e9, fabs(r_dev(&rps[i])/1e9), x_dev(&rps[i], det), y_dev(&rps[i], det)); //fprintf(fh, "%4i %4i %4i 0.0 - 0.0 1 %7.2f %7.2f %s\n", // h, k, l, write_fs, write_ss, p->name); } fclose(fh); STATUS("Wrote %i pairs to %s\n", n, filename); } static int cmpd2(const void *av, const void *bv) { struct reflpeak *a, *b; a = (struct reflpeak *)av; b = (struct reflpeak *)bv; if ( fabs(r_dev(a)) < fabs(r_dev(b)) ) return -1; return 1; } static int check_outlier_transition(struct reflpeak *rps, int n, struct detector *det) { int i; if ( n < 3 ) return n; qsort(rps, n, sizeof(struct reflpeak), cmpd2); //write_pairs("pairs-before-outlier.lst", rps, n, det); for ( i=1; ifeatures); i++ ) { struct imagefeature *f; double h, k, l, hd, kd, ld; Reflection *refl; /* Assume all image "features" are genuine peaks */ f = image_get_feature(image->features, i); if ( f == NULL ) continue; /* Decimal and fractional Miller indices of nearest reciprocal * lattice point */ hd = f->rx * ax + f->ry * ay + f->rz * az; kd = f->rx * bx + f->ry * by + f->rz * bz; ld = f->rx * cx + f->ry * cy + f->rz * cz; h = lrint(hd); k = lrint(kd); l = lrint(ld); /* Don't pair with 000, because that can cause trouble later */ if ( (h==0) && (k==0) && (l==0) ) continue; refl = reflection_new(h, k, l); if ( refl == NULL ) { ERROR("Failed to create reflection\n"); return 0; } add_refl_to_list(refl, all_reflist); set_symmetric_indices(refl, h, k, l); /* It doesn't matter if the actual predicted location * doesn't fall on this panel. We're only interested * in how far away it is from the peak location. * The predicted position and excitation errors will be * filled in by update_predictions(). */ set_panel(refl, f->p); rps[n].refl = refl; rps[n].peak = f; rps[n].panel = f->p; n++; } /* Get the excitation errors and detector positions for the candidate * reflections */ crystal_set_reflections(cr, all_reflist); update_predictions(cr); /* Pass over the peaks again, keeping only the ones which look like * good pairings */ for ( i=0; ifs, 2.0) + pow(ss - rps[i].peak->ss, 2.0); if ( pd > 10.0 * 10.0 ) continue; /* FIXME Hardcoded distance */ rps[n_acc] = rps[i]; rps[n_acc].refl = reflection_new(h, k, l); copy_data(rps[n_acc].refl, refl); n_acc++; } reflist_free(all_reflist); /* Sort the pairings by excitation error and look for a transition * between good pairings and outliers */ n_final = check_outlier_transition(rps, n_acc, image->det); /* Add the final accepted reflections to the caller's list */ if ( reflist != NULL ) { for ( i=0; ifeatures) * sizeof(struct reflpeak)); if ( rps == NULL ) return 1; reflist = reflist_new(); n_acc = pair_peaks(image, cr, reflist, rps); if ( n_acc < 3 ) { free(rps); reflist_free(reflist); return 1; } crystal_set_reflections(cr, reflist); update_predictions(cr); crystal_set_reflections(cr, NULL); qsort(rps, n_acc, sizeof(struct reflpeak), cmpd2); n = (n_acc-1) - n_acc/50; if ( n < 2 ) n = 2; /* n_acc is always >= 2 */ crystal_set_profile_radius(cr, fabs(r_dev(&rps[n]))); reflist_free(reflist); free(rps); return 0; } static void update_detector(struct detector *det, double xoffs, double yoffs, double coffs) { int i; for ( i=0; in_panels; i++ ) { struct panel *p = &det->panels[i]; p->cnx += xoffs * p->res; p->cny += yoffs * p->res; p->clen += coffs; } } static int iterate(struct reflpeak *rps, int n, UnitCell *cell, struct image *image, double *total_x, double *total_y, double *total_z) { int i; gsl_matrix *M; gsl_vector *v; gsl_vector *shifts; double asx, asy, asz; double bsx, bsy, bsz; double csx, csy, csz; /* Number of parameters to refine */ M = gsl_matrix_calloc(num_params, num_params); v = gsl_vector_calloc(num_params); for ( i=0; i k ) continue; M_c = w * gradients[g] * gradients[k]; M_curr = gsl_matrix_get(M, k, g); gsl_matrix_set(M, k, g, M_curr + M_c); gsl_matrix_set(M, g, k, M_curr + M_c); } v_c = w * r_dev(&rps[i]); v_c *= -gradients[k]; v_curr = gsl_vector_get(v, k); gsl_vector_set(v, k, v_curr + v_c); } /* Positional x terms */ for ( k=0; k k ) continue; M_c = gradients[g] * gradients[k]; M_curr = gsl_matrix_get(M, k, g); gsl_matrix_set(M, k, g, M_curr + M_c); gsl_matrix_set(M, g, k, M_curr + M_c); } v_c = x_dev(&rps[i], image->det); v_c *= -gradients[k]; v_curr = gsl_vector_get(v, k); gsl_vector_set(v, k, v_curr + v_c); } /* Positional y terms */ for ( k=0; k k ) continue; M_c = gradients[g] * gradients[k]; M_curr = gsl_matrix_get(M, k, g); gsl_matrix_set(M, k, g, M_curr + M_c); gsl_matrix_set(M, g, k, M_curr + M_c); } v_c = y_dev(&rps[i], image->det); v_c *= -gradients[k]; v_curr = gsl_vector_get(v, k); gsl_vector_set(v, k, v_curr + v_c); } } int k; for ( k=0; kdet, gsl_vector_get(shifts, 9), gsl_vector_get(shifts, 10), gsl_vector_get(shifts, 11)); *total_x += gsl_vector_get(shifts, 9); *total_y += gsl_vector_get(shifts, 10); *total_z += gsl_vector_get(shifts, 11); cell_set_reciprocal(cell, asx, asy, asz, bsx, bsy, bsz, csx, csy, csz); gsl_vector_free(shifts); gsl_matrix_free(M); gsl_vector_free(v); return 0; } static double pred_residual(struct reflpeak *rps, int n, struct detector *det) { int i; double res = 0.0; double r; r = 0.0; for ( i=0; ifeatures) * sizeof(struct reflpeak)); if ( rps == NULL ) return 1; reflist = reflist_new(); n = pair_peaks(image, cr, reflist, rps); if ( n < 10 ) { free(rps); reflist_free(reflist); return 1; } crystal_set_reflections(cr, reflist); /* Normalise the intensities to max 1 */ max_I = -INFINITY; for ( i=0; iintensity; if ( cur_I > max_I ) max_I = cur_I; } if ( max_I <= 0.0 ) { ERROR("All peaks negative?\n"); free(rps); return 1; } for ( i=0; iintensity / max_I; } //STATUS("Initial residual = %e\n", pred_residual(rps, n, image->det)); /* Refine */ for ( i=0; idet)); } //STATUS("Final residual = %e\n", pred_residual(rps, n, image->det)); snprintf(tmp, 255, "predict_refine/final_residual = %e", pred_residual(rps, n, image->det)); crystal_add_notes(cr, tmp); crystal_set_det_shift(cr, total_x, total_y); crystal_set_reflections(cr, NULL); reflist_free(reflist); n = pair_peaks(image, cr, NULL, rps); free_rps_noreflist(rps, n); if ( n < 10 ) { return 1; } return 0; }