Libcootapi’s Documentation
Libcootapi is a new (as of 2022) interface to the functionality of Coot. It is designed to be (is intended to be) a clean, consistent and easy-to-use interface primarily targeting WebAssembly/JavaScript and secondarily, Python.
The previous API for Coot still exists and is available for use with Python and
Guile/Scheme (it provides the coot
module). That interface is considerably more
extensive than this one (consisting of several thousand API functions), but is
significantly less easy to use as it embeds OpenGL and GTK libraries (and all of their
dependencies, of course).
The Molecules Container
Coot’s molecules are referred to by index and should only be accessed via member
functions of the molecules_container_t
class.
The functions in this class for the most part return simple types. There are some functions
that return a more complex type, such as coot::simple_mesh_t
or coot::validation_information_t
.
Note: The private types, functions, attributes and members are listed here, but they are not, for the most part, useful for exporting. Which is not to rule out that there may be something there that could usefully be declared as public.
-
class molecules_container_t
- #include <molecules-container.hh>
the container of molecules. The class for all libcootapi functions.
Basic Utilities
-
bool make_backups_flag
the backup-enable state (raw public if needed/preferred)
-
inline void set_make_backups(bool state)
Allow the user to disable/enable backups
- Parameters:
state – is True to mean that it is enabled. The default is True.
-
inline bool get_make_backups() const
Get the state of the backups
- Returns:
the backup-enabled state
-
std::string file_name_to_string(const std::string &file_name) const
File name to string
- Parameters:
file_name – is the name of the given file
- Returns:
the string of the contents of the given file-name.
-
inline unsigned int get_number_of_molecules() const
Get the number of molecules (either map or model)
- Returns:
the number of molecules
-
void create_empty_molecules(unsigned int n_empty)
Add a number of empty molecules to the internal vector/list of molecules
Note this is not like STL
reserve
as it will increase the molecule index of the next added molecule byn_empty
.- Parameters:
n_empty – the number of empty molecules to create
-
inline void set_imol_refinement_map(int i)
Set the map used for refinement and fitting
- Parameters:
i – the map molecule index used for refinement and fitting
-
inline void set_map_weight(float w)
Set the map weight
- Parameters:
w – the map weight to be used for refinement, e.g. 50.0
-
inline float get_map_weight() const
Get the map weight
- Returns:
the map weight
-
coot::atom_spec_t atom_cid_to_atom_spec(int imol, const std::string &cid) const
Convert atom cid string to a coot atom specifier
- Parameters:
imol – is the model molecule index
cid – is the atom selection CID e.g “//A/15/OH” (atom OH in residue 15 of chain A)
- Returns:
the atom spec.,
spec.empty()
is true on failure.
-
coot::residue_spec_t residue_cid_to_residue_spec(int imol, const std::string &cid) const
Convert residue cid string to a coot residue specifier
- Parameters:
imol – is the model molecule index
cid – is the atom selection CID e.g “//A/15” (all the atoms in residue 15 of chain A)
- Returns:
the residues spec.,
spec.empty()
is true on failure.
-
inline void set_show_timings(bool s)
Set the show_timings flag
Various (not all) functions in this class can calculate how long they took to run. Setting this will write the time to taken (in milliseconds) to stdout.
- Parameters:
s – is True to mean that it is enabled. The default is True.
-
moorhen::header_info_t get_header_info(int imol) const
Get header info
(the header info is sparce at the moment)
- Parameters:
imol – is the model molecule index
- Returns:
an object with header info.
-
int get_imol_enc_any() const
Get imol_enc_any (enc: encoded)
imol_enc_any refers to the molecule number for dictionary that can be used with any molecule
- Returns:
the value of
imol_enc_any
Coordinates Utils
-
int read_coordinates(const std::string &file_name)
Read a coordinates file (mmcif or PDB)
- Parameters:
file_name – is the name of the coordinates file
- Returns:
the new molecule index on success and -1 on failure
-
int read_pdb(const std::string &file_name)
Read a PDB file (or mmcif coordinates file, despite the name)
It does the same job as
read_coordinates
but has (perhaps) a more familiar name- Parameters:
file_name – is the name of the coordinates file
- Returns:
the new molecule index on success and -1 on failure
-
void print_secondary_structure_info(int imol) const
Print the secondary structure information
- Parameters:
imol – is the model molecule index
-
void replace_molecule_by_model_from_file(int imol, const std::string &pdb_file_name)
Read a PDB file (or mmcif coordinates file, despite the name) to replace the current molecule
This will only work if the molecules is already a model molecule
- Parameters:
imol – is the model molecule index
file_name – is the name of the coordinates file
-
std::vector<int> split_multi_model_molecule(int imol)
Split an NMR model into multiple models
- Parameters:
imol – is the model molecule index
- Returns:
the vector/list of new molecule indices
-
int make_ensemble(const std::string &model_molecules_list)
Make a multi-model molecule given the input molecules
- Parameters:
model_molecules_list – is a colon-separated list of molecules, (e.g. “2:3:4”)
- Returns:
the new molecule index, -1 if no models were found in the model molecules list
-
std::string molecule_to_PDB_string(int imol) const
Get the molecule as a PDB string
- Parameters:
imol – is the model molecule index
- Returns:
the model molecule imol as a string. Return emtpy string on error
-
std::string molecule_to_mmCIF_string(int imol) const
Get the molecule as an mmCIF string
- Parameters:
imol – is the model molecule index
- Returns:
the model molecule imol as a string. Return emtpy string on error
-
std::pair<int, std::string> get_active_atom(float x, float y, float z, const std::string &displayed_model_molecules_list) const
Get the active atom given the screen centre
- Parameters:
x – is x position of the centre of the screen
y – is y position of the centre of the screen
z – is z position of the centre of the screen
displayed_model_molecules_list – is a colon-separated list of molecules, (e.g. “2:3:4”)
- Returns:
the molecule index and the atom cid. On failure (no molecules with atoms in them) then return -1 and a blank string.
-
int import_cif_dictionary(const std::string &cif_file_name, int imol_enc)
Import a dictionary cif
- Parameters:
imol_enc – is used to specify to which molecule this dictionary should apply. Use IMOL_ENC_ANY to mean “it applies to all molecules.” IMOL_ENC_ANY = -999999
- Returns:
1 on success and 0 on failure
-
std::string get_cif_file_name(const std::string &comp_id, int imol_enc) const
Get the cif file name
- Parameters:
comp_id – is the 3-letter code for the residue/ligand, e.g. “ALA” for alanine
imol_enc – is the molecule index for the residue type/compound_id
- Returns:
the dictionary read for the given residue type, return an empty string on failure
-
std::string get_cif_restraints_as_string(const std::string &comp_id, int imol_enc) const
Get the cif restraints as a string
- Parameters:
comp_id – is the 3-letter code for the residue/ligand, e.g. “ALA” for alanine
imol_enc – is the molecule index for the residue type/compound_id
- Returns:
a string that is the contents of a dictionary cif file
-
bool copy_dictionary(const std::string &monomer_name, int imol_current, int imol_new)
Copy the dictionary that is specific for
imol_current
so that it can be used with a new molecule- Parameters:
monomer_name – is the 3 letter code of the monomer in the dictionary, e.g. “ALA” for alanine
imol_current – is the model molecule index with the dictionary to be copied from
imol_new – is the model molecule index the dictionary will be copied into
-
int get_monomer(const std::string &monomer_name)
Get a monomer
- Parameters:
monomer_name – is the 3-letter code of the monomer in the dictionary, e.g. “ALA” for alanine
- Returns:
the new molecule index on success and -1 on failure
-
int get_monomer_from_dictionary(const std::string &comp_id, int imol, bool idealised_flag)
Get a monomer for a particular molecule
- Parameters:
comp_id – is the 3-letter code for the residue/ligand, e.g. “ALA” for alanine
imol – is the model molecule index, use -999999 (IMOL_ENC_ANY) if no molecule-specific dictionary is needed
idealised_flag – means that the coordinates have been minimised with a molecular modelling minimisation algo, usually the value is True
- Returns:
the new molecule index on success and -1 on failure
-
int get_monomer_and_position_at(const std::string &comp_id, int imol, float x, float y, float z)
Get monomer and place it at the given position for a particular molecule
- Parameters:
comp_id – is the 3-letter code for the residue/ligand, e.g. “ALA” for alanine
imol – is the model molecule index, use -999999 (IMOL_ENC_ANY) if no molecule-specific dictionary is needed
x – is the x value of the target position
y – is the y value of the target position
z – is the z value of the target position
- Returns:
the new molecule index on success and -1 on failure
-
std::map<std::string, std::string> dictionary_atom_name_map(const std::string &comp_id_1, int imol_1, const std::string &comp_id_2, int imol_2)
Match atom between 2 dictionaries
- Parameters:
comp_id_1 – is the 3-letter code for the residue/ligand in the first model, e.g. “ALA” for alanine
imol_1 – is the model molecule index of the first model
comp_id_2 – is the 3-letter code for the residue/ligand in the second model, e.g. “ALA” for alanine
imol_2 – is the model molecule index of the second model
- Returns:
the atom name match on superposing the atoms of the given dictionaries
-
std::vector<std::string> get_groups_for_monomers(const std::vector<std::string> &residue_names) const
Get the groups for a vector/list of monomers
e.g. “NON POLYMER”, “PEPTIDE”, etc
- Parameters:
residue_names – is a list of residue names, e.g. [“ALA”, “TRP”]
- Returns:
the group for the given list of residue names
-
std::string get_group_for_monomer(const std::string &residue_name) const
Get the group for a particular monomer
e.g. “NON POLYMER”, “PEPTIDE”, etc
- Parameters:
residue_name – is the is the 3-letter code for the residue, e.g. “ALA” for alanine
- Returns:
the group for the given residue name
-
std::string get_hb_type(const std::string &compound_id, int imol_enc, const std::string &atom_name) const
Get the hydrogen bond type of a particular atom in a given residue type
- Parameters:
compound_id – is the 3-letter code for the residue/ligand in the first model, e.g. “TYR” for tyrosine
imol_enc – is the molecule index for the residue type/compound_id
atom_name – is the name of the atom, e.g. “OH”
- Returns:
the hb_type for the given atom. On failure return an empty string. Valid types are: “HB_UNASSIGNED” ,”HB_NEITHER”, “HB_DONOR”, “HB_ACCEPTOR”, “HB_BOTH”, “HB_HYDROGEN”.
-
std::vector<std::pair<std::string, std::string>> get_gphl_chem_comp_info(const std::string &compound_id, int imol_enc)
Get the GPhL extra restraint information (from the input cif file)
- Parameters:
compound_id – is the 3-letter code for the residue/ligand in the first model, e.g. “TYR” for tyrosine
imol_enc – is the molecule index for the residue type/compound_id
- Returns:
a vector/list of string pairs that were part of a gphl_chem_comp_info. return an empty vector/list on failure to find any such info.
-
std::vector<std::pair<std::string, std::string>> get_acedrg_atom_types(const std::string &compound_id, int imol_enc) const
Get a list of atom names and their associated AceDRG atom types
- Parameters:
compound_id – is the 3-letter code for the residue/ligand in the first model, e.g. “TYR” for tyrosine
imol_enc – is the molecule index for the residue type/compound_id
- Returns:
a list of atom names and their associated AceDRG atom types, return an empty list on failure (e.g. when atoms types are not in the dictionary)
-
coot::acedrg_types_for_residue_t get_acedrg_atom_types_for_ligand(int imol, const std::string &residue_cid) const
Get AceDRG atom types for ligand bonds
- Parameters:
imol – is the model molecule index
residue_cid – is the atom selection CID e.g “//A/15” (all the atoms in residue 15 of chain A)
- Returns:
a
coot::acedrg_types_for_residue_t
- which contains a vector/list of bond descriptions.
-
void write_png(const std::string &compound_id, int imol, const std::string &file_name) const
Write a PNG for the given compound_id.
Currently this function does nothing (drawing is done with the not-allowed cairo)
- Parameters:
compound_id – is the 3-letter code for the residue/ligand in the first model, e.g. “TYR” for tyrosine
imol – is the model molecule index
-
int write_coordinates(int imol, const std::string &file_name) const
Write the coordinates to the given file name
- Parameters:
imol – is the model molecule index
file_name – is the name of the new model
- Returns:
1 on success and 0 on failure
-
void set_draw_missing_residue_loops(bool state)
Set the state for drawing missing loops
By default missing loops are drawn. This function allows missing loops to not be drawn. Sometimes that can clarify the representation. This is a lightweight function that sets a flag that is used by subsequent calls to get_bonds_mesh()
- Parameters:
state – is True to mean it is enabled
-
coot::simple_mesh_t get_bonds_mesh(int imol, const std::string &mode, bool against_a_dark_background, float bond_width, float atom_radius_to_bond_width_ratio, int smoothness_factor)
Get the bonds mesh
- Parameters:
mode – is “COLOUR-BY-CHAIN-AND-DICTIONARY”, “CA+LIGANDS” or “VDW-BALLS”
against_a_dark_background – allows the bond colours to be relevant for the background. When the background is dark, the colours should (as a rule) be bright and pastelly. When the background is light/white, the colour are darker and more saturated.
smoothness_factor – controls the number of triangles used to make the bond cylinders and spheres for the atoms - it rises in powers of 4. 1 is the smallest smoothness_factor, 2 looks nice and 3 is best.
bond_width – is the bond width in Angstroms. 0.12 is a reasonable default value.
atom_radius_to_bond_width_ratio – allows the representation of “ball and stick”. To do so use a value between 1.5 and 3.0. The ratio for “liquorice” representation is 1.0.
- Returns:
a
simple_mesh_t
-
coot::instanced_mesh_t get_bonds_mesh_instanced(int imol, const std::string &mode, bool against_a_dark_background, float bond_width, float atom_radius_to_bond_width_ratio, int smoothness_factor)
Get the instanced bonds mesh.
- Parameters:
mode – is “COLOUR-BY-CHAIN-AND-DICTIONARY” - more modes to follow
against_a_dark_background – allows the bond colours to be relevant for the background. When the background is dark, the colours should (as a rule) be bright and pastelly. When the background is light/white, the colour are darker and more saturated.
smoothness_factor – controls the number of triangles used to make the bond cylinders and spheres for the atoms - it rises in powers of 4. 1 is the smallest smoothness_factor, 2 looks nice and 3 is best. Instancing may mean that smoothness factor 3 should be used by default.
bond_width – is the bond width in Angstroms. 0.12 is a reasonable default value.
atom_radius_to_bond_width_ratio – allows the representation of “ball and stick”. To do so use a value between 1.5 and 3.0. The ratio for “liquorice” representation is 1.0.
- Returns:
a
instanced_mesh_t
-
coot::instanced_mesh_t get_bonds_mesh_for_selection_instanced(int imol, const std::string &atom_selection_cid, const std::string &mode, bool against_a_dark_background, float bond_width, float atom_radius_to_bond_width_ratio, int smoothness_factor)
As
get_bonds_mesh_instanced
above, but only return the bonds for the atom selection. Typically one would call this with a wider bond width than one would use for standards atoms (all molecule)- Parameters:
atom_selection_cid – e.g. “//A/15” (all the atoms in residue 15 of chain A)
- Returns:
a
instanced_mesh_t
-
coot::instanced_mesh_t get_goodsell_style_mesh_instanced(int imol, float colour_wheel_rotation_step, float saturation, float goodselliness)
Get the Goodsell style mesh
- Parameters:
colour_wheel_rotation_step' – the amount, in degrees, that the colour wheel advances between different chains. 97 degrees is a reasonable starting value
saturation' – a number between 0 and 1, where 0 is grey and 1 is “lego-like” colour scheme. 0.5 is a nice middle value
goodselliness – is the degree to which the C-atoms are desaturated, 0.3 is a reasonable value
- Returns:
a
instanced_mesh_t
-
void export_map_molecule_as_gltf(int imol, float pos_x, float pos_y, float pos_z, float radius, float contour_level, const std::string &file_name)
Export map molecule as glTF file
glTF files can be imported into Blender or other 3D graphics applications
- Parameters:
imol – is the model molecule index
x – x of the center of the screen
y – y of the center of the screen
z – z of the center of the screen
radius – e.g. 12.0 for X-ray map and 100.0 for Cryo-EM map
contour_level – e.g. 1.5 sd for X-ray, 5.0 sd for cryo-EM
file_name – extension should be .glb
-
void export_model_molecule_as_gltf(int imol, const std::string &selection_cid, const std::string &mode, bool against_a_dark_background, float bonds_width, float atom_radius_to_bond_width_ratio, int smoothness_factor, bool draw_hydrogen_atoms_flag, bool draw_missing_residue_loops, const std::string &file_name)
Export model molecule as glTF file
glTF files can be imported into Blender or other 3D graphics applications
Same parameters as the
get_bonds_mesh
function.draw_hydrogen_atoms_flag
anddraw_missing_residue_loops
are typically False. This API will change - we want to specify surfaces and ribbons too.
-
void export_molecular_representation_as_gltf(int imol, const std::string &atom_selection_cid, const std::string &colour_scheme, const std::string &style, int secondary_structure_usage_flag, const std::string &file_name)
Export molecular representation as glTF file
glTF files can be imported into Blender
- Parameters:
imol – is the model molecule index
atom_selection_cid – : e.g “//A/15” (all the atoms in residue 15 of chain A)
colour_scheme – is one of “colorRampChainsScheme”, “colorBySecondaryScheme”, “Chain”
style – “Ribbon” or “MolecularSurface”
secondary_structure_usage_flag – 0 (USE_HEADER), 1 (DONT_USE) or 2 (CALC_SECONDARY_STRUCTURE)
-
void export_chemical_features_as_gltf(int imol, const std::string &cid, const std::string &file_name) const
export chemical features for the specified residue
-
std::vector<glm::vec4> get_colour_table(int imol, bool against_a_dark_background) const
Get colour table (for testing)
- Parameters:
imol – is the model molecule index
against_a_dark_background – allows the bond colours to be relevant for the background.
- Returns:
the colour table
-
void set_colour_wheel_rotation_base(int imol, float r)
Set the colour wheel rotation base for the specified molecule (in degrees)
- Parameters:
imol – is the model molecule index
r – is the rotation angle in degrees
-
void set_base_colour_for_bonds(int imol, float r, float g, float b)
Set the base colour - to be used as a base for colour wheel rotation
RGB colour codes, e.g. green is r:0, g: 255, b:0
- Parameters:
imol – is the model molecule index
-
void add_to_non_drawn_bonds(int imol, const std::string &atom_selection_cid)
Add an atom selection cid for atoms and bonds not to be drawn
- Parameters:
imol – is the model molecule index
atom_selection_cid – e.g “//A/15” (all the atoms in residue 15 of chain A)
-
void clear_non_drawn_bonds(int imol)
Clear the set of non-drawn atoms (so that they can be displayed again)
- Parameters:
imol – is the model molecule index
-
void print_non_drawn_bonds(int imol) const
Print non-drawn bonds
- Parameters:
imol – is the model molecule index
-
void set_user_defined_bond_colours(int imol, const std::map<unsigned int, std::array<float, 4>> &colour_map)
User-defined colour-index to colour.
-
void set_user_defined_atom_colour_by_selection(int imol, const std::vector<std::pair<std::string, unsigned int>> &indexed_residues_cids, bool colour_applies_to_non_carbon_atoms_also)
Set the user-defined residue selections (CIDs) to colour index
- Parameters:
imol – is the model molecule index
-
void add_colour_rule(int imol, const std::string &selection_cid, const std::string &colour)
Add a colour rule for M2T representations.
-
void add_colour_rules_multi(int imol, const std::string &selections_and_colours_combo_string)
Add multiple colour rules
- Parameters:
selections_and_colours_combo_string – e.g. “//A/1^#cc0000|//A/2^#cb0002|//A/3^#c00007”, where “|” is the separator for each rule and “^” is the separator for the selection string and the colour string
-
void delete_colour_rules(int imol)
Delete the colour rules for the given molecule
- Parameters:
imol – is the model molecule index
-
std::vector<std::pair<std::string, std::string>> get_colour_rules(int imol) const
Get the colour rules
- Parameters:
imol – is the model molecule index
-
void print_colour_rules(int imol) const
Print the colour rules
- Parameters:
imol – is the model molecule index
-
void set_use_bespoke_carbon_atom_colour(int imol, bool state)
Use bespoke carbon atom colour
- Parameters:
imol – is the model molecule index
-
void set_bespoke_carbon_atom_colour(int imol, const coot::colour_t &col)
Set bespoke carbon atom colour
- Parameters:
imol – is the model molecule index
-
void M2T_updateFloatParameter(int imol, const std::string ¶m_name, float value)
Update float parameter for MoleculesToTriangles molecular mesh.
-
void M2T_updateIntParameter(int imol, const std::string ¶m_name, int value)
Update int parameter for MoleculesToTriangles molecular mesh.
-
coot::simple_mesh_t get_molecular_representation_mesh(int imol, const std::string &cid, const std::string &colour_scheme, const std::string &style, int secondary_structure_usage_flag)
Get ribbon and molecular surface representation
- Parameters:
imol – is the model molecule index
cid – is the atom selection CID e.g “//A/15” (all the atoms in residue 15 of chain A)
colour_scheme – should be one of “colorRampChainsScheme”, “colorBySecondaryScheme”, “Chain”
style – “Ribbon” or “MolecularSurface”
secondary_structure_usage_flag – 0 (USE_HEADER), 1 (DONT_USE) or 2 (CALC_SECONDARY_STRUCTURE).
- Returns:
a
simple_mesh_t
-
coot::simple_mesh_t get_gaussian_surface(int imol, float sigma, float contour_level, float box_radius, float grid_scale, float b_factor) const
Get a Gaussian surface representation
- Parameters:
imol – is the model molecule index
sigma – default 4.4
contour_level – default 4.0
box_radius – default 5.0
grid_scale – default 0.7
b_factor – default 100.0 (use 0.0 for no FFT-B-factor smoothing)
- Returns:
a
simple_mesh_t
composed of a number of Gaussian surfaces (one for each chain)
-
coot::simple_mesh_t get_chemical_features_mesh(int imol, const std::string &cid) const
Get chemical features for the specified residue
- Parameters:
imol – is the model molecule index
cid – is the atom selection CID e.g “//A/15” (all the atoms in residue 15 of chain A)
- Returns:
a
simple_mesh_t
-
std::vector<double> get_residue_CA_position(int imol, const std::string &cid) const
get the residue CA position
- Returns:
a vector. The length of the vector is 0 on failure, otherwise it is the x,y,z values
-
std::vector<double> get_residue_average_position(int imol, const std::string &cid) const
get the avarge residue position
- Returns:
a vector. The length of the vector is 0 on failure, otherwise it is the x,y,z values
-
std::vector<double> get_residue_sidechain_average_position(int imol, const std::string &cid) const
get the avarge residue side-chain position
- Returns:
a vector. The length of the vector is 0 on failure, otherwise it is the x,y,z values
-
unsigned int get_number_of_atoms(int imol) const
Get number of atoms
- Parameters:
imol – is the model molecule index
- Returns:
the number of atoms in the specified model, or 0 on error
-
float get_molecule_diameter(int imol) const
Get molecule diameter
- Parameters:
imol – is the model molecule index
- Returns:
an estimate of the diameter of the model molecule (-1 on failure)
-
int get_number_of_hydrogen_atoms(int imol) const
Get number of hydrogen atoms
- Parameters:
imol – is the model molecule index
- Returns:
the number of hydrogen atoms in the specified model, or -1 on error
-
std::vector<std::string> get_chains_in_model(int imol) const
Get the chain IDs in the given molecule
- Parameters:
imol – is the model molecule index
- Returns:
vector/list of chain-ids for the given molecule
Get the chains that are related by NCS or molecular symmetry
- Parameters:
imol – is the model molecule index
- Returns:
a vector/list of vector/list of chain ids, e.g. [[A,C], [B,D]] (for hemoglobin).
-
std::vector<std::pair<coot::residue_spec_t, std::string>> get_single_letter_codes_for_chain(int imol, const std::string &chain_id) const
Get the single letter codes for the residues in the specified chain
- Parameters:
imol – is the model molecule index
chain_id – e.g. “A” for chain A
- Returns:
vector/list of single letter codes - in a pair with the given residue spec
-
std::vector<std::string> get_residue_names_with_no_dictionary(int imol) const
Get a list of residues that don’t have a dictionary
- Parameters:
imol – is the model molecule index
- Returns:
a list of residue names (compound_ids) for which there is no dictionary in the geometry store
-
std::string get_residue_name(int imol, const std::string &chain_id, int res_no, const std::string &ins_code) const
Get residue name
- Parameters:
imol – is the model molecule index
chain_id – e.g. “A” for chain A
res_no – is the residue number, e.g. 12
ins_code – is the insertion code, e.g. “A”
- Returns:
the residue name, return a blank string on residue not found.
-
std::vector<coot::residue_spec_t> residues_with_missing_atoms(int imol)
Get residues with missing atoms
- Parameters:
imol – is the model molecule index
- Returns:
an object that has information about residues without dictionaries and residues with missing atom in the the specified molecule
-
std::vector<coot::residue_spec_t> get_residues_near_residue(int imol, const std::string &residue_cid, float dist) const
Get a list of residues specs that have atoms within distance of the atoms of the specified residue
- Parameters:
imol – is the model molecule index
residue_cid – is the atom selection CID e.g “//A/15” (all the atoms in residue 15 of chain A)
dist – is the distance in Angstrom
- Returns:
a list of residue specs
-
superpose_results_t SSM_superpose(int imol_ref, const std::string &chain_id_ref, int imol_mov, const std::string &chain_id_mov)
Superposition (using SSM)
The specified chain of the moving molecule is superposed onto the chain in the reference molecule (if possible).
- Parameters:
imol_ref – the reference model molecule index
chain_id_ref – the chain ID for the reference chain
imol_mov – the moving model molecule index
chain_id_mov – the chain ID for the moving chain
-
void add_lsq_superpose_match(const std::string &chain_id_ref, int res_no_ref_start, int res_no_ref_end, const std::string &chain_id_mov, int res_no_mov_start, int res_no_mov_end, int match_type)
Superpose using LSQ - setup the matches
- Parameters:
chain_id_ref – the chain ID for the reference chain
res_no_ref_start – the starting residue number in the reference chain
res_no_ref_end – the ending residue number in the reference chain
chain_id_mov – the chain ID for the moving chain
res_no_mov_start – the starting residue number in the moving chain
res_no_mov_end – the ending residue number in the moving chain
match_type – 0: all, 1: main, 2: CAs, 3: N, CA, C
-
void add_lsq_superpose_atom_match(const std::string &chain_id_ref, int res_no_ref, const std::string &atom_name_ref, const std::string &chain_id_mov, int res_no_mov, const std::string &atom_name_mov)
-
void clear_lsq_matches()
Clear any existing lsq matchers.
-
void lsq_superpose(int imol_ref, int imol_mov)
Apply the superposition using LSQ
- Parameters:
imol_ref – the reference model molecule index
imol_mov – the moving model molecule index
-
int transform_map_using_lsq_matrix(int imol_map, lsq_results_t lsq_matrix, float x, float y, float z, float radius)
transform a map and create a new map
- Returns:
the molecule index of the new map, -1 for failure
-
lsq_results_t get_lsq_matrix(int imol_ref, int imol_mov, bool summary_to_screen) const
Get LSQ matrix
don’t apply it to the coordinates
- Parameters:
imol_ref – the reference model molecule index
imol_mov – the moving model molecule index
summary_to_screen – if True, write a summary of the statistics to the output
- Returns:
the transformation matrix in a simple class
-
coot::symmetry_info_t get_symmetry(int imol, float symmetry_search_radius, float centre_x, float centre_y, float centre_z) const
Get symmetry
now comes in a simple container that also includes the cell
-
::api::cell_t get_cell(int imol) const
Get the cell
Check that
is_set
is true before use.- Parameters:
imol – is the model molecule index
- Returns:
a
cell_t
-
coot::util::map_molecule_centre_info_t get_map_molecule_centre(int imol) const
Get the middle of the “molecule blob” in cryo-EM reconstruction maps
- Parameters:
imol – is the map molecule index
- Returns:
-
int undo(int imol)
Undo
- Parameters:
imol – is the model molecule index
- Returns:
1 on successful undo, return 0 on failure
-
int redo(int imol)
Redo
- Parameters:
imol – is the model molecule index
- Returns:
1 on successful redo, return 0 on failure
Testing functions
-
ltj_stats_t long_term_job_stats
long term job
-
bool interrupt_long_term_job
-
void testing_start_long_term_job(unsigned int n_seconds)
Testing function
start a long-term job.
- Parameters:
n_seconds – is the number of seconds, if is 0, then run forever (or until interrupted)
-
void testing_stop_long_term_job()
Testing function
stop the long-term job runnning
-
inline ltj_stats_t testing_interrogate_long_term_job()
Testing function
get the stats for the long-term job
-
inline double get_contouring_time() const
Testing function
get the time for contouring in milliseconds
-
void set_max_number_of_threads(unsigned int n_threads)
Testing function
set the maximum number of threads for both the thread pool and the vector of threads
- Parameters:
n_threads – is the number of threads
-
void set_max_number_of_threads_in_thread_pool(unsigned int n_threads)
Testing function
Deprecated name for the “set_max_number_of_threads()” function
-
double test_the_threading(int n_threads)
Testing function
get the time to run a test function in milliseconds
- Parameters:
n_threads – is the number of threads
-
double test_launching_threads(unsigned int n_threads_per_batch, unsigned int n_batches) const
Testing function
- Parameters:
n_threads_per_batch – is the number of threads per batch
n_batches – is the number batches
- Returns:
the time per batch in microseconds
-
double test_thread_pool_threads(unsigned int n_threads)
Testing function
- Parameters:
n_threads – is the number of threads
- Returns:
time in microseconds
-
int mmcif_tests(bool last_test_only)
Testing function
a test for mmdb/gemmi/mmcif functionality
- Parameters:
last_test_only – is True to mean that only that last test should be run. The default is False. This is useful to set to True while a test is being developed.
- Returns:
the success status: 1 means that all the tests passed.
Generic Utils
-
std::string get_molecule_name(int imol) const
Get the molecule name
- Parameters:
imol – is the model molecule index
- Returns:
the name of the molecule
-
void set_molecule_name(int imol, const std::string &new_name)
Set the molecule name
- Parameters:
imol – is the model or map molecule index
new_name – is the new name of the model or map
-
void display_molecule_names_table() const
Debugging function: display the table of molecule and names.
-
bool is_valid_model_molecule(int imol) const
Check if the model index is valid
e.g. if the molecule is a map you will have an invalid model
- Parameters:
imol – is the model molecule index
- Returns:
True or False
-
bool is_valid_map_molecule(int imol_map) const
Check if the map index is valid
e.g. if the map is a model you will have an invalid map
- Parameters:
imol_map – is the map molecule index
- Returns:
True or False
-
bool is_a_difference_map(int imol_map) const
Check if it the map is a difference map
- Parameters:
imol_map – is the map molecule index
- Returns:
True or False
-
int new_molecule(const std::string &name)
Create an empty molecule
- Returns:
the index of the new molecule
-
int close_molecule(int imol)
Close the molecule (and delete dynamically allocated memory)
- Parameters:
imol – is the model molecule index
- Returns:
1 on successful closure and 0 on failure to close
-
void end_delete_closed_molecules()
Delete the most recent/last closed molecule in the molecule vector, until the first non-closed molecule is found (working from the end)
-
void pop_back()
Delete the most recent/last molecule in the molecule vector.
-
void clear()
Delete all molecules.
-
std::vector<double> get_eigenvalues(int imol, const std::string &chain_id, int res_no, const std::string &ins_code)
Get the eigenvalues of the specified residue
- Parameters:
imol – is the model molecule index
chain_id – e.g. “A” for chain A
res_no – is the residue number, e.g. 12
ins_code – is the insertion code, e.g. “A”
- Returns:
the eigenvalues of the atoms in the specified residue
-
coot::simple_mesh_t test_origin_cube() const
Get a simple test mesh
- Returns:
the mesh of a unit solid cube at the origin
-
inline void fill_rotamer_probability_tables()
Fill the rotamer probability tables (currently not ARG and LYS)
-
void accept_rotamer_probability_tables_compressed_data(const std::string &data_stream)
Access to a compressed file that contains the rotamer probabilities
libcootapi will fill the rotamer probabilities tables from this compressed data stream (placeholder only)
Backup and Saving
-
inline bool contains_unsaved_models() const
Check if there are unsaved changes for this model
e.g. as yet not written to disk
- Returns:
a flag of unsaved models state - e.g. if any of them are unsaved, then this returns True.
-
inline void save_unsaved_model_changes()
Save the unsaved model - this function has not yet been written!
Geometry and Dictionaries
-
void geometry_init_standard()
Read the standard list of residues.
-
std::vector<std::string> non_standard_residue_types_in_model(int imol) const
Get a list of non-standard residues in the given molecule
- Parameters:
imol – is the model molecule index
- Returns:
a vector/list of non-standard residues
Map Utils
-
inline float get_map_sampling_rate()
Get map sampling rate
- Returns:
the map sampling rate, the default is 1.8
-
inline void set_map_sampling_rate(float msr)
Set the map sampling rate
Higher numbers mean smoother maps, but they take longer to generate, longer to transfer, longer to parse and longer to draw
- Parameters:
msr – is the map sampling rate to set, the default is 1.8
-
int read_mtz(const std::string &file_name, const std::string &f, const std::string &phi, const std::string &weight, bool use_weight, bool is_a_difference_map)
Read the given mtz file
- Parameters:
file_name – is the name of the MTZ file
f – F column, “FWT”
phi – phi column, “PHWT”
weight – weight column, “W”
use_weight – flag for weights usage, False
is_a_difference_map – False
- Returns:
the new molecule number or -1 on failure
-
int replace_map_by_mtz_from_file(int imol, const std::string &file_name, const std::string &f, const std::string &phi, const std::string &weight, bool use_weight)
Replace map by mtz
- Parameters:
imol – is the map molecule index
f – F column, “FWT”
phi – phi column, “PHWT”
weight – weight column, “W”
use_weight – flag for weights usage, False
-
std::vector<auto_read_mtz_info_t> auto_read_mtz(const std::string &file_name)
Auto read the given mtz file
- Parameters:
file_name – is the name of the MTZ file
- Returns:
a vector of the maps created from reading the file
-
int read_ccp4_map(const std::string &file_name, bool is_a_difference_map)
Read a CCP4 (or MRC) map
There is currently a size limit of 1000 pixels per edge.
- Parameters:
file_name – is the name of the map file
is_a_difference_map – False
- Returns:
the new molecule number or -1 on failure
-
int write_map(int imol, const std::string &file_name) const
Write a map
- Parameters:
imol – is the map molecule index
file_name – is the name of the new map file
- Returns:
1 on a successful write, return 0 on failure.
-
float get_map_mean(int imol) const
Get map mean
- Parameters:
imol – is the map molecule index
- Returns:
the mean of the map or -1 if imol is not a valid map molecule index
-
float get_map_rmsd_approx(int imol_map) const
Get map rmsd approx
the function is approximate because the epsilon factor is not taken into account
- Parameters:
imol_map – is the map molecule index
- Returns:
the map rmsd. -1 is returned if imol_map is not a valid map molecule index.
-
coot::molecule_t::histogram_info_t get_map_histogram(int imol, unsigned int n_bins, float zoom_factor) const
Get map histogram
- Parameters:
imol – is the map molecule index
n_bins – is the number of bins - 200 is a reasonable default.
zoom_factor – (reduces the range by the given factor) centred around the median (typically 1.0 but usefully can vary until ~20.0).
- Returns:
the map histogram
-
float get_suggested_initial_contour_level(int imol) const
- Parameters:
imol – is the map molecule index
- Returns:
the suggested initial contour level. Return -1 on not-a-map
-
bool is_EM_map(int imol) const
Check if a map is an EM map or not
- Parameters:
imol – is the map molecule index
- Returns:
the “EM” status of this molecule. Return false on not-a-map.
-
int sharpen_blur_map(int imol_map, float b_factor, bool in_place_flag)
Create a new map that is blurred/sharpened
- Parameters:
imol – is the map molecule index
b_factor – e.g. 100.0
in_place_flag – True if you want to replace the current map, False if you want to create a new map
- Returns:
the molecule index of the new map or -1 on failure or if
in_place_flag
was true.
-
int sharpen_blur_map_with_resample(int imol_map, float b_factor, float resample_factor, bool in_place_flag)
Create a new map that is blurred/sharpened and resampling
Note that resampling can be slow, a resample_factor of 1.5 is about the limit of the trade of of prettiness for speed.
- Parameters:
imol_map – is the map molecule index
b_factor – e.g. 100.0
resample_factor – e.g. 1.4
in_place_flag – True if you want to replace the current map, False if you want to create a new map
- Returns:
the molecule index of the new map or -1 on failure or if
in_place_flag
was true.
-
int mask_map_by_atom_selection(int imol_coords, int imol_map, const std::string &cid, float atom_radius, bool invert_flag)
Mask map by atom selection
(note the argument order is reversed compared to the coot api).
- Parameters:
imol_coords – is the model molecule index
imol_map – is the map molecule index
cid – is the atom selection CID e.g “//A/15” (all the atoms in residue 15 of chain A)
atom_radius – is the atom radius. Use a negative number to mean “default”.
invert_flag – changes the parts of the map that are masked, so to highlight the density for a ligand one would pass the cid for the ligand and invert_flag as true, so that the parts of the map that are not the ligand are suppressed.
- Returns:
the index of the new map or -1 on failure
-
std::vector<int> partition_map_by_chain(int imol_map, int imol_model)
Partition the input map
Each voxel in the map is assigned to the chain to which it is nearest. Unlike masking, the generated maps are not restricted to be “close” to the atoms in the atom selection.
e.g. maskChains for ChimeraX - JiangLab
- Parameters:
imol_map – is the map molecule index
imol_model – is the model molecule index
- Returns:
a vector/list of the molecules indices of the newly created maps
-
int make_mask(int imol_map_ref, int imol_model, const std::string &atom_selection_cid, float radius)
Make a masked map
- Parameters:
imol_map_ref – is the map molecule index
imol_model – is the model molecule index
atom_selection_cid – is the atom selection CID e.g “//A/15” (all the atoms in residue 15 of chain A)
radius – is the radius of the map, e.g. 12.0 for X-ray map and 100.0 for Cryo-EM map
- Returns:
the index of the newly created mask. Return -1 on failure.
-
int flip_hand(int imol_map)
Generate a new map which is the hand-flipped version of the input map
- Parameters:
imol_map – is the map molecule index
- Returns:
the molecule index of the new map, or -1 on failure.
-
std::vector<int> make_masked_maps_split_by_chain(int imol, int imol_map)
Make a vector/list of maps that are split by chain-id of the input imol
- Parameters:
imol – is the model molecule index
imol_map – is the map molecule index
- Returns:
a vector/list of the map molecule indices.
-
void set_map_colour(int imol, float r, float g, float b)
Set the map colour
The next time a map mesh is requested, it will have this colour. This does not affect the colour of the difference maps.
RGB colour codes, e.g. green is r:0, g: 255, b:0
- Parameters:
imol – is the model molecule index
-
void set_map_is_contoured_with_thread_pool(bool state)
Set the state of the mode of the threading in map contouring
- Parameters:
state – is True to mean that it is enabled. The default is True.
-
coot::simple_mesh_t get_map_contours_mesh(int imol, double position_x, double position_y, double position_z, float radius, float contour_level)
Get the mesh for the map contours
This function is not const because the internal state of a
coot_molecule_t
is changed.- Parameters:
imol – is the model molecule index
position_x – is the x coordinate of the target position
position_y – is the y coordinate of the target position
position_z – is the z coordinate of the target position
radius – is the radius of the map, e.g. 12.0 for X-ray map and 100.0 for Cryo-EM map
contour_level – e.g. 1.5 sd for X-ray, 5.0 sd for cryo-EM
- Returns:
a
simple_mesh_t
for the map contours of the specified map
-
coot::simple_mesh_t get_map_contours_mesh_using_other_map_for_colours(int imol_ref, int imol_map_for_colouring, double position_x, double position_y, double position_z, float radius, float contour_level, float other_map_for_colouring_min_value, float other_map_for_colouring_max_value, bool invert_colour_ramp)
Get the mesh for the map contours using another map for colouring
- Parameters:
imol_ref – is the reference map index
imol_map_for_coloring – is the target map index
position_x – is the x coordinate of the target position
position_y – is the y coordinate of the target position
position_z – is the z coordinate of the target position
radius – is the radius of the map, e.g. 12.0 for X-ray map and 100.0 for Cryo-EM map
contour_level – e.g. 1.5 sd for X-ray, 5.0 sd for cryo-EM
other_map_for_colouring_min_value – e.g. -1.0 in the case of correlation map
other_map_for_colouring_max_value – e.g. 1.0 in the case of correlation map
invert_colour_ramp – e.g. red to blue rather than blue to red
- Returns:
a
simple_mesh_t
for the map contours of the specified map
-
void set_map_colour_saturation(int imol, float s)
Set the map saturation
- Parameters:
s – is the map saturation, e.g. a number between 0 and 1, where 0 is grey and 1 is “lego-like” colour scheme. 0.5 is a nice middle value
-
inline coot::util::sfcalc_genmap_stats_t get_latest_sfcalc_stats() const
Get the latest sfcalc stats
- Returns:
a sfcalc_genmap_stats_t object
-
r_factor_stats get_r_factor_stats()
Get the R-factors
- Returns:
an object
r_factor_stats
-
std::string r_factor_stats_as_string(const r_factor_stats &rfs) const
Get the R factor stats as a string
- Parameters:
rfs – is the name of the string
- Returns:
a string with the R-factor stats
-
int average_map(const std::string &imol_maps, std::vector<float> &scales)
Get the average map
This function does no normalization of the scales, presuming that they are pre-normalized.
- Parameters:
imol_maps – is a colon-separated list of map indices e.g. “2:3:4”
scales – is the list of weights corresponding to the list of maps. The number of scales factors should match the number of maps
- Returns:
the index of the new map, or -1 on failure.
-
bool regen_map(int imol_map, const std::string &imol_maps, const std::vector<float> &scales)
This function does no normalisation of the scales, presuming that they are pre-normalized.
The number of maps in imol_maps should match the size of the scale vector. If not, nothing will happen and False will be returned
- Parameters:
imol_map – is the map molecule index
imol_maps – is a colon-separated list of map indices e.g. “2:3:4”
scales – is the list of weights corresponding to the list of maps
- Returns:
the success status
Coordinates Modelling
-
int auto_fit_rotamer(int imol, const std::string &chain_id, int res_no, const std::string &ins_code, const std::string &alt_conf, int imol_map)
Auto-fit rotamer
- Parameters:
imol – is the model molecule index
chain_id – e.g. “A” for chain A
res_no – is the residue number, e.g. 12
ins_code – is the insertion code, e.g. “A”
alt_conf – is the alternate conformation, e.g. “A” or “B”
imol_map – is the map molecule index
- Returns:
1 on successful modification, return 0 on failure
-
coot::molecule_t::rotamer_change_info_t change_to_next_rotamer(int imol, const std::string &residue_cid, const std::string &alt_conf)
Change to the next rotamer (rotamer cycling is implicit if needed)
- Parameters:
imol – is the model molecule index
residue_cid – is the atom selection CID e.g “//A/15” (all the atoms in residue 15 of chain A)
alt_conf – is the alternate conformation, e.g. “A” or “B”
- Returns:
the change information.
-
coot::molecule_t::rotamer_change_info_t change_to_previous_rotamer(int imol, const std::string &residue_cid, const std::string &alt_conf)
Change to the previous rotamer (rotamer cycling is implicit if needed)
- Parameters:
imol – is the model molecule index
residue_cid – is the atom selection CID e.g “//A/15” (all the atoms in residue 15 of chain A)
alt_conf – is the alternate conformation, e.g. “A” or “B”
- Returns:
the change information.
-
coot::molecule_t::rotamer_change_info_t change_to_first_rotamer(int imol, const std::string &residue_cid, const std::string &alt_conf)
Change to the first (0th) rotamer
- Parameters:
imol – is the model molecule index
residue_cid – is the atom selection CID e.g “//A/15” (all the atoms in residue 15 of chain A)
alt_conf – is the alternate conformation, e.g. “A” or “B”
- Returns:
the change information.
-
std::pair<int, unsigned int> delete_using_cid(int imol, const std::string &cid, const std::string &scope)
Delete item
- Parameters:
imol – is the model molecule index
cid – is the selection CID e.g “//A/15” (all the atoms in residue 15 of chain A)
scope – is one of the strings: [“ATOM”, “WATER”, “RESIDUE”,” CHAIN”,” MOLECULE”, “LITERAL”]
- Returns:
1 on successful deletion, return 0 on failure
-
std::pair<int, unsigned int> delete_atom(int imol, const std::string &chain_id, int res_no, const std::string &ins_code, const std::string &atom_name, const std::string &alt_conf)
Delete atom
- Parameters:
imol – is the model molecule index
chain_id – e.g. “A” for chain A
res_no – is the residue number, e.g. 12
ins_code – is the insertion code, e.g. “A”
atom_name – is the name of the atom, e.g. “OH”
alt_conf – is the alternate conformation, e.g. “A” or “B”
- Returns:
1 on successful deletion, return 0 on failure to delete.
-
std::pair<int, unsigned int> delete_atom_using_cid(int imol, const std::string &cid)
Delete atom using cid
- Parameters:
imol – is the model molecule index
cid – is the atom selection CID e.g “//A/15/OH” (atom OH in residue 15 of chain A)
- Returns:
1 on successful deletion, return 0 on failure to delete.
-
std::pair<int, unsigned int> delete_residue(int imol, const std::string &chain_id, int res_no, const std::string &ins_code)
Delete residue
- Parameters:
imol – is the model molecule index
chain_id – e.g. “A” for chain A
res_no – is the residue number, e.g. 12
ins_code – is the insertion code, e.g. “A”
- Returns:
1 on successful deletion, return 0 on failure to delete.
-
std::pair<int, unsigned int> delete_residue_using_cid(int imol, const std::string &cid)
Delete residue using cid
- Parameters:
imol – is the model molecule index
cid – is the residue selection CID e.g “//A/15” (all the atoms in residue 15 of chain A)
- Returns:
1 on successful deletion, return 0 on failure to delete.
-
std::pair<int, unsigned int> delete_residue_atoms_with_alt_conf(int imol, const std::string &chain_id, int res_no, const std::string &ins_code, const std::string &alt_conf)
Delete residue atoms using alt_conf
- Parameters:
imol – is the model molecule index
chain_id – e.g. “A” for chain A
res_no – is the residue number, e.g. 12
ins_code – is the insertion code, e.g. “A”
alt_conf – is the alternate conformation, e.g. “A” or “B”
- Returns:
1 on successful deletion, return 0 on failure to delete.
-
std::pair<int, unsigned int> delete_residue_atoms_using_cid(int imol, const std::string &cid)
Delete residue atoms using cid
This is the same as
delete_atom_using_cid
. It will be deleted in the future- Returns:
1 on successful deletion, return 0 on failure to delete.
-
std::pair<int, unsigned int> delete_side_chain(int imol, const std::string &chain_id, int res_no, const std::string &ins_code)
Delete side chain
- Parameters:
imol – is the model molecule index
chain_id – e.g. “A” for chain A
res_no – is the residue number, e.g. 12
ins_code – is the insertion code, e.g. “A”
- Returns:
1 on successful deletion, return 0 on failure to delete.
-
std::pair<int, unsigned int> delete_side_chain_using_cid(int imol, const std::string &cid)
Delete side chain using cid
- Parameters:
imol – is the model molecule index
cid – is the residue selection CID e.g “//A/15” (all the atoms in residue 15 of chain A)
- Returns:
1 on successful deletion, return 0 on failure to delete.
-
std::pair<int, unsigned int> delete_chain_using_cid(int imol, const std::string &cid)
Delete chain using chain cid
- Parameters:
imol – is the model molecule index
cid – is the selection CID e.g “//A” (chain A), “//*” (all chains)
- Returns:
1 on successful deletion, return 0 on failure to delete.
-
std::pair<int, unsigned int> delete_literal_using_cid(int imol, const std::string &cid)
Delete the atoms specified in the cid selection
- Parameters:
imol – is the model molecule index
cid – is the selection CID e.g “//A/15/OH” (atom OH in residue 15)
- Returns:
1 on successful deletion, return 0 on failure to delete.
-
std::pair<int, std::string> add_terminal_residue_directly(int imol, const std::string &chain_id, int res_no, const std::string &ins_code)
Add a residue onto the end of the chain by fitting to density
- Parameters:
imol – is the model molecule index
chain_id – e.g. “A” for chain A
res_no – is the residue number, e.g. 12
ins_code – is the insertion code, e.g. “A”
- Returns:
1 on success.
-
int add_terminal_residue_directly_using_cid(int imol, const std::string &cid)
Add a residue onto the end of the chain by fitting to density using cid
- Parameters:
imol – is the model molecule index
cid – is the selection CID e.g “//A/15/OH” (atom OH in residue 15)
-
int add_terminal_residue_directly_using_bucca_ml_growing_using_cid(int imol, const std::string &cid)
Add a residue onto the end of the chain by fitting to density using Buccaneer building and cid
- Parameters:
imol – is the model molecule index
cid – is the atom selection CID e.g “//A/15/OH” (atom OH in residue 15)
-
int add_terminal_residue_directly_using_bucca_ml_growing(int imol, const coot::residue_spec_t &spec)
Add a residue onto the end of the chain by fitting to density using Buccaneer building
- Parameters:
imol – is the model molecule index
spec – is the residue specifier, residue_spec_t(“A”, 10, “”)
-
inline void set_add_waters_water_to_protein_distance_lim_min(float d)
Parameter for
add_waters
function- Parameters:
d – is the min distance, default 2.4
-
inline void set_add_waters_water_to_protein_distance_lim_max(float d)
Parameter for
add_waters
function- Parameters:
d – is the max distance, default 3.4
-
inline void set_add_waters_variance_limit(float d)
Parameter for
add_waters
function- Parameters:
d – is the variance limit, default is 0.1
-
inline void set_add_waters_sigma_cutoff(float d)
Parameter for
add_waters
function- Parameters:
d – is the sigma cutoff, default is 1.75
-
int add_waters(int imol_model, int imol_map)
Add waters
- Parameters:
imol – is the model molecule index
imol_map – is the map molecule index
- Returns:
the number of waters added on a success, -1 on failure.
-
int add_hydrogen_atoms(int imol_model)
Add hydrogen atoms
- Parameters:
imol_model – is the model molecule index
- Returns:
1 on success, 0 on failure.
-
int delete_hydrogen_atoms(int imol_model)
Delete hydrogen atoms
- Parameters:
imol_model – is the model molecule index
- Returns:
1 on a successful deletion, 0 on failure.
-
int add_alternative_conformation(int imol_model, const std::string &cid)
Add an alternative conformation for the specified residue
- Parameters:
imol_model – is the model molecule index
cid – is the selection CID e.g “//A/15” (residue 15 in chain A)
- Returns:
1 on a successful addition, 0 on failure.
-
int fill_partial_residue(int imol, const std::string &chain_id, int res_no, const std::string &ins_code)
Fill the specified residue
- Parameters:
imol – is the model molecule index
chain_id – e.g. “A” for chain A
res_no – is the residue number, e.g. 12
ins_code – is the insertion code, e.g. “A”
- Returns:
1 on a successful fill, 0 on failure.
-
int fill_partial_residue_using_cid(int imol, const std::string &cid)
Fill the specified residue using cid
- Parameters:
imol – is the model molecule index
cid – is the selection CID e.g “//A/15” (residue 15 in chain A)
- Returns:
1 on a successful fill, 0 on failure.
-
int fill_partial_residues(int imol)
Fill all the the partially-filled residues in the molecule
- Parameters:
imol – is the model molecule index
- Returns:
1 on a successful fill, 0 on failure.
-
int flip_peptide(int imol, const coot::atom_spec_t &atom_spec, const std::string &alt_conf)
Flip peptide
- Parameters:
imol – is the model molecule index
atom_spec – is the atom specifier, atom_spec_t(“A”, 10, “”, “CA”, “”)
alt_conf – is the alternate conformation, e.g. “A” or “B”
- Returns:
1 on a successful flip
-
int flip_peptide_using_cid(int imol, const std::string &atom_cid, const std::string &alt_conf)
Flip peptide using cid
- Parameters:
imol – is the model molecule index
atom_cid – is the atom selection CID e.g “//A/15/OH” (atom OH in residue 15 of chain A)
alt_conf – is the alternate conformation, e.g. “A” or “B”
- Returns:
1 on a successful flip
-
void eigen_flip_ligand(int imol, const std::string &chain_id, int res_no, const std::string &ins_code)
Eigen-flip the specified ligand
- Parameters:
imol – is the model molecule index
chain_id – e.g. “A” for chain A
res_no – is the residue number, e.g. 12
ins_code – is the insertion code, e.g. “A”
-
void eigen_flip_ligand_using_cid(int imol, const std::string &residue_cid)
Eigen-flip ligand using cid
- Parameters:
imol – is the model molecule index
residue_cid – is the residue selection CID e.g “//A/15” (residue 15 of chain A)
-
int mutate(int imol, const std::string &cid, const std::string &new_residue_type)
Mutate residue
- Parameters:
imol – is the model molecule index
residue_cid – is the residue selection CID e.g “//A/15” (residue 15 of chain A)
new_residue_type – is the 3-letter code of the new residue, e.g. “TYR” for tyrosine
- Returns:
1 on a successful move, 0 on failure.
-
int side_chain_180(int imol, const std::string &atom_cid)
Rotate last chi angle of the side chain by 180 degrees
- Parameters:
imol – is the model molecule index
atom_cid – is the atom selection CID e.g “//A/15/OH” (atom OH of residue 15 of chain A)
- Returns:
1 on a successful move, 0 on failure.
-
std::string jed_flip(int imol, const std::string &atom_cid, bool invert_selection)
JED-Flip the ligand (or residue) at the specified atom
- Parameters:
imol – is the model molecule index
atom_cid – is the residue selection CID e.g “//A/15” (residue 15 of chain A)
invert_selection – is True if you want to use the larger fragment
- Returns:
a non-blank message if there is a problem
-
int move_molecule_to_new_centre(int imol, float x, float y, float z)
Move the molecule to the given centre
- Parameters:
imol – is the model molecule index
x – is the x coordinate of the new centre of the screen
y – is the y coordinate of the new centre of the screen
z – is the z coordinate of the new centre of the screen
- Returns:
1 on a successful move, 0 on failure.
-
void multiply_residue_temperature_factors(int imol, const std::string &cid, float factor)
Interactive B-factor refinement
- Parameters:
imol – is the model molecule index
cid – is the selection CID e.g “//A/15” (residue 15 of chain A)
factor – might typically be 0.9 or 1.1
-
coot::Cartesian get_molecule_centre(int imol) const
Get molecule centre
- Parameters:
imol – is the model molecule index
- Returns:
the molecule centre
-
int copy_fragment_using_cid(int imol, const std::string &multi_cid)
Copy a fragment given the multi_cid selection string
- Parameters:
imol – is the model molecule index
multi_cids – is a “||”-separated list of residues CIDs, e.g. “//A/12-52||//A/14-15||/B/56-66”
- Returns:
the new molecule number (or -1 on no atoms selected)
-
int copy_fragment_for_refinement_using_cid(int imol, const std::string &multi_cid)
Copy a fragment given the multi_cid selection string for refinement
Use this in preference to
copy_fragment_using_cid
when copying a molecule fragment to make a molten zone for refinement. That is because this version quietly also copies the residues near the residues of the selection, so that those residues can be used for links and non-bonded contact restraints.- Parameters:
imol – is the model molecule index
multi_cids – is a “||”-separated list of residues CIDs, e.g. “//A/12-52||//A/14-15||/B/56-66”
- Returns:
the new molecule number (or -1 on no atoms selected)
-
int copy_fragment_using_residue_range(int imol, const std::string &chain_id, int res_no_start, int res_no_end)
Copy a residue-range fragment
- Parameters:
imol – is the model molecule index
chain_id – e.g. “A”
res_no_start – the starting residue number
res_no_ref_end – the ending residue number
- Returns:
the new molecule number (or -1 on no atoms selected)
-
int apply_transformation_to_atom_selection(int imol, const std::string &atoms_selection_cid, int n_atoms, float m00, float m01, float m02, float m10, float m11, float m12, float m20, float m21, float m22, float c0, float c1, float c2, float t0, float t1, float t2)
Apply transformation to atom selection in the given molecule
- Returns:
the number of atoms moved.
-
int new_positions_for_residue_atoms(int imol, const std::string &residue_cid, std::vector<coot::api::moved_atom_t> &moved_atoms)
Update the positions of the atoms in the residue
- Parameters:
imol – is the model molecule index
residue_cid – is the residue selection CID e.g “//A/15” (residue 15 of chain A)
moved_atoms – is a list of the atoms moved in the specified residue, e.g. moved_atom_t(“CA”, 1, 2, 3)
-
int new_positions_for_atoms_in_residues(int imol, const std::vector<coot::api::moved_residue_t> &moved_residues)
Update the positions of the atoms in the residues
- Parameters:
imol – is the model molecule index
moved_residue – is a list of the residues with moved atoms, e.g. moved_residue_t(“A”, 10, “”)
-
std::pair<int, std::vector<merge_molecule_results_info_t>> merge_molecules(int imol, const std::string &list_of_other_molecules)
Merge molecules
- Parameters:
imol – is the model molecule index
list_of_other_molecules – is a colon-separated list of molecules, e.g. “2:3:4”
- Returns:
the first is a flag set to 1 if a merge occurred (and 0 if it did not) the second is a vector of merge results, i.e. if you merged a ligand, what is the new residue spec of the ligand, and if you merged a (polymer) chain, what is the new chain-id of that chain.
-
int cis_trans_convert(int imol, const std::string &atom_cid)
Convert a cis peptide to a trans or vice versa
- Parameters:
imol – is the model molecule index
atom_cid – is the atom selection CID e.g “//A/15/OH” (atom OH residue 15 of chain A)
- Returns:
1 on a successful conversion.
-
void replace_residue(int imol, const std::string &residue_cid, const std::string &new_residue_type, int imol_enc)
Replace a residue
This has a different meaning of “replace” to replace_fragment(). In this function the atoms are not merely moved/”slotted in to place”, but the residue type is changed - new atoms are introduce and others are deleted (typically).
Change the type of a residue (for example, “TYR” to “CYS”) The algorithm will superpose the mainchain CA, C and N and try to set matching torsion to the angle that they were in the reference structure.
- Parameters:
imol – is the model molecule index
residue_cid – is the residue selection CID e.g “//A/15” (residue 15 of chain A)
new_residue_type – is the 3-letter code of the new residue, e.g “CYS”
imol_enc – is the molecule index for the residue type/compound_id
-
int replace_fragment(int imol_base, int imol_reference, const std::string &atom_selection)
Replace a fragment
- Parameters:
imol_base – is the base model index
imol_reference – is the reference model index
atom_selection – is the selection CID e.g “//A/15-17” (residue 15, 16 and 17 of chain A)
- Returns:
the success status
-
int rigid_body_fit(int imol, const std::string &multi_cid, int imol_map)
Rigid-body fitting
- Parameters:
imol – is the model molecule index
multi_cids – is a “||”-separated list of residues CIDs, e.g. “//A/12-52||//A/14-15||/B/56-66”
imol_map – is the map molecule index
-
int rotate_around_bond(int imol, const std::string &residue_cid, const std::string &atom_name_1, const std::string &atom_name_2, const std::string &atom_name_3, const std::string &atom_name_4, double torsion_angle)
Rotate atoms around torsion
the bond is presumed to be between atom-2 and atom-3. Atom-1 and atom-4 are used to define the absolute torsion angle.
- Returns:
status 1 if successful, 0 if not.
-
std::pair<int, std::string> change_chain_id(int imol, const std::string &from_chain_id, const std::string &to_chain_id, bool use_resno_range, int start_resno, int end_resno)
Change the chain ID
- Parameters:
imol – is the model molecule index
from_chain_id – e.g. “A”
to_chain_id – e.g. “C”
use_resno_range – use residue number range, typically True
start_resno – the starting residue number of the range
end_resno – the ending residue number of the range
- Returns:
-1 on a conflict, 1 on good, 0 on did nothing, return also an information/error message
-
void associate_sequence(int imol, const std::string &name_or_chain_id, const std::string &sequence)
Associate a sequence with a molecule
- Parameters:
imol – is the model molecule index
name_or_chain_id – e.g. “A”
sequence – is the model sequence
-
void assign_sequence(int imol_model, int imol_map)
Assign a sequence to a molecule
Often one might copy out a fragment from a more complete molecule (and then copy it back after the sequence has been added). This runs
backrub_rotamer()
on the newly assigned residues- Parameters:
imol – is the model molecule index
imol_map – is the map molecule index
Coordinates Refinement
-
int refine_residues_using_atom_cid(int imol, const std::string &cid, const std::string &mode, int n_cycles)
Refine the residues using cid
- Parameters:
imol – is the model molecule index
cid – is the selection CID e.g “//A/15” (residue 15 of chain A)
mode – is the mode of real space refinement e.g. “SINGLE”, “TRIPLE”, “QUINTUPLE”, “HEPTUPLE”, “SPHERE”, “BIG_SPHERE”, “CHAIN”, “ALL”
n_cycles – is the number of refinement cycles
- Returns:
a value of 1 if the refinement was performed and 0 if it was not.
-
int refine_residues(int imol, const std::string &chain_id, int res_no, const std::string &ins_code, const std::string &alt_conf, const std::string &mode, int n_cycles)
Refine the residues
- Parameters:
imol – is the model molecule index
chain_id – e.g. “A” for chain A
res_no – is the residue number, e.g. 12
ins_code – is the insertion code, e.g. “A”
alt_conf – is the alternate conformation, e.g. “A” or “B”
mode – is the mode of real space refinement e.g. “SINGLE”, “TRIPLE”, “QUINTUPLE”, “HEPTUPLE”, “SPHERE”, “BIG_SPHERE”, “CHAIN”, “ALL”
n_cycles – is the number of refinement cycles
- Returns:
a value of 1 if the refinement was performed and 0 if it was not.
-
int refine_residue_range(int imol, const std::string &chain_id, int res_no_start, int res_no_end, int n_cycles)
Refine residue range
- Parameters:
imol – is the model molecule index
chain_id – e.g. “A” for chain A
res_no_start – the starting residue number
res_no_ref_end – the ending residue number
n_cycles – is the number of refinement cycles
- Returns:
a value of 1 if the refinement was performed and 0 if it was not.
-
std::pair<int, coot::instanced_mesh_t> minimize_energy(int imol, const std::string &atom_selection_cid, int n_cycles, bool do_rama_plot_restraints, float rama_plot_weight, bool do_torsion_restraints, float torsion_weight, bool refinement_is_quiet)
Minimise/optimise the geometry of the specified residue(s)
The use of “energy” should not be taken literally here
- Parameters:
imol – is the model molecule index
atom_selection_cid – is the selection CID e.g. “//A/15” (residue 15 of chain A)
n_cycles – is the number of refinement cycles. If you pass n_cycles = 100 (or some such) then you can get the mesh for the partially optimized ligand/residues
do_rama_plot_restraints – is the flag for the usage of Ramachandran plot restraints
rama_plot_weight – is the flag to set the Ramachandran plot restraints weight
do_torsion_restraints – is the flag for the usage of torsion restraints
torsion_weight – is the flag to set the torsion restraints weight
refinement_is_quiet – is used to reduce the amount of diagnostic text written to the output
- Returns:
the success status 1 if the minimization was performed and 0 if it was not.
-
float minimize(int imol, const std::string &atom_selection_cid, int n_cycles, bool do_rama_plot_restraints, float rama_plot_weight, bool do_torsion_restraints, float torsion_weight, bool refinement_is_quiet)
- Parameters:
imol – is the model molecule index
atom_selection_cid – is the selection CID e.g. “//A/15” (residue 15 of chain A)
n_cycles – is the number of refinement cycles. If you pass n_cycles = 100 (or some such) then you can get the mesh for the partially optimized ligand/residues
do_rama_plot_restraints – is the flag for the usage of Ramachandran plot restraints
rama_plot_weight – is the flag to set the Ramachandran plot restraints weight
do_torsion_restraints – is the flag for the usage of torsion restraints
torsion_weight – is the flag to set the torsion restraints weight
refinement_is_quiet – is used to reduce the amount of diagnostic text written to the output
- Returns:
the function value at termination
-
void fix_atom_selection_during_refinement(int imol, const std::string &atom_selection_cid)
Fix atoms during refinement
Does nothing at the moment
- Parameters:
imol – is the model molecule index
atom_selection_cid – is the selection CID e.g “//A/15/OH” (atom OH of residue 15 of chain A)
-
void add_target_position_restraint(int imol, const std::string &atom_cid, float pos_x, float pos_y, float pos_z)
Add or update restraint (if it has a pull restraint already)
- Parameters:
imol – is the model molecule index
atom_cid – is the selection CID e.g “//A/15/OH” (atom OH of residue 15 of chain A)
pos_x – is the x coordinate of the target position of the specified atom
pos_y – is the y coordinate of the target position of the specified atom
pos_z – is the z coordinate of the target position of the specified atom
-
void clear_target_position_restraint(int imol, const std::string &atom_cid)
Clear target_position restraint
- Parameters:
imol – is the model molecule index
atom_cid – is the selection CID e.g “//A/15/OH” (atom OH of residue 15 of chain A)
-
void turn_off_when_close_target_position_restraint(int imol)
Clear target_position restraint if it is (or they are) close to their target position
- Parameters:
imol – is the model molecule index
-
inline void set_use_rama_plot_restraints(bool state)
Turn on or off rama restraints
- Parameters:
state – is True to mean that it is enabled
-
inline bool get_use_rama_plot_restraints() const
Get the state of the rama plot restraints usage in refinement
- Returns:
the state
-
inline void set_rama_plot_restraints_weight(float f)
Set the Ramachandran plot restraints weight
- Parameters:
f – is the weight to set, default 1.0
-
inline float get_rama_plot_restraints_weight() const
Get the Ramachandran plot restraints weight
- Returns:
the Ramachandran plot restraints weight
-
inline void set_use_torsion_restraints(bool state)
Turn on or off torsion restraints
- Parameters:
state – is True to mean that it is enabled
-
inline bool get_use_torsion_restraints() const
Get the state of the rama plot restraints usage in refinement
- Returns:
the state
-
inline void set_torsion_restraints_weight(float f)
Set the torsiont restraints weight
- Parameters:
f – is the weight to set, default value is 1.0
-
inline float get_torsion_restraints_weight() const
Get the torsion restraints weight
- Returns:
the torsion restraints weight
-
void init_refinement_of_molecule_as_fragment_based_on_reference(int imol_frag, int imol_ref, int imol_map)
Initialise the refinement of (all of) molecule imol_frag
- Parameters:
imol_frag – is the model molecule index of the fragment
imol_ref – is the model molecule index of the reference
imol_map – is the map molecule index
-
std::pair<int, coot::instanced_mesh_t> refine(int imol, int n_cycles)
Run some cycles of refinement and return a mesh
That way we can see the molecule animate as it refines
- Parameters:
imol – is the model molecule index
n_cycles – is the number of refinement cycles
- Returns:
a pair: the first of which is the status of the refinement: GSL_CONTINUE, GSL_SUCCESS, GSL_ENOPROG (no progress). i.e. don’t call this function again unless the status is GSL_CONTINUE (-2); The second is a coot::instanced_mesh_t
-
coot::instanced_mesh_t add_target_position_restraint_and_refine(int imol, const std::string &atom_cid, float pos_x, float pos_y, float pos_z, int n_cycles)
Create a new position for the given atom and create a new bonds mesh based on that
This is currently “heavyweight” as the bonds mesh is calculated from scratch (it is not (yet) merely a distortion of an internally-stored mesh).
- Parameters:
imol – is the model molecule index
atom_cid – is the selection CID e.g “//A/15/OH” (atom OH of residue 15 of chain A)
pos_x – is the x coordinate of the target position of the specified atom
pos_y – is the y coordinate of the target position of the specified atom
pos_z – is the z coordinate of the target position of the specified atom
n_cycles – specifies the number of refinement cyles to run after the target position of the atom has been applied. If n_cycles is -1 then, no cycles are done and the mesh is bonds merely calculated.
- Returns:
a
instanced_mesh_t
-
void clear_target_position_restraints(int imol)
Clear any and all drag-atom target position restraints
- Parameters:
imol – is the model molecule index
-
void clear_refinement(int imol)
Call this after molecule refinement has finished (say when the molecule molecule is accepted into the original molecule)
- Parameters:
imol – is the model molecule index
-
inline void set_refinement_is_verbose(bool state)
For debugging the refinement - write out some diagnositics - some might be useful.
API change 20240226 - this function now takes a boolean argument
- Parameters:
state – is True to mean that it is enabled
-
inline void set_refinement_geman_mcclure_alpha(float a)
Set the refinement Geman-McClure alpha
- Parameters:
a – is the Geman-McClure alpha, e.g. 0.01
-
inline float get_geman_mcclure_alpha() const
Get the refinement Geman-McClure alpha
- Returns:
the Geman-McClure alpha
-
int generate_self_restraints(int imol, float local_dist_max)
Generate GM self restraints for the whole molecule
- Parameters:
imol – is the model molecule index
local_dist_max – is the maximum distance, e.g. 4.6
-
void generate_chain_self_restraints(int imol, float local_dist_max, const std::string &chain_id)
Generate GM self restraints for the given chain
- Parameters:
imol – is the model molecule index
local_dist_max – is the maximum distance, e.g. 4.6
chain_id – e.g. “A” for chain A
-
void generate_local_self_restraints(int imol, float local_dist_max, const std::string &residue_cids)
Generate GM self restraints for the given residues
- Parameters:
imol – is the model molecule index
local_dist_max – is the maximum distance, e.g. 4.6
residue_cids – is a “||”-separated list of residues, e.g. “//A/12||//A/14||//B/56”
-
void add_parallel_plane_restraint(int imol, const std::string &residue_cid_1, const std::string &residue_cid_2)
Generate parallel plane restraints (for RNA and DNA)
- Parameters:
imol – is the model molecule index
residue_cid_1 – is the selection CID e.g “//A/15” (residue 15 of chain A)
residue_cid_2 – is the selection CID e.g “//A/17” (residue 17 of chain A)
-
coot::instanced_mesh_t get_extra_restraints_mesh(int imol, int mode)
Get the mesh for extra restraints
- Parameters:
imol – is the model molecule index
mode – is currently unused
-
void read_extra_restraints(int imol, const std::string &file_name)
Read extra restraints (e.g. from ProSMART)
- Parameters:
imol – is the model molecule index
-
void clear_extra_restraints(int imol)
Clear the extra restraints
- Parameters:
imol – is the model molecule index
Coordinates Validation
-
coot::simple_mesh_t get_rotamer_dodecs(int imol)
Get the rotamer dodecs for the model
- Parameters:
imol – is the model molecule index
- Returns:
a
simple_mesh_t
-
coot::instanced_mesh_t get_rotamer_dodecs_instanced(int imol)
Get the rotamer dodecs for the model instanced
- Parameters:
imol – is the model molecule index
- Returns:
an
instanced_mesh_t
-
coot::simple_mesh_t get_ramachandran_validation_markup_mesh(int imol) const
Get the Ramachandran validation markup mesh
20221126-PE: the function was renamed from ramachandran_validation_markup_mesh().
- Parameters:
imol – is the model molecule index
- Returns:
a
simple_mesh_t
-
std::vector<coot::phi_psi_prob_t> ramachandran_validation(int imol) const
Get the data for Ramachandran validation, which importantly contains probability information
- Parameters:
imol – is the model molecule index
- Returns:
a vector/list of
phi_psi_prob_t
-
coot::instanced_mesh_t contact_dots_for_ligand(int imol, const std::string &cid, unsigned int smoothness_factor) const
Contact dots for ligand
- Parameters:
imol – is the model molecule index
cid – is the selection CID e.g “//A/15” (residue 15 of chain A)
smoothness_factor – is 1, 2 or 3 (3 is the most smooth). Recently added (20230202)
- Returns:
the instanced mesh for the specified ligand
-
coot::instanced_mesh_t all_molecule_contact_dots(int imol, unsigned int smoothness_factor) const
Contact dots for the whole molecule/model
- Parameters:
imol – is the model molecule index
smoothness_factor – is 1, 2 or 3 (3 is the most smooth). Recently added (20230202)
- Returns:
the instanced mesh for the specified molecule.
-
coot::simple::molecule_t get_simple_molecule(int imol, const std::string &residue_cid, bool draw_hydrogen_atoms_flag)
Get a simple molecule
- Parameters:
imol – is the model molecule index
residue_cid – is the selection CID e.g “//A/15” (residue 15 of chain A)
draw_hydrogen_atoms_flag – is the flag for drawing H atoms
- Returns:
a simple::molecule_t for the specified residue.
-
generic_3d_lines_bonds_box_t make_exportable_environment_bond_box(int imol, coot::residue_spec_t &spec)
- Parameters:
spec – is the residue specifier, e.g. residue_spec_t(“A”, 10, “”)
- Returns:
a vector of lines for non-bonded contacts and hydrogen bonds
-
std::vector<moorhen::h_bond> get_h_bonds(int imol, const std::string &cid_str, bool mcdonald_and_thornton_mode) const
Get hydrogen bonds
- Parameters:
imol – is the model molecule index
cid – is the selection CID e.g “//A/15” (residue 15 of chain A)
mcdonald_and_thornton_mode – turns on the McDonald & Thornton algorithm - using explicit hydrogen atoms
- Returns:
a vector of hydrogen bonds around the specified residue (typically a ligand)
-
coot::simple_mesh_t get_mesh_for_ligand_validation_vs_dictionary(int imol, const std::string &ligand_cid)
Get the mesh for ligand validation vs dictionary, coloured by badness
Greater then 3 standard deviations is fully red. Less than 0.5 standard deviations is fully green
- Parameters:
imol – is the model molecule index
ligand_cid – is the ligand selection CID e.g “//A/15” (ligand 15 of chain A)
-
std::vector<coot::geometry_distortion_info_container_t> get_ligand_validation_vs_dictionary(int imol, const std::string &ligand_cid, bool include_non_bonded_contacts)
Ligand validation
- Parameters:
imol – is the model molecule index
ligand_cid – is the ligand selection CID e.g “//A/15” (ligand 15 of chain A)
include_non_bonded_contacts – is the flag to include non bonded contacts
- Returns:
a vector/list of interesting geometry
-
std::pair<int, double> get_ligand_distortion(int imol, const std::string &ligand_cid, bool include_non_bonded_contacts)
Get ligand distortion
a more simple interface to the above
- Returns:
a pair: the first is the status (1 for OK, 0 for fail)
-
bool match_ligand_torsions(int imol_ligand, int imol_ref, const std::string &chain_id_ref, int resno_ref)
Match ligand torsions
- Parameters:
imol_ligand – is the ligand molecule index
imol_ref – is the reference model molecule index
chain_id_ref – is the reference chain, e.g. “A”
resno_ref – is the reference residue number, e.g. 12
- Returns:
the success status
-
bool match_ligand_position(int imol_ligand, int imol_ref, const std::string &chain_id_ref, int resno_ref)
Match ligand positions
i.e. do a least-squares superposition of the atoms that match in the graphs of the two specified ligands - typically one would use this function after matching ligand torsions.
- Parameters:
imol_ligand – is the ligand molecule index
imol_ref – is the reference model molecule index
chain_id_ref – is the reference chain, e.g. “A”
resno_ref – is the reference residue number, e.g. 12
- Returns:
the success status
-
bool match_ligand_torsions_and_position(int imol_ligand, int imol_ref, const std::string &chain_id_ref, int resno_ref)
Match ligand torsions and positions
- Parameters:
imol_ligand – is the ligand molecule index
imol_ref – is the reference model molecule index
chain_id_ref – is the reference chain, e.g. “A”
resno_ref – is the reference residue number, e.g. 12
- Returns:
the success status.
-
bool match_ligand_torsions_and_position_using_cid(int imol_ligand, int imol_ref, const std::string &cid)
Match ligand torsions and positions using cid
- Parameters:
imol_ligand – is the ligand molecule index
imol_ref – is the reference model molecule index
cid – is the selection CID e.g “//A/15” (residue 15 of chain A)
-
coot::atom_overlaps_dots_container_t get_overlap_dots(int imol)
- Parameters:
imol – is the model molecule index
-
coot::atom_overlaps_dots_container_t get_overlap_dots_for_ligand(int imol, const std::string &cid_ligand)
This function not const because it can dynamically add dictionaries
- Parameters:
imol – is the model molecule index
cid_ligand – is the ligand selection CID e.g “//A/15” (ligand 15 of chain A)
Coordinates and Map Validation
-
coot::validation_information_t density_fit_analysis(int imol_model, int imol_map) const
Density fit validation information
- Parameters:
imol_model – is the model molecule index
imol_map – is the map molecule index
- Returns:
an object
validation_information_t
-
double get_sum_density_for_atoms_in_residue(int imol, const std::string &cid, const std::vector<std::string> &atom_names, int imol_map)
- Returns:
the sum of the density of the given atoms in the specified CID return -1001 on failure to find the residue or any atoms in the residue or if imol_map is not a map
-
coot::validation_information_t density_correlation_analysis(int imol_model, int imol_map) const
Get the density correlation validation information
- Parameters:
imol_model – is the model molecule index
imol_map – is the map molecule index
- Returns:
an object
validation_information_t
-
coot::validation_information_t rotamer_analysis(int imol_model) const
Get the rotamer validation information
- Parameters:
imol_model – is the model molecule index
- Returns:
an object
validation_information_t
-
coot::validation_information_t ramachandran_analysis(int imol_model) const
Get the ramachandran validation information (formatted for a graph, not 3D)
- Parameters:
imol_model – is the model molecule index
- Returns:
an object
validation_information_t
-
coot::validation_information_t ramachandran_analysis_for_chain(int imol_model, const std::string &chain_id) const
Get the ramachandran validation information (formatted for a graph, not 3D) for a given chain in a given molecule
This function does not exist yet (20230127-PE)
- Parameters:
imol_model – is the model molecule index
chain_id – e.g. “A”
- Returns:
an object
validation_information_t
-
coot::validation_information_t peptide_omega_analysis(int imol_model) const
Peptide omega validation information
- Parameters:
imol_model – is the model molecule index
- Returns:
an object
validation_information_t
-
float get_median_temperature_factor(int imol) const
Get the median temperature factor for the model
- Parameters:
imol – is the model molecule index
- Returns:
a negative number on failure
-
std::vector<coot::molecule_t::interesting_place_t> get_interesting_places(int imol, const std::string &mode) const
Get interesting places
This function does not work yet
- Returns:
a vector/list of
validation_information_t
-
std::vector<coot::molecule_t::interesting_place_t> difference_map_peaks(int imol_map, int imol_protein, float n_rmsd) const
Get difference map peaks
- Parameters:
imol_map – is the map molecule index
imol_protein – is the model molecule index
n_rmsd – number of sd, e.g. 4.8
- Returns:
a vector/list of
validation_information_t
-
std::vector<coot::molecule_t::interesting_place_t> pepflips_using_difference_map(int imol_coords, int imol_difference_map, float n_sigma) const
Get pepflips based on the difference map
- Parameters:
imol_coords – is the model molecule index
imol_difference_map – is the difference map molecule index
n_sigma – number of sd, e.g. 4.8
- Returns:
a vector/list of
validation_information_t
-
std::vector<coot::molecule_t::interesting_place_t> unmodelled_blobs(int imol_model, int imol_map) const
Unmodelled blobs
- Parameters:
imol_model – is the model molecule index
imol_map – is the map molecule index
- Returns:
a vector/list of
validation_information_t
-
std::vector<coot::atom_spec_t> find_water_baddies(int imol_model, int imol_map, float b_factor_lim, float outlier_sigma_level, float min_dist, float max_dist, bool ignore_part_occ_contact_flag, bool ignore_zero_occ_flag)
Check waters, using implicit logical OR
- Parameters:
imol_model – is the model molecule index
imol_map – is the map molecule index
b_factor_lim – typical value is 60.0
outlier_sigma_level – typical value is 0.8
min_dist – typical value is 2.3
max_dist – typical value is 3.5
ignore_part_occ_contact_flag – typical value is False
ignore_zero_occ_flag – typical value is False
- Returns:
a vector/list of atom specifiers
-
std::vector<std::pair<double, double>> fourier_shell_correlation(int imol_map_1, int imol_map_2) const
Fourier Shell Correlation (FSC) between maps
- Parameters:
imol_map_1 – is the first map molecule index
imol_map_2 – is the second map molecule index
- Returns:
a vector/list or pairs of graph points (resolution, correlation). The resolution is in inverse Angstroms squared. An empty list is returned on failure
-
int make_power_scaled_map(int imol_ref, int imol_map_for_scaling)
Make a FSC-scaled map
- Parameters:
imol_ref – is the reference map molecule index
imol_map_for_scaling – is the second map molecule index
- Returns:
the molecule index of the new map
-
coot::validation_information_t get_q_score(int imol_model, int imol_map) const
Get the Q Score (Pintilie et al.)
- Parameters:
imol_model – is the model molecule index
imol_map – is the map molecule index
- Returns:
a
validation_information_t
object
-
coot::validation_information_t get_q_score_for_cid(int imol_model, const std::string &cid, int imol_map) const
Get the Pintile et al. Q Score for a particular residue (typically a ligand)
- Parameters:
cid – If the
cid
matches more than one residue the score will be returned for all of the residues covered in thecid
. Typically, of course thecid
will be something like “//A/301”.- Returns:
a
validation_information_t
object
Rail Points!
-
int calculate_new_rail_points()
Calling this adds to the rail_points history. Make this pairs when we add model scoring
- Returns:
the new rail points (since last modification)
-
int rail_points_total() const
The total rail points
- Returns:
the sum of all rail points accumulated since the maps were connected.
Updating Maps
-
void associate_data_mtz_file_with_map(int imol, const std::string &data_mtz_file_name, const std::string &f_col, const std::string &sigf_col, const std::string &free_r_col)
Associate a data mtz file with a molecule
This function is called before calling “connect_updating_maps()”
- Parameters:
imol – is the map molecule index
data_mtz_file_name – is the name of the mtz file
f_col – is the F column, e.g. “FOBS”
sigf_col – e.g. “SIGFOBS”
free_r_col – e.g. “RFREE”
-
int connect_updating_maps(int imol_model, int imol_with_data_info_attached, int imol_map_2fofc, int imol_map_fofc)
Connect updating maps
Reset the rail_points (calls “reset_the_rail_points()”), updates the maps (using internal/clipper SFC). Update your contour lines meshes after calling this function.
- Parameters:
imol_model – is the model molecule index
imol_with_data_info_attached – is the map index with the data have been attached by the previous function (associate_data_mtz_file_with_map)
imol_map_2fofc – is the map molecule index of the 2FO-FC map
imol_map_fofc – is the map molecule index of the FO-FC map
- Returns:
1 if the connection was successful
-
void sfcalc_genmap(int imol_model, int imol_map_with_data_attached, int imol_updating_difference_map)
Calculate SF and re-generate maps
This is a low-level function - generally one would use the updating maps method rather than this
- Parameters:
imol_model – is the model molecule index
imol_map_with_data_attached – is the map index with the data have been attached by the previous function (associate_data_mtz_file_with_map)
imol_updating_difference_map – is the index of the difference map that you want to update when the model updates
-
coot::util::sfcalc_genmap_stats_t sfcalc_genmaps_using_bulk_solvent(int imol_model, int imol_2fofc_map, int imol_updating_difference_map, int imol_map_with_data_attached)
Calculate SF and re-generate maps using bulk solvent
This is a low-level function. Call this function after connecting maps for updating maps to set the initial R-factor and store the initial map flatness.
- Parameters:
imol_model – is the model molecule index
imol_2fofc_map – is the map molecule index of the 2FO-FC map
imol_updating_difference_map – is the index of the difference map that you want to update when the model updates
imol_map_with_data_attached – is the map index with the data have been attached by the previous function (associate_data_mtz_file_with_map)
- Returns:
a class of interesting statistics. On failure to calculate SFs and generate the maps the returned r_factor in the returned stats will be set to -1.
-
bool shift_field_b_factor_refinement(int imol, int imol_with_data_attached)
Shift-field B-factor refinement
This function presumes that the Fobs, sigFobs and RFree data have been filled in the
imol_map_with_data_attached
molecule- Parameters:
imol – is the model molecule index
imol_map_with_data_attached – is the map index with the data have been attached by the previous function (associate_data_mtz_file_with_map)
- Returns:
success status
-
float get_density_at_position(int imol_map, float x, float y, float z) const
Get density at position
- Parameters:
imol_map – is the map molecule index
x – is the x coordinate of the target position
y – is the y coordinate of the target position
z – is the z coordinate of the target position
- Returns:
the density value
-
std::vector<std::pair<clipper::Coord_orth, float>> get_diff_diff_map_peaks(int imol_diff_map, float screen_centre_x, float screen_centre_y, float screen_centre_z) const
- Parameters:
imol_diff_map – is the map molecule index of the difference map
screen_centre_x – is the position x of the center of the screen
screen_centre_y – is the position y of the center of the screen
screen_centre_z – is the position z of the center of the screen
- Returns:
a vector of the position where the difference map has been flattened. The associated float value is the amount that the map has been flattened.
-
std::string get_data_set_file_name(int imol) const
Get the stored data set file name
- Parameters:
imol – is the model molecule index
Go to Blob
-
std::pair<bool, clipper::Coord_orth> go_to_blob(float x1, float y1, float z1, float x2, float y2, float z2, float contour_level)
Given a point on the front clipping plane (x1, y1, z1) and a point on the back clipping plane (x2, y2, z2) this function searches imol_refinement_map (if set) to find a the centre of a blob above the contour level. Blobs at the “front” are selected in preference to blobs at the back. If no blob is found, then the first of the pair is false. If it is found, then the second is the centre of the blob.
In future, this function should/will be provided with a list of displayed maps and their contour levels - but for now, it uses (only) imol_refinement_map.
Use a string to pass the map information (map index and contour level), something like “0 0.45:1 1.2:2 0.1”
- Parameters:
x1 – is the x point of the front clipping plane
y1 – is the y point of the front clipping plane
z1 – is the z point of the front clipping plane
x2 – is the x point of the back clipping plane
y2 – is the y point of the back clipping plane
z2 – is the z point of the back clipping plane
Ligand Functions
-
std::vector<int> fit_ligand_right_here(int imol_protein, int imol_map, int imol_ligand, float x, float y, float z, float n_rmsd, bool use_conformers, unsigned int n_conformers)
Fit the ligand at specified position. You can expect this to take about 20 seconds. For trivial (i.e non-flexible) ligands you should instead use the jiggle-fit algorithm, which takes a fraction of a second. (That is the algorithm used for “Add Other Solvent Molecules” in Coot.)
- Parameters:
imol_protein – is the model molecule index
imol_map – is the map molecule index
imol_ligand – is the ligand molecule index
x – is the x position of the blob
y – is the y position of the blob
z – is the z position of the blob
n_rmsd – number of sd, e.g. 4.8
use_conformers – is True for flexible ligands
n_conformers – set the number of conformers
- Returns:
a vector/list of indices of molecules for the best fitting ligands to this blob.
-
std::vector<fit_ligand_info_t> fit_ligand(int imol_protein, int imol_map, int imol_ligand, float n_rmsd, bool use_conformers, unsigned int n_conformers)
Fit ligand
- Parameters:
imol_protein – is the model molecule index
imol_map – is the map molecule index
imol_ligand – is the ligand molecule index
n_rmsd – number of sd, e.g. 4.8
use_conformers – is True for flexible ligands
n_conformers – set the number of conformers
- Returns:
a vector/list of interesting information about the fitted ligands
-
std::vector<fit_ligand_info_t> fit_ligand_multi_ligand(int imol_protein, int imol_map, const std::string &multi_ligand_molecule_number_list, float n_rmsd, bool use_conformers, unsigned int n_conformers)
Fit multiple ligands (place-holder)
- Parameters:
imol_protein – is the model molecule index
imol_map – is the map molecule index
multi_ligand_molecule_number_list – is a colon-separated list of molecules, e.g. “2:3:4”
n_rmsd – is number of sd, e.g. 4.8
use_conformers – is True for flexible ligands
n_conformers – set the number of conformers
- Returns:
an empty vector (at the moment)
-
float fit_to_map_by_random_jiggle(int imol, const coot::residue_spec_t &res_spec, int n_trials, float translation_scale_factor)
Jiggle-Fit Ligand
- Parameters:
imol – is the model molecule index
res_spec – is the residue specifier, e.g. residue_spec_t(“A”, 10, “”)
n_trials – is the number of trials, if n_trials is 0, then a sensible default value will be used.
translation_scale_factor – if is negative then a sensible default value will be used.
- Returns:
a value less than -99.9 on failure to fit.
-
float fit_to_map_by_random_jiggle_using_cid(int imol, const std::string &cid, int n_trials, float translation_scale_factor)
Jiggle-Fit Ligand using cid
- Parameters:
imol – is the model molecule index
cid – is the selection CID, e.g “//A/15” (ligand 15 of chain A)
n_trials – is the number of trials, if n_trials is 0, then a sensible default value will be used.
translation_scale_factor – if is negative then a sensible default value will be used.
- Returns:
a value less than -99.9 on failure to fit.
-
float fit_to_map_by_random_jiggle_with_blur_using_cid(int imol, int imol_map, const std::string &cid, float b_factor, int n_trials, float translation_scale_factor)
Jiggle-Fit an atom selection, typically a whole molecule or a chain
- Parameters:
imol – is the model molecule index
cid – is the selection CID, e.g. “//A” (chain A)
b_factor – e.g. 100.0
n_trials – is the number of trials, if n_trials is 0, then a sensible default value will be used.
translation_scale_factor – if is negative then a sensible default value will be used.
- Returns:
a value less than -99.9 on failure to fit.
-
std::string get_svg_for_residue_type(int imol, const std::string &comp_id, bool use_rdkit_svg, bool dark_background_flag)
Get svg for residue type
It won’t work unless the dictionary for that ligand has been imported. The output renderings are not very good at the moment.
- Parameters:
imol – is the model molecule index, except for unusual cases, it will be IMOL_ENC_ANY (-999999)
comp_id – is the 3-letter code for the residue/ligand, e.g. “ALA” for alanine
use_rdkit_svg – is the flag for using the rdkit svg renderer
dark_background_flag – returns a representation suitable for rendering on a dark background
- Returns:
the string for the SVG representation.
-
int add_compound(int imol, const std::string &tlc, int imol_dict, int imol_map, float x, float y, float z)
This function is for adding compounds/molecules like buffer agents and precipitants or anions and cations. e.g. those ligands that can be positioned without need for internal torsion angle manipulation.
- Parameters:
imol – is the model molecule index
tlc – is the 3-letter-code/compound-id
imol_dict – is the molecule to which the ligand is attached (if any). Typically this will be IMOL_ENC_ANY (-999999).
imol_map – is the map molecule index
x – is the x position
y – is the y position
z – is the z position
- Returns:
the success status, 1 for good, 0 for not good.
-
std::vector<coot::residue_spec_t> get_non_standard_residues_in_molecule(int imol) const
Get non-standard residues in a model
- Parameters:
imol – is the model molecule index
- Returns:
a vector/list of residue specifiers - the residue name is encoded in the
string_user_data
data item of the residue specifier
-
std::vector<int> get_dictionary_conformers(const std::string &comp_id, int imol_enc, bool remove_internal_clash_conformers)
Get the conformers that can be generated by variation around rotatable bonds as described in the dictionary.
Torsions that are marked as “const” are excluded from the variation, as are pyranose ring torsions and torsions that rotate hydrogen atoms
- Parameters:
comp_id – is the 3-letter code for the residue/ligand, e.g. “ALA” for alanine
imol_enc – is the molecule index for the residue type/compound_id
remove_internal_clash_conformers – is the flag for removing internal clash
- Returns:
a vector/list of indices of the new molecules
-
texture_as_floats_t get_map_section_texture(int imol, int section_id, int axis, float data_value_for_bottom, float data_value_for_top) const
The new arguments,
data_value_for_bottom
,data_value_for_top
should be pre-calculated (don’t calculate them for every call to this function).- Parameters:
imol – is the map molecule index
section_id – e.g. 2
axis – e.g. 0 for X-axis, 1 for Y-axis, 2 for Z-axis
- Returns:
a texture_as_floats_t object for the given section. On failure, the image_data vector is empty.
-
int get_number_of_map_sections(int imol_map, int axis_id) const
- Parameters:
imol_map – is the map molecule index
axis_id – is 0 for X-axis, 1 for Y-axis, 2 for Z-axis
- Returns:
the number of sections in the map along the given axis, -1 on failure.
Other Features
-
coot::simple_mesh_t make_mesh_from_gltf_file(const std::string &file_name)
- Parameters:
file_name – is the glTF file name
- Returns:
a
simple_mesh_t
from the given file.
-
coot::simple_mesh_t get_octahemisphere(unsigned int n_divisions) const
Get octahemisphere
- Parameters:
n_divisions – is a number divisible by 2, at least 4 (typically 16)
- Returns:
a unit-vector end-cap octohemisphere mesh
-
std::string pae_png(const std::string &pae_file_name) const
Predicted alignment error (AlphaFold)
- Returns:
a string of a png
Functions for Blender Interface
-
void make_mesh_for_map_contours_for_blender(int imol, float x, float y, float z, float level, float radius)
Function for Blender interface.
-
void make_mesh_for_bonds_for_blender(int imol, const std::string &mode, bool against_a_dark_background, float bond_width, float atom_radius_to_bond_width_ratio, int smoothness_factor)
Function for Blender interface.
-
void make_mesh_for_molecular_representation_for_blender(int imol, const std::string &cid, const std::string &colour_scheme, const std::string &style, int secondary_structure_usage_flag)
Function for Blender interface
Make an (internal) mesh
This function doesn’t return a value, instead it stores a
blender_mesh_t
blender_mesh in this model. One then (shortly later) uses get_triangles_for_blender(imol) (etc) to import this mesh into blender.@modifies internal state to fill the internal
blender_mesh
object
-
void make_mesh_for_gaussian_surface_for_blender(int imol, float sigma, float contour_level, float box_radius, float grid_scale, float b_factor)
Function for Blender interface.
-
void make_mesh_for_goodsell_style_for_blender(int imol, float colour_wheel_rotation_step, float saturation, float goodselliness)
blender
Function for Blender interface
-
std::vector<float> get_colour_table_for_blender(int imol)
Function for Blender interface.
-
std::vector<float> get_vertices_for_blender(int imol)
Function for Blender interface.
-
std::vector<int> get_triangles_for_blender(int imol)
Function for Blender interface.
Public Functions
-
inline explicit molecules_container_t(bool verbose = true)
the one and only constructor
Real space refinement.
the
n_cycles
parameter allows partial refinement - so for an animated representation one would call this with a small number (10, 20, 100?) and call it again if the refine status is still yet to reach completion GSL_CONTINUE (-2). And then make a call to get the bonds mesh (or other molecular representation). If n_cycles is negative, this means “refine to completion.”- Returns:
success/progress status
-
~molecules_container_t()
-
inline void set_use_gemmi(bool state)
Set the state of using GEMMI for coordinates parsing
- Parameters:
state – is True to mean that it is enabled. The default is True.
-
inline bool get_use_gemmi()
Get the state of using GEMMI for coordinates parsing.
Public Members
-
int imol_refinement_map
the refinement map - direct access. When refinement is performed, this is the map that will be used. Many (indeed most) of thesee functions explicity take a map. If the map is not known by the calling function then this map can be used as the map molecule index
-
int imol_difference_map
the difference map - direct access
I am not sure that this is needed - or will ever be.
-
bool use_gemmi
Private Types
Private Functions
-
std::string adjust_refinement_residue_name(const std::string &resname) const
- Parameters:
resname – is the 3 letter code for the residue, e.g. “ALA” for alanine
- Returns:
the state of having found restraints
Private Members
-
std::vector<coot::molecule_t> molecules
-
ramachandrans_container_t ramachandrans_container
-
bool draw_missing_residue_loops_flag
-
std::vector<rail_points_t> rail_point_history
-
float map_weight
-
float geman_mcclure_alpha
-
bool use_rama_plot_restraints
-
float rama_plot_restraints_weight
-
bool use_torsion_restraints
-
float torsion_restraints_weight
-
ctpl::thread_pool thread_pool
-
bool show_timings
-
bool continue_threaded_refinement_loop
-
bool refinement_is_quiet
-
int cif_dictionary_read_number
-
bool refinement_immediate_replacement_flag = true
-
int imol_moving_atoms
-
short int moving_atoms_asc_type
Private Static Attributes
-
static std::atomic<bool> on_going_updating_map_lock
-
static clipper::Xmap<float> *dummy_xmap
-
class auto_read_mtz_info_t
- #include <molecules-container.hh>
class for the information about columns extracted from auto-reading the given mtz file
Public Functions
-
inline auto_read_mtz_info_t()
-
inline auto_read_mtz_info_t(int index, const std::string &F_in, const std::string &phi_in)
-
inline void set_fobs_sigfobs(const std::string &f, const std::string &s)
Public Members
-
int idx
molecule index
-
std::string F
F column.
-
std::string phi
phi column
-
std::string w
weights column
-
bool weights_used
flag for weights usage
-
std::string F_obs
F_obs column. There were not available if the return value is empty.
-
std::string sigF_obs
sigF_obs column
-
std::string Rfree
R-Free column. There were not available if the return value is empty.
-
inline auto_read_mtz_info_t()
-
class fit_ligand_info_t
- #include <molecules-container.hh>
Ligand Fitting
- Return:
a vector of indices of molecules for the best fitting ligands each of the “possible ligand” blobs.
-
class ltj_stats_t
Public Functions
-
inline ltj_stats_t()
-
inline void update_timer()
This function is called by the long-term job, udating the timer and count.
-
inline double time_difference()
This function is called by the interrogation function - and is to help the uer know how the job is going.
20230127-PE A nice graph of the change of the function value seems like a good idea
-
inline ltj_stats_t()
-
class r_factor_stats
-
class rail_points_t
Public Functions
-
inline explicit rail_points_t(float rmsd)
-
inline rail_points_t(float rmsd_diff_map_current, const rail_points_t &rail_points_prev)
-
inline int rail_points_delta(const rail_points_t &prev)
Public Members
-
int model_rail_points_delta
-
int map_rail_points_delta
-
float rmsd_of_difference_map
Public Static Functions
-
static inline int total(const std::vector<rail_points_t> &rail_point_history)
-
inline explicit rail_points_t(float rmsd)
-
bool make_backups_flag
Mesh Objects
Functions
-
std::ostream &operator<<(std::ostream &s, const g_triangle &t)
-
class g_triangle
- #include <g_triangle.hh>
g_triangle is the container for the vertices of a triangles (indexing into the vertices vector).
Subclassed by g_triangle_with_colour_index
Public Functions
-
inline g_triangle(const unsigned int &a0, const unsigned int &a1, const unsigned int &a2)
constructor
-
inline g_triangle()
-
inline unsigned int &operator[](const unsigned int &i)
int colour_index;
-
inline const unsigned int &operator[](const unsigned int &i) const
-
inline void rebase(const unsigned int &idx_base)
use
rebase()
when adding more vertices and triangles into a mesh.
-
inline void reverse_winding()
in-place reverse the winding of this triangle
Public Members
-
unsigned int point_id[3]
Friends
-
friend std::ostream &operator<<(std::ostream &s, const g_triangle &t)
-
inline g_triangle(const unsigned int &a0, const unsigned int &a1, const unsigned int &a2)
-
class g_triangle_with_colour_index : public g_triangle
- #include <g_triangle.hh>
a triangle with a colour index - helpful, I think when transfering colour/material information to Blender.
Public Functions
-
inline g_triangle_with_colour_index(const unsigned int &a0, const unsigned int &a1, const unsigned int &a2)
Public Members
-
int colour_index
-
inline g_triangle_with_colour_index(const unsigned int &a0, const unsigned int &a1, const unsigned int &a2)
-
namespace coot
-
namespace api
-
class vertex_with_rotation_translation
- #include <vertex.hh>
20230108-PE copied from generic-vertex.hh (then edited).
vertex with rotation and translation (e.g. for oriented bonds)
Public Functions
-
inline vertex_with_rotation_translation(const glm::vec3 &p, const glm::vec3 &n, const glm::vec4 &c)
constructor
-
inline vertex_with_rotation_translation(const vnc_vertex &v, const glm::vec3 &atom_position, float scale)
constructor
-
inline vertex_with_rotation_translation()
constructor
-
inline vertex_with_rotation_translation(const glm::vec3 &p, const glm::vec3 &n, const glm::vec4 &c)
-
class vn_vertex
- #include <vertex.hh>
a vertex with (only) postion and normal. Useful for instancing, perhaps
-
class vnc_vertex
- #include <vertex.hh>
a vertex with postion, normal and colour
-
class vertex_with_rotation_translation
-
namespace api
-
namespace coot
-
class simple_mesh_t
- #include <simple-mesh.hh>
The basic mesh for transfering mesh geometry and colours.
Public Functions
-
inline simple_mesh_t()
constructor (for vectors)
-
inline explicit simple_mesh_t(const std::string &name_in)
constructor with name
-
inline simple_mesh_t(const std::vector<api::vnc_vertex> &vertices_in, const std::vector<g_triangle> &triangles_in)
-
inline explicit simple_mesh_t(const cylinder &cyl)
-
void translate(const glm::vec3 &t)
-
void add_submesh(const simple_mesh_t &submesh)
utilty function
-
void fill_colour_map()
if the colour map is empty then go through the vector of vertices findling colours and putting them into a colour table. This is for Blender - where the colour are assigned to a Material, and a Material is assigned to a face
-
void export_to_gltf(const std::string &file_name, bool use_binary_format) const
export to glTF
-
inline void set_name(const std::string &n)
-
void scale(float scale_factor)
-
void change_colour(const glm::vec4 &c)
-
void clear()
-
std::string vandt() const
Public Members
-
int status
check this status before using a
simple_mesh_t
. 1 is good, 0 is bad (0 is set when when we get a bad_alloc)
-
std::vector<api::vnc_vertex> vertices
vertices
-
std::vector<g_triangle> triangles
vertex index triples
-
std::string name
mesh name
-
std::map<int, glm::vec4> colour_index_to_colour_map
Public Static Functions
-
static simple_mesh_t make_sphere()
-
inline simple_mesh_t()
-
class simple_mesh_t
-
namespace coot
Functions
-
simple_mesh_t instanced_mesh_to_simple_mesh(const instanced_mesh_t &im)
-
class instanced_geometry_t
- #include <instancing.hh>
instancing container for vertices, triangles and instancing data
Public Functions
-
inline instanced_geometry_t()
-
inline explicit instanced_geometry_t(const std::string &n)
-
inline instanced_geometry_t(const std::vector<api::vn_vertex> &v, const std::vector<g_triangle> &t)
Public Members
-
std::vector<g_triangle> triangles
triangle indices
-
std::string name
-
std::vector<instancing_data_type_A_t> instancing_data_A
a vector of type A instancing
-
std::vector<instancing_data_type_B_t> instancing_data_B
a vector of type B instancing
-
inline instanced_geometry_t()
-
class instanced_mesh_t
- #include <instancing.hh>
a simple container for instancing containers - for multiple geometries that are instanced. (e.g. balls and sticks).
Public Functions
-
inline instanced_mesh_t()
-
inline void add(const instanced_geometry_t &ig)
-
inline void clear()
clear
-
void export_to_glTF(const std::string &file_name, bool use_binary_format) const
Public Members
-
std::vector<instanced_geometry_t> geom
-
simple_mesh_t markup
cis-peptide markup can’t be drawn instanced
-
inline instanced_mesh_t()
-
class instancing_data_type_A_t
- #include <instancing.hh>
class for A type instancing data - this does not contain an orientation matrix
-
class instancing_data_type_B_t
- #include <instancing.hh>
class for B type instancing data - this does contain an orientation matrix
-
simple_mesh_t instanced_mesh_to_simple_mesh(const instanced_mesh_t &im)
Atom and Residue Specifiers
-
namespace coot
Functions
-
bool compare_atom_specs_user_float(const atom_spec_t &a1, const atom_spec_t &a2)
-
bool compare_atom_specs_user_float_in_pair(const std::pair<atom_spec_t, std::string> &a, const std::pair<atom_spec_t, std::string> &b)
-
std::pair<atom_spec_t, atom_spec_t> link_atoms(mmdb::Link *link, mmdb::Model *model_p = 0)
-
std::pair<atom_spec_t, atom_spec_t> link_atoms(mmdb::LinkR *link, mmdb::Model *model_p = 0)
-
class atom_spec_t
Public Functions
-
inline atom_spec_t()
-
inline atom_spec_t(const std::string &chain_in, int resno_in, const std::string &insertion_code_in, const std::string &atom_name_in, const std::string &alt_conf_in)
-
inline explicit atom_spec_t(mmdb::Atom *at)
-
inline atom_spec_t(mmdb::Atom *at, const std::string &user_data_string)
-
inline bool empty() const
-
inline void selectatoms(mmdb::Manager *mol, int SelHnd)
-
bool matches_spec(mmdb::Atom *atom) const
-
std::string format() const
-
std::string label() const
-
std::string label(const std::string &residue_name) const
-
std::string simple_label(const std::string &residue_name = "") const
-
inline bool is_same(const atom_spec_t &matcher) const
-
mmdb::Atom *get_atom(mmdb::Manager *mol) const
-
inline atom_spec_t()
-
class residue_spec_t
Public Functions
-
inline explicit residue_spec_t(int r)
-
inline residue_spec_t(const std::string &chain_in, int r)
-
inline residue_spec_t(int model_number_in, const std::string &chain_in, int r, const std::string &ins_code_in)
-
inline residue_spec_t(const std::string &chain_in, int r, const std::string &ins_code_in)
-
inline explicit residue_spec_t(mmdb::Residue *res)
-
inline explicit residue_spec_t(const atom_spec_t &atom_spec)
-
inline residue_spec_t()
-
inline bool unset_p() const
-
inline bool empty() const
-
inline residue_spec_t next() const
-
inline residue_spec_t previous() const
-
std::string format() const
-
std::string label() const
-
std::string label(const std::string &residue_name) const
-
mmdb::Residue *get_residue(mmdb::Manager *mol) const
-
int select_atoms(mmdb::Manager *mol, int selhnd, mmdb::SELECTION_KEY selection_key)
-
inline explicit residue_spec_t(int r)
-
bool compare_atom_specs_user_float(const atom_spec_t &a1, const atom_spec_t &a2)
Validation Information
-
namespace coot
Enums
-
class chain_validation_information_t
Public Functions
-
inline chain_validation_information_t()
-
inline explicit chain_validation_information_t(const std::string &chain_id_in)
-
inline void add(const residue_validation_information_t &rvi)
-
inline void add_residue_validation_information(const residue_validation_information_t &rvi)
-
inline chain_validation_information_t()
-
class validation_information_min_max_t
-
class validation_information_t
Public Functions
-
inline validation_information_t()
-
inline validation_information_t(graph_data_type gdt, const validation_information_min_max_t &min_max_in)
-
inline unsigned int get_index_for_chain(const std::string &chain_id)
-
inline void add_residue_validation_information(const residue_validation_information_t &rvi, const std::string &chain_id)
-
inline void add(const chain_validation_information_t &cvi)
-
inline bool empty() const
return true of there are no data
-
inline void set_min_max()
set the min and max for the graph - internal use.
Public Members
-
std::string name
-
std::vector<chain_validation_information_t> cviv
-
enum graph_data_type type
-
inline validation_information_t()
-
class chain_validation_information_t
-
namespace coot
-
class residue_validation_information_t
Public Functions
-
inline residue_validation_information_t()
-
inline residue_validation_information_t(const coot::residue_spec_t &rs, const coot::atom_spec_t &atom_spec_in, double d, const std::string &l)
Public Members
-
residue_spec_t residue_spec
-
atom_spec_t atom_spec
-
double function_value
-
std::string label
-
inline residue_validation_information_t()
-
class residue_validation_information_t
-
namespace coot
Enums
-
namespace util
-
class correlation_parts_t
-
class density_correlation_stats_info_t
- #include <coot-density-stats.hh>
density correlation state - between map and model typically.
Public Functions
-
inline density_correlation_stats_info_t()
constructor
-
inline density_correlation_stats_info_t(double n_in, double sum_xy_in, double sum_sqrd_x_in, double sum_sqrd_y_in, double sum_x_in, double sum_y_in)
constructor
-
inline double var_x() const
-
inline double var_y() const
-
inline void add(const double &x, const double &y)
-
inline double correlation() const
the correlation
-
inline correlation_parts_t correlation_by_parts() const
the correlation by parts
-
inline bool operator<(const density_correlation_stats_info_t &dcsi) const
-
inline density_correlation_stats_info_t()
-
class correlation_parts_t
-
namespace util
Superposition
-
class superpose_results_t
- #include <superpose-results.hh>
superposition results
Public Members
-
std::string superpose_info
the JSON string with the superposition statistics
-
std::pair<std::string, std::string> alignment
the superposition alignment
-
std::vector<coot::validation_information_t> alignment_info_vec
post-alignment distances between aligned residues encoded as validation information for the reference sequence
-
std::vector<std::pair<coot::residue_validation_information_t, coot::residue_validation_information_t>> aligned_pairs
-
std::string superpose_info
Molecule Internals
-
namespace coot
Enums
Values:
-
enumerator RESIDUE_NUMBER_UNSET
-
enumerator RESIDUE_NUMBER_UNSET
-
class molecule_t
Public Functions
-
inline molecule_t(const std::string &name_in, int mol_no_in)
constructor
-
inline molecule_t(const std::string &name_in, int mol_no_in, short int is_em_map)
constructor, when we know we are giving it an em map
-
inline molecule_t(const std::string &name_in, int mol_no_in, const clipper::Xmap<float> &xmap_in, bool is_em_map_flag)
constructor
-
inline explicit molecule_t(atom_selection_container_t asc, int imol_no_in, const std::string &name_in)
constructor
-
float get_median_temperature_factor() const
-
int close_yourself()
-
inline bool is_closed() const
-
inline void set_really_do_backups(bool state)
-
void add_neighbor_residues_for_refinement_help(mmdb::Manager *mol)
-
void fill_fobs_sigfobs()
-
inline clipper::HKL_data<clipper::data32::F_sigF> *get_original_fobs_sigfobs() const
-
inline clipper::HKL_data<clipper::data32::Flag> *get_original_rfree_flags() const
-
int sfcalc_genmap(const clipper::HKL_data<clipper::data32::F_sigF> &fobs, const clipper::HKL_data<clipper::data32::Flag> &free, clipper::Xmap<float> *xmap_p)
-
util::sfcalc_genmap_stats_t sfcalc_genmaps_using_bulk_solvent(const clipper::HKL_data<clipper::data32::F_sigF> &fobs, const clipper::HKL_data<clipper::data32::Flag> &free, clipper::Xmap<float> *xmap_2fofc_p, clipper::Xmap<float> *xmap_fofc_p)
-
std::pair<std::string, std::string> make_import_datanames(const std::string &fcol, const std::string &phi_col, const std::string &weight_col, int use_weights) const
-
inline std::string Refmac_fobs_col() const
-
inline std::string Refmac_sigfobs_col() const
-
inline std::string Refmac_mtz_filename() const
-
void associate_data_mtz_file_with_map(const std::string &data_mtz_file_name, const std::string &f_col, const std::string &sigf_col, const std::string &r_free_col)
-
util::map_molecule_centre_info_t get_map_molecule_centre() const
-
void replace_molecule_by_model_from_file(const std::string &pdb_file_name)
-
inline std::string get_name() const
-
inline void set_molecule_name(const std::string &n)
-
inline int get_molecule_index() const
-
bool is_valid_model_molecule() const
-
bool is_valid_map_molecule() const
-
unsigned int get_number_of_atoms() const
-
int get_number_of_hydrogen_atoms() const
-
float get_molecule_diameter() const
-
mmdb::Residue *cid_to_residue(const std::string &cid) const
-
std::vector<mmdb::Residue*> cid_to_residues(const std::string &cid) const
-
mmdb::Atom *cid_to_atom(const std::string &cid) const
-
std::pair<bool, residue_spec_t> cid_to_residue_spec(const std::string &cid) const
-
std::pair<bool, atom_spec_t> cid_to_atom_spec(const std::string &cid) const
-
std::vector<std::string> get_residue_names_with_no_dictionary(const protein_geometry &geom) const
-
int insert_waters_into_molecule(const minimol::molecule &water_mol)
-
void make_bonds(protein_geometry *geom, rotamer_probability_tables *rot_prob_tables_p, bool draw_hydrogen_atoms_flag, bool draw_missing_loops_flag)
-
std::vector<glm::vec4> make_colour_table(bool against_a_dark_background) const
useful for debugging, perhaps
-
std::vector<glm::vec4> make_colour_table_for_goodsell_style(float colour_wheel_rotation_step, float saturation, float goodselliness) const
-
void print_colour_table(const std::string &debugging_label) const
-
void print_secondary_structure_info() const
-
mmdb::Atom *get_atom(const atom_spec_t &atom_spec) const
-
mmdb::Residue *get_residue(const residue_spec_t &residue_spec) const
-
mmdb::Residue *get_residue(const std::string &residue_cid) const
-
std::string get_residue_name(const residue_spec_t &residue_spec) const
-
inline bool have_unsaved_changes() const
-
int undo()
-
int redo()
-
int write_coordinates(const std::string &file_name) const
-
std::string molecule_to_PDB_string() const
- Returns:
a model molecule imol as a string. Return emtpy string on error
-
std::string molecule_to_mmCIF_string() const
- Returns:
a model molecule imol as a string. Return emtpy string on error
-
std::vector<atom_spec_t> get_fixed_atoms() const
-
std::vector<std::string> chains_in_model() const
-
std::vector<std::pair<residue_spec_t, std::string>> get_single_letter_codes_for_chain(const std::string &chain_id) const
-
residue_spec_t get_residue_closest_to(mmdb::Manager *mol, const clipper::Coord_orth &co) const
-
std::vector<std::string> get_chain_ids() const
Get the chains that are related by NCS:
-
std::vector<double> get_residue_CA_position(const std::string &cid) const
get the residue CA position
- Returns:
a vector. The length of the vector is 0 on failure, otherwise it is the x,y,z values
-
std::vector<double> get_residue_average_position(const std::string &cid) const
get the avarge residue position
- Returns:
a vector. The length of the vector is 0 on failure, otherwise it is the x,y,z values
-
std::vector<double> get_residue_sidechain_average_position(const std::string &cid) const
get the avarge residue side-chain position
- Returns:
a vector. The length of the vector is 0 on failure, otherwise it is the x,y,z values
-
simple_mesh_t get_bonds_mesh(const std::string &mode, protein_geometry *geom, bool against_a_dark_background, float bonds_width, float atom_radius_to_bond_width_ratio, int smoothness_factor, bool draw_hydrogen_atoms_flag, bool draw_missing_residue_loops)
-
simple_mesh_t get_goodsell_style_mesh(protein_geometry *geom_p, float colour_wheel_rotation_step, float saturation, float goodselliness)
-
instanced_mesh_t get_bonds_mesh_instanced(const std::string &mode, protein_geometry *geom, bool against_a_dark_background, float bonds_width, float atom_radius_to_bond_width_ratio, int smoothness_factor, bool draw_hydrogen_atoms_flag, bool draw_missing_residue_loops)
-
instanced_mesh_t get_bonds_mesh_for_selection_instanced(const std::string &mode, const std::string &selection_cid, protein_geometry *geom, bool against_a_dark_background, float bonds_width, float atom_radius_to_bond_width_ratio, int smoothness_factor, bool draw_hydrogen_atoms_flag, bool draw_missing_residue_loops)
-
instanced_mesh_t get_goodsell_style_mesh_instanced(protein_geometry *geom_p, float colour_wheel_rotation_step, float saturation, float goodselliness)
-
void set_user_defined_bond_colours(const std::map<unsigned int, std::array<float, 4>> &colour_map)
user-defined colour-index to colour (internallly, this converts the
colour_map
to the above vector of colour holders, so it’s probably a good idea if the colour (index) keys are less than 200 or so.
-
void set_user_defined_atom_colour_by_selections(const std::vector<std::pair<std::string, unsigned int>> &indexed_residues_cids, bool colour_applies_to_non_carbon_atoms_also, mmdb::Manager *mol)
user-defined atom selection to colour index.
-
void store_user_defined_atom_colour_selections(const std::vector<std::pair<std::string, unsigned int>> &indexed_residues_cids, bool colour_applies_to_non_carbon_atoms_also)
-
void apply_user_defined_atom_colour_selections(const std::vector<std::pair<std::string, unsigned int>> &indexed_residues_cids, bool colour_applies_to_non_carbon_atoms_also, mmdb::Manager *mol)
-
void set_colour_wheel_rotation_base(float r)
set the colour wheel rotation base for the specified molecule
-
void set_base_colour_for_bonds(float r, float g, float b)
set the base colour - to be used as a base for colour wheel rotation
-
void add_to_non_drawn_bonds(const std::string &atom_selection_cid)
-
inline void clear_non_drawn_bonds()
-
void print_non_drawn_bonds() const
-
void fill_default_colour_rules()
-
void add_colour_rule(const std::string &selection, const std::string &colour_name)
Add a colour rule: eg. (“//A”, “red”)
-
void delete_colour_rules()
delete all the colour rules
-
void print_colour_rules() const
-
std::vector<std::pair<std::string, std::string>> get_colour_rules() const
get the colour rules. Preferentially return the user-defined colour rules.
- Returns:
If there are no user-defined colour rules, then return the stand-in rules
-
void M2T_updateFloatParameter(const std::string ¶m_name, float value)
Update float parameter for MoleculesToTriangles molecular mesh.
-
void M2T_updateFloatParameter(const std::string ¶m_name, int value)
Update int parameter for MoleculesToTriangles molecular mesh.
-
void print_M2T_FloatParameters() const
-
void print_M2T_IntParameters() const
-
void M2T_updateIntParameter(const std::string ¶m_name, int value)
Update int parameter for MoleculesToTriangles molecular mesh.
-
simple_mesh_t get_molecular_representation_mesh(const std::string &cid, const std::string &colour_scheme, const std::string &style, int secondaryStructureUsageFlag) const
-
simple_mesh_t get_gaussian_surface(float sigma, float contour_level, float box_radius, float grid_scale, float fft_b_factor) const
-
simple_mesh_t get_chemical_features_mesh(const std::string &cid, const protein_geometry &geom) const
-
inline bool hydrogen_atom_should_be_drawn() const
-
inline void set_use_bespoke_carbon_atom_colour(bool state)
-
void export_map_molecule_as_gltf(clipper::Coord_orth &position, float radius, float contour_level, const std::string &file_name)
export map molecule as glTF
-
void export_model_molecule_as_gltf(const std::string &mode, const std::string &selection_cid, protein_geometry *geom, bool against_a_dark_background, float bonds_width, float atom_radius_to_bond_width_ratio, int smoothness_factor, bool draw_hydrogen_atoms_flag, bool draw_missing_residue_loops, const std::string &file_name)
export model molecule as glTF - this is the bonds and atoms API
-
void export_molecular_representation_as_gltf(const std::string &atom_selection_cid, const std::string &colour_scheme, const std::string &style, int secondary_structure_usage_flag, const std::string &file_name)
-
void export_chemical_features_as_gltf(const std::string &cid, const protein_geometry &geom, const std::string &file_name) const
-
inline void set_show_symmetry(bool f)
-
inline bool get_show_symmetry()
-
void transform_by(mmdb::mat44 SSMAlign_TMatrix)
-
void transform_by(const clipper::RTop_orth &rtop, mmdb::Residue *res)
-
void transform_by(const clipper::RTop_orth &rtop)
-
std::vector<std::string> non_standard_residue_types_in_model() const
-
std::vector<phi_psi_prob_t> ramachandran_validation(const ramachandrans_container_t &rc) const
-
simple_mesh_t get_rotamer_dodecs(protein_geometry *geom_p, rotamer_probability_tables *rpt)
-
instanced_mesh_t get_rotamer_dodecs_instanced(protein_geometry *geom_p, rotamer_probability_tables *rpt)
-
omega_distortion_info_container_t peptide_omega_analysis(const protein_geometry &geom, const std::string &chain_id, bool mark_cis_peptides_as_bad_flag) const
-
std::vector<residue_spec_t> get_non_standard_residues_in_molecule() const
-
instanced_mesh_t contact_dots_for_ligand(const std::string &cid, const protein_geometry &geom, unsigned int num_subdivisions) const
- Returns:
the instanced mesh for the specified ligand
-
instanced_mesh_t all_molecule_contact_dots(const coot::protein_geometry &geom, unsigned int num_subdivisions) const
- Returns:
the instanced mesh for the specified molecule
-
generic_3d_lines_bonds_box_t make_exportable_environment_bond_box(coot::residue_spec_t &spec, coot::protein_geometry &geom) const
-
simple::molecule_t get_simple_molecule(int imol, const std::string &residue_cid, const bool draw_hydrogen_atoms_flag, coot::protein_geometry *geom_p)
we pass the imol because we use that to look up the residue type in the dictionary annoyingly, we pass a non-const pointer to the protein-geometry because that is what is passed in the Bond_lines_container. Think about changin that one day.
-
simple::molecule_t get_simple_molecule(int imol, mmdb::Residue *residue_p, bool draw_hydrogen_atoms_flag, coot::protein_geometry *geom_p)
-
coot::simple_mesh_t get_mesh_for_ligand_validation_vs_dictionary(const std::string &ligand_cid, coot::protein_geometry &geom, ctpl::thread_pool &static_thread_pool)
get the mesh for ligand validation vs dictionary, coloured by badness. greater then 3 standard deviations is fully red. Less than 0.5 standard deviations is fully green.
-
std::vector<coot::geometry_distortion_info_container_t> geometric_distortions_from_mol(const std::string &ligand_cid, bool with_nbcs, coot::protein_geometry &geom, ctpl::thread_pool &static_thread_pool)
-
std::pair<int, double> simple_geometric_distortions_from_mol(const std::string &ligand_cid, bool with_nbcs, coot::protein_geometry &geom, ctpl::thread_pool &static_thread_pool)
-
coot::instanced_mesh_t get_extra_restraints_mesh(int mode) const
-
std::vector<coot::residue_spec_t> residues_near_residue(const std::string &residue_cid, float dist) const
- Returns:
a list of residues specs that have atoms within dist of the atoms of the specified residue
-
std::vector<plain_atom_overlap_t> get_overlaps(protein_geometry *geom_p)
not const because it can dynamically add dictionaries
-
std::vector<plain_atom_overlap_t> get_overlaps_for_ligand(const std::string &cid_ligand, protein_geometry *geom_p)
not const because it can dynamically add dictionaries
-
coot::atom_overlaps_dots_container_t get_overlap_dots(protein_geometry *geom_p)
not const because it can dynamically add dictionaries
-
coot::atom_overlaps_dots_container_t get_overlap_dots_for_ligand(const std::string &cid_ligand, protein_geometry *geom_p)
not const because it can dynamically add dictionaries
-
int flip_peptide(const atom_spec_t &rs, const std::string &alt_conf)
-
int auto_fit_rotamer(const std::string &chain_id, int res_no, const std::string &ins_code, const std::string &alt_conf, const clipper::Xmap<float> &xmap, const coot::protein_geometry &pg)
-
std::pair<bool, float> backrub_rotamer(mmdb::Residue *residue_p, const clipper::Xmap<float> &xmap, const coot::protein_geometry &pg)
-
std::pair<bool, float> backrub_rotamer(const std::string &chain_id, int res_no, const std::string &ins_code, const std::string &alt_conf, const clipper::Xmap<float> &xmap, const coot::protein_geometry &pg)
-
int delete_atoms(const std::vector<atom_spec_t> &atoms)
-
int delete_atom(atom_spec_t &atom_spec)
-
int delete_residue(residue_spec_t &residue_spec)
-
int delete_residue_atoms_with_alt_conf(coot::residue_spec_t &residue_spec, const std::string &alt_conf)
-
int delete_chain_using_atom_cid(const std::string &cid)
-
int delete_literal_using_cid(const std::string &cid)
-
std::pair<int, std::string> add_terminal_residue_directly(const residue_spec_t &spec, const std::string &new_res_type, const protein_geometry &geom, const clipper::Xmap<float> &xmap, ctpl::thread_pool &static_thread_pool)
-
int add_compound(const dictionary_residue_restraints_t &monomer_restraints, const Cartesian &position, const clipper::Xmap<float> &xmap, float map_rmsd)
-
int add_alternative_conformation(const std::string &cid)
add an alternative conformation for the specified residue
-
int fill_partial_residue(const residue_spec_t &res_spec, const std::string &alt_conf, const clipper::Xmap<float> &xmap, const protein_geometry &geom)
add atoms to a partially-filled side chaain
-
int fill_partial_residues(const clipper::Xmap<float> &xmap, protein_geometry *geom)
add atoms to a partially-filled side chaain
-
int mutate(const residue_spec_t &spec, const std::string &new_res_type)
-
int side_chain_180(const residue_spec_t &residue_spec, const std::string &alt_conf, coot::protein_geometry *geom_p)
-
int delete_side_chain(const residue_spec_t &residue_spec)
-
std::string jed_flip(coot::residue_spec_t &spec, const std::string &atom_name, const std::string &alt_conf, bool invert_selection, protein_geometry *geom)
-
std::string jed_flip_internal(coot::atom_tree_t &tree, const std::vector<coot::dict_torsion_restraint_t> &interesting_torsions, const std::string &atom_name, bool invert_selection)
-
std::string jed_flip_internal(coot::atom_tree_t &tree, const dict_torsion_restraint_t &torsion, const std::string &atom_name, bool invert_selection)
-
int match_torsions(mmdb::Residue *res_ref, const std::vector<coot::dict_torsion_restraint_t> &tr_ligand, const coot::protein_geometry &geom)
-
coot::minimol::molecule eigen_flip_residue(const residue_spec_t &residue_spec)
-
int apply_transformation_to_atom_selection(const std::string &atom_selection_cid, int n_atoms_in_selection, clipper::Coord_orth &rotation_centre, clipper::RTop_orth &rtop)
-
void multiply_residue_temperature_factors(const std::string &cid, float factor)
Interactive B-factor refinement (fun). “factor” might typically be say 0.9 or 1.1
-
int add_hydrogen_atoms(protein_geometry *geom)
- Returns:
1 on a successful additions, 0 on failure.
-
int delete_hydrogen_atoms()
- Returns:
1 on a successful additions, 0 on failure.
-
std::pair<int, std::string> change_chain_id(const std::string &from_chain_id, const std::string &to_chain_id, bool use_resno_range, int start_resno, int end_resno)
-
std::pair<int, std::string> change_chain_id_with_residue_range(const std::string &from_chain_id, const std::string &to_chain_id, int start_resno, int end_resno)
-
void change_chain_id_with_residue_range_helper_insert_or_add(mmdb::Chain *to_chain_p, mmdb::Residue *new_residue)
-
int new_positions_for_residue_atoms(const std::string &residue_cid, const std::vector<api::moved_atom_t> &moved_atoms)
set new positions for the atoms in the specified residue
-
int new_positions_for_atoms_in_residues(const std::vector<api::moved_residue_t> &moved_residues)
set new positions for the atoms of the specified residues
-
int new_positions_for_residue_atoms(mmdb::Residue *residue_p, const std::vector<api::moved_atom_t> &moved_atoms, bool do_backup)
not for wrapping (should be private). We don’t want this function to backup if the backup happens in the calling function (i.e. new_positions_for_atoms_in_residues).
-
int merge_molecules(const std::vector<mmdb::Manager*> &mols)
merge molecules - copy the atom of mols into this molecule
- Returns:
the number of atoms added.
-
float fit_to_map_by_random_jiggle(const residue_spec_t &res_spec, const clipper::Xmap<float> &xmap, float map_rmsd, int n_trials, float translation_scale_factor)
-
float fit_to_map_by_random_jiggle_using_atom_selection(const std::string &cid, const clipper::Xmap<float> &xmap, float map_rmsd, int n_trials, float translation_scale_factor)
My ligands don’t jiggle-jiggle…
Hey, what do you know, they actually do.
-
int cis_trans_conversion(const std::string &atom_cid, mmdb::Manager *standard_residues_mol)
-
int replace_residue(const std::string &residue_cid, const std::string &new_residue_type, int imol_enc, const protein_geometry &geom)
-
int mutate_by_overlap(mmdb::Residue *residue_p, const dictionary_residue_restraints_t &restraints)
-
int replace_fragment(atom_selection_container_t asc)
- Returns:
the success status
-
int replace_fragment(mmdb::Manager *mol_ref, int old_atom_index_handle, int SelHnd)
-
rotamer_change_info_t change_to_next_rotamer(const coot::residue_spec_t &res_spec, const std::string &alt_conf, const coot::protein_geometry &pg)
change rotamers
-
rotamer_change_info_t change_to_previous_rotamer(const coot::residue_spec_t &res_spec, const std::string &alt_conf, const coot::protein_geometry &pg)
-
rotamer_change_info_t change_to_first_rotamer(const coot::residue_spec_t &res_spec, const std::string &alt_conf, const coot::protein_geometry &pg)
-
rotamer_change_info_t change_rotamer_number(const coot::residue_spec_t &res_spec, const std::string &alt_conf, int rotamer_change_direction, const coot::protein_geometry &pg)
-
void associate_sequence_with_molecule(const std::string &chain_id, const std::string &sequence)
-
void assign_sequence(const clipper::Xmap<float> &xmap, const coot::protein_geometry &geom)
try to fit all of the sequences to all of the chains
-
bool is_het_residue(mmdb::Residue *residue_p) const
-
std::pair<short int, int> next_residue_number_in_chain(mmdb::Chain *w, bool new_res_no_by_hundreds = false) const
-
mmdb::Residue *copy_and_add_residue_to_chain(mmdb::Chain *this_model_chain, mmdb::Residue *add_model_residue, bool new_resno_by_hundreds_flag = true)
-
void copy_and_add_chain_residues_to_chain(mmdb::Chain *new_chain, mmdb::Chain *this_molecule_chain)
-
std::vector<std::string> map_chains_to_new_chains(const std::vector<std::string> &adding_model_chains, const std::vector<std::string> &this_model_chains) const
-
std::string suggest_new_chain_id(const std::string ¤t_chain_id) const
-
std::pair<bool, std::vector<std::string>> try_add_by_consolidation(mmdb::Manager *adding_mol)
-
bool merge_molecules_just_one_residue_homogeneous(atom_selection_container_t molecule_to_add)
-
bool merge_molecules_just_one_residue_at_given_spec(atom_selection_container_t molecule_to_add, residue_spec_t target_spec)
-
std::pair<bool, coot::residue_spec_t> merge_ligand_to_near_chain(mmdb::Manager *mol)
-
std::pair<int, std::vector<merge_molecule_results_info_t>> merge_molecules(const std::vector<atom_selection_container_t> &add_molecules)
-
void read_extra_restraints(const std::string &file_name)
read extra restraints (e.g. from ProSMART)
-
std::vector<mmdb::Residue*> select_residues(const residue_spec_t &spec, const std::string &mode) const
refinement tool
-
std::vector<mmdb::Residue*> select_residues(const std::string &chain_id, int resno_start, int resno_end) const
resno_start and resno_end are inclusive
-
std::vector<mmdb::Residue*> select_residues(const std::string &multi_cid, const std::string &mode) const
select residues given a multi-cid
-
int refine_direct(std::vector<mmdb::Residue*> rv, const std::string &alt_loc, const clipper::Xmap<float> &xmap, float map_weight, int n_cycles, const coot::protein_geometry &geom, bool do_rama_plot_restraints, float rama_plot_weight, bool do_torsion_restraints, float torsion_weight, bool refinement_is_quiet)
real space refinement
-
int minimize(const std::string &atom_selection_cid, int n_cycles, bool do_rama_plot_restraints, float rama_plot_weight, bool do_torsion_restraints, float torsion_weight, bool refinement_is_quiet, coot::protein_geometry *geom_p)
-
bool shiftfield_b_factor_refinement(const clipper::HKL_data<clipper::data32::F_sigF> &F_sigF, const clipper::HKL_data<clipper::data32::Flag> &free_flag)
-
void fix_atom_selection_during_refinement(const std::string &atom_selection_cid)
-
void init_all_molecule_refinement(int imol_ref_mol, coot::protein_geometry &geom, const clipper::Xmap<float> &xmap, float map_weight, ctpl::thread_pool *thread_pool)
-
void add_target_position_restraint(const std::string &atom_cid, float pos_x, float pos_y, float pos_z)
-
void turn_off_when_close_target_position_restraint()
-
instanced_mesh_t add_target_position_restraint_and_refine(const std::string &atom_cid, float pos_x, float pos_y, float pos_z, int n_cyles, coot::protein_geometry *geom_p)
-
void clear_target_position_restraint(const std::string &atom_cid)
clear
-
int refine_using_last_restraints(int n_steps)
refine (again).
- Returns:
the status of the refinement: GSL_CONTINUE, GSL_SUCCESS, GSL_ENOPROG (no progress). i.e. don’t call thus function again unless the status is GSL_CONTINUE (-2);
-
inline restraints_container_t *get_last_restraints()
-
void clear_target_position_restraints()
clear any and all drag-atom target position restraints
-
void clear_refinement()
call this after molecule refinement has finished (say when the molecule molecule is accepted into the original molecule)
-
void generate_chain_self_restraints(float local_dist_max, const std::string &chain_id, const coot::protein_geometry &geom)
-
void generate_local_self_restraints(float local_dist_max, const std::vector<coot::residue_spec_t> &residue_specs, const coot::protein_geometry &geom)
-
void generate_local_self_restraints(float local_dist_max, const std::string &multi_selection_cid, const coot::protein_geometry &geom)
-
void generate_local_self_restraints(int selHnd, float local_dist_max, const coot::protein_geometry &geom)
-
void add_parallel_plane_restraint(coot::residue_spec_t spec_1, coot::residue_spec_t spec_2)
-
std::vector<std::string> nucelotide_residue_name_to_base_atom_names(const std::string &rn) const
-
std::vector<std::string> residue_name_to_plane_atom_names(const std::string &rn) const
-
void clear_extra_restraints()
-
int rigid_body_fit(const std::string &mult_cids, const clipper::Xmap<float> &xmap)
-
int rotate_around_bond(const std::string &residue_cid, const std::string &alt_conf, coot::atom_name_quad quad, double torsion_angle, protein_geometry &geom)
-
bool is_EM_map() const
-
float get_density_at_position(const clipper::Coord_orth &pos) const
-
float get_map_mean() const
-
float get_map_rmsd_approx() const
-
int write_map(const std::string &file_name) const
-
void set_map_is_difference_map(bool flag)
-
bool is_difference_map_p() const
-
inline void set_updating_maps_diff_diff_map_peaks(const std::vector<std::pair<clipper::Coord_orth, float>> &v)
-
std::vector<std::pair<clipper::Coord_orth, float>> get_updating_maps_diff_diff_map_peaks(const clipper::Coord_orth &screen_centre) const
does the peaks-move operation.
-
float get_suggested_initial_contour_level() const
- Returns:
the suggested initial contour level. Return -1 on not-a-map
-
simple_mesh_t get_map_contours_mesh(clipper::Coord_orth position, float radius, float contour_level, bool use_thread_pool, ctpl::thread_pool *thread_pool_p)
-
simple_mesh_t get_map_contours_mesh_using_other_map_for_colours(const clipper::Coord_orth &position, float radius, float contour_level, const clipper::Xmap<float> &xmap)
-
histogram_info_t get_map_histogram(unsigned int n_bins, float zoom_factor) const
- Returns:
the map histogram The caller should select the number of bins - 200 is a reasonable default. The caller should also set the zoom factor (which reduces the range by the given factor) centred around the median (typically 1.0 but usefully can vary until ~20.0).
-
void set_map_colour(colour_holder holder)
-
inline void set_map_colour_saturation(float s)
-
void set_other_map_for_colouring_min_max(float min_v, float max_v)
Set the limit for the colour range for the values from the other map. If the other map were, for example, a map of correlation values, then you’d pass -1.0 and 1.0.
-
inline void set_other_map_for_colouring_invert_colour_ramp(bool state)
-
double sum_density_for_atoms_in_residue(const std::string &cid, const std::vector<std::string> &atom_names, const clipper::Xmap<float> &xmap) const
-
std::vector<interesting_place_t> difference_map_peaks(mmdb::Manager *mol, float n_rmsd) const
-
texture_as_floats_t get_map_section_texture(int section_index, int axis, float data_value_for_bottom, float data_value_for_top) const
-
int get_number_of_map_sections(int axis_id) const
- Returns:
the number of section in the map along the give axis. (0 for X-axis, 1 for y-axis, 2 for Z-axis). return -1 on failure.
-
std::vector<float> get_vertices_for_blender() const
-
std::vector<int> get_triangles_for_blender() const
-
std::vector<float> get_colour_table_for_blender() const
-
void make_mesh_for_bonds_for_blender(const std::string &mode, protein_geometry *geom, bool against_a_dark_background, float bond_width, float atom_radius_to_bond_width_ratio, int smoothness_factor)
-
void make_mesh_for_molecular_representation_for_blender(const std::string &cid, const std::string &colour_scheme, const std::string &style, int secondary_structure_usage_flag)
Make an (internal) mesh
this function doesn’t return a value, instead it stores a
blender_mesh_t
blender_mesh in this model@modifies internal state to fill the internal
blender_mesh
object
-
void make_mesh_for_goodsell_style_for_blender(protein_geometry *geom_p, float colour_wheel_rotation_step, float saturation, float goodselliness)
-
void make_mesh_for_gaussian_surface_for_blender(float sigma, float contour_level, float box_radius, float grid_scale, float b_factor)
Public Members
-
atom_selection_container_t atom_sel
-
float default_temperature_factor_for_new_atoms
-
std::vector<std::pair<bool, mmdb::Residue*>> neighbouring_residues
-
std::string refmac_fobs_col
-
std::string refmac_sigfobs_col
-
std::string refmac_mtz_filename
-
std::string refmac_r_free_col
-
bool refmac_r_free_flag_sensible
-
clipper::Xmap<float> xmap
-
std::map<unsigned int, colour_holder> user_defined_bond_colours
-
std::vector<std::pair<std::string, unsigned int>> indexed_user_defined_colour_selection_cids
-
bool indexed_user_defined_colour_selection_cids_apply_to_non_carbon_atoms_also
-
colour_holder base_colour_for_bonds
-
std::set<int> no_bonds_to_these_atom_indices
-
std::vector<std::pair<std::string, std::string>> colour_rules
If any colour rule has been set for this molecule, then we will use these (and that his its internal colour-by-chain colouring scheme).
the
colour_rules
is a vector of things like: (“//A”, “red”)
-
std::vector<std::pair<std::string, float>> M2T_float_params
-
std::vector<std::pair<std::string, int>> M2T_int_params
-
std::vector<std::pair<mmdb::Atom*, clipper::Coord_orth>> atoms_with_position_restraints
-
clipper::Xmap<float> updating_maps_previous_difference_map
-
std::vector<std::pair<clipper::Coord_orth, float>> updating_maps_diff_diff_map_peaks
-
blender_mesh_t blender_mesh
Public Static Attributes
-
static std::atomic<bool> draw_vector_sets_lock
Private Functions
-
void makebonds(protein_geometry *geom, rotamer_probability_tables *rotamer_tables_p, const std::set<int> &no_bonds_to_these_atoms, bool draw_hydrogen_atoms_flag, bool draw_missing_loops_flag)
-
void make_bonds_type_checked(protein_geometry *geom, const char *s = 0)
-
void make_bonds_type_checked(protein_geometry *geom, rotamer_probability_tables *rot_prob_tables_p, bool draw_hydrogen_atoms_flag, bool draw_missing_loops_flag, const char *s = 0)
-
inline api_bond_colour_t get_bonds_box_type() const
-
void make_colour_by_chain_bonds(protein_geometry *geom, const std::set<int> &no_bonds_to_these_atoms, bool change_c_only_flag, bool goodsell_mode, bool draw_hydrogen_atoms_flag, bool draw_missing_loops_flag, bool do_rota_markup = false, rotamer_probability_tables *rotamer_tables_p = nullptr, bool force_rebonding = true)
-
void make_ca_bonds()
-
glm::vec4 get_bond_colour_by_colour_wheel_position(int icol, api_bond_colour_t bonds_box_type) const
-
void update_map_triangles_using_thread_pool(float radius, Cartesian centre, float contour_level, ctpl::thread_pool *thread_pool_p)
-
short int is_em_map_cached_state()
-
inline bool has_xmap() const
-
glm::vec4 position_to_colour_using_other_map(const clipper::Coord_orth &position, const clipper::Xmap<float> &other_map_for_colouring) const
-
glm::vec4 fraction_to_colour(float f) const
-
void clear_draw_vecs()
-
void clear_diff_map_draw_vecs()
-
void replace_coords(const atom_selection_container_t &asc, bool change_altconf_occs_flag, bool replace_coords_with_zero_occ_flag)
-
bool movable_atom(mmdb::Atom *mol_atom, bool replace_coords_with_zero_occ_flag) const
-
bool moving_atom_matches(mmdb::Atom *at, int this_mol_index_maybe) const
-
void adjust_occupancy_other_residue_atoms(mmdb::Atom *at, mmdb::Residue *residue, short int force_sum_1_flag)
-
int full_atom_spec_to_atom_index(const atom_spec_t &atom_spec) const
-
int full_atom_spec_to_atom_index(const std::string &chain, int resno, const std::string &insertion_code, const std::string &atom_name, const std::string &alt_conf) const
-
std::string name_for_display_manager() const
-
std::string dotted_chopped_name() const
-
std::string make_backup(const std::string &modification_type)
-
void save_history_file_name(const std::string &file)
-
void restore_from_backup(int mod_index, const std::string &cwd)
-
std::pair<int, mmdb::Residue*> find_serial_number_for_insert(int seqnum_for_new, const std::string &ins_code_for_new, const std::string &chain_id) const
-
void remove_TER_internal(mmdb::Residue *res_p)
-
void remove_TER_on_last_residue(mmdb::Chain *chain_p)
-
std::pair<bool, std::string> unused_chain_id() const
-
int append_to_molecule(const minimol::molecule &water_mol)
-
glm::vec4 colour_holder_to_glm(const colour_holder &ch) const
-
void setup_cylinder_clashes(instanced_mesh_t &im, const atom_overlaps_dots_container_t &c, float ball_size, unsigned int num_subdivisions, const std::string &molecule_name_stub) const
modify the im
-
void setup_dots(instanced_mesh_t &im, const atom_overlaps_dots_container_t &c, float ball_size, unsigned int num_subdivisions, const std::string &molecule_name_stub) const
modify the im
-
std::pair<int, std::string> write_shelx_ins_file(const std::string &filename) const
-
int read_shelx_ins_file(const std::string &filename)
-
int add_shelx_string_to_molecule(const std::string &str)
-
inline bool is_from_shelx_ins() const
-
void trim_atom_label_table()
-
void delete_ghost_selections()
-
void update_symmetry()
-
void delete_any_link_containing_residue(const residue_spec_t &res_spec)
this doesn’t do a backup - the calling function is in charge of that
-
void delete_link(mmdb::Link *link, mmdb::Model *model_p)
-
bool sanity_check_atoms(mmdb::Manager *mol) const
-
int set_residue_to_rotamer_move_atoms(mmdb::Residue *res, mmdb::Residue *moving_res)
-
float fit_to_map_by_random_jiggle(mmdb::PPAtom atom_selection, int n_atoms, const clipper::Xmap<float> &xmap, float map_sigma, int n_trials, float jiggle_scale_factor, bool use_biased_density_scoring, std::vector<mmdb::Chain*> chains_for_moving)
-
minimol::molecule rigid_body_fit(const minimol::molecule &mol_in, const clipper::Xmap<float> &xmap, float map_sigma) const
-
std::vector<coot::geometry_distortion_info_container_t> geometric_distortions_from_mol(const atom_selection_container_t &asc, bool with_nbcs, coot::protein_geometry &geom, ctpl::thread_pool &static_thread_pool)
-
inline void init()
Private Members
-
modification_info_t modification_info
-
bool use_gemmi
-
int imol_no
-
bool is_closed_flag
-
int ligand_flip_number
-
std::string name
-
bool is_from_shelx_ins_flag
-
ShelxIns shelxins
-
std::map<residue_spec_t, int> current_rotamer_map
-
bool really_do_backups
-
api_bond_colour_t bonds_box_type
-
graphical_bonds_container bonds_box
-
float bonds_colour_map_rotation
-
bool use_bespoke_grey_colour_for_carbon_atoms
-
short int is_em_map_cached_flag
-
ghost_molecule_display_t map_ghost_info
-
bool xmap_is_diff_map
-
colour_holder map_colour
-
float other_map_for_colouring_min_value
-
float other_map_for_colouring_max_value
-
bool radial_map_colour_invert_flag
-
float radial_map_colour_saturation
-
bool original_fphis_filled
-
bool original_fobs_sigfobs_filled
-
bool original_fobs_sigfobs_fill_tried_and_failed
-
clipper::HKL_data<clipper::datatypes::F_phi<float>> *original_fphis_p
-
clipper::HKL_data<clipper::datatypes::F_sigF<float>> *original_fobs_sigfobs_p
-
clipper::HKL_data<clipper::data32::Flag> *original_r_free_flags_p
-
std::vector<density_contour_triangles_container_t> draw_vector_sets
-
std::vector<density_contour_triangles_container_t> draw_diff_map_vector_sets
-
std::vector<std::pair<int, TRIANGLE>> map_triangle_centres
-
std::vector<std::string> history_filename_vec
-
std::string save_time_string
-
std::vector<atom_spec_t> fixed_atom_specs
-
bool show_symmetry
Private Static Functions
-
static std::string file_to_string(const std::string &fn)
-
class difference_map_peaks_info_t
- #include <coot-molecule.hh>
difference maps peaks class
Public Functions
-
inline difference_map_peaks_info_t(const clipper::Coord_orth &p, float ph)
-
inline difference_map_peaks_info_t(const clipper::Coord_orth &p, float ph)
-
class histogram_info_t
- #include <coot-molecule.hh>
map histogram class
-
class interesting_place_t
- #include <coot-molecule.hh>
The container class for an interesting place.
This documentation doesn’t work and I don’t know why.
Public Functions
-
inline interesting_place_t()
constructor
-
inline interesting_place_t(const std::string &feature_type, const residue_spec_t &rs, const clipper::Coord_orth &pt, const std::string &bl)
constructor
-
inline interesting_place_t(const std::string &feature_type, const clipper::Coord_orth &pt, const std::string &button_label)
constructor
-
inline void set_feature_value(const float &f)
internal to libcootapi function to set the values
-
inline void set_badness_value(const float &b)
-
inline interesting_place_t()
-
class modification_info_t
Public Functions
-
inline modification_info_t()
-
inline modification_info_t(const std::string &mol_name_for_backup, bool is_mmcif)
-
bool have_unsaved_changes() const
-
inline void set_molecule_name(const std::string &molecule_name, bool is_mmcif)
-
inline std::string get_index_string(int idx) const
-
std::string get_backup_file_name_from_index(int idx) const
-
std::string make_backup(mmdb::Manager *mol, const std::string &modification_info_string)
- Returns:
a string that when non-empty is the error message
-
mmdb::Manager *undo(mmdb::Manager *mol)
- Returns:
non-null on success
-
mmdb::Manager *redo()
- Returns:
success status
-
inline modification_info_t()
-
class molecule_save_info_t
Public Functions
-
inline molecule_save_info_t()
-
inline void new_modification(const std::string &mod_string)
-
inline void made_a_save()
-
inline bool have_unsaved_changes() const
-
inline std::string index_string(int idx) const
-
inline void set_modification_index(int idx)
-
inline int get_previous_modification_index() const
-
inline int get_next_modification_index() const
-
inline molecule_save_info_t()
-
inline molecule_t(const std::string &name_in, int mol_no_in)
Hydrogen Bonds
-
namespace moorhen
-
class h_bond
Public Functions
-
inline h_bond()
Public Members
-
h_bond_atom hb_hydrogen
-
h_bond_atom donor
-
h_bond_atom acceptor
-
h_bond_atom donor_neigh
-
h_bond_atom acceptor_neigh
-
double angle_1
-
double angle_2
-
double angle_3
-
double dist
-
bool ligand_atom_is_donor
-
bool hydrogen_is_ligand_atom
-
bool bond_has_hydrogen_flag
-
inline h_bond()
-
class h_bond
3D lines
-
class generic_3d_lines_bonds_box_t
- #include <generic-3d-lines.hh>
This class is for generic line set information. When used for environment distances, the indexing is important/useful - index 0 is for non-bonded contacts and index 1 is for hydrogen bonds (or other interesting bonds).
Public Functions
-
inline generic_3d_lines_bonds_box_t()
-
inline explicit generic_3d_lines_bonds_box_t(const graphical_bonds_container &gbc)
the bonds gnerator makes a graphical_bonds_container and libcootapi uses that to convert to this simple container
Public Members
-
std::vector<std::vector<coot::CartesianPair>> line_segments
a vector of line position pairs
-
inline generic_3d_lines_bonds_box_t()
Map Utils
-
namespace coot
-
namespace util
Typedefs
-
typedef std::pair<double, double> phitheta
Functions
-
clipper::RTop_orth make_rtop_orth_from(mmdb::mat44 mat)
-
float density_at_point(const clipper::Xmap<float> &map_in, const clipper::Coord_orth &co)
-
float density_at_point_by_cubic_interp(const clipper::NXmap<float> &map_in, const clipper::Coord_map &cm)
-
float density_at_point_by_linear_interpolation(const clipper::Xmap<float> &map_in, const clipper::Coord_orth &co)
-
float density_at_point_by_nearest_grid(const clipper::Xmap<float> &map_in, const clipper::Coord_orth &co)
-
float density_at_map_point(const clipper::Xmap<float> &map_in, const clipper::Coord_map &cg)
-
clipper::Grad_orth<double> gradient_at_point(const clipper::Xmap<float> &map_in, const clipper::Coord_orth &co)
-
std::pair<float, float> mean_and_variance(const clipper::Xmap<float> &xmap)
-
density_stats_info_t density_around_point(const clipper::Coord_orth &point, const clipper::Xmap<float> &xmap, float d)
-
bool map_fill_from_mtz(clipper::Xmap<float> *xmap, std::string mtz_file_name, std::string f_col, std::string phi_col, std::string weight_col, short int use_weights, float sampling_rate = 1.5)
-
bool map_fill_from_mtz(clipper::Xmap<float> *xmap, std::string mtz_file_name, std::string f_col, std::string phi_col, std::string weight_col, short int use_weights, float reso_limit_high, short int use_reso_limit_high, float sampling_rate = 1.5)
-
void filter_by_resolution(clipper::HKL_data<clipper::datatypes::F_phi<float>> *fphidata, const float &reso_low, const float &reso_high)
-
float max_gridding(const clipper::Xmap<float> &xmap)
-
float map_score(mmdb::PPAtom atom_selection, int n_selected_atoms, const clipper::Xmap<float> &xmap, short int with_atomic_weighting)
-
float map_score(std::vector<mmdb::Atom*> atoms, const clipper::Xmap<float> &xmap)
-
float map_score_atom(mmdb::Atom *atom, const clipper::Xmap<float> &xmap)
-
float map_score_by_residue_specs(mmdb::Manager *mol, const std::vector<residue_spec_t> &res_specs, const clipper::Xmap<float> &xmap, bool main_chain_only_flag = false)
-
clipper::Xmap<float> sharpen_blur_map(const clipper::Xmap<float> &xmap_in, float b_factor)
-
void sharpen_blur_map(clipper::Xmap<float> *xmap, float b_factor)
-
clipper::Xmap<float> sharpen_blur_map_with_resample(const clipper::Xmap<float> &xmap_in, float b_factor, float resample_factor)
-
clipper::Xmap<float> sharpen_blur_map_with_reduced_sampling(const clipper::Xmap<float> &xmap_in, float b_factor, float resample_factor)
-
void multi_sharpen_blur_map(const clipper::Xmap<float> &xmap_in, const std::vector<float> &b_factors, std::vector<clipper::Xmap<float>> *maps_p)
-
clipper::Xmap<float> sharpen_map(const clipper::Xmap<float> &xmap_in, float sharpen_factor)
-
map_molecule_centre_info_t map_molecule_centre(const clipper::Xmap<float> &xmap)
-
map_molecule_centre_info_t map_molecule_recentre_from_position(const clipper::Xmap<float> &xmap, const clipper::Coord_orth ¤t_centre)
-
std::vector<amplitude_vs_resolution_point> amplitude_vs_resolution(const clipper::Xmap<float> &xmap_in, int n_bins = -1)
-
float b_factor(const std::vector<amplitude_vs_resolution_point> &fsqrd_data, std::pair<bool, float> reso_low_invresolsq = std::pair<bool, float>(false, -1), std::pair<bool, float> reso_higy_invresolsq = std::pair<bool, float>(false, -1))
-
clipper::Xmap<float> transform_map(const clipper::Xmap<float> &xmap_in, const clipper::Spacegroup &new_space_group, const clipper::Cell &new_cell, const clipper::RTop_orth &rtop, const clipper::Coord_orth &about_pt, float box_size)
-
clipper::Grid_sampling suggested_grid_sampling(const clipper::Grid_sampling &orig_sampling, const clipper::Cell &orig_cell, const clipper::Spacegroup &new_space_group, const clipper::Cell &new_cell)
-
clipper::Xmap<float> laplacian_transform(const clipper::Xmap<float> &xmap_in)
-
std::vector<float> density_map_points_in_sphere(clipper::Coord_orth pt, float radius, const clipper::Xmap<float> &xmap_in)
-
clipper::Xmap<float> calc_atom_map(mmdb::Manager *mol, int atom_selection_handle, const clipper::Cell &cell, const clipper::Spacegroup &space_group, const clipper::Grid_sampling &sampling)
-
clipper::Xmap<float> mask_map(const clipper::Xmap<float> &xmap_in, const std::vector<mmdb::Residue*> &neighbs)
-
clipper::Xmap<float> make_map_mask(const clipper::Spacegroup &space_group, const clipper::Cell &cell, const clipper::Grid_sampling &grid_sampling, mmdb::Manager *mol, int atom_selection_handle, float radius, float smooth)
-
float map_to_model_correlation(mmdb::Manager *mol, const std::vector<residue_spec_t> &specs_for_correl, const std::vector<residue_spec_t> &specs_for_masking_neighbs, unsigned short int atom_mask_mode, float atom_radius, const clipper::Xmap<float> &xmap_from_sfs)
-
density_correlation_stats_info_t map_to_model_correlation_stats(mmdb::Manager *mol, const std::vector<residue_spec_t> &specs_for_correl, const std::vector<residue_spec_t> &specs_for_masking_neighbs, unsigned short int atom_mask_mode, float atom_radius, const clipper::Xmap<float> &xmap_from_sfs, map_stats_t map_stats_flag)
-
std::vector<std::pair<residue_spec_t, float>> map_to_model_correlation_per_residue(mmdb::Manager *mol, const std::vector<residue_spec_t> &specs, unsigned short int atom_mask_mode, float atom_radius, const clipper::Xmap<float> &xmap_from_sfs)
-
std::map<coot::residue_spec_t, density_stats_info_t> map_to_model_correlation_stats_per_residue(mmdb::Manager *mol, const std::vector<residue_spec_t> &specs, unsigned short int atom_mask_mode, float atom_radius, const clipper::Xmap<float> &xmap)
-
std::pair<std::map<coot::residue_spec_t, density_correlation_stats_info_t>, std::map<coot::residue_spec_t, density_correlation_stats_info_t>> map_to_model_correlation_stats_per_residue_run(mmdb::Manager *mol, const std::string &chain_id, const clipper::Xmap<float> &xmap, unsigned int n_residues_per_run, bool exclude_CON, float atom_mask_radius = 2.8, float NOC_mask_radius = 1.8)
-
std::pair<clipper::Coord_frac, clipper::Coord_frac> find_struct_fragment_coord_fracs_v2(const std::pair<clipper::Coord_orth, clipper::Coord_orth> &selection_extents, const clipper::Cell &cell)
-
std::vector<std::pair<double, double>> qq_plot_for_map_over_model(mmdb::Manager *mol, const std::vector<coot::residue_spec_t> &specs, const std::vector<coot::residue_spec_t> &nb_residues, int atom_mask_mode, const clipper::Xmap<float> &xmap)
-
std::pair<clipper::Xmap<float>, float> difference_map(const clipper::Xmap<float> &xmap_in_1, const clipper::Xmap<float> &xmap_in_2, float map_scale)
-
clipper::Xmap<float> average_map(const std::vector<std::pair<clipper::Xmap<float>, float>> &maps_and_scales_vec)
-
clipper::Xmap<float> variance_map(const std::vector<std::pair<clipper::Xmap<float>, float>> &maps_and_scales_vec)
-
void regen_weighted_map(clipper::Xmap<float> *xmap_in, const std::vector<std::pair<clipper::Xmap<float>*, float>> &maps_and_scales_vec)
-
std::pair<float, float> spin_search(const clipper::Xmap<float> &xmap, mmdb::Residue *res, coot::torsion tors)
-
clipper::Xmap<float> reinterp_map_fine_gridding(const clipper::Xmap<float> &xmap)
-
clipper::Xmap<float> reinterp_map(const clipper::Xmap<float> &xmap_in, float sampling_multiplier)
-
clipper::Xmap<float> reinterp_map(const clipper::Xmap<float> &xmap_in, const clipper::Xmap<float> &reference_xmap)
-
void reverse_map(clipper::Xmap<float> *xmap_p)
-
std::vector<std::pair<std::string, clipper::Xmap<float>>> partition_map_by_chain(const clipper::Xmap<float> &xmap, mmdb::Manager *mol, std::string *state_string_p)
-
bool is_EM_map(const clipper::Xmap<float> &xmap)
-
float average_of_sample_map_at_sphere_points(clipper::Coord_orth ¢re, float radius, const std::vector<phitheta> &phi_thetas, clipper::Xmap<float> &xmap)
-
std::vector<std::pair<clipper::Resolution, double>> fsc(const clipper::Xmap<float> &xmap_1, const clipper::Xmap<float> &xmap_2)
-
clipper::Xmap<float> power_scale(const clipper::Xmap<float> &xmap_ref, const clipper::Xmap<float> &xmap_for_scaling)
-
void compare_structure_factors(const clipper::Xmap<float> &xmap_1, const clipper::Xmap<float> &xmap_2)
-
void flip_hand(clipper::Xmap<float> *xmap_p)
-
clipper::Xmap<float> analyse_map_point_density_change(const std::vector<std::pair<clipper::Xmap<float>*, float>> &xmaps, const clipper::Xmap<float> &xmap_for_mask)
-
clipper::Xmap<float> zero_dose_extrapolation(const std::vector<std::pair<clipper::Xmap<float>*, float>> &xmaps, const clipper::Xmap<float> &xmap_for_mask)
-
clipper::Xmap<float> real_space_zero_dose_extrapolation(const std::vector<clipper::Xmap<float>*> &xmaps, const clipper::Xmap<float> &xmap_for_mask)
-
class backrub_residue_triple_t : public coot::util::residue_triple_t
Public Functions
-
inline backrub_residue_triple_t(mmdb::Residue *this_residue_in, mmdb::Residue *prev_residue_in, mmdb::Residue *next_residue_in, const std::string &alt_conf_in)
-
void trim_this_residue_atoms()
-
void trim_prev_residue_atoms()
-
void trim_next_residue_atoms()
-
void trim_residue_atoms_generic(mmdb::Residue *residue_p, std::vector<std::string> keep_atom_vector, bool use_keep_atom_vector)
-
inline backrub_residue_triple_t(mmdb::Residue *this_residue_in, mmdb::Residue *prev_residue_in, mmdb::Residue *next_residue_in, const std::string &alt_conf_in)
-
class map_fragment_info_t
Public Functions
-
map_fragment_info_t(const clipper::Xmap<float> &xmap, const clipper::Coord_orth ¢re, float radius, bool centre_at_origin = false)
-
void unshift(clipper::Xmap<float> *xmap_p, const clipper::Coord_orth ¢re)
-
void simple_origin_shift(const clipper::Xmap<float> &ip_xmap, const clipper::Coord_orth ¢re, float radius)
-
clipper::Grid_map make_grid_map(const clipper::Xmap<float> &input_xmap, const clipper::Coord_orth ¢re) const
-
map_fragment_info_t(const clipper::Xmap<float> &xmap, const clipper::Coord_orth ¢re, float radius, bool centre_at_origin = false)
-
class map_molecule_centre_info_t
- #include <coot-map-utils.hh>
map molecule centre
Public Functions
-
inline map_molecule_centre_info_t()
-
inline map_molecule_centre_info_t()
-
class map_ref_triple_t
Public Functions
-
inline map_ref_triple_t(const double &d_in, const clipper::Xmap<float>::Map_reference_coord &iw_in, const float &den_in)
-
inline map_ref_triple_t()
-
inline bool operator<(const map_ref_triple_t &mrt) const
-
inline map_ref_triple_t(const double &d_in, const clipper::Xmap<float>::Map_reference_coord &iw_in, const float &den_in)
-
class map_stats_holder_helper_t
-
class residue_triple_t
Subclassed by coot::util::backrub_residue_triple_t
Public Functions
-
inline residue_triple_t()
-
inline residue_triple_t(mmdb::Residue *this_residue_in, mmdb::Residue *prev_residue_in, mmdb::Residue *next_residue_in, const std::string &alt_conf_in)
-
inline ~residue_triple_t()
-
inline residue_triple_t deep_copy()
-
inline residue_triple_t()
-
class segment_map
Public Functions
-
inline segment_map()
-
std::pair<int, clipper::Xmap<int>> segment(const clipper::Xmap<float> &xmap_in, float low_level)
-
std::pair<int, clipper::Xmap<int>> segment_emsley_flood(const clipper::Xmap<float> &xmap_in, float low_level)
-
std::pair<int, clipper::Xmap<int>> segment(const clipper::Xmap<float> &xmap_in, float low_level, float b_factor_increment, int n_rounds)
Private Functions
-
void flood_fill_segmented_map(clipper::Xmap<std::pair<bool, int>> *segmented_map, const clipper::Xmap<float> &xmap, const clipper::Coord_grid &seed_point, int from_val, int to_val)
-
std::vector<clipper::Coord_grid> path_to_peak(const clipper::Coord_grid &start_point, const clipper::Xmap<float> &xmap_new)
-
int find_biggest_segment(const std::map<int, std::vector<clipper::Coord_grid>> &segment_id_map, const std::map<int, int> &segment_id_counter_map) const
-
int find_smallest_segment(const std::map<int, std::vector<clipper::Coord_grid>> &segment_id_map, const std::map<int, int> &segment_id_counter_map) const
-
void resegment_watershed_points(clipper::Xmap<int> *xmap_int, const clipper::Xmap<float> &xmap) const
-
inline segment_map()
-
class simple_residue_triple_t
-
class soi_variance
Public Functions
-
inline soi_variance(const clipper::Xmap<float> &xmap_in)
-
void proc(float solvent_content_frac)
Public Static Functions
-
static bool mri_var_pair_sorter(const std::pair<clipper::Xmap_base::Map_reference_index, float> &p1, const std::pair<clipper::Xmap_base::Map_reference_index, float> &p2)
Private Functions
-
clipper::Xmap<float> make_variance_map() const
-
clipper::Xmap<float> solvent_treatment_map() const
-
clipper::Xmap<float> protein_treatment_map() const
Private Members
-
const clipper::Xmap<float> &xmap
Private Static Functions
-
static void apply_variance_values(clipper::Xmap<float> &variance_map, const clipper::Xmap<float> &xmap, const std::vector<clipper::Coord_grid> &soi_gps, const std::vector<clipper::Xmap_base::Map_reference_index> &grid_indices)
-
inline soi_variance(const clipper::Xmap<float> &xmap_in)
-
typedef std::pair<double, double> phitheta
-
namespace util
Colour
-
namespace coot
-
-
class colour_t
Public Functions
-
inline colour_t()
-
inline colour_t(float r, float g, float b)
-
inline void set(float r, float g, float b)
-
inline float &operator[](const unsigned int &idx)
-
inline const float &operator[](const unsigned int &idx) const
-
void rotate(float f)
-
inline void brighter(float f)
-
inline glm::vec4 to_glm() const
-
inline colour_holder to_colour_holder() const
Public Members
-
std::vector<float> col
-
inline colour_t()
-
class colour_t
Cartesian Coordinates
Enums
-
namespace coot
Functions
-
std::ostream &operator<<(std::ostream&, CartesianPair)
-
std::ofstream &operator<<(std::ofstream&, CartesianPair)
-
short int is_an_in_triangle(surface_face_data face, const Cartesian &b, const Cartesian &c)
this should be a boolean
-
class Cartesian
Public Functions
-
inline float get_x() const
get x (disapproved)
-
inline float get_y() const
get y (disapproved)
-
inline float get_z() const
get z (disapproved)
-
inline const float &x() const
get x
-
inline const float &y() const
get y
-
inline const float &z() const
get z
-
inline Cartesian(float xi, float yi, float zi)
-
Cartesian()
-
inline explicit Cartesian(const clipper::Coord_orth &pt)
-
inline void set_them(float a, float b, float c)
-
float amplitude(void) const
-
inline float amplitude_squared(void) const
-
inline float length(void) const
-
short int normalize()
-
void operator*=(float scale)
-
void operator/=(float scale)
-
void invert_z(void)
-
inline void unit_vector_yourself()
Public Static Functions
-
static double DihedralAngle(const Cartesian &a, const Cartesian &b, const Cartesian &c, const Cartesian &d)
-
static Cartesian GetCartFrom3Carts(const Cartesian &Atom1, double blength, const Cartesian &Atom2, double angle1, const Cartesian &Atom3, double angle2, int chiral = 0)
-
static Cartesian GetCartFrom3Carts_intermediate(const Cartesian &Atom1, const Cartesian &Atom2, const Cartesian &Atom3, double blength, double angle1, double angle2, int chiral = 0)
Friends
-
friend surface_face_data on_a_face(const Cartesian &a, const Cartesian &b)
-
inline float get_x() const
-
class CartesianPair
- #include <Cartesian.h>
Cartesian pair.
Public Functions
-
CartesianPair(void)
-
float amplitude() const
get the length of the vector between the coordinates
Friends
-
friend std::ostream &operator<<(std::ostream&, CartesianPair)
-
friend std::ofstream &operator<<(std::ofstream&, CartesianPair)
-
CartesianPair(void)
-
std::ostream &operator<<(std::ostream&, CartesianPair)
Statistics Container
-
namespace coot
-
namespace stats
Functions
-
double get_kolmogorov_smirnov_vs_normal(const std::vector<double> &v1, const double &reference_mean, const double &reference_sd)
-
class pnorm
Public Functions
-
inline pnorm()
-
double erf(const double &z) const
-
inline double get(const double &x) const
Private Functions
-
inline void init()
-
inline pnorm()
-
class single
Public Functions
-
inline single()
-
inline explicit single(const std::vector<double> &v_in)
-
inline unsigned int size() const
-
inline bool empty() const
-
inline void add(const double &a)
-
inline double mean() const
-
inline double variance() const
-
inline double skew() const
-
inline double kurtosis() const
-
inline std::pair<double, double> median_and_iqr() const
-
inline double get_ith_highest(unsigned int idx) const
-
inline double get_ith_lowest(unsigned int idx) const
Public Members
-
std::vector<double> v
-
inline single()
-
double get_kolmogorov_smirnov_vs_normal(const std::vector<double> &v1, const double &reference_mean, const double &reference_sd)
-
namespace stats
Symmetry
Functions
-
std::ostream &operator<<(std::ostream &s, const coot::coot_mat44 &m)
-
std::string to_string(const std::pair<symm_trans_t, Cell_Translation> &sts)
-
mmdb::PPAtom translated_atoms(atom_selection_container_t AtomSel, symm_trans_t symm_trans)
-
coot::Cartesian translate_atom(atom_selection_container_t AtomSel, int i, symm_trans_t symm_trans)
-
coot::Cartesian translate_atom_with_pre_shift(atom_selection_container_t AtomSel, int i, const std::pair<symm_trans_t, Cell_Translation> &symm_trans)
-
int set_mmdb_cell_and_symm(atom_selection_container_t asc, std::pair<std::vector<float>, std::string> cell_spgr)
-
class symm_trans_t
- #include <mmdb-crystal.h>
A class for symmetry operators.
Public Functions
-
inline symm_trans_t(int n, int x, int y, int z)
-
inline explicit symm_trans_t(int idx)
-
inline symm_trans_t()
-
inline int isym() const
the operator index (from SYMINFO)
-
inline int x() const
the x-shift
-
inline int y() const
the y-shift
-
inline int z() const
the z-shift
-
inline void add_shift(int xs, int ys, int zs)
-
bool is_identity()
- Returns:
true if this symmetry operator is the identity matrix
-
std::string str(bool expanded_flag) const
symmetry as string
-
void as_mat44(mmdb::mat44 *mat, mmdb::Manager *mol)
fill mat using mol
-
inline void fill_mat(mmdb::Manager *mol)
fill m
Friends
-
friend std::ostream &operator<<(std::ostream &s, const symm_trans_t &st)
-
inline symm_trans_t(int n, int x, int y, int z)
-
class Cell_Translation
- #include <mmdb-crystal.h>
class for the cell translation for generation/application of symmery-related molecules
Public Functions
-
inline Cell_Translation()
constructor
-
Cell_Translation(int a, int b, int c)
constructor
-
inline Cell_Translation inv() const
inversion
Friends
-
friend std::ostream &operator<<(std::ostream &s, Cell_Translation ct)
-
inline Cell_Translation()
-
class molecule_extents_t
Public Functions
-
molecule_extents_t(atom_selection_container_t, float expansion_size)
-
~molecule_extents_t()
-
Cell_Translation coord_to_unit_cell_translations(coot::Cartesian point, atom_selection_container_t AtomSel) const
-
std::vector<std::pair<symm_trans_t, Cell_Translation>> which_boxes(coot::Cartesian point, atom_selection_container_t AtomSel, int shift_search_size = 1) const
-
std::vector<std::pair<int, symm_trans_t>> which_strict_ncs(const coot::Cartesian ¢re_pt, atom_selection_container_t &AtomSel, const std::vector<coot::coot_mat44> &strict_ncs_mats, const Cell_Translation &c_t) const
-
coot::trans_selection_t trans_sel_o(mmdb::Manager *mol, const symm_trans_t &symm_trans) const
-
mmdb::PPAtom trans_sel(mmdb::Cryst *my_cryst, symm_trans_t symm_trans) const
-
mmdb::PPAtom trans_sel(mmdb::Manager *mol, const symm_trans_t &symm_trans) const
-
mmdb::PPAtom trans_sel(mmdb::Manager *mol, mmdb::mat44 my_mat, int x_shift, int y_shift, int z_shift) const
Private Functions
-
void shift_matrix(mmdb::Manager *mol, mmdb::mat44 my_matt, int x_shift, int y_shift, int z_shift, mmdb::mat44 new_matrix) const
Private Members
-
mmdb::PPAtom extents_selection
-
float expansion_size_
-
Cell_Translation atom_sel_cell_trans
Friends
-
friend std::ostream &operator<<(std::ostream &s, const molecule_extents_t &e)
-
molecule_extents_t(atom_selection_container_t, float expansion_size)
Structure Factors
-
namespace coot
-
namespace util
Functions
-
void sfcalc_genmap(mmdb::Manager *mol, const clipper::HKL_data<clipper::data32::F_sigF> &fobs, const clipper::HKL_data<clipper::data32::Flag> &free, clipper::Xmap<float> *xmap_p)
-
sfcalc_genmap_stats_t sfcalc_genmaps_using_bulk_solvent(mmdb::Manager *mol, const clipper::HKL_data<clipper::data32::F_sigF> &fobs, const clipper::HKL_data<clipper::data32::Flag> &free, const clipper::Cell &cell_for_fobs, clipper::Xmap<float> *xmap_2fofc_p, clipper::Xmap<float> *xmap_fofc_p)
-
class sfcalc_genmap_stats_t
Public Functions
-
inline sfcalc_genmap_stats_t(float r_factor_in, float free_r_factor, float bulk_solvent_volume, float bulk_correction, unsigned int n_splines, const loc_table_t &loc_table)
-
inline sfcalc_genmap_stats_t()
Public Members
-
float r_factor
-
float free_r_factor
-
float bulk_solvent_volume
-
float bulk_correction
-
unsigned int n_splines
-
loc_table_t loc_table
-
class loc_table_t
Public Functions
-
inline loc_table_t()
-
inline void add(const loc_table_item_t &item)
-
inline std::size_t size() const
Public Members
-
std::vector<loc_table_item_t> items
-
inline loc_table_t()
-
inline sfcalc_genmap_stats_t(float r_factor_in, float free_r_factor, float bulk_solvent_volume, float bulk_correction, unsigned int n_splines, const loc_table_t &loc_table)
-
void sfcalc_genmap(mmdb::Manager *mol, const clipper::HKL_data<clipper::data32::F_sigF> &fobs, const clipper::HKL_data<clipper::data32::Flag> &free, clipper::Xmap<float> *xmap_p)
-
namespace util