aboutsummaryrefslogtreecommitdiff
path: root/src/calibrate_detector.c
blob: 831fcec1e2cdd2a08287c71682f359596aea0744 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
/*
 * calibrate_detector.c
 *
 * Attempt to refine detector geometry
 *
 * (c) 2006-2011 Thomas White <taw@physics.org>
 *
 * Part of CrystFEL - crystallography with a FEL
 *
 */


#ifdef HAVE_CONFIG_H
#include <config.h>
#endif

#include <stdarg.h>
#include <stdlib.h>
#include <stdio.h>
#include <string.h>
#include <unistd.h>
#include <getopt.h>
#include <fenv.h>
#include <math.h>
#include <gsl/gsl_linalg.h>

#include "utils.h"
#include "image.h"
#include "detector.h"
#include "index.h"
#include "hdf5-file.h"
#include "stream.h"
#include "peaks.h"



static void show_help(const char *s)
{
	printf("Syntax: %s [options] -i <file.h5>\n\n", s);
	printf(
"Stream-based optimisation of detector geometry.\n"
"\n"
"  -h, --help                 Display this help message.\n"
"  -g. --geometry=<file>      Get detector geometry from file.\n"
"  -i, --input=<file>         Input filename.\n"
"  -m, --method=<method>      The calibration method.\n"
"               xy            Determine panel shifts in plane of detector\n"
"  -o, --output=<file>        Output results here"
"  -n, --npeaks=<number>      Don't refine unless this many peaks observed\n"
"                             (in the whole stream, not a single shot)\n"
"\n");
}



