aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
-rw-r--r--libcrystfel/src/image.h7
-rw-r--r--libcrystfel/src/index.c27
-rw-r--r--libcrystfel/src/indexers/asdf.c11
-rw-r--r--libcrystfel/src/indexers/dirax.c7
-rw-r--r--libcrystfel/src/indexers/felix.c11
-rw-r--r--libcrystfel/src/indexers/mosflm.c13
-rw-r--r--libcrystfel/src/indexers/taketwo.c13
-rw-r--r--libcrystfel/src/indexers/xds.c13
-rw-r--r--libcrystfel/src/indexers/xgandalf.c12
-rw-r--r--libcrystfel/src/integration.c30
-rw-r--r--libcrystfel/src/predict-refine.c11
11 files changed, 85 insertions, 70 deletions
diff --git a/libcrystfel/src/image.h b/libcrystfel/src/image.h
index 91a0afc4..fb33f2d9 100644
--- a/libcrystfel/src/image.h
+++ b/libcrystfel/src/image.h
@@ -68,13 +68,6 @@ struct imagefeature {
int pn; /**< Panel number */
double intensity; /**< Intensity */
- /** \name Reciprocal space coordinates (m^-1) of this feature */
- /** @{ */
- double rx;
- double ry;
- double rz;
- /** @} */
-
const char *name; /**< Text name, e.g. "5,3,-1" */
};
diff --git a/libcrystfel/src/index.c b/libcrystfel/src/index.c
index 457fd7fe..130b4262 100644
--- a/libcrystfel/src/index.c
+++ b/libcrystfel/src/index.c
@@ -563,32 +563,6 @@ void cleanup_indexing(IndexingPrivate *ipriv)
}
-void map_all_peaks(struct image *image)
-{
- int i, n;
-
- n = image_feature_count(image->features);
-
- /* Map positions to 3D */
- for ( i=0; i<n; i++ ) {
-
- struct imagefeature *f;
- double r[3];
-
- f = image_get_feature(image->features, i);
- if ( f == NULL ) continue;
-
- detgeom_transform_coords(&image->detgeom->panels[f->pn],
- f->fs, f->ss,
- image->lambda, r);
- f->rx = r[0];
- f->ry = r[1];
- f->rz = r[2];
-
- }
-}
-
-
/* Return 0 for cell OK, 1 for cell incorrect */
static int check_cell(IndexingFlags flags, Crystal *cr, UnitCell *target,
double *tolerance)
@@ -953,7 +927,6 @@ void index_pattern_3(struct image *image, IndexingPrivate *ipriv, int *ping,
if ( ipriv == NULL ) return;
- map_all_peaks(image);
image->crystals = NULL;
image->n_crystals = 0;
diff --git a/libcrystfel/src/indexers/asdf.c b/libcrystfel/src/indexers/asdf.c
index 806501cc..68b714a9 100644
--- a/libcrystfel/src/indexers/asdf.c
+++ b/libcrystfel/src/indexers/asdf.c
@@ -1116,14 +1116,19 @@ int run_asdf(struct image *image, void *ipriv)
for ( i=0; i<n; i++ ) {
struct imagefeature *f;
+ double r[3];
f = image_get_feature(image->features, i);
if ( f == NULL ) continue;
+ detgeom_transform_coords(&image->detgeom->panels[f->pn],
+ f->fs, f->ss, image->lambda,
+ r);
+
reflections[N_reflections] = gsl_vector_alloc(3);
- gsl_vector_set(reflections[N_reflections], 0, f->rx/1e10);
- gsl_vector_set(reflections[N_reflections], 1, f->ry/1e10);
- gsl_vector_set(reflections[N_reflections], 2, f->rz/1e10);
+ gsl_vector_set(reflections[N_reflections], 0, r[0]/1e10);
+ gsl_vector_set(reflections[N_reflections], 1, r[1]/1e10);
+ gsl_vector_set(reflections[N_reflections], 2, r[2]/1e10);
N_reflections++;
}
diff --git a/libcrystfel/src/indexers/dirax.c b/libcrystfel/src/indexers/dirax.c
index 24be87ba..a21dd7e2 100644
--- a/libcrystfel/src/indexers/dirax.c
+++ b/libcrystfel/src/indexers/dirax.c
@@ -490,12 +490,17 @@ static void write_drx(struct image *image)
for ( i=0; i<image_feature_count(image->features); i++ ) {
struct imagefeature *f;
+ double r[3];
f = image_get_feature(image->features, i);
if ( f == NULL ) continue;
+ detgeom_transform_coords(&image->detgeom->panels[f->pn],
+ f->fs, f->ss, image->lambda,
+ r);
+
fprintf(fh, "%10f %10f %10f %8f\n",
- f->rx/1e10, f->ry/1e10, f->rz/1e10, 1.0);
+ r[0]/1e10, r[1]/1e10, r[2]/1e10, 1.0);
}
fclose(fh);
diff --git a/libcrystfel/src/indexers/felix.c b/libcrystfel/src/indexers/felix.c
index b94227d0..0db5ecbb 100644
--- a/libcrystfel/src/indexers/felix.c
+++ b/libcrystfel/src/indexers/felix.c
@@ -364,14 +364,19 @@ static void write_gve(struct image *image, struct felix_private *gp)
for ( i=0; i<image_feature_count(image->features); i++ ) {
struct imagefeature *f;
+ double r[3];
f = image_get_feature(image->features, i);
if ( f == NULL ) continue;
+ detgeom_transform_coords(&image->detgeom->panels[f->pn],
+ f->fs, f->ss, image->lambda,
+ r);
+
fprintf(fh, "%.6f %.6f %.6f 0 0 %.6f %.6f %.6f 0\n",
- f->rz/1e10, f->rx/1e10, f->ry/1e10,
- modulus(f->rx, f->ry, f->rz)/1e10, /* dstar */
- rad2deg(atan2(f->ry, f->rx)), 0.0); /* eta, omega */
+ r[2]/1e10, r[0]/1e10, r[1]/1e10,
+ modulus(r[0], r[1], r[2])/1e10, /* dstar */
+ rad2deg(atan2(r[1], r[0])), 0.0); /* eta, omega */
}
fclose(fh);
diff --git a/libcrystfel/src/indexers/mosflm.c b/libcrystfel/src/indexers/mosflm.c
index bacd345f..4c1d2906 100644
--- a/libcrystfel/src/indexers/mosflm.c
+++ b/libcrystfel/src/indexers/mosflm.c
@@ -357,16 +357,21 @@ static void write_spt(struct image *image, const char *filename)
struct imagefeature *f;
double ttx, tty, x, y;
+ double r[3];
f = image_get_feature(image->features, i);
if ( f == NULL ) continue;
+ detgeom_transform_coords(&image->detgeom->panels[f->pn],
+ f->fs, f->ss, image->lambda,
+ r);
+
ttx = angle_between_2d(0.0, 1.0,
- f->rx, 1.0/image->lambda + f->rz);
+ r[0], 1.0/image->lambda + r[2]);
tty = angle_between_2d(0.0, 1.0,
- f->ry, 1.0/image->lambda + f->rz);
- if ( f->rx < 0.0 ) ttx *= -1.0;
- if ( f->ry < 0.0 ) tty *= -1.0;
+ r[1], 1.0/image->lambda + r[2]);
+ if ( r[0] < 0.0 ) ttx *= -1.0;
+ if ( r[1] < 0.0 ) tty *= -1.0;
x = -tan(ttx)*FAKE_CLEN;
y = tan(tty)*FAKE_CLEN;
diff --git a/libcrystfel/src/indexers/taketwo.c b/libcrystfel/src/indexers/taketwo.c
index 4243adcb..18845ce5 100644
--- a/libcrystfel/src/indexers/taketwo.c
+++ b/libcrystfel/src/indexers/taketwo.c
@@ -2157,11 +2157,18 @@ int taketwo_index(struct image *image, void *priv)
rlps = malloc((image_feature_count(image->features)+1)*sizeof(struct rvec));
for ( i=0; i<image_feature_count(image->features); i++ ) {
+
+ double r[3];
struct imagefeature *pk = image_get_feature(image->features, i);
if ( pk == NULL ) continue;
- rlps[n_rlps].u = pk->rx;
- rlps[n_rlps].v = pk->ry;
- rlps[n_rlps].w = pk->rz;
+
+ detgeom_transform_coords(&image->detgeom->panels[pk->pn],
+ pk->fs, pk->ss, image->lambda,
+ r);
+
+ rlps[n_rlps].u = r[0];
+ rlps[n_rlps].v = r[1];
+ rlps[n_rlps].w = r[2];
n_rlps++;
}
rlps[n_rlps].u = 0.0;
diff --git a/libcrystfel/src/indexers/xds.c b/libcrystfel/src/indexers/xds.c
index 02610d5e..76c657c6 100644
--- a/libcrystfel/src/indexers/xds.c
+++ b/libcrystfel/src/indexers/xds.c
@@ -213,17 +213,22 @@ static void write_spot(struct image *image)
{
struct imagefeature *f;
double ttx, tty, x, y;
+ double r[3];
f = image_get_feature(image->features, i);
if ( f == NULL ) continue;
if ( f->intensity <= 0 ) continue;
+ detgeom_transform_coords(&image->detgeom->panels[f->pn],
+ f->fs, f->ss, image->lambda,
+ r);
+
ttx = angle_between_2d(0.0, 1.0,
- f->rx, 1.0/image->lambda + f->rz);
+ r[0], 1.0/image->lambda + r[2]);
tty = angle_between_2d(0.0, 1.0,
- f->ry, 1.0/image->lambda + f->rz);
- if ( f->rx < 0.0 ) ttx *= -1.0;
- if ( f->ry < 0.0 ) tty *= -1.0;
+ r[1], 1.0/image->lambda + r[2]);
+ if ( r[0] < 0.0 ) ttx *= -1.0;
+ if ( r[1] < 0.0 ) tty *= -1.0;
x = tan(ttx)*FAKE_CLEN;
y = tan(tty)*FAKE_CLEN;
diff --git a/libcrystfel/src/indexers/xgandalf.c b/libcrystfel/src/indexers/xgandalf.c
index a0ce768d..a4773e18 100644
--- a/libcrystfel/src/indexers/xgandalf.c
+++ b/libcrystfel/src/indexers/xgandalf.c
@@ -92,14 +92,20 @@ int run_xgandalf(struct image *image, void *ipriv)
reciprocalPeaks_1_per_A->peakCount = 0;
for ( i = 0; i < peakCountMax && i < MAX_PEAK_COUNT_FOR_INDEXER; i++) {
struct imagefeature *f;
+ double r[3];
+
f = image_get_feature(image->features, i);
if (f == NULL) {
continue;
}
- reciprocalPeaks_1_per_A->coordinates_x[reciprocalPeaks_1_per_A->peakCount] = f->rx * 1e-10;
- reciprocalPeaks_1_per_A->coordinates_y[reciprocalPeaks_1_per_A->peakCount] = f->ry * 1e-10;
- reciprocalPeaks_1_per_A->coordinates_z[reciprocalPeaks_1_per_A->peakCount] = f->rz * 1e-10;
+ detgeom_transform_coords(&image->detgeom->panels[f->pn],
+ f->fs, f->ss, image->lambda,
+ r);
+
+ reciprocalPeaks_1_per_A->coordinates_x[reciprocalPeaks_1_per_A->peakCount] = r[0] * 1e-10;
+ reciprocalPeaks_1_per_A->coordinates_y[reciprocalPeaks_1_per_A->peakCount] = r[1] * 1e-10;
+ reciprocalPeaks_1_per_A->coordinates_z[reciprocalPeaks_1_per_A->peakCount] = r[2] * 1e-10;
reciprocalPeaks_1_per_A->peakCount++;
}
diff --git a/libcrystfel/src/integration.c b/libcrystfel/src/integration.c
index 4a66363b..c15d66ba 100644
--- a/libcrystfel/src/integration.c
+++ b/libcrystfel/src/integration.c
@@ -1533,7 +1533,7 @@ int integrate_rings_once(Reflection *refl,
}
-static double estimate_resolution(UnitCell *cell, ImageFeatureList *flist)
+static double estimate_resolution(UnitCell *cell, struct image *image)
{
int i;
const double min_dist = 0.25;
@@ -1549,27 +1549,33 @@ static double estimate_resolution(UnitCell *cell, ImageFeatureList *flist)
return INFINITY;
}
- for ( i=0; i<image_feature_count(flist); i++ ) {
+ for ( i=0; i<image_feature_count(image->features); i++ ) {
struct imagefeature *f;
double h, k, l, hd, kd, ld;
-
- /* Assume all image "features" are genuine peaks */
- f = image_get_feature(flist, i);
- if ( f == NULL ) continue;
-
double ax, ay, az;
double bx, by, bz;
double cx, cy, cz;
+ double r[3];
+
+ /* Assume all image "features" are genuine peaks */
+ f = image_get_feature(image->features, i);
+ if ( f == NULL ) continue;
cell_get_cartesian(cell,
- &ax, &ay, &az, &bx, &by, &bz, &cx, &cy, &cz);
+ &ax, &ay, &az,
+ &bx, &by, &bz,
+ &cx, &cy, &cz);
+
+ detgeom_transform_coords(&image->detgeom->panels[f->pn],
+ f->fs, f->ss, image->lambda,
+ r);
/* 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;
+ hd = r[0] * ax + r[1] * ay + r[2] * az;
+ kd = r[0] * bx + r[1] * by + r[2] * bz;
+ ld = r[0] * cx + r[1] * cy + r[2] * cz;
h = lrint(hd);
k = lrint(kd);
l = lrint(ld);
@@ -1680,7 +1686,7 @@ void integrate_all_5(struct image *image, IntegrationMethod meth,
}
res = estimate_resolution(crystal_get_cell(image->crystals[i]),
- image->features);
+ image);
crystal_set_resolution_limit(image->crystals[i], res);
list = predict_to_res(image->crystals[i], res+push_res);
diff --git a/libcrystfel/src/predict-refine.c b/libcrystfel/src/predict-refine.c
index 43a52cae..1c9e8504 100644
--- a/libcrystfel/src/predict-refine.c
+++ b/libcrystfel/src/predict-refine.c
@@ -191,16 +191,21 @@ static int pair_peaks(struct image *image, Crystal *cr,
struct imagefeature *f;
double h, k, l, hd, kd, ld;
Reflection *refl;
+ double r[3];
/* Assume all image "features" are genuine peaks */
f = image_get_feature(image->features, i);
if ( f == NULL ) continue;
+ detgeom_transform_coords(&image->detgeom->panels[f->pn],
+ f->fs, f->ss, image->lambda,
+ r);
+
/* 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;
+ hd = r[0] * ax + r[1] * ay + r[2] * az;
+ kd = r[0] * bx + r[1] * by + r[2] * bz;
+ ld = r[0] * cx + r[1] * cy + r[2] * cz;
h = lrint(hd);
k = lrint(kd);
l = lrint(ld);