Class Avogadro::Core::Molecule#

class Molecule#

The Molecule class represents a chemical molecule.

Subclassed by Molecule

Unnamed Group

const CustomElementMap &customElementMap() const#

A map of custom element atomic numbers to string identifiers. These ids can be used to override the generic custom element names returned by the Elements class, and should be somewhat meaningful to the user.

Note

Custom element atomic numbers lie between CustomElementMin and CustomElementMax.

void setCustomElementMap(const CustomElementMap &map)#

A map of custom element atomic numbers to string identifiers. These ids can be used to override the generic custom element names returned by the Elements class, and should be somewhat meaningful to the user.

Note

Custom element atomic numbers lie between CustomElementMin and CustomElementMax.

Unnamed Group

virtual BondType addBond(Index atom1, Index atom2, unsigned char order = 1)#

Create a new bond in the molecule.

Parameters:
  • atom1 – The first atom in the bond.

  • atom2 – The second atom in the bond.

  • order – The bond order.

Returns:

The new bond object. Will be invalid if atom1 or atom2 does not exist.

virtual BondType addBond(const AtomType &atom1, const AtomType &atom2, unsigned char order = 1)#

Create a new bond in the molecule.

Parameters:
  • atom1 – The first atom in the bond.

  • atom2 – The second atom in the bond.

  • order – The bond order.

Returns:

The new bond object. Will be invalid if atom1 or atom2 does not exist.

Unnamed Group

virtual bool removeBond(Index atom1, Index atom2)#

Remove the specified bond.

Parameters:
  • atom1 – One atom in the bond.

  • atom2 – The other atom in the bond.

Returns:

True on success, false if the bond was not found. This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

virtual bool removeBond(const AtomType &atom1, const AtomType &atom2)#

Remove the specified bond.

Parameters:
  • atom1 – One atom in the bond.

  • atom2 – The other atom in the bond.

Returns:

True on success, false if the bond was not found. This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

Unnamed Group

Array<BondType> bonds(const AtomType &a)#

Get all bonds to a.

Returns:

A vector of bonds to the supplied atom a.

Array<BondType> bonds(Index a)#

Get all bonds to a.

Returns:

A vector of bonds to the supplied atom a.

Array<const BondType*> bonds(Index a) const#

Get all bonds to a.

Returns:

A vector of bonds to the supplied atom a.

Unnamed Group

void setUnitCell(UnitCell *uc)#

The unit cell for this molecule. May be nullptr for non-periodic structures.

inline UnitCell *unitCell()#

The unit cell for this molecule. May be nullptr for non-periodic structures.

inline const UnitCell *unitCell() const#

The unit cell for this molecule. May be nullptr for non-periodic structures.

Unnamed Group

inline void setHallNumber(unsigned short hallNumber)#

The space group for this molecule. It is updated after every space group operation.

inline unsigned short hallNumber() const#

The space group for this molecule. It is updated after every space group operation.

Public Types

typedef Atom AtomType#

Typedef for Atom class.

typedef Bond BondType#

Typedef for Bond class.

typedef std::map<unsigned char, std::string> CustomElementMap#

Type for custom element map.

typedef std::bitset<element_count> ElementMask#

Type for element masks (e.g., does this molecule contain certain elements)

Public Functions

Molecule()#

Creates a new, empty molecule.

Molecule(const Molecule &other)#

Copy constructor

Molecule(Molecule &&other) noexcept#

Move constructor

Molecule &operator=(const Molecule &other)#

Assignment operator

Molecule &operator=(Molecule &&other) noexcept#

Move assignment operator

virtual ~Molecule()#

Destroys the molecule object.

void readProperties(const Molecule &other)#

Adds the properties from the supplied molecule to this molecule. Does not otherwise modify atoms / bonds / residues, etc.

void setData(const std::string &name, const Variant &value)#

Sets the data value with name to value.

Variant data(const std::string &name) const#
Returns:

the data value for name.

bool hasData(const std::string &name) const#
Returns:

