/* * facetron.c * * Profile fitting for coherent nanocrystallography * * (c) 2006-2010 Thomas White * * Part of CrystFEL - crystallography with a FEL * */ #ifdef HAVE_CONFIG_H #include #endif #include #include #include #include #include #include #include #include #include #include #include #include "utils.h" #include "hdf5-file.h" #include "symmetry.h" #include "reflections.h" #include "stream.h" #include "geometry.h" #include "peaks.h" #include "thread-pool.h" #include "beam-parameters.h" /* Maximum number of iterations of NLSq to do for each image per macrocycle. */ #define MAX_CYCLES (100) /* Refineable parameters */ enum { REF_SCALE, NUM_PARAMS }; static void show_help(const char *s) { printf("Syntax: %s [options]\n\n", s); printf( "Post-refinement and profile fitting for coherent nanocrystallography.\n" "\n" " -h, --help Display this help message.\n" "\n" " -i, --input= Specify the name of the input 'stream'.\n" " (must be a file, not e.g. stdin)\n" " -o, --output= Output filename. Default: facetron.hkl.\n" " -g. --geometry= Get detector geometry from file.\n" " -b, --beam= Get beam parameters from file (provides initial\n" " values for parameters, and nominal wavelengths\n" " if no per-shot value is found in an HDF5 file.\n" " -x, --prefix=

Prefix filenames from input file with

