![]() |
Bsoft 2.1.4
Bernard's software package
|
Library routines used for atomic coordinates. More...
#include "rwmolecule.h"
#include "rwatomprop.h"
#include "rwresprop.h"
#include "mol_util.h"
#include "seq_util.h"
#include "linked_list.h"
#include "Matrix.h"
#include "utilities.h"
#include <map>
Functions | |
int | molgroup_Bfactors (Bmolgroup *molgroup) |
Calculates the average Bfactors for residues in molecules. More... | |
int | molgroup_print_sequence (Bmolgroup *molgroup) |
Prints the sequence of a molecule as a single letter code string. More... | |
int | molgroup_set_radius (Bmolgroup *molgroup, double radius) |
Shifts all coordinates to a given radius. More... | |
int | molgroup_rename (Bmolgroup *molgroup, char first_name) |
Renames molecules. More... | |
int | molgroup_residue_renumber (Bmolgroup *molgroup, int first) |
Renumbers residues. More... | |
int | molgroup_atom_renumber (Bmolgroup *molgroup, int first) |
Renumbers all atoms, disregarding molecule distinction. More... | |
int | molgroup_coor_reset_occupancy (Bmolgroup *molgroup, int range_first, int range_last, double occupancy) |
Sets occupancies to new values.. More... | |
double | molgroup_weight_from_atoms (Bmolgroup *molgroup) |
Calculates the weight of each molecule. More... | |
double | molgroup_weight_from_sequence (Bmolgroup *molgroup) |
Calculates the weight of each molecule. More... | |
Vector3< double > | molgroup_center_of_mass (Bmolgroup *molgroup) |
Calculates the center of mass of a molecule group. More... | |
Vector3< double > | molgroup_selected_center_of_mass (Bmolgroup *molgroup) |
Calculates the center of mass of a molecule group. More... | |
Vector3< double > | molgroup_show_center_of_mass (Bmolgroup *molgroup) |
Calculates the center of mass of a molecule group. More... | |
Vector3< double > | mol_center_of_mass (Bmolecule *mol) |
Calculates the center of mass of a molecule. More... | |
Vector3< double > | mol_show_center_of_mass (Bmolecule *mol) |
Calculates the center of mass of a molecule. More... | |
Vector3< double > | molgroup_principal_axes (Bmolgroup *molgroup, Vector3< double > *eigenvec) |
Calculates the principal axes of a molecule group. More... | |
Vector3< double > | mol_principal_axes (Bmolecule *mol, Vector3< double > *eigenvec) |
Calculates the principal axes of a molecule. More... | |
double | molgroup_sphericity (Bmolgroup *molgroup, double da) |
Calculates how close the coordinate set is to a spherical shape. More... | |
double | molgroup_density (Bmolgroup *molgroup) |
Calculates the number and mass density of a molecular system. More... | |
double | molgroup_volume (Bmolgroup *molgroup, Bstring ¶mfile, int wrap) |
Calculates the volume of a molecular system. More... | |
int | molgroup_composition (Bmolgroup *molgroup, Bstring ¶mfile) |
Calculates the elemental composition of a molecular system. More... | |
JSvalue | molgroup_elements (Bmolgroup *molgroup) |
Calculates the elemental composition. More... | |
JSvalue | molgroup_elements (Bmolgroup *molgroup, Bstring ¶mfile) |
Calculates the elemental composition from residue compositions. More... | |
JSvalue | protein_composition_default (double mass) |
Default protein composition adjusted by mass. More... | |
int | molgroup_radial_density (Bmolgroup *molgroup, double interval, double cutoff, int wrap) |
Calculates a radial distribution function. More... | |
Latom ** | molgroup_atom_mesh_lists (Bmolgroup *molgroup, Vector3< int > size, Vector3< double > sampling) |
Generates lists of atoms based on a grid. More... | |
vector< Batom * > | atom_get_array (Bmolgroup *molgroup) |
Generates an array of pointers to atom structures. More... | |
Variables | |
int | verbose |
Library routines used for atomic coordinates.
Generates an array of pointers to atom structures.
*molgroup | molecule group structure to be modified. |
Calculates the center of mass of a molecule.
*mol | molecule structure. |
Calculates the principal axes of a molecule.
*mol | molecule structure. |
*eigenvec | eigen vectors (can be NULL). |
Calculates the center of mass of a molecule.
*mol | molecule structure. |
Latom ** molgroup_atom_mesh_lists | ( | Bmolgroup * | molgroup, |
Vector3< int > | size, | ||
Vector3< double > | sampling | ||
) |
Generates lists of atoms based on a grid.
The coordinates must be positive to fit into a grid starting at {0,0,0}.
*molgroup | molecule group structure to be modified. |
size | size of grid. |
sampling | spacing in each dimension. |
int molgroup_atom_renumber | ( | Bmolgroup * | molgroup, |
int | first | ||
) |
Renumbers all atoms, disregarding molecule distinction.
*molgroup | molecule group structure. |
first | number of first atom. |
int molgroup_Bfactors | ( | Bmolgroup * | molgroup | ) |
Calculates the average Bfactors for residues in molecules.
*molgroup | molecule group structure. |
Calculates the center of mass of a molecule group.
*molgroup | molecule group structure. |
Calculates the elemental composition of a molecular system.
*molgroup | molecule group structure. |
¶mfile | atomic parameter file. |
int molgroup_coor_reset_occupancy | ( | Bmolgroup * | molgroup, |
int | range_first, | ||
int | range_last, | ||
double | occupancy | ||
) |
Sets occupancies to new values..
*molgroup | molecule group structure. |
range_first | first residue in range to change. |
range_last | last residue in range to change. |
occupancy | value to set occupancy to. |
double molgroup_density | ( | Bmolgroup * | molgroup | ) |
Calculates the number and mass density of a molecular system.
*molgroup | molecule group structure. |
Requirement: The box size must be correctly specified in the molecular group structure. The number of molecules, residues and atoms are calculated and the density for each number reported. The total mass is calculated and the mass density reported as Da/A^3 and g/ml.
Calculates the elemental composition.
*molgroup | molecule group. |
Calculates the elemental composition from residue compositions.
*molgroup | molecule group. |
¶mfile | residue parameter file. |
Calculates the principal axes of a molecule group.
*molgroup | molecule group structure. |
*eigenvec | eigen vectors (can be NULL). |
int molgroup_print_sequence | ( | Bmolgroup * | molgroup | ) |
Prints the sequence of a molecule as a single letter code string.
*molgroup | molecule group structure. |
int molgroup_radial_density | ( | Bmolgroup * | molgroup, |
double | interval, | ||
double | cutoff, | ||
int | wrap | ||
) |
Calculates a radial distribution function.
*molgroup | molecule group. |
interval | interval between bins. |
cutoff | distance cutoff. |
wrap | flag to use periodic boundaries. |
The atoms are not distinguished by type or properties.
int molgroup_rename | ( | Bmolgroup * | molgroup, |
char | first_name | ||
) |
Renames molecules.
*molgroup | molecule group structure. |
first_name | letter of first molecule. |
Letters are assigned to the molecules in sequence starting from the given first letter. If it extends beyond 'Z', it restarts at 'A'.
int molgroup_residue_renumber | ( | Bmolgroup * | molgroup, |
int | first | ||
) |
Renumbers residues.
*molgroup | molecule group structure. |
first | number of first residue. |
Calculates the center of mass of a molecule group.
*molgroup | molecule group structure. |
int molgroup_set_radius | ( | Bmolgroup * | molgroup, |
double | radius | ||
) |
Shifts all coordinates to a given radius.
*molgroup | molecule group structure. |
radius | spherical radius. |
Calculates the center of mass of a molecule group.
*molgroup | molecule group structure. |
double molgroup_sphericity | ( | Bmolgroup * | molgroup, |
double | da | ||
) |
Calculates how close the coordinate set is to a spherical shape.
*molgroup | molecule group structure. |
da | angular step size. (in radians) |
For each solid angle, the algorithm finds the coordinates with the maximum distance from the center-of-mass.
Calculates the volume of a molecular system.
*molgroup | molecule group structure. |
¶mfile | atomic parameter file. |
wrap | flag to turn on wrapping. |
The volume occupied by all the atoms is estimated from the footprint under the Van der Waals volume of each atom.
double molgroup_weight_from_atoms | ( | Bmolgroup * | molgroup | ) |
Calculates the weight of each molecule.
*molgroup | molecule group structure. |
double molgroup_weight_from_sequence | ( | Bmolgroup * | molgroup | ) |
Calculates the weight of each molecule.
*molgroup | molecule group structure. |
JSvalue protein_composition_default | ( | double | mass | ) |
Default protein composition adjusted by mass.
mass | molecular weight. |
|
extern |