cctbx.xray.structure module

class cctbx.xray.structure.conservative_pair_proxies(structure, bond_sym_table, conserve_angles)

Bases: object

class cctbx.xray.structure.meaningful_site_cart_differences(xs1, xs2)

Bases: object

Differences between the Cartesian coordinates of corresponding sites in two structures, cancelling continuous origin shifts if any.

This is especially useful to compare a refined structure to a reference structure as the former may have drifted along a continuous shift direction during refinement, therefore spoiling a naive comparison of corresponding sites.

max_absolute()
show()
class cctbx.xray.structure.scattering_type_registry_params(custom_dict=None, d_min=None, table=None, types_without_a_scattering_contribution=None)

Bases: object

class cctbx.xray.structure.structure(special_position_settings=None, scatterers=None, site_symmetry_table=None, non_unit_occupancy_implies_min_distance_sym_equiv_zero=False, scattering_type_registry=None, crystal_symmetry=None, wavelength=None)

Bases: special_position_settings

A class to describe and handle information related to a crystal structure.

It offers various methods to modify the crystallographic information contained.

Important members are:

  • .special_position_settings (base class)

  • .scatterers

  • .site_symmetry

  • .crystal_symmetry

add_scatterer(scatterer, site_symmetry_ops=None, insert_at_index=None)
add_scatterers(scatterers, site_symmetry_table=None, non_unit_occupancy_implies_min_distance_sym_equiv_zero=False)
adjust_occupancy(occ_max, occ_min, selection=None)

Adjust site occupancy factor for selected sites to be between occ_min and occ_max.

Parameters:
  • occ_max (float) – maximal site occupancy factor

  • occ_min (float) – minimal site occupancy factor

  • selection (boolean[]) – an array of bools to select scatterers to be adjusted

Returns:

none

adjust_u_iso()
all_selection()

Get a selector array for all scatterers of the structure.

Returns:

an array to select all scatterers of the structure

Return type:

boolean[]

apply_rigid_body_shift(rot, trans, selection=None, recompute_site_symmetries=True)
apply_rigid_body_shift_obj(sites_cart, sites_frac, rot, trans, selection, unit_cell, atomic_weights)
apply_shift(shift, recompute_site_symmetries=False)
apply_special_position_ops_d_target_d_site(d_target_d_site)
apply_symmetry_sites()
apply_symmetry_u_stars()
as_cif_block(covariance_matrix=None, cell_covariance_matrix=None, format='mmCIF')
as_cif_simple(out=None, data_name='global', format='mmCIF')
as_emma_model()
as_pdb_file(remark=None, remarks=[], fractional_coordinates=False, resname=None, connect=None)
as_py_code(indent='')

eval(self.as_py_code()) is similar (but sometimes not identical) to self and meant to enable quick formatting of self as Python code for unit tests.

asu_mappings(buffer_thickness, asu_is_inside_epsilon=None)
asymmetric_unit_in_p1()
atomic_weights()
b_iso_min_max_mean()

Get the minimal, maximal and mean isotropic Debye-Waller/temperature/B factors of all scatterers in this structure.

Returns:

minimal b_iso, maximal b_iso, mean b_iso

Return type:

float, float, float

b_iso_or_b_equiv()
by_index_selection(selected_scatterers)

Get a selector array for scatterers with specified index from the structure. For example you can select scatterers 1,4 and 5 by passing (1,4,5) as argument.

Parameters:

selected_scatterers (list(int)) – list of scatterers to select

Returns:

an array to select the desired scatterers of the structure

Return type:

flex.bool[]

center_of_mass(atomic_weights=None)
change_basis(cb_op)
change_hand()
closest_distances(sites_frac, distance_cutoff, use_selection=None)
concatenate(other)
concatenate_inplace(other)
conservative_pair_proxies(bond_sym_table, conserve_angles)
convert_to_anisotropic(selection=None)
convert_to_isotropic(selection=None)
coordinate_degrees_of_freedom_counts(selection=None)
crystal_density()

Get the value of the diffraction-determined density for the crystal, suitable for the CIF item _exptl_crystal_density_diffrn

Density values are calculated from the crystal cell and contents. The units are megagrams per cubic metre (=grams per cubic centimetre).

Equivalent to:

1.66042 * _chemical_formula_weight * _cell_formula_units_Z / _cell_volume

Returns:

chemical density in megagrams per cubic metre (=grams per cm^3)

Return type:

float

crystal_symmetry()

Get crystal symmetry of the structure

Returns:

a new crystal symmetry object

Return type:

