Namespace Scine::Utils::Geometry

namespace Geometry

Functionalities working with an entire geometry (PositionCollection).

Functions

PositionCollection translatePositions(const PositionCollection &positions, const Eigen::Ref<Eigen::RowVector3d> &translation)

Translates a set of postions by a given displacement.

Return

PositionCollection Returns the translated positions.

Parameters
  • positions: The original positions.

  • translation: The displacement to be added.

void translatePositionsInPlace(PositionCollection &positions, const Eigen::Ref<Eigen::RowVector3d> &translation)

Translates a set of positions by a given displacement.

The translation happens in-place.

Parameters
  • positions: The initial positions, will be transformed in-place.

  • translation: The displacement to be added.

int getIndexOfClosestAtom(const PositionCollection &positions, const Position &targetPosition, double squaredDistanceConsideredZero = -1.0)

Get the index of closest position (atom) to a given position in space.

Return

int The index of the closest atom to the given target.

Parameters
  • positions: A set of positions to be traversed.

  • targetPosition: The target position.

  • distanceConsideredZero: Squared distance between two positions resulting in them being considered equal and therefore skipping this position as it is already present in the given PositionCollection. The default is set to a negative value, so that all positions are considered as possible candidates for the closest one.

int getIndexOfAtomInStructure(const AtomCollection &structure, const Atom &atom, double squaredDistanceConsideredZero = 1e-4)

Get the index of a given atom in a given molecular structure.

Return

int The index of the given atom in the given molecular structure.

Parameters
  • structure: The structure to be traversed.

  • atom: The target atom that is tried to be found in the given structure.

  • distanceConsideredZero: Squared distance between two atoms resulting in them being considered equal.

Exceptions
  • std::runtime_error: Function throws error if the atom is not found in the given structure.

Eigen::MatrixXd positionVectorToMatrix(const Eigen::VectorXd &v)

Transforms a 3N-dimensional vector {x0, y0, z0, x1, y1, z1, …} to a Nx3 matrix.

Return

Eigen::MatrixXd Returns the final Nx3 matrix.

Parameters
  • v: The positions in vector form.

Eigen::VectorXd positionMatrixToVector(const Eigen::MatrixXd &m)

Transforms a Nx3 matrix to a 3N-dimensional vector {x0, y0, z0, x1, y1, z1, …}.

Return

Eigen::VectorXd Returns the final vector.

Parameters
  • m: The positions in matrix form.

void alignPositions(const PositionCollection &reference, PositionCollection &positions)

Rotate and translate positions so that it is as close as possible to referencePositions.

Parameters
  • reference: The reference positions.

  • positions: The positions to be aligned, will be transformed in place.

std::vector<double> getMasses(const ElementTypeCollection &elements)

Get a vector of all masses (in a.u.).

Return

std::vector<double> Returns the masses listed in a vector.

Parameters
  • elements: A collection of elements.

Position getCenterOfMass(const PositionCollection &positions, const std::vector<double> &masses)

Get the center of mass.

Return

Position Returns the center of mass (COM).

Parameters
  • positions: The positions.

  • masses: The masses (sorted according to the positions).

Position getCenterOfMass(const AtomCollection &structure)

Get the center of mass.

Return

Position Returns the center of mass (COM).

Parameters
  • structure: The structure (positions and masses are relevant).

Position getAveragePosition(const PositionCollection &positions)

Get the average Position.

(The average Position is identical to the center of mass if all masses are identical)

Return

Position Returns the average position.

Parameters
  • positions: The positions.

Eigen::Matrix3d calculateInertiaTensor(const PositionCollection &positions, const std::vector<double> &masses, const Position &centerOfMass)

Calculates the inertia tensor.

Return

Eigen::Matrix3d Returns the inertia tensor.

Parameters
  • positions: The positions.

  • masses: The masses (sorted according to the positions).

  • centerOfMass: The center of mass.

PrincipalMomentsOfInertia calculatePrincipalMoments(const PositionCollection &positions, const std::vector<double> &masses, const Position &centerOfMass)

Calculates the principal moments of inertia.

Return

PrincipalMomentsOfInertia Returns the principal moments of inertia.

Parameters
  • positions: The positions.

  • masses: The masses (sorted according to the positions).

  • centerOfMass: The center of mass.

Eigen::MatrixXd calculateTranslationAndRotationModes(const PositionCollection &positions, const ElementTypeCollection &elements)

Calculated the cartesian modes corresponding to translations and rotations of the entire system.

Return

Eigen::MatrixXd The rotation and translation modes. Translation modes (x,y,z) first-third column, rotation modes 3-final column (depending on the geometry)

Parameters
  • positions: The positions of all atoms.

  • elements: The ElemenetTypes of all atoms.

Eigen::MatrixXd calculateRotTransFreeTransformMatrix(const PositionCollection &positions, const ElementTypeCollection &elements, bool massWeighted = false)

Generates the matrix removing rotation and translation modes from the given geometry if applied.

Return

Eigen::MatrixXd The transformation matrix (applied as X^T*H*X to the Hessian).

Parameters
  • positions: The positions of all atoms.

  • elements: The ElemenetTypes of all atoms.

  • massWeighted: True if Hessian to be transformed is mass-weighted

class PrincipalMomentsOfInertia
#include <Geometry.h>

The principal moments of inertia stored.