19#include "math_eigen_impl.h"
54 for (i = 0; i < size(); i++) {
68 if (data.size() > 0) {
95 data.assign(data.size(), T(0.0));
98 inline size_t size()
const
103 inline void resize(
size_t const n)
113 inline T & operator [] (
size_t const i) {
117 inline T
const & operator [] (
size_t const i)
const {
121 inline static void check_sizes(vector1d<T>
const &v1, vector1d<T>
const &v2)
123 if (v1.size() != v2.size()) {
124 cvm::error(
"Error: trying to perform an operation between vectors of different sizes, "+
129 inline void operator += (vector1d<T>
const &v)
131 check_sizes(*
this, v);
133 for (i = 0; i < this->size(); i++) {
138 inline void operator -= (vector1d<T>
const &v)
140 check_sizes(*
this, v);
142 for (i = 0; i < this->size(); i++) {
150 for (i = 0; i < this->size(); i++) {
158 for (i = 0; i < this->size(); i++) {
163 inline friend vector1d<T> operator + (vector1d<T>
const &v1,
164 vector1d<T>
const &v2)
166 check_sizes(v1.size(), v2.size());
167 vector1d<T> result(v1.size());
169 for (i = 0; i < v1.size(); i++) {
170 result[i] = v1[i] + v2[i];
175 inline friend vector1d<T> operator - (vector1d<T>
const &v1,
176 vector1d<T>
const &v2)
178 check_sizes(v1.size(), v2.size());
179 vector1d<T> result(v1.size());
181 for (i = 0; i < v1.size(); i++) {
182 result[i] = v1[i] - v2[i];
187 inline friend vector1d<T> operator * (vector1d<T>
const &v,
cvm::real a)
189 vector1d<T> result(v.size());
191 for (i = 0; i < v.size(); i++) {
192 result[i] = v[i] * a;
197 inline friend vector1d<T> operator * (
cvm::real a, vector1d<T>
const &v)
202 inline friend vector1d<T> operator / (vector1d<T>
const &v,
cvm::real a)
204 vector1d<T> result(v.size());
206 for (i = 0; i < v.size(); i++) {
207 result[i] = v[i] / a;
215 check_sizes(v1.size(), v2.size());
218 for (i = 0; i < v1.size(); i++) {
219 prod += v1[i] * v2[i];
229 for (i = 0; i < this->size(); i++) {
230 result += (*this)[i] * (*this)[i];
244 for (i = 0; i < this->size(); i++) {
245 result += (*this)[i];
253 if ((i2 < i1) || (i2 >= this->size())) {
254 cvm::error(
"Error: trying to slice a vector using incorrect boundaries.\n");
258 for (i = 0; i < (i2 - i1); i++) {
259 result[i] = (*this)[i1+i];
268 if ((i2 < i1) || (i2 >= this->size())) {
269 cvm::error(
"Error: trying to slice a vector using incorrect boundaries.\n");
272 for (i = 0; i < (i2 - i1); i++) {
273 (*this)[i1+i] = v[i];
281 return real_width*(this->size()) + 3*(this->size()-1) + 4;
284 inline friend std::istream & operator >> (std::istream &is,
287 if (v.size() == 0)
return is;
288 std::streampos
const start_pos = is.tellg();
290 if ( !(is >> sep) || !(sep ==
'(') ) {
292 is.seekg(start_pos, std::ios::beg);
293 is.setstate(std::ios::failbit);
297 while ( (is >> v[count]) &&
298 (count < (v.size()-1) ? ((is >> sep) && (sep ==
',')) :
true) ) {
299 if (++count == v.size())
break;
301 if (count < v.size()) {
303 is.seekg(start_pos, std::ios::beg);
304 is.setstate(std::ios::failbit);
309 inline friend std::ostream & operator << (std::ostream &os,
312 std::streamsize
const w = os.width();
313 std::streamsize
const p = os.precision();
318 for (i = 0; i < v.size()-1; i++) {
319 os.width(w); os.precision(p);
322 os.width(w); os.precision(p);
323 os << v[v.size()-1] <<
" )";
327 inline std::string to_simple_string()
const
329 if (this->size() == 0)
return std::string(
"");
330 std::ostringstream os;
331 os.setf(std::ios::scientific, std::ios::floatfield);
335 for (i = 1; i < this->size(); i++) {
336 os <<
" " << (*this)[i];
341 inline int from_simple_string(std::string
const &s)
343 std::stringstream stream(s);
346 while ((stream >> (*
this)[i]) && (i < this->size())) {
349 if (i < this->size()) {
350 return COLVARS_ERROR;
354 while (stream >> input) {
355 if ((i % 100) == 0) {
356 data.reserve(data.size()+100);
358 data.resize(data.size()+1);
386 inline row(T *
const row_data,
size_t const inner_length)
387 : data(row_data), length(inner_length)
389 inline T & operator [] (
size_t const j) {
392 inline T
const & operator [] (
size_t const j)
const {
401 if (v.size() != length) {
402 return cvm::error(
"Error: setting a matrix row from a vector of "
403 "incompatible size.\n", COLVARS_BUG_ERROR);
405 for (
size_t i = 0; i < length; i++) data[i] = v[i];
411 std::vector<row> rows;
412 std::vector<T *> pointers;
417 inline void resize(
size_t const ol,
size_t const il)
419 if ((ol > 0) && (il > 0)) {
421 if (data.size() > 0) {
424 std::vector<T> new_data(ol * il);
425 for (i = 0; i < outer_length; i++) {
426 for (j = 0; j < inner_length; j++) {
427 new_data[il*i+j] = data[inner_length*i+j];
430 data.resize(ol * il);
434 data.resize(ol * il);
440 if (data.size() > 0) {
444 rows.reserve(outer_length);
446 pointers.reserve(outer_length);
447 for (i = 0; i < outer_length; i++) {
448 rows.push_back(
row(&(data[0])+inner_length*i, inner_length));
449 pointers.push_back(&(data[0])+inner_length*i);
468 data.assign(data.size(), T(0.0));
471 inline size_t size()
const
478 : outer_length(0), inner_length(0)
483 inline matrix2d(
size_t const ol,
size_t const il)
484 : outer_length(ol), inner_length(il)
486 this->
resize(outer_length, inner_length);
492 : outer_length(m.outer_length), inner_length(m.inner_length)
495 this->
resize(outer_length, inner_length);
517 inline row & operator [] (
size_t const i)
521 inline row
const & operator [] (
size_t const i)
const
529 if ((outer_length != m.outer_length) || (inner_length != m.inner_length)){
531 outer_length = m.outer_length;
532 inner_length = m.inner_length;
533 this->
resize(outer_length, inner_length);
541 if (rows.size() > 0) {
542 return &(pointers[0]);
550 if ((m1.outer_length != m2.outer_length) ||
551 (m1.inner_length != m2.inner_length)) {
552 cvm::error(
"Error: trying to perform an operation between "
553 "matrices of different sizes, "+
561 inline void operator += (matrix2d<T>
const &m)
563 check_sizes(*
this, m);
565 for (i = 0; i < data.size(); i++) {
566 data[i] += m.data[i];
570 inline void operator -= (matrix2d<T>
const &m)
572 check_sizes(*
this, m);
574 for (i = 0; i < data.size(); i++) {
575 data[i] -= m.data[i];
582 for (i = 0; i < data.size(); i++) {
590 for (i = 0; i < data.size(); i++) {
595 inline friend matrix2d<T> operator + (matrix2d<T>
const &m1,
596 matrix2d<T>
const &m2)
599 matrix2d<T> result(m1.outer_length, m1.inner_length);
601 for (i = 0; i < m1.data.size(); i++) {
602 result.data[i] = m1.data[i] + m2.data[i];
607 inline friend matrix2d<T> operator - (matrix2d<T>
const &m1,
608 matrix2d<T>
const &m2)
611 matrix2d<T> result(m1.outer_length, m1.inner_length);
613 for (i = 0; i < m1.data.size(); i++) {
614 result.data[i] = m1.data[i] - m1.data[i];
619 inline friend matrix2d<T> operator * (matrix2d<T>
const &m,
cvm::real a)
621 matrix2d<T> result(m.outer_length, m.inner_length);
623 for (i = 0; i < m.data.size(); i++) {
624 result.data[i] = m.data[i] * a;
629 inline friend matrix2d<T> operator * (
cvm::real a, matrix2d<T>
const &m)
634 inline friend matrix2d<T> operator / (matrix2d<T>
const &m,
cvm::real a)
636 matrix2d<T> result(m.outer_length, m.inner_length);
638 for (i = 0; i < m.data.size(); i++) {
639 result.data[i] = m.data[i] * a;
649 if (m.outer_length != v.size()) {
650 cvm::error(
"Error: trying to multiply a vector and a matrix "
651 "of incompatible sizes, "+
657 for (i = 0; i < m.inner_length; i++) {
658 for (k = 0; k < m.outer_length; k++) {
659 result[i] += m[k][i] * v[k];
670 std::streamsize
const w = os.width();
671 std::streamsize
const p = os.precision();
676 for (i = 0; i < m.outer_length; i++) {
679 for (j = 0; j < m.inner_length-1; j++) {
682 os << m[i][j] <<
" , ";
686 os << m[i][m.inner_length-1] <<
" )";
693 inline std::string to_simple_string()
const
695 if (this->size() == 0)
return std::string(
"");
696 std::ostringstream os;
697 os.setf(std::ios::scientific, std::ios::floatfield);
701 for (i = 1; i < data.size(); i++) {
702 os <<
" " << data[i];
707 inline int from_simple_string(std::string
const &s)
709 std::stringstream stream(s);
711 while ((i < data.size()) && (stream >> data[i])) {
714 if (i < data.size()) {
715 return COLVARS_ERROR;
748 set(v[0], v[1], v[2]);
772 return (i == 0) ? x : (i == 1) ? y : (i == 2) ? z : x;
777 return (i == 0) ? x : (i == 1) ? y : (i == 2) ? z : x;
810 inline void operator /= (
cvm::real const& v)
819 return (x*x + y*y + z*z);
833 static inline size_t output_width(
size_t real_width)
835 return 3*real_width + 10;
843 -v1.x*v2.z + v2.x*v1.z,
844 v1.x*v2.y - v2.x*v1.y);
855 return cvm::rvector(v1.x + v2.x, v1.y + v2.y, v1.z + v2.z);
860 return cvm::rvector(v1.x - v2.x, v1.y - v2.y, v1.z - v2.z);
867 return v1.x * v2.x + v1.y * v2.y + v1.z * v2.z;
885 std::string to_simple_string()
const;
886 int from_simple_string(std::string
const &s);
896 cvm::real xx, xy, xz, yx, yy, yz, zx, zy, zz;
923 xx = xy = xz = yx = yy = yz = zx = zy = zz = 0.0;
930 ( xx * (yy*zz - zy*yz))
931 - (yx * (xy*zz - zy*xz))
932 + (zx * (xy*yz - yy*xz));
946 m.yx*r.x + m.yy*r.y + m.yz*r.z,
947 m.zx*r.x + m.zy*r.y + m.zz*r.z);
963 : q0(qv[0]), q1(qv[1]), q2(qv[2]), q3(qv[3])
971 : q0(q0i), q1(q1i), q2(q2i), q3(q3i)
975 : q0(v[0]), q1(v[1]), q2(v[2]), q3(v[3])
987 q0 = q1 = q2 = q3 = 0.0;
993 return 4*real_width + 13;
996 std::string to_simple_string()
const;
997 int from_simple_string(std::string
const &s);
1016 cvm::error(
"Error: incorrect quaternion component.\n");
1033 cvm::error(
"Error: trying to access a quaternion "
1034 "component which is not between 0 and 3.\n");
1052 return q0*q0 + q1*q1 + q2*q2 + q3*q3;
1069 q0 *= a; q1 *= a; q2 *= a; q3 *= a;
1074 q0 /= a; q1 /= a; q2 /= a; q3 /= a;
1077 inline void set_positive()
1079 if (q0 > 0.0)
return;
1088 q0+=h.q0; q1+=h.q1; q2+=h.q2; q3+=h.q3;
1092 q0-=h.q0; q1-=h.q1; q2-=h.q2; q3-=h.q3;
1120 h.q0*q.q1 + h.q1*q.q0 + h.q2*q.q3 - h.q3*q.q2,
1121 h.q0*q.q2 + h.q2*q.q0 + h.q3*q.q1 - h.q1*q.q3,
1122 h.q0*q.q3 + h.q3*q.q0 + h.q1*q.q2 - h.q2*q.q1);
1163 R.xx = q0*q0 + q1*q1 - q2*q2 - q3*q3;
1164 R.yy = q0*q0 - q1*q1 + q2*q2 - q3*q3;
1165 R.zz = q0*q0 - q1*q1 - q2*q2 + q3*q3;
1167 R.xy = 2.0 * (q1*q2 - q0*q3);
1168 R.xz = 2.0 * (q0*q2 + q1*q3);
1170 R.yx = 2.0 * (q0*q3 + q1*q2);
1171 R.yz = 2.0 * (q2*q3 - q0*q1);
1173 R.zx = 2.0 * (q1*q3 - q0*q2);
1174 R.zy = 2.0 * (q0*q1 + q2*q3);
1211 return std::array<cvm::real, 4>{{
1212 2.0 * ( q0 * C[0][0] - q3 * C[0][1] + q2 * C[0][2] +
1213 q3 * C[1][0] + q0 * C[1][1] - q1 * C[1][2] +
1214 -q2 * C[2][0] + q1 * C[2][1] + q0 * C[2][2]),
1215 2.0 * ( q1 * C[0][0] + q2 * C[0][1] + q3 * C[0][2] +
1216 q2 * C[1][0] - q1 * C[1][1] - q0 * C[1][2] +
1217 q3 * C[2][0] + q0 * C[2][1] - q1 * C[2][2]),
1218 2.0 * (-q2 * C[0][0] + q1 * C[0][1] + q0 * C[0][2] +
1219 q1 * C[1][0] + q2 * C[1][1] + q3 * C[1][2] +
1220 -q0 * C[2][0] + q3 * C[2][1] - q2 * C[2][2]),
1221 2.0 * (-q3 * C[0][0] - q0 * C[0][1] + q1 * C[0][2] +
1222 q0 * C[1][0] - q3 * C[1][1] + q2 * C[1][2] +
1223 q1 * C[2][0] + q2 * C[2][1] + q3 * C[2][2])
1232 return 2.0*iprod*iprod - 1.0;
1240 cvm::real const cos_omega = this->q0*Q2.q0 + this->q1*Q2.q1 +
1241 this->q2*Q2.q2 + this->q3*Q2.q3;
1244 ( (cos_omega < -1.0) ? -1.0 : cos_omega) );
1247 if (cos_omega > 0.0)
1248 return omega * omega;
1250 return (PI-omega) * (PI-omega);
1257 cvm::real const cos_omega = this->q0*Q2.q0 + this->q1*Q2.q1 + this->q2*Q2.q2 + this->q3*Q2.q3;
1259 ( (cos_omega < -1.0) ? -1.0 : cos_omega) );
1268 grad1((-1.0)*sin_omega*Q2.q0 + cos_omega*(this->q0-cos_omega*Q2.q0)/sin_omega,
1269 (-1.0)*sin_omega*Q2.q1 + cos_omega*(this->q1-cos_omega*Q2.q1)/sin_omega,
1270 (-1.0)*sin_omega*Q2.q2 + cos_omega*(this->q2-cos_omega*Q2.q2)/sin_omega,
1271 (-1.0)*sin_omega*Q2.q3 + cos_omega*(this->q3-cos_omega*Q2.q3)/sin_omega);
1273 if (cos_omega > 0.0) {
1274 return 2.0*omega*grad1;
1276 return -2.0*(PI-omega)*grad1;
1284 cvm::real const cos_omega = this->q0*Q2.q0 + this->q1*Q2.q1 +
1285 this->q2*Q2.q2 + this->q3*Q2.q3;
1286 if (cos_omega < 0.0) Q2 *= -1.0;
1293 cvm::real const prod = this->q0*Q2.q0 + this->q1*Q2.q1 +
1294 this->q2*Q2.q2 + this->q3*Q2.q3;
1301#ifndef COLVARS_LAMMPS
1303int diagonalize_matrix(
cvm::real m[4][4],
1347 const std::vector<cvm::real> &pos1,
1348 const std::vector<cvm::real> &pos2,
1349 const size_t num_atoms_pos1,
1350 const size_t num_atoms_pos2);
1363 std::vector<atom_pos>
const &pos2);
1364 void calc_optimal_rotation_soa(
1365 std::vector<cvm::real>
const &pos1,
1366 std::vector<cvm::real>
const &pos2,
1367 const size_t num_atoms_pos1,
1368 const size_t num_atoms_pos2);
1409 while (alpha > 180.0) alpha -= 360;
1410 while (alpha < -180.0) alpha += 360;
1424 (180.0/PI) * 2.0 * (1.0 / (1.0 + (iprod*iprod)/(
q.q0*
q.q0)));
1427 dspindx * ((1.0/
q.q0) * axis.x),
1428 dspindx * ((1.0/
q.q0) * axis.y),
1429 dspindx * ((1.0/
q.q0) * axis.z));
1434 return cvm::quaternion((180.0/PI) * 2.0 * ((-1.0)/iprod), 0.0, 0.0, 0.0);
1444 (180.0/PI) * 2.0 *
cvm::atan2(axis * q_vec,
q.q0);
1447 cvm::real const cos_theta_2 = ( (cos_spin_2 != 0.0) ?
1448 (
q.q0 / cos_spin_2) :
1451 return 2.0 * (cos_theta_2*cos_theta_2) - 1.0;
1465 (4.0 *
q.q0 / (cos_spin_2*cos_spin_2)) *
1466 (1.0 - (iprod*iprod)/(
q.q0*
q.q0) / (1.0 + (iprod*iprod)/(
q.q0*
q.q0)));
1469 (4.0 *
q.q0 / (cos_spin_2*cos_spin_2) *
1470 (iprod/
q.q0) / (1.0 + (iprod*iprod)/(
q.q0*
q.q0)));
1473 d_cos_theta_dqn * axis.x,
1474 d_cos_theta_dqn * axis.y,
1475 d_cos_theta_dqn * axis.z);
1479 (4.0 / (cos_spin_2*cos_spin_2 * iprod));
1482 d_cos_theta_dqn * axis.x,
1483 d_cos_theta_dqn * axis.y,
1484 d_cos_theta_dqn * axis.z);
1502 std::vector<cvm::atom_pos>
const &pos2);
Definition: colvartypes.h:382
Arbitrary size array (two dimensions) suitable for linear algebra operations (i.e....
Definition: colvartypes.h:373
matrix2d(matrix2d< T > const &m)
Copy constructor.
Definition: colvartypes.h:491
T ** c_array()
Return the 2-d C array.
Definition: colvartypes.h:540
void clear()
Deallocation routine.
Definition: colvartypes.h:460
~matrix2d()
Destructor.
Definition: colvartypes.h:501
std::vector< T > const & data_array() const
Return a reference to the data.
Definition: colvartypes.h:512
matrix2d< T > & operator=(matrix2d< T > const &m)
Assignment.
Definition: colvartypes.h:527
matrix2d()
Default constructor.
Definition: colvartypes.h:477
std::vector< T > & data_array()
Return a reference to the data.
Definition: colvartypes.h:506
void resize(size_t const ol, size_t const il)
Allocation routine, used by all constructors.
Definition: colvartypes.h:417
friend std::ostream & operator<<(std::ostream &os, matrix2d< T > const &m)
Formatted output.
Definition: colvartypes.h:667
void reset()
Set all elements to zero.
Definition: colvartypes.h:466
1-dimensional vector of real numbers with four components and a quaternion algebra
Definition: colvartypes.h:955
cvm::real & operator[](int i)
Access the quaternion as a 4-d array (return a reference)
Definition: colvartypes.h:1005
cvm::real cosine(cvm::quaternion const &q) const
Return the cosine between the orientation frame associated to this quaternion and another.
Definition: colvartypes.h:1229
cvm::rvector get_vector() const
Return the vector component.
Definition: colvartypes.h:1096
friend std::istream & operator>>(std::istream &is, cvm::quaternion &q)
Formatted input operator.
Definition: colvartypes.cpp:123
cvm::quaternion rotate(cvm::quaternion const &Q2) const
Rotate Q2 through this quaternion (put it in the rotated reference frame)
Definition: colvartypes.h:1152
cvm::real norm() const
Norm of the quaternion.
Definition: colvartypes.h:1056
quaternion(cvm::real const qv[4])
Constructor component by component.
Definition: colvartypes.h:962
std::array< cvm::real, 4 > derivative_element_wise_product_sum(const cvm::real(&C)[3][3]) const
Calculate the sums of element-wise products of a given matrix with respect to dR/dq0,...
Definition: colvartypes.h:1210
quaternion()
Default constructor.
Definition: colvartypes.h:979
void reset()
Set all components to zero (null quaternion)
Definition: colvartypes.h:985
friend std::ostream & operator<<(std::ostream &os, cvm::quaternion const &q)
Formatted output operator.
Definition: colvartypes.cpp:104
cvm::rvector rotate(cvm::rvector const &v) const
Rotate v through this quaternion (put it in the rotated reference frame)
Definition: colvartypes.h:1144
cvm::quaternion dist2_grad(cvm::quaternion const &Q2) const
Definition: colvartypes.h:1255
void match(cvm::quaternion &Q2) const
Choose the closest between Q2 and -Q2 and save it back. Not required for dist2() and dist2_grad()
Definition: colvartypes.h:1282
static size_t output_width(size_t real_width)
Tell the number of characters required to print a quaternion, given that of a real number.
Definition: colvartypes.h:991
cvm::real norm2() const
Square norm of the quaternion.
Definition: colvartypes.h:1050
cvm::rmatrix rotation_matrix() const
Return the 3x3 matrix associated to this quaternion.
Definition: colvartypes.h:1159
cvm::quaternion conjugate() const
Return the conjugate quaternion.
Definition: colvartypes.h:1062
cvm::real inner(cvm::quaternion const &Q2) const
Inner product (as a 4-d vector) with Q2; requires match() if the largest overlap is looked for.
Definition: colvartypes.h:1291
quaternion(cvm::real q0i, cvm::real q1i, cvm::real q2i, cvm::real q3i)
Constructor component by component.
Definition: colvartypes.h:967
friend cvm::quaternion operator*(cvm::quaternion const &h, cvm::quaternion const &q)
Provides the quaternion product. NOTE: for the inner product use: h.inner (q);
Definition: colvartypes.h:1116
cvm::real dist2(cvm::quaternion const &Q2) const
Square distance from another quaternion on the 4-dimensional unit sphere: returns the square of the a...
Definition: colvartypes.h:1238
2-dimensional array of real numbers with three components along each dimension (works with colvarmodu...
Definition: colvartypes.h:892
rmatrix(cvm::real xxi, cvm::real xyi, cvm::real xzi, cvm::real yxi, cvm::real yyi, cvm::real yzi, cvm::real zxi, cvm::real zyi, cvm::real zzi)
Constructor component by component.
Definition: colvartypes.h:905
cvm::real determinant() const
Return the determinant.
Definition: colvartypes.h:927
rmatrix()
Default constructor.
Definition: colvartypes.h:899
A rotation between two sets of coordinates (for the moment a wrapper for colvarmodule::quaternion)
Definition: colvartypes.h:1313
~rotation()
Destructor.
Definition: colvartypes.cpp:198
cvm::rmatrix C
Correlation matrix C (3, 3)
Definition: colvartypes.h:1316
rotation()
Default constructor.
Definition: colvartypes.cpp:160
cvm::real cos_theta(cvm::rvector const &axis) const
Return the projection of the orientation vector onto a predefined axis.
Definition: colvartypes.h:1440
cvm::real S_eigval[4]
Eigenvalues of S.
Definition: colvartypes.h:1322
cvm::real S_backup[4][4]
Used for debugging gradients.
Definition: colvartypes.h:1328
cvm::rvector rotate(cvm::rvector const &v) const
Return the rotated vector.
Definition: colvartypes.h:1386
int init()
Initialize member data.
Definition: colvartypes.cpp:151
cvm::real S[4][4]
Overlap matrix S (4, 4)
Definition: colvartypes.h:1319
cvm::quaternion dspin_angle_dq(cvm::rvector const &axis) const
Return the derivative of the spin angle with respect to the quaternion.
Definition: colvartypes.h:1416
cvm::quaternion dcos_theta_dq(cvm::rvector const &axis) const
Return the derivative of the tilt wrt the quaternion.
Definition: colvartypes.h:1455
cvm::real spin_angle(cvm::rvector const &axis) const
Return the spin angle (in degrees) with respect to the provided axis (which MUST be normalized)
Definition: colvartypes.h:1405
cvm::real S_eigvec[4][4]
Eigenvectors of S.
Definition: colvartypes.h:1325
void calc_optimal_rotation(std::vector< atom_pos > const &pos1, std::vector< atom_pos > const &pos2)
Calculate the optimal rotation and store the corresponding eigenvalue and eigenvector in the argument...
Definition: colvartypes.cpp:379
static cvm::real crossing_threshold
Threshold for the eigenvalue crossing test.
Definition: colvartypes.h:1491
cvm::quaternion q_old
Previous value of the rotation (used to warn the user when the structure changes too much,...
Definition: colvartypes.h:1498
static bool monitor_crossings
Whether to test for eigenvalue crossing.
Definition: colvartypes.h:1489
bool b_debug_gradients
Perform gradient tests.
Definition: colvartypes.h:1332
void build_correlation_matrix(std::vector< cvm::atom_pos > const &pos1, std::vector< cvm::atom_pos > const &pos2)
Build the correlation matrix C (used by calc_optimal_rotation())
Definition: colvartypes.cpp:209
cvm::rmatrix matrix() const
Return the associated 3x3 matrix.
Definition: colvartypes.h:1398
void * jacobi
Pointer to instance of Jacobi solver.
Definition: colvartypes.h:1511
void compute_overlap_matrix()
Compute the overlap matrix S (used by calc_optimal_rotation())
Definition: colvartypes.cpp:229
void debug_gradients(cvm::rotation &rot, const std::vector< cvm::real > &pos1, const std::vector< cvm::real > &pos2, const size_t num_atoms_pos1, const size_t num_atoms_pos2)
Function for debugging gradients.
Definition: colvartypes.cpp:291
void calc_optimal_rotation_impl()
Actual implementation of calc_optimal_rotation (and called by it)
Definition: colvartypes.cpp:429
cvm::rotation inverse() const
Return the inverse of this rotation.
Definition: colvartypes.h:1392
cvm::quaternion q
The rotation itself (implemented as a quaternion)
Definition: colvartypes.h:1335
vector of real numbers with three components
Definition: colvartypes.h:724
cvm::real & operator[](int i)
Access cartesian components by index.
Definition: colvartypes.h:771
void set(cvm::real x_i, cvm::real y_i, cvm::real z_i)
Assign all components.
Definition: colvartypes.h:763
void reset()
Set all components to zero.
Definition: colvartypes.h:736
void set(cvm::real value)
Set all components to a scalar.
Definition: colvartypes.h:757
friend cvm::real operator*(cvm::rvector const &v1, cvm::rvector const &v2)
Inner (dot) product.
Definition: colvartypes.h:864
Arbitrary size array (one dimensions) suitable for linear algebra operations (i.e....
Definition: colvartypes.h:34
void reset()
Set all elements to zero.
Definition: colvartypes.h:93
size_t output_width(size_t real_width) const
Formatted output.
Definition: colvartypes.h:279
cvm::real norm2() const
Squared norm.
Definition: colvartypes.h:225
void sliceassign(size_t const i1, size_t const i2, vector1d< T > const &v)
Assign a vector to a slice of this vector.
Definition: colvartypes.h:265
std::vector< T > & data_array()
Return a reference to the data.
Definition: colvartypes.h:76
vector1d & operator=(const vector1d &)=default
Explicit Copy assignement.
vector1d(const vector1d &)=default
Explicit Copy constructor.
vector1d(size_t const n, T const *t)
Constructor from C array.
Definition: colvartypes.h:49
vector1d(size_t const n=0)
Default constructor.
Definition: colvartypes.h:42
T * c_array()
Return a pointer to the data location.
Definition: colvartypes.h:66
vector1d< T > const slice(size_t const i1, size_t const i2) const
Slicing.
Definition: colvartypes.h:251
std::vector< T > const & data_array() const
Return a reference to the data.
Definition: colvartypes.h:82
double real
Defining an abstract real number allows to switch precision.
Definition: colvarmodule.h:139
static real cos(real const &x)
Reimplemented to work around MS compiler issues.
Definition: colvarmodule.h:189
static int error(std::string const &message, int code=-1)
Print a message to the main log and set global error code.
Definition: colvarmodule.cpp:2096
static real atan2(real const &x, real const &y)
Reimplemented to work around MS compiler issues.
Definition: colvarmodule.h:234
static size_t const cv_prec
Number of digits to represent a collective variables value(s)
Definition: colvarmodule.h:744
static real sqrt(real const &x)
Reimplemented to work around MS compiler issues.
Definition: colvarmodule.h:177
static real acos(real const &x)
Reimplemented to work around compiler issues; return hard-coded values for boundary conditions.
Definition: colvarmodule.h:218
static real fabs(real const &x)
Reimplemented to work around MS compiler issues.
Definition: colvarmodule.h:171
static real sin(real const &x)
Reimplemented to work around MS compiler issues.
Definition: colvarmodule.h:183
static std::string to_str(char const *s)
Convert to string for output purposes.
Definition: colvarmodule.cpp:2446
Collective variables main module.
Helper class for calculating the derivative of rotation.
Definition: colvar_rotation_derivative.h:41