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/prefered)

inline void set_make_backups(bool state)

Allow the user to disable/enable backups (state is true for “enable”). The default is true.

inline bool get_make_backups() const
Returns

the backup-enabled state

std::string file_name_to_string(const std::string &file_name) const
Returns

the string of the contents of the given file-name.

inline unsigned int get_number_of_molecules() const
Returns

the number of molecules

void create_empty_molecules(unsigned int n_empty)

The adds a number of empty molecules to the internal vector of molecules Note that his is not like reserve as it will increase the molecule index of the next added molecule by n_empty.

inline void set_imol_refinement_map(int i)

set the map used for refinement and fitting

inline void set_map_weight(float w)

set the map weight

inline float get_map_weight() const
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. The test for these failing is spec.empty()

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.

Returns

the residues spec. spec.empty() is true on failure.

inline void set_show_timings(bool s)

this 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. The default is true.

inline coot::protein_geometry &get_geom()

Generic Utils

std::string get_molecule_name(int imol) const
Returns

the name of the molecule

void display_molecule_names_table() const

debugging function: display the table of molecule and names

bool is_valid_model_molecule(int imol) const
Returns

is this a valid model?

bool is_valid_map_molecule(int imol_map) const
Returns

is this a valid map?

bool is_a_difference_map(int imol_map) const
Returns

is this a difference map?

int close_molecule(int imol)

close the molecule (and delete dynamically allocated memory)

Returns

1 on successful closure and 0 on failure to close

void end_delete_closed_molecules()
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)
Returns

the eigenvalues of the atoms in the specified residue

coot::simple_mesh_t test_origin_cube() const
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)

the caller has 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
Returns

a flag of unsaved models state - i.e. 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 stardard list of residues

std::vector<std::string> non_standard_residue_types_in_model(int imol) const
Returns

a vector of non-standard residues (so that they can be used for auxiliary dictionary import)

Coordinates Utils

int read_coordinates(const std::string &file_name)

read a coordinates file (mmcif or PDB)

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.

Returns

the new molecule index on success and -1 on failure

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

std::vector<int> split_multi_model_molecule(int imol)

split an NMR model into multiple models - all in MODEL 1.

Returns

the vector of new molecule indices.

int make_ensemble(const std::string &model_molecule_list)

make a multi-model molecule given the input molecules 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
Returns

the model molecule imol as a string. Return emtpy string on error

std::string molecule_to_mmCIF_string(int imol) const
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

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, say) then return -1 and a blank string.

int import_cif_dictionary(const std::string &cif_file_name, int imol_enc)

import a dictionary cif - 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
Returns

the dictionary read for the give residue type, return an empty string on failure to lookup the residue type

std::string get_cif_restraints_as_string(const std::string &comp_id, int imol_enc) const
Returns

a string that is the contents of a dictionary cif file

int get_monomer(const std::string &monomer_name)

get a monomer

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 - use -999999 (IMOL_ENC_ANY) if no molecule-specific dictionary is needed.

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 - use -999999 if no molecule-specific dictionary is needed

Returns

the new molecule index on success and -1 on failure

std::vector<std::string> get_groups_for_monomers(const std::vector<std::string> &residue_names) const
Returns

the group for the given list of residue names.

std::string get_group_for_monomer(const std::string &residue_name) const
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
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)
Returns

a vector of string pairs that were part of a gphl_chem_comp_info. return an empty vector on failure to find any such info.

void write_png(const std::string &compound_id, int imol, const std::string &file_name) const

write a PNG for the given compound_id. imol can be IMOL_ENC_ANY Currently this function does nothing (drawing is done with the not-allowed cairo)

int write_coordinates(int imol, const std::string &file_name) const

write the coordinate to the give file name

Returns

1 on success and 0 on failure

void set_draw_missing_residue_loops(bool state)

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().

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.

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 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 (but maybe is slower to transfer) 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 (say) 1.5 and 3.0. The ratio for “liquorice” representation is 1.0 (of course).

Returns

a coot::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.

The arguments are as above:

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 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 (say) 1.5 and 3.0. The ratio for “liquorice” representation is 1.0 (of course). 1.7 or 1.8 looks nice.

Returns

a coot::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 above, but only return the bonds for the atom selection. Typically one would call this with a wider bond_with than one would use for standards atoms (all molecule)

Returns

a coot::instanced_mesh_t

coot::instanced_mesh_t get_goodsell_style_mesh_instanced(int imol, float colour_wheel_rotation_step, float saturation, float goodselliness)
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

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 - This API will change - we want to specify surfaces and ribbons too.

void export_molecular_represenation_as_gltf(int imol, const std::string &atom_selection_cid, const std::string &colour_scheme, const std::string &style, const std::string &file_name)
std::vector<glm::vec4> get_colour_table(int imol, bool against_a_dark_background) const

return the colur table (for testing)

void set_colour_wheel_rotation_base(int imol, float r)

