dials.algorithms.indexing¶
dials.algorithms.indexing.indexer¶
dials.algorithms.indexing.stills_indexer¶
- class dials.algorithms.indexing.stills_indexer.StillsIndexer(reflections, experiments, params=None)[source]¶
Bases:
dials.algorithms.indexing.indexer.Indexer
Class for indexing stills
- class dials.algorithms.indexing.stills_indexer.StillsIndexerBasisVectorSearch(reflections, experiments, params=None)[source]¶
Bases:
dials.algorithms.indexing.stills_indexer.StillsIndexer
,dials.algorithms.indexing.lattice_search.BasisVectorSearch
- class dials.algorithms.indexing.stills_indexer.StillsIndexerKnownOrientation(reflections, experiments, params, known_orientations)[source]¶
Bases:
dials.algorithms.indexing.known_orientation.IndexerKnownOrientation
,dials.algorithms.indexing.stills_indexer.StillsIndexer
- class dials.algorithms.indexing.stills_indexer.StillsIndexerLatticeSearch(reflections, experiments, params=None)[source]¶
Bases:
dials.algorithms.indexing.stills_indexer.StillsIndexer
,dials.algorithms.indexing.lattice_search.LatticeSearch
dials.algorithms.indexing.basis_vector_search¶
- class dials.algorithms.indexing.basis_vector_search.FFT1D(max_cell, params=None, *args, **kwargs)[source]¶
Bases:
dials.algorithms.indexing.basis_vector_search.strategy.Strategy
Basis vector search using 1D FFTs in reciprocal space.
A set of dimensionless radial unit vectors, typically ~7000 in total, is chosen so that they are roughly evenly spaced in solid angle over a hemisphere. The reciprocal space displacements of the measured spot centroids are then projected onto each of these radial vectors in turn (that is, we calculate the scalar product of each displacement with each unit vector). A 1D FFT of the linear density of projected spot positions is performed along each direction. Aggregating the results of all the transforms, the three shortest non-collinear wave vectors with the greatest spectral weight correspond to the basis vectors of the direct lattice.
- See:
Steller, I., Bolotovsky, R. & Rossmann, M. G. (1997). J. Appl. Cryst. 30, 1036-1040. Sauter, N. K., Grosse-Kunstleve, R. W. & Adams, P. D. (2004). J. Appl. Cryst. 37, 399-409.
- find_basis_vectors(reciprocal_lattice_vectors)[source]¶
Find a list of likely basis vectors.
- Parameters
reciprocal_lattice_vectors (scitbx.array_family.flex.vec3_double) – The list of reciprocal lattice vectors to search for periodicity.
- Returns
A tuple containing the list of basis vectors and a flex.bool array identifying which reflections were used in indexing.
- phil_help = 'Search for the basis vectors of the direct lattice by performing a series of 1D FFTs along various directions in reciprocal space. This has a lower memory requirement than a single 3D FFT (the fft3d method). This method may also be more appropriate than a 3D FFT if the reflections are from narrow wedges of rotation data or from stills data.'¶
- phil_scope = <libtbx.phil.scope object>¶
- class dials.algorithms.indexing.basis_vector_search.FFT3D(max_cell, min_cell=3, params=None, *args, **kwargs)[source]¶
Bases:
dials.algorithms.indexing.basis_vector_search.strategy.Strategy
Basis vector search using a 3D FFT in reciprocal space.
Reciprocal space is sampled as a 3D Cartesian grid, aligned with the basis of the laboratory frame. The reciprocal-space positions of the centroids of measured spots are ascribed a value of 1, with the rest of the grid assigned a value of 0. A 3D FFT is performed and the three shortest non-collinear reciprocal spatial wave vectors with appreciable spectral weight correspond to the basis vectors of the real space lattice.
Because this procedure requires a sampling of all of reciprocal space, up to the d* value of the measured spot with the highest resolution, it can be more memory intensive than alternative approaches. To mitigate this, the 3D FFT will sometimes be curtailed to a region of reciprocal space below a certain resolution, and higher-resolution spots will be ignored.
- See:
Bricogne, G. (1986). Proceedings of the EEC Cooperative Workshop on Position-Sensitive Detector Software (Phase III), p. 28. Paris: LURE. Campbell, J. W. (1998). J. Appl. Cryst. 31, 407-413.
- __init__(max_cell, min_cell=3, params=None, *args, **kwargs)[source]¶
Construct an FFT3D object.
- Parameters
max_cell (float) – An estimate of the maximum cell dimension of the primitive cell.
n_points (int) – The size of the fft3d grid.
d_min (float) – The high resolution limit in Angstrom for spots to include in the initial indexing. If Auto then calculated as d_min = 5 * max_cell/n_points.
b_iso (float) – Apply an isotropic b_factor weight to the points when doing the FFT. If Auto then calculated as b_iso = -4 * d_min ** 2 * math.log(0.05).
rmsd_cutoff (float) – RMSD cutoff applied to the transformed real-space map prior to performing the peak search.
peak_volume_cutoff (float) – Only include peaks that are larger than this fraction of the volume of the largest peak in the transformed real-space map.
min_cell (float) – A conservative lower bound on the minimum possible primitive unit cell dimension.
- find_basis_vectors(reciprocal_lattice_vectors)[source]¶
Find a list of likely basis vectors.
- Parameters
reciprocal_lattice_vectors (scitbx.array_family.flex.vec3_double) – The list of reciprocal lattice vectors to search for periodicity.
- Returns
A tuple containing the list of basis vectors and a flex.bool array identifying which reflections were used in indexing.
- phil_help = 'Search for the basis vectors of the direct lattice by performing a 3D FFT in reciprocal space of the density of found spots. Since this can be quite memory-intensive, the data used for indexing may automatically be constrained to just the lower resolution spots.'¶
- phil_scope = <libtbx.phil.scope object>¶
- class dials.algorithms.indexing.basis_vector_search.RealSpaceGridSearch(max_cell, target_unit_cell, params=None, *args, **kwargs)[source]¶
Bases:
dials.algorithms.indexing.basis_vector_search.strategy.Strategy
Basis vector search using a real space grid search.
Search strategy to index found spots based on known unit cell parameters. It is often useful for difficult cases of narrow-wedge rotation data or stills data, especially where there is diffraction from multiple crystals.
A set of dimensionless radial unit vectors, typically ~7000 in total, is chosen so that they are roughly evenly spaced in solid angle over a hemisphere. For each direction, each of the three known unit cell vectors is aligned with the unit vector and is scored according to how well it accords with the periodicity in that direction of the reconstructed reciprocal space positions of the observed spot centroids. Examining the highest-scoring combinations, any basis vectors in orientations that are nearly collinear with a shorter basis vector are eliminated. The highest-scoring remaining combinations are selected as the basis of the direct lattice.
- See:
Gildea, R. J., Waterman, D. G., Parkhurst, J. M., Axford, D., Sutton, G., Stuart, D. I., Sauter, N. K., Evans, G. & Winter, G. (2014). Acta Cryst. D70, 2652-2666.
- __init__(max_cell, target_unit_cell, params=None, *args, **kwargs)[source]¶
Construct a real_space_grid_search object.
- Parameters
max_cell (float) – An estimate of the maximum cell dimension of the primitive cell.
target_unit_cell (cctbx.uctbx.unit_cell) – The target unit cell.
- static compute_functional(vector, reciprocal_lattice_vectors)[source]¶
Compute the functional for a single direction vector.
- Parameters
vector (tuple) – The vector at which to compute the functional.
reciprocal_lattice_vectors (scitbx.array_family.flex.vec3_double) – The list of reciprocal lattice vectors.
- Returns
The functional for the given vector.
- find_basis_vectors(reciprocal_lattice_vectors)[source]¶
Find a list of likely basis vectors.
- Parameters
reciprocal_lattice_vectors (scitbx.array_family.flex.vec3_double) – The list of reciprocal lattice vectors to search for periodicity.
- phil_help = 'Index the found spots by testing a known unit cell in various orientations until the best match is found. This strategy is often useful for difficult cases of narrow-wedge rotation data or stills data, especially where there is diffraction from multiple crystals.'¶
- phil_scope = <libtbx.phil.scope object>¶
- score_vectors(reciprocal_lattice_vectors)[source]¶
Compute the functional for the given directions.
- Parameters
directions – An iterable of the search directions.
reciprocal_lattice_vectors (scitbx.array_family.flex.vec3_double) – The list of reciprocal lattice vectors.
- Returns
A tuple containing the list of search vectors and their scores.
- property search_directions¶
Generator of the search directions (i.e. vectors with length 1).
- property search_vectors¶
Generator of the search vectors.
The lengths of the vectors correspond to the target unit cell dimensions.
- class dials.algorithms.indexing.basis_vector_search.Strategy(max_cell, params=None, *args, **kwargs)[source]¶
Bases:
object
A base class for basis vector search strategies.
- __init__(max_cell, params=None, *args, **kwargs)[source]¶
Construct the strategy.
- Parameters
max_cell (float) – An estimate of the maximum cell dimension of the primitive cell.
- find_basis_vectors(reciprocal_lattice_vectors)[source]¶
Find a list of likely basis vectors.
- Parameters
reciprocal_lattice_vectors (scitbx.array_family.flex.vec3_double) – The list of reciprocal lattice vectors to search for periodicity.
- Returns
A tuple containing the list of basis vectors and a flex.bool array identifying which reflections were used in indexing.
- phil_help = None¶
- phil_scope = None¶
- dials.algorithms.indexing.basis_vector_search.combinations.candidate_orientation_matrices(basis_vectors, max_combinations=None)[source]¶
- dials.algorithms.indexing.basis_vector_search.combinations.filter_known_symmetry(crystal_models, target_symmetry, relative_length_tolerance=0.1, absolute_angle_tolerance=5, max_delta=5)[source]¶
Filter crystal models for known symmetry.
- Parameters
crystal_models (list) – A list of
dxtbx.model.Crystal
objects.target_symmetry (cctbx.crystal.symmetry) – The target symmetry for filtering.
relative_length_tolerance (float) – Relative tolerance for unit cell lengths in unit cell comparison (default value is 0.1).
absolute_angle_tolerance (float) – Angular tolerance (in degrees) in unit cell comparison (default value is 5).
max_delta (float) – Maximum allowed Le Page delta used in searching for basis vector combinations that are consistent with the given symmetry (default value is 5).
- dials.algorithms.indexing.basis_vector_search.combinations.filter_similar_orientations(crystal_models, other_crystal_models, minimum_angular_separation=5)[source]¶
- class dials.algorithms.indexing.basis_vector_search.optimise.BasisVectorMinimiser(reciprocal_lattice_points, vector, lbfgs_termination_params=None, lbfgs_core_params=<scitbx.lbfgs.core_parameters object>)[source]¶
Bases:
object
dials.algorithms.indexing.lattice_search¶
- class dials.algorithms.indexing.lattice_search.LowResSpotMatch(target_symmetry_primitive, max_lattices, params=None, *args, **kwargs)[source]¶
Bases:
dials.algorithms.indexing.lattice_search.strategy.Strategy
A lattice search strategy matching low resolution spots to candidate indices.
The match is based on resolution and reciprocal space distance between observed spots. A prior unit cell and space group are required and solutions are assessed by matching the low resolution spots against candidate reflection positions predicted from the known cell. This lattice search strategy is a special case designed to work for electron diffraction still images, in which one typically only collects reflections from the zero-order Laue zone. In principle, it is not limited to this type of data, but probably works best with narrow wedges, good initial geometry and a small beam-stop shadow so that a good number of low-order reflections are collected.
- __init__(target_symmetry_primitive, max_lattices, params=None, *args, **kwargs)[source]¶
Construct a LowResSpotMatch object.
- Parameters
target_symmetry_primitive (cctbx.crystal.symmetry) – The target crystal symmetry and unit cell
max_lattices (int) – The maximum number of lattice models to find
- find_crystal_models(reflections, experiments)[source]¶
Find a list of candidate crystal models.
- Parameters
reflections (dials.array_family.flex.reflection_table) – The found spots centroids and associated data
experiments (dxtbx.model.experiment_list.ExperimentList) – The experimental geometry models
- phil_help = 'A lattice search strategy that matches low resolution spots to candidate indices based on a known unit cell and space group. Designed primarily for electron diffraction still images.'¶
- phil_scope = <libtbx.phil.scope object>¶
- class dials.algorithms.indexing.lattice_search.Strategy(params=None, *args, **kwargs)[source]¶
Bases:
object
A base class for lattice search strategies.
- __init__(params=None, *args, **kwargs)[source]¶
Construct the strategy.
- Parameters
params – an extracted PHIL scope containing the parameters
- find_crystal_models(reflections, experiments)[source]¶
Find a list of likely crystal models.
- Parameters
reflections (dials.array_family.flex.reflection_table) – The found spots centroids and associated data
experiments (dxtbx.model.experiment_list.ExperimentList) – The experimental geometry models
- Returns
A list of candidate crystal models.
- phil_help = None¶
- phil_scope = None¶
dials.algorithms.indexing.model_evaluation¶
- class dials.algorithms.indexing.model_evaluation.ModelRankFilter(check_doubled_cell=True, likelihood_cutoff=0.8, volume_cutoff=1.25, n_indexed_cutoff=0.9)[source]¶
- class dials.algorithms.indexing.model_evaluation.ModelRankWeighted(power=2, volume_weight=1, n_indexed_weight=1, rmsd_weight=1)[source]¶
- class dials.algorithms.indexing.model_evaluation.Result(model_likelihood, crystal, rmsds, n_indexed, fraction_indexed, hkl_offset)¶
Bases:
tuple
- crystal¶
Alias for field number 1
- fraction_indexed¶
Alias for field number 4
- hkl_offset¶
Alias for field number 5
- model_likelihood¶
Alias for field number 0
- n_indexed¶
Alias for field number 3
- rmsds¶
Alias for field number 2
dials.algorithms.indexing.max_cell¶
- dials.algorithms.indexing.max_cell.find_max_cell(reflections, max_cell_multiplier=1.3, step_size=45, nearest_neighbor_percentile=None, histogram_binning='linear', nn_per_bin=5, max_height_fraction=0.25, filter_ice=True, filter_overlaps=True, overlaps_border=0)[source]¶
- class dials.algorithms.indexing.nearest_neighbor.NeighborAnalysis(reflections, step_size=45, tolerance=1.5, max_height_fraction=0.25, percentile=None, histogram_binning='linear', nn_per_bin=5)[source]¶
Bases:
object
- class dials.algorithms.indexing.assign_indices.AssignIndicesGlobal(tolerance=0.3)[source]¶
Bases:
dials.algorithms.indexing.assign_indices.AssignIndicesStrategy
- class dials.algorithms.indexing.assign_indices.AssignIndicesLocal(d_min=None, epsilon=0.05, delta=8, l_min=0.8, nearest_neighbours=20)[source]¶
Bases:
dials.algorithms.indexing.assign_indices.AssignIndicesStrategy
- class dials.algorithms.indexing.assign_indices.AssignIndicesStrategy(d_min=None)[source]¶
Bases:
object
- dials.algorithms.indexing.compare_orientation_matrices.difference_rotation_matrix_axis_angle(crystal_a, crystal_b, target_angle=0)[source]¶
- dials.algorithms.indexing.compare_orientation_matrices.rotation_matrix_differences(crystal_models, miller_indices=None, comparison='pairwise')[source]¶
- class dials.algorithms.indexing.symmetry.SymmetryHandler(unit_cell=None, space_group=None, max_delta=5)[source]¶
Bases:
object
- apply_symmetry(crystal_model)[source]¶
Apply symmetry constraints to a crystal model.
Returns the crystal model (with symmetry constraints applied) in the same setting as provided as input. The cb_op returned by the method is that necessary to transform that model to the user-provided target symmetry.
- Parameters
crystal_model (dxtbx.model.Crystal) – The input crystal model to which to apply symmetry constraints.
Returns: (dxtbx.model.Crystal, cctbx.sgtbx.change_of_basis_op): The crystal model with symmetry constraints applied, and the change_of_basis_op that transforms the returned model to the user-specified target symmetry.