cctbx.miller
index
/net/chevy/raid1/nat/src/cctbx_project/cctbx/miller/__init__.py

# -*- coding: utf-8 -*-

 
Package Contents
       
display
reindexing
tst_map_to_asu_isym
tst_reindexing

 
Classes
       
__builtin__.object
array_info
binned_data
merge_equivalents
normalised_amplitudes
cctbx.crystal.symmetry(__builtin__.object)
set
array
cctbx.maptbx.crystal_gridding(__builtin__.object)
fft_map
crystal_symmetry_is_compatible_with_symmetry_from_file
cctbx_miller_ext.binner(cctbx_miller_ext.binning)
binner

 
class array(set)
    Extension of the set class with addition of data and sigmas flex arrays.
 
 
Method resolution order:
array
set
cctbx.crystal.symmetry
__builtin__.object

Methods defined here:
__abs__(self)
Return a copy of the array with data replaced by absolute values, i.e.
complex arrays will be converted to amplitudes.
__add__(self, other)
__getitem__(self, slice_object)
__imul__(self, other)
__init__(self, miller_set, data=None, sigmas=None)
__iter__(self)
__itruediv__(self, other)
__mul__(self, other)
__truediv__(self, other)
This requires from __future__ import division
adopt_set(self, other, assert_is_similar_symmetry=True)
amplitude_normalisations(self, asu_contents, wilson_plot=None)
Overriden version of set.amplitude_normalisation which computes
the Wilson parameters from the array data if wilson_plot is None.
amplitude_quasi_normalisations(self, d_star_power=1)
A miller.array whose data N(h) are the normalisations to convert
between locally normalised E's and F's:
E(h) = F(h) / N(h)
 
self features the F's, which are then binned with the current binner
and N(h) is the average of F's in the bin h belongs to.
amplitudes(self)
For a complex array, return array of absolute values.
analyze_intensity_statistics(self, d_min=2.5, log=None)
Detect translational pseudosymmetry and twinning, using methods in
Xtriage.  Returns a mmtbx.scaling.twin_analyses.twin_law_interpretation
object.  (Requires mmtbx to be configured to be functional.)
anomalous_completeness(self, use_binning=False, d_min_tolerance=1e-06, d_max=None, d_min=None, relative_to_complete_set=True)
Return the percent of acenric reflections with both h,k,l and -h,-k,-l
observed (only meaningful for amplitude and intensity arrays).  By default
this is calculated relative to the complete set.
anomalous_differences(self)
Returns an array object with DANO (i.e. F(+) - F(-)) as data, and
optionally SIGDANO as sigmas.
anomalous_signal(self, use_binning=False)
Get the anomalous signal according to this formula:
 
.. math::
   \sqrt{\dfrac{<||F(+)|-|F(-)||^2>}{\frac{1}{2} (<|F(+)|>^2 + <|F(-)|>^2)}}
 
