occ
Loading...
Searching...
No Matches
occ::core::Molecule Class Reference

Storage class for relevant information of a Molecule. More...

#include <molecule.h>

Public Types

enum  Origin { Cartesian , Centroid , CenterOfMass }
 An enum to specify the origin used in some calculations on the Molecule object. More...
 

Public Member Functions

 Molecule ()
 
 Molecule (const IVec &nums, const Mat3N &pos)
 Construct a Molecule from a vector of atomic numbers and positions.
 
 Molecule (const std::vector< core::Atom > &atoms)
 Construct a Molecule from a vector of Atom objects.
 
 Molecule (const std::vector< Element > &els, const std::vector< std::array< double, 3 > > &pos)
 Construct a Molecule from a vector of Element objects, and std::array<double, 3> positions.
 
template<typename N , typename D >
 Molecule (const std::vector< N > &nums, const std::vector< std::array< D, 3 > > &pos)
 Construct a Molecule from a vector of atomic numbers a vector arrays of positions.
 
size_t size () const
 The number of atoms in this molecule.
 
void set_name (const std::string &name)
 Set the name for this molecule.
 
const std::string & name () const
 Get the name for this molecule.
 
Vec interatomic_distances () const
 A vector representing the compressed distance matrix for this Molecule.
 
const auto & elements () const
 The Element objects for the atoms in this Molecule.
 
const Mat3N & positions () const
 The positions of the atoms in this Molecule.
 
const IVec & atomic_numbers () const
 The atomic numbers for the atoms in this Molecule.
 
Vec vdw_radii () const
 The van der Waals radii for the atoms in this Molecule.
 
Vec covalent_radii () const
 The van der Waals radii for the atoms in this Molecule.
 
Vec atomic_masses () const
 The atomic masses of the atoms in this Molecule.
 
double molar_mass () const
 The total molecular mass.
 
std::vector< core::Atomatoms () const
 Convert the atoms represented by this molecule to a std:vector<Atom>
 
void set_cell_shift (const IVec3 &shift, bool update_atoms=true)
 Set the unit cell offset for this Molecule (default 000)
 
const IVec3cell_shift () const
 Get the unit cell offset for this Molecule (default 000)
 
Vec3 centroid () const
 The geometric centre of this molecule.
 
Vec3 center_of_mass () const
 The mass-weighted geometric centre of this molecule.
 
Mat3 inertia_tensor () const
 The tensor of inertia (in matrix form).
 
Vec3 principal_moments_of_inertia () const
 The principal moments of inertia, derived from the intertia tensor.
 
Vec3 rotational_constants () const
 The rotational constants of this Molecule in GHz.
 
double rotational_free_energy (double T) const
 The rotational component of the entropic free energy of this Molecule (rigid molecule, ideal gas approximation), units are kJ/mol.
 
double translational_free_energy (double T) const
 The translational component of the entropic free energy of this Molecule (rigid molecule, ideal gas approximation), units are kJ/mol.
 
std::tuple< size_t, size_t, double > nearest_atom (const Molecule &rhs) const
 The nearest atom-atom pair between atoms in this Molecule and atoms in another.
 
void add_bond (size_t l, size_t r)
 Manually add a bond connection between two atoms in this Molecule.
 
void set_bonds (const std::vector< std::pair< size_t, size_t > > &bonds)
 Set the known bond connections in this Molecule.
 
const std::vector< std::pair< size_t, size_t > > & bonds () const
 Get the list of bonds in this Molecule.
 
int charge () const
 The net charge of this Molecule.
 
void set_charge (int c)
 Set the net charge of this Molecule.
 
int multiplicity () const
 Get the spin multiplicity of this molecule.
 
void set_multiplicity (int m)
 Set the spin multiplicity of this molecule.
 
int num_electrons () const
 The number of electrons in this molecule.
 
int asymmetric_molecule_idx () const
 Get the index into Crystal::symmetry_unique_molecules for this Molecule.
 
void set_asymmetric_molecule_idx (size_t idx)
 Set the index into Crystal::symmetry_unique_molecules for this Molecule.
 
int unit_cell_molecule_idx () const
 Get the index into Crystal::unit_cell_molecules for this Molecule.
 
void set_unit_cell_molecule_idx (size_t idx)
 Set the index into Crystal::unit_cell_molecules for this Molecule.
 