set the colour wheel rotation base for the specified molecule (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

void add_to_non_drawn_bonds(int imol, const std::string &atom_selection_cid)

add a atom selection cid for atoms and bonds not to be drawn

void clear_non_drawn_bonds(int imol)

clear the set of non-drawn atoms (so that they can be displayed again)

void print_non_drawn_bonds(int imol) const
void set_user_defined_bond_colours(int imol, const std::map<unsigned int, std::array<float, 3>> &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

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, combined like the following “//A/1^#cc0000|//A/2^#cb0002|//A/3^#c00007” i.e. “|” 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

std::vector<std::pair<std::string, std::string>> get_colour_rules(int imol) const

get the colour rules

void print_colour_rules(int imol) const

print the colour rules

void set_use_bespoke_carbon_atom_colour(int imol, bool state)

use bespoke carbon atom colour

void set_bespoke_carbon_atom_colour(int imol, const coot::colour_t &col)

set bespoke carbon atom colour

void M2T_updateFloatParameter(int imol, const std::string &param_name, float value)

Update float parameter for MoleculesToTriangles molecular mesh.

void M2T_updateIntParameter(int imol, const std::string &param_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)

get ribbon and surface representation

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

These values seem to give a reasonable quite smooth surface:

sigma = 4.4

contour_level = 4.0

box_radius = 5.0

grid_scale = 0.7

b_factor = 100.0 (use 0.0 for no FFT-B-factor smoothing)

Returns

a simple mesh 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 feaatures for the specified residue

unsigned int get_number_of_atoms(int imol) const
Returns

the number of atoms in the specified model, or 0 on error

float get_molecule_diameter(int imol) const
Returns

an estimate of the diameter of the model molecule (-1 on failure)

int get_number_of_hydrogen_atoms(int imol) const
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
Returns

vector of chain-ids for the given molecule

Get the chains that are related by NCS or molecular symmetry:

Returns

a vector of vector 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
Returns

vector 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
Returns

a list of residue that don’t have a dictionary

std::string get_residue_name(int imol, const std::string &chain_id, int res_no, const std::string &ins_code) const

get residue name

Returns

the residue name, return a blank string on residue not found.

std::vector<coot::residue_spec_t> residues_with_missing_atoms(int imol)
Returns

an object that has information about residues without dictionaries and residues with missing atom in the the specified molecule

coot::util::missing_atom_info missing_atoms_info_raw(int imol)

Ths function is not const because missing_atoms() takes a non-const pointer to the geometry.

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
Returns

a list of residues specs that have atoms within dist of the atoms of the specified residue

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). There is some alignment screen output that would be better added to the return value.

coot::symmetry_info_t get_symmetry(int imol, float symmetry_search_radius, float centre_x, float centre_y, float centre_z) const

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.

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

Returns

a coot::util::map_molecule_centre_info_t.

int undo(int imol)

undo

Returns

1 on successful undo, return 0 on failure

int redo(int imol)

redo

Returns

1 on successful redo, return 0 on failure

Map Utils

float map_sampling_rate
Return

the map sampling rate (default is 1.8)

inline void set_map_sampling_rate(float msr)

set the map sampling rate (default is 1.8). Higher numbers mean smoother maps, but they take longer to generate, longer to transfer, longer to parse and longer to draw

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.

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

std::vector<auto_read_mtz_info_t> auto_read_mtz(const std::string &file_name)

Read the given 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)
Returns

the new molecule number or -1 on failure

int write_map(int imol, const std::string &file_name) const

write a map. This function was be renamed from writeMap

Returns

1 on a successful write, return 0 on failure.

float get_map_rmsd_approx(int imol_map) const
Returns

the map rmsd (epsilon testing is not used). -1 is returned if imol_map is not a map molecule index.

coot::molecule_t::histogram_info_t get_map_histogram(int imol, 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).

float get_suggested_initial_contour_level(int imol) const
Returns

the suggested initial contour level. Return -1 on not-a-map

bool is_EM_map(int imol) const
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

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.

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).

atom_radius is the atom radius (funnily enough). Use a negative number to mean “default”.

the 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

int flip_hand(int imol_map)

generate a new map which is the hand-flipped version of the input map.

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 of maps that are split by chain-id of the input imol

Returns

a vector 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.

void set_map_is_contoured_with_thread_pool(bool state)

set the state of the mode of the threading in map contouring

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.

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

void set_map_colour_saturation(int imol, float s)

set the map saturation

inline coot::util::sfcalc_genmap_stats_t get_latest_sfcalc_stats() const
r_factor_stats get_r_factor_stats()
std::string r_factor_stats_as_string(const r_factor_stats &rfs) const

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

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)

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 next rotamer (rotamer cycling is implicit if needed)

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

std::pair<int, unsigned int> delete_using_cid(int imol, const std::string &cid, const std::string &scope)

delete item

where scope is one of the strings: [“ATOM”,”WATER”,”RESIDUE”,”CHAIN”,”MOLECULE”, “LITERAL”]

Returns

1 on successful modification, 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

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 atom cid

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

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

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

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

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

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

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.

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

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

Returns

a first of 1 on success. Return a useful message in second if the addition did not work

int add_terminal_residue_directly_using_cid(int imol, const std::string &cid)

the cid is for an atom. This used to return a pair, but I removed it so that I could compile the binding.

Returns

an status.

int add_terminal_residue_directly_using_bucca_ml_growing_using_cid(int imol, const std::string &cid)

the cid is for an atom. buccaneer building

int add_terminal_residue_directly_using_bucca_ml_growing(int imol, const coot::residue_spec_t &spec)

buccaneer building, called by the above

int add_waters(int imol_model, int imol_map)

add waters, updating imol_model (of course)

Returns

the number of waters added on a success, -1 on failure.

int add_hydrogen_atoms(int imol_model)

add hydrogen atoms, updating imol_model (of course)

Returns

1 on success, 0 on failure.

int delete_hydrogen_atoms(int imol_model)

delete hydrogen atoms, updating imol_model (of course)

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

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

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

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

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

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 an atom CID

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 ligand

void eigen_flip_ligand_using_cid(int imol, const std::string &residue_cid)