.\n" " --basename Remove the directory parts of the filenames.\n" " --no-check-prefix Don't attempt to correct the --prefix.\n" " -y, --symmetry= Merge according to symmetry .\n" " -n, --iterations= Run cycles of post-refinement.\n" "\n" " -j Run analyses in parallel.\n"); } struct refine_args { const char *sym; ReflItemList *obs; double *i_full; struct image *image; FILE *graph; }; /* Return the gradient of parameter 'k' given the current status of 'image'. */ static double gradient(struct image *image, int k, struct cpeak spot, double I_partial) { switch ( k ) { case REF_SCALE : return I_partial; } ERROR("No gradient defined for parameter %i\n", k); abort(); } /* Apply the given shift to the 'k'th parameter of 'image'. */ static void apply_shift(struct image *image, int k, double shift) { switch ( k ) { case REF_SCALE : image->osf += shift; break; default : ERROR("No shift defined for parameter %i\n", k); abort(); } } static double mean_partial_dev(struct image *image, struct cpeak *spots, int n, const char *sym, double *i_full, FILE *graph) { int h; double delta_I = 0.0; for ( h=0; hosf; get_asymm(hind, kind, lind, &ha, &ka, &la, sym); I_full = lookup_intensity(i_full, ha, ka, la); delta_I += I_partial - spots[h].p * I_full; if ( graph != NULL ) { fprintf(graph, "%3i %3i %3i %5.2f (at %5.2f,%5.2f)" " %5.2f %5.2f\n", hind, kind, lind, I_partial/spots[h].p, xc, yc, spots[h].p, I_partial / I_full); } } return delta_I / (double)n; } static double iterate(struct image *image, double *i_full, const char *sym, struct cpeak *spots, int n) { gsl_matrix *M; gsl_vector *v; gsl_vector *shifts; int h, param; M = gsl_matrix_alloc(NUM_PARAMS, NUM_PARAMS); v = gsl_vector_alloc(NUM_PARAMS); /* Construct the equations, one per reflection in this image */ for ( h=0; hosf; get_asymm(hind, kind, lind, &ha, &ka, &la, sym); I_full = lookup_intensity(i_full, ha, ka, la); delta_I = I_partial - spots[h].p * I_full; for ( k=0; kosf); free(spots); spots = find_intersections(image, image->indexed_cell, &n, 0); return mean_partial_dev(image, spots, n, sym, i_full, NULL); } static void refine_image(int mytask, void *tasks) { struct refine_args *all_args = tasks; struct refine_args *pargs = &all_args[mytask]; struct image *image = pargs->image; double nominal_photon_energy = pargs->image->beam->photon_energy; struct hdfile *hdfile; struct cpeak *spots; int n, i; double dev, last_dev; hdfile = hdfile_open(image->filename); if ( hdfile == NULL ) { ERROR("Couldn't open '%s'\n", image->filename); return; } else if ( hdfile_set_image(hdfile, "/data/data0") ) { ERROR("Couldn't select path\n"); hdfile_close(hdfile); return; } if ( hdf5_read(hdfile, pargs->image, 0, nominal_photon_energy) ) { ERROR("Couldn't read '%s'\n", image->filename); hdfile_close(hdfile); return; } spots = find_intersections(image, image->indexed_cell, &n, 0); dev = +INFINITY; i = 0; do { last_dev = dev; dev = iterate(image, pargs->i_full, pargs->sym, spots, n); STATUS("Iteration %2i: mean dev = %5.2f\n", i, dev); i++; } while ( (fabs(last_dev - dev) > 1.0) || (i == MAX_CYCLES) ); mean_partial_dev(image, spots, n, pargs->sym, pargs->i_full, pargs->graph); free(image->data); if ( image->flags != NULL ) free(image->flags); hdfile_close(hdfile); /* Muppet proofing */ image->data = NULL; image->flags = NULL; } struct integrate_args { const char *sym; ReflItemList *obs; double *i_full; unsigned int *cts; pthread_mutex_t *list_lock; struct image *image; }; static void integrate_image(int mytask, void *tasks) { struct integrate_args *all_args = tasks; struct integrate_args *pargs = &all_args[mytask]; struct cpeak *spots; int j, n; struct hdfile *hdfile; struct image *image = pargs->image; double nominal_photon_energy = pargs->image->beam->photon_energy; hdfile = hdfile_open(image->filename); if ( hdfile == NULL ) { ERROR("Couldn't open '%s'\n", image->filename); return; } else if ( hdfile_set_image(hdfile, "/data/data0") ) { ERROR("Couldn't select path\n"); hdfile_close(hdfile); return; } if ( hdf5_read(hdfile, pargs->image, 0, nominal_photon_energy) ) { ERROR("Couldn't read '%s'\n", image->filename); hdfile_close(hdfile); return; } /* Figure out which spots should appear in this pattern */ spots = find_intersections(image, image->indexed_cell, &n, 0); /* For each reflection, estimate the partiality */ for ( j=0; josf; i_full_est = i_partial / spots[j].p; get_asymm(h, k, l, &ha, &ka, &la, pargs->sym); pthread_mutex_lock(pargs->list_lock); integrate_intensity(pargs->i_full, ha, ka, la, i_full_est); integrate_count(pargs->cts, ha, ka, la, 1); if ( !find_item(pargs->obs, ha, ka, la) ) { add_item(pargs->obs, ha, ka, la); } pthread_mutex_unlock(pargs->list_lock); } free(image->data); if ( image->flags != NULL ) free(image->flags); hdfile_close(hdfile); free(spots); /* Muppet proofing */ image->data = NULL; image->flags = NULL; } static void refine_all(struct image *images, int n_total_patterns, struct detector *det, const char *sym, ReflItemList *obs, double *i_full, int nthreads, FILE *graph) { struct refine_args *tasks; int i; tasks = malloc(n_total_patterns * sizeof(struct refine_args)); for ( i=0; ih, it->k, it->l); total /= lookup_count(cts, it->h, it->k, it->l); set_intensity(i_full, it->h, it->k, it->l, total); } } int main(int argc, char *argv[]) { int c; char *infile = NULL; char *outfile = NULL; char *geomfile = NULL; char *prefix = NULL; char *sym = NULL; FILE *fh; int nthreads = 1; int config_basename = 0; int config_checkprefix = 1; struct detector *det; double *i_full; unsigned int *cts; ReflItemList *obs; int i; int n_total_patterns; struct image *images; int n_iter = 10; struct beam_params *beam = NULL; /* Long options */ const struct option longopts[] = { {"help", 0, NULL, 'h'}, {"input", 1, NULL, 'i'}, {"output", 1, NULL, 'o'}, {"geometry", 1, NULL, 'g'}, {"beam", 1, NULL, 'b'}, {"prefix", 1, NULL, 'x'}, {"basename", 0, &config_basename, 1}, {"no-check-prefix", 0, &config_checkprefix, 0}, {"symmetry", 1, NULL, 'y'}, {"iterations", 1, NULL, 'n'}, {0, 0, NULL, 0} }; /* Short options */ while ((c = getopt_long(argc, argv, "hi:g:x:j:y:o:b:", longopts, NULL)) != -1) { switch (c) { case 'h' : show_help(argv[0]); return 0; case 'i' : infile = strdup(optarg); break; case 'g' : geomfile = strdup(optarg); break; case 'x' : prefix = strdup(optarg); break; case 'j' : nthreads = atoi(optarg); break; case 'y' : sym = strdup(optarg); break; case 'o' : outfile = strdup(optarg); break; case 'n' : n_iter = atoi(optarg); break; case 'b' : beam = get_beam_parameters(optarg); if ( beam == NULL ) { ERROR("Failed to load beam parameters" " from '%s'\n", optarg); return 1; } break; case 0 : break; default : return 1; } } /* Sanitise input filename and open */ if ( infile == NULL ) { infile = strdup("-"); } if ( strcmp(infile, "-") == 0 ) { fh = stdin; } else { fh = fopen(infile, "r"); } if ( fh == NULL ) { ERROR("Failed to open input file '%s'\n", infile); return 1; } free(infile); /* Sanitise output filename */ if ( outfile == NULL ) { outfile = strdup("facetron.hkl"); } /* Sanitise prefix */ if ( prefix == NULL ) { prefix = strdup(""); } else { if ( config_checkprefix ) { prefix = check_prefix(prefix); } } if ( sym == NULL ) sym = strdup("1"); /* Get detector geometry */ det = get_detector_geometry(geomfile); if ( det == NULL ) { ERROR("Failed to read detector geometry from '%s'\n", geomfile); return 1; } free(geomfile); if ( beam == NULL ) { ERROR("You must provide a beam parameters file.\n"); return 1; } /* Prepare for iteration */ i_full = new_list_intensity(); obs = new_items(); n_total_patterns = count_patterns(fh); STATUS("There are %i patterns to process\n", n_total_patterns); images = malloc(n_total_patterns * sizeof(struct image)); if ( images == NULL ) { ERROR("Couldn't allocate memory for images.\n"); return 1; } /* Fill in what we know about the images so far */ rewind(fh); for ( i=0; idivergence; images[i].bw = beam->bandwidth; images[i].det = det; images[i].beam = beam; images[i].osf = 1.0; /* Muppet proofing */ images[i].data = NULL; images[i].flags = NULL; free(filename); progress_bar(i, n_total_patterns-1, "Loading pattern data"); } fclose(fh); free(prefix); cts = new_list_count(); /* Make initial estimates */ estimate_full(images, n_total_patterns, det, sym, obs, i_full, cts, nthreads); /* Iterate */ for ( i=0; ipanels); free(det); free(beam); free(cts); for ( i=0; i