void set_asymmetric_unit_transformation (const Mat3 &rot, const Vec3 &trans)
 Set the transformation from the corresponding molecule in the asymmetric unit of a Crystal to this geometry.
 
std::pair< Mat3, Vec3asymmetric_unit_transformation () const
 Get the transformation from the corresponding molecule in the asymmetric unit of a Crystal to this geometry.
 
void set_unit_cell_idx (const IVec &idx)
 Set the unit cell atom indices for all atoms in this Molecule.
 
void set_unit_cell_shift (const IMat3N &shifts)
 Set the unit cell atom shifts for all atoms in this Molecule.
 
void set_asymmetric_unit_idx (const IVec &idx)
 Set the asymmetric unit atom indices for all atoms in this Molecule.
 
void set_asymmetric_unit_symop (const IVec &symop)
 Set the associated SymmetryOperation for all all atoms in from their asymmetric unit counterpart (encoded as an int).
 
const auto & unit_cell_idx () const
 Get the unit cell atom indices for all atoms in this Molecule.
 
const auto & unit_cell_shift () const
 Get the unit cell atom shifts for all atoms in this Molecule.
 
const auto & asymmetric_unit_idx () const
 Get the asymmetric unit atom indices for all atoms in this Molecule.
 
const auto & asymmetric_unit_symop () const
 Get the set of SymmetryOperation for all atoms in this Molecule (i.e.
 
void rotate (const Eigen::Affine3d &r, Origin o=Cartesian)
 Rotate this Molecule by the given rotation, about the specified origin.
 
void rotate (const Mat3 &r, Origin o=Cartesian)
 Rotate this Molecule by the given rotation, about the specified origin.
 
void rotate (const Mat3 &r, const Vec3 &o)
 Rotate this Molecule by the given rotation, about the specified origin.
 
void transform (const Mat4 &t, Origin o=Cartesian)
 Transform this Molecule by the given homogeneous transformation matrix, with rotation performed about the specified origin.
 
void transform (const Mat4 &t, const Vec3 &o)
 Transform this Molecule by the given homogeneous transformation matrix, with rotation performed about the specified origin.
 
void translate (const Vec3 &t)
 Translate this Molecule by the given vector.
 
Molecule rotated (const Eigen::Affine3d &r, Origin o=Cartesian) const
 A copy of this Molecule transformed by the given rotation, about the specified origin.
 
Molecule rotated (const Mat3 &r, Origin o=Cartesian) const
 A copy of this Molecule transformed by the given rotation, about the specified origin.
 
Molecule rotated (const Mat3 &r, const Vec3 &o) const
 A copy of this Molecule transformed by the given rotation, about the specified origin.
 
Molecule transformed (const Mat4 &t, Origin o=Cartesian) const
 A copy of this Molecule transformed by the given homogeneous transformation matrix, with rotation performed about the specified origin.
 
Molecule transformed (const Mat4 &t, const Vec3 &o) const
 A copy of this Molecule transformed by the given homogeneous transformation matrix, with rotation performed about the specified origin.
 
Molecule translated (const Vec3 &t) const
 A copy of this Molecule translated by the given vector.
 
bool is_comparable_to (const Molecule &rhs) const
 Determine whether this Molecule and another can be sensibly compared.
 
bool is_equivalent_to (const Molecule &rhs) const
 Determine whether this Molecule and another are equivalent.
 
void set_partial_charges (const Vec &v)
 
const auto & partial_charges ()
 
Vec esp_partial_charges (const occ::Mat3N &positions) const
 

Detailed Description

Storage class for relevant information of a Molecule.

The role of the Molecule class is to store basic information about atoms which are bonded, and to facilitate convenient calculations of properties of these atoms.

Member Enumeration Documentation

◆ Origin

An enum to specify the origin used in some calculations on the Molecule object.

Enumerator
Cartesian 

The Cartesian origin i.e.

(0, 0, 0) in R3

Centroid 

The molecular centroid i.e.

average position of atoms (ignoring mass)

CenterOfMass 

The centre of mass i.e.

the weighted average positon

Constructor & Destructor Documentation

◆ Molecule() [1/5]

occ::core::Molecule::Molecule ( )
inlineexplicit

◆ Molecule() [2/5]

occ::core::Molecule::Molecule ( const IVec &  nums,
const Mat3N &  pos 
)

Construct a Molecule from a vector of atomic numbers and positions.

Parameters
numsVector of length N atomic numbers
posMatrix (3, N) of atomic positions (Angstroms)

◆ Molecule() [3/5]

occ::core::Molecule::Molecule ( const std::vector< core::Atom > &  atoms)

Construct a Molecule from a vector of Atom objects.

Parameters
atomsstd::vector<Atom> of length N atoms.

Will convert from the expected Bohr unit of the Atom objects internally to Angstroms.

◆ Molecule() [4/5]

occ::core::Molecule::Molecule ( const std::vector< Element > &  els,
const std::vector< std::array< double, 3 > > &  pos 
)

Construct a Molecule from a vector of Element objects, and std::array<double, 3> positions.

Parameters
elsstd::vector<Element> of length N Element objects.
posstd::vector<std::array<double, 3> > of length N of atomic positions.

Positions are assumed to be in Angstroms.

◆ Molecule() [5/5]

template<typename N , typename D >
occ::core::Molecule::Molecule ( const std::vector< N > &  nums,
const std::vector< std::array< D, 3 > > &  pos 
)
inline

Construct a Molecule from a vector of atomic numbers a vector arrays of positions.

Parameters
numsstd::vector<N> of atomic numbers, convertible to int
posstd::vector<std::array<D, 3> > of atomic positions, convertible to double.

Positions are assumed to be in Angstroms.

Member Function Documentation

◆ add_bond()

void occ::core::Molecule::add_bond ( size_t  l,
size_t  r 
)
inline

Manually add a bond connection between two atoms in this Molecule.

Parameters
lthe first index of an atom into this Molecule.
rthe second index of an atom into this Molecule.

Internally, the set of bonds in this Molecule will be updated.

◆ asymmetric_molecule_idx()

int occ::core::Molecule::asymmetric_molecule_idx ( ) const
inline

Get the index into Crystal::symmetry_unique_molecules for this Molecule.

Returns
an integer representing the index (default is -1)
Warning
This will return an invalid index (-1) by default, and should be checked or ensured to have a value before use.

◆ asymmetric_unit_idx()

const auto & occ::core::Molecule::asymmetric_unit_idx ( ) const
inline

Get the asymmetric unit atom indices for all atoms in this Molecule.

Returns
a const ref to the internal IVec containing the relevant indices into the set of asymmetric unit atoms for the Crystal associated with this Molecule.

See also Molecule::asymmetric_molecule_idx if you wish to get the mapping for the molecule index rather than atom indices.

◆ asymmetric_unit_symop()

const auto & occ::core::Molecule::asymmetric_unit_symop ( ) const
inline

Get the set of SymmetryOperation for all atoms in this Molecule (i.e.

mapping from their corresponding asymmetric unit atom to their position here).

Returns
a const ref to the internal IVec containing the relevant integer encoded SymmeteryOperations for the mapping from asymmetric unit atoms for the Crystal associated with this Molecule.

See also Molecule::asymmetric_unit_transformation if you wish to get the mapping for the molecule index rather than atom indices.

◆ asymmetric_unit_transformation()

std::pair< Mat3, Vec3 > occ::core::Molecule::asymmetric_unit_transformation ( ) const
inline

Get the transformation from the corresponding molecule in the asymmetric unit of a Crystal to this geometry.

Returns
a std::pair<Mat3, Vec3> representing the rotation part of the transform. and the translation part of the transform.

By default the rotation is the identity matrix, and the translation is the zero vector.

See also Molecule::asymmetric_molecule_idx

◆ atomic_masses()

Vec occ::core::Molecule::atomic_masses ( ) const

The atomic masses of the atoms in this Molecule.

Returns
a newly constructed Vec containing the masses.

Creates a vector, they are not stored internally to the Molecule object.

◆ atomic_numbers()

const IVec & occ::core::Molecule::atomic_numbers ( ) const
inline

The atomic numbers for the atoms in this Molecule.

Returns
a const ref to the internal vector of atomic numbers.

◆ atoms()

std::vector< core::Atom > occ::core::Molecule::atoms ( ) const

Convert the atoms represented by this molecule to a std:vector<Atom>

Returns
a vector of Atom objects (positions in Bohr).

◆ bonds()

const std::vector< std::pair< size_t, size_t > > & occ::core::Molecule::bonds ( ) const
inline

Get the list of bonds in this Molecule.

Returns
a std::vector<std::pair<size_t, size_t>> i.e. list of the bond edges that have been set in this Molecule.

◆ cell_shift()

const IVec3 & occ::core::Molecule::cell_shift ( ) const

Get the unit cell offset for this Molecule (default 000)

Returns
the CellShift stored in this Molecule

See Molecule::set_cell_shift

◆ center_of_mass()

Vec3 occ::core::Molecule::center_of_mass ( ) const

The mass-weighted geometric centre of this molecule.

Returns
a Vec3 representing the average atomic position (including masses)

See Molecule::centroid if you don't need to incorporate masses

◆ centroid()

Vec3 occ::core::Molecule::centroid ( ) const

The geometric centre of this molecule.

Returns
a Vec3 representing the average atomic position (ignoring masses)

See Molecule::center_of_mass if you need to incorporate masses

◆ charge()

int occ::core::Molecule::charge ( ) const
inline

The net charge of this Molecule.

Returns
an integer representing the net charge. Default is neutral (0)

◆ covalent_radii()

Vec occ::core::Molecule::covalent_radii ( ) const

The van der Waals radii for the atoms in this Molecule.

Returns
a newly constructed Vec containing van der Waals radii

Creates a vector, they are not stored internally to the Molecule object.

◆ elements()

const auto & occ::core::Molecule::elements ( ) const
inline

The Element objects for the atoms in this Molecule.

Returns
a const std::vector<Element>& to the internal vector of elements

◆ esp_partial_charges()

Vec occ::core::Molecule::esp_partial_charges ( const occ::Mat3N positions) const

◆ inertia_tensor()

Mat3 occ::core::Molecule::inertia_tensor ( ) const

The tensor of inertia (in matrix form).

Returns
a Mat3 representing the distribution of mass about the Molecule::center_of_mass of this Molecule.

◆ interatomic_distances()

Vec occ::core::Molecule::interatomic_distances ( ) const

A vector representing the compressed distance matrix for this Molecule.

Returns
distances a Vec object representing the distances

Since the NxN distance matrix for all atoms in this molecule is by definition symmetric, and the diagonals would be 0, this will return a N * (N - 1) / 2 length vector where the distance for atom i to atom j can be calculated in the following loop:

for (size_t i = 0; i < N; i++) {
for (size_t j = i + 1; j < N; j++) {
// distance
}
}

◆ is_comparable_to()

bool occ::core::Molecule::is_comparable_to ( const Molecule rhs) const

Determine whether this Molecule and another can be sensibly compared.

Returns
true if they are comparable Molecules per the definition, false otherwise.

The working definition here of whether two molecules are comparable is whether they have the exact same IVec of atomic numbers, in the same order (i.e. the same chemical composition, and the order of atoms is not permuted.

◆ is_equivalent_to()

bool occ::core::Molecule::is_equivalent_to ( const Molecule rhs) const

Determine whether this Molecule and another are equivalent.

Parameters
rhsThe other Molecule object to check equivalence
Returns
true if they are comparable Molecules per the definition, false otherwise.

The working definition here of whether two molecules are equivalent if a) they are comparable (i.e. all atomic numbers are the same) and b) all interatomic distances are the same.