int main(int argc, char *argv[])
{


	char c;
	struct image image;
	struct panel p;
	char *filename = NULL;
	char *geometry = NULL;
	char *method = NULL;
	char *outputfile = NULL;
	FILE *fh = NULL;
	FILE *outfh = NULL;
	int nFeatures;
	int i;
	int fail;
	int nChunks;
	int minpeaks=0;	

	/* Long options */
	const struct option longopts[] = {
		{"help",               0, NULL,               'h'},
		{"input",              1, NULL,               'i'},
		{"geometry",           1, NULL,               'g'},
		{"method",             1, NULL,               'm'},
		{"output",             1, NULL,               'o'},
		{"npeaks",             1, NULL,               'n'},
		{0, 0, NULL, 0}
	};

	/* Short options */
	while ((c = getopt_long(argc, argv, "hi:g:m:o:n:", longopts, NULL)) != -1) {

		switch (c) {
		case 'h' :
			show_help(argv[0]);
			return 0;
		case 0 :
			break;
		case 'i' :
			filename = strdup(optarg);
			break;
		case 'g' :
			geometry = strdup(optarg);
			break;
		case 'm' :
			method = strdup(optarg);
			break;
		case 'o' :
			outputfile = strdup(optarg);
			break;
		case 'n' :
			minpeaks = atoi(optarg);
			break;
		default :
			return 1;
		}

	}

	if ( filename == NULL ) {
		ERROR("You must specify the input filename with -i\n");
		return 1;
	}
	fh = fopen(filename,"r");
   if ( fh == NULL ) {
		ERROR("Problem opening file\n");
		return 1;
	}
	printf("Read stream file: %s\n",filename);

	if ( geometry == NULL ) {
		ERROR("You need to specify a geometry file with --geometry\n");
		return 1;
	}	

	image.det = get_detector_geometry(geometry);
	if ( image.det == NULL ) {
		ERROR("Failed to read detector geometry from %s\n", geometry);
		return 1;
	}
	printf("Read geometry file: %s\n",geometry);
	image.width = image.det->max_fs;
	image.height = image.det->max_ss;
	free(geometry);

	if ( !(outputfile == NULL) ) {
		printf("Writing result to file: %s\n",outputfile);
      outfh = fopen(outputfile,"w");
   } else {
		ERROR("You did not specify an output file.\n");
		return 1;
	}
	if ( minpeaks == 0 ) {
		minpeaks = 1;
		printf("You did not specify minimum number of peaks."
             " Setting default value of %d\n",minpeaks);
	}
	
	if ( method == NULL ) {
		printf("You did not specify a refinement method-- setting default.\n");
		method = strdup("xy");
	}


	if ( !strcmp(method,"xy" ) ) {

		printf("Using refinement method %s\n",method);

		double * weightedSumFS, * weightedSumSS;
		double * summedWeights;
		double * meanShiftFS, * meanShiftSS;
		int * peaksFound;
		double cnx, cny;
		double xsh, ysh;

		weightedSumFS = (double *)calloc(sizeof(double), image.det->n_panels);
		weightedSumSS = (double *)calloc(sizeof(double), image.det->n_panels);
		summedWeights = (double *)calloc(sizeof(double), image.det->n_panels);
		peaksFound = (int *)calloc(sizeof(int), image.det->n_panels);
		meanShiftFS = (double *)calloc(sizeof(double), image.det->n_panels);
		meanShiftSS = (double *)calloc(sizeof(double), image.det->n_panels);

		/* initialize arrays (is there a standard function to do this?)  */
		int pi;
		for (pi=0; pi<image.det->n_panels; pi++) {
			weightedSumFS[pi] = 0;
	      weightedSumSS[pi] = 0;
   	   summedWeights[pi] = 0;
     		peaksFound[pi] = 0;
      	meanShiftFS[pi] = 0;
      	meanShiftSS[pi] = 0;
		}


		fesetround(1);  /* Round towards nearest */

		nChunks = 0;
		while (1) {

		   // check for peaks in next file
			fail = read_chunk(fh, &image);
		   if ( fail == 1 ) {
		      // FIXME:should check if this is EOF, or broken file handle
		      break;
		   }
			nChunks += 1;


			// move on to the next chunk if no peaks found
			if ( image.features == NULL ) {
				continue;
			}

			// now loop through peaks to determine mean panel shift
			nFeatures = image_feature_count(image.features);

			for (i=0; i<nFeatures; i++) {

				struct panel * p;
				struct imagefeature * thisFeature;			
				double ax, ay, az;
				double bx, by, bz;
				double cx, cy, cz;
				double hd, kd, ld;  /* Indices with decimal places */
				signed int h, k, l;
				struct rvec q;
				double twotheta;
				double fs, ss;   /* observed peaks */
				double pfs, pss; /* predicted peaks */
				double dfs, dss; /* observed - predicted */
				int pi;
				double thisWeight;

	
				/* if we find a feature, determine peak position */
				thisFeature = image_get_feature(image.features,i);	
				if ( thisFeature == NULL ) {
					continue;
				}

				fs = thisFeature->fs;
				ss = thisFeature->ss;

				p = find_panel(image.det, fs, ss);
				if ( p == NULL ) {
					continue;
				}
				if ( p->no_index ) continue;

				/* now determine the predicted peak position */

				/* scattering vector of this peak */
				q = get_q(&image, fs, ss, &twotheta, 1.0/image.lambda);

				/* miller indices of nearest Bragg reflection */
				cell_get_cartesian(image.indexed_cell, &ax, &ay, &az,
				                                       &bx, &by, &bz,
				                                       &cx, &cy, &cz);
			
				hd = q.u * ax + q.v * ay + q.w * az;
				kd = q.u * bx + q.v * by + q.w * bz;
				ld = q.u * cx + q.v * cy + q.w * cz;
				
				h = lrint(hd);
				k = lrint(kd);
				l = lrint(ld);

				/* now get scattering vector for reflectin [hkl]
				 * this means solving the equation U*q = h */

				double U[] = {ax, ay, az,
				              bx, by, bz,
				              cx, cy, cz};
     
				double hvec[] = {h,k,l};
     
				gsl_matrix_view m 
					= gsl_matrix_view_array (U, 3, 3);
     
				gsl_vector_view b
					= gsl_vector_view_array (hvec, 3);
     
				gsl_vector *x = gsl_vector_alloc (3);
       
				int s;
     
				gsl_permutation * perm = gsl_permutation_alloc (3);
				gsl_linalg_LU_decomp (&m.matrix, perm, &s);
				gsl_linalg_LU_solve (&m.matrix, perm, &b.vector, x);
   

				/* outgoing wavevector for [hkl] */
				double gx = x->data[0];
				double gy = x->data[1];
				double gz = x->data[2];

				double kk;
				double xd, yd, cl;
				double plx, ply;

				kk = 1/image.lambda;
				const double den = kk + gz;

				/* Camera length for this panel */
				cl = p->clen;

				/* Coordinates of peak relative to central beam, in m */
				xd = cl * gx / den;
				yd = cl * gy / den;
				
				/* Convert to pixels */
				xd *= p->res;
				yd *= p->res;
				
				/* Convert to relative to the panel corner */
				plx = xd - p->cnx;
				ply = yd - p->cny;
				
				pfs = p->xfs*plx + p->yfs*ply;
				pss = p->xss*plx + p->yss*ply;
				
				pfs += p->min_fs;
				pss += p->min_ss;

				/* Now, is this on this panel? */
				if ( fs < p->min_fs ) continue;
				if ( fs > p->max_fs ) continue;
				if ( ss < p->min_ss ) continue;
				if ( ss > p->max_ss ) continue;

				/* Finally, we have the shift in position of this peak */				
				dfs = pfs - fs;
				dss = pss - ss;

				/* Add this shift to the weighted sum over shifts */
				pi = find_panel_number(image.det,fs,ss);
				thisWeight = 1; // FIXME: use real weighting some day
				weightedSumFS[pi] += thisWeight*dfs;
				weightedSumSS[pi] += thisWeight*dss;
				summedWeights[pi] += thisWeight;
				peaksFound[pi] += 1;

			} /* end loop through image features */
			
		} /* end loop through stream chunks */
 


		/* calculate weighted average shift in peak positions */
		for (pi=0; pi < image.det->n_panels; pi++) {
			meanShiftFS[pi] = weightedSumFS[pi]/summedWeights[pi];
			meanShiftSS[pi] = weightedSumSS[pi]/summedWeights[pi];
		}
		
		/* now generate a new geometry file */
		for (pi=0; pi < image.det->n_panels; pi++) {
	
			p = image.det->panels[pi];

			xsh = 0;
			ysh = 0;
	
			if ( peaksFound[pi] >= minpeaks ) {

				//printf("meanShift: %f %f\n",meanShiftFS[pi],meanShiftSS[pi]);

				/* convert shifts from raw coords to lab frame */
				xsh = meanShiftFS[pi]*p.fsx + meanShiftSS[pi]*p.ssx;
				ysh = meanShiftFS[pi]*p.fsy + meanShiftSS[pi]*p.ssy;
		
				/* add shifts to original panel corner locations */
				cnx = p.cnx + xsh;
				cny = p.cny + ysh;
	
				//printf("new panel shift: %f %f\n",cnx,cny);

	
			} else {

				/* not refined: use original coordinates */
				cnx = p.cnx;
				cny = p.cny;

			}

			printf("panel %s, # peaks: %10d, mean shifts: %f %f\n",p.name, peaksFound[pi],xsh,ysh);

			//FIXME: there should be a function in geometry.c to write 
			//       these values to text file, since it will be useful on 
			//       other places as well (e.g. writing the geometry to 
			//       the data stream)
			fprintf(outfh,"%s/min_fs = %d\n",p.name,p.min_fs);
			fprintf(outfh,"%s/min_ss = %d\n",p.name,p.min_ss);
			fprintf(outfh,"%s/max_fs = %d\n",p.name,p.max_fs);
	      fprintf(outfh,"%s/max_ss = %d\n",p.name,p.max_ss);
			fprintf(outfh,"%s/badrow_direction = %C\n",p.name,p.badrow);
			fprintf(outfh,"%s/res = %g\n",p.name,p.res);
			fprintf(outfh,"%s/peak_sep = %g\n",p.name,p.peak_sep);
			fprintf(outfh,"%s/clen = %s\n",p.name,p.clen_from);
			//FIXME: the following is sketchy, but it will work for now.  we need
			//       to generalise the parser in detector.c
			char coord;
			char sign;
			if (p.fsx != 0){
				if (p.fsx>0){sign='+';}else{sign='-';}
				coord = 'x';
			} else {
				if (p.fsy>0){sign='+';}else{sign='-';}
				coord = 'y';
			}
			fprintf(outfh,"%s/fs = %C%C\n",p.name, sign, coord);
			if (p.ssx != 0){
            if (p.ssx>0){sign='+';}else{sign='-';}
            coord = 'x';
         } else {
            if (p.ssy>0){sign='+';}else{sign='-';}
            coord = 'y';
         }
			fprintf(outfh,"%s/ss = %C%C\n",p.name, sign, coord);
			fprintf(outfh,"%s/corner_x = %g\n",p.name,cnx);
	      fprintf(outfh,"%s/corner_y = %g\n",p.name,cny);
			if ( peaksFound[pi] < minpeaks ) {
				fprintf(outfh,"%s/no_index = %d\n",p.name,1);
			} else {
				fprintf(outfh,"%s/no_index = %d\n",p.name,p.no_index);
			}
			fprintf(outfh,"\n\n");
	
		}

	} else {

		printf("Refinement method %s not recognized\n",method);
   	return 1;

	}
	
	if ( !(outputfile == NULL) ) {
      fclose(outfh);
   }

	

	printf("Done!\n");

	return 0;
}