#include <BaseMolecule.h>
Inheritance diagram for BaseMolecule:
Public Types | |
enum | dataset_flag { NODATA = 0x0000, INSERTION = 0x0001, OCCUPANCY = 0x0002, BFACTOR = 0x0004, MASS = 0x0008, CHARGE = 0x0010, RADIUS = 0x0020, ALTLOC = 0x0040, ATOMICNUMBER = 0x0080, BONDS = 0x0100, BONDORDERS = 0x0200, ANGLES = 0x0400, CTERMS = 0x0800, BONDTYPES = 0x1000, ANGLETYPES = 0x2000, ATOMFLAGS = 0x4000 } |
flags indicating what data was loaded or user-specified vs. guessed. More... | |
Public Methods | |
void | set_dataset_flag (int flag) |
void | unset_dataset_flag (int flag) |
int | test_dataset_flag (int flag) |
float * | radius () |
various standard per-atom fields. More... | |
float * | mass () |
float * | charge () |
float * | beta () |
float * | occupancy () |
unsigned char * | flags () |
void | set_radii_changed () |
void | get_radii_minmax (float &min, float &max) |
float * | bondorders () |
number of electron pairs, also fractional. More... | |
void | setbondorder (int atom, int bond, float order) |
float | getbondorder (int atom, int bond) |
int * | bondtypes () |
per atom/bond bond type list. More... | |
void | setbondtype (int atom, int bond, int type) |
set new bond type. More... | |
int | getbondtype (int atom, int bond) |
quest bond type. More... | |
int | has_structure () const |
has structure information already been loaded for this molecule? More... | |
void | clear_bonds (void) |
clear the entire bond list for all atoms. More... | |
int | count_bonds (void) |
return the number of unique bonds in the molecule. More... | |
BaseMolecule (int) | |
constructor takes molecule ID. More... | |
virtual | ~BaseMolecule (void) |
destructor. More... | |
int | init_atoms (int n) |
Try to set the number of atoms to n. n must be positive. May be called more than once with the same n. Return success. More... | |
void | find_bonds_from_timestep () |
compute molecule's bonds using distance bond search from 1st timestep. More... | |
void | find_unique_bonds_from_timestep () |
int | add_atoms (int natoms, const char *name, const char *type, int atomicnumber, const char *resname, int resid, const char *chainid, const char *segname, const char *insertion=(char *)"", const char *altloc="") |
add a new atom; return it's index. More... | |
int | add_bond (int, int, float, int, int=ATOMNORMAL) |
add a new bond from a to b; return total number of bonds added so far. More... | |
int | add_bond_dupcheck (int, int, float, int) |
add a bond after checking for duplicates. More... | |
void | clear_angles (void) |
clear list of angles and types. More... | |
int | num_angles () |
count angle list entries. More... | |
int | add_angle (int a, int b, int c, int type=-1) |
add an angle definition to existing list with optional angletype. checks for duplicates and returns angle index number. More... | |
int | set_angletype (int a, int type) |
set new angle type. More... | |
int | get_angletype (int a) |
query angle type. More... | |
void | clear_dihedrals (void) |
clear list of dihedrals and types. More... | |
int | num_dihedrals () |
count dihedral list entries. More... | |
int | add_dihedral (int a, int b, int c, int d, int type=-1) |
append a dihedral definition to existing list with optional dihedraltype. checks for duplicates. return dihedral index number. More... | |
int | set_dihedraltype (int d, int type) |
set new dihedral type. return type index number. More... | |
int | get_dihedraltype (int d) |
query dihedral type for number. More... | |
void | clear_impropers (void) |
clear list of impropers and types. More... | |
int | num_impropers () |
count improper list entries. More... | |
int | add_improper (int a, int b, int c, int d, int type=-1) |
append a improper definition to existing list with optional impropertype. checks for duplicates. return improper index number. More... | |
int | set_impropertype (int i, int type) |
set new improper type. returns type index number. More... | |
int | get_impropertype (int i) |
query improper type for number;. More... | |
void | clear_cterms () |
clear list of improper definitions. More... | |
int | num_cterms () |
count number of crossterm list entries. More... | |
void | add_cterm (int a, int b, int c, int d, int e, int f, int g, int h) |
add cross term definition. More... | |
void | analyze (void) |
find higher level constructs given the atom/bond information. More... | |
void | add_instance (Matrix4 &inst) |
int | num_instances (void) |
void | clear_instances (void) |
int | id (void) const |
return id number of this molecule. More... | |
const char * | molname () const |
return molecule name. More... | |
MolAtom * | atom (int n) |
return Nth atom. More... | |
Residue * | residue (int) |
return Nth residue. More... | |
Fragment * | fragment (int) |
return Nth fragment. More... | |
Residue * | atom_residue (int) |
return residue atom is in. More... | |
Fragment * | atom_fragment (int) |
return fragment atom is in. More... | |
void | add_volume_data (const char *name, const float *o, const float *xa, const float *ya, const float *za, int x, int y, int z, float *voldata, float *gradient=NULL, float *variance=NULL) |
add volumetric data to a molecule. More... | |
void | add_volume_data (const char *name, const double *o, const double *xa, const double *ya, const double *za, int x, int y, int z, float *voldata) |
void | remove_volume_data (int idx) |
remove a volumetric data object from a molecule. More... | |
int | num_volume_data () |
return number of volume data sets in molecule. More... | |
const VolumetricData * | get_volume_data (int) |
return requested data set. More... | |
VolumetricData * | modify_volume_data (int) |
return requested data set. More... | |
void | compute_volume_gradient (VolumetricData *) |
compute negated normalized volume gradient map. More... | |
int | find_atom_in_residue (int atomnameindex, int residue) |
find first occurance of an atom name in the residue, returns -3 not found. More... | |
int | find_atom_in_residue (const char *atomname, int residue) |
find first occurance of an atom name in the residue, returns -3 not found. More... | |
float | default_charge (const char *) |
return 'default' charge, mass, occupancy value and beta value. More... | |
float | default_mass (const char *) |
return 'default' charge, mass, occupancy value and beta value. More... | |
float | default_radius (const char *) |
return 'default' charge, mass, occupancy value and beta value. More... | |
float | default_occup (void) |
return 'default' charge, mass, occupancy value and beta value. More... | |
float | default_beta (void) |
return 'default' charge, mass, occupancy value and beta value. More... | |
Public Attributes | |
int | nAtoms |
number of atoms. More... | |
int | nResidues |
number of residues. More... | |
int | nWaters |
number of waters. More... | |
int | nSegments |
number of segments. More... | |
int | nFragments |
total number of fragments. More... | |
int | nProteinFragments |
number of protein fragments. More... | |
int | nNucleicFragments |
number of nucleic fragments. More... | |
NameList< int > | atomNames |
list of unique atom names. More... | |
NameList< int > | atomTypes |
list of unique atom types. More... | |
NameList< int > | resNames |
list of unique residue names. More... | |
NameList< int > | chainNames |
list of unique chain names. More... | |
NameList< int > | segNames |
list of unique segment names. More... | |
NameList< int > | altlocNames |
list of alternate location names. More... | |
ResizeArray< Residue * > | residueList |
residue connectivity list. More... | |
ResizeArray< Fragment * > | fragList |
list of connected residues. More... | |
ResizeArray< Fragment * > | pfragList |
list of connected protein residues this is a single chain from N to C. More... | |
ResizeArray< int > | pfragCyclic |
flag indicating cyclic fragment. More... | |
ResizeArray< Fragment * > | nfragList |
ditto for nuc; from 5' to 3'. More... | |
ResizeArray< int > | nfragCyclic |
flag indicating cyclic fragment. More... | |
ResizeArray< Matrix4 > | instances |
list of instance transform matrices. More... | |
ResizeArray< int > | angles |
list of angles. More... | |
ResizeArray< int > | dihedrals |
list of dihedrals. More... | |
ResizeArray< int > | impropers |
list of impropers. More... | |
ResizeArray< int > | cterms |
list of cross-terms. More... | |
NameList< float * > | extraflt |
optional floating point data. More... | |
NameList< int * > | extraint |
optional integer data. More... | |
NameList< unsigned char * > | extraflg |
optional bitwise flags. More... | |
ResizeArray< int > | angleTypes |
type of angle. More... | |
ResizeArray< int > | dihedralTypes |
type of dihedral. More... | |
ResizeArray< int > | improperTypes |
type of improper. More... | |
NameList< int > | bondTypeNames |
list of unique bond type names. More... | |
NameList< int > | angleTypeNames |
list of unique angle type names. More... | |
NameList< int > | dihedralTypeNames |
list of unique dihedral type names. More... | |
NameList< int > | improperTypeNames |
list of unique improper type names. More... | |
QMData * | qm_data |
int | datasetflags |
int | radii_minmax_need_update |
methods for querying the min/max atom radii, used by the the QuickSurf representation. More... | |
float | radii_min |
float | radii_max |
Protected Attributes | |
char * | moleculename |
name of the molcule. More... | |
int | need_find_bonds |
whether to compute bonds from the first timestep. More... | |
ResizeArray< VolumetricData * > | volumeList |
array of volume data sets. More... |
Definition at line 61 of file BaseMolecule.h.
|
flags indicating what data was loaded or user-specified vs. guessed.
Definition at line 140 of file BaseMolecule.h. |
|
constructor takes molecule ID.
Definition at line 50 of file BaseMolecule.C. References datasetflags, moleculename, nAtoms, need_find_bonds, nfragList, nFragments, nNucleicFragments, NODATA, nProteinFragments, nResidues, nSegments, NULL, nWaters, pfragList, qm_data, radii_max, radii_min, and radii_minmax_need_update. |
|
destructor.
Definition at line 84 of file BaseMolecule.C. References NameList< unsigned char * >::data, NameList< int * >::data, NameList< float * >::data, extraflg, extraflt, extraint, fragList, nfragList, NameList< unsigned char * >::num, NameList< int * >::num, NameList< float * >::num, ResizeArray< VolumetricData * >::num, ResizeArray< Fragment * >::num, ResizeArray< Residue * >::num, pfragList, qm_data, residueList, and volumeList. |
|
add an angle definition to existing list with optional angletype. checks for duplicates and returns angle index number.
Definition at line 384 of file BaseMolecule.C. References angles, ResizeArray< int >::append3, n, num_angles, and set_angletype. Referenced by VMDApp::molecule_from_selection_list, molinfo_set, MolFilePlugin::read_optional_structure, MolFilePlugin::read_structure, and topo_add_angle. |
|
add a new atom; return it's index.
Definition at line 210 of file BaseMolecule.C. References NameList< int >::add_name, MolAtom::altlocindex, altlocNames, atom, MolAtom::atomicnumber, atomNames, atomTypes, MolAtom::chainindex, chainNames, MolAtom::init, MolAtom::nameindex, nAtoms, NameList< int >::num, MolAtom::resnameindex, resNames, MolAtom::segnameindex, segNames, and MolAtom::typeindex. Referenced by VMDApp::molecule_from_selection_list, VMDApp::molecule_new, and MolFilePlugin::read_structure. |
|
add a new bond from a to b; return total number of bonds added so far.
Definition at line 305 of file BaseMolecule.C. References MolAtom::add_bond, atom, MolAtom::bonds, MAXBONDERRORS, nAtoms, setbondorder, and setbondtype. Referenced by add_bond_dupcheck, VMDApp::molecule_from_selection_list, PickModeAddBond::pick_molecule_end, MolFilePlugin::read_optional_structure, MolFilePlugin::read_structure, and vmd_bond_search. |
|
add a bond after checking for duplicates.
Definition at line 363 of file BaseMolecule.C. References add_bond, atom, MolAtom::bonds, MolAtom::bondTo, and nAtoms. Referenced by MolFilePlugin::read_structure, topo_add_bond, and vmd_bond_search. |
|
add cross term definition.
Definition at line 450 of file BaseMolecule.h. References ResizeArray::append. Referenced by VMDApp::molecule_from_selection_list, molinfo_set, MolFilePlugin::read_optional_structure, and MolFilePlugin::read_structure. |
|
append a dihedral definition to existing list with optional dihedraltype. checks for duplicates. return dihedral index number.
Definition at line 423 of file BaseMolecule.C. References ResizeArray< int >::append4, dihedrals, n, num_dihedrals, and set_dihedraltype. Referenced by VMDApp::molecule_from_selection_list, molinfo_set, MolFilePlugin::read_optional_structure, MolFilePlugin::read_structure, and topo_add_dihed. |
|
append a improper definition to existing list with optional impropertype. checks for duplicates. return improper index number.
Definition at line 461 of file BaseMolecule.C. References ResizeArray< int >::append4, impropers, n, num_impropers, and set_impropertype. Referenced by VMDApp::molecule_from_selection_list, molinfo_set, MolFilePlugin::read_optional_structure, MolFilePlugin::read_structure, and topo_add_improp. |
|
Definition at line 463 of file BaseMolecule.h. References ResizeArray::append. Referenced by VMDApp::molecule_add_instance. |
|
Definition at line 2588 of file BaseMolecule.C. References ResizeArray< VolumetricData * >::append, VolumetricData::compute_volume_gradient, VolumetricData::datarange, VolumetricData::gridsize, VolumetricData::name, volumeList, and z. |
|
add volumetric data to a molecule.
Definition at line 2531 of file BaseMolecule.C. References ResizeArray< VolumetricData * >::append, VolumetricData::compute_volume_gradient, VolumetricData::datarange, VolumetricData::gridsize, VolumetricData::name, VolumetricData::set_volume_gradient, volumeList, and z. Referenced by VMDApp::molecule_add_volumetric, and MolFilePlugin::read_volumetric. |
|
find higher level constructs given the atom/bond information.
Definition at line 726 of file BaseMolecule.C. References angleTypeNames, atom, ATOMHYDROGEN, atomNames, ATOMNORMAL, MolAtom::atomType, bondTypeNames, ResizeArray< int >::clear, ResizeArray< Fragment * >::clear, ResizeArray< Residue * >::clear, count_bonds, dihedralTypeNames, fragList, improperTypeNames, IS_HYDROGEN, NameList< int >::name, MolAtom::nameindex, nAtoms, need_find_bonds, nfragCyclic, nfragList, nFragments, nNucleicFragments, nProteinFragments, nResidues, nSegments, NULL, ResizeArray< Fragment * >::num, NameList< int >::num, num_angles, num_cterms, num_dihedrals, num_impropers, nWaters, pfragCyclic, pfragList, residueList, MolAtom::residueType, RESNUCLEIC, RESPROTEIN, wkf_timer_create, wkf_timer_destroy, wkf_timer_start, wkf_timer_timenow, and wkf_timerhandle. Referenced by VMDApp::molecule_from_selection_list, VMDApp::molecule_load, VMDApp::molecule_reanalyze, and MolFilePlugin::read_optional_structure. |
|
|
return fragment atom is in.
Definition at line 590 of file BaseMolecule.C. References atom, fragment, Residue::fragment, n, NULL, residue, and MolAtom::uniq_resid. Referenced by PickModeForceFragment::pick_molecule_move, PickModeForceFragment::pick_molecule_start, PickModeMoveFragment::rotate, and PickModeMoveFragment::translate. |
|
return residue atom is in.
Definition at line 579 of file BaseMolecule.C. References atom, n, NULL, residue, and MolAtom::uniq_resid. Referenced by PickModeForceResidue::pick_molecule_move, PickModeForceResidue::pick_molecule_start, PickModeMoveResidue::rotate, and PickModeMoveResidue::translate. |
|
Definition at line 173 of file BaseMolecule.h. References NameList::data. Referenced by GeometryMol::atom_formatted_name, AtomColor::find, colvarproxy_vmd::load_atoms, colvarproxy_vmd::load_coords, VMDApp::molecule_from_selection_list, VMDApp::molecule_new, MolFilePlugin::read_optional_structure, MolFilePlugin::read_structure, and MolFilePlugin::write_structure. |
|
number of electron pairs, also fractional.
Definition at line 211 of file BaseMolecule.h. References NameList::data. |
|
Definition at line 216 of file BaseMolecule.h. References NameList::data. |
|
|
clear list of angles and types.
Definition at line 397 of file BaseMolecule.h. References ResizeArray::clear. Referenced by molinfo_set, MolFilePlugin::read_optional_structure, and MolFilePlugin::read_structure. |
|
clear the entire bond list for all atoms.
Definition at line 714 of file BaseMolecule.C. References MolAtom::bonds, and nAtoms. Referenced by DrawMolecule::recalc_bonds, and topo_del_all_bonds. |
|
clear list of improper definitions.
Definition at line 444 of file BaseMolecule.h. References ResizeArray::clear. Referenced by molinfo_set. |
|
clear list of dihedrals and types.
Definition at line 413 of file BaseMolecule.h. References ResizeArray::clear. Referenced by molinfo_set, MolFilePlugin::read_optional_structure, and MolFilePlugin::read_structure. |
|
clear list of impropers and types.
Definition at line 429 of file BaseMolecule.h. References ResizeArray::clear. Referenced by molinfo_set, MolFilePlugin::read_optional_structure, and MolFilePlugin::read_structure. |
|
Definition at line 465 of file BaseMolecule.h. References ResizeArray::clear. Referenced by VMDApp::molecule_delete_all_instances. |
|
compute negated normalized volume gradient map.
|
|
return the number of unique bonds in the molecule.
Definition at line 696 of file BaseMolecule.C. References MolAtom::bonds, MolAtom::bondTo, and nAtoms. Referenced by analyze, DrawMolecule::recalc_bonds, and topo_del_all_bonds. |
|
return 'default' charge, mass, occupancy value and beta value.
Definition at line 504 of file BaseMolecule.h. Referenced by VMDApp::molecule_new, and MolFilePlugin::read_structure. |
|
return 'default' charge, mass, occupancy value and beta value.
Definition at line 687 of file BaseMolecule.C. Referenced by VMDApp::molecule_new, and MolFilePlugin::read_structure. |
|
return 'default' charge, mass, occupancy value and beta value.
Definition at line 625 of file BaseMolecule.C. Referenced by VMDApp::molecule_new, and MolFilePlugin::read_structure. |
|
return 'default' charge, mass, occupancy value and beta value.
Definition at line 503 of file BaseMolecule.h. Referenced by VMDApp::molecule_new, and MolFilePlugin::read_structure. |
|
return 'default' charge, mass, occupancy value and beta value.
Definition at line 600 of file BaseMolecule.C. Referenced by VMDApp::molecule_new, and MolFilePlugin::read_structure. |
|
find first occurance of an atom name in the residue, returns -3 not found.
Definition at line 1465 of file BaseMolecule.C. References atomNames, find_atom_in_residue, residue, and NameList< int >::typecode. |
|
find first occurance of an atom name in the residue, returns -3 not found.
Definition at line 486 of file BaseMolecule.h. References atom, Residue::atoms, atoms, MolAtom::nameindex, ResizeArray::num, and num. Referenced by find_atom_in_residue. |
|
compute molecule's bonds using distance bond search from 1st timestep.
Definition at line 379 of file BaseMolecule.h. References need_find_bonds. Referenced by MolFilePlugin::read_structure. |
|
Definition at line 380 of file BaseMolecule.h. References need_find_bonds. Referenced by MolFilePlugin::read_structure. |
|
Definition at line 177 of file BaseMolecule.h. References NameList::data. |
|
return Nth fragment.
Definition at line 572 of file BaseMolecule.C. Referenced by atom_fragment, and compoundmapper. |
|
query angle type.
Definition at line 415 of file BaseMolecule.C. References angleTypes, and ResizeArray< int >::num. Referenced by molinfo_get, topo_get_angle, and MolFilePlugin::write_structure. |
|
query dihedral type for number.
Definition at line 454 of file BaseMolecule.C. References dihedralTypes, and ResizeArray< int >::num. Referenced by molinfo_get, topo_get_dihed, and MolFilePlugin::write_structure. |
|
query improper type for number;.
Definition at line 492 of file BaseMolecule.C. References improperTypes, and ResizeArray< int >::num. Referenced by molinfo_get, topo_get_impro, and MolFilePlugin::write_structure. |
|
Definition at line 186 of file BaseMolecule.h. References NameList::data, minmax_1fv_aligned, nAtoms, NULL, radii_max, radii_min, and radii_minmax_need_update. Referenced by calc_density_bounds, QuickSurf::calc_surf, measure_sasa, measure_sasa_perresidue, and measure_sasa_thread. |
|
return requested data set.
Definition at line 2637 of file BaseMolecule.C. References id, and modify_volume_data. Referenced by atomsel_gridindex_array, atomsel_interp_volume_array, atomsel_volume_array, colvarproxy_vmd::compute_volmap, mdff_cc, SaveTrajectoryFltkMenu::molchooser_activate_selection, py_mol_get_volumetric, segment_volume, vmd_measure_volinterior, vmd_volmap_compare, GraphicsFltkMenu::volindex_update, and MolFilePlugin::write_volumetric. |
|
Definition at line 521 of file BaseMolecule.C. References atom, NameList< float * >::data, extraflt, MAXATOMBONDS, and NULL. Referenced by access_tcl_atomsel, assign_atoms, topo_get_bond, and MolFilePlugin::write_structure. |
|
quest bond type.
Definition at line 553 of file BaseMolecule.C. References atom, NameList< int * >::data, extraint, MAXATOMBONDS, and NULL. Referenced by access_tcl_atomsel, topo_get_bond, and MolFilePlugin::write_structure. |
|
has structure information already been loaded for this molecule?
Definition at line 225 of file BaseMolecule.h. Referenced by VMDApp::molecule_load, and MolFilePlugin::write_structure. |
|
|
Try to set the number of atoms to n. n must be positive. May be called more than once with the same n. Return success.
Definition at line 141 of file BaseMolecule.C. References NameList< unsigned char * >::add_name, NameList< int * >::add_name, NameList< float * >::add_name, NameList< unsigned char * >::data, NameList< int * >::data, NameList< float * >::data, data, extraflg, extraflt, extraint, n, nAtoms, NULL, NameList< unsigned char * >::num, NameList< int * >::num, and NameList< float * >::num. Referenced by VMDApp::molecule_from_selection_list, VMDApp::molecule_load, VMDApp::molecule_new, and MolFilePlugin::read_structure. |
|
|
return requested data set.
Definition at line 2631 of file BaseMolecule.C. References id, NULL, ResizeArray< VolumetricData * >::num, and volumeList. Referenced by density_add, density_average, density_binmask, density_clamp, density_com, density_correlate, density_crop, density_downsample, density_histogram, density_info, density_mdff_potential, density_move, density_moveto, density_multiply, density_range, density_sadd, density_save, density_sigma, density_smooth, density_smult, density_subtract, density_supersample, density_trim, AtomColor::find, fit, get_volume_data, mask, text_cmd_mol, GraphicsFltkMenu::update_molchooser, and vmd_measure_volinterior. |
|
return molecule name.
Definition at line 472 of file BaseMolecule.h. References moleculename. Referenced by GeometryMol::atom_formatted_name, VMDApp::molecule_from_selection_list, VMDApp::molecule_name, molinfo_get, and MolBrowser::update. |
|
count angle list entries.
Definition at line 400 of file BaseMolecule.h. References ResizeArray::num. Referenced by add_angle, analyze, VMDApp::molecule_from_selection_list, molinfo_get, set_angletype, topo_del_all_angles, topo_del_angle, topo_get_angle, and MolFilePlugin::write_structure. |
|
count number of crossterm list entries.
Definition at line 447 of file BaseMolecule.h. References ResizeArray::num. Referenced by analyze, VMDApp::molecule_from_selection_list, molinfo_get, and MolFilePlugin::write_structure. |
|
count dihedral list entries.
Definition at line 416 of file BaseMolecule.h. References ResizeArray::num. Referenced by add_dihedral, analyze, VMDApp::molecule_from_selection_list, molinfo_get, set_dihedraltype, topo_del_all_dihed, topo_del_dihed, topo_get_dihed, and MolFilePlugin::write_structure. |
|
count improper list entries.
Definition at line 432 of file BaseMolecule.h. References ResizeArray::num. Referenced by add_improper, analyze, VMDApp::molecule_from_selection_list, molinfo_get, set_impropertype, topo_del_all_impropers, topo_del_improper, topo_get_impro, and MolFilePlugin::write_structure. |
|
Definition at line 464 of file BaseMolecule.h. References ResizeArray::num. Referenced by VMDApp::molecule_num_instances. |
|
return number of volume data sets in molecule.
Definition at line 2627 of file BaseMolecule.C. References ResizeArray< VolumetricData * >::num, and volumeList. Referenced by VolMapCreateILS::add_map_to_molecule, colvarproxy_vmd::check_volmap_by_id, AtomColor::find, SaveTrajectoryFltkMenu::molchooser_activate_selection, VMDApp::molecule_load, VMDApp::molecule_savetrajectory, molinfo_get, py_mol_del_volumetric, py_mol_get_volumetric, py_mol_num_volumetric, text_cmd_mol, MolBrowser::update, GraphicsFltkMenu::update_molchooser, vmd_volmap_compare, GraphicsFltkMenu::volindex_update, VolMapCreateILS::write_map, and MolFilePlugin::write_volumetric. |
|
Definition at line 174 of file BaseMolecule.h. References NameList::data. Referenced by GeometryMol::atom_formatted_name, AtomColor::find, colvarproxy_vmd::load_atoms, colvarproxy_vmd::load_coords, VMDApp::molecule_from_selection_list, VMDApp::molecule_new, MolFilePlugin::read_optional_structure, MolFilePlugin::read_structure, and MolFilePlugin::write_structure. |
|
various standard per-atom fields.
Definition at line 170 of file BaseMolecule.h. References NameList::data. Referenced by calc_cc, mdff_cc, mdff_sim, VMDApp::molecule_from_selection_list, VMDApp::molecule_new, MolFilePlugin::read_optional_structure, MolFilePlugin::read_structure, vmd_bond_search, vmd_measure_volinterior, and MolFilePlugin::write_structure. |
|
remove a volumetric data object from a molecule.
Definition at line 2619 of file BaseMolecule.C. References ResizeArray< VolumetricData * >::num, ResizeArray< VolumetricData * >::remove, and volumeList. Referenced by py_mol_del_volumetric, and text_cmd_mol. |
|
return Nth residue.
Definition at line 566 of file BaseMolecule.C. References n, and residueList. Referenced by atom_fragment, atom_residue, colvarproxy_vmd::check_atom_id, compoundmapper, AtomColor::find, find_atom_in_residue, PickModeForceFragment::pick_molecule_move, PickModeForceFragment::pick_molecule_start, PickModeMoveFragment::rotate, and PickModeMoveFragment::translate. |
|
set new angle type.
Definition at line 397 of file BaseMolecule.C. References ANGLETYPES, angleTypes, ResizeArray< int >::appendN, ResizeArray< int >::num, num_angles, and set_dataset_flag. Referenced by add_angle. |
|
Definition at line 159 of file BaseMolecule.h. References datasetflags. Referenced by access_tcl_atomsel, VMDApp::molecule_set_dataset_flag, molinfo_set, PickModeAddBond::pick_molecule_end, MolFilePlugin::read_optional_structure, MolFilePlugin::read_structure, set_angletype, set_dihedraltype, set_impropertype, topo_add_angle, topo_add_bond, topo_add_dihed, topo_add_improp, topo_del_angle, topo_del_bond, topo_del_dihed, and topo_del_improper. |
|
set new dihedral type. return type index number.
Definition at line 436 of file BaseMolecule.C. References ANGLETYPES, ResizeArray< int >::appendN, dihedralTypes, ResizeArray< int >::num, num_dihedrals, and set_dataset_flag. Referenced by add_dihedral. |
|
set new improper type. returns type index number.
Definition at line 474 of file BaseMolecule.C. References ANGLETYPES, ResizeArray< int >::appendN, improperTypes, ResizeArray< int >::num, num_impropers, and set_dataset_flag. Referenced by add_improper. |
|
Definition at line 184 of file BaseMolecule.h. References radii_minmax_need_update. Referenced by MolFilePlugin::read_optional_structure, and MolFilePlugin::read_structure. |
|
Definition at line 501 of file BaseMolecule.C. References NameList< float * >::add_name, atom, NameList< float * >::data, extraflt, MAXATOMBONDS, nAtoms, and NULL. Referenced by access_tcl_atomsel, and add_bond. |
|
set new bond type.
Definition at line 533 of file BaseMolecule.C. References NameList< int * >::add_name, atom, NameList< int * >::data, extraint, MAXATOMBONDS, nAtoms, and NULL. Referenced by access_tcl_atomsel, and add_bond. |
|
Definition at line 165 of file BaseMolecule.h. References datasetflags. Referenced by MolFilePlugin::write_structure. |
|
Definition at line 162 of file BaseMolecule.h. References datasetflags. Referenced by VMDApp::molecule_set_dataset_flag. |
|
list of alternate location names.
Definition at line 79 of file BaseMolecule.h. Referenced by add_atoms, MoleculeList::add_color_names, GeometryMol::atom_formatted_name, AtomColor::find, VMDApp::molecule_from_selection_list, vmd_bond_search, and MolFilePlugin::write_structure. |
|
list of angles.
Definition at line 113 of file BaseMolecule.h. Referenced by add_angle, VMDApp::molecule_from_selection_list, molinfo_get, topo_del_all_angles, topo_del_angle, topo_get_angle, and MolFilePlugin::write_structure. |
|
list of unique angle type names.
Definition at line 129 of file BaseMolecule.h. Referenced by analyze, molinfo_get, molinfo_set, MolFilePlugin::read_optional_structure, MolFilePlugin::read_structure, topo_add_angle, topo_angletypes, topo_get_angle, and MolFilePlugin::write_structure. |
|
type of angle.
Definition at line 125 of file BaseMolecule.h. Referenced by get_angletype, molinfo_get, set_angletype, topo_del_all_angles, and topo_del_angle. |
|
list of unique atom names.
Definition at line 74 of file BaseMolecule.h. Referenced by add_atoms, MoleculeList::add_color_names, analyze, GeometryMol::atom_formatted_name, GeometryMol::atom_short_name, atomsel_name, colvarproxy_vmd::check_atom_id, AtomColor::find, find_atom_in_residue, VMDApp::molecule_from_selection_list, print_atom_info, vmd_bond_search, write_ss_input_pdb, and MolFilePlugin::write_structure. |
|
list of unique atom types.
Definition at line 75 of file BaseMolecule.h. Referenced by add_atoms, MoleculeList::add_color_names, GeometryMol::atom_formatted_name, atomsel_type, AtomColor::find, VMDApp::molecule_from_selection_list, print_atom_info, and MolFilePlugin::write_structure. |
|
list of unique bond type names.
Definition at line 128 of file BaseMolecule.h. Referenced by access_tcl_atomsel, analyze, MolFilePlugin::read_optional_structure, MolFilePlugin::read_structure, topo_add_bond, topo_bondtypes, topo_get_bond, and MolFilePlugin::write_structure. |
|
list of unique chain names.
Definition at line 77 of file BaseMolecule.h. Referenced by add_atoms, MoleculeList::add_color_names, GeometryMol::atom_formatted_name, AtomColor::find, VMDApp::molecule_from_selection_list, print_atom_info, write_ss_input_pdb, and MolFilePlugin::write_structure. |
|
list of cross-terms.
Definition at line 116 of file BaseMolecule.h. Referenced by VMDApp::molecule_from_selection_list, molinfo_get, and MolFilePlugin::write_structure. |
|
Definition at line 158 of file BaseMolecule.h. Referenced by BaseMolecule, VMDApp::molecule_from_selection_list, set_dataset_flag, test_dataset_flag, and unset_dataset_flag. |
|
list of dihedrals.
Definition at line 114 of file BaseMolecule.h. Referenced by add_dihedral, VMDApp::molecule_from_selection_list, molinfo_get, topo_del_all_dihed, topo_del_dihed, topo_get_dihed, and MolFilePlugin::write_structure. |
|
list of unique dihedral type names.
Definition at line 130 of file BaseMolecule.h. Referenced by analyze, molinfo_get, molinfo_set, MolFilePlugin::read_optional_structure, MolFilePlugin::read_structure, topo_add_dihed, topo_dihetypes, topo_get_dihed, and MolFilePlugin::write_structure. |
|
type of dihedral.
Definition at line 126 of file BaseMolecule.h. Referenced by get_dihedraltype, molinfo_get, set_dihedraltype, topo_del_all_dihed, and topo_del_dihed. |
|
optional bitwise flags.
Definition at line 122 of file BaseMolecule.h. Referenced by init_atoms, and ~BaseMolecule. |
|
optional floating point data.
Definition at line 120 of file BaseMolecule.h. Referenced by build_xyzr_from_sel, calc_density_bounds, VolMapCreate::calculate_max_radius, VolMapCreateCoulombPotentialMSM::compute_frame, VolMapCreateCoulombPotential::compute_frame, VolMapCreateDistance::compute_frame, VolMapCreateOccupancy::compute_frame, VolMapCreateDensity::compute_frame, getbondorder, init_atoms, measure_sasa_thread, minmax, sasa, sasa_perresidue, setbondorder, topo_del_bond, vmd_measure_minmax, vmd_measure_sasa, vmd_measure_sasaperresidue, vmd_volmap_ils, and ~BaseMolecule. |
|
optional integer data.
Definition at line 121 of file BaseMolecule.h. Referenced by getbondtype, init_atoms, setbondtype, topo_del_bond, and ~BaseMolecule. |
|
list of connected residues.
Definition at line 82 of file BaseMolecule.h. Referenced by analyze, fragment, and ~BaseMolecule. |
|
list of impropers.
Definition at line 115 of file BaseMolecule.h. Referenced by add_improper, VMDApp::molecule_from_selection_list, molinfo_get, topo_del_all_impropers, topo_del_improper, topo_get_impro, and MolFilePlugin::write_structure. |
|
list of unique improper type names.
Definition at line 131 of file BaseMolecule.h. Referenced by analyze, molinfo_get, molinfo_set, MolFilePlugin::read_optional_structure, MolFilePlugin::read_structure, topo_add_improp, topo_get_impro, topo_imptypes, and MolFilePlugin::write_structure. |
|
type of improper.
Definition at line 127 of file BaseMolecule.h. Referenced by get_impropertype, molinfo_get, set_impropertype, topo_del_all_impropers, and topo_del_improper. |
|
list of instance transform matrices.
Definition at line 109 of file BaseMolecule.h. |
|
name of the molcule.
Definition at line 356 of file BaseMolecule.h. Referenced by BaseMolecule, Molecule::Molecule, molname, Molecule::rename, and Molecule::~Molecule. |
|
|
whether to compute bonds from the first timestep.
Definition at line 357 of file BaseMolecule.h. Referenced by analyze, DrawMolecule::append_frame, BaseMolecule, DrawMolecule::DrawMolecule, find_bonds_from_timestep, and find_unique_bonds_from_timestep. |
|
flag indicating cyclic fragment.
Definition at line 94 of file BaseMolecule.h. Referenced by analyze. |
|
ditto for nuc; from 5' to 3'.
Definition at line 93 of file BaseMolecule.h. Referenced by analyze, BaseMolecule, and ~BaseMolecule. |
|
total number of fragments.
Definition at line 70 of file BaseMolecule.h. Referenced by analyze, BaseMolecule, compoundmapper, AtomColor::find, and VMDApp::molecule_from_selection_list. |
|
number of nucleic fragments.
Definition at line 72 of file BaseMolecule.h. Referenced by analyze, BaseMolecule, and VMDApp::molecule_from_selection_list. |
|
number of protein fragments.
Definition at line 71 of file BaseMolecule.h. Referenced by analyze, BaseMolecule, and VMDApp::molecule_from_selection_list. |
|
number of residues.
Definition at line 67 of file BaseMolecule.h. Referenced by analyze, BaseMolecule, colvarproxy_vmd::check_atom_id, compoundmapper, and VMDApp::molecule_from_selection_list. |
|
number of segments.
Definition at line 69 of file BaseMolecule.h. Referenced by analyze, and BaseMolecule. |
|
number of waters.
Definition at line 68 of file BaseMolecule.h. Referenced by analyze, and BaseMolecule. |
|
flag indicating cyclic fragment.
Definition at line 86 of file BaseMolecule.h. Referenced by analyze. |
|
list of connected protein residues this is a single chain from N to C.
Definition at line 84 of file BaseMolecule.h. Referenced by analyze, BaseMolecule, and ~BaseMolecule. |
|
Definition at line 134 of file BaseMolecule.h. Referenced by BaseMolecule, VMDApp::molecule_orblocalize, molinfo_get, MolFilePlugin::read_qm_data, GraphicsFltkRepOrbital::regen_excitationlist, GraphicsFltkRepOrbital::regen_orbitallist, GraphicsFltkRepOrbital::regen_wavefunctypes, and ~BaseMolecule. |
|
Definition at line 183 of file BaseMolecule.h. Referenced by BaseMolecule, and get_radii_minmax. |
|
Definition at line 182 of file BaseMolecule.h. Referenced by BaseMolecule, and get_radii_minmax. |
|
methods for querying the min/max atom radii, used by the the QuickSurf representation.
Definition at line 181 of file BaseMolecule.h. Referenced by BaseMolecule, get_radii_minmax, and set_radii_changed. |
|
residue connectivity list.
Definition at line 81 of file BaseMolecule.h. Referenced by analyze, QuickSurf::calc_surf, read_dssp_record, read_stride_record, residue, and ~BaseMolecule. |
|
list of unique residue names.
Definition at line 76 of file BaseMolecule.h. Referenced by add_atoms, MoleculeList::add_color_names, GeometryMol::atom_formatted_name, GeometryMol::atom_short_name, AtomColor::find, VMDApp::molecule_from_selection_list, print_atom_info, write_ss_input_pdb, and MolFilePlugin::write_structure. |
|
list of unique segment names.
Definition at line 78 of file BaseMolecule.h. Referenced by add_atoms, MoleculeList::add_color_names, GeometryMol::atom_formatted_name, colvarproxy_vmd::check_atom_id, AtomColor::find, VMDApp::molecule_from_selection_list, print_atom_info, and MolFilePlugin::write_structure. |
|
array of volume data sets.
Definition at line 525 of file BaseMolecule.h. Referenced by add_volume_data, DrawMolecule::cov, modify_volume_data, num_volume_data, remove_volume_data, DrawMolecule::scale_factor, and ~BaseMolecule. |