:param use_binning: If 'True' the anomalous signal will be calculated for     each bin of the data set individually
:type use_binning: boolean
:returns: the anomalous signal
:rtype: float or cctbx.miller.binned_data
apply_debye_waller_factors(self, u_iso=None, b_iso=None, u_cart=None, b_cart=None, u_cif=None, u_star=None, apply_to_sigmas=True, exp_arg_limit=50, truncate_exp_arg=False)
apply_scaling(self, target_max=None, factor=None)
apply_shelxl_extinction_correction(self, x, wavelength)
arg(self, deg=False)
as_amplitude_array(self, algorithm='xtal_3_7')
Convert the array to simple amplitudes if not already in that format.
Only valid for complex (i.e. F,PHI), intensity, or amplitude arrays.
as_cif_block(self, array_type)
as_cif_simple(self, array_type, out=None, data_name='global')
as_double(self)
as_intensity_array(self, algorithm='simple')
Convert the array to intensities if not already in that format.  Only valid
for complex (F,PHI), amplitude, or intensity arrays.
as_mtz_dataset(self, column_root_label, column_types=None, label_decorator=None, title=None, crystal_name='crystal', project_name='project', dataset_name='dataset', wavelength=None)
as_non_anomalous_array(self)
as_phases_phs(self, out, scale_amplitudes=True, phases=None, phases_deg=None, figures_of_merit=None)
as_xray_observations(self, scale_indices=None, twin_fractions=None, twin_components=None)
average_bijvoet_mates(self)
Given an anomalous array, merge the anomalous pairs and return the
non-anomalous average.
bijvoet_ratios(self, obs_type='intensity', measurable_only=True, cutoff=3.0)
cc_anom(self, *args, **kwds)
cc_one_half(self, use_binning=False, n_trials=1, anomalous_flag=False)
Calculate the correlation between two randomly assigned pools of unmerged
data ("CC 1/2").  If desired the mean over multiple trials can be taken.
See Karplus PA & Diederichs K (2012) Science 336:1030-3 for motivation.
This method assumes that the reflections still have the original indices
and maps them to the ASU first; the underlying compute_cc_one_half
function skips this method.
change_basis(self, cb_op, deg=None)
complete_array(self, d_min_tolerance=1e-06, d_min=None, d_max=None, new_data_value=-1, new_sigmas_value=-1)
concatenate(self, other, assert_is_similar_symmetry=True)
conjugate(self)
copy(self)
correlation(self, other, use_binning=False, assert_is_similar_symmetry=True)
count_and_fraction_in_bins(self, data_value_to_count, count_not_equal=False)
crystal_symmetry_is_compatible_with_symmetry_from_file(self, unit_cell_relative_length_tolerance=0.02, unit_cell_absolute_angle_tolerance=3.0, working_point_group=None)
customized_copy(self, miller_set=<class libtbx.utils.Keep>, data=<class libtbx.utils.Keep>, sigmas=<class libtbx.utils.Keep>, crystal_symmetry=<class libtbx.utils.Keep>, indices=<class libtbx.utils.Keep>, anomalous_flag=<class libtbx.utils.Keep>, unit_cell=<class libtbx.utils.Keep>, space_group_info=<class libtbx.utils.Keep>, observation_type=<class libtbx.utils.Keep>)
data(self)
deep_copy(self)
detwin_data(self, twin_law, alpha)
Detwin data using a known twin fraction, returning an array with the same
original data type.
direct_summation_at_point(self, site_frac, sigma=None)
disagreeable_reflections(self, f_calc_sq, n_reflections=20)
discard_sigmas(self)
double_step_filtration(self, complete_set=None, vol_cutoff_plus_percent=5.0, vol_cutoff_minus_percent=5.0, resolution_factor=0.25, scale_to=None)
eliminate_sys_absent(self, integral_only=False, log=None, prefix='')
ellipsoidal_resolutions_and_indices_by_sigma(self, sigma_cutoff=3)
ellipsoidal_truncation_by_sigma(self, sigma_cutoff=3)
enforce_positive_amplitudes(self, i_sig_level=-4.0)
Takes in an intensity array (including negatives) and spits out amplitudes.
The basic assumption is that
P(Itrue) \propto exp(-(Itrue-Iobs)**2/(2*s))
where Itrue>=0 (positivity constraint on error free amplitudes).
For amplitudes, this results in
P(Ftrue) \propto 2 Ftrue exp( -(Ftrue**2-Iobs)**2/(2s) )
A Gaussian approximation is fitted to the Mode of this distribution.
An analytical solution exists and is implemented below.
This method does not require any Wilson statistics assumptions.
expand_to_p1(self, phase_deg=None, return_iselection=False)
export_as_cns_hkl(self, file_object, file_name=None, info=[], array_names=None, r_free_flags=None)
export_as_scalepack_unmerged(self, file_object=None, file_name=None, batch_numbers=None, spindle_flags=None, scale_intensities_for_scalepack_merge=False)
export_as_shelx_hklf(self, file_object=None, normalise_if_format_overflow=False)
f_as_f_sq(self, algorithm='simple')
Convert amplitudes (and associated sigmas, if present) to intensities.
f_obs_f_calc_fan_outlier_selection(self, f_calc, offset_low=0.05, offset_high=0.1, also_return_x_and_y=False)
Preconditions (not checked explicitly):
self is amplitude array,
f_calc is complex array or amplitude array.
f_obs_minus_f_calc(self, f_obs_factor, f_calc)
f_sq_as_f(self, algorithm='xtal_3_7', tolerance=1e-06)
fft_map(self, resolution_factor=0.3333333333333333, d_min=None, grid_step=None, crystal_gridding=None, symmetry_flags=None, mandatory_factors=None, max_prime=5, assert_shannon_sampling=True, f_000=None)
french_wilson(self, **kwds)
generate_bijvoet_mates(self)
Given a non-anomalous array, expand to generate anomalous pairs.
half_dataset_anomalous_correlation(self, use_binning=False)
Calculate the correlation of the anomalous differences of two randomly
assigned half-datasets (starting from unmerged data).
has_twinning(self, d_min=2.5)
Convenience method for identifying twinned data.  Note that this is
hugely inefficient if any other Xtriage analyses are planned, since it
discards the other results.  Requires mmtbx.
hemisphere_acentrics(self, plus_or_minus)
hemispheres_acentrics(self)
hoppe_gassmann_modification(self, mean_scale, n_iterations, resolution_factor=0.25, d_min=None)
i_over_sig_i(self, use_binning=False, return_fail=None)
<I/sigma_I>
info(self)
Return the associated info object, or None if undefined.
intensities(self)
is_bool_array(self)
is_complex_array(self)
is_hendrickson_lattman_array(self)
is_integer_array(self)
is_real_array(self)
is_string_array(self)
is_unmerged_intensity_array(self)
Determine whether the array contains unmerged experimental observations
or not.  In some files only the centric reflections will appear to be
unmerged, so we specifically check the acentrics (if present).
is_xray_amplitude_array(self)
is_xray_data_array(self)
is_xray_intensity_array(self)
is_xray_reconstructed_amplitude_array(self)
local_standard_deviation_map(self, radius, mean_solvent_density=0, resolution_factor=0.3333333333333333, d_min=None, grid_step=None, symmetry_flags=None, mandatory_factors=None, max_prime=5, assert_shannon_sampling=True, f_000=None)
map_correlation(self, other)
map_to_asu(self, deg=None)
matching_set(self, other, data_substitute, sigmas_substitute=None, assert_is_similar_symmetry=True)
mean(self, use_binning=False, use_multiplicities=False, squared=False, rms=False)
mean_of_intensity_divided_by_epsilon(self, use_binning=False, return_fail=None)
<I/epsilon>
mean_of_squared_sigma_divided_by_epsilon(self, use_binning=False, return_fail=None)
<sigma^2/epsilon>
mean_phase_error(self, phase_source)
mean_sq(self, use_binning=False, use_multiplicities=False)
mean_weighted_phase_error(self, phase_source)
measurability(self, use_binning=False, cutoff=3.0, return_fail=None)
Fraction of reflections for which
(:math:`\dfrac{|\Delta I|}{\sigma_{dI}}` > cutoff and
:math:`min(\dfrac{I_{+}}{\sigma_{+}},\dfrac{I_{-}}{\sigma_{-}})` > cutoff
merge_equivalents(self, algorithm='gaussian', incompatible_flags_replacement=None, use_internal_variance=True)
Given a non-unique array, merge the symmetry-related reflections (keeping
anomalous flag).
 
:returns: a merge_equivalents object, from which the merged array may     be extracted by calling the array() method.
min_f_over_sigma(self, return_none_if_zero_sigmas=False)
norm(self)
normalised_amplitudes(self, asu_contents, wilson_plot=None)
observation_type(self)
patterson_map(self, resolution_factor=0.3333333333333333, d_min=None, symmetry_flags=None, mandatory_factors=None, max_prime=5, assert_shannon_sampling=True, f_000=None, sharpening=False, origin_peak_removal=False)
patterson_symmetry(self)
permute_d_range(self, d_max, d_min)
Randomly shuffle reflections within a given resolution range.  Used for
control refinements to validate the information content of a dataset.
phase_entropy(self, exponentiate=False, return_binned_data=False, return_mean=False)
Get phase entropy as measured in terms of an base-360 entropy
(base-2 for centrics).
 
An entropy of 0, indicates that the phase uncertainity is as low as possible
An entropy of 1 however, indicates that the uncertainty is maximal:
all phases are equally likely!
 
:param return_binned_data: if 'True' you receive a binned object rather     then a raw array
:type return_binned_data: boolean
:param exponentiate: whether or not to exponentiate the entropy. This will     return a phase uncertainty in degrees (or the 'alphabet size')
:type exponentiate: boolean
phase_integrals(self, n_steps=None, integrator=None)
phase_transfer(self, phase_source, epsilon=1e-10, deg=False, phase_integrator_n_steps=None)
Combines phases of phase_source with self's data if real (keeping
the sign of self's data) or with self's amplitudes if complex.
 
Centric reflections are forced to be compatible with the phase restrictions.
 
phase_source can be a miller.array or a plain flex array.
 
epsilon is only used when phase_source is a complex array. If both the
real and the imaginary part of phase_source[i] < epsilon the phase is
assumed to be 0.
 
deg is only used if phase_source is an array of doubles.
deg=True indicates that the phases are given in degrees,
deg=False indicates phases are given in radians.
 
phase_integrator_n_steps is only used if phase_source is an
array of Hendrickson-Lattman coefficients. The centroid
phases are determined on the fly using the given step size.
phased_translation_function_coeff(self, phase_source, f_calc, fom=None)
phases(self, deg=False)
For a complex array, return the array of its phases (in radians by default).
quasi_normalize_structure_factors(self, d_star_power=1)
quasi_normalized_as_normalized(self)
r1_factor(self, other, scale_factor=None, assume_index_matching=False, use_binning=False)
Get the R1 factor according to this formula
 
.. math::
   R1 = \dfrac{\sum{||F| - k|F'||}}{\sum{|F|}}
 
where F is data() and F' is other.data() and
k is the factor to put F' on the same scale as F
r_free_flags_accumulation(self)
randomize_phases(self)
remove_cone(self, fraction_percent, vertex=(0, 0, 0), axis_point_1=(0, 0, 0), axis_point_2=(0, 0, 1), negate=False)
remove_patterson_origin_peak(self)
rms(self, use_binning=False, use_multiplicities=False)
rms_filter(self, cutoff_factor, use_binning=False, use_multiplicities=False, negate=False)
scale_factor(self, f_calc, weights=None, cutoff_factor=None, use_binning=False)
The analytical expression for the least squares scale factor.
 
K = sum(w * yo * yc) / sum(w * yc^2)
 
If the optional cutoff_factor argument is provided, only the reflections
whose magnitudes are greater than cutoff_factor * max(yo) will be included
in the calculation.
second_moment(self, use_binning=False)
<data^2>/(<data>)^2
second_moment_of_intensities(self, use_binning=False)
<I^2>/(<I>)^2 (2.0 for untwinned, 1.5 for twinned data)
select(self, selection, negate=False, anomalous_flag=None)
select_indices(self, indices, map_indices_to_asu=False, negate=False)
select_sys_absent(self, integral_only=False)
set(self, crystal_symmetry=<class libtbx.utils.Keep>, indices=<class libtbx.utils.Keep>, anomalous_flag=<class libtbx.utils.Keep>, unit_cell=<class libtbx.utils.Keep>, space_group_info=<class libtbx.utils.Keep>)
set_info(self, info)
set_observation_type(self, observation_type)
set_observation_type_xray_amplitude(self)
set_observation_type_xray_intensity(self)
set_sigmas(self, sigmas)
shelxl_extinction_correction(self, x, wavelength)
Extinction parameter x, where Fc is multiplied by:
  k[1 + 0.001 x Fc^2 wavelength^3 / sin(2theta)]^(-1/4)
 
See SHELX-97 manual, page 7-7 for more information.
 
Note: The scale factor, k, is not applied nor calculated by
      this function. The scale factor should be calculated
      and applied ***AFTER*** the application of the extinction
      corrections.
show_array(self, f=None, prefix='', deg=None)
Listing of Miller indices and data
show_comprehensive_summary(self, f=None, prefix='')
show_disagreeable_reflections(self, f_calc_sq, n_reflections=20, out=None)
show_mean_data_over_sigma_along_a_b_c_star(self)
show_r_free_flags_info(self, n_bins=10, binner_range='used', out=None, prefix='')
show_summary(self, f=None, prefix='')
sigma_filter(self, cutoff_factor, negate=False)
Return a copy of the array filtered to remove reflections whose value is
less than cutoff_factor*sigma (or the reverse, if negate=True).
sigmas(self)
sigmas_are_sensible(self, critical_ratio=0.75, epsilon=1e-06)
size(self)
sort_permutation(self, by_value='resolution', reverse=False)
statistical_mean(self, use_binning=0)
sum(self, use_binning=False, use_multiplicities=False, squared=False)
sum_sq(self, use_binning=False, use_multiplicities=False)
symmetry_agreement_factor(self, op, assert_is_similar_symmetry=True)
The factor phi_{sym} quantifying whether complex structure factors
are invariant under the given symmetry operator, as used in Superflip.
Ref: J. Appl. Cryst. (2008). 41, 975-984
twin_data(self, twin_law, alpha)
Apply a twin law to the data, returning an array of the same original type.
value_at_index(self, hkl)
Extract the value of the array for the specified reflection h,k,l
wilson_plot(self, use_binning=False)
<F^2>
wilson_ratio(self, use_binning=False)
(<F>)^2/<F^2> (0.785 for untwinned, 0.885 for twinned data)

Class methods defined here:
from_cif(cls, file_object=None, file_path=None, data_block_name=None) from __builtin__.type

Methods inherited from set:
all_selection(self)
anomalous_flag(self)
array(self, data=None, sigmas=None)
Create an array object, given data and/or sigma arrays of identical
dimensions to the indices array.
 
:param data: a flex array (any format) or None
:param sigmas: a flex array (any format, but almost always double) or None
as_non_anomalous_set(self)
at_first_index(self, ary, miller_index)
Returns the element `ary` coresponding to the `miller_index` if
`miller_index exists, otherwise returns None.
 
:param miller_index: Miller index as a 3-tuple
:type miller_index: tuple
:param ary: any array (e.g. data(), sigmas())
:type ary: sequence (list, array, etc)
:returns: type of contents of `ary`, or None
auto_anomalous(self, min_n_bijvoet_pairs=None, min_fraction_bijvoet_pairs=None)
binner(self)
centric_flags(self)
Generate a boolean Miller array flagging centric reflections.
clear_binner(self)
combine(self, other, scale=True, scale_for_lones=1)
common_set(self, other, assert_is_similar_symmetry=True)
Match the indices in the current set and another set, and return a set
(or array) containing only those reflections present in both.  Assumes that
both sets are already in the asymmetric unit (ASU).
common_sets(self, other, assert_is_similar_symmetry=True, assert_no_singles=False)
Like common_set(other), but returns a tuple containing matching copies of
both sets (or arrays).
complete_set(self, d_min_tolerance=1e-06, d_min=None, d_max=None, max_index=None)
Generate the complete set of Miller indices expected for the current
symmetry, excepting systematic absences.
 
:param d_min_tolerance: tolerance factor for d_min (avoid precision errors)
:param d_min: High-resolution limit (default = d_min of current set)
:param d_max: Low-resolution limit (default = d_max of current set)
complete_with(self, other, scale=False, replace_phases=False)
complete_with_bin_average(self, reflections_per_bin=100)
completeness(self, use_binning=False, d_min_tolerance=1e-06, return_fail=None, d_max=None)
crystal_gridding(self, resolution_factor=0.3333333333333333, d_min=None, grid_step=None, symmetry_flags=None, mandatory_factors=None, max_prime=5, assert_shannon_sampling=True)
crystal_symmetry(self)
Get crystal symmetry of the miller set
 
:returns: a new crystal.symmetry object
:rtype: cctbx.crystal.symmetry
d_max_min(self)
Low- and high-resolution limits.
d_min(self)
High-resolution limit.
d_min_along_a_b_c_star(self)
Returns the effective resolution limits along the reciprocal space axes.
d_spacings(self)
Generate a double Miller array containing the resolution d of each
index.
d_star_cubed(self)
d_star_sq(self)
data_at_first_index(self, miller_index)
Returns the value of data of the first index matching
`miller_index`. If the `miller_index` is not found in `self`,
then returns ``None``.
 
:param miller_index: Miller index as a 3-tuple
:type miller_index: tuple
:returns: int, float, complex, None -- data value or None
debye_waller_factors(self, u_iso=None, b_iso=None, u_cart=None, b_cart=None, u_cif=None, u_star=None, exp_arg_limit=50, truncate_exp_arg=False)
delete_index(self, hkl)
Remove all reflections with the specified Miller index.
delete_indices(self, other)
Delete multiple reflections, as specified by the Miller indices of
another set.
epsilons(self)
expand_to_p1_iselection(self, build_iselection=True)
f_obs_minus_xray_structure_f_calc(self, f_obs_factor, xray_structure, structure_factor_algorithm=None, cos_sin_table=False, quality_factor=None, u_base=None, b_base=None, wing_cutoff=None, exp_table_one_over_step_size=None)
first_index(self, miller_index)
Returns the first index of the item  matching
`miller_index`. If the `miller_index` is not found in `self`,
then returns ``None``.
 
:param miller_index: Miller index as a 3-tuple
:type miller_index: tuple
:returns: int, None -- index of first occurrence of
          `miller_index` or None
generate_r_free_flags(self, fraction=0.1, max_free=2000, lattice_symmetry_max_delta=5.0, use_lattice_symmetry=False, use_dataman_shells=False, n_shells=20, format='cns')
Create an array of R-free flags for the current set, keeping anomalous
pairs together.  Requires that the set already be unique under symmetry,
and generally assumes that the set is in the ASU.
 
:param fraction: fraction of reflections to flag for the test set
:param max_free: limit on size of test set, overrides fraction
:param lattice_symmetry_max_delta: limit on lattice symmetry calculation
:param use_lattice_symmetry: given the current symmetry, determine the     highest possible lattice symmetry and generate flags for that symmetry,     then expand to the current (lower) symmetry if necessary.  This is almost     always a good idea.
:param use_dataman_shells: generate flags in thin resolution shells to     avoid bias due to non-crystallographic symmetry.
:param n_shells: number of resolution shells if use_dataman_shells=True
:param format: convention of the resulting flags.  'cns' will return a     boolean array (True = free), 'ccp4' will return an integer array from     0 to X (0 = free, X dependent on fraction), 'shelx' will return an     integer array with values 1 (work) or -1 (free).
 
:returns: a boolean or integer Miller array, depending on format.
generate_r_free_flags_basic(self, fraction=0.1, max_free=2000, use_dataman_shells=False, n_shells=20, format='cns')
generate_r_free_flags_on_lattice_symmetry(self, fraction=0.1, max_free=2000, max_delta=5.0, return_integer_array=False, n_partitions=None, use_dataman_shells=False, n_shells=20, format='cns')
index_span(self)
indices(self)
is_in_asu(self)
is_unique_set_under_symmetry(self)
log_binning(self, n_reflections_in_lowest_resolution_bin=100, eps=0.0001, max_number_of_bins=30, min_reflections_in_bin=50)
lone_set(self, other, assert_is_similar_symmetry=True)
Match the indices in the current set and another set, and return a set
(or array) containing reflections which are present only in the current
set.  Assumes that both sets are already in the asymmetric unit.
lone_sets(self, other, assert_is_similar_symmetry=True)
Like lone_set(other), but returns a tuple containing the reflections
unique to each set (or array).
match_bijvoet_mates(self, assert_is_unique_set_under_symmetry=True)
match_indices(self, other, assert_is_similar_symmetry=True)
miller_indices_as_pdb_file(self, file_name=None, expand_to_p1=False)
Write out Miller indices as pseudo-waters for visualization.  Note that
this treats the indices as literal coordinates (times a scale factor),
and the distances between points will not correspond to the distances
in reciprocal space.
 
See cctbx/miller/display.py and crys3d/hklview for an alternative (but
less lightweight) approach.
min_max_d_star_sq(self)
min_max_indices(self)
minimum_wavelength_based_on_d_min(self, tolerance=0.01)
multiplicities(self)
multiscale(self, other, reflections_per_bin=None)
n_bijvoet_pairs(self)
random_phases_compatible_with_phase_restrictions(self, deg=False)
randomize_amplitude_and_phase(self, amplitude_error, phase_error_deg, selection=None, random_seed=None)
reflection_intensity_symmetry(self)
remove_systematic_absences(self, negate=False)
resolution_filter(self, d_max=0, d_min=0, negate=0)
resolution_filter_selection(self, d_max=None, d_min=None)
resolution_range(self)
scale(self, other, resolution_dependent=False)
select_acentric(self)
select_centric(self)
setup_binner(self, d_max=0, d_min=0, auto_binning=False, reflections_per_bin=0, n_bins=0)
setup_binner_counting_sorted(self, d_max=0, d_min=0, reflections_per_bin=None, d_tolerance=1e-10)
setup_binner_d_star_sq_step(self, auto_binning=True, d_max=None, d_min=None, d_star_sq_step=None)
show_completeness(self, reflections_per_bin=500, out=None)
sigma_at_first_index(self, miller_index)
Returns the value of sigmas of the first index matching
`miller_index`. If the `miller_index` is not found in `self`,
then returns ``None``.
 
:param miller_index: Miller index as a 3-tuple
:type miller_index: tuple
:returns: int, float, complex, None -- sigmas value or None
sin_theta_over_lambda_sq(self)
slice(self, axis=None, axis_index=None, slice_index=None, slice_start=None, slice_end=None)
sort(self, by_value='resolution', reverse=False)
Reorder reflections by resolution or Miller index.
 
:param by_value: 'resolution' or 'packed_indices'
structure_factors_from_map(self, map, in_place_fft=False, use_scale=False, anomalous_flag=None, use_sg=False)
structure_factors_from_scatterers(self, xray_structure, algorithm=None, cos_sin_table=False, grid_resolution_factor=0.3333333333333333, quality_factor=None, u_base=None, b_base=None, wing_cutoff=None, exp_table_one_over_step_size=None)
sys_absent_flags(self, integral_only=False)
Generate a boolean Miller array flagging those reflections which are
systematically absent under the current symmetry.
two_theta(self, wavelength, deg=False)
Generate a double Miller array containing the scattering angle of each
index.
unique_under_symmetry(self)
unique_under_symmetry_selection(self)
use_binner_of(self, other)
use_binning(self, binning)
use_binning_of(self, other)

Methods inherited from cctbx.crystal.symmetry:
as_py_code(self, indent='')
as_reference_setting(self)
asu_mappings(self, buffer_thickness, asu_is_inside_epsilon=None)
average_b_cart(self, b_cart)
average_u_cart(self, u_cart)
best_cell(self, angular_tolerance=None)
build_miller_set(self, anomalous_flag, d_min, d_max=None)
cell_equivalent_p1(self)
change_of_basis_op_to_best_cell(self, angular_tolerance=None, best_monoclinic_beta=True)
change_of_basis_op_to_inverse_hand(self)
change_of_basis_op_to_minimum_cell(self)
change_of_basis_op_to_niggli_cell(self, relative_epsilon=None, iteration_limit=None)
change_of_basis_op_to_primitive_setting(self)
change_of_basis_op_to_reference_setting(self)
direct_space_asu(self)
gridding(self, d_min=None, resolution_factor=None, step=None, symmetry_flags=None, mandatory_factors=None, max_prime=5, assert_shannon_sampling=True)
inverse_hand(self)
is_compatible_unit_cell(self)
is_patterson_symmetry(self)
is_similar_symmetry(self, other, relative_length_tolerance=0.01, absolute_angle_tolerance=1.0)
join_symmetry(self, other_symmetry, force=False, raise_sorry_if_incompatible_unit_cell=False)
miller_set(self, indices, anomalous_flag)
minimum_cell(self)
niggli_cell(self, relative_epsilon=None, iteration_limit=None)
primitive_setting(self)
space_group(self)
space_group_info(self)
special_position_settings(self, min_distance_sym_equiv=0.5, u_star_tolerance=0, assert_min_distance_sym_equiv=True)
subtract_continuous_allowed_origin_shifts(self, translation_cart)
unit_cell(self)

Data descriptors inherited from cctbx.crystal.symmetry:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class array_info(__builtin__.object)
    Container for metadata associated with a Miller array, including labels
read from a data file.
 
  Methods defined here:
__init__(self, source=None, source_type=None, history=None, labels=None, merged=False, systematic_absences_eliminated=False, crystal_symmetry_from_file=None, type_hints_from_file=None, wavelength=None)
__setstate__(self, state)
__str__(self)
as_string_part_2(self)
customized_copy(self, source=<class libtbx.utils.Keep>, source_type=<class libtbx.utils.Keep>, history=<class libtbx.utils.Keep>, labels=<class libtbx.utils.Keep>, merged=<class libtbx.utils.Keep>, systematic_absences_eliminated=<class libtbx.utils.Keep>, crystal_symmetry_from_file=<class libtbx.utils.Keep>, type_hints_from_file=<class libtbx.utils.Keep>, wavelength=<class libtbx.utils.Keep>)
label_string(self)
A combined representation of the data labels extracted from the input file.
This is generally how downstream programs will identify and select Miller
arrays.

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class binned_data(__builtin__.object)
     Methods defined here:
__init__(self, binner, data, data_fmt=None)
as_simple_table(self, data_label, data_fmt=None)
Export table rows for display elsewhere.  (Used in Xtriage)
show(self, data_fmt=None, show_bin_number=True, show_bin_range=True, show_d_range=None, show_counts=True, show_unused=True, bin_range_as='d', wavelength=None, f=None, prefix='')

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class binner(cctbx_miller_ext.binner)
    
Method resolution order:
binner
cctbx_miller_ext.binner
cctbx_miller_ext.binning
Boost.Python.instance
__builtin__.object

Methods defined here:
__getinitargs__(self)
__init__(self, binning, miller_set)
as_simple_table(self, data, data_label, data_fmt=None, show_bin_number=False, show_unused=False)
Export table rows for display elsewhere.
bin_legend(self, i_bin, show_bin_number=True, show_bin_range=True, show_d_range=None, show_counts=True, bin_range_as='d', wavelength=None, join_string=True)
counts_complete(self, include_centric=True, include_acentric=True, d_min_tolerance=1e-06)
counts_given(self)
n_bin_d_too_large(self)
n_bin_d_too_large_or_small(self)
n_bin_d_too_small(self)
show_data(self, data, data_fmt=None, show_bin_number=True, show_bin_range=True, show_d_range=None, show_counts=True, show_unused=True, bin_range_as='d', wavelength=None, f=None, prefix='')
show_summary(self, show_bin_number=True, show_bin_range=True, show_d_range=None, show_counts=True, bin_range_as='d', wavelength=None, f=None, prefix='')

Methods inherited from cctbx_miller_ext.binner:
__reduce__ = (...)
array_indices(...)
array_indices( (binner)arg1, (int)arg2) -> size_t :
 
    C++ signature :
        scitbx::af::shared<unsigned long> array_indices(cctbx::miller::binner {lvalue},unsigned long)
bin_centers(...)
bin_centers( (binner)arg1, (float)arg2) -> double :
 
    C++ signature :
        scitbx::af::shared<double> bin_centers(cctbx::miller::binner {lvalue},double)
bin_indices(...)
bin_indices( (binner)arg1) -> size_t :
 
    C++ signature :
        scitbx::af::shared<unsigned long> bin_indices(cctbx::miller::binner {lvalue})
count(...)
count( (binner)arg1, (int)arg2) -> int :
 
    C++ signature :
        unsigned long count(cctbx::miller::binner {lvalue},unsigned long)
counts(...)
counts( (binner)arg1) -> size_t :
 
    C++ signature :
        scitbx::af::shared<unsigned long> counts(cctbx::miller::binner {lvalue})
interpolate(...)
interpolate( (binner)arg1, (object)arg2, (float)arg3) -> double :
 
    C++ signature :
        scitbx::af::shared<double> interpolate(cctbx::miller::binner {lvalue},scitbx::af::const_ref<double, scitbx::af::trivial_accessor>,double)
miller_indices(...)
miller_indices( (binner)arg1) -> miller_index :
 
    C++ signature :
        scitbx::af::shared<cctbx::miller::index<int> > miller_indices(cctbx::miller::binner {lvalue})
selection(...)
selection( (binner)arg1, (int)arg2) -> bool :
 
    C++ signature :
        scitbx::af::shared<bool> selection(cctbx::miller::binner {lvalue},unsigned long)

Data and other attributes inherited from cctbx_miller_ext.binner:
__safe_for_unpickling__ = True

Methods inherited from cctbx_miller_ext.binning:
bin_d_min(...)
bin_d_min( (binning)arg1, (int)arg2) -> float :
 
    C++ signature :
        double bin_d_min(cctbx::miller::binning {lvalue},unsigned long)
bin_d_range(...)
bin_d_range( (binning)arg1, (int)arg2) -> tuple :
 
    C++ signature :
        scitbx::af::tiny<double, 2ul> bin_d_range(cctbx::miller::binning {lvalue},unsigned long)
d_max(...)
d_max( (binning)arg1) -> float :
 
    C++ signature :
        double d_max(cctbx::miller::binning {lvalue})
d_min(...)
d_min( (binning)arg1) -> float :
 
    C++ signature :
        double d_min(cctbx::miller::binning {lvalue})
get_i_bin(...)
get_i_bin( (binning)arg1, (float)arg2) -> int :
 
    C++ signature :
        unsigned long get_i_bin(cctbx::miller::binning {lvalue},double)
 
get_i_bin( (binning)arg1, (object)arg2) -> int :
 
    C++ signature :
        unsigned long get_i_bin(cctbx::miller::binning {lvalue},cctbx::miller::index<int>)
i_bin_d_too_large(...)
i_bin_d_too_large( (binning)arg1) -> int :
 
    C++ signature :
        unsigned long i_bin_d_too_large(cctbx::miller::binning {lvalue})
i_bin_d_too_small(...)
i_bin_d_too_small( (binning)arg1) -> int :
 
    C++ signature :
        unsigned long i_bin_d_too_small(cctbx::miller::binning {lvalue})
limits(...)
limits( (binning)arg1) -> double :
 
    C++ signature :
        scitbx::af::shared<double> limits(cctbx::miller::binning {lvalue})
n_bins_all(...)
n_bins_all( (binning)arg1) -> int :
 
    C++ signature :
        unsigned long n_bins_all(cctbx::miller::binning {lvalue})
n_bins_used(...)
n_bins_used( (binning)arg1) -> int :
 
    C++ signature :
        unsigned long n_bins_used(cctbx::miller::binning {lvalue})
range_all(...)
range_all( (binning)arg1) -> object :
 
    C++ signature :
        boost::python::api::object range_all(cctbx::miller::binning)
range_used(...)
range_used( (binning)arg1) -> object :
 
    C++ signature :
        boost::python::api::object range_used(cctbx::miller::binning)
unit_cell(...)
unit_cell( (binning)arg1) -> unit_cell :
 
    C++ signature :
        cctbx::uctbx::unit_cell unit_cell(cctbx::miller::binning {lvalue})

Data descriptors inherited from Boost.Python.instance:
__dict__
__weakref__

Data and other attributes inherited from Boost.Python.instance:
__new__ = <built-in method __new__ of Boost.Python.class object>
T.__new__(S, ...) -> a new object with type S, a subtype of T

 
class crystal_symmetry_is_compatible_with_symmetry_from_file
     Methods defined here:
__init__(self, miller_array, unit_cell_relative_length_tolerance=0.02, unit_cell_absolute_angle_tolerance=3.0, working_point_group=None)
format_error_message(self, data_description)

 
class fft_map(cctbx.maptbx.crystal_gridding)
    
Method resolution order:
fft_map
cctbx.maptbx.crystal_gridding
__builtin__.object

Methods defined here:
__init__(self, crystal_gridding, fourier_coefficients, f_000=None)
anomalous_flag(self)
apply_fourier_scaling(self)
apply_scaling(self, scale)
apply_sigma_scaling(self)
apply_volume_scaling(self)
as_ccp4_map(self, file_name, gridding_first=None, gridding_last=None, labels=['cctbx.miller.fft_map'])
as_xplor_map = cctbx_miller_fft_map_as_xplor_map(self, file_name, title_lines=['cctbx.miller.fft_map'], gridding_first=None, gridding_last=None, average=None, standard_deviation=None)
complex_map(self)
peak_search(self, parameters=None, verify_symmetry=True)
real_map(self, direct_access=True)
real_map_unpadded(self, in_place=True)
statistics(self)

Methods inherited from cctbx.maptbx.crystal_gridding:
change_space_group(self, space_group_info)
crystal_symmetry(self)
d_min(self)
mandatory_factors(self)
max_prime(self)
n_grid_points(self)
n_real(self)
resolution_factor(self)
space_group(self)
space_group_info(self)
symmetry_flags(self)
tags(self)
unit_cell(self)

Data descriptors inherited from cctbx.maptbx.crystal_gridding:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class merge_equivalents(__builtin__.object)
     Methods defined here:
__init__(self, miller_array, algorithm='gaussian', incompatible_flags_replacement=None, use_internal_variance=True)
array(self)
inconsistent_equivalents(self)
r_int(self)
r_linear(self)
R-linear = sum(abs(data - mean(data))) / sum(abs(data))
r_meas(self)
r_merge(self)
r_pim(self)
r_sigma(self)
r_square(self)
R-square = sum((data - mean(data))**2) / sum(data**2)
redundancies(self)
show_summary(self, n_bins=10, out=None, prefix='')

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class normalised_amplitudes(__builtin__.object)
    E-values and related statistics
 
  Methods defined here:
__init__(self, miller_array, asu_contents, wilson_plot=None)
array(self)
mean_e_sq_minus_1(self)
percent_e_sq_gt_2(self)

Data descriptors defined here:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
class set(cctbx.crystal.symmetry)
    Basic class for handling sets of Miller indices (h,k,l), including sorting
and matching functions, symmetry handling, generation of R-free flags, and
extraction of associated statistics.  Does not actually contain data, but
this can be added using the array(...) method.
 
 
Method resolution order:
set
cctbx.crystal.symmetry
__builtin__.object

Methods defined here:
__getitem__(self, slice_object)
__init__(self, crystal_symmetry, indices, anomalous_flag=None)
all_selection(self)
amplitude_normalisations(self, asu_contents, wilson_plot)
A miller.array whose data N(h) are the normalisations to convert
between E's and F's:
E(h) = F(h) / N(h)
The argument wilson_plot shall feature attributes
- wilson_intensity_scale_factor
- wilson_b
anomalous_flag(self)
array(self, data=None, sigmas=None)
Create an array object, given data and/or sigma arrays of identical
dimensions to the indices array.
 
:param data: a flex array (any format) or None
:param sigmas: a flex array (any format, but almost always double) or None
as_non_anomalous_set(self)
at_first_index(self, ary, miller_index)
Returns the element `ary` coresponding to the `miller_index` if
`miller_index exists, otherwise returns None.
 
:param miller_index: Miller index as a 3-tuple
:type miller_index: tuple
:param ary: any array (e.g. data(), sigmas())
:type ary: sequence (list, array, etc)
:returns: type of contents of `ary`, or None
auto_anomalous(self, min_n_bijvoet_pairs=None, min_fraction_bijvoet_pairs=None)
binner(self)
centric_flags(self)
Generate a boolean Miller array flagging centric reflections.
change_basis(self, cb_op)
Get a transformation of the miller set with a new basis specified by cb_op
 
:param cb_op: object describing the desired transformation of the basis
:type cb_op: string or sgtbx.change_of_basis_operator
 
:returns: a new miller set with the new basis
:rtype: cctbx.miller.set
clear_binner(self)
combine(self, other, scale=True, scale_for_lones=1)
common_set(self, other, assert_is_similar_symmetry=True)
Match the indices in the current set and another set, and return a set
(or array) containing only those reflections present in both.  Assumes that
both sets are already in the asymmetric unit (ASU).
common_sets(self, other, assert_is_similar_symmetry=True, assert_no_singles=False)
Like common_set(other), but returns a tuple containing matching copies of
both sets (or arrays).
complete_set(self, d_min_tolerance=1e-06, d_min=None, d_max=None, max_index=None)
Generate the complete set of Miller indices expected for the current
symmetry, excepting systematic absences.
 
:param d_min_tolerance: tolerance factor for d_min (avoid precision errors)
:param d_min: High-resolution limit (default = d_min of current set)
:param d_max: Low-resolution limit (default = d_max of current set)
complete_with(self, other, scale=False, replace_phases=False)
complete_with_bin_average(self, reflections_per_bin=100)
completeness(self, use_binning=False, d_min_tolerance=1e-06, return_fail=None, d_max=None)
concatenate(self, other, assert_is_similar_symmetry=True)
Combine two Miller sets.  Both must have the same anomalous flag, and
similar symmetry is also assumed.
copy(self)
crystal_gridding(self, resolution_factor=0.3333333333333333, d_min=None, grid_step=None, symmetry_flags=None, mandatory_factors=None, max_prime=5, assert_shannon_sampling=True)
crystal_symmetry(self)
Get crystal symmetry of the miller set
 
:returns: a new crystal.symmetry object
:rtype: cctbx.crystal.symmetry
customized_copy(self, crystal_symmetry=<class libtbx.utils.Keep>, indices=<class libtbx.utils.Keep>, anomalous_flag=<class libtbx.utils.Keep>, unit_cell=<class libtbx.utils.Keep>, space_group_info=<class libtbx.utils.Keep>)
Create a copy of the set, optionally changing the symmetry, indices,
and/or anomalous flag (default = keep all unmodified).
d_max_min(self)
Low- and high-resolution limits.
d_min(self)
High-resolution limit.
d_min_along_a_b_c_star(self)
Returns the effective resolution limits along the reciprocal space axes.
d_spacings(self)
Generate a double Miller array containing the resolution d of each
index.
d_star_cubed(self)
d_star_sq(self)
data_at_first_index(self, miller_index)
Returns the value of data of the first index matching
`miller_index`. If the `miller_index` is not found in `self`,
then returns ``None``.
 
:param miller_index: Miller index as a 3-tuple
:type miller_index: tuple
:returns: int, float, complex, None -- data value or None
debye_waller_factors(self, u_iso=None, b_iso=None, u_cart=None, b_cart=None, u_cif=None, u_star=None, exp_arg_limit=50, truncate_exp_arg=False)
deep_copy(self)
delete_index(self, hkl)
Remove all reflections with the specified Miller index.
delete_indices(self, other)
Delete multiple reflections, as specified by the Miller indices of
another set.
epsilons(self)
expand_to_p1(self, return_iselection=False)
Get a transformation of the miller set to spacegroup P1
 
:returns: a new set of parameters (symmetry, miller indices, anomalous_flag) in spacegroup P1
:rtype: set(cctbx.crystal.symmetry, cctbx.miller.indices, boolean)
expand_to_p1_iselection(self, build_iselection=True)
f_obs_minus_xray_structure_f_calc(self, f_obs_factor, xray_structure, structure_factor_algorithm=None, cos_sin_table=False, quality_factor=None, u_base=None, b_base=None, wing_cutoff=None, exp_table_one_over_step_size=None)
first_index(self, miller_index)
Returns the first index of the item  matching
`miller_index`. If the `miller_index` is not found in `self`,
then returns ``None``.
 
:param miller_index: Miller index as a 3-tuple
:type miller_index: tuple
:returns: int, None -- index of first occurrence of
          `miller_index` or None
generate_r_free_flags(self, fraction=0.1, max_free=2000, lattice_symmetry_max_delta=5.0, use_lattice_symmetry=False, use_dataman_shells=False, n_shells=20, format='cns')
Create an array of R-free flags for the current set, keeping anomalous
pairs together.  Requires that the set already be unique under symmetry,
and generally assumes that the set is in the ASU.
 
:param fraction: fraction of reflections to flag for the test set
:param max_free: limit on size of test set, overrides fraction
:param lattice_symmetry_max_delta: limit on lattice symmetry calculation
:param use_lattice_symmetry: given the current symmetry, determine the     highest possible lattice symmetry and generate flags for that symmetry,     then expand to the current (lower) symmetry if necessary.  This is almost     always a good idea.
:param use_dataman_shells: generate flags in thin resolution shells to     avoid bias due to non-crystallographic symmetry.
:param n_shells: number of resolution shells if use_dataman_shells=True
:param format: convention of the resulting flags.  'cns' will return a     boolean array (True = free), 'ccp4' will return an integer array from     0 to X (0 = free, X dependent on fraction), 'shelx' will return an     integer array with values 1 (work) or -1 (free).
 
:returns: a boolean or integer Miller array, depending on format.
generate_r_free_flags_basic(self, fraction=0.1, max_free=2000, use_dataman_shells=False, n_shells=20, format='cns')
generate_r_free_flags_on_lattice_symmetry(self, fraction=0.1, max_free=2000, max_delta=5.0, return_integer_array=False, n_partitions=None, use_dataman_shells=False, n_shells=20, format='cns')
index_span(self)
indices(self)
is_in_asu(self)
is_unique_set_under_symmetry(self)
log_binning(self, n_reflections_in_lowest_resolution_bin=100, eps=0.0001, max_number_of_bins=30, min_reflections_in_bin=50)
lone_set(self, other, assert_is_similar_symmetry=True)
Match the indices in the current set and another set, and return a set
(or array) containing reflections which are present only in the current
set.  Assumes that both sets are already in the asymmetric unit.
lone_sets(self, other, assert_is_similar_symmetry=True)
Like lone_set(other), but returns a tuple containing the reflections
unique to each set (or array).
map_to_asu(self)
Convert all indices to lie within the canonical asymmetric unit for the
current space group (while preserving anomalous flag).  Required for many
downstream steps.
match_bijvoet_mates(self, assert_is_unique_set_under_symmetry=True)
match_indices(self, other, assert_is_similar_symmetry=True)
miller_indices_as_pdb_file(self, file_name=None, expand_to_p1=False)
Write out Miller indices as pseudo-waters for visualization.  Note that
this treats the indices as literal coordinates (times a scale factor),
and the distances between points will not correspond to the distances
in reciprocal space.
 
See cctbx/miller/display.py and crys3d/hklview for an alternative (but
less lightweight) approach.
min_max_d_star_sq(self)
min_max_indices(self)
minimum_wavelength_based_on_d_min(self, tolerance=0.01)
multiplicities(self)
multiscale(self, other, reflections_per_bin=None)
n_bijvoet_pairs(self)
patterson_symmetry(self)
random_phases_compatible_with_phase_restrictions(self, deg=False)
randomize_amplitude_and_phase(self, amplitude_error, phase_error_deg, selection=None, random_seed=None)
reflection_intensity_symmetry(self)
remove_systematic_absences(self, negate=False)
resolution_filter(self, d_max=0, d_min=0, negate=0)
resolution_filter_selection(self, d_max=None, d_min=None)
resolution_range(self)
scale(self, other, resolution_dependent=False)
select(self, selection, negate=False, anomalous_flag=None)
select_acentric(self)
select_centric(self)
setup_binner(self, d_max=0, d_min=0, auto_binning=False, reflections_per_bin=0, n_bins=0)
setup_binner_counting_sorted(self, d_max=0, d_min=0, reflections_per_bin=None, d_tolerance=1e-10)
setup_binner_d_star_sq_step(self, auto_binning=True, d_max=None, d_min=None, d_star_sq_step=None)
show_completeness(self, reflections_per_bin=500, out=None)
show_comprehensive_summary(self, f=None, prefix='')
Comprehensive Miller set or array summary
show_summary(self, f=None, prefix='')
Minimal Miller set summary
sigma_at_first_index(self, miller_index)
Returns the value of sigmas of the first index matching
`miller_index`. If the `miller_index` is not found in `self`,
then returns ``None``.
 
:param miller_index: Miller index as a 3-tuple
:type miller_index: tuple
:returns: int, float, complex, None -- sigmas value or None
sin_theta_over_lambda_sq(self)
size(self)
slice(self, axis=None, axis_index=None, slice_index=None, slice_start=None, slice_end=None)
sort(self, by_value='resolution', reverse=False)
Reorder reflections by resolution or Miller index.
 
:param by_value: 'resolution' or 'packed_indices'
sort_permutation(self, by_value='resolution', reverse=False)
structure_factors_from_map(self, map, in_place_fft=False, use_scale=False, anomalous_flag=None, use_sg=False)
structure_factors_from_scatterers(self, xray_structure, algorithm=None, cos_sin_table=False, grid_resolution_factor=0.3333333333333333, quality_factor=None, u_base=None, b_base=None, wing_cutoff=None, exp_table_one_over_step_size=None)
sys_absent_flags(self, integral_only=False)
Generate a boolean Miller array flagging those reflections which are
systematically absent under the current symmetry.
two_theta(self, wavelength, deg=False)
Generate a double Miller array containing the scattering angle of each
index.
unique_under_symmetry(self)
unique_under_symmetry_selection(self)
use_binner_of(self, other)
use_binning(self, binning)
use_binning_of(self, other)

Methods inherited from cctbx.crystal.symmetry:
as_cif_block(self)
as_py_code(self, indent='')
as_reference_setting(self)
asu_mappings(self, buffer_thickness, asu_is_inside_epsilon=None)
average_b_cart(self, b_cart)
average_u_cart(self, u_cart)
best_cell(self, angular_tolerance=None)
build_miller_set(self, anomalous_flag, d_min, d_max=None)
cell_equivalent_p1(self)
change_of_basis_op_to_best_cell(self, angular_tolerance=None, best_monoclinic_beta=True)
change_of_basis_op_to_inverse_hand(self)
change_of_basis_op_to_minimum_cell(self)
change_of_basis_op_to_niggli_cell(self, relative_epsilon=None, iteration_limit=None)
change_of_basis_op_to_primitive_setting(self)
change_of_basis_op_to_reference_setting(self)
direct_space_asu(self)
gridding(self, d_min=None, resolution_factor=None, step=None, symmetry_flags=None, mandatory_factors=None, max_prime=5, assert_shannon_sampling=True)
inverse_hand(self)
is_compatible_unit_cell(self)
is_patterson_symmetry(self)
is_similar_symmetry(self, other, relative_length_tolerance=0.01, absolute_angle_tolerance=1.0)
join_symmetry(self, other_symmetry, force=False, raise_sorry_if_incompatible_unit_cell=False)
miller_set(self, indices, anomalous_flag)
minimum_cell(self)
niggli_cell(self, relative_epsilon=None, iteration_limit=None)
primitive_setting(self)
space_group(self)
space_group_info(self)
special_position_settings(self, min_distance_sym_equiv=0.5, u_star_tolerance=0, assert_min_distance_sym_equiv=True)
subtract_continuous_allowed_origin_shifts(self, translation_cart)
unit_cell(self)

Data descriptors inherited from cctbx.crystal.symmetry:
__dict__
dictionary for instance variables (if defined)
__weakref__
list of weak references to the object (if defined)

 
Functions
       
as_hendrickson_lattman(...)
as_hendrickson_lattman( (bool)centric_flag, (complex)phase_integral, (float)max_figure_of_merit) -> tuple :
 
    C++ signature :
        cctbx::hendrickson_lattman<double> as_hendrickson_lattman(bool,std::complex<double>,double)
build_set(crystal_symmetry, anomalous_flag, d_min=None, d_max=None, max_index=None)
compute_cc_one_half(unmerged, n_trials=1)
Implementation of array.cc_one_half, assuming that the reflections are
already in the ASU.  Because the implementation uses random numbers, the
function has the option to calculate the mean over multiple trials.
is_unique_set_under_symmetry(...)
is_unique_set_under_symmetry( (space_group_type)space_group_type, (bool)anomalous_flag, (miller_index)miller_indices) -> bool :
 
    C++ signature :
        bool is_unique_set_under_symmetry(cctbx::sgtbx::space_group_type,bool,scitbx::af::const_ref<cctbx::miller::index<int>, scitbx::af::trivial_accessor>)
make_lookup_dict(indices)
map_to_asu(...)
map_to_asu( (space_group_type)arg1, (bool)arg2, (miller_index)arg3) -> None :
 
    C++ signature :
        void map_to_asu(cctbx::sgtbx::space_group_type,bool,scitbx::af::ref<cctbx::miller::index<int>, scitbx::af::trivial_accessor>)
 
map_to_asu( (space_group_type)arg1, (bool)arg2, (miller_index)arg3, (object)arg4) -> None :
 
    C++ signature :
        void map_to_asu(cctbx::sgtbx::space_group_type,bool,scitbx::af::ref<cctbx::miller::index<int>, scitbx::af::trivial_accessor>,scitbx::af::ref<double, scitbx::af::trivial_accessor>)
 
map_to_asu( (space_group_type)arg1, (bool)arg2, (miller_index)arg3, (object)arg4, (bool)arg5) -> None :
 
    C++ signature :
        void map_to_asu(cctbx::sgtbx::space_group_type,bool,scitbx::af::ref<cctbx::miller::index<int>, scitbx::af::trivial_accessor>,scitbx::af::ref<double, scitbx::af::trivial_accessor>,bool)
 
map_to_asu( (space_group_type)arg1, (bool)arg2, (miller_index)arg3, (complex_double)arg4) -> None :
 
    C++ signature :
        void map_to_asu(cctbx::sgtbx::space_group_type,bool,scitbx::af::ref<cctbx::miller::index<int>, scitbx::af::trivial_accessor>,scitbx::af::ref<std::complex<double>, scitbx::af::trivial_accessor>)
 
map_to_asu( (space_group_type)arg1, (bool)arg2, (miller_index)arg3, (hendrickson_lattman)arg4) -> None :
 
    C++ signature :
        void map_to_asu(cctbx::sgtbx::space_group_type,bool,scitbx::af::ref<cctbx::miller::index<int>, scitbx::af::trivial_accessor>,scitbx::af::ref<cctbx::hendrickson_lattman<double>, scitbx::af::trivial_accessor>)
map_to_asu_isym(...)
map_to_asu_isym( (space_group_type)arg1, (bool)arg2, (miller_index)arg3, (int)arg4) -> None :
 
    C++ signature :
        void map_to_asu_isym(cctbx::sgtbx::space_group_type,bool,scitbx::af::ref<cctbx::miller::index<int>, scitbx::af::trivial_accessor>,scitbx::af::ref<int, scitbx::af::trivial_accessor>)
multi_slice(...)
multi_slice( (miller_index)indices, (int)slice_axis, (int)slice_start, (int)slice_end) -> bool :
 
    C++ signature :
        scitbx::af::shared<bool> multi_slice(scitbx::af::const_ref<cctbx::miller::index<int>, scitbx::af::trivial_accessor>,unsigned int,int,int)
patterson_map(crystal_gridding, f_patt, f_000=None, sharpening=False, origin_peak_removal=False)
phase_transfer(...)
phase_transfer( (space_group)arg1, (miller_index)arg2, (object)arg3, (complex_double)arg4, (float)arg5) -> complex_double :
 
    C++ signature :
        scitbx::af::shared<std::complex<double> > phase_transfer(cctbx::sgtbx::space_group,scitbx::af::const_ref<cctbx::miller::index<int>, scitbx::af::trivial_accessor>,scitbx::af::const_ref<double, scitbx::af::trivial_accessor>,scitbx::af::const_ref<std::complex<double>, scitbx::af::trivial_accessor>,double)
 
phase_transfer( (space_group)arg1, (miller_index)arg2, (object)arg3, (object)arg4, (bool)arg5) -> complex_double :
 
    C++ signature :
        scitbx::af::shared<std::complex<double> > phase_transfer(cctbx::sgtbx::space_group,scitbx::af::const_ref<cctbx::miller::index<int>, scitbx::af::trivial_accessor>,scitbx::af::const_ref<double, scitbx::af::trivial_accessor>,scitbx::af::const_ref<double, scitbx::af::trivial_accessor>,bool)
randomize_amplitude_and_phase(...)
randomize_amplitude_and_phase( (complex_double)data, (bool)selection, (float)amplitude_error, (float)phase_error_deg, (int)random_seed) -> complex_double :
 
    C++ signature :
        scitbx::af::shared<std::complex<double> > randomize_amplitude_and_phase(scitbx::af::const_ref<std::complex<double>, scitbx::af::trivial_accessor>,scitbx::af::const_ref<bool, scitbx::af::trivial_accessor>,double,double,int)
raw_array_summary(array)
simple_slice(...)
simple_slice( (miller_index)indices, (int)slice_axis, (int)slice_index) -> bool :
 
    C++ signature :
        scitbx::af::shared<bool> simple_slice(scitbx::af::const_ref<cctbx::miller::index<int>, scitbx::af::trivial_accessor>,unsigned int,int)
statistical_mean(...)
statistical_mean( (space_group)arg1, (bool)arg2, (miller_index)arg3, (object)arg4) -> float :
 
    C++ signature :
        double statistical_mean(cctbx::sgtbx::space_group,bool,scitbx::af::const_ref<cctbx::miller::index<int>, scitbx::af::trivial_accessor>,scitbx::af::const_ref<double, scitbx::af::trivial_accessor>)
union_of_sets(miller_sets)
unique_under_symmetry_selection(...)
unique_under_symmetry_selection( (space_group_type)space_group_type, (bool)anomalous_flag, (miller_index)miller_indices) -> size_t :
 
    C++ signature :
        scitbx::af::shared<unsigned long> unique_under_symmetry_selection(cctbx::sgtbx::space_group_type,bool,scitbx::af::const_ref<cctbx::miller::index<int>, scitbx::af::trivial_accessor>)

 
Data
        Auto = <libtbx.AutoType object>
division = _Feature((2, 2, 0, 'alpha', 2), (3, 0, 0, 'alpha', 0), 8192)
fp_eps_double = 2.220446049250313e-16
generate_r_free_params_str = ' fraction = 0.1\n .type=float\n .shor... .short_caption = Number of resolution shells\n'