◆ molar_mass()

double occ::core::Molecule::molar_mass ( ) const

The total molecular mass.

Returns
a double representing the sum of atomic masses in this molecule

◆ multiplicity()

int occ::core::Molecule::multiplicity ( ) const
inline

Get the spin multiplicity of this molecule.

Returns
an integer representing the spin multiplicity of this molecule (default 1)

◆ name()

const std::string & occ::core::Molecule::name ( ) const
inline

Get the name for this molecule.

Returns
name std::string representing the name/identifier for this Molecule.

See Molecule::set_name

◆ nearest_atom()

std::tuple< size_t, size_t, double > occ::core::Molecule::nearest_atom ( const Molecule rhs) const

The nearest atom-atom pair between atoms in this Molecule and atoms in another.

Parameters
rhsThe other Molecule object to find distances.
Returns
a std::tuple<size_t, size_t, double> consisting of the index into atoms in this Molecule, the index into the atoms of the Molecule rhs, and the distance in Angstroms.

See Molecule::nearest_atom_distance if you only need the distance

◆ num_electrons()

int occ::core::Molecule::num_electrons ( ) const
inline

The number of electrons in this molecule.

Returns
an integer representing total number of electrons.

Molecule::charge is incorporated, so the result should be the sum of atomic numbers - net charge (as electrons have -ve charge).