true if the molecule has data with the given key, false otherwise.

void setDataMap(const VariantMap &map)#

Set the molecule’s variant data to the entries in map.

const VariantMap &dataMap() const#
Returns:

the molecule’s variant data.

VariantMap &dataMap()#

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

MatrixX spectra(const std::string &name) const#
Returns:

a specific spectra entry

void setSpectra(const std::string &name, const MatrixX &value)#

Sets the spectra value with name to value.

std::set<std::string> spectraTypes() const#
Returns:

the list of available spectra

void setPartialCharges(const std::string &type, const MatrixX &value)#

Sets atomic partial charges with type to value.

MatrixX partialCharges(const std::string &type) const#
Returns:

the atomic partial charges of type type

std::set<std::string> partialChargeTypes() const#
Returns:

the types of partial charges available stored with this molecule.

Array<AtomHybridization> &hybridizations()#
Returns:

a vector of hybridizations for the atoms in the molecule.

const Array<AtomHybridization> &hybridizations() const#

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

inline AtomHybridization hybridization(Index atomId) const#

Get the hybridization for the requested atom.

Parameters:

atomId – The index of the atom.

Returns:

The hybridization of the atom indexed at atomId, or 0 if atomId is invalid.

inline bool setHybridizations(const Core::Array<AtomHybridization> &hybs)#

Replace the current array of hybridizations.

Parameters:

hybs – The new hybridization array. Must be of length atomCount().

Returns:

True on success, false otherwise.

inline bool setHybridization(Index atomId, AtomHybridization hybridization)#

Set the hybridization of a single atom.

Parameters:
  • atomId – The index of the atom to modify.

  • hybridization – The new hybridization.

Returns:

True on success, false otherwise.

Array<signed char> &formalCharges()#
Returns:

a vector of formal charges for the atoms in the molecule.

const Array<signed char> &formalCharges() const#

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

signed char totalCharge() const#

Get the total charge on the molecule. The method will first check to see if a total charge has been set. If not, it will calculate the total charge from the formal charges (if set). If neither has been set, it will assume the total charge is zero.

Returns:

The total charge of the molecule.

char totalSpinMultiplicity() const#

Get the total spin multiplicity of the molecule. The method will first check to see if a total spin has been set. If not, it will either suggest a singlet if an even number of electrons are present, or a doublet if an odd number of electrons are present.

Returns:

The total spin multiplicity of the molecule.

inline signed char formalCharge(Index atomId) const#

Get the formal charge for the requested atom.

Parameters:

atomId – The index of the atom.

Returns:

The formal charge of the atom indexed at atomId, or 0 if atomId is invalid.

inline bool setFormalCharges(const Core::Array<signed char> &charges)#

Replace the current array of formal charges.

Parameters:

charges – The new formal charge array. Must be of length atomCount().

Returns:

True on success, false otherwise.

inline bool setFormalCharge(Index atomId, signed char charge)#

Set the formal charge of a single atom.

Parameters:
  • atomId – The index of the atom to modify.

  • charge – The new formal charge.

Returns:

True on success, false otherwise.

Array<Vector3ub> &colors()#
Returns:

a vector of colors for the atoms in the moleucle.

const Array<Vector3ub> &colors() const#

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

inline Vector3ub color(Index atomId) const#

Get the color for the requested atom.

Parameters:

atomId – The index of the atom.

Returns:

The color of the atom indexed at atomId, or (0,0,0) if atomId is invalid. If no color is set for the given atomId, the default color for the atomic number of the atomId is returned.

inline bool setColors(const Core::Array<Vector3ub> &colors)#

Replace the current array of colors.

Parameters:

colors – The new color array. Must be of length atomCount().

Returns:

True on success, false otherwise.

inline bool setColor(Index atomId, Vector3ub color)#

Set the color of a single atom.

Parameters:
  • atomId – The index of the atom to modify.

  • color – The new color.

Returns:

True on success, false otherwise.