eigen-flip ligand using CID

int mutate(int imol, const std::string &cid, const std::string &new_residue_type)

mutate residue

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

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.

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

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 (fun). “factor” might typically be say 0.9 or 1.1

coot::Cartesian get_molecule_centre(int imol) const

get molecule centre

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.

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 - 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. `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

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::molecule_t::moved_atom_t> &moved_atoms)

update the positions of the atoms in the residue

int new_positions_for_atoms_in_residues(int imol, const std::vector<coot::molecule_t::moved_residue_t> &moved_residues)

update the positions of the atoms in the residues

std::pair<int, std::vector<merge_molecule_results_info_t>> merge_molecules(int imol, const std::string &list_of_other_molecules)

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.

std::pair<int, std::vector<merge_molecule_results_info_t>> merge_molecules(int imol, std::vector<mmdb::Manager*> mols)

this is called by the above function and is useful for other non-api functions (such as add_compound()).

int cis_trans_convert(int imol, const std::string &atom_cid)

Convert a cis peptide to a trans or vice versa.

Returns

1 on a successful conversion.

int replace_fragment(int imol_base, int imol_reference, const std::string &atom_selection)

replace a fragment

i.e. replace the atoms of imol_base by those of the atom selection atom_selection in imol_reference (imol_base is the molecule that is modified).

Returns

the success status

int rigid_body_fit(int imol, const std::string &multi_cid, int imol_map)

Rigid-body fitting

`multi_cids” is a “||”-separated list of residues CIDs, e.g. “//A/12-52||//A/14-15||/B/56-66”

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

Returns

-1 on a conflict 1 on good. 0 on did nothing return also an information/error message

Coordinates Refinement

int refine_residues_using_atom_cid(int imol, const std::string &cid, const std::string &mode, int n_cycles)

refine the residues

mode is one of {SINGLE, TRIPLE, QUINTUPLE, HEPTUPLE, SPHERE, BIG_SPHERE, CHAIN, ALL};

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

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

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)
void fix_atom_selection_during_refinement(int imol, const std::string &atom_selection_cid)

fix atoms during refinement. Does nothing at the moment.

void add_target_position_restraint(int imol, const std::string &atom_cid, float pos_x, float pos_y, float pos_z)

add or update (if it has a pull restraint already)

void clear_target_position_restraint(int imol, const std::string &atom_cid)

clear target_position restraint

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

inline void set_use_rama_plot_restraints(bool state)

turn on or off rama restraints

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

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

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 Ramachandran plot restraints weight

inline float get_torsion_restraints_weight() const

get the Ramachandran plot restraints weight

Returns

the Ramachandran plot 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

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

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 thus 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). 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 coot::instanced_mesh_t

void clear_target_position_restraints(int imol)

clear any and all drag-atom target position restraints

void clear_refinement(int imol)

call this after molecule refinement has finished (say when the molecule molecule is accepted into the original molecule)

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

inline void set_refinement_geman_mcclure_alpha(float a)

set the refinement Geman-McClure alpha

inline float get_geman_mcclure_alpha() const

set the refinement Geman-McClure alpha

int generate_self_restraints(int imol, float local_dist_max)

generate GM self restraints for the whole molecule

Returns

nothing useful.

void generate_chain_self_restraints(int imol, float local_dist_max, const std::string &chain_id)

generate GM self restraints for the given chain

void generate_local_self_restraints(int imol, float local_dist_max, const std::string &residue_cids)

