18#include "math_eigen_impl.h"
24#define PI 3.14159265358979323846
57 for (i = 0; i < size(); i++) {
71 if (data.size() > 0) {
98 data.assign(data.size(), T(0.0));
101 inline size_t size()
const
106 inline void resize(
size_t const n)
116 inline T & operator [] (
size_t const i) {
120 inline T
const & operator [] (
size_t const i)
const {
124 inline static void check_sizes(vector1d<T>
const &v1, vector1d<T>
const &v2)
126 if (v1.size() != v2.size()) {
127 cvm::error(
"Error: trying to perform an operation between vectors of different sizes, "+
132 inline void operator += (vector1d<T>
const &v)
134 check_sizes(*
this, v);
136 for (i = 0; i < this->size(); i++) {
141 inline void operator -= (vector1d<T>
const &v)
143 check_sizes(*
this, v);
145 for (i = 0; i < this->size(); i++) {
153 for (i = 0; i < this->size(); i++) {
161 for (i = 0; i < this->size(); i++) {
166 inline friend vector1d<T> operator + (vector1d<T>
const &v1,
167 vector1d<T>
const &v2)
169 check_sizes(v1.size(), v2.size());
170 vector1d<T> result(v1.size());
172 for (i = 0; i < v1.size(); i++) {
173 result[i] = v1[i] + v2[i];
178 inline friend vector1d<T> operator - (vector1d<T>
const &v1,
179 vector1d<T>
const &v2)
181 check_sizes(v1.size(), v2.size());
182 vector1d<T> result(v1.size());
184 for (i = 0; i < v1.size(); i++) {
185 result[i] = v1[i] - v2[i];
190 inline friend vector1d<T> operator * (vector1d<T>
const &v,
cvm::real a)
192 vector1d<T> result(v.size());
194 for (i = 0; i < v.size(); i++) {
195 result[i] = v[i] * a;
200 inline friend vector1d<T> operator * (
cvm::real a, vector1d<T>
const &v)
205 inline friend vector1d<T> operator / (vector1d<T>
const &v,
cvm::real a)
207 vector1d<T> result(v.size());
209 for (i = 0; i < v.size(); i++) {
210 result[i] = v[i] / a;
218 check_sizes(v1.size(), v2.size());
221 for (i = 0; i < v1.size(); i++) {
222 prod += v1[i] * v2[i];
232 for (i = 0; i < this->size(); i++) {
233 result += (*this)[i] * (*this)[i];
247 for (i = 0; i < this->size(); i++) {
248 result += (*this)[i];
256 if ((i2 < i1) || (i2 >= this->size())) {
257 cvm::error(
"Error: trying to slice a vector using incorrect boundaries.\n");
261 for (i = 0; i < (i2 - i1); i++) {
262 result[i] = (*this)[i1+i];
271 if ((i2 < i1) || (i2 >= this->size())) {
272 cvm::error(
"Error: trying to slice a vector using incorrect boundaries.\n");
275 for (i = 0; i < (i2 - i1); i++) {
276 (*this)[i1+i] = v[i];
284 return real_width*(this->size()) + 3*(this->size()-1) + 4;
287 inline friend std::istream & operator >> (std::istream &is,
290 if (v.size() == 0)
return is;
291 std::streampos
const start_pos = is.tellg();
293 if ( !(is >> sep) || !(sep ==
'(') ) {
295 is.seekg(start_pos, std::ios::beg);
296 is.setstate(std::ios::failbit);
300 while ( (is >> v[count]) &&
301 (count < (v.size()-1) ? ((is >> sep) && (sep ==
',')) :
true) ) {
302 if (++count == v.size())
break;
304 if (count < v.size()) {
306 is.seekg(start_pos, std::ios::beg);
307 is.setstate(std::ios::failbit);
312 inline friend std::ostream & operator << (std::ostream &os,
315 std::streamsize
const w = os.width();
316 std::streamsize
const p = os.precision();
321 for (i = 0; i < v.size()-1; i++) {
322 os.width(w); os.precision(p);
325 os.width(w); os.precision(p);
326 os << v[v.size()-1] <<
" )";
330 inline std::string to_simple_string()
const
332 if (this->size() == 0)
return std::string(
"");
333 std::ostringstream os;
334 os.setf(std::ios::scientific, std::ios::floatfield);
338 for (i = 1; i < this->size(); i++) {
339 os <<
" " << (*this)[i];
344 inline int from_simple_string(std::string
const &s)
346 std::stringstream stream(s);
349 while ((stream >> (*
this)[i]) && (i < this->size())) {
352 if (i < this->size()) {
353 return COLVARS_ERROR;
357 while (stream >> input) {
358 if ((i % 100) == 0) {
359 data.reserve(data.size()+100);
361 data.resize(data.size()+1);
389 inline row(T *
const row_data,
size_t const inner_length)
390 : data(row_data), length(inner_length)
392 inline T & operator [] (
size_t const j) {
395 inline T
const & operator [] (
size_t const j)
const {
404 if (v.size() != length) {
405 return cvm::error(
"Error: setting a matrix row from a vector of "
406 "incompatible size.\n", COLVARS_BUG_ERROR);
408 for (
size_t i = 0; i < length; i++) data[i] = v[i];
414 std::vector<row> rows;
415 std::vector<T *> pointers;
420 inline void resize(
size_t const ol,
size_t const il)
422 if ((ol > 0) && (il > 0)) {
424 if (data.size() > 0) {
427 std::vector<T> new_data(ol * il);
428 for (i = 0; i < outer_length; i++) {
429 for (j = 0; j < inner_length; j++) {
430 new_data[il*i+j] = data[inner_length*i+j];
433 data.resize(ol * il);
437 data.resize(ol * il);
443 if (data.size() > 0) {
447 rows.reserve(outer_length);
449 pointers.reserve(outer_length);
450 for (i = 0; i < outer_length; i++) {
451 rows.push_back(
row(&(data[0])+inner_length*i, inner_length));
452 pointers.push_back(&(data[0])+inner_length*i);
471 data.assign(data.size(), T(0.0));
474 inline size_t size()
const
481 : outer_length(0), inner_length(0)
486 inline matrix2d(
size_t const ol,
size_t const il)
487 : outer_length(ol), inner_length(il)
489 this->
resize(outer_length, inner_length);
495 : outer_length(m.outer_length), inner_length(m.inner_length)
498 this->
resize(outer_length, inner_length);
520 inline row & operator [] (
size_t const i)
524 inline row
const & operator [] (
size_t const i)
const
532 if ((outer_length != m.outer_length) || (inner_length != m.inner_length)){
534 outer_length = m.outer_length;
535 inner_length = m.inner_length;
536 this->
resize(outer_length, inner_length);
544 if (rows.size() > 0) {
545 return &(pointers[0]);
553 if ((m1.outer_length != m2.outer_length) ||
554 (m1.inner_length != m2.inner_length)) {
555 cvm::error(
"Error: trying to perform an operation between "
556 "matrices of different sizes, "+
564 inline void operator += (matrix2d<T>
const &m)
566 check_sizes(*
this, m);
568 for (i = 0; i < data.size(); i++) {
569 data[i] += m.data[i];
573 inline void operator -= (matrix2d<T>
const &m)
575 check_sizes(*
this, m);
577 for (i = 0; i < data.size(); i++) {
578 data[i] -= m.data[i];
585 for (i = 0; i < data.size(); i++) {
593 for (i = 0; i < data.size(); i++) {
598 inline friend matrix2d<T> operator + (matrix2d<T>
const &m1,
599 matrix2d<T>
const &m2)
602 matrix2d<T> result(m1.outer_length, m1.inner_length);
604 for (i = 0; i < m1.data.size(); i++) {
605 result.data[i] = m1.data[i] + m2.data[i];
610 inline friend matrix2d<T> operator - (matrix2d<T>
const &m1,
611 matrix2d<T>
const &m2)
614 matrix2d<T> result(m1.outer_length, m1.inner_length);
616 for (i = 0; i < m1.data.size(); i++) {
617 result.data[i] = m1.data[i] - m1.data[i];
622 inline friend matrix2d<T> operator * (matrix2d<T>
const &m,
cvm::real a)
624 matrix2d<T> result(m.outer_length, m.inner_length);
626 for (i = 0; i < m.data.size(); i++) {
627 result.data[i] = m.data[i] * a;
632 inline friend matrix2d<T> operator * (
cvm::real a, matrix2d<T>
const &m)
637 inline friend matrix2d<T> operator / (matrix2d<T>
const &m,
cvm::real a)
639 matrix2d<T> result(m.outer_length, m.inner_length);
641 for (i = 0; i < m.data.size(); i++) {
642 result.data[i] = m.data[i] * a;
652 if (m.outer_length != v.size()) {
653 cvm::error(
"Error: trying to multiply a vector and a matrix "
654 "of incompatible sizes, "+
660 for (i = 0; i < m.inner_length; i++) {
661 for (k = 0; k < m.outer_length; k++) {
662 result[i] += m[k][i] * v[k];
673 std::streamsize
const w = os.width();
674 std::streamsize
const p = os.precision();
679 for (i = 0; i < m.outer_length; i++) {
682 for (j = 0; j < m.inner_length-1; j++) {
685 os << m[i][j] <<
" , ";
689 os << m[i][m.inner_length-1] <<
" )";
696 inline std::string to_simple_string()
const
698 if (this->size() == 0)
return std::string(
"");
699 std::ostringstream os;
700 os.setf(std::ios::scientific, std::ios::floatfield);
704 for (i = 1; i < data.size(); i++) {
705 os <<
" " << data[i];
710 inline int from_simple_string(std::string
const &s)
712 std::stringstream stream(s);
714 while ((i < data.size()) && (stream >> data[i])) {
717 if (i < data.size()) {
718 return COLVARS_ERROR;
751 set(v[0], v[1], v[2]);
775 return (i == 0) ? x : (i == 1) ? y : (i == 2) ? z : x;
780 return (i == 0) ? x : (i == 1) ? y : (i == 2) ? z : x;
813 inline void operator /= (
cvm::real const& v)
822 return (x*x + y*y + z*z);
836 static inline size_t output_width(
size_t real_width)
838 return 3*real_width + 10;
846 -v1.x*v2.z + v2.x*v1.z,
847 v1.x*v2.y - v2.x*v1.y);
858 return cvm::rvector(v1.x + v2.x, v1.y + v2.y, v1.z + v2.z);
863 return cvm::rvector(v1.x - v2.x, v1.y - v2.y, v1.z - v2.z);
870 return v1.x * v2.x + v1.y * v2.y + v1.z * v2.z;
888 std::string to_simple_string()
const;
889 int from_simple_string(std::string
const &s);
899 cvm::real xx, xy, xz, yx, yy, yz, zx, zy, zz;
926 xx = xy = xz = yx = yy = yz = zx = zy = zz = 0.0;
933 ( xx * (yy*zz - zy*yz))
934 - (yx * (xy*zz - zy*xz))
935 + (zx * (xy*yz - yy*xz));
949 m.yx*r.x + m.yy*r.y + m.yz*r.z,
950 m.zx*r.x + m.zy*r.y + m.zz*r.z);
966 : q0(0.0), q1(x), q2(y), q3(z)
971 : q0(qv[0]), q1(qv[1]), q2(qv[2]), q3(qv[3])
979 : q0(q0i), q1(q1i), q2(q2i), q3(q3i)
983 : q0(v[0]), q1(v[1]), q2(v[2]), q3(v[3])
1015 q0 = q1 = q2 = q3 = value;
1035 return 4*real_width + 13;
1038 std::string to_simple_string()
const;
1039 int from_simple_string(std::string
const &s);
1058 cvm::error(
"Error: incorrect quaternion component.\n");
1075 cvm::error(
"Error: trying to access a quaternion "
1076 "component which is not between 0 and 3.\n");
1094 return q0*q0 + q1*q1 + q2*q2 + q3*q3;
1111 q0 *= a; q1 *= a; q2 *= a; q3 *= a;
1116 q0 /= a; q1 /= a; q2 /= a; q3 /= a;
1119 inline void set_positive()
1121 if (q0 > 0.0)
return;
1130 q0+=h.q0; q1+=h.q1; q2+=h.q2; q3+=h.q3;
1134 q0-=h.q0; q1-=h.q1; q2-=h.q2; q3-=h.q3;
1162 h.q0*q.q1 + h.q1*q.q0 + h.q2*q.q3 - h.q3*q.q2,
1163 h.q0*q.q2 + h.q2*q.q0 + h.q3*q.q1 - h.q1*q.q3,
1164 h.q0*q.q3 + h.q3*q.q0 + h.q1*q.q2 - h.q2*q.q1);
1205 R.xx = q0*q0 + q1*q1 - q2*q2 - q3*q3;
1206 R.yy = q0*q0 - q1*q1 + q2*q2 - q3*q3;
1207 R.zz = q0*q0 - q1*q1 - q2*q2 + q3*q3;
1209 R.xy = 2.0 * (q1*q2 - q0*q3);
1210 R.xz = 2.0 * (q0*q2 + q1*q3);
1212 R.yx = 2.0 * (q0*q3 + q1*q2);
1213 R.yz = 2.0 * (q2*q3 - q0*q1);
1215 R.zx = 2.0 * (q1*q3 - q0*q2);
1216 R.zy = 2.0 * (q0*q1 + q2*q3);
1262 return cvm::quaternion(2.0 * (vec.x * ( q0 * pos.x - q3 * pos.y + q2 * pos.z) +
1263 vec.y * ( q3 * pos.x + q0 * pos.y - q1 * pos.z) +
1264 vec.z * (-q2 * pos.x + q1 * pos.y + q0 * pos.z)),
1265 2.0 * (vec.x * ( q1 * pos.x + q2 * pos.y + q3 * pos.z) +
1266 vec.y * ( q2 * pos.x - q1 * pos.y - q0 * pos.z) +
1267 vec.z * ( q3 * pos.x + q0 * pos.y - q1 * pos.z)),
1268 2.0 * (vec.x * (-q2 * pos.x + q1 * pos.y + q0 * pos.z) +
1269 vec.y * ( q1 * pos.x + q2 * pos.y + q3 * pos.z) +
1270 vec.z * (-q0 * pos.x + q3 * pos.y - q2 * pos.z)),
1271 2.0 * (vec.x * (-q3 * pos.x - q0 * pos.y + q1 * pos.z) +
1272 vec.y * ( q0 * pos.x - q3 * pos.y + q2 * pos.z) +
1273 vec.z * ( q1 * pos.x + q2 * pos.y + q3 * pos.z)));
1282 return 2.0*iprod*iprod - 1.0;
1290 cvm::real const cos_omega = this->q0*Q2.q0 + this->q1*Q2.q1 +
1291 this->q2*Q2.q2 + this->q3*Q2.q3;
1294 ( (cos_omega < -1.0) ? -1.0 : cos_omega) );
1297 if (cos_omega > 0.0)
1298 return omega * omega;
1300 return (PI-omega) * (PI-omega);
1307 cvm::real const cos_omega = this->q0*Q2.q0 + this->q1*Q2.q1 + this->q2*Q2.q2 + this->q3*Q2.q3;
1309 ( (cos_omega < -1.0) ? -1.0 : cos_omega) );
1318 grad1((-1.0)*sin_omega*Q2.q0 + cos_omega*(this->q0-cos_omega*Q2.q0)/sin_omega,
1319 (-1.0)*sin_omega*Q2.q1 + cos_omega*(this->q1-cos_omega*Q2.q1)/sin_omega,
1320 (-1.0)*sin_omega*Q2.q2 + cos_omega*(this->q2-cos_omega*Q2.q2)/sin_omega,
1321 (-1.0)*sin_omega*Q2.q3 + cos_omega*(this->q3-cos_omega*Q2.q3)/sin_omega);
1323 if (cos_omega > 0.0) {
1324 return 2.0*omega*grad1;
1326 return -2.0*(PI-omega)*grad1;
1334 cvm::real const cos_omega = this->q0*Q2.q0 + this->q1*Q2.q1 +
1335 this->q2*Q2.q2 + this->q3*Q2.q3;
1336 if (cos_omega < 0.0) Q2 *= -1.0;
1343 cvm::real const prod = this->q0*Q2.q0 + this->q1*Q2.q1 +
1344 this->q2*Q2.q2 + this->q3*Q2.q3;
1351#ifndef COLVARS_LAMMPS
1353int diagonalize_matrix(
cvm::real m[4][4],
1387 template <
typename T1,
typename T2>
1390 template<
typename T1,
typename T2>
1393 const std::vector<T1> &pos1,
1394 const std::vector<T2> &pos2);
1407 std::vector<atom_pos>
const &pos2);
1409 std::vector<atom_pos>
const &pos2);
1450 while (alpha > 180.0) alpha -= 360;
1451 while (alpha < -180.0) alpha += 360;
1465 (180.0/PI) * 2.0 * (1.0 / (1.0 + (iprod*iprod)/(
q.q0*
q.q0)));
1468 dspindx * ((1.0/
q.q0) * axis.x),
1469 dspindx * ((1.0/
q.q0) * axis.y),
1470 dspindx * ((1.0/
q.q0) * axis.z));
1475 return cvm::quaternion((180.0/PI) * 2.0 * ((-1.0)/iprod), 0.0, 0.0, 0.0);
1485 (180.0/PI) * 2.0 *
cvm::atan2(axis * q_vec,
q.q0);
1488 cvm::real const cos_theta_2 = ( (cos_spin_2 != 0.0) ?
1489 (
q.q0 / cos_spin_2) :
1492 return 2.0 * (cos_theta_2*cos_theta_2) - 1.0;
1506 (4.0 *
q.q0 / (cos_spin_2*cos_spin_2)) *
1507 (1.0 - (iprod*iprod)/(
q.q0*
q.q0) / (1.0 + (iprod*iprod)/(
q.q0*
q.q0)));
1510 (4.0 *
q.q0 / (cos_spin_2*cos_spin_2) *
1511 (iprod/
q.q0) / (1.0 + (iprod*iprod)/(
q.q0*
q.q0)));
1514 d_cos_theta_dqn * axis.x,
1515 d_cos_theta_dqn * axis.y,
1516 d_cos_theta_dqn * axis.z);
1520 (4.0 / (cos_spin_2*cos_spin_2 * iprod));
1523 d_cos_theta_dqn * axis.x,
1524 d_cos_theta_dqn * axis.y,
1525 d_cos_theta_dqn * axis.z);
1543 std::vector<cvm::atom_pos>
const &pos2);
1545 std::vector<cvm::atom_pos>
const &pos2);
Definition: colvartypes.h:385
Arbitrary size array (two dimensions) suitable for linear algebra operations (i.e....
Definition: colvartypes.h:376
matrix2d(matrix2d< T > const &m)
Copy constructor.
Definition: colvartypes.h:494
T ** c_array()
Return the 2-d C array.
Definition: colvartypes.h:543
void clear()
Deallocation routine.
Definition: colvartypes.h:463
~matrix2d()
Destructor.
Definition: colvartypes.h:504
std::vector< T > const & data_array() const
Return a reference to the data.
Definition: colvartypes.h:515
matrix2d< T > & operator=(matrix2d< T > const &m)
Assignment.
Definition: colvartypes.h:530
matrix2d()
Default constructor.
Definition: colvartypes.h:480
std::vector< T > & data_array()
Return a reference to the data.
Definition: colvartypes.h:509
void resize(size_t const ol, size_t const il)
Allocation routine, used by all constructors.
Definition: colvartypes.h:420
friend std::ostream & operator<<(std::ostream &os, matrix2d< T > const &m)
Formatted output.
Definition: colvartypes.h:670
void reset()
Set all elements to zero.
Definition: colvartypes.h:469
1-dimensional vector of real numbers with four components and a quaternion algebra
Definition: colvartypes.h:958
cvm::real & operator[](int i)
Access the quaternion as a 4-d array (return a reference)
Definition: colvartypes.h:1047
void set_from_euler_angles(cvm::real phi_in, cvm::real theta_in, cvm::real psi_in)
Definition: colvartypes.h:989
cvm::real cosine(cvm::quaternion const &q) const
Return the cosine between the orientation frame associated to this quaternion and another.
Definition: colvartypes.h:1279
cvm::rvector get_vector() const
Return the vector component.
Definition: colvartypes.h:1138
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:1194
cvm::real norm() const
Norm of the quaternion.
Definition: colvartypes.h:1098
quaternion(cvm::real const qv[4])
Constructor component by component.
Definition: colvartypes.h:970
quaternion()
Default constructor.
Definition: colvartypes.h:1007
void reset()
Set all components to zero (null quaternion)
Definition: colvartypes.h:1019
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:1186
cvm::quaternion dist2_grad(cvm::quaternion const &Q2) const
Definition: colvartypes.h:1305
cvm::quaternion position_derivative_inner(cvm::rvector const &pos, cvm::rvector const &vec) const
Multiply the given vector by the derivative of the given (rotated) position with respect to the quate...
Definition: colvartypes.h:1260
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:1332
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:1033
cvm::real norm2() const
Square norm of the quaternion.
Definition: colvartypes.h:1092
cvm::rmatrix rotation_matrix() const
Return the 3x3 matrix associated to this quaternion.
Definition: colvartypes.h:1201
cvm::quaternion conjugate() const
Return the conjugate quaternion.
Definition: colvartypes.h:1104
void set(cvm::real value)
Set all components to a scalar.
Definition: colvartypes.h:1013
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:1341
quaternion(cvm::real q0i, cvm::real q1i, cvm::real q2i, cvm::real q3i)
Constructor component by component.
Definition: colvartypes.h:975
quaternion(cvm::real x, cvm::real y, cvm::real z)
Constructor from a 3-d vector.
Definition: colvartypes.h:965
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:1158
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:1288
void reset_rotation()
Set the q0 component to 1 and the others to 0 (quaternion representing no rotation)
Definition: colvartypes.h:1026
2-dimensional array of real numbers with three components along each dimension (works with colvarmodu...
Definition: colvartypes.h:895
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:908
cvm::real determinant() const
Return the determinant.
Definition: colvartypes.h:930
rmatrix()
Default constructor.
Definition: colvartypes.h:902
A rotation between two sets of coordinates (for the moment a wrapper for colvarmodule::quaternion)
Definition: colvartypes.h:1363
~rotation()
Destructor.
Definition: colvartypes.cpp:198
cvm::rmatrix C
Correlation matrix C (3, 3)
Definition: colvartypes.h:1366
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:1481
cvm::real S_eigval[4]
Eigenvalues of S.
Definition: colvartypes.h:1372
cvm::real S_backup[4][4]
Used for debugging gradients.
Definition: colvartypes.h:1378
cvm::rvector rotate(cvm::rvector const &v) const
Return the rotated vector.
Definition: colvartypes.h:1427
int init()
Initialize member data.
Definition: colvartypes.cpp:151
cvm::real S[4][4]
Overlap matrix S (4, 4)
Definition: colvartypes.h:1369
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:1457
friend void debug_gradients(cvm::rotation &rot, const std::vector< T1 > &pos1, const std::vector< T2 > &pos2)
Function for debugging gradients (allow using either std::vector<cvm::atom_pos> or std::vector<cvm::a...
Definition: colvar_rotation_derivative.h:527
cvm::quaternion dcos_theta_dq(cvm::rvector const &axis) const
Return the derivative of the tilt wrt the quaternion.
Definition: colvartypes.h:1496
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:1446
cvm::real S_eigvec[4][4]
Eigenvectors of S.
Definition: colvartypes.h:1375
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...
static cvm::real crossing_threshold
Threshold for the eigenvalue crossing test.
Definition: colvartypes.h:1532
cvm::quaternion q_old
Previous value of the rotation (used to warn the user when the structure changes too much,...
Definition: colvartypes.h:1539
static bool monitor_crossings
Whether to test for eigenvalue crossing.
Definition: colvartypes.h:1530
bool b_debug_gradients
Perform gradient tests.
Definition: colvartypes.h:1382
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:1439
void * jacobi
Pointer to instance of Jacobi solver.
Definition: colvartypes.h:1554
void compute_overlap_matrix()
Compute the overlap matrix S (used by calc_optimal_rotation())
Definition: colvartypes.cpp:248
void calc_optimal_rotation_impl()
Actual implementation of calc_optimal_rotation (and called by it)
Definition: colvartypes.cpp:340
cvm::rotation inverse() const
Return the inverse of this rotation.
Definition: colvartypes.h:1433
cvm::quaternion q
The rotation itself (implemented as a quaternion)
Definition: colvartypes.h:1385
vector of real numbers with three components
Definition: colvartypes.h:727
cvm::real & operator[](int i)
Access cartesian components by index.
Definition: colvartypes.h:774
void set(cvm::real x_i, cvm::real y_i, cvm::real z_i)
Assign all components.
Definition: colvartypes.h:766
void reset()
Set all components to zero.
Definition: colvartypes.h:739
void set(cvm::real value)
Set all components to a scalar.
Definition: colvartypes.h:760
friend cvm::real operator*(cvm::rvector const &v1, cvm::rvector const &v2)
Inner (dot) product.
Definition: colvartypes.h:867
Arbitrary size array (one dimensions) suitable for linear algebra operations (i.e....
Definition: colvartypes.h:37
void reset()
Set all elements to zero.
Definition: colvartypes.h:96
size_t output_width(size_t real_width) const
Formatted output.
Definition: colvartypes.h:282
cvm::real norm2() const
Squared norm.
Definition: colvartypes.h:228
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:268
std::vector< T > & data_array()
Return a reference to the data.
Definition: colvartypes.h:79
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:52
vector1d(size_t const n=0)
Default constructor.
Definition: colvartypes.h:45
T * c_array()
Return a pointer to the data location.
Definition: colvartypes.h:69
vector1d< T > const slice(size_t const i1, size_t const i2) const
Slicing.
Definition: colvartypes.h:254
std::vector< T > const & data_array() const
Return a reference to the data.
Definition: colvartypes.h:85
double real
Defining an abstract real number allows to switch precision.
Definition: colvarmodule.h:95
static real cos(real const &x)
Reimplemented to work around MS compiler issues.
Definition: colvarmodule.h:145
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:2046
static real atan2(real const &x, real const &y)
Reimplemented to work around MS compiler issues.
Definition: colvarmodule.h:163
static size_t const cv_prec
Number of digits to represent a collective variables value(s)
Definition: colvarmodule.h:664
static real sqrt(real const &x)
Reimplemented to work around MS compiler issues.
Definition: colvarmodule.h:133
static real acos(real const &x)
Reimplemented to work around MS compiler issues.
Definition: colvarmodule.h:157
static real fabs(real const &x)
Reimplemented to work around MS compiler issues.
Definition: colvarmodule.h:127
static real sin(real const &x)
Reimplemented to work around MS compiler issues.
Definition: colvarmodule.h:139
static std::string to_str(char const *s)
Convert to string for output purposes.
Definition: colvarmodule.cpp:2392
Collective variables main module.
Helper class for calculating the derivative of rotation.
Definition: colvar_rotation_derivative.h:59