Class Scine::Utils::QuaternionFit

class QuaternionFit

Compares two sets of positions and calculates translation and rotation to go from one to the other.

References:

Public Functions

template<typename DerivedA, typename DerivedB>
QuaternionFit(const Eigen::DenseBase<DerivedA> &refMat, const Eigen::DenseBase<DerivedB> &fitMat, const Eigen::VectorXd &weights, bool improperRotationIsAllowed = false)

Construct a new Quaternion Fit object.

Parameters
  • refMat: The positions of the reference data.

  • fitMat: The positions of the data to be fitted.

  • weights: The individual weights for each point in the set of data.

  • improperRotationIsAllowed: If true allows the algorithm to invert the geometry to generate a better fit.

template<typename DerivedA, typename DerivedB>
QuaternionFit(const Eigen::DenseBase<DerivedA> &refMat, const Eigen::DenseBase<DerivedB> &fitMat, bool improperRotationIsAllowed = false)

Construct a new Quaternion Fit object.

Parameters
  • Rref: The positions of the reference data.

  • Rfit: The positions of the data to be fitted.

  • improperRotationIsAllowed: If true allows the algorithm to invert the geometry to generate a better fit.

template<typename DerivedA, typename DerivedB>
QuaternionFit(const Eigen::DenseBase<DerivedA> &refMat, const Eigen::DenseBase<DerivedB> &fitMat, const ElementTypeCollection &elements, bool improperRotationIsAllowed = false)

Construct a new Quaternion Fit object.

If this contructor is used, the element masses will be used as weights for the fit.

Parameters
  • Rref: The atom positions of the reference geometry.

  • Rfit: The atom positions of the geometry to be fitted.

  • elements: The elements for the given atoms coordinates, these will be used as weights.

  • improperRotationIsAllowed: If true allows the algorithm to invert the geometry to generate a better fit.

Eigen::Matrix3d getRotationMatrix() const

Getter for the reverse of the applied rotation.

Return

Eigen::Matrix3d The rotation from the reference structure to the second structure before fitting.

Eigen::Vector3d getTransVector() const

Getter for the reverse of the applied translation.

Return

Eigen::Vector3d The translation from the reference structure to the second structure before fitting.

Eigen::Matrix<double, Eigen::Dynamic, 3, Eigen::RowMajor> getFittedData() const

Getter for the fitted data as Eigen3 object.

Return

Eigen::Matrix3d The fitted data.

double getRotRMSD() const

Getter for the RMSD not using any weights that might be stored.

Return

double The RMSD.

double getRMSD() const

Getter for the RMSD due to differences in rotation only.

Return

double The RMSD.

double getWeightedRMSD(const Eigen::VectorXd &weights) const

Getter for the RMSD using the given weights.

Return

double The RMSD.

Parameters
  • weights: The individaul weights for each point in the set of data.

double getWeightedRMSD() const

Getter for the RMSD using the internal weights given/implied in the constructor.

Return

double The RMSD.