cctbx.crystal.symmetry

cubic_unit_cell_around_centered_scatterers(buffer_size)
customized_copy(crystal_symmetry=Keep, unit_cell=Keep, space_group_info=Keep, non_unit_occupancy_implies_min_distance_sym_equiv_zero=Keep, wavelength=Keep)
deep_copy_scatterers()

Create a deep copy of the structure with all scatterers

Returns:

a new cctbx.xray.structure object

Return type:

cctbx.xray.structure

difference_vectors_cart(other)
discard_scattering_type_registry()
distances(other, selection=None)

Calculates pairwise distances between the atoms of this structure and another structure with the same number of scatterers.

Parameters:
  • other (cctbx.xray.structure) – the other structure

  • selection (boolean[]) – an array of bools to select scatterers to be taken into calculation

Returns:

an array of distances for the selected scatterers

Return type:

float[]

element_selection(*elements)

Get a selector array for scatterers of specified element type(s) of the structure.

Parameters:

elements (list(string) or set(string) or tuple(string)) – tuple of element symbols to select

Returns:

an array to select the desired scatterers of the structure

Return type:

boolean[]

erase_scatterers()

Remove all scatterers from structure

Returns:

none

expand_to_p1(append_number_to_labels=False, sites_mod_positive=False)

Get the current structure expanded into spacegroup P1. This turns all symmetry induced scatterers into independent individual scatterers. The expanded structure may have sites with negative or > 1.0 coordinates. Use ‘.sites_mod_positive()’ or ‘.sites_mod_short()’ on the result, or alternatively set sites_mod_positive to ‘True’ in case you want to change this behaviour.

