iotbx.pdb¶
API documentation¶
Input methods and objects¶
- iotbx.pdb.input(file_name=None, source_info=<class 'iotbx.pdb.Please_pass_string_or_None'>, lines=None, pdb_id=None, raise_sorry_if_format_error=False)¶
Main input method for both PDB and mmCIF files; will automatically determine the actual format and return the appropriate data type.
- Parameters:
file_name (path to PDB or mmCIF file)
source_info (string describing source of input (e.g. file name))
lines (flex.std_string array of input lines)
pdb_id (PDB ID to automatically retrieve from local mirror)
raise_sorry_if_format_error (re-raise any low-level parser errors as a) – libtbx.utils.Sorry exception instance for clean user feedback
- Returns:
An object representing the result of parsing, including an array of atom
objects; the actual class will differ depending on the input format. Much of
the API will be the same in either case.
- class iotbx.pdb.ext.input¶
Bases:
instance
This class parses PDB format, including non-ATOM records. Atom objects will be created as part of the parsing, but the full PDB hierarchy object requires calling the construct_hierarchy() method.
- as_pdb_string(crystal_symmetry=<libtbx.AutoType object>, cryst1_z=<libtbx.AutoType object>, write_scale_records=True, append_end=False, atom_hetatm=True, sigatm=True, anisou=True, siguij=True, cstringio=None, link_records=<libtbx.AutoType object>, return_cstringio=<libtbx.AutoType object>)¶
Generate standard PDB format. Will use built-in crystal symmetry if available.
- atoms((input)arg1) af_shared_atom : ¶
- C++ signature :
scitbx::af::shared<iotbx::pdb::hierarchy::atom> atoms(iotbx::pdb::input {lvalue})
- atoms_with_labels((input)arg1) af_shared_atom_with_labels : ¶
- C++ signature :
scitbx::af::shared<iotbx::pdb::hierarchy::atom_with_labels> atoms_with_labels(iotbx::pdb::input {lvalue})
- bookkeeping_section((input)arg1) std_string : ¶
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > bookkeeping_section(iotbx::pdb::input {lvalue})
- break_indices((input)arg1) size_t : ¶
- C++ signature :
scitbx::af::shared<unsigned long> break_indices(iotbx::pdb::input {lvalue})
- chain_indices((input)arg1) stl_vector_unsigned : ¶
- C++ signature :
scitbx::af::shared<std::vector<unsigned int, std::allocator<unsigned int> > > chain_indices(iotbx::pdb::input {lvalue})
- connectivity_annotation_section((input)arg1) std_string : ¶
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > connectivity_annotation_section(iotbx::pdb::input {lvalue})
- connectivity_section((input)arg1) std_string : ¶
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > connectivity_section(iotbx::pdb::input {lvalue})
- construct_hierarchy((input)arg1[, (object)residue_group_post_processing=True[, (object)set_atom_i_seq=True[, (object)sort_atoms=True]]]) root : ¶
- C++ signature :
iotbx::pdb::hierarchy::root construct_hierarchy(iotbx::pdb::input {lvalue} [,bool=True [,bool=True [,bool=True]]])
- construct_hierarchy_BIOMT_expanded(sort_atoms=True)¶
- construct_hierarchy_MTRIX_expanded(sort_atoms=True)¶
- crystal_symmetry(crystal_symmetry=None, weak_symmetry=False)¶
- crystal_symmetry_from_cryst1()¶
- crystallographic_section((input)arg1) std_string : ¶
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > crystallographic_section(iotbx::pdb::input {lvalue})
- deposition_date(us_style=True)¶
Placeholder to match mmCIF functionality. Probably could parse REVDAT.
- extract_LINK_records()¶
Collect link records from PDB file
- extract_authors()¶
- extract_connectivity()¶
Parse CONECT records and extract the indices of bonded atoms. Returns a scitbx.array_family.shared.stl_set_unsigned object corresponding to the atoms array, with each element being the list of indices of bonded atoms (if any). If no CONECT records are found, returns None.
Note that the ordering of atoms may be altered by construct_hierarchy(), so this method should probably be called after the hierarchy is made.
- extract_cryst1_z_columns()¶
- extract_f_model_core_constants()¶
- extract_header_year()¶
- extract_remark_iii_records(iii)¶
- extract_secondary_structure(log=None)¶
- extract_tls_params(hierarchy)¶
- extract_wavelength(first_only=True)¶
- file_type()¶
- get_experiment_type()¶
- get_link_records()¶
- get_matthews_coeff()¶
- get_program_name()¶
- get_r_rfree_sigma(file_name=None)¶
- get_restraints_used()¶
- get_solvent_content()¶
- heterogen_section((input)arg1) std_string : ¶
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > heterogen_section(iotbx::pdb::input {lvalue})
- miscellaneous_features_section((input)arg1) std_string : ¶
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > miscellaneous_features_section(iotbx::pdb::input {lvalue})
- model_atom_counts((input)arg1) size_t : ¶
- C++ signature :
scitbx::af::shared<unsigned long> model_atom_counts(iotbx::pdb::input {lvalue})
- model_ids((input)arg1) std_string : ¶
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > model_ids(iotbx::pdb::input {lvalue})
- model_indices((input)arg1) size_t : ¶
- C++ signature :
scitbx::af::shared<unsigned long> model_indices(iotbx::pdb::input {lvalue})
- primary_structure_section((input)arg1) std_string : ¶
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > primary_structure_section(iotbx::pdb::input {lvalue})
- process_BIOMT_records()¶
- process_MTRIX_records()¶
- record_type_counts((input)arg1) dict : ¶
- C++ signature :
boost::python::dict record_type_counts(iotbx::pdb::input)
- remark_section((input)arg1) std_string : ¶
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > remark_section(iotbx::pdb::input {lvalue})
- resolution()¶
- scale_matrix()¶
- secondary_structure_section((input)arg1) std_string : ¶
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > secondary_structure_section(iotbx::pdb::input {lvalue})
- sequence_from_SEQRES()¶
- source_info((input)arg1) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > source_info(iotbx::pdb::input {lvalue})
- special_position_settings(special_position_settings=None, weak_symmetry=False, min_distance_sym_equiv=0.5, u_star_tolerance=0)¶
- ter_indices((input)arg1) size_t : ¶
- C++ signature :
scitbx::af::shared<unsigned long> ter_indices(iotbx::pdb::input {lvalue})
- title_section((input)arg1) std_string : ¶
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > title_section(iotbx::pdb::input {lvalue})
- unknown_section((input)arg1) std_string : ¶
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > unknown_section(iotbx::pdb::input {lvalue})
- used_amber_restraints()¶
- used_cdl_restraints()¶
- used_omega_cdl_restraints()¶
- write_pdb_file(file_name, open_append=False, crystal_symmetry=<libtbx.AutoType object>, cryst1_z=<libtbx.AutoType object>, write_scale_records=True, append_end=False, atom_hetatm=True, sigatm=True, anisou=True, siguij=True)¶
- xray_structure_simple(crystal_symmetry=None, weak_symmetry=False, cryst1_substitution_buffer_layer=None, unit_cube_pseudo_crystal=False, fractional_coordinates=False, use_scale_matrix_if_available=True, min_distance_sym_equiv=0.5, non_unit_occupancy_implies_min_distance_sym_equiv_zero=True, scattering_type_exact=False, enable_scattering_type_unknown=False, atom_names_scattering_type_const=['PEAK', 'SITE'])¶
Create a single cctbx.xray.structure object from the atom records, using only the first model found.
- xray_structures_simple(one_structure_for_each_model=True, crystal_symmetry=None, weak_symmetry=False, cryst1_substitution_buffer_layer=None, unit_cube_pseudo_crystal=False, fractional_coordinates=False, min_distance_sym_equiv=0.5, non_unit_occupancy_implies_min_distance_sym_equiv_zero=True, use_scale_matrix_if_available=True, scattering_type_exact=False, enable_scattering_type_unknown=False, atom_names_scattering_type_const=['PEAK', 'SITE'])¶
Create a list of cctbx.xray.structure objects, one per model in the input file. Note that for most single-model structures (i.e. nearly all crystal structures), this will be a single-item list.
- class iotbx.pdb.mmcif.cif_input(file_name=None, cif_object=None, source_info=<class 'iotbx.pdb.Please_pass_string_or_None'>, lines=None, pdb_id=None, raise_sorry_if_format_error=False)¶
Bases:
pdb_input_mixin
- atoms()¶
- atoms_with_labels()¶
- connectivity_annotation_section()¶
- connectivity_section()¶
- construct_hierarchy(set_atom_i_seq=True, sort_atoms=True)¶
- crystal_symmetry(crystal_symmetry=None, weak_symmetry=False)¶
- crystal_symmetry_from_cryst1()¶
- crystallographic_section()¶
- deposition_date(us_style=True)¶
Placeholder to match mmCIF functionality. Probably could parse REVDAT.
- extract_cryst1_z_columns()¶
- extract_f_model_core_constants()¶
- extract_header_year()¶
- extract_secondary_structure(log=None)¶
- extract_tls_params(hierarchy)¶
- extract_wavelength(first_only=True)¶
- file_type()¶
- get_experiment_type()¶
- get_matthews_coeff()¶
- get_program_name()¶
- get_r_rfree_sigma(file_name=None)¶
- get_restraints_used()¶
- get_solvent_content()¶
- heterogen_section()¶
- model_ids()¶
- model_indices()¶
- process_BIOMT_records()¶
- process_MTRIX_records()¶
Read MTRIX records from a pdb file
Arguments:¶
Returns:¶
- resultobject containing information on all NCS operations
iotbx.mtrix_biomt.container
- remark_section()¶
- resolution()¶
- scale_matrix()¶
- source_info()¶
- ter_indices()¶
- title_section()¶
- used_amber_restraints()¶
- used_cdl_restraints()¶
- used_omega_cdl_restraints()¶
PDB hierarchy objects¶
- class iotbx.pdb.hierarchy.root¶
Bases:
instance
Root node of the PDB hierarchy object. This is returned by the method construct_hierarchy() of the PDB/mmCIF input objects, but it may also be created programatically. Note that it does not contain any reference to crystal symmetry or source scattering information, meaning that in practice it must often be tracked alongside an equivalent cctbx.xray.structure object. Pickling is supported, simply by writing out and reading back the PDB-format representation of the hierarchy.
Examples
>>> hierarchy = iotbx.pdb.hierarchy.root()
- adopt_xray_structure(xray_structure)¶
Apply the current (refined) atomic parameters from the cctbx.xray.structure object to the atoms in the PDB hierarchy.
- altloc_indices((root)arg1) object : ¶
- C++ signature :
boost::python::api::object altloc_indices(iotbx::pdb::hierarchy::root)
- altlocs_present(skip_blank=True)¶
- append_model((root)arg1, (model)model) None : ¶
- C++ signature :
void append_model(iotbx::pdb::hierarchy::root {lvalue},iotbx::pdb::hierarchy::model {lvalue})
- apply_atom_selection(atom_selection)¶
Apply atom selection string and return deep copy with selected atoms
- apply_rotation_translation(rot_matrices, trans_vectors)¶
LIMITATION: ANISOU records in resulting hierarchy will be invalid!!!
- as_cif_block(crystal_symmetry=None, coordinate_precision=5, occupancy_precision=3, b_iso_precision=5, u_aniso_precision=5, segid_as_auth_segid=False, output_break_records=False)¶
- as_dict_of_chain_id_resseq_as_int_residue_names()¶
- as_forward_compatible_hierarchy(conversion_info=None)¶
Convert a standard hierarchy to a forward_compatible_hierarchy
:params conversion_info
- :returns pdb_hierarchy with chain ID and residue names converted to
strings compatible with PDB formatting. Returned hierarchy is a deep_copy and contains the attribute _conversion_info with the conversion_info used
- Typical use: running a method with cmd_text (text commands) and supplying
a hierarchy (or string from it). Convert the hierarchy and the cmd_text, run the method, convert the results back:
# Convert the hierarchy to forward compatible ph_fc = ph.as_forward_compatible_hierarchy()
# Convert any commands. Can be done one word at a time also cmd_text_fc = ph_fc.convert_multi_word_text_to_forward_compatible(cmd_text)
# Run the method with converted hierarchy and commands result = do_something(ph = ph_fc, command_text = cmd_text_fc)
# Convert back any resulting hierarchy new_ph = result.ph.forward_compatible_hierarchy_as_standard(
conversion_info = ph_fc.conversion_info())
# Convert back any words in the results that referred to converted # items. Keys are chain_id and resname new_result_items = [] for result_item,key in zip(results.text_words, results.text_keys):
- new_result_item = ph_fc.conversion_info(). get_full_text_from_forward_compatible_pdb_text(key = key,
forward_compatible_pdb_text = result_item)
new_result_items.append(new_result_item)
- as_forward_compatible_string(**kw)¶
- Create a forward_compatible PDB string from a hierarchy and
throw away the conversion information.
One-way conversion useful for creating a file that is in PDB format.
- Params **kw**kw:
any params suitable for as_pdb_string()
:returns text string
- as_list_of_residue_names()¶
- as_mmcif_string(crystal_symmetry=None, data_block_name=None, segid_as_auth_segid=False, output_break_records=False)¶
- as_model_manager(crystal_symmetry, unit_cell_crystal_symmetry=None, shift_cart=None)¶
Returns simple version of model object based on this hierarchy Expects but does not require crystal_symmetry. Optional unit_cell_crystal_symmetry and shift_cart.
Note that if crystal_symmetry is not supplied, this uses a text representation of the hierarchy so that values for xyz, occ, b, and crystal_symmetry are all rounded.
- as_pdb_input(crystal_symmetry=None, segid_as_auth_segid=True)¶
Generate corresponding pdb.input object. Note that this uses a text representation of the hierarchy so that values for xyz, occ, b, and crystal_symmetry are all rounded.
- as_pdb_or_mmcif_string(target_format=None, segid_as_auth_segid=True, remark_section=None, **kw)¶
- Shortcut for pdb_or_mmcif_string_info with write_file=False, returning
only the string representing this hierarchy. The string may be in PDB or mmCIF format, with target_format used if it is feasible.
Method to allow shifting from general writing as pdb to writing as mmcif, with the change in two places (here and model.py) Use default of segid_as_auth_segid=True here (different than
as_mmcif_string())
- Parameters:
target_format – desired output format, pdb or mmcif
segid_as_auth_segid – use the segid in hierarchy as the auth_segid in mmcif output
remark_section – if supplied and format is pdb, add this text
**kw –
any keywords suitable for as_pdb_string() and as_mmcif_string()
:returns text string representing this hierarchy
- as_pdb_string(crystal_symmetry=None, cryst1_z=None, write_scale_records=True, append_end=False, interleaved_conf=0, atoms_reset_serial_first_value=None, atom_hetatm=True, sigatm=True, anisou=True, siguij=True, output_break_records=True, force_write=False, cstringio=None, return_cstringio=<libtbx.AutoType object>)¶
Generate complete PDB-format string representation. External crystal symmetry is strongly recommended if this is being output to a file.
- Parameters:
crystal_symmetry – cctbx.crystal.symmetry object or equivalent (such as an xray.structure object or Miller array)
write_scale_records – write fractional scaling records (SCALE) if crystal symmetry is provided
anisou – write ANISOU records for anisotropic atoms
sigatm – write SIGATM records if applicable
siguij – write SIGUIJ records if applicable
force_write – write even if it does not fit in pdb format
- Returns:
Python str
- as_sequence(substitute_unknown='X', substitute_unknown_na='N', ignore_all_unknown=None, as_string=False, only_one_model=False)¶
Uses chain.as_sequence() for all chains and returns the catenation :param substitute_unknown: character to use for unrecognized 3-letter codes :param substitute_unknown_na: character to use for unrecognized na codes :param ignore_all_unknown: set substitute_unknown and substitute_unknown_na to ‘’ :param as_string: return string (default is to return list) :param only_one_model: Only use the first model if more than one
- as_str(prefix='', level_id=None, level_id_exception=<class 'ValueError'>)¶
Alias for show().
- atom_groups()¶
Iterate over all atom groups (by model, then chain, then residue group)
- atom_selection_cache(special_position_settings=None)¶
- atoms((root)arg1[, (object)interleaved_conf=0]) af_shared_atom : ¶
- C++ signature :
scitbx::af::shared<iotbx::pdb::hierarchy::atom> atoms(iotbx::pdb::hierarchy::root {lvalue} [,int=0])
- atoms_reset_serial((root)arg1[, (object)interleaved_conf=0[, (object)first_value=1]]) None : ¶
- C++ signature :
void atoms_reset_serial(iotbx::pdb::hierarchy::root {lvalue} [,int=0 [,int=1]])
- atoms_size((root)arg1) int : ¶
- C++ signature :
unsigned int atoms_size(iotbx::pdb::hierarchy::root {lvalue})
- atoms_with_i_seq_mismatch((root)arg1) af_shared_atom : ¶
- C++ signature :
scitbx::af::shared<iotbx::pdb::hierarchy::atom> atoms_with_i_seq_mismatch(iotbx::pdb::hierarchy::root {lvalue})
- atoms_with_labels()¶
Generator for atom_with_labels objects, presented in the same order as the array returned by the atoms() method.
- chain_ids(unique_only=False)¶
Get list of chain IDS, return unique set if unique_only=True
- chain_type()¶
Inspect residue names and counts to determine what chain types are present If only one chain type, return it. Otherwise return None
- chain_types()¶
Inspect residue names and counts to determine what chain types are present
- chains()¶
Iterate over all chains in all models.
- chunk_selections(residues_per_chunk)¶
- clear_label_asym_id_lookups()¶
Make sure we have fresh lookups in case the hierarchy was modified since they were calculated.
- composition()¶
- contains_break_records()¶
- contains_dna(oc=None)¶
Inspect residue names and counts to determine if any of them are DNA.
- contains_hetero()¶
- contains_nucleic_acid(min_content=0, oc=None)¶
Inspect residue names and counts to determine if enough of them are RNA or DNA.
- contains_protein(min_content=0, oc=None)¶
Inspect residue names and counts to determine if enough of them are protein.
- contains_rna(oc=None)¶
Inspect residue names and counts to determine if any of them are RNA.
- conversion_info()¶
Get the conversion info for this forward_compatible hierarchy
- convert_met_to_semet()¶
- convert_multi_word_text_to_forward_compatible(text)¶
Use conversion info to convert words in a text string to forward-compatible equivalents :params text: text to convert :returns modified text
- convert_semet_to_met()¶
- de_deuterate()¶
Remove all D atoms and replace with H. Keep only H at hydrogen/deuterium sites. Changes hierarchy in place.
- deep_copy((root)arg1) root : ¶
- C++ signature :
iotbx::pdb::hierarchy::root deep_copy(iotbx::pdb::hierarchy::root {lvalue})
- distance_based_simple_two_way_bond_sets(fallback_expected_bond_length=1.4, fallback_search_max_distance=2.5)¶
- exchangeable_hd_selections()¶
- expand_to_p1(crystal_symmetry, exclude_self=False)¶
- extract_xray_structure(crystal_symmetry=None, min_distance_sym_equiv=None)¶
Generate the equivalent cctbx.xray.structure object. If the crystal symmetry is not provided, this will be placed in a P1 box. In practice it is usually best to keep the original xray structure object around, but this method is helpful in corner cases.
- find_model_index((root)arg1, (model)model[, (object)must_be_present=False]) int : ¶
- C++ signature :
long find_model_index(iotbx::pdb::hierarchy::root {lvalue},iotbx::pdb::hierarchy::model [,bool=False])
- first_chain_id()¶
Get first chain ID
- first_resseq_as_int(chain_id=None)¶
Return residue number of first residue in specified chain, as integer. If chain not specified, first residue in hierarchy.
- fits_in_pdb_format((root)arg1[, (object)use_hybrid36=True]) bool : ¶
- C++ signature :
bool fits_in_pdb_format(iotbx::pdb::hierarchy::root {lvalue} [,bool=True])
- flip_symmetric_amino_acids()¶
- format_correction_for_H(verbose=False)¶
- format_fasta(substitute_unknown='X', substitute_unknown_na='N', ignore_all_unknown=None, as_string=False)¶
uses format_fasta for all chains and returns catenation :param substitute_unknown: character to use for unrecognized 3-letter codes :param substitute_unknown_na: character to use for unrecognized na codes :param ignore_all_unknown: set substitute_unknown and substitute_unknown_na to ‘’ :param as_string: return string (default is to return list of lines)
- forward_compatible_hierarchy_as_standard(conversion_info=None)¶
- Convert a forward_compatible_hierarchy to a standard one.
Inverse of as_forward_compatible_hierarchy. Restores chain IDs and residue names using conversion_info object
- Params:
conversion_info: optional conversion_info object specifying conversion to be applied
:returns pdb_hierarchy with original (standard) chain ID and residue names
- get_atom_selection_cache((root)arg1, (object)arg2) None : ¶
- C++ signature :
void get_atom_selection_cache(iotbx::pdb::hierarchy::root,boost::python::api::object)
- get_auth_asym_id(chain, segid_as_auth_segid=False)¶
- get_auth_asym_id_iseq(iseq)¶
- get_auth_seq_id(rg)¶
- get_auth_seq_id_iseq(iseq)¶
- get_conformer_indices()¶
- get_label_alt_id_atom(atom)¶
- get_label_alt_id_iseq(iseq)¶
- get_label_asym_id(residue_group)¶
- get_label_asym_id_iseq(iseq)¶
- get_label_seq_id(atom_group)¶
- get_label_seq_id_iseq(iseq)¶
- get_overall_counts((root)arg1, (object)arg2) None : ¶
- C++ signature :
void get_overall_counts(iotbx::pdb::hierarchy::root,boost::python::api::object)
- get_peptide_c_alpha_selection()¶
Extract atom selection (flex.size_t) for protein C-alpha atoms.
- guess_chemical_elements(check_pseudo=False, convert_atom_names_to_uppercase=True, allow_incorrect_spacing=None)¶
Attempt to guess chemical elements for all atoms in hierarchy where this is not set. Normally used only just after reading in PDB-formatted files that do not have elements specified
- has_icodes()¶
- property info¶
- insert_model((root)arg1, (object)i, (model)model) None : ¶
- C++ signature :
void insert_model(iotbx::pdb::hierarchy::root {lvalue},long,iotbx::pdb::hierarchy::model {lvalue})
- is_ca_only()¶
Determine if hierarchy consists only from CA atoms. Upgrade options:
implement threshold for cases where several residues are present in full;
figure out how to deal with HETATM records of the same chain.
Ignore possible incorrect alignment of atom names.
- is_forward_compatible_hierarchy()¶
Determine if this is a forward_compatible hierarchy
- is_hierarchy_altloc_consistent(verbose=False)¶
- is_similar_hierarchy((root)arg1, (root)other) bool : ¶
- C++ signature :
bool is_similar_hierarchy(iotbx::pdb::hierarchy::root {lvalue},iotbx::pdb::hierarchy::root)
- last_resseq_as_int(chain_id=None)¶
Return residue number of last residue in specified chain, as integer. If chain not specified, last residue in hierarchy.
- memory_id((root)arg1) int : ¶
- C++ signature :
unsigned long memory_id(iotbx::pdb::hierarchy::root {lvalue})
- merge_atoms_at_end_to_residues()¶
Transfered from qrefine for merging single H/D atoms from the end of the PDB input to the correct residue object
- models((root)arg1) list : ¶
- C++ signature :
boost::python::list models(iotbx::pdb::hierarchy::root)
- models_size((root)arg1) int : ¶
- C++ signature :
unsigned int models_size(iotbx::pdb::hierarchy::root {lvalue})
- occupancy_counts()¶
- occupancy_groups_simple(common_residue_name_class_only=None, always_group_adjacent=True, ignore_hydrogens=True)¶
- only_atom()¶
- only_atom_group()¶
- only_chain()¶
- only_conformer()¶
- only_model()¶
- only_residue()¶
- only_residue_group()¶
- overall_counts(only_one_model=None)¶
Calculate basic statistics for contents of the PDB hierarchy, including number of residues of each type.
- Parameters:
only_one_model – return results for first model only
- Returns:
iotbx.pdb.hierarchy.overall_counts object
- pdb_or_mmcif_string_info(target_format=None, target_filename=None, data_manager=None, overwrite=True, segid_as_auth_segid=True, write_file=False, remark_section=None, **kw)¶
NOTE: Normally use instead either as_pdb_or_mmcif_string write_pdb_or_mmcif_file.
Method to allow shifting from general writing as pdb to writing as mmcif, with the change in two places (here and model.py) Use default of segid_as_auth_segid=True here (different than
as_mmcif_string())
- Parameters:
target_format – desired output format, pdb or mmcif
target_filename – desired output file name, to be modified to match the output format
data_manager – data_manager to write files
overwrite – parameter to set overwrite=True in data_manager if True
segid_as_auth_segid – use the segid in hierarchy as the auth_segid in mmcif output
write_file – Write the string to the target file
remark_section – if supplied and format is pdb, add this text
**kw –
any keywords suitable for as_pdb_string() and as_mmcif_string()
- :returns group_args object with attributes
pdb_string, file_name (the actual file name used) and is_mmcif
- pre_allocate_models((root)arg1, (object)number_of_additional_models) None : ¶
- C++ signature :
void pre_allocate_models(iotbx::pdb::hierarchy::root {lvalue},unsigned int)
- remove_alt_confs(always_keep_one_conformer, altloc_to_keep=None, keep_occupancy=False)¶
- remove_atoms(fraction)¶
- remove_hd(reset_i_seq=False)¶
Remove all hydrogen/deuterium atoms in-place. Returns the number of atoms deleted.
- remove_hetero()¶
- remove_incomplete_main_chain_protein(required_atom_names=['CA', 'N', 'C', 'O'])¶
- remove_model((root)arg1, (object)i) None : ¶
- C++ signature :
void remove_model(iotbx::pdb::hierarchy::root {lvalue},long)
remove_model( (root)arg1, (model)model) -> None :
- C++ signature :
void remove_model(iotbx::pdb::hierarchy::root {lvalue},iotbx::pdb::hierarchy::model {lvalue})
- remove_residue_groups_with_atoms_on_special_positions_selective(crystal_symmetry)¶
- remove_ter_or_break()¶
- rename_chain_id(old_id, new_id)¶
- reset_atom_i_seqs((root)arg1) int : ¶
- C++ signature :
unsigned long reset_atom_i_seqs(iotbx::pdb::hierarchy::root {lvalue})
- reset_i_seq_if_necessary()¶
- residue_groups()¶
Iterate over all residue groups (by model and then chain)
- round_occupancies_in_place(ndigits)¶
Round occupancies of those alternative conformations that cannot be rounded properly to the sum of 1 by standard round procedure in the output. The rest occupancies left intact.
- Parameters:
ndigits (int) – number of significant digits after the dot
- select((root)arg1, (bool)atom_selection[, (object)copy_atoms=False]) root : ¶
- C++ signature :
iotbx::pdb::hierarchy::root select(iotbx::pdb::hierarchy::root {lvalue},scitbx::af::const_ref<bool, scitbx::af::trivial_accessor> [,bool=False])
select( (root)arg1, (size_t)atom_selection [, (object)copy_atoms=False]) -> root :
- C++ signature :
iotbx::pdb::hierarchy::root select(iotbx::pdb::hierarchy::root {lvalue},scitbx::af::const_ref<unsigned long, scitbx::af::trivial_accessor> [,bool=False])
- set_atomic_charge(iselection, charge)¶
- shift_to_origin(crystal_symmetry)¶
- show(out=None, prefix='', level_id=None, level_id_exception=<class 'ValueError'>)¶
Display a summary of hierarchy contents.
- sort_atoms_in_place((root)arg1) None : ¶
- C++ signature :
void sort_atoms_in_place(iotbx::pdb::hierarchy::root {lvalue})
- sort_chains_by_id()¶
- transfer_chains_from_other(other)¶
- truncate_to_poly(atom_names_set={})¶
- truncate_to_poly_ala()¶
- truncate_to_poly_gly()¶
- write_mmcif_file(file_name, crystal_symmetry=None, data_block_name=None, segid_as_auth_segid=False, output_break_records=False)¶
- write_pdb_file(file_name, open_append=False, crystal_symmetry=None, cryst1_z=None, write_scale_records=True, append_end=False, interleaved_conf=0, atoms_reset_serial_first_value=None, atom_hetatm=True, sigatm=True, anisou=True, siguij=True, link_records=None)¶
- write_pdb_or_mmcif_file(target_filename, target_format=None, data_manager=None, overwrite=True, segid_as_auth_segid=True, remark_section=None, **kw)¶
- Shortcut for pdb_or_mmcif_string_info with write_file=True, returning
only the name of the file that is written. The file may be written in PDB or mmCIF format, with target_format used if feasible.
Method to allow shifting from general writing as pdb to writing as mmcif, with the change in two places (here and model.py) Use default of segid_as_auth_segid=True here (different than
as_mmcif_string())
- Parameters:
target_format – desired output format, pdb or mmcif
target_filename – desired output file name, to be modified to match the output format
data_manager – data_manager to write files
overwrite – parameter to set overwrite=True in data_manager if True
segid_as_auth_segid – use the segid in hierarchy as the auth_segid in mmcif output
remark_section – if supplied and format is pdb, add this text
**kw –
any keywords suitable for as_pdb_string() and as_mmcif_string()
:returns name of file that is written
- class iotbx.pdb.hierarchy.model¶
Bases:
instance
Class representing MODEL blocks in a PDB file (or equivalent mmCIF). There will always be at least one of these in a hierarchy root extracted from a PDB file even if no MODEL records are present.
Example
>>> hierarchy = iotbx.pdb.hierarchy.root() >>> model = iotbx.pdb.hierarchy.model(id="1") >>> hierarchy.append_model(model) >>> model = hierarchy.only_model()
- append_chain((model)arg1, (chain)chain) None : ¶
- C++ signature :
void append_chain(iotbx::pdb::hierarchy::model {lvalue},iotbx::pdb::hierarchy::chain {lvalue})
- atom_groups()¶
- atoms((model)arg1[, (object)interleaved_conf=0]) af_shared_atom : ¶
- C++ signature :
scitbx::af::shared<iotbx::pdb::hierarchy::atom> atoms(iotbx::pdb::hierarchy::model {lvalue} [,int=0])
- atoms_size((model)arg1) int : ¶
- C++ signature :
unsigned int atoms_size(iotbx::pdb::hierarchy::model {lvalue})
- chains((model)arg1) list : ¶
- C++ signature :
boost::python::list chains(iotbx::pdb::hierarchy::model)
- chains_size((model)arg1) int : ¶
- C++ signature :
unsigned int chains_size(iotbx::pdb::hierarchy::model {lvalue})
- detached_copy((model)arg1) model : ¶
- C++ signature :
iotbx::pdb::hierarchy::model detached_copy(iotbx::pdb::hierarchy::model {lvalue})
- find_chain_index((model)arg1, (chain)chain[, (object)must_be_present=False]) int : ¶
- C++ signature :
long find_chain_index(iotbx::pdb::hierarchy::model {lvalue},iotbx::pdb::hierarchy::chain [,bool=False])
- property id¶
- insert_chain((model)arg1, (object)i, (chain)chain) None : ¶
- C++ signature :
void insert_chain(iotbx::pdb::hierarchy::model {lvalue},long,iotbx::pdb::hierarchy::chain {lvalue})
- is_ca_only()¶
Determine if model consists only from CA atoms. Upgrade options:
implement threshold for cases where several residues are present in full;
figure out how to deal with HETATM records of the same chain.
Ignore possible incorrect alignment of atom names.
- is_identical_hierarchy((model)arg1, (model)other) bool : ¶
- C++ signature :
bool is_identical_hierarchy(iotbx::pdb::hierarchy::model {lvalue},iotbx::pdb::hierarchy::model)
- is_similar_hierarchy((model)arg1, (model)other) bool : ¶
- C++ signature :
bool is_similar_hierarchy(iotbx::pdb::hierarchy::model {lvalue},iotbx::pdb::hierarchy::model)
- memory_id((model)arg1) int : ¶
- C++ signature :
unsigned long memory_id(iotbx::pdb::hierarchy::model {lvalue})
- only_atom()¶
- only_atom_group()¶
- only_chain()¶
- only_conformer()¶
- only_residue()¶
- only_residue_group()¶
- parent((model)arg1[, (object)optional=True]) object : ¶
- C++ signature :
boost::python::api::object parent(iotbx::pdb::hierarchy::model [,bool=True])
- pre_allocate_chains((model)arg1, (object)number_of_additional_chains) None : ¶
- C++ signature :
void pre_allocate_chains(iotbx::pdb::hierarchy::model {lvalue},unsigned int)
- remove_chain((model)arg1, (object)i) None : ¶
- C++ signature :
void remove_chain(iotbx::pdb::hierarchy::model {lvalue},long)
remove_chain( (model)arg1, (chain)chain) -> None :
- C++ signature :
void remove_chain(iotbx::pdb::hierarchy::model {lvalue},iotbx::pdb::hierarchy::chain {lvalue})
- residue_groups()¶
- transfer_chains_from_other((model)arg1, (model)other) None : ¶
- C++ signature :
void transfer_chains_from_other(iotbx::pdb::hierarchy::model {lvalue},iotbx::pdb::hierarchy::model {lvalue})
- class iotbx.pdb.hierarchy.chain¶
Bases:
instance
Class representing a continuous chain of atoms, as defined by the combination of chain ID field and TER records (or the chain index in mmCIF format). Note that this does not necessarily correspond to a covalently linked entity, as it may be used to group various heteroatoms (including water), but chemically distinct protein or nucleic acid chains will typically be grouped into exactly one chain object apiece.
- append_residue_group((chain)arg1, (residue_group)residue_group) None : ¶
- C++ signature :
void append_residue_group(iotbx::pdb::hierarchy::chain {lvalue},iotbx::pdb::hierarchy::residue_group {lvalue})
- as_dict_of_resseq_as_int_residue_names()¶
- as_dict_of_resseq_residue_names(strip_resseq=True)¶
- as_list_of_residue_names()¶
- as_new_hierarchy()¶
- as_padded_sequence(missing_char='X', skip_insertions=False, pad=True, substitute_unknown='X', pad_at_start=True, ignore_hetatm=False)¶
Extract protein or nucleic acid sequence, taking residue numbering into account so that apparent gaps will be filled with substitute characters.
- as_sequence(substitute_unknown='X', substitute_unknown_na='N', ignore_all_unknown=None, as_string=False)¶
Naively extract single-character protein or nucleic acid sequence, without accounting for residue numbering.
- Parameters:
substitute_unknown – character to use for unrecognized 3-letter codes
substitute_unknown_na – character to use for unrecognized na codes
ignore_all_unknown – set substitute_unknown and substitute_unknown_na to ‘’
as_string – return string (default is to return list)
- atom_groups()¶
- atoms((chain)arg1[, (object)interleaved_conf=0]) af_shared_atom : ¶
- C++ signature :
scitbx::af::shared<iotbx::pdb::hierarchy::atom> atoms(iotbx::pdb::hierarchy::chain {lvalue} [,int=0])
- atoms_size((chain)arg1) int : ¶
- C++ signature :
unsigned int atoms_size(iotbx::pdb::hierarchy::chain {lvalue})
- conformers((chain)arg1) object : ¶
- C++ signature :
boost::python::api::object conformers(iotbx::pdb::hierarchy::chain)
- detached_copy((chain)arg1) chain : ¶
- C++ signature :
iotbx::pdb::hierarchy::chain detached_copy(iotbx::pdb::hierarchy::chain {lvalue})
- find_pure_altloc_ranges((chain)arg1[, (str)common_residue_name_class_only=None]) tiny_size_t_2 : ¶
- C++ signature :
scitbx::af::shared<scitbx::af::tiny<unsigned long, 2ul> > find_pure_altloc_ranges(iotbx::pdb::hierarchy::chain {lvalue} [,char const*=None])
- find_residue_group_index((chain)arg1, (residue_group)residue_group[, (object)must_be_present=False]) int : ¶
- C++ signature :
long find_residue_group_index(iotbx::pdb::hierarchy::chain {lvalue},iotbx::pdb::hierarchy::residue_group [,bool=False])
- format_fasta(max_line_length=79, substitute_unknown='X', substitute_unknown_na='N', ignore_all_unknown=None, as_string=False)¶
Format this chain as Fasta :param max_line_length: length of lines in formatted output :param substitute_unknown: character to use for unrecognized 3-letter codes :param substitute_unknown_na: character to use for unrecognized na codes :param ignore_all_unknown: set substitute_unknown and substitute_unknown_na to ‘’ :param as_string: return string (default is to return list of lines)
- get_residue_ids(skip_insertions=False, pad=True, pad_at_start=True, ignore_hetatm=False)¶
- get_residue_names_and_classes()¶
Extract the residue names and counts of each residue type (protein, nucleic acid, etc) within the chain.
- Returns:
a tuple containing a list of residue names, and a dictionary of residue type frequencies.
- get_residue_names_padded(skip_insertions=False, pad=True, pad_at_start=True, ignore_hetatm=False)¶
- property id¶
- insert_residue_group((chain)arg1, (object)i, (residue_group)residue_group) None : ¶
- C++ signature :
void insert_residue_group(iotbx::pdb::hierarchy::chain {lvalue},long,iotbx::pdb::hierarchy::residue_group {lvalue})
- is_ca_only()¶
Determine if chain consists only from CA atoms. Upgrade options:
implement threshold for cases where several residues are present in full;
figure out how to deal with HETATM records of the same chain.
Ignore possible incorrect alignment of atom names.
- is_identical_hierarchy((chain)arg1, (chain)other) bool : ¶
- C++ signature :
bool is_identical_hierarchy(iotbx::pdb::hierarchy::chain {lvalue},iotbx::pdb::hierarchy::chain)
- is_na(min_content=0.8, ignore_water=True)¶
Determine whether the chain represents a nucleic acid polymer, based on the frequency of base names. Very slow due to usage of residue_name_plus_atom_names_interpreter in get_residue_names_and_classes (majority of the processing is unnecessary)
- is_protein(min_content=0.8, ignore_water=True)¶
Determine whether the chain represents an amino acid polymer, based on the frequency of residue names. Very slow due to usage of residue_name_plus_atom_names_interpreter in get_residue_names_and_classes (majority of the processing is unnecessary)
- is_similar_hierarchy((chain)arg1, (chain)other) bool : ¶
- C++ signature :
bool is_similar_hierarchy(iotbx::pdb::hierarchy::chain {lvalue},iotbx::pdb::hierarchy::chain)
- memory_id((chain)arg1) int : ¶
- C++ signature :
unsigned long memory_id(iotbx::pdb::hierarchy::chain {lvalue})
- merge_disconnected_residue_groups_with_pure_altloc((chain)arg1) size_t : ¶
- C++ signature :
scitbx::af::shared<unsigned long> merge_disconnected_residue_groups_with_pure_altloc(iotbx::pdb::hierarchy::chain {lvalue})
- merge_residue_groups((chain)arg1, (residue_group)primary, (residue_group)secondary) None : ¶
- C++ signature :
void merge_residue_groups(iotbx::pdb::hierarchy::chain {lvalue},iotbx::pdb::hierarchy::residue_group {lvalue},iotbx::pdb::hierarchy::residue_group {lvalue})
- occupancy_groups_simple(common_residue_name_class_only=None, always_group_adjacent=True)¶
- only_atom()¶
- only_atom_group()¶
- only_conformer()¶
- only_residue()¶
- only_residue_group()¶
- parent((chain)arg1[, (object)optional=True]) object : ¶
- C++ signature :
boost::python::api::object parent(iotbx::pdb::hierarchy::chain [,bool=True])
- pre_allocate_residue_groups((chain)arg1, (object)number_of_additional_residue_groups) None : ¶
- C++ signature :
void pre_allocate_residue_groups(iotbx::pdb::hierarchy::chain {lvalue},unsigned int)
- remove_residue_group((chain)arg1, (object)i) None : ¶
- C++ signature :
void remove_residue_group(iotbx::pdb::hierarchy::chain {lvalue},long)
remove_residue_group( (chain)arg1, (residue_group)residue_group) -> None :
- C++ signature :
void remove_residue_group(iotbx::pdb::hierarchy::chain {lvalue},iotbx::pdb::hierarchy::residue_group {lvalue})
- residue_groups((chain)arg1) list : ¶
- C++ signature :
boost::python::list residue_groups(iotbx::pdb::hierarchy::chain)
- residue_groups_size((chain)arg1) int : ¶
- C++ signature :
unsigned int residue_groups_size(iotbx::pdb::hierarchy::chain {lvalue})
- residues()¶
- class iotbx.pdb.hierarchy.residue_group¶
Bases:
instance
- append_atom_group((residue_group)arg1, (atom_group)atom_group) None : ¶
- C++ signature :
void append_atom_group(iotbx::pdb::hierarchy::residue_group {lvalue},iotbx::pdb::hierarchy::atom_group {lvalue})
- atom_groups((residue_group)arg1) list : ¶
- C++ signature :
boost::python::list atom_groups(iotbx::pdb::hierarchy::residue_group)
- atom_groups_size((residue_group)arg1) int : ¶
- C++ signature :
unsigned int atom_groups_size(iotbx::pdb::hierarchy::residue_group {lvalue})
- atoms((residue_group)arg1[, (object)interleaved_conf=0]) af_shared_atom : ¶
- C++ signature :
scitbx::af::shared<iotbx::pdb::hierarchy::atom> atoms(iotbx::pdb::hierarchy::residue_group {lvalue} [,int=0])
- atoms_size((residue_group)arg1) int : ¶
- C++ signature :
unsigned int atoms_size(iotbx::pdb::hierarchy::residue_group {lvalue})
- conformers((residue_group)arg1) object : ¶
- C++ signature :
boost::python::api::object conformers(iotbx::pdb::hierarchy::residue_group)
- detached_copy((residue_group)arg1) residue_group : ¶
- C++ signature :
iotbx::pdb::hierarchy::residue_group detached_copy(iotbx::pdb::hierarchy::residue_group {lvalue})
- edit_blank_altloc((residue_group)arg1) tuple : ¶
- C++ signature :
scitbx::af::tiny<unsigned int, 2ul> edit_blank_altloc(iotbx::pdb::hierarchy::residue_group {lvalue})
- find_atom_group_index((residue_group)arg1, (atom_group)atom_group[, (object)must_be_present=False]) int : ¶
- C++ signature :
long find_atom_group_index(iotbx::pdb::hierarchy::residue_group {lvalue},iotbx::pdb::hierarchy::atom_group [,bool=False])
- have_conformers((residue_group)arg1) bool : ¶
- C++ signature :
bool have_conformers(iotbx::pdb::hierarchy::residue_group {lvalue})
- property icode¶
- id_str()¶
- insert_atom_group((residue_group)arg1, (object)i, (atom_group)atom_group) None : ¶
- C++ signature :
void insert_atom_group(iotbx::pdb::hierarchy::residue_group {lvalue},long,iotbx::pdb::hierarchy::atom_group {lvalue})
- is_identical_hierarchy((residue_group)arg1, (residue_group)other) bool : ¶
- C++ signature :
bool is_identical_hierarchy(iotbx::pdb::hierarchy::residue_group {lvalue},iotbx::pdb::hierarchy::residue_group)
- is_similar_hierarchy((residue_group)arg1, (residue_group)other) bool : ¶
- C++ signature :
bool is_similar_hierarchy(iotbx::pdb::hierarchy::residue_group {lvalue},iotbx::pdb::hierarchy::residue_group)
- property link_to_previous¶
- memory_id((residue_group)arg1) int : ¶
- C++ signature :
unsigned long memory_id(iotbx::pdb::hierarchy::residue_group {lvalue})
- merge_atom_groups((residue_group)arg1, (atom_group)primary, (atom_group)secondary) None : ¶
- C++ signature :
void merge_atom_groups(iotbx::pdb::hierarchy::residue_group {lvalue},iotbx::pdb::hierarchy::atom_group {lvalue},iotbx::pdb::hierarchy::atom_group {lvalue})
- move_blank_altloc_atom_groups_to_front((residue_group)arg1) int : ¶
- C++ signature :
unsigned int move_blank_altloc_atom_groups_to_front(iotbx::pdb::hierarchy::residue_group {lvalue})
- only_atom()¶
- only_atom_group()¶
- parent((residue_group)arg1[, (object)optional=True]) object : ¶
- C++ signature :
boost::python::api::object parent(iotbx::pdb::hierarchy::residue_group [,bool=True])
- pre_allocate_atom_groups((residue_group)arg1, (object)number_of_additional_atom_groups) None : ¶
- C++ signature :
void pre_allocate_atom_groups(iotbx::pdb::hierarchy::residue_group {lvalue},unsigned int)
- remove_atom_group((residue_group)arg1, (object)i) None : ¶
- C++ signature :
void remove_atom_group(iotbx::pdb::hierarchy::residue_group {lvalue},long)
remove_atom_group( (residue_group)arg1, (atom_group)atom_group) -> None :
- C++ signature :
void remove_atom_group(iotbx::pdb::hierarchy::residue_group {lvalue},iotbx::pdb::hierarchy::atom_group {lvalue})
- resid((residue_group)arg1) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > resid(iotbx::pdb::hierarchy::residue_group {lvalue})
- property resseq¶
- resseq_as_int((residue_group)arg1) int : ¶
- C++ signature :
int resseq_as_int(iotbx::pdb::hierarchy::residue_group {lvalue})
- unique_resnames((residue_group)arg1) std_string : ¶
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > unique_resnames(iotbx::pdb::hierarchy::residue_group {lvalue})
- class iotbx.pdb.hierarchy.atom_group¶
Bases:
instance
- property altloc¶
- append_atom((atom_group)arg1, (atom)atom) None : ¶
- C++ signature :
void append_atom(iotbx::pdb::hierarchy::atom_group {lvalue},iotbx::pdb::hierarchy::atom {lvalue})
- append_atom_with_other_parent((atom_group)arg1, (atom)atom) None : ¶
- C++ signature :
void append_atom_with_other_parent(iotbx::pdb::hierarchy::atom_group {lvalue},iotbx::pdb::hierarchy::atom)
- atoms((atom_group)arg1) af_shared_atom : ¶
- C++ signature :
scitbx::af::shared<iotbx::pdb::hierarchy::atom> atoms(iotbx::pdb::hierarchy::atom_group)
- atoms_size((atom_group)arg1) int : ¶
- C++ signature :
unsigned int atoms_size(iotbx::pdb::hierarchy::atom_group {lvalue})
- confid((atom_group)arg1) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > confid(iotbx::pdb::hierarchy::atom_group {lvalue})
- detached_copy((atom_group)arg1) atom_group : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom_group detached_copy(iotbx::pdb::hierarchy::atom_group {lvalue})
- find_atom_index((atom_group)arg1, (atom)atom[, (object)must_be_present=False]) int : ¶
- C++ signature :
long find_atom_index(iotbx::pdb::hierarchy::atom_group {lvalue},iotbx::pdb::hierarchy::atom [,bool=False])
- get_atom((atom_group)arg1, (str)name) atom : ¶
- C++ signature :
boost::optional<iotbx::pdb::hierarchy::atom> get_atom(iotbx::pdb::hierarchy::atom_group {lvalue},char const*)
- id_str(suppress_segid=None)¶
- insert_atom((atom_group)arg1, (object)i, (atom)atom) None : ¶
- C++ signature :
void insert_atom(iotbx::pdb::hierarchy::atom_group {lvalue},long,iotbx::pdb::hierarchy::atom {lvalue})
- memory_id((atom_group)arg1) int : ¶
- C++ signature :
unsigned long memory_id(iotbx::pdb::hierarchy::atom_group {lvalue})
- occupancy(raise_error_if_non_uniform=False)¶
Calculate the mean occupancy for atoms in this group, with option of raising ValueError if they differ.
- only_atom()¶
- parent((atom_group)arg1[, (object)optional=True]) object : ¶
- C++ signature :
boost::python::api::object parent(iotbx::pdb::hierarchy::atom_group [,bool=True])
- pre_allocate_atoms((atom_group)arg1, (object)number_of_additional_atoms) None : ¶
- C++ signature :
void pre_allocate_atoms(iotbx::pdb::hierarchy::atom_group {lvalue},unsigned int)
- remove_atom((atom_group)arg1, (object)i) None : ¶
- C++ signature :
void remove_atom(iotbx::pdb::hierarchy::atom_group {lvalue},long)
remove_atom( (atom_group)arg1, (atom)atom) -> None :
- C++ signature :
void remove_atom(iotbx::pdb::hierarchy::atom_group {lvalue},iotbx::pdb::hierarchy::atom {lvalue})
- property resname¶
- sort_atoms_in_place((atom_group)arg1) None : ¶
- C++ signature :
void sort_atoms_in_place(iotbx::pdb::hierarchy::atom_group {lvalue})
- class iotbx.pdb.hierarchy.atom¶
Bases:
instance
The basic unit of the PDB hierarchy (or the PDB input object in general), representing a single point scatterer corresponding to an ATOM or HETATM record in PDB format (plus associated ANISOU or related records if present). Note that this does not directly store attributes of higher-level entities whose identity is also recorded in ATOM records, such as the chain ID or residue name. These may be retrieved either by walking up the hierarchy starting with atom.parent(), or by calling atom.fetch_labels().
- angle((atom)arg1, (object)atom_1_xyz, (object)atom_3_xyz[, (object)deg=False]) object : ¶
- C++ signature :
boost::optional<double> angle(iotbx::pdb::hierarchy::atom {lvalue},scitbx::vec3<double>,scitbx::vec3<double> [,bool=False])
angle( (atom)arg1, (atom)atom_1, (atom)atom_3 [, (object)deg=False]) -> object :
- C++ signature :
boost::optional<double> angle(iotbx::pdb::hierarchy::atom {lvalue},iotbx::pdb::hierarchy::atom,iotbx::pdb::hierarchy::atom [,bool=False])
- property b¶
- chain()¶
Convenience method for fetching the chain object associated with this atom (or None of not defined).
- property charge¶
- charge_as_int()¶
Extract the atomic charge from the (string) charge field.
- Returns:
Python int, defaulting to zero
- charge_tidy((atom)arg1[, (object)strip=False]) object : ¶
- C++ signature :
boost::optional<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > charge_tidy(iotbx::pdb::hierarchy::atom {lvalue} [,bool=False])
- static data_offsets() dict : ¶
- C++ signature :
boost::python::dict data_offsets()
- detached_copy((atom)arg1) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom detached_copy(iotbx::pdb::hierarchy::atom {lvalue})
- determine_chemical_element_simple((atom)arg1) object : ¶
- C++ signature :
boost::optional<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > determine_chemical_element_simple(iotbx::pdb::hierarchy::atom {lvalue})
- distance((atom)arg1, (object)other_xyz) float : ¶
- C++ signature :
double distance(iotbx::pdb::hierarchy::atom {lvalue},scitbx::vec3<double>)
distance( (atom)arg1, (atom)other) -> float :
- C++ signature :
double distance(iotbx::pdb::hierarchy::atom {lvalue},iotbx::pdb::hierarchy::atom)
- property element¶
- element_is_hydrogen((atom)arg1) bool : ¶
- C++ signature :
bool element_is_hydrogen(iotbx::pdb::hierarchy::atom {lvalue})
- element_is_ion((atom)arg1) bool : ¶
- C++ signature :
bool element_is_ion(iotbx::pdb::hierarchy::atom {lvalue})
- element_is_negative_ion((atom)arg1) bool : ¶
- C++ signature :
bool element_is_negative_ion(iotbx::pdb::hierarchy::atom {lvalue})
- element_is_positive_ion((atom)arg1) bool : ¶
- C++ signature :
bool element_is_positive_ion(iotbx::pdb::hierarchy::atom {lvalue})
- property fdp¶
- fetch_labels((atom)arg1) atom_with_labels : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom_with_labels fetch_labels(iotbx::pdb::hierarchy::atom {lvalue})
- format_anisou_record((atom)arg1) object : ¶
- C++ signature :
boost::python::api::object format_anisou_record(iotbx::pdb::hierarchy::atom)
- format_atom_record((atom)self[, (str)replace_floats_with=None]) str : ¶
- C++ signature :
boost::python::str format_atom_record(iotbx::pdb::hierarchy::atom [,char const*=None])
- format_atom_record_group((atom)self[, (object)atom_hetatm=True[, (object)sigatm=True[, (object)anisou=True[, (object)siguij=True]]]]) object : ¶
- C++ signature :
boost::python::api::object format_atom_record_group(iotbx::pdb::hierarchy::atom [,bool=True [,bool=True [,bool=True [,bool=True]]]])
- format_sigatm_record((atom)arg1) object : ¶
- C++ signature :
boost::python::api::object format_sigatm_record(iotbx::pdb::hierarchy::atom)
- format_siguij_record((atom)arg1) object : ¶
- C++ signature :
boost::python::api::object format_siguij_record(iotbx::pdb::hierarchy::atom)
- property fp¶
- static has_siguij() bool : ¶
- C++ signature :
bool has_siguij()
- property hetero¶
- property i_seq¶
- id_str((atom)arg1[, (object)pdbres=False[, (object)suppress_segid=False]]) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > id_str(iotbx::pdb::hierarchy::atom {lvalue} [,bool=False [,bool=False]])
- is_in_same_conformer_as(other)¶
Indicate whether two atoms are part of the same conformer and thus are capable of interacting directly, as defined by the parent atom_group and model object(s).
- memory_id((atom)arg1) int : ¶
- C++ signature :
unsigned long memory_id(iotbx::pdb::hierarchy::atom {lvalue})
- property name¶
- property occ¶
- parent((atom)arg1[, (object)optional=True]) object : ¶
- C++ signature :
boost::python::api::object parent(iotbx::pdb::hierarchy::atom [,bool=True])
- pdb_element_charge_columns((atom)arg1) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > pdb_element_charge_columns(iotbx::pdb::hierarchy::atom {lvalue})
- pdb_label_columns((atom)arg1) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > pdb_label_columns(iotbx::pdb::hierarchy::atom {lvalue})
- quote((atom)arg1[, (object)full=False]) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > quote(iotbx::pdb::hierarchy::atom {lvalue} [,bool=False])
- property segid¶
- property serial¶
- serial_as_int((atom)arg1) int : ¶
- C++ signature :
int serial_as_int(iotbx::pdb::hierarchy::atom {lvalue})
- set_b((atom)arg1, (object)new_b) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom set_b(iotbx::pdb::hierarchy::atom,double)
- set_charge((atom)arg1, (str)new_charge) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom {lvalue} set_charge(iotbx::pdb::hierarchy::atom {lvalue},char const*)
- set_chemical_element_simple_if_necessary((atom)arg1[, (object)tidy_existing=True]) bool : ¶
- C++ signature :
bool set_chemical_element_simple_if_necessary(iotbx::pdb::hierarchy::atom {lvalue} [,bool=True])
- set_element((atom)arg1, (str)new_element) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom {lvalue} set_element(iotbx::pdb::hierarchy::atom {lvalue},char const*)
- set_element_and_charge_from_scattering_type_if_necessary(scattering_type)¶
- set_fdp((atom)arg1, (object)new_fdp) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom set_fdp(iotbx::pdb::hierarchy::atom,double)
- set_fp((atom)arg1, (object)new_fp) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom set_fp(iotbx::pdb::hierarchy::atom,double)
- set_hetero((atom)arg1, (object)new_hetero) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom set_hetero(iotbx::pdb::hierarchy::atom,bool)
- set_name((atom)arg1, (str)new_name) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom {lvalue} set_name(iotbx::pdb::hierarchy::atom {lvalue},char const*)
- set_occ((atom)arg1, (object)new_occ) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom set_occ(iotbx::pdb::hierarchy::atom,double)
- set_segid((atom)arg1, (str)new_segid) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom {lvalue} set_segid(iotbx::pdb::hierarchy::atom {lvalue},char const*)
- set_serial((atom)arg1, (object)new_serial) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom {lvalue} set_serial(iotbx::pdb::hierarchy::atom {lvalue},boost::python::api::object)
- set_sigb((atom)arg1, (object)new_sigb) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom set_sigb(iotbx::pdb::hierarchy::atom,double)
- set_sigocc((atom)arg1, (object)new_sigocc) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom set_sigocc(iotbx::pdb::hierarchy::atom,double)
- set_sigxyz((atom)arg1, (object)new_sigxyz) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom set_sigxyz(iotbx::pdb::hierarchy::atom,scitbx::vec3<double>)
- set_uij((atom)arg1, (object)new_uij) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom set_uij(iotbx::pdb::hierarchy::atom,scitbx::sym_mat3<double>)
- set_xyz((atom)arg1, (object)new_xyz) atom : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom set_xyz(iotbx::pdb::hierarchy::atom,scitbx::vec3<double>)
- property sigb¶
- property sigocc¶
- siguij_erase((atom)arg1) None : ¶
- C++ signature :
void siguij_erase(iotbx::pdb::hierarchy::atom {lvalue})
- siguij_is_defined((atom)arg1) bool : ¶
- C++ signature :
bool siguij_is_defined(iotbx::pdb::hierarchy::atom {lvalue})
- property sigxyz¶
- static sizeof_data() int : ¶
- C++ signature :
unsigned long sizeof_data()
- property tmp¶
- property uij¶
- uij_erase((atom)arg1) None : ¶
- C++ signature :
void uij_erase(iotbx::pdb::hierarchy::atom {lvalue})
- uij_is_defined((atom)arg1) bool : ¶
- C++ signature :
bool uij_is_defined(iotbx::pdb::hierarchy::atom {lvalue})
- property xyz¶
- class iotbx.pdb.hierarchy.atom_with_labels¶
Bases:
atom
Stand-in for atom object, which explicitly records the attributes normally reserved for parent classes such as residue name, chain ID, etc.
- property altloc¶
- property chain_id¶
- detached_copy((atom_with_labels)arg1) atom_with_labels : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom_with_labels detached_copy(iotbx::pdb::hierarchy::atom_with_labels {lvalue})
- fetch_labels((atom)arg1) atom_with_labels : ¶
- C++ signature :
iotbx::pdb::hierarchy::atom_with_labels fetch_labels(iotbx::pdb::hierarchy::atom {lvalue})
- format_anisou_record((atom_with_labels)arg1) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > format_anisou_record(iotbx::pdb::hierarchy::atom_with_labels {lvalue})
- format_atom_record((atom_with_labels)arg1[, (str)replace_floats_with=None]) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > format_atom_record(iotbx::pdb::hierarchy::atom_with_labels {lvalue} [,char const*=None])
- format_atom_record_group((atom_with_labels)arg1[, (object)atom_hetatm=True[, (object)sigatm=True[, (object)anisou=True[, (object)siguij=True]]]]) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > format_atom_record_group(iotbx::pdb::hierarchy::atom_with_labels {lvalue} [,bool=True [,bool=True [,bool=True [,bool=True]]]])
- format_sigatm_record((atom_with_labels)arg1) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > format_sigatm_record(iotbx::pdb::hierarchy::atom_with_labels {lvalue})
- format_siguij_record((atom_with_labels)arg1) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > format_siguij_record(iotbx::pdb::hierarchy::atom_with_labels {lvalue})
- property icode¶
- id_str((atom_with_labels)arg1[, (object)pdbres=False[, (object)suppress_segid=False]]) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > id_str(iotbx::pdb::hierarchy::atom_with_labels {lvalue} [,bool=False [,bool=False]])
- property is_first_after_break¶
- property is_first_in_chain¶
- property model_id¶
- quote((atom_with_labels)arg1[, (object)full=False]) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > quote(iotbx::pdb::hierarchy::atom_with_labels {lvalue} [,bool=False])
- resid((atom_with_labels)arg1) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > resid(iotbx::pdb::hierarchy::atom_with_labels {lvalue})
- property resname¶
- property resseq¶
- resseq_as_int((atom_with_labels)arg1) int : ¶
- C++ signature :
int resseq_as_int(iotbx::pdb::hierarchy::atom_with_labels {lvalue})
- serial_as_int((atom_with_labels)arg1) int : ¶
- C++ signature :
int serial_as_int(iotbx::pdb::hierarchy::atom_with_labels {lvalue})
- class iotbx.pdb.hierarchy.conformer¶
Bases:
instance
Alternate view into a chain object, grouping sequential residues with equivalent altlocs. As a general rule it is preferrable to iterate over chain.residue_groups() instead.
- property altloc¶
- as_padded_sequence(missing_char='X', skip_insertions=False, pad=True, substitute_unknown='X', pad_at_start=True)¶
- as_sec_str_sequence(helix_sele, sheet_sele, missing_char='X', pad=True, pad_at_start=True)¶
- as_sequence(substitute_unknown='X')¶
- atoms((conformer)arg1) af_shared_atom : ¶
- C++ signature :
scitbx::af::shared<iotbx::pdb::hierarchy::atom> atoms(iotbx::pdb::hierarchy::conformer {lvalue})
- atoms_size((conformer)arg1) int : ¶
- C++ signature :
unsigned int atoms_size(iotbx::pdb::hierarchy::conformer {lvalue})
- format_fasta(max_line_length=79)¶
- get_residue_ids(skip_insertions=False, pad=True, pad_at_start=True)¶
- get_residue_names_and_classes()¶
- get_residue_names_padded(skip_insertions=False, pad=True, pad_at_start=True)¶
- is_na(min_content=0.8)¶
- is_protein(min_content=0.8)¶
- memory_id((conformer)arg1) int : ¶
- C++ signature :
unsigned long memory_id(iotbx::pdb::hierarchy::conformer {lvalue})
- only_atom()¶
- only_residue()¶
- parent((conformer)arg1[, (object)optional=True]) object : ¶
- C++ signature :
boost::python::api::object parent(iotbx::pdb::hierarchy::conformer [,bool=True])
- residues((conformer)arg1) list : ¶
- C++ signature :
boost::python::list residues(iotbx::pdb::hierarchy::conformer)
- residues_size((conformer)arg1) int : ¶
- C++ signature :
unsigned int residues_size(iotbx::pdb::hierarchy::conformer {lvalue})
- class iotbx.pdb.hierarchy.residue¶
Bases:
instance
- atoms((residue)arg1) af_shared_atom : ¶
- C++ signature :
scitbx::af::shared<iotbx::pdb::hierarchy::atom> atoms(iotbx::pdb::hierarchy::residue)
- atoms_size((residue)arg1) int : ¶
- C++ signature :
unsigned int atoms_size(iotbx::pdb::hierarchy::residue {lvalue})
- find_atom_by((residue)arg1, (str)name) atom : ¶
- C++ signature :
boost::optional<iotbx::pdb::hierarchy::atom> find_atom_by(iotbx::pdb::hierarchy::residue {lvalue},char const*)
- property icode¶
- id_str((residue)arg1[, (object)suppress_segid=0]) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > id_str(iotbx::pdb::hierarchy::residue {lvalue} [,int=0])
- property is_pure_main_conf¶
- property link_to_previous¶
- memory_id((residue)arg1) int : ¶
- C++ signature :
unsigned long memory_id(iotbx::pdb::hierarchy::residue {lvalue})
- only_atom()¶
- parent((residue)arg1[, (object)optional=True]) object : ¶
- C++ signature :
boost::python::api::object parent(iotbx::pdb::hierarchy::residue [,bool=True])
- resid((residue)arg1) str : ¶
- C++ signature :
std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > resid(iotbx::pdb::hierarchy::residue {lvalue})
- residue_name_plus_atom_names_interpreter(translate_cns_dna_rna_residue_names=None, return_mon_lib_dna_name=False)¶
- property resname¶
- property resseq¶
- resseq_as_int((residue)arg1) int : ¶
- C++ signature :
int resseq_as_int(iotbx::pdb::hierarchy::residue {lvalue})
- root((residue)arg1) object : ¶
- C++ signature :
boost::python::api::object root(iotbx::pdb::hierarchy::residue)
- standalone_copy()¶
Bases:
instance
- C++ signature :
void append(scitbx::af::shared<iotbx::pdb::hierarchy::atom> {lvalue},iotbx::pdb::hierarchy::atom)
- C++ signature :
boost::python::dict build_dict(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor> [,bool=False [,bool=False [,bool=False [,bool=True]]]])
- C++ signature :
void clear(scitbx::af::shared<iotbx::pdb::hierarchy::atom> {lvalue})
- C++ signature :
scitbx::af::shared<iotbx::pdb::hierarchy::atom> deep_copy(scitbx::af::shared<iotbx::pdb::hierarchy::atom> {lvalue})
- C++ signature :
void extend(scitbx::af::shared<iotbx::pdb::hierarchy::atom> {lvalue},scitbx::af::shared<iotbx::pdb::hierarchy::atom>)
- C++ signature :
scitbx::af::shared<double> extract_b(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > extract_element(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor> [,bool=False])
- C++ signature :
scitbx::af::shared<double> extract_fdp(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<double> extract_fp(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<unsigned long> extract_hetero(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<unsigned long> extract_i_seq(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > extract_name(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<double> extract_occ(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > extract_segid(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > > extract_serial(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<double> extract_sigb(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<double> extract_sigocc(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<scitbx::vec3<double> > extract_sigxyz(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<unsigned long> extract_tmp_as_size_t(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<scitbx::sym_mat3<double> > extract_uij(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<scitbx::vec3<double> > extract_xyz(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
void insert(scitbx::af::shared<iotbx::pdb::hierarchy::atom> {lvalue},long,iotbx::pdb::hierarchy::atom)
- C++ signature :
void reserve(scitbx::af::shared<iotbx::pdb::hierarchy::atom> {lvalue},unsigned long)
- C++ signature :
void reset_i_seq(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
void reset_serial(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor> [,int=1])
- C++ signature :
std::auto_ptr<iotbx::pdb::hierarchy::atoms::atom_tmp_sentinel> reset_tmp(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor> [,int=0 [,int=1]])
- C++ signature :
std::auto_ptr<iotbx::pdb::hierarchy::atoms::atom_tmp_sentinel> reset_tmp_for_occupancy_groups_simple(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::shared<iotbx::pdb::hierarchy::atom> select(scitbx::af::shared<iotbx::pdb::hierarchy::atom>,scitbx::af::const_ref<bool, scitbx::af::trivial_accessor>)
select( (af_shared_atom)self, (object)indices [, (object)reverse=False]) -> af_shared_atom :
- C++ signature :
scitbx::af::shared<iotbx::pdb::hierarchy::atom> select(scitbx::af::shared<iotbx::pdb::hierarchy::atom>,scitbx::af::const_ref<unsigned int, scitbx::af::trivial_accessor> [,bool=False])
select( (af_shared_atom)self, (size_t)indices [, (object)reverse=False]) -> af_shared_atom :
- C++ signature :
scitbx::af::shared<iotbx::pdb::hierarchy::atom> select(scitbx::af::shared<iotbx::pdb::hierarchy::atom>,scitbx::af::const_ref<unsigned long, scitbx::af::trivial_accessor> [,bool=False])
- C++ signature :
void set_adps_from_scatterers(scitbx::af::const_ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>,scitbx::af::const_ref<cctbx::xray::scatterer<double, std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> >, std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> > >, scitbx::af::trivial_accessor>,cctbx::uctbx::unit_cell)
- C++ signature :
scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor> set_b(scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>,scitbx::af::const_ref<double, scitbx::af::trivial_accessor>)
- C++ signature :
unsigned long set_chemical_element_simple_if_necessary(scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor> [,bool=True])
- C++ signature :
scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor> set_fdp(scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>,scitbx::af::const_ref<double, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor> set_fp(scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>,scitbx::af::const_ref<double, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor> set_occ(scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>,scitbx::af::const_ref<double, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor> set_sigb(scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>,scitbx::af::const_ref<double, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor> set_sigocc(scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>,scitbx::af::const_ref<double, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor> set_sigxyz(scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>,scitbx::af::const_ref<scitbx::vec3<double>, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor> set_uij(scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>,scitbx::af::const_ref<scitbx::sym_mat3<double>, scitbx::af::trivial_accessor>)
- C++ signature :
scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor> set_xyz(scitbx::af::ref<iotbx::pdb::hierarchy::atom, scitbx::af::trivial_accessor>,scitbx::af::const_ref<scitbx::vec3<double>, scitbx::af::trivial_accessor>)
- C++ signature :
unsigned long size(scitbx::af::shared<iotbx::pdb::hierarchy::atom> {lvalue})
Accessory classes and functions¶
- class iotbx.pdb.hierarchy.overall_counts¶
Bases:
object
- as_str(prefix='', residue_groups_max_show=10, duplicate_atom_labels_max_show=10)¶
- errors()¶
- errors_and_warnings()¶
- get_n_residues_of_classes(classes)¶
- raise_chains_with_mix_of_proper_and_improper_alt_conf_if_necessary()¶
- raise_duplicate_atom_labels_if_necessary(max_show=10)¶
- raise_improper_alt_conf_if_necessary()¶
- raise_residue_groups_with_multiple_resnames_using_same_altloc_if_necessary(max_show=10)¶
- show(out=None, prefix='', flag_errors=True, flag_warnings=True, residue_groups_max_show=10, duplicate_atom_labels_max_show=10)¶
- show_chains_with_mix_of_proper_and_improper_alt_conf(out=None, prefix='')¶
- show_consecutive_residue_groups_with_same_resid(out=None, prefix='', max_show=10)¶
- show_duplicate_atom_labels(out=None, prefix='', max_show=10)¶
- show_improper_alt_conf(out=None, prefix='')¶
- show_residue_groups_with_multiple_resnames_using_same_altloc(out=None, prefix='', max_show=10)¶
- warnings()¶