inline bool setLayer(Index atomId, size_t layer)#
inline size_t layer(Index atomId) const#
const Array<Vector2> &atomPositions2d() const#
Returns:

a vector of 2d atom positions for the atoms in the molecule.

Array<Vector2> &atomPositions2d()#

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

inline Vector2 atomPosition2d(Index atomId) const#

Get the 2D position of a single atom.

Parameters:

atomId – The index of the atom.

Returns:

The position of the atom, or Vector3::Zero() if no position information has been set.

inline bool setAtomPositions2d(const Core::Array<Vector2> &pos)#

Replace the current array of 2D atomic coordinates.

Parameters:

pos – The new coordinate array. Must be of length atomCount().

Returns:

True on success, false otherwise.

inline bool setAtomPosition2d(Index atomId, const Vector2 &pos)#

Set the 2D position of a single atom.

Parameters:
  • atomId – The index of the atom to modify.

  • pos – The new position of the atom.

Returns:

True on success, false otherwise.

const Array<Vector3> &atomPositions3d() const#
Returns:

a vector of 3d atom positions for the atoms in the molecule.

Array<Vector3> &atomPositions3d()#

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

inline Vector3 atomPosition3d(Index atomId) const#

Get the 3D position of a single atom.

Parameters:

atomId – The index of the atom.

Returns:

The position of the atom, or Vector3::Zero() if no position information has been set.

inline bool setAtomPositions3d(const Core::Array<Vector3> &pos)#

Replace the current array of 3D atomic coordinates.

Parameters:

pos – The new coordinate array. Must be of length atomCount().

Returns:

True on success, false otherwise.

inline bool setAtomPosition3d(Index atomId, const Vector3 &pos)#

Set the 3D position of a single atom.

Parameters:
  • atomId – The index of the atom to modify.

  • pos – The new position of the atom.

Returns:

True on success, false otherwise.

inline std::string atomLabel(Index atomId) const#
Parameters:

atomId – The index of the atom.

Returns:

Any custom label for the requested atom.

inline bool setAtomLabel(Index atomId, const std::string &label)#

Set the custom label of a single atom.

Parameters:
  • atomId – The index of the atom to modify.

  • label – The new label of the atom.

Returns:

True on success, false otherwise.

inline const Core::Array<std::string> atomLabels() const#
inline bool setAtomLabels(const Core::Array<std::string> &label)#

Set all the atom labels in the molecule.

Parameters:

label – The new label array. Must be of length atomCount().

Returns:

True on success, false otherwise.

inline void setAtomSelected(Index atomId, bool selected)#

Set whether the specified atom is selected or not.

inline bool atomSelected(Index atomId) const#

Query whether the supplied atom index has been selected.

inline bool isSelectionEmpty() const#
Returns:

whether the selection is empty or not

inline const ElementMask elements() const#
Returns:

the elements currently in this molecule

virtual AtomType addAtom(unsigned char atomicNumber)#

Adds an atom to the molecule.

AtomType addAtom(unsigned char atomicNumber, Vector3 position3d)#
virtual bool removeAtom(Index index)#

Remove the specified atom from the molecule.

Parameters:

index – The index of the atom to be removed.

Returns:

True on success, false if the atom was not found.

virtual bool removeAtom(const AtomType &atom)#

Remove the specified atom from the molecule.

Parameters:

atom – The atom to be removed.

Returns:

True on success, false if the atom was not found. This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

virtual void clearAtoms()#

Remove all atoms from the molecule.

AtomType atom(Index index) const#
Returns:

the atom at index in the molecule.

virtual bool removeBond(Index index)#

Remove the specified bond.

Parameters:

index – The index of the bond to be removed.

Returns:

True on success, false if the bond was not found.

virtual bool removeBond(const BondType &bond)#

Remove the specified bond.

Parameters:

bond – The bond to be removed.

Returns:

True on success, false if the bond was not found. This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

virtual void clearBonds()#

Remove all bonds from the molecule.