◆ partial_charges()

const auto & occ::core::Molecule::partial_charges ( )
inline

◆ positions()

const Mat3N & occ::core::Molecule::positions ( ) const
inline

The positions of the atoms in this Molecule.

Returns
a const ref to the internal matrix of positions (Angstroms).

◆ principal_moments_of_inertia()

Vec3 occ::core::Molecule::principal_moments_of_inertia ( ) const

The principal moments of inertia, derived from the intertia tensor.

Returns
a Vec3 representing the three moments of inertia. Eigenvalues calculated using Eigen::SelfAdjointEigenSolver.

◆ rotate() [1/3]

void occ::core::Molecule::rotate ( const Eigen::Affine3d &  r,
Origin  o = Cartesian 
)

Rotate this Molecule by the given rotation, about the specified origin.

Parameters
rThe rotation.
oThe desired Origin about which to rotate.

All rotations will be in the Cartesian axis frame (though not necessarily about the Cartesian origin)

See Molecule::rotated if you wish to return a copy of this Molecule.

◆ rotate() [2/3]

void occ::core::Molecule::rotate ( const Mat3 r,
const Vec3 o 
)

Rotate this Molecule by the given rotation, about the specified origin.

Overload for any Mat3.

Parameters
rThe rotation.
oThe desired position about which to rotate (Angstroms).

