2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2018, 2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #ifndef GMX_MATH_VEC_H
38 #define GMX_MATH_VEC_H
40 /*! \brief Mathematical operations on (deprecated) rvec and matrix classes
42 * \todo Remove functions as Rvec replaces rvec and BasicMatrix3x3 replaces matrix
44 * \ingroup module_math
47 collection of in-line ready operations:
50 void rvec_add(const rvec a,const rvec b,rvec c) c = a + b
51 void dvec_add(const dvec a,const dvec b,dvec c) c = a + b
52 void rvec_inc(rvec a,const rvec b) a += b
53 void dvec_inc(dvec a,const dvec b) a += b
54 void ivec_inc(ivec a,const ivec b) a += b
55 void rvec_sub(const rvec a,const rvec b,rvec c) c = a - b
56 void dvec_sub(const dvec a,const dvec b,dvec c) c = a - b
57 void rvec_dec(rvec a,rvec b) a -= b
58 void copy_rvec(const rvec a,rvec b) b = a (reals)
59 void copy_dvec(const dvec a,dvec b) b = a (reals)
60 void copy_rvec_to_dvec(const rvec a,dvec b) b = a (reals)
61 void copy_dvec_to_rvec(const dvec a,rvec b) b = a (reals)
62 void copy_ivec(const ivec a,ivec b) b = a (integers)
63 void ivec_sub(const ivec a,const ivec b,ivec c) c = a - b
64 void svmul(real a,rvec v1,rvec v2) v2 = a * v1
65 void dsvmul(double a,dvec v1,dvec v2) v2 = a * v1
66 void clear_rvec(rvec a) a = 0
67 void clear_dvec(dvec a) a = 0
68 void clear_ivec(rvec a) a = 0
69 void clear_rvecs(int n,rvec v[])
70 real iprod(rvec a,rvec b) = a . b (inner product)
71 double diprod(dvec a,dvec b) = a . b (inner product)
72 real norm2(rvec a) = | a |^2 ( = x*y*z )
73 double dnorm2(dvec a) = | a |^2 ( = x*y*z )
74 real norm(rvec a) = | a |
75 double dnorm(dvec a) = | a |
76 void cprod(rvec a,rvec b,rvec c) c = a x b (cross product)
77 void dcprod(dvec a,dvec b,dvec c) c = a x b (cross product)
78 void dprod(rvec a,rvec b,rvec c) c = a * b (direct product)
79 real cos_angle(rvec a,rvec b)
80 real distance2(rvec v1, rvec v2) = | v2 - v1 |^2
81 void unitv(rvec src,rvec dest) dest = src / |src|
83 matrix (3x3) operations:
84 ! indicates that dest should not be the same as a, b or src
85 the _ur0 varieties work on matrices that have only zeros
86 in the upper right part, such as box matrices, these varieties
87 could produce less rounding errors, not due to the operations themselves,
88 but because the compiler can easier recombine the operations
89 void copy_mat(matrix a,matrix b) b = a
90 void clear_mat(matrix a) a = 0
91 void mmul(matrix a,matrix b,matrix dest) ! dest = a . b
92 void mmul_ur0(matrix a,matrix b,matrix dest) dest = a . b
93 void transpose(matrix src,matrix dest) ! dest = src*
94 void tmmul(matrix a,matrix b,matrix dest) ! dest = a* . b
95 void mtmul(matrix a,matrix b,matrix dest) ! dest = a . b*
96 real det(matrix a) = det(a)
97 void m_add(matrix a,matrix b,matrix dest) dest = a + b
98 void m_sub(matrix a,matrix b,matrix dest) dest = a - b
99 void msmul(matrix m1,real r1,matrix dest) dest = r1 * m1
100 void mvmul(matrix a,rvec src,rvec dest) ! dest = a . src
101 void mvmul_ur0(matrix a,rvec src,rvec dest) dest = a . src
102 void tmvmul_ur0(matrix a,rvec src,rvec dest) dest = a* . src
103 real trace(matrix m) = trace(m)
108 #include <type_traits>
110 #include "gromacs/math/functions.h"
111 #include "gromacs/math/vectypes.h"
112 #include "gromacs/utility/real.h"
114 static inline void rvec_add(const rvec a, const rvec b, rvec c)
127 static inline void ivec_add(const ivec a, const ivec b, ivec c)
140 static inline void rvec_inc(rvec a, const rvec b)
153 static inline void dvec_inc(dvec a, const dvec b)
166 static inline void rvec_sub(const rvec a, const rvec b, rvec c)
179 static inline void dvec_sub(const dvec a, const dvec b, dvec c)
192 static inline void rvec_dec(rvec a, const rvec b)
205 static inline void copy_rvec(const rvec a, rvec b)
212 static inline void copy_rvec_to_dvec(const rvec a, dvec b)
219 static inline void copy_dvec_to_rvec(const dvec a, rvec b)
221 b[XX] = static_cast<real>(a[XX]);
222 b[YY] = static_cast<real>(a[YY]);
223 b[ZZ] = static_cast<real>(a[ZZ]);
226 static inline void copy_rvecn(const rvec* a, rvec* b, int startn, int endn)
229 for (i = startn; i < endn; i++)
237 static inline void copy_dvec(const dvec a, dvec b)
244 static inline void copy_ivec(const ivec a, ivec b)
251 static inline void ivec_sub(const ivec a, const ivec b, ivec c)
264 static inline void copy_mat(const matrix a, matrix b)
266 copy_rvec(a[XX], b[XX]);
267 copy_rvec(a[YY], b[YY]);
268 copy_rvec(a[ZZ], b[ZZ]);
271 static inline void svmul(real a, const rvec v1, rvec v2)
278 static inline void dsvmul(double a, const dvec v1, dvec v2)
285 static inline real distance2(const rvec v1, const rvec v2)
287 return gmx::square(v2[XX] - v1[XX]) + gmx::square(v2[YY] - v1[YY]) + gmx::square(v2[ZZ] - v1[ZZ]);
290 static inline void clear_rvec(rvec a)
292 /* The ibm compiler has problems with inlining this
293 * when we use a const real variable
300 static inline void clear_dvec(dvec a)
302 /* The ibm compiler has problems with inlining this
303 * when we use a const real variable
310 static inline void clear_ivec(ivec a)
317 static inline void clear_rvecs(int n, rvec v[])
321 for (i = 0; (i < n); i++)
327 static inline void clear_mat(matrix a)
329 const real nul = 0.0;
331 a[XX][XX] = a[XX][YY] = a[XX][ZZ] = nul;
332 a[YY][XX] = a[YY][YY] = a[YY][ZZ] = nul;
333 a[ZZ][XX] = a[ZZ][YY] = a[ZZ][ZZ] = nul;
336 static inline real iprod(const rvec a, const rvec b)
338 return (a[XX] * b[XX] + a[YY] * b[YY] + a[ZZ] * b[ZZ]);
341 static inline double diprod(const dvec a, const dvec b)
343 return (a[XX] * b[XX] + a[YY] * b[YY] + a[ZZ] * b[ZZ]);
346 static inline real norm2(const rvec a)
348 return a[XX] * a[XX] + a[YY] * a[YY] + a[ZZ] * a[ZZ];
351 static inline double dnorm2(const dvec a)
353 return a[XX] * a[XX] + a[YY] * a[YY] + a[ZZ] * a[ZZ];
357 * As dnorm() uses sqrt() (which is slow) _only_ use it if you are sure you
358 * don't need 1/dnorm(), otherwise use dnorm2()*dinvnorm(). */
359 static inline double dnorm(const dvec a)
361 return std::sqrt(diprod(a, a));
365 * As norm() uses sqrt() (which is slow) _only_ use it if you are sure you
366 * don't need 1/norm(), otherwise use norm2()*invnorm(). */
367 static inline real norm(const rvec a)
369 return std::sqrt(iprod(a, a));
373 * Do _not_ use these routines to calculate the angle between two vectors
374 * as acos(cos_angle(u,v)). While it might seem obvious, the acos function
375 * is very flat close to -1 and 1, which will lead to accuracy-loss.
376 * Instead, use the new gmx_angle() function directly.
378 static inline real cos_angle(const rvec a, const rvec b)
381 * ax*bx + ay*by + az*bz
382 * cos-vec (a,b) = ---------------------
387 double aa, bb, ip, ipa, ipb, ipab; /* For accuracy these must be double! */
389 ip = ipa = ipb = 0.0;
390 for (m = 0; (m < DIM); m++) /* 18 */
401 cosval = static_cast<real>(ip * gmx::invsqrt(ipab)); /* 7 */
420 static inline void cprod(const rvec a, const rvec b, rvec c)
422 c[XX] = a[YY] * b[ZZ] - a[ZZ] * b[YY];
423 c[YY] = a[ZZ] * b[XX] - a[XX] * b[ZZ];
424 c[ZZ] = a[XX] * b[YY] - a[YY] * b[XX];
427 static inline void dcprod(const dvec a, const dvec b, dvec c)
429 c[XX] = a[YY] * b[ZZ] - a[ZZ] * b[YY];
430 c[YY] = a[ZZ] * b[XX] - a[XX] * b[ZZ];
431 c[ZZ] = a[XX] * b[YY] - a[YY] * b[XX];
434 /* This routine calculates the angle between a & b without any loss of accuracy close to 0/PI.
435 * If you only need cos(theta), use the cos_angle() routines to save a few cycles.
436 * This routine is faster than it might appear, since atan2 is accelerated on many CPUs (e.g. x86).
438 static inline real gmx_angle(const rvec a, const rvec b)
448 return std::atan2(wlen, s);
451 static inline double gmx_angle_between_dvecs(const dvec a, const dvec b)
461 return std::atan2(wlen, s);
464 static inline void mmul_ur0(const matrix a, const matrix b, matrix dest)
466 dest[XX][XX] = a[XX][XX] * b[XX][XX];
469 dest[YY][XX] = a[YY][XX] * b[XX][XX] + a[YY][YY] * b[YY][XX];
470 dest[YY][YY] = a[YY][YY] * b[YY][YY];
472 dest[ZZ][XX] = a[ZZ][XX] * b[XX][XX] + a[ZZ][YY] * b[YY][XX] + a[ZZ][ZZ] * b[ZZ][XX];
473 dest[ZZ][YY] = a[ZZ][YY] * b[YY][YY] + a[ZZ][ZZ] * b[ZZ][YY];
474 dest[ZZ][ZZ] = a[ZZ][ZZ] * b[ZZ][ZZ];
477 static inline void mmul(const matrix a, const matrix b, matrix dest)
479 dest[XX][XX] = a[XX][XX] * b[XX][XX] + a[XX][YY] * b[YY][XX] + a[XX][ZZ] * b[ZZ][XX];
480 dest[YY][XX] = a[YY][XX] * b[XX][XX] + a[YY][YY] * b[YY][XX] + a[YY][ZZ] * b[ZZ][XX];
481 dest[ZZ][XX] = a[ZZ][XX] * b[XX][XX] + a[ZZ][YY] * b[YY][XX] + a[ZZ][ZZ] * b[ZZ][XX];
482 dest[XX][YY] = a[XX][XX] * b[XX][YY] + a[XX][YY] * b[YY][YY] + a[XX][ZZ] * b[ZZ][YY];
483 dest[YY][YY] = a[YY][XX] * b[XX][YY] + a[YY][YY] * b[YY][YY] + a[YY][ZZ] * b[ZZ][YY];
484 dest[ZZ][YY] = a[ZZ][XX] * b[XX][YY] + a[ZZ][YY] * b[YY][YY] + a[ZZ][ZZ] * b[ZZ][YY];
485 dest[XX][ZZ] = a[XX][XX] * b[XX][ZZ] + a[XX][YY] * b[YY][ZZ] + a[XX][ZZ] * b[ZZ][ZZ];
486 dest[YY][ZZ] = a[YY][XX] * b[XX][ZZ] + a[YY][YY] * b[YY][ZZ] + a[YY][ZZ] * b[ZZ][ZZ];
487 dest[ZZ][ZZ] = a[ZZ][XX] * b[XX][ZZ] + a[ZZ][YY] * b[YY][ZZ] + a[ZZ][ZZ] * b[ZZ][ZZ];
490 static inline void transpose(const matrix src, matrix dest)
492 dest[XX][XX] = src[XX][XX];
493 dest[YY][XX] = src[XX][YY];
494 dest[ZZ][XX] = src[XX][ZZ];
495 dest[XX][YY] = src[YY][XX];
496 dest[YY][YY] = src[YY][YY];
497 dest[ZZ][YY] = src[YY][ZZ];
498 dest[XX][ZZ] = src[ZZ][XX];
499 dest[YY][ZZ] = src[ZZ][YY];
500 dest[ZZ][ZZ] = src[ZZ][ZZ];
503 static inline void tmmul(const matrix a, const matrix b, matrix dest)
505 /* Computes dest=mmul(transpose(a),b,dest) - used in do_pr_pcoupl */
506 dest[XX][XX] = a[XX][XX] * b[XX][XX] + a[YY][XX] * b[YY][XX] + a[ZZ][XX] * b[ZZ][XX];
507 dest[XX][YY] = a[XX][XX] * b[XX][YY] + a[YY][XX] * b[YY][YY] + a[ZZ][XX] * b[ZZ][YY];
508 dest[XX][ZZ] = a[XX][XX] * b[XX][ZZ] + a[YY][XX] * b[YY][ZZ] + a[ZZ][XX] * b[ZZ][ZZ];
509 dest[YY][XX] = a[XX][YY] * b[XX][XX] + a[YY][YY] * b[YY][XX] + a[ZZ][YY] * b[ZZ][XX];
510 dest[YY][YY] = a[XX][YY] * b[XX][YY] + a[YY][YY] * b[YY][YY] + a[ZZ][YY] * b[ZZ][YY];
511 dest[YY][ZZ] = a[XX][YY] * b[XX][ZZ] + a[YY][YY] * b[YY][ZZ] + a[ZZ][YY] * b[ZZ][ZZ];
512 dest[ZZ][XX] = a[XX][ZZ] * b[XX][XX] + a[YY][ZZ] * b[YY][XX] + a[ZZ][ZZ] * b[ZZ][XX];
513 dest[ZZ][YY] = a[XX][ZZ] * b[XX][YY] + a[YY][ZZ] * b[YY][YY] + a[ZZ][ZZ] * b[ZZ][YY];
514 dest[ZZ][ZZ] = a[XX][ZZ] * b[XX][ZZ] + a[YY][ZZ] * b[YY][ZZ] + a[ZZ][ZZ] * b[ZZ][ZZ];
517 static inline void mtmul(const matrix a, const matrix b, matrix dest)
519 /* Computes dest=mmul(a,transpose(b),dest) - used in do_pr_pcoupl */
520 dest[XX][XX] = a[XX][XX] * b[XX][XX] + a[XX][YY] * b[XX][YY] + a[XX][ZZ] * b[XX][ZZ];
521 dest[XX][YY] = a[XX][XX] * b[YY][XX] + a[XX][YY] * b[YY][YY] + a[XX][ZZ] * b[YY][ZZ];
522 dest[XX][ZZ] = a[XX][XX] * b[ZZ][XX] + a[XX][YY] * b[ZZ][YY] + a[XX][ZZ] * b[ZZ][ZZ];
523 dest[YY][XX] = a[YY][XX] * b[XX][XX] + a[YY][YY] * b[XX][YY] + a[YY][ZZ] * b[XX][ZZ];
524 dest[YY][YY] = a[YY][XX] * b[YY][XX] + a[YY][YY] * b[YY][YY] + a[YY][ZZ] * b[YY][ZZ];
525 dest[YY][ZZ] = a[YY][XX] * b[ZZ][XX] + a[YY][YY] * b[ZZ][YY] + a[YY][ZZ] * b[ZZ][ZZ];
526 dest[ZZ][XX] = a[ZZ][XX] * b[XX][XX] + a[ZZ][YY] * b[XX][YY] + a[ZZ][ZZ] * b[XX][ZZ];
527 dest[ZZ][YY] = a[ZZ][XX] * b[YY][XX] + a[ZZ][YY] * b[YY][YY] + a[ZZ][ZZ] * b[YY][ZZ];
528 dest[ZZ][ZZ] = a[ZZ][XX] * b[ZZ][XX] + a[ZZ][YY] * b[ZZ][YY] + a[ZZ][ZZ] * b[ZZ][ZZ];
531 static inline real det(const matrix a)
533 return (a[XX][XX] * (a[YY][YY] * a[ZZ][ZZ] - a[ZZ][YY] * a[YY][ZZ])
534 - a[YY][XX] * (a[XX][YY] * a[ZZ][ZZ] - a[ZZ][YY] * a[XX][ZZ])
535 + a[ZZ][XX] * (a[XX][YY] * a[YY][ZZ] - a[YY][YY] * a[XX][ZZ]));
539 static inline void m_add(const matrix a, const matrix b, matrix dest)
541 dest[XX][XX] = a[XX][XX] + b[XX][XX];
542 dest[XX][YY] = a[XX][YY] + b[XX][YY];
543 dest[XX][ZZ] = a[XX][ZZ] + b[XX][ZZ];
544 dest[YY][XX] = a[YY][XX] + b[YY][XX];
545 dest[YY][YY] = a[YY][YY] + b[YY][YY];
546 dest[YY][ZZ] = a[YY][ZZ] + b[YY][ZZ];
547 dest[ZZ][XX] = a[ZZ][XX] + b[ZZ][XX];
548 dest[ZZ][YY] = a[ZZ][YY] + b[ZZ][YY];
549 dest[ZZ][ZZ] = a[ZZ][ZZ] + b[ZZ][ZZ];
552 static inline void m_sub(const matrix a, const matrix b, matrix dest)
554 dest[XX][XX] = a[XX][XX] - b[XX][XX];
555 dest[XX][YY] = a[XX][YY] - b[XX][YY];
556 dest[XX][ZZ] = a[XX][ZZ] - b[XX][ZZ];
557 dest[YY][XX] = a[YY][XX] - b[YY][XX];
558 dest[YY][YY] = a[YY][YY] - b[YY][YY];
559 dest[YY][ZZ] = a[YY][ZZ] - b[YY][ZZ];
560 dest[ZZ][XX] = a[ZZ][XX] - b[ZZ][XX];
561 dest[ZZ][YY] = a[ZZ][YY] - b[ZZ][YY];
562 dest[ZZ][ZZ] = a[ZZ][ZZ] - b[ZZ][ZZ];
565 static inline void msmul(const matrix m1, real r1, matrix dest)
567 dest[XX][XX] = r1 * m1[XX][XX];
568 dest[XX][YY] = r1 * m1[XX][YY];
569 dest[XX][ZZ] = r1 * m1[XX][ZZ];
570 dest[YY][XX] = r1 * m1[YY][XX];
571 dest[YY][YY] = r1 * m1[YY][YY];
572 dest[YY][ZZ] = r1 * m1[YY][ZZ];
573 dest[ZZ][XX] = r1 * m1[ZZ][XX];
574 dest[ZZ][YY] = r1 * m1[ZZ][YY];
575 dest[ZZ][ZZ] = r1 * m1[ZZ][ZZ];
578 static inline void mvmul(const matrix a, const rvec src, rvec dest)
580 dest[XX] = a[XX][XX] * src[XX] + a[XX][YY] * src[YY] + a[XX][ZZ] * src[ZZ];
581 dest[YY] = a[YY][XX] * src[XX] + a[YY][YY] * src[YY] + a[YY][ZZ] * src[ZZ];
582 dest[ZZ] = a[ZZ][XX] * src[XX] + a[ZZ][YY] * src[YY] + a[ZZ][ZZ] * src[ZZ];
586 static inline void mvmul_ur0(const matrix a, const rvec src, rvec dest)
588 dest[ZZ] = a[ZZ][XX] * src[XX] + a[ZZ][YY] * src[YY] + a[ZZ][ZZ] * src[ZZ];
589 dest[YY] = a[YY][XX] * src[XX] + a[YY][YY] * src[YY];
590 dest[XX] = a[XX][XX] * src[XX];
593 static inline void tmvmul_ur0(const matrix a, const rvec src, rvec dest)
595 dest[XX] = a[XX][XX] * src[XX] + a[YY][XX] * src[YY] + a[ZZ][XX] * src[ZZ];
596 dest[YY] = a[YY][YY] * src[YY] + a[ZZ][YY] * src[ZZ];
597 dest[ZZ] = a[ZZ][ZZ] * src[ZZ];
600 static inline void unitv(const rvec src, rvec dest)
604 linv = gmx::invsqrt(norm2(src));
605 dest[XX] = linv * src[XX];
606 dest[YY] = linv * src[YY];
607 dest[ZZ] = linv * src[ZZ];
610 static inline real trace(const matrix m)
612 return (m[XX][XX] + m[YY][YY] + m[ZZ][ZZ]);
618 * \brief Forward operations on C Array style vectors to C implementations.
620 * Since vec.h and vectypes.h independently declare `norm` and `norm2` in
621 * different namespaces, code that includes both headers but does not specify
622 * the namespace from which to use `norm` and `norm2` cannot properly resolve
623 * overloads without the following helper templates.
624 * \tparam T array element type (e.g. real, int, etc.)
625 * \param v address of first vector element
626 * \return magnitude or squared magnitude of vector
630 std::remove_const_t<T> norm(T* v)
635 std::remove_const_t<T> norm2(T* v)