BondType bond(Index index) const#
Returns:

the bond at index in the molecule.

BondType bond(const AtomType &a, const AtomType &b) const#
Returns:

the bond between atoms a and b.

BondType bond(Index atomId1, Index atomId2) const#
Returns:

the bond between atomId1 and atomId2.

inline std::string bondLabel(Index bondIndex) const#
Parameters:

bondIndex – The index of the bond.

Returns:

Any custom label for the requested bond.

inline bool setBondLabel(Index bondIndex, const std::string &label)#

Set the custom label of a single bond.

Parameters:
  • bondIndex – The index of the bond to modify.

  • label – The new label of the bond.

Returns:

True on success, false otherwise.

inline const Core::Array<std::string> bondLabels() const#
inline bool setBondLabels(const Core::Array<std::string> &label)#

Set all the atom labels in the molecule.

Parameters:

label – The new label array. Must be of length atomCount().

Returns:

True on success, false otherwise.

Mesh *addMesh()#

Add a mesh to the molecule.

Returns:

The mesh object added to the molecule.

Mesh *mesh(Index index)#
const Mesh *mesh(Index index) const#
inline Index meshCount() const#
void clearMeshes()#
Cube *addCube()#

Add a cube to the molecule.

Returns:

The cube object added to the molecule.

Cube *cube(Index index)#
const Cube *cube(Index index) const#
inline Index cubeCount() const#
void clearCubes()#
inline std::vector<Cube*> cubes()#

Get the cubes vector set (if present) for the molecule.

Returns:

The cube vector for the molecule

inline const std::vector<Cube*> cubes() const#
std::string formula(const std::string &delimiter = "", int showCountsOver = 1) const#
Parameters:
  • delimiter – Delimiter to insert between tokens, defaults to none.

  • showCountsOver – Show atom counts above this (defaults to 1).

Returns:

the chemical formula of the molecule.

double mass() const#
Returns:

The mass of the molecule obtained by summing constituent atomic masses.

Vector3 centerOfGeometry() const#
Returns:

The center of geometry of the molecule obtained by summing the coordinates of the atoms.

Vector3 centerOfMass() const#
Returns:

The center of mass of the molecule obtained by summing the coordinates of the atoms weighted by mass.

double radius() const#
Returns:

The minimum radius of a sphere centered on centerOfGeometry containing all the centers of the atoms.

std::pair<Vector3, Vector3> bestFitPlane() const#
Returns:

The (centroid, normal vector) pair of the best-fit plane of the atoms of the molecule.

inline void setBasisSet(BasisSet *basis)#

Set the basis set for the molecule, note that the molecule takes ownership of the object.

inline BasisSet *basisSet()#
Returns:

the basis set (if present) for the molecule.

inline const BasisSet *basisSet() const#
Array<double> vibrationFrequencies() const#
void setVibrationFrequencies(const Array<double> &freq)#
Array<double> vibrationIRIntensities() const#
void setVibrationIRIntensities(const Array<double> &intensities)#
Array<double> vibrationRamanIntensities() const#
void setVibrationRamanIntensities(const Array<double> &intensities)#
Array<Vector3> vibrationLx(int mode) const#
void setVibrationLx(const Array<Array<Vector3>> &lx)#
void perceiveBondsSimple(const double tolerance = 0.45, const double minDistance = 0.32)#

Perceives bonds in the molecule based on the 3D coordinates of the atoms. atoms are considered bonded if within the sum of radii plus a small tolerance.

Parameters:
  • tolerance – The calculation tolerance.

  • minDistance – = atoms closer than the square of this are ignored

void perceiveBondsFromResidueData()#

Perceives bonds in the molecule based on preset residue data.

Use this if you have residue data available (e.g., reading PDB or MMTF files) Otherwise consider

See also

perceiveBondsSimple and

See also

perceiveBondOrders

void perceiveBondOrders()#
void perceiveSubstitutedCations()#

