mmtbx.ions.identify module¶
Examines a structure for metal ions. Can iterate over all atoms, examining their density and chemical environment to determine if they are correctly identified, or if there are better candidate ions that they can be replaced with.
See mmtbx.ions.build to actually modify the structure, the code in this module handles ion identification, but only prints out messages to a log.
- class mmtbx.ions.identify.AtomProperties(i_seq, manager)¶
Bases:
object
Collect physical attributes of an atom, including B, occupancy, and map statistics, and track those which are at odds with its chemical identity.
- ANOM_PEAK = 7¶
- BAD_COORD_ATOM = 17¶
- BAD_COORD_RESIDUE = 19¶
- BAD_FPP = 18¶
- BAD_GEOMETRY = 9¶
- BAD_HALIDE = 21¶
- BAD_VALENCES = 12¶
- BAD_VECTORS = 11¶
- CLOSE_CONTACT = 24¶
- COORDING_GEOMETRY = 23¶
- FOFC_HOLE = 6¶
- FOFC_PEAK = 5¶
- HIGH_2FOFC = 22¶
- HIGH_B = 1¶
- HIGH_OCC = 3¶
- LIKE_COORD = 16¶
- LOW_B = 0¶
- LOW_OCC = 2¶
- NO_2FOFC_PEAK = 4¶
- NO_ANOM_PEAK = 8¶
- NO_GEOMETRY = 10¶
- TOO_FEW_COORD = 14¶
- TOO_FEW_NON_WATERS = 13¶
- TOO_MANY_COORD = 15¶
- VERY_BAD_VALENCES = 20¶
- atom_weight(manager)¶
Evaluates whether factors indicate that the atom is lighter, heavier, or isoelectric to what it is currently identified as.
- Parameters:
manager (mmtbx.ions.identify.manager)
- Returns:
-1 if lighter, 0 if isoelectronic, and 1 if heavier.
- Return type:
- check_fpp_ratio(ion_params, wavelength, fpp_ratio_min=0.3, fpp_ratio_max=1.05)¶
Compare the refined and theoretical f’’ values if available.
- Parameters:
ion_params (mmtbx.ions.metal_parameters)
wavelength (float)
fpp_ratio_min (float, optional)
fpp_ratio_max (float, optional)
- Returns:
f’’ / f’’_expected
- Return type:
- check_ion_environment(ion_params, wavelength=None, require_valence=True)¶
Checks whether or not the specified ion satisfies the metal-coordination parameters, specified by ion_params, such as valence sum, geometry, etc. The criteria used here are quite strict, but many of the analyses are saved for later if we want to use looser critera.
- Parameters:
ion_params (mmtbx.ions.metal_parameters)
wavelength (float, optional)
require_valence (bool, optional)
- error_strs = {0: 'Abnormally low b-factor', 1: 'Abnormally high b-factor', 2: 'Abnormally low occupancy', 3: 'Abnormally high occupancy', 4: 'No 2mFo-DFc map peak', 5: 'Peak in mFo-DFc map', 6: 'Negative peak in dFo-mFc map', 7: 'Peak in the anomalous map', 8: 'No peak in the anomalous map', 9: 'Unexpected geometry for coordinating atoms', 10: 'No distinct geometry for coordinating atoms', 11: 'VECSUM above cutoff', 12: 'BVS above or below cutoff', 13: 'Too few non-water coordinating atoms', 14: 'Too few coordinating atoms', 15: 'Too many coordinating atoms', 16: 'Like charge coordinating like', 17: 'Disallowed coordinating atom', 18: "Bad refined f'' value", 19: 'Disallowed coordinating residue', 20: 'BVS far above or below cutoff', 21: 'Bad halide site', 22: 'Unexpectedly high 2mFo-DFc value', 23: 'No distinct geometry and coordinating another atom with distinct geometry', 24: 'Close contact to oxygen atom'}¶
- get_atom_type(params, aggressive=False)¶
Checks the atom characteristics against what we would expect for a water. Updates self with any inaccuracies noticed (Surpringly low b-factor, high occupancy, etc). Note that HEAVY_ION does not necessarily rule out NA/MG, as these often have mFo-DFc peaks at high resolution.
- has_compatible_ligands(identity)¶
Indicates whether the coordinating atoms are of the allowed type (e.g. no N or S atoms coordinating CA, etc.) and residue (e.g. Ser is not an appropriate ligand for ZN).
- Parameters:
identity (mmtbx.ions.metal_parameters)
- Return type:
- identity(ion=None)¶
Covers an atom into a string representing its element and charge.
- Parameters:
ion (iotbx.pdb.hierarchy.atom, optional)
- Return type:
- is_compatible_anomalous_scattering(ion_params)¶
- is_compatible_site(ion_params, require_anom=True, ignore_valence=False)¶
More minimal criteria for determining whether a site is chemically compatible, allowing for incomplete coordination shells.
- Parameters:
ion_params (mmtbx.ions.metal_parameters)
require_anom (bool, optional)
ignore_valence (bool, optional)
- Return type:
- is_correctly_identified(identity=None)¶
Calculates whether factors indicate that the atom was correctly identified.
- Parameters:
identity (mmtbx.ions.metal_parameters)
- Return type:
- looks_like_nucleotide_phosphate_site(min_phosphate_oxygen_atoms=2, distance_cutoff=2.5)¶
Decide whether the atom is coordinating phosphate oxygens from a nucleotide, based on common atom names.
- number_of_atoms_within_radius(distance_cutoff)¶
Counts the number of coordinating atoms within a given radius.
- Parameters:
float
- Return type:
- number_of_backbone_oxygens(distance_cutoff=3.0)¶
Counts the number of backbone oxygens coordinating a site.
- number_of_favored_ligand_residues(ion_params, distance=3.0, exclude_atoms=())¶
Counts the number of preferred residues coordinating the atom. Used for approximate detection of transition-metal binding sites.
- Parameters:
ion_params (mmtbx.ions.metal_parameters)
distance (float, optional)
exclude_atoms (tuple of ...)
- Return type:
- show_ion_results(identity=None, out=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>, confirmed=False, valence_used=True)¶
Show statistics for a proposed element identity.
- Parameters:
identity (mmtbx.ions.metal_parameters)
out (file, optional)
confirmed (bool, optional)
valence_used (bool, optional)
- show_properties(identity, out=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>)¶
Show atomic properties that are independent of the suspected identity.
- Parameters:
identity (mmtbx.ions.metal_parameters)
out (file, optional)
- mmtbx.ions.identify.create_manager(pdb_hierarchy, geometry_restraints_manager, fmodel, wavelength, params, resolution_factor=0.25, nproc=<libtbx.AutoType object>, verbose=False, log=None, manager_class=None)¶
Wrapper around mmtbx.ions.identify.manager init method. Retrieves the connectivity and xray_structure from fmodel automatically.
- Parameters:
pdb_hierarchy (iotbx.pdb.hierarchy.root)
geometry_restraints_manager (cctbx.geometry_restraints.manager.manager)
fmodel (mmtbx.f_model.manager)
wavelength (float)
params (libtbx.phil.scope_extract)
resolution_factor (float, optional)
nproc (int, optional)
verbose (bool, optional)
log (file, optional)
manager_class (class, optional)
- Return type:
- mmtbx.ions.identify.find_anomalous_scatterers(*args, **kwds)¶
Wrapper for corresponding method in phaser.substructure, if phaser is available and configured.
- mmtbx.ions.identify.ion_master_phil()¶
- class mmtbx.ions.identify.ion_result(**kwds)¶
Bases:
atom
Class for validation results for an ion site. This includes information about both the atom modeled at that site as well as its chemical and scattering environment.
- altloc¶
- atom_selection¶
- b_iso¶
- chain_id¶
- chem_env¶
- element¶
- icode¶
- model_id¶
- name¶
- occupancy¶
- outlier¶
- resname¶
- resseq¶
- scatter_env¶
- score¶
- segid¶
- show(out=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>, prefix='')¶
- show_brief(out=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>, prefix='')¶
- symop¶
- xyz¶
- class mmtbx.ions.identify.manager(fmodel, pdb_hierarchy, xray_structure, params=None, wavelength=None, connectivity=None, nproc=1, verbose=False, log=None)¶
Bases:
object
Wrapper object for extracting local environment and adding or modifying scatterers.
- analyze_substructure(log=None, verbose=True)¶
Given a list of AX pseudo-atoms placed by Phaser, finds the nearest real atoms, and if possible determines their equivalence.
- Parameters:
log (file, optional)
verbose (bool, optional)
- analyze_water(i_seq, debug=True, candidates=<libtbx.AutoType object>, no_final=False, out=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>)¶
Examines the environment around a single atom to determine if it is actually a misidentified metal ion.
no_final is used internally, when candidates is not Auto, to try all of the Auto without selecting any as a final choice when none of the elements in candidates appear to probable ion identities.
- analyze_waters(out=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>, debug=True, candidates=<libtbx.AutoType object>)¶
Iterates through all of the waters in a model, examining the maps and local environment to check their identity and suggest likely ions where appropriate.
- find_nearby_atoms(i_seq, far_distance_cutoff=3.0, near_distance_cutoff=1.5, filter_by_bonding=True, filter_by_two_fofc=True)¶
Given site in the structure, return a list of contacts, optionally filtering by connectivity and 2mFo-DFc map level.
- get_b_iso(i_seq)¶
Calculates the isotropic b-factor of a site.
- get_fp(i_seq)¶
Retrieve the refined f’ for a site.
- get_fpp(i_seq)¶
Retrieve the refined f’’ for a site. Because this can come from either the built-in anomalous refinement or Phaser, it is handled differently than the f’ retrieval.
- get_initial_b_iso()¶
Calculates the isotropic b-factor to assign to newly labeled ions during the building process. Uses either the mean b-factor of waters or of all atoms if the former is unavailable.
- Return type:
- get_map(map_type)¶
Creates the real-space version of a given map type.
- Parameters:
map_type (str)
- Return type:
- get_map_gaussian_fit(map_type, i_seq)¶
Retrieves the two parameters of a gaussian function, fit to a given map at a site.
- get_map_sphere_variance(i_seq)¶
Calculate the density levels for points on a sphere around a given point. This will give us some indication whether the density falls off in the expected manner around a single atom.
- get_strict_valence_flag()¶
Checks whether the resolution of a model is high enough to use more strict thresholds for ion valences.
- Return type:
- guess_b_iso_real(i_seq)¶
Guess an approximate B_iso for an atom by averaging the values for coordinating atoms. This should partially compensate for the deflation of the B-factor due to incorrect scattering type, and consequently make the occupancy refinement more accurate.
- guess_molecular_weight(i_seq)¶
Guesses the molecular weight of a site by scaling the signal from the mFo map by the average signal from carbon atoms.
- looks_like_halide_ion(i_seq, element='CL', assume_hydrogens_all_missing=<libtbx.AutoType object>)¶
Given a site, analyze the nearby atoms (taking geometry into account) to guess whether it might be a halide ion. Among other things, halides will often be coordinated by amide hydrogens, and not in close contact with negatively charged atoms. Note that this procedure has a very high false positive rate when used by itself, so additional information about the electron count (map, occ, b) is absolutely essential.
- map_stats(i_seq)¶
Given a site in the structure, find the signal of the 2FoFc, FoFc, and anomalous map (When available).
- Parameters:
i_seq (int)
- Returns:
Object with .two_fofc, .fofc, and .anom properties for the associated signals in each map.
- Return type:
- principal_axes_of_inertia(i_seq)¶
Extracts the map grid points around a site, and calculates the axes of inertia (using the density values as weights). This is used to calculate the sphericity of the blob of density around a site.
- Parameters:
i_seq (int)
- Return type:
…
- refine_anomalous_substructure(log)¶
Run simple “substructure completion” implemented using CCTBX tools (the command-line equivalent is mmtbx.refine_anomalous_substructure). A faster alternative to Phaser, but not clear how effective this is at the moment.
- Parameters:
log (file)
- show_current_scattering_statistics(out=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>)¶
Prints out information about an entire model’s scattering statistics, such as mean map heights and b-factors for carbons and waters.
- Parameters:
out (file, optional)
- update_maps()¶
Generate new maps for the current structure, including anomalous map if data are anomalous.
If Phaser is installed and anomalous data are available, the anomalous log-likelihood gradient (LLG) map can be used instead of the conventional anomalous difference map. This is only really useful if the anomalous scattering of existing atoms is modeled (and ideally, refined).
- update_structure(pdb_hierarchy, xray_structure, connectivity=None, log=None, refine_if_necessary=True)¶
Set the current atomic data: PDB hierarchy, Xray structure, and simple connectivity list.
- Parameters:
pdb_hierarchy (iotbx.pdb.hierarchy.root)
xray_structure (cctbx.xray.structure.structure)
connectivity
log (file, optional)
refine_if_necessary (bool, optional)
- validate_ion(i_seq, out=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>, debug=True)¶
Examines one site in the model and determines if it was correctly modelled, returning a boolean indicating correctness.
- validate_ions(out=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>, debug=True, segid=None)¶
Iterate over all the ions built into the model by this module and double check their correctness. Looks for tell-tale signs such as negative peaks in the mFo-DFc and anomalous difference maps, abnormal b-factor, and disallowed chemical environments.
Prints out a table of bad ions and their information and returns a list of their sites’ i_seqs.
- water_selection()¶
Fetchs the selection for all waters in the model.
- Return type:
- class mmtbx.ions.identify.water_result(atom_props, filtered_candidates, matching_candidates, rejected_candidates, nuc_phosphate_site, atom_type, looks_like_halide, ambiguous_valence_cutoff, valence_used, final_choice, wavelength=None, no_final=False)¶
Bases:
object
Container for storing the results of manager.analyze_water for later display and retrieval.