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
|
/*
* hrs-scaling.c
*
* Intensity scaling using generalised HRS target function
*
* (c) 2006-2010 Thomas White <taw@physics.org>
*
* Part of CrystFEL - crystallography with a FEL
*
*/
#ifdef HAVE_CONFIG_H
#include <config.h>
#endif
#include <stdlib.h>
#include <assert.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_vector.h>
#include <gsl/gsl_linalg.h>
#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; i<r; i++ ) {
STATUS("[ ");
for ( j=0; j<r; j++ ) {
STATUS("%+9.3e ", gsl_matrix_get(M, i, j));
}
STATUS("][ a%2i ] = [ %+9.3e ]\n", i, gsl_vector_get(v, i));
}
}
static void sum_GI(struct image *images, int n, const char *sym,
signed int hat, signed int kat, signed int lat,
double *sigma_GI, double *sigma_Gsq)
{
int k;
*sigma_GI = 0.0;
*sigma_Gsq = 0.0;
for ( k=0; k<n; k++ ) {
int hi;
struct image *image = &images[k];
struct cpeak *spots = images->cpeaks;
int found = 0;
for ( hi=0; hi<image->n_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; find<image->n_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 crossed = 0;
int n_ref;
M = gsl_matrix_calloc(n-1, n-1);
v = gsl_vector_calloc(n-1);
n_ref = num_items(obs);
for ( l=0; l<n; l++ ) { /* "Equation number": one equation per frame */
int m; /* Frame (scale factor) number */
int h;
int lcomp;
double vc_tot = 0.0;
struct image *imagel = &images[l];
if ( l == crossed ) continue;
/* Determine the "solution" vector component */
for ( h=0; h<n_ref; h++ ) {
double sigma_GI, sigma_Gsq;
double vc;
double Ihl;
struct refl_item *it = get_item(obs, h);
sum_GI(images, n, sym, it->h, 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;
}
lcomp = l;
if ( l > crossed ) lcomp--;
gsl_vector_set(v, lcomp, vc_tot);
/* Now fill in the matrix components */
for ( m=0; m<n; m++ ) {
double mc_tot = 0.0;
int mcomp;
struct image *imagem = &images[m];
if ( m > l ) continue; /* Matrix is symmetric */
if ( m == crossed ) continue;
for ( h=0; h<n_ref; h++ ) {
double mc = 0.0;
double Ihl, Ihm;
struct refl_item *it = get_item(obs, h);
double sigma_GI, sigma_Gsq;
sum_GI(images, n, sym, it->h, 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;
}
mcomp = m;
if ( m > crossed ) mcomp--;
gsl_matrix_set(M, lcomp, mcomp, mc_tot);
gsl_matrix_set(M, mcomp, lcomp, mc_tot);
}
}
show_matrix_eqn(M, v, n-1);
shifts = gsl_vector_alloc(n-1);
gsl_linalg_HH_solve(M, v, shifts);
max_shift = 0.0;
for ( l=0; l<n-1; l++ ) {
double shift = gsl_vector_get(shifts, l);
int limg;
limg = l;
if ( limg >= crossed ) limg++;
images[limg].osf += shift;
if ( fabs(shift) > 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; i<num_items(obs); i++ ) {
signed int h, k, l;
struct refl_item *it = get_item(obs, i);
double num = 0.0;
double den = 0.0;
int m;
get_asymm(it->h, it->k, it->l, &h, &k, &l, sym);
/* For each frame */
for ( m=0; m<n; m++ ) {
double G;
int a;
G = images[m].osf;
/* For each peak */
for ( a=0; a<images[m].n_cpeaks; a++ ) {
signed int ha, ka, la;
if ( !images[m].cpeaks[a].valid ) continue;
if ( images[m].cpeaks[a].p < 0.1 ) continue;
/* Correct reflection? */
get_asymm(images[m].cpeaks[a].h,
images[m].cpeaks[a].k,
images[m].cpeaks[a].l,
&ha, &ka, &la, sym);
if ( ha != h ) continue;
if ( ka != k ) continue;
if ( la != l ) continue;
num += images[m].cpeaks[a].intensity
* images[m].cpeaks[a].p * G;
den += pow(images[m].cpeaks[a].p, 2.0)
* pow(G, 2.0);
}
}
set_intensity(I_full, h, k, l, num/den);
}
return I_full;
}
static void normalise_osfs(struct image *images, int n)
{
int m;
double tot = 0.0;
double mean;
for ( m=0; m<n; m++ ) {
tot += images[m].osf;
}
mean = tot / (double)n;
for ( m=0; m<n; m++ ) {
images[m].osf /= mean;
}
}
/* Scale the stack of images */
double *scale_intensities(struct image *images, int n, const char *sym,
ReflItemList *obs)
{
int m;
double *I_full;
int i;
double max_shift;
/* Start with all scale factors equal */
for ( m=0; m<n; m++ ) {
double tot = 0.0;
int j;
for ( j=0; j<images[m].n_cpeaks; j++ ) {
tot += images[m].cpeaks[j].intensity
/ images[m].cpeaks[j].p;
}
images[m].osf = tot;
}
normalise_osfs(images, n);
/* Iterate */
i = 0;
do {
max_shift = iterate_scale(images, n, obs, sym);
STATUS("Iteration %2i: max shift = %5.2f\n", i, max_shift);
i++;
normalise_osfs(images, n);
} while ( (max_shift > 0.01) && (i < MAX_CYCLES) );
for ( m=0; m<n; m++ ) {
images[m].osf /= images[0].osf;
}
I_full = lsq_intensities(images, n, obs, sym);
return I_full;
}
|