Perceives all-carbon-substituted onium ions of nitrogen, oxygen, phosphorus, sulfur, arsenic and selenium.

int coordinate3dCount()#
bool setCoordinate3d(int coord)#
Array<Vector3> coordinate3d(int index) const#
bool setCoordinate3d(const Array<Vector3> &coords, int index)#
void clearCoordinate3d()#

Clear coordinate sets (except the default set)

bool setTimeStep(double timestep, int index)#

Timestep property is used when molecular dynamics trajectories are read

double timeStep(int index, bool &status)#
const Array<Vector3> &forceVectors() const#
Returns:

a vector of forces for the atoms in the molecule.

Array<Vector3> &forceVectors()#

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

inline Vector3 forceVector(Index atomId) const#

Get the force of a single atom.

Parameters:

atomId – The index of the atom.

Returns:

The force vector of the atom, or Vector3::Zero() if no force information has been set.

inline bool setForceVectors(const Core::Array<Vector3> &forces)#

Replace the current array of force vectors.

Parameters:

forces – The new coordinate array. Must be of length atomCount().

Returns:

True on success, false otherwise.

inline bool setForceVector(Index atomId, const Vector3 &force)#

Set the 3D position of a single atom.

Parameters:
  • atomId – The index of the atom to modify.

  • force – The new position of the atom.

Returns:

True on success, false otherwise.

Residue &addResidue(std::string &name, Index &number, char &id)#
void addResidue(Residue &residue)#
Residue &residue(Index index)#
inline Array<Residue> &residues()#
inline const Array<Residue> &residues() const#
Index residueCount() const#
Returns:

The number of residues in the molecule.

inline Index atomCount() const#
Returns:

The number of atoms in the molecule.

Index atomCount(unsigned char atomicNumber) const#

Get the number of atoms in the molecule that match atomicNumber.

Parameters:

atomicNumber – The atomic number to match on.

Returns:

The number of atoms with the supplied atomic number.

inline Index bondCount() const#
Returns:

the number of bonds in the molecule.

inline std::pair<Index, Index> bondPair(Index bondId) const#

Get the set of bonded atoms corresponding to bondId.

Parameters:

bondId – The index of the requested bond.

Returns:

The bonded atom pair, represented as a pair of atom indices.

bool setBondPairs(const Array<std::pair<Index, Index>> &pairs)#

Replace the current array of bonded atoms.

Note

The bonded atoms are represented as a pair of bond indices.

Note

If needed, the elements in pairs will be modified to ensure that the first atom index is less than the second.

Parameters:

pairs – The array.

Returns:

True on success, false on failure.

inline const Array<std::pair<Index, Index>> &bondPairs() const#
Returns:

a vector of pairs of atom indices of the bonds in the molecule.

inline const Array<unsigned char> &bondOrders() const#
Returns:

a vector of the bond orders for the bonds in the molecule.

inline const Graph &graph() const#
Returns:

the graph for the molecule.

inline const Array<unsigned char> &atomicNumbers() const#
Returns:

a vector of atomic numbers for the atoms in the molecule.

inline unsigned char atomicNumber(Index atomId) const#

Get the atomic number for the requested atom.

Parameters:

atomId – The index of the atom.

Returns:

The atomic number of the atom indexed at atomId, or Avogadro::InvalidElement if atomId is invalid.

bool setBondOrders(const Array<unsigned char> &orders)#

Replace the current array of bond orders.

Parameters:

orders – The new array.

Returns:

True on success, false on failure.

bool setBondOrder(Index bondId, unsigned char order)#

Set the order of a bond in the molecule.

Parameters:
  • bondId – The bond’s index.

  • order – The new order of the bond.

Returns:

True on success, false on failure.

bool hasCustomElements() const#

Note

Custom element atomic numbers lie between CustomElementMin and CustomElementMax.

Returns:

True if custom elements exist in the molecule.

bool setBondPair(Index bondId, const std::pair<Index, Index> &pair)#

Set the bonded atoms for a bond.