All rotations will be in the Cartesian axis frame (though not necessarily about the Cartesian origin).

Warning
The provided Mat3 is not checked to be a valid rotation matrix.

See Molecule::rotated if you wish to return a copy of this Molecule.

◆ rotate() [3/3]

void occ::core::Molecule::rotate ( const Mat3 r,
Origin  o = Cartesian 
)

Rotate this Molecule by the given rotation, about the specified origin.

Overload for any Mat3.

Parameters
rThe rotation.
oThe desired Origin about which to rotate.

All rotations will be in the Cartesian axis frame (though not necessarily about the Cartesian origin).

Warning
The provided Mat3 is not checked to be a valid rotation matrix.

See Molecule::rotated if you wish to return a copy of this Molecule.

◆ rotated() [1/3]

Molecule occ::core::Molecule::rotated ( const Eigen::Affine3d &  r,
Origin  o = Cartesian 
) const

A copy of this Molecule transformed by the given rotation, about the specified origin.

Parameters
rThe rotation.
oThe desired Origin about which to rotate.
Returns
a rotated copy of this Molecule.

All rotations will be in the Cartesian axis frame (though not necessarily about the Cartesian origin)

See Molecule::rotate if you wish to modify this Molecule.

◆ rotated() [2/3]

Molecule occ::core::Molecule::rotated ( const Mat3 r,
const Vec3 o 
) const

A copy of this Molecule transformed by the given rotation, about the specified origin.

Parameters
rThe rotation.
oThe desired position about which to rotate (Angstroms).
Returns
a rotated copy of this Molecule.

All rotations will be in the Cartesian axis frame (though not necessarily about the Cartesian origin)

Warning
The provided Mat3 is not checked to be a valid rotation matrix.

See Molecule::rotate if you wish to modify this Molecule.

◆ rotated() [3/3]

Molecule occ::core::Molecule::rotated ( const Mat3 r,
Origin  o = Cartesian 
) const

A copy of this Molecule transformed by the given rotation, about the specified origin.

Parameters
rThe rotation.
oThe desired Origin about which to rotate.
Returns
a rotated copy of this Molecule.

All rotations will be in the Cartesian axis frame (though not necessarily about the Cartesian origin)

Warning
The provided Mat3 is not checked to be a valid rotation matrix.

See Molecule::rotate if you wish to modify this Molecule.

◆ rotational_constants()

Vec3 occ::core::Molecule::rotational_constants ( ) const

The rotational constants of this Molecule in GHz.

Returns
a Vec3 representing the three rotational constants. Calculated directly from Molecule::principal_moments_of_inertia.

◆ rotational_free_energy()

double occ::core::Molecule::rotational_free_energy ( double  T) const

The rotational component of the entropic free energy of this Molecule (rigid molecule, ideal gas approximation), units are kJ/mol.

Parameters
Tthe temperature (K).
Returns
a double representing the rotational free energy for this Molecule at the given temperature in kJ/mol

◆ set_asymmetric_molecule_idx()

void occ::core::Molecule::set_asymmetric_molecule_idx ( size_t  idx)
inline