generate GM self restraints for the given residues. `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)

coot::instanced_mesh_t get_extra_restraints_mesh(int imol, int mode)

get the mesh for extra restraints currently mode is unused.

void read_extra_restraints(int imol, const std::string &file_name)

read extra restraints (e.g. from ProSMART)

void clear_extra_restraints(int imol)

clear the extra restraints

Coordinates Validation

coot::simple_mesh_t get_rotamer_dodecs(int imol)

get the rotamer dodecs for the model, not const because it regenerates the bonds.

Returns

a coot::simple_mesh_t

coot::instanced_mesh_t get_rotamer_dodecs_instanced(int imol)

get the rotamer dodecs for the model, not const because it regenerates the bonds.

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().

Returns

a coot::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

Returns

a vector of phi_psi_prob_t

coot::instanced_mesh_t contact_dots_for_ligand(int imol, const std::string &cid, unsigned int smoothness_factor) const

Recently (20230202) the smoothness factor has been added as an extra argument smoothness_factor is 1, 2 or 3 (3 is the most smooth).

Returns

the instanced mesh for the specified ligand

coot::instanced_mesh_t all_molecule_contact_dots(int imol, unsigned int smoothness_factor) const

Recently (20230202) the smoothness factor has been added as an extra argument smoothness_factor is 1, 2 or 3 (3 is the most smooth).

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)
Returns

a simple::molecule_t for the specified residue. this function is not const because we pass a pointer to the protein_geometry geom.

generic_3d_lines_bonds_box_t make_exportable_environment_bond_box(int imol, coot::residue_spec_t &spec)
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

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.

bool match_ligand_torsions(int imol_ligand, int imol_ref, const std::string &chain_id_ref, int resno_ref)

match ligand torsions - return 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 - return 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

bool match_ligand_torsions_and_position_using_cid(int imol_ligand, int imol_ref, const std::string &cid)

match ligand torsions and positions, different api

coot::atom_overlaps_dots_container_t get_overlap_dots(int imol)

not const because it can dynamically add dictionaries

coot::atom_overlaps_dots_container_t get_overlap_dots_for_ligand(int imol, const std::string &cid_ligand)

not const because it can dynamically add dictionaries

std::vector<coot::plain_atom_overlap_t> get_overlaps(int imol)

not const because it can dynamically add dictionaries

std::vector<coot::plain_atom_overlap_t> get_overlaps_for_ligand(int imol, const std::string &cid_ligand)

not const because it can dynamically add dictionaries

Coordinates and Map Validation

coot::validation_information_t density_fit_analysis(int imol_model, int imol_map) const

density fit validation information

Returns

a coot::validation_information_t

coot::validation_information_t density_correlation_analysis(int imol_model, int imol_map) const

density correlation validation information

Returns

a coot::validation_information_t

coot::validation_information_t rotamer_analysis(int imol_model) const

rotamer validation information

Returns

a coot::validation_information_t

coot::validation_information_t ramachandran_analysis(int imol_model) const

ramachandran validation information (formatted for a graph, not 3d)

Returns

a coot::validation_information_t

coot::validation_information_t ramachandran_analysis_for_chain(int imol_model, const std::string &chain_id) const

ramachandran validation information (formatted for a graph, not 3d) for a given chain in a given molecule 20230127-PE This function does not exist yet.

Returns

a coot::validation_information_t

coot::validation_information_t peptide_omega_analysis(int imol_model) const

peptide omega validation information

Returns

a coot::validation_information_t

std::vector<coot::molecule_t::interesting_place_t> get_interesting_places(int imol, const std::string &mode) const

get interesting places (does not work yet)

Returns

a vector of coot::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

Returns

a vector of coot::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

Returns

a vector of coot::validation_information_t

std::vector<coot::molecule_t::interesting_place_t> unmodelled_blobs(int imol_model, int imol_map) const

unmodelled blobs

Returns

a vector of coot::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, implicit OR

typical values for b_factor_lim is 60.0 typical values for outlier_sigma_level is 0.8 typical values for min_dist is 2.3 typical values for max_dist is 3.5

Returns

a vector of atom specifiers Use the string_user_data of the spec for the button label

std::vector<std::pair<double, double>> fourier_shell_correlation(int imol_map_1, int imol_map_2) const

Calculate the MMRRCC for the residues in the chain Multi Masked Residue Range Corellation Coefficient calculate the MMRRCC for the residues in the chain Multi Masked Residue Range Corellation Coefficient Fourier Shell Correlation (FSC) between maps

Returns

a vector or pairs of graph points (resolution, correlation). The resolution is in inverse Angstroms squared. An empty list is returned on failure

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

call this before calling connect_updating_maps().

int connect_updating_maps(int imol_model, int imol_with_data_info_attached, int imol_map_2fofc, int imol_map_fofc)

reset the rail_points (calls reset_the_rail_points()), updates the maps (using internal/clipper SFC) so, update your contour lines meshes after calling this function.

Returns

1 if the connection was successful.

void sfcalc_genmap(int imol_model, int imol_map_with_data_attached, int imol_updating_difference_map)

sfcalc and re-generate maps. This is a low-level function - generally one would use the updating maps method rather than this

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)

sfcalc and re-generate maps (low-level function). This functions uses bulk solvent.

Call this function after connecting maps for updating maps to set the initial R-factor and store the initial map flatness.

On failure to calculate SFS and generate the maps the returned r_factor in the returned stats will be set to -1.

Returns

a class of interesting statistics.

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.

Returns

success status

float get_density_at_position(int imol_map, float x, float y, float z) const

get density at position

Returns

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

This is a light-weight fetch, the values have already been computed, here were are merely copying them.

Returns

a vector the position where the differenc map has been flattened. The associated float value is the ammount that the map has been flattened.

std::string get_data_set_file_name(int imol) const

the stored data set file name

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 (obviously) the centre of the blob.

20221022-PE: 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”

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)

Ligand Fitting

I am not yet clear what extra cut-offs and flags need to be added here. 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.)

Returns

a vector 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

Returns

a vector 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 ligands (place-holder) multi_ligand_molecule_number_list is a colon-separated list of molecules, e.g. “2:3:4”

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” if n_trials is 0, then a sensible default value will be used. if translation_scale_factor 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” with a different interface - one that can use any atom selection (instead of just a ligand). As above, if n_trials is 0, then a sensible default value will be used. if translation_scale_factor 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 As above, if n_trials is 0, then a sensible default value will be used. if translation_scale_factor 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)

This is a ligand function, not really a ligand-fitting function.

It won’t work unless the dictionary for that ligand has been imported. The output renderings are not very good at the moment.

Except for unusual cases, imol will be IMOL_ENC_ANY (-666666)

dark_background_flag returns a representation suitable for rendering on a dark background (funnily enough).

This function is not const because it caches the svgs if it can.

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. i.e. those ligands that can be positioned without need for internal torsion angle manipulation.

tlc is the three-letter-code/compound-id

imol_dict is the molecule to which the ligand is attached (if any). Typically this will be IMOL_ENC_ANY (-666666).

imol_map is the molecule number of the map that will be used for fitting.

Returns

the success status, 1 or good, 0 for not good.

std::vector<coot::residue_spec_t> get_non_standard_residues_in_molecule(int imol) const
Returns

a vector of residue specifiers for the ligand residues - the residue name is encoded in the string_user_data data item of the residue specifier

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)

start a long-term job.

if n_seconds is 0, then run forever (or until interrupted)

void testing_stop_long_term_job()

stop the long-term job runnning (testing function)

inline ltj_stats_t testing_interrogate_long_term_job()

get the stats for the long-term job (testing function)

inline double get_contouring_time() const

get the time for conntouring in milliseconds

void set_max_number_of_threads(unsigned int n_threads)

set the maximum number of threads for both the thread pool and the vector of threads

void set_max_number_of_threads_in_thread_pool(unsigned int n_threads)

deprecated name for the above function

double test_the_threading(int n_threads)

get the time to run a test function in milliseconds

double test_launching_threads(unsigned int n_threads_per_batch, unsigned int n_batches) const
Returns

the time per batch in microseconds

double test_thread_pool_threads(unsigned int n_threads) const
Returns

time in microsections

inline coot::protein_geometry &get_geometry()
void make_mesh_for_map_contours_for_blender(int imol, float x, float y, float z, float level, float radius)
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)
void make_mesh_for_molecular_representation_for_blender(int imol, const std::string &cid, const std::string &colour_scheme, const std::string &style)
void make_mesh_for_gaussian_surface_for_blender(int imol, float sigma, float contour_level, float box_radius, float grid_scale, float b_factor)
void make_mesh_for_goodsell_style_for_blender(int imol, float colour_wheel_rotation_step, float saturation, float goodselliness)
std::vector<float> get_colour_table_for_blender(int imol)
std::vector<float> get_vertices_for_blender(int imol)
std::vector<int> get_triangles_for_blender(int imol)

Python functions

enum mesh_mode_t

Values:

enumerator UNKNOWN
enumerator SINGLE_COLOUR
enumerator MULTI_COLOUR
PyObject *simple_mesh_to_pythonic_mesh(const coot::simple_mesh_t &mesh, int mesh_mode)
PyObject *get_pythonic_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)
PyObject *get_pythonic_map_mesh(int imol, float x, float y, float z, float radius, float contour_level)
PyObject *get_pythonic_molecular_representation_mesh(int imol, const std::string &atom_selection, const std::string &colour_sheme, const std::string &style)
PyObject *get_pythonic_gaussian_surface_mesh(int imol, float sigma, float contour_level, float box_radius, float grid_scale, float fft_b_factor)
PyObject *get_pythonic_simple_molecule(int imol, const std::string &cid, bool include_hydrogen_atoms_flag)

make a “proper” simple molecule python class one day.

0: atom-name

1: atom-element

2: position (a list of 3 floats)

3: formal charge (an integer)

4: aromaticity flag (boolean)

Returns

a pair - the first of which (index 0) is the list of atoms, the second (index 1) is the list of bonds. An atom is a list:

Public Functions

inline explicit molecules_container_t(bool verbose = true)

the one and only constructor

~molecules_container_t()
inline void set_use_gemmi(bool state)

Set the state of using gemmi for coordinates parsing. The default is false.

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

enum moving_atoms_asc_t

Values:

enumerator NEW_COORDS_UNSET
enumerator NEW_COORDS_ADD
enumerator NEW_COORDS_REPLACE
enumerator NEW_COORDS_REPLACE_CHANGE_ALTCONF
enumerator NEW_COORDS_INSERT
enumerator NEW_COORDS_INSERT_CHANGE_ALTCONF

Private Functions

void set_updating_maps_need_an_update(int imol)
void update_updating_maps(int imol)

update the updating maps without generating a mesh

std::string adjust_refinement_residue_name(const std::string &resname) const
bool make_last_restraints(const std::vector<std::pair<bool, mmdb::Residue*>> &local_residues, const std::vector<mmdb::Link> &links, const coot::protein_geometry &geom, mmdb::Manager *mol_for_residue_selection, const std::vector<coot::atom_spec_t> &fixed_atom_specs, coot::restraint_usage_Flags flags, bool use_map_flag, const clipper::Xmap<float> *xmap_p)
coot::refinement_results_t refine_residues_vec(int imol, const std::vector<mmdb::Residue*> &residues, const std::string &alt_conf, mmdb::Manager *mol)
int find_serial_number_for_insert(int seqnum_new, const std::string &ins_code_for_new, mmdb::Chain *chain_p) const
atom_selection_container_t make_moving_atoms_asc(mmdb::Manager *residues_mol, const std::vector<mmdb::Residue*> &residues) const
std::pair<int, std::vector<std::string>> check_dictionary_for_residue_restraints(int imol, mmdb::PResidue *SelResidues, int nSelResidues)
std::pair<int, std::vector<std::string>> check_dictionary_for_residue_restraints(int imol, const std::vector<mmdb::Residue*> &residues)
std::pair<mmdb::Manager*, std::vector<mmdb::Residue*>> create_mmdbmanager_from_res_vector(const std::vector<mmdb::Residue*> &residues, int imol, mmdb::Manager *mol_in, std::string alt_conf)
coot::refinement_results_t generate_molecule_and_refine(int imol, const std::vector<mmdb::Residue*> &residues, const std::string &alt_conf, mmdb::Manager *mol, bool use_map_flag = true)
std::vector<std::pair<mmdb::Residue*, std::vector<coot::dict_torsion_restraint_t>>> make_rotamer_torsions(const std::vector<std::pair<bool, mmdb::Residue*>> &local_residues) const
int refine_direct(int imol, std::vector<mmdb::Residue*> rv, const std::string &alt_loc, int n_cycles)

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

double phi_psi_probability(const coot::util::phi_psi_t &phi_psi, const ramachandrans_container_t &rc) const
void read_standard_residues()

read the standard protein, RNA, and DNA dictionaries.

int install_model(const coot::molecule_t &m)
coot::graph_match_info_t overlap_ligands_internal(int imol_ligand, int imol_ref, const std::string &chain_id_ref, int resno_ref, bool apply_rtop_flag)
superpose_results_t superpose_with_atom_selection(atom_selection_container_t asc_ref, atom_selection_container_t asc_mov, int imol_mov, std::string moving_mol_name, std::string referennce_mol_name, bool move_copy_of_imol2_flag)
int valid_labels(const std::string &mtz_file_name, const std::string &f_col, const std::string &phi_col, const std::string &weight_col, int use_weights) const
inline void init()
void debug() const

Private Members

std::vector<coot::molecule_t> molecules
coot::protein_geometry geom
coot::rotamer_probability_tables rot_prob_tables
ramachandrans_container_t ramachandrans_container
bool draw_missing_residue_loops_flag
std::vector<rail_points_t> rail_point_history
updating_maps_info_f updating_maps_info
coot::util::sfcalc_genmap_stats_t latest_sfcalc_stats
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
bool show_timings
coot::restraints_container_t *last_restraints
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
bool particles_have_been_shown_already_for_this_round_flag
std::map<svg_store_key_t, std::string> ligand_svg_store
atom_selection_container_t standard_residues_asc
bool map_is_contoured_using_thread_pool_flag
double contouring_time

Private Static Functions

static void thread_for_refinement_loop_threaded()
static void get_restraints_lock(const std::string &calling_function_name)
static void release_restraints_lock(const std::string &calling_function_name)
static void all_atom_pulls_off()
static void atom_pull_off(const coot::atom_spec_t &spec)
static void atom_pulls_off(const std::vector<coot::atom_spec_t> &specs)

Private Static Attributes

static std::atomic<bool> on_going_updating_map_lock
static clipper::Xmap<float> *dummy_xmap
static ctpl::thread_pool static_thread_pool
static std::atomic<bool> restraints_lock
static std::string restraints_locking_function_name
static std::vector<atom_pull_info_t> atom_pulls
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 avaliable if the return value is empty.

std::string sigF_obs

sigF_obs column

std::string Rfree

R-Free column. There were not avaliable if the return value is empty.

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.

Public Functions

inline fit_ligand_info_t(int i, int c, int l)
inline fit_ligand_info_t()
inline float get_fitting_score() const
Returns

the fitting score

inline float get_cluster_volume() const

Public Members

int imol
int cluster_idx
int ligand_idx
float fitting_score
float cluster_volume
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

Public Members

unsigned int count
float function_value
std::chrono::time_point<std::chrono::high_resolution_clock> timer_start
std::chrono::time_point<std::chrono::high_resolution_clock> timer
class r_factor_stats

Public Members

float r_factor
float free_r_factor
int rail_points_total
int rail_points_new
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)
class updating_maps_info_f

Public Functions

inline updating_maps_info_f()

Public Members

bool maps_need_an_update
int imol_model
int imol_2fofc
int imol_fofc
int imol_with_data_info_attached

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)
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
namespace coot
namespace api
class vn_vertex
#include <vertex.hh>

a vertex with (only) postion and normal. Useful for instancing, perhaps

Public Functions

inline vn_vertex(const glm::vec3 &pos_in, const glm::vec3 &norm_in)
inline vn_vertex()

Public Members

glm::vec3 pos
glm::vec3 normal
class vnc_vertex
#include <vertex.hh>

a vertex with postion, normal and colour

Public Functions

inline vnc_vertex(const glm::vec3 &pos_in, const glm::vec3 &norm_in, const glm::vec4 &col_in)
inline explicit vnc_vertex(const vn_vertex &vn)
inline vnc_vertex(const vn_vertex &vn, const glm::vec4 &c)
inline vnc_vertex()

Public Members

glm::vec3 pos
glm::vec3 normal
glm::vec4 color
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

Public Members

glm::mat3 model_rotation_matrix
glm::vec3 model_translation
glm::vec3 pos
glm::vec3 normal
glm::vec4 colour
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()
namespace coot

Functions

simple_mesh_t instanced_mesh_to_simple_mesh(const instanced_mesh_t &im)
class instancing_data_type_A_t
#include <instancing.hh>

class for A type instancing data - this does not contain an orientation matrix

Public Functions

inline instancing_data_type_A_t(const glm::vec3 &position_in, const glm::vec4 &colour_in, const glm::vec3 &size_in)
inline instancing_data_type_A_t()

Public Members

glm::vec3 position
glm::vec4 colour
glm::vec3 size
class instancing_data_type_B_t
#include <instancing.hh>

class for B type instancing data - this does contain an orientation matrix

Public Functions

inline instancing_data_type_B_t(const glm::vec3 &position_in, const glm::vec4 &colour_in, const glm::vec3 &size_in, glm::mat4 &ori)
inline instancing_data_type_B_t()

Public Members

glm::vec3 position
glm::vec4 colour
glm::vec3 size
glm::mat4 orientation

the orientation matrix rotates the vector away from “z is up”

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<api::vn_vertex> vertices

vertices (containing positions and normals)

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

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

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)
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

Public Members

std::string chain_id
int res_no
std::string ins_code
std::string atom_name
std::string alt_conf
int int_user_data
float float_user_data
std::string string_user_data
int model_number
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)

Public Members

int model_number
std::string chain_id
int res_no
std::string ins_code
int int_user_data
float float_user_data
std::string string_user_data

Validation Information

namespace coot

Enums

enum graph_data_type

Values:

enumerator UNSET
enumerator DENSITY
enumerator DISTORTION
enumerator ENERGY
enumerator PROBABILITY
enumerator CORRELATION
enumerator LOG_PROBABILITY
enumerator TORSION_ANGLE
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_residue_validation_information(const residue_validation_information_t &rvi)

Public Members

std::string chain_id
std::vector<residue_validation_information_t> rviv
class validation_information_min_max_t

Public Functions

inline validation_information_min_max_t()
inline validation_information_min_max_t(const double &min_in, const double &max_in)

Public Members

bool is_set
double min
double max
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 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.

inline stats::single get_stats() const

get the stats for the data in this validation container

Public Members

std::string name
validation_information_min_max_t min_max
std::vector<chain_validation_information_t> cviv
enum graph_data_type type
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
namespace coot

Enums

enum map_stats_t

Values:

enumerator SIMPLE
enumerator WITH_KOLMOGOROV_SMIRNOV_DIFFERENCE_MAP_TEST
namespace util
class density_stats_info_t

Public Functions

inline density_stats_info_t()
inline void add(const double &v)
inline void add(const double &v, const double &weight)
inline std::pair<double, double> mean_and_variance() const

Public Members

double n
double sum_sq
double sum
double sum_weight
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 bool operator<(const density_correlation_stats_info_t &dcsi) const

Public Members

double n
double sum_xy
double sum_sqrd_x
double sum_sqrd_y
double sum_x
double sum_y
std::vector<double> density_values

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

Molecule Internals

Warning

doxygenfile: Cannot find file “coot_molecule.hh

Hydrogen Bonds

namespace moorhen
class h_bond_atom

Public Functions

inline h_bond_atom()

Public Members

int serial
float x
float y
float z
float charge
float occ
float b_iso
std::string element
std::string name
int model
std::string chain
int res_no
std::string residue_name
std::string altLoc
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

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

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 &current_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)
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)
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)
bool is_EM_map(const clipper::Xmap<float> &xmap)
std::vector<phitheta> make_phi_thetas(unsigned int n_pts)
float average_of_sample_map_at_sphere_points(clipper::Coord_orth &centre, 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)
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)
class map_molecule_centre_info_t
#include <coot-map-utils.hh>

map molecule centre

Public Functions

inline map_molecule_centre_info_t()

Public Members

bool success

success flag

clipper::Coord_orth updated_centre

new centre

float suggested_contour_level

suggested contour level

float suggested_radius

the suggested radius

double sum_of_densities

sum of densities - for whatever use that may be.

class map_stats_holder_helper_t

Public Functions

inline map_stats_holder_helper_t()
inline void add_xy(const double &x, const double &y)

Public Members

double sum_x
double sum_x_squared
double sum_y
double sum_y_squared
double sum_xy
int n
class map_fragment_info_t

Public Functions

map_fragment_info_t(const clipper::Xmap<float> &xmap, const clipper::Coord_orth &centre, float radius, bool centre_at_origin = false)
void unshift(clipper::Xmap<float> *xmap_p, const clipper::Coord_orth &centre)
void simple_origin_shift(const clipper::Xmap<float> &ip_xmap, const clipper::Coord_orth &centre, float radius)
clipper::Grid_map make_grid_map(const clipper::Xmap<float> &input_xmap, const clipper::Coord_orth &centre) const

Public Members

clipper::Xmap<float> xmap
clipper::Coord_grid offset

Private Functions

void init(const clipper::Xmap<float> &xmap, const clipper::Coord_orth &centre, float radius)
void init_making_map_centred_at_origin(const clipper::Xmap<float> &xmap, const clipper::Coord_orth &centre, float radius)

Private Members

float box_radius_a_internal
float box_radius_b_internal
float box_radius_c_internal
class simple_residue_triple_t

Public Functions

inline simple_residue_triple_t()
inline simple_residue_triple_t(mmdb::Residue *this_residue_in, mmdb::Residue *prev_residue_in, mmdb::Residue *next_residue_in, const std::string &alt_conf_in)

Public Members

mmdb::Residue *this_residue
mmdb::Residue *next_residue
mmdb::Residue *prev_residue
std::string alt_conf
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()

Public Members

mmdb::Residue *this_residue
mmdb::Residue *next_residue
mmdb::Residue *prev_residue
std::string alt_conf
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)
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

Public Members

double dist_sq
clipper::Xmap<float>::Map_reference_coord iw
float density
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 Types

enum [anonymous]

Values:

enumerator UNASSIGNED
enumerator TOO_LOW

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

Private Static Functions

static bool compare_density_values_map_refs(const std::pair<clipper::Xmap_base::Map_reference_index, float> &v1, const std::pair<clipper::Xmap_base::Map_reference_index, float> &v2)
static bool sort_segment_vec(const std::pair<int, int> &a, const std::pair<int, int> &b)
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)
namespace api
class cell_t
#include <api-cell.hh>

POD cell type for moorhen.

the angles are in radians.

Public Functions

inline cell_t()
inline cell_t(float a, float b, float c, float alpha, float beta, float gamma)

Public Members

float a

a

float b

b

float c

c

float alpha

alpha

float beta

beta

float gamma

gamma

bool is_set

is_set

Colour

namespace coot

Functions

std::ostream &operator<<(std::ostream &s, colour_t col)
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 average(const colour_t &other)
inline void brighter(float f)
inline glm::vec4 to_glm() const
inline colour_holder to_colour_holder() const

Public Members

std::vector<float> col

Private Functions

inline void init(float r, float g, float b)
std::vector<float> convert_to_hsv() const
void convert_from_hsv(const std::vector<float> &hsv)

Cartesian Coordinates

Enums

enum surface_face_data

Values:

enumerator NO_FACE
enumerator X_FACE
enumerator Y_FACE
enumerator Z_FACE
namespace coot

Functions

std::ostream &operator<<(std::ostream&, Cartesian pt)
std::ofstream &operator<<(std::ofstream&, Cartesian pt)
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

double cos_angle_btwn_vecs(const Cartesian &a, const Cartesian &b)
float dot_product(const Cartesian &a, const Cartesian &b)
Cartesian cross_product(const Cartesian &a, const Cartesian &b)
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()
inline Cartesian operator+(const Cartesian &in1) const
inline Cartesian operator-(const Cartesian &in1) const
Cartesian operator*(const float &f) const
void operator+=(const Cartesian&)
void operator-=(const Cartesian&)
void operator*=(float scale)
void operator/=(float scale)
Cartesian by_scalar(float scale)
void invert_z(void)
inline Cartesian &operator=(const Cartesian &in)
inline void unit_vector_yourself()
inline Cartesian unit() const
Cartesian rotate_about_vector(const coot::Cartesian &direction, const coot::Cartesian &origin, double angle) const
float distance_to_line(const Cartesian &front, const Cartesian &back) const
int within_box(const Cartesian &front, const Cartesian &back) const
Cartesian mid_point(const Cartesian &other) const
std::vector<Cartesian> third_points(const Cartesian &other) const

Public Static Functions

static float LineLength(const Cartesian &a, const Cartesian &b)
static double Angle(const Cartesian &a, const Cartesian &b, const Cartesian &c)
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)
static Cartesian position_by_torsion(const Cartesian &Atom_1, const Cartesian &Atom_2, const Cartesian &Atom_3, float theta_2, float torsion, float dist)
static Cartesian CrossProduct(const Cartesian &Atom_1, const Cartesian &Atom_2)
static float lengthsq(const Cartesian &c1, const Cartesian &c2)

Private Members

float x_
float y_
float z_

Friends

friend double cos_angle_btwn_vecs(const Cartesian &a, const Cartesian &b)
friend float dot_product(const Cartesian &a, const Cartesian &b)
friend Cartesian cross_product(const Cartesian &a, const Cartesian &b)
friend surface_face_data on_a_face(const Cartesian &a, const Cartesian &b)
friend std::ostream &operator<<(std::ostream&, Cartesian pt)
friend std::ofstream &operator<<(std::ofstream&, Cartesian pt)
class CartesianPair
#include <Cartesian.h>

Cartesian pair.

Public Functions

CartesianPair(const Cartesian &start_, const Cartesian &finish_)

constructor

CartesianPair(void)
inline void extentsBox(Cartesian centre, float dist)
inline const Cartesian &getStart() const

get the start Cartesian

inline const Cartesian &getFinish() const

get the end Cartesian

float amplitude() const

get the length of the vector between the coordinates

Private Members

Cartesian start
Cartesian finish

Friends

friend std::ostream &operator<<(std::ostream&, CartesianPair)
friend std::ofstream &operator<<(std::ofstream&, CartesianPair)
class CartesianPairInfo

Public Functions

inline CartesianPairInfo()

Public Members

CartesianPair *data
int size

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 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 void add(const single &s)
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
class pnorm

Public Functions

inline pnorm()
double erf(const double &z) const
inline double get(const double &x) const

Private Functions

inline void init()

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

Public Members

std::string symm_as_string

symmetry as symm

mmdb::mat44 mat

Private Members

int symm_no
int x_shift_
int y_shift_
int z_shift_

Friends

friend std::ostream &operator<<(std::ostream &s, const symm_trans_t &st)
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

Public Members

int us

unit cell steps

int vs
int ws

Friends

friend std::ostream &operator<<(std::ostream &s, Cell_Translation ct)
class molecule_extents_t

Public Functions

molecule_extents_t(atom_selection_container_t, float expansion_size)
~molecule_extents_t()
coot::Cartesian get_front() const
coot::Cartesian get_back() const
coot::Cartesian get_left() const
coot::Cartesian get_right() const
coot::Cartesian get_top() const
coot::Cartesian get_bottom() const
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 &centre_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
bool point_is_in_box(const coot::Cartesian &point, mmdb::PPAtom TransSel) 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

coot::Cartesian front
coot::Cartesian back
coot::Cartesian left
coot::Cartesian right
coot::Cartesian top
coot::Cartesian bottom
coot::Cartesian centre
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)
class SymmMatrix

Public Functions

SymmMatrix()
SymmMatrix(double **in_mat)
double **getMat() const
void add_unit_shift(int x, int y, int z)

Private Members

double mat[4][4]

Friends

friend std::ostream &operator<<(std::ostream&, SymmMatrix)
namespace coot
class coot_v4

Public Members

std::vector<float> v4
class coot_mat44

Public Functions

inline coot_mat44()
inline bool is_close_to_unit_matrix() const

Public Members

std::vector<coot_v4> m
class trans_selection_t

Public Functions

bool point_is_in_box(const coot::Cartesian &point) const

Public Members

Cartesian front
Cartesian back
Cartesian left
Cartesian right
Cartesian top
Cartesian bottom

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
class loc_table_item_t

Public Functions

inline loc_table_item_t(float i, float s, float loc)

Public Members

float invresolsq
float scale
float lack_of_closure

Indices and tables