Note

If needed, pair will be modified to ensure that the first atom index is less than the second.

Parameters:
  • bondId – The bond to modify.

  • pair – The new bond pair.

Returns:

True on success, false otherwise.

unsigned char bondOrder(Index bondId) const#

Get the order of a bond.

Parameters:

bondId – The id of the bond.

Returns:

The bond order.

bool setAtomicNumbers(const Core::Array<unsigned char> &nums)#

Replace the current array of atomic numbers.

Parameters:

nums – The new atomic number array. Must be of length atomCount().

Returns:

True on success, false otherwise.

bool setAtomicNumber(Index atomId, unsigned char atomicNumber)#

Set the atomic number of a single atom.

Parameters:
  • atomId – The index of the atom to modify.

  • atomicNumber – The new atomic number.

Returns:

True on success, false otherwise.

void setFrozenAtom(Index atomId, bool frozen)#

Freeze or unfreeze an atom for optimization

bool frozenAtom(Index atomId) const#

Get the frozen status of an atom

void setFrozenAtomAxis(Index atomId, int axis, bool frozen)#

Freeze or unfreeze X, Y, or Z coordinate of an atom for optimization

Parameters:
  • atomId – The index of the atom to modify.

  • axis – The axis to freeze (0, 1, or 2 for X, Y, or Z)

  • frozen – True to freeze, false to unfreeze

inline Eigen::VectorXd frozenAtomMask() const#
std::map<unsigned char, size_t> composition() const#
Returns:

a map of components and count.

Array<std::pair<Index, Index>> getAtomBonds(Index index) const#
Array<unsigned char> getAtomOrders(Index index) const#
bool removeBonds(Index atom)#
void addBonds(const Array<std::pair<Index, Index>> &bonds, const Array<unsigned char> &orders)#
void swapBond(Index a, Index b)#
void swapAtom(Index a, Index b)#
std::list<Index> getAtomsAtLayer(size_t layer)#
Layer &layer()#
const Layer &layer() const#
void boundingBox(Vector3 &boxMin, Vector3 &boxMax, const double radius = 1.0) const#

Calculte and return bounding box of the whole molecule or selected atoms only.

Parameters:
  • boxMin – [out] the minimum corner (first end of the box diagonal)

  • boxMax – [out] the maximum corner (second end of the box diagonal)

  • radius – [in] radius of a single sphere

Public Static Functions

static std::pair<Vector3, Vector3> bestFitPlane(const Array<Vector3> &pos)#
Returns:

The normal vector of the best-fit plane of some specific atoms.

static inline std::pair<Index, Index> makeBondPair(const Index &a, const Index &b)#

Protected Attributes

VariantMap m_data#
std::map<std::string, MatrixX> m_partialCharges#

Sets of atomic partial charges.

std::map<std::string, MatrixX> m_spectra#

Sets of spectra.

CustomElementMap m_customElementMap#
ElementMask m_elements#

Which elements this molecule contains (e.g., for force fields)

Array<Vector2> m_positions2d#
Array<Vector3> m_positions3d#
Array<std::string> m_atomLabels#
Array<std::string> m_bondLabels#
Array<Array<Vector3>> m_coordinates3d#

Store conformers/trajectories.

Array<double> m_timesteps#
Array<AtomHybridization> m_hybridizations#
Array<signed char> m_formalCharges#
Array<Vector3> m_forceVectors#
Array<Vector3ub> m_colors#
Array<double> m_vibrationFrequencies#
Array<double> m_vibrationIRIntensities#
Array<double> m_vibrationRamanIntensities#
Array<Array<Vector3>> m_vibrationLx#
std::vector<bool> m_selectedAtoms#
std::vector<Mesh*> m_meshes#
std::vector<Cube*> m_cubes#
BasisSet *m_basisSet#
UnitCell *m_unitCell#
Array<Residue> m_residues#
unsigned short m_hallNumber = 0#
Eigen::VectorXd m_frozenAtomMask#