Parameters:
  • append_number_to_labels (boolean) – If set to ‘True’ scatterers generated from symmetry will be labelled with a numerical suffix

  • sites_mod_positive (boolean) – If set to ‘True’ xyz coordinates of the scatterers will be kept inside [0,1[

Returns:

a new instance of the structure expanded into P1

Return type:

cctbx.xray.structure

extract_u_cart_plus_u_iso()
extract_u_iso_or_u_equiv()
f_000(include_inelastic_part=False)

Get the effective number of electrons in the crystal unit cell contributing to F(000), suitable for the CIF item _exptl_crystal_F_000.

According to the CIF definition, this item may contain dispersion contributions.

Parameters:

include_inelastic_part (boolean) – If ‘True’ contributions due to dispersion are included in F(000).

Returns:

F(000)

Return type:

float

classmethod from_cif(file_object=None, file_path=None, data_block_name=None)
classmethod from_shelx(*args, **kwds)
grads_and_curvs_target_simple(miller_indices, da_db, daa_dbb_dab)
guess_scattering_type_is_a_mixture_of_xray_and_neutron()
guess_scattering_type_neutron()
hd_selection()

Get a selector array for all hydrogen and deuterium scatterers of the structure.

Returns:

an array to select all H and D scatterers of the structure

Return type:

boolean[]

intersection_of_scatterers(i_seq, j_seq)

For a pair of scatterers, calculates their overlap given the coordinates and displacement parameters (using adptbx.intersection).

is_positive_definite_u(u_cart_tolerance=None)
is_similar(other, eps=1e-06)

Check if similar. Can be endlessly fortified as needed.

make_scatterer_labels_shelx_compatible_in_place()
max_distance(other, selection=None)

Calculates the maximum pairwise distance between the atoms of this structure and another structure with the same number of scatterers.

Parameters:
  • other (cctbx.xray.structure) – the other structure

  • selection (boolean[]) – an array of bools to select scatterers to be taken into calculation

Returns:

the maximum distance of two corresponding scatterers out of the selected scatterers

Return type:

float

mean_distance(other, selection=None)

Calculates the arithmetic mean pairwise distance between the atoms of this structure and another structure with the same number of scatterers.

Parameters:
  • other (cctbx.xray.structure) – the other structure

  • selection (boolean[]) – an array of bools to select scatterers to be taken into calculation

Returns:

the mean pairwise distance of the selected scatterers

Return type:

float

mean_scattering_density()
min_distance(other, selection=None)

Calculates the minimum pairwise distance between the atoms of this structure and another structure with the same number of scatterers.

Parameters:
  • other (cctbx.xray.structure) – the other structure

  • selection (boolean[]) – an array of bools to select scatterers to be taken into calculation

Returns:

the minimum distance of two corresponding scatterers out of the selected scatterers

Return type:

float

min_u_cart_eigenvalue()
n_grad_u_aniso()
n_grad_u_iso()
n_parameters(considering_site_symmetry_constraints=False)
n_undefined_multiplicities()
non_unit_occupancy_implies_min_distance_sym_equiv_zero()
orthorhombic_unit_cell_around_centered_scatterers(buffer_size)
pair_asu_table(distance_cutoff=None, asu_mappings_buffer_thickness=None, asu_is_inside_epsilon=None, min_cubicle_edge=5)
pair_sym_table_show(pair_sym_table, is_full_connectivity=False, out=None)
pair_sym_table_show_distances(pair_sym_table, show_cartesian=False, skip_j_seq_less_than_i_seq=False, skip_sym_equiv=False, out=None)
parameter_map()
principal_axes_of_inertia(atomic_weights=None)
random_remove_sites_selection(fraction)
random_shift_sites(max_shift_cart=0.2)
re_apply_symmetry(i_scatterer)
replace_scatterers(scatterers, site_symmetry_table='existing')
replace_sites_cart(new_sites, selection=None)
replace_sites_frac(new_sites, selection=None)
rms_difference(other)
scale_adp(factor, selection=None)

Scale the atomic displacement parameters of the selected scatterers of the structure with the specified factor. If no selection is given, all scatterers will be handled as if selected.

Parameters:
  • factor (float) – scale factor to apply to the adps of the selected scatterers

  • selection (boolean[]) – an array of bools to select the scatterers to have their adps scaled

Returns:

none

scale_adps(scale_factor)
scatterer_flags()
scatterers()

Get all scatterers of the structure

Returns:

a reference to an array of cctbx.xray.scatterer

Return type:

cctbx.xray.scatterer[]

scattering_dictionary_as_string()
scattering_type_registry(custom_dict=None, d_min=None, table=None, types_without_a_scattering_contribution=None, explicitly_allow_mixing=False)
scattering_types()
scattering_types_counts_and_occupancy_sums()
select(selection, negate=False)
select_inplace(selection)
selection_within(radius, selection)
set_b_iso(value=None, values=None, selection=None)

Set isotropic Debye-Waller/temperature/B factors with automatic conversion to u_iso

Parameters:
  • value (double) – a single double value to set all b_iso of selected scatterers to

  • values (double[]) – an array of double values to set all b_iso of selected scatterers to

  • selection (boolean[]) – an array of bools to select scatterers to be updated with new b_iso values

Returns:

the modified base object

Return type:

cctbx.xray.structure

set_custom_inelastic_form_factors(table, set_use_fp_fdp=True, source='custom')

Expects a dictionary of tuples like ‘C’ : (fp, fdp). If an element is not in the dictionary, the fp and fdp are reset to 0 and the use_fp_fdp is set to false.

set_fdps(value, selection=None)
set_fps(value, selection=None)
set_inelastic_form_factors(photon, table, set_use_fp_fdp=True)
set_non_unit_occupancy_implies_min_distance_sym_equiv_zero(value)
set_occupancies(value, selection=None)
set_scatterer_flags(scatterer_flags)
set_sites_cart(sites_cart)

Set the the cartesian coordinates of all sites of the structure to ‘sites_cart’.

Parameters:

sites_cart (scitbx.array_family.flex.vec3_double) – a list of the cartesian coordinates for all scatterers

Returns:

none

set_sites_frac(sites_frac)

Set the the fractional coordinates of all sites of the structure to ‘sites_frac’.

Parameters:

sites_frac (scitbx.array_family.flex.vec3_double) – a list of the fractional coordinates for all scatterers

Returns:

none

set_u_cart(u_cart, selection=None)
set_u_iso(value=None, values=None, selection=None)

Set isotropic mean thermic displacements of scatterers

Parameters:
  • value (double) – a single double value to set all u_iso of selected scatterers to

  • values (double[]) – an array of double values to set all u_iso of selected scatterers to

  • selection (boolean[]) – an array of bools to select scatterers to be updated with new u_iso values

Returns:

the modified base object

Return type:

cctbx.xray.structure

shake_adp(b_max=None, b_min=None, spread=10.0, aniso_spread=0.1, keep_anisotropic=False, random_u_cart_scale=1.0, selection=None)
shake_adp_if_all_equal(b_iso_tolerance=0.1)
shake_fdps(selection=None)
shake_fps(selection=None)
shake_occupancies(selection=None)
shake_sites_in_place(rms_difference=None, mean_distance=None, selection=None, allow_all_fixed=False, random_double=None)

Shake the coordinates of the selected scatterers in this structure.

Parameters:
  • rms_difference (float) – radial mean square displacement (>=0) to apply to selected scatterers

  • mean_distance (float) – a mean distance shift (>=0) to apply to selected scatterers

  • selection (boolean[]) – an array of bools to select scatterers to be shaken

  • allow_all_fixed (boolean) – if set to ‘True’ shaking a structure with all scatterers on fixed special positions will not cause an error

  • random_double (float[]) – “random” numbers to use for displacements

Returns:

‘True’ if at least one scatterer was moved, ‘False’ otherwise

Return type:

boolean

shift_occupancies(q_shift, selection=None)
shift_sites_in_place(shift_length, mersenne_twister=None)

Shifts the coordinates of all scatterers in this structure.

Parameters:
  • shift_length (float) – the distance to shift each scatterer with

  • mersenne_twister (flex.mersenne_twister) – a mersenne twister to use as entropy source

Returns:

none

shift_us(u_shift=None, b_shift=None, selection=None)
show_angles(distance_cutoff=None, asu_mappings_buffer_thickness=None, asu_is_inside_epsilon=None, pair_asu_table=None, keep_pair_asu_table=False, out=None)
show_dihedral_angles(distance_cutoff=None, max_d=1.7, max_angle=170, asu_mappings_buffer_thickness=None, asu_is_inside_epsilon=None, pair_asu_table=None, keep_pair_asu_table=False, out=None)
show_distances(distance_cutoff=None, asu_mappings_buffer_thickness=None, asu_is_inside_epsilon=None, min_cubicle_edge=5, pair_asu_table=None, show_cartesian=False, keep_pair_asu_table=False, out=None)
show_scatterer_flags_summary(out=None)
show_scatterers(f=None, special_positions_only=False)
show_special_position_shifts(sites_frac_original=None, sites_cart_original=None, out=None, prefix='')
show_summary(f=None, prefix='')
show_u_statistics(text='', out=None, use_hydrogens=False)
site_symmetry_table()
sites_cart()

Get the the cartesian coordinates of all sites of the structure.

Returns:

a list of the sites of the structure in cartesian coordinates

Return type:

scitbx.array_family.flex.vec3_double

sites_frac()
sites_mod_positive()

Get the current structure converted into a structure with x,y,z of all scatterers in the interval [0,1[

Returns:

the same instance of the structure with only posive coordinates of its scatterers

Return type:

cctbx.xray.structure

sites_mod_short()

Get the current structure converted into a structure with short coordinates vectors of all scatterers

Returns:

the same instance of the structure with only short coordinates vectors of its scatterers

Return type:

cctbx.xray.structure

sort(by_value='occupancy', reverse=False)
special_position_indices()
structure_factors(anomalous_flag=None, d_min=None, algorithm=None, cos_sin_table=False, quality_factor=None, u_base=None, b_base=None, wing_cutoff=None)

Calculate structure factors for the current scatterers using either direct summation or FFT method; by default the appropriate method will be guessed automatically.

Parameters:
  • anomalous_flag (bool) – toggles whether the returned structure factors are anomalous or not

  • d_min (float) – resolution cutoff (required)

  • algorithm (str or None) – specifies ‘direct’ or ‘fft’, or if None (default), the algorithm will be chosen automatically

  • cos_sin_table (bool) – if True, uses interpolation of values from a pre-calculated lookup table for trigonometric function calls in the structure factor calculations in preference to calling the system libraries

  • quality_factor (float) – determines accuracy of sampled density in fft method

  • u_base (float) – additional smearing B-factor for scatterers which have too small Bs so they fall between the grid nodes

  • b_base (float) – same as u_base (but 8*pi^2 larger)

  • wing_cutoff (float) – is how far away from atomic center you sample density around atom

Returns:

a custom Python object (exact type depends on method used), from which f_calc() may be called

Return type:

derived from cctbx.xray.structure_factors.manager.managed_calculation_base

switch_to_neutron_scattering_dictionary()
tidy_us(u_min=1e-06, u_max=12.665021303812669, anisotropy_min=0.25)

Clean up atomic displacements so they fall within a sensible range (this is especially important when writing out PDB format).

translate(x=0, y=0, z=0)

Translates all scatterers of this structure by x,y,z.

Parameters:
  • x (float) – x component of the translation vector

  • y (float) – y component of the translation vector

  • z (float) – z component of the translation vector

Returns:

a new translated copy of the structure

Return type:

cctbx.xray.structure

truncate_at_pdb_format_precision()
unit_cell_content(omit=None)

The content of the unit cell as a chemical formula

use_u_aniso()
use_u_iso()