Set the index into Crystal::symmetry_unique_molecules for this Molecule.

Parameters
idxan integer representing the index.
Warning
No check is performed to ensure this is a sensible value.

◆ set_asymmetric_unit_idx()

void occ::core::Molecule::set_asymmetric_unit_idx ( const IVec &  idx)
inline

Set the asymmetric unit atom indices for all atoms in this Molecule.

Parameters
idxan IVec containing the relevant indices into the set of asymmetric unit atoms for the Crystal associated with this Molecule.

See also Molecule::set_asymmetric_molecule_idx if you wish to set the molecule index rather than atom indices.

◆ set_asymmetric_unit_symop()

void occ::core::Molecule::set_asymmetric_unit_symop ( const IVec &  symop)

Set the associated SymmetryOperation for all all atoms in from their asymmetric unit counterpart (encoded as an int).

Parameters
symopan IVec containing the integer encoded SymmetryOperation from the corresponding asymmetric unit atoms in the Crystal for all atoms in.

See also Molecule::set_asymmetric_unit_transformation if you wish to set the transformation from the molecule. And note that these are not the same thing.

◆ set_asymmetric_unit_transformation()

void occ::core::Molecule::set_asymmetric_unit_transformation ( const Mat3 rot,
const Vec3 trans 
)
inline

Set the transformation from the corresponding molecule in the asymmetric unit of a Crystal to this geometry.

Parameters
rotA Mat3 representing the rotation part of the transform.
transA Vec3 representing the translation part of the transform.

By default the rotation is the identity matrix, and the translation is the zero vector. Molecule::set_asymmetric_molecule_index should also be set for this to make sense. It is the responsibility of the caller to ensure these values correspond.

◆ set_bonds()

void occ::core::Molecule::set_bonds ( const std::vector< std::pair< size_t, size_t > > &  bonds)
inline

Set the known bond connections in this Molecule.

Parameters
bondsa std::vector<std::pair<size_t, size_t>> i.e. list of the bond edges to set in this Molecule.
Warning
Existing bonds will be removed, this method will overwrite bonds already stored.

◆ set_cell_shift()

void occ::core::Molecule::set_cell_shift ( const IVec3 shift,
bool  update_atoms = true 
)

Set the unit cell offset for this Molecule (default 000)

Parameters
shiftthe desired CellShift

See Molecule::cell_shift

◆ set_charge()

void occ::core::Molecule::set_charge ( int  c)
inline

Set the net charge of this Molecule.

Parameters
can integer representing the net charge.

No checks are performed to ensure this charge is sensible.

◆ set_multiplicity()

void occ::core::Molecule::set_multiplicity ( int  m)
inline

Set the spin multiplicity of this molecule.

Parameters
man integer representing the spin multiplicity of this molecule.

No checks are performed to ensure this is a sensible value.

◆ set_name()

void occ::core::Molecule::set_name ( const std::string &  name)
inline

Set the name for this molecule.

Parameters
namestd::string representing the name/identifier for this Molecule.

See Molecule::name

◆ set_partial_charges()

void occ::core::Molecule::set_partial_charges ( const Vec &  v)
inline

◆ set_unit_cell_idx()

void occ::core::Molecule::set_unit_cell_idx ( const IVec &  idx)
inline

Set the unit cell atom indices for all atoms in this Molecule.

Parameters
idxan IVec containing the relevant indices into the set of unit cell atoms for the Crystal associated with this Molecule.

See also Molecule::set_unit_cell_molecule_idx if you wish to set the molecule index rather than atom indices.

◆ set_unit_cell_molecule_idx()

void occ::core::Molecule::set_unit_cell_molecule_idx ( size_t  idx)
inline

Set the index into Crystal::unit_cell_molecules for this Molecule.

Parameters
idxan integer representing the index.
Warning
No check is performed to ensure this is a sensible value.

◆ set_unit_cell_shift()

void occ::core::Molecule::set_unit_cell_shift ( const IMat3N shifts)
inline

Set the unit cell atom shifts for all atoms in this Molecule.

Parameters
shiftsan IMat3N containing the relevant shifts to generate the atoms of this molecule from the unit cell atoms of the Crystal associated with this Molecule.

See also Molecule::set_cell_shift if you wish to set the unit cell molecule shift rather than atom indices.

◆ size()

size_t occ::core::Molecule::size ( ) const
inline

The number of atoms in this molecule.

