absOrientation.hpp
Go to the documentation of this file.
1 //-----------------------------------------------------------------------------------
2 // d-SEAMS - Deferred Structural Elucidation Analysis for Molecular Simulations
3 //
4 // Copyright (c) 2018--present d-SEAMS core team
5 //
6 // This program is free software: you can redistribute it and/or modify
7 // it under the terms of the MIT License as published by
8 // the Open Source Initiative.
9 //
10 // A copy of the MIT License is included in the LICENSE file of this repository.
11 // You should have received a copy of the MIT License along with this program.
12 // If not, see <https://opensource.org/licenses/MIT>.
13 //-----------------------------------------------------------------------------------
14 
15 #ifndef __ABSORIENTATION_H_
16 #define __ABSORIENTATION_H_
17 
18 #include <algorithm>
19 #include <array>
20 #include <fstream>
21 #include <iostream>
22 #include <iterator>
23 #include <math.h>
24 #include <memory>
25 #include <sstream>
26 #include <string>
27 #include <sys/stat.h>
28 #include <vector>
29 
30 // External
31 #include <Spectra/SymEigsShiftSolver.h>
32 #include <Spectra/SymEigsSolver.h>
33 #include <Eigen/Core>
34 #include <Eigen/Dense>
35 
36 #include <mol_sys.hpp>
37 #include <ring.hpp>
38 #include <seams_input.hpp>
39 #include <seams_output.hpp>
40 
41 // Inspired by AStar_Dual_Tree_HandPose by jsupancic
42 
43 namespace absor {
44 
48 int hornAbsOrientation(const Eigen::MatrixXd &refPoints,
49  const Eigen::MatrixXd &targetPoints,
50  std::vector<double> *quat, double *rmsd,
51  std::vector<double> *rmsdList, double *scale);
52 
55 Eigen::MatrixXd calcMatrixS(const Eigen::MatrixXd &centeredRefPnts,
56  const Eigen::MatrixXd &centeredTargetPnts, int nop,
57  int dim);
58 
61 Eigen::MatrixXd calcMatrixN(const Eigen::MatrixXd &S);
62 
64 Eigen::MatrixXd centerWRTcentroid(const Eigen::MatrixXd &pointSet);
65 
67 double calcScaleFactor(const Eigen::MatrixXd &rightSys,
68  const Eigen::MatrixXd &leftSys, int n);
69 
71 Eigen::MatrixXd quat2RotMatrix(const Eigen::VectorXd &quat);
72 
74 double getRMSD(const Eigen::MatrixXd &centeredRefPnts,
75  const Eigen::MatrixXd &centeredTargetPnts,
76  const Eigen::VectorXd &quat, std::vector<double> *rmsdList,
77  int nop, double scale);
78 
79 } // namespace absor
80 
81 #endif // __ABSORIENTATION_H_
The main molecular system handler.
double getRMSD(const Eigen::MatrixXd &centeredRefPnts, const Eigen::MatrixXd &centeredTargetPnts, const Eigen::VectorXd &quat, std::vector< double > *rmsdList, int nop, double scale)
Calculate the RMSD.
Eigen::MatrixXd calcMatrixS(const Eigen::MatrixXd &centeredRefPnts, const Eigen::MatrixXd &centeredTargetPnts, int nop, int dim)
double calcScaleFactor(const Eigen::MatrixXd &rightSys, const Eigen::MatrixXd &leftSys, int n)
Calculate the scale factor from the centered left and right point sets.
Eigen::MatrixXd centerWRTcentroid(const Eigen::MatrixXd &pointSet)
Center a point set wrt the centroid.
int hornAbsOrientation(const Eigen::MatrixXd &refPoints, const Eigen::MatrixXd &targetPoints, std::vector< double > *quat, double *rmsd, std::vector< double > *rmsdList, double *scale)
Get the absolute orientation using Horn's algorithm (with quaternions)
Eigen::MatrixXd calcMatrixN(const Eigen::MatrixXd &S)
Eigen::MatrixXd quat2RotMatrix(const Eigen::VectorXd &quat)
Get a rotation matrix from a unit quaternion.
File containing common functions used by bulk and confined topological network critera.
File for functions that read in files).