Returns
size size_t representing the number of atoms in this Molecule.

Calculated based on the internal IVec of atomic numbers.

◆ transform() [1/2]

void occ::core::Molecule::transform ( const Mat4 t,
const Vec3 o 
)

Transform this Molecule by the given homogeneous transformation matrix, with rotation performed about the specified origin.

Parameters
tThe homogeneous transformation matrix.
oThe desired position about which to rotate (Angstroms).

All rotations will be in the Cartesian axis frame (though not necessarily about the Cartesian origin)

See Molecule::transformed if you wish to return a copy of this Molecule.

◆ transform() [2/2]

void occ::core::Molecule::transform ( const Mat4 t,
Origin  o = Cartesian 
)

Transform this Molecule by the given homogeneous transformation matrix, with rotation performed about the specified origin.

Parameters
tThe homogeneous transformation matrix.
oThe desired Origin about which to rotate.

All rotations will be in the Cartesian axis frame (though not necessarily about the Cartesian origin)

See Molecule::transformed if you wish to return a copy of this Molecule.

◆ transformed() [1/2]

Molecule occ::core::Molecule::transformed ( const Mat4 t,
const Vec3 o 
) const

A copy of this Molecule transformed by the given homogeneous transformation matrix, with rotation performed about the specified origin.

Parameters
tThe homogeneous transformation matrix.
oThe desired position about which to rotate (Angstroms).
Returns
transformed copy of this Molecule.

All rotations will be in the Cartesian axis frame (though not necessarily about the Cartesian origin)

See Molecule::transforme if you wish to modify this Molecule.

◆ transformed() [2/2]

Molecule occ::core::Molecule::transformed ( const Mat4 t,
Origin  o = Cartesian 
) const

A copy of this Molecule transformed by the given homogeneous transformation matrix, with rotation performed about the specified origin.

Parameters
tThe homogeneous transformation matrix.
oThe desired origin about which to rotate.
Returns
transformed copy of this Molecule.

All rotations will be in the Cartesian axis frame (though not necessarily about the Cartesian origin)

See Molecule::transform if you wish to modify this Molecule.

◆ translate()

void occ::core::Molecule::translate ( const Vec3 t)

Translate this Molecule by the given vector.

Parameters
ttranslation vector (Angstroms).

See Molecule::translated if you wish to return a copy of this Molecule.

◆ translated()

Molecule occ::core::Molecule::translated ( const Vec3 t) const

A copy of this Molecule translated by the given vector.

Parameters
ttranslation vector (Angstroms).
Returns
a translated copy of this Molecule.

See Molecule::translate if you wish to modify this Molecule.

◆ translational_free_energy()

double occ::core::Molecule::translational_free_energy ( double  T) const

The translational component of the entropic free energy of this Molecule (rigid molecule, ideal gas approximation), units are kJ/mol.

Parameters
Tthe temperature (K).
Returns
a double representing the translational component of the free energy for this Molecule at the given temperature in kJ/mol

◆ unit_cell_idx()

const auto & occ::core::Molecule::unit_cell_idx ( ) const
inline

Get the unit cell atom indices for all atoms in this Molecule.

Returns
a const ref to the internal IVec containing the relevant indices into the set of unit cell atoms for the Crystal associated with this Molecule.

See also Molecule::unit_cell_molecule_idx if you wish to get the mapping for the molecule index rather than atom indices.

◆ unit_cell_molecule_idx()

int occ::core::Molecule::unit_cell_molecule_idx ( ) const
inline

Get the index into Crystal::unit_cell_molecules for this Molecule.

Returns
an integer representing the index (default is -1)
Warning
This will return an invalid index (-1) by default, and should be checked or ensured to have a value before use.

◆ unit_cell_shift()

const auto & occ::core::Molecule::unit_cell_shift ( ) const
inline

Get the unit cell atom shifts for all atoms in this Molecule.

Returns
an IMat3N containing the relevant shifts to generate the atoms of this molecule from the unit cell atoms of the Crystal associated with this Molecule.

See also Molecule::cell_shift if you wish to get the unit cell molecule shift rather than atom indices.

◆ vdw_radii()

Vec occ::core::Molecule::vdw_radii ( ) const

The van der Waals radii for the atoms in this Molecule.

Returns
a newly constructed Vec containing van der Waals radii

Creates a vector, they are not stored internally to the Molecule object.


The documentation for this class was generated from the following file: