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) 2012,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gromacs/linearalgebra/nrjac.h"
46 #include "gromacs/math/functions.h"
47 #include "gromacs/math/utilities.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/smalloc.h"
52 real calc_similar_ind(gmx_bool bRho, int nind, const int* index, const real mass[], rvec x[], rvec xp[])
55 real m, tm, xs, xd, rs, rd;
60 for (j = 0; j < nind; j++)
72 for (d = 0; d < DIM; d++)
74 xd = x[i][d] - xp[i][d];
75 rd += m * gmx::square(xd);
78 xs = x[i][d] + xp[i][d];
79 rs += m * gmx::square(xs);
85 return 2 * std::sqrt(rd / rs);
89 return std::sqrt(rd / tm);
93 real rmsdev_ind(int nind, int index[], real mass[], rvec x[], rvec xp[])
95 return calc_similar_ind(FALSE, nind, index, mass, x, xp);
98 real rmsdev(int natoms, real mass[], rvec x[], rvec xp[])
100 return calc_similar_ind(FALSE, natoms, nullptr, mass, x, xp);
103 real rhodev_ind(int nind, int index[], real mass[], rvec x[], rvec xp[])
105 return calc_similar_ind(TRUE, nind, index, mass, x, xp);
108 real rhodev(int natoms, real mass[], rvec x[], rvec xp[])
110 return calc_similar_ind(TRUE, natoms, nullptr, mass, x, xp);
113 void calc_fit_R(int ndim, int natoms, const real* w_rls, const rvec* xp, rvec* x, matrix R)
115 int c, r, n, j, i, irot, s;
116 double **omega, **om;
117 double d[2 * DIM], xnr, xpc;
123 if (ndim != 3 && ndim != 2)
125 gmx_fatal(FARGS, "calc_fit_R called with ndim=%d instead of 3 or 2", ndim);
128 snew(omega, 2 * ndim);
130 for (i = 0; i < 2 * ndim; i++)
132 snew(omega[i], 2 * ndim);
133 snew(om[i], 2 * ndim);
136 for (i = 0; i < 2 * ndim; i++)
139 for (j = 0; j < 2 * ndim; j++)
146 /*calculate the matrix U*/
148 for (n = 0; (n < natoms); n++)
150 if ((mn = w_rls[n]) != 0.0)
152 for (c = 0; (c < ndim); c++)
155 for (r = 0; (r < ndim); r++)
158 u[c][r] += mn * xnr * xpc;
165 /*omega is symmetric -> omega==omega' */
166 for (r = 0; r < 2 * ndim; r++)
168 for (c = 0; c <= r; c++)
170 if (r >= ndim && c < ndim)
172 omega[r][c] = u[r - ndim][c];
173 omega[c][r] = u[r - ndim][c];
183 /*determine h and k*/
184 jacobi(omega, 2 * ndim, d, om, &irot);
185 /*real **omega = input matrix a[0..n-1][0..n-1] must be symmetric
186 * int natoms = number of rows and columns
187 * real NULL = d[0]..d[n-1] are the eigenvalues of a[][]
188 * real **v = v[0..n-1][0..n-1] contains the vectors in columns
189 * int *irot = number of jacobi rotations
192 if (debug && irot == 0)
194 fprintf(debug, "IROT=0\n");
197 index = 0; /* For the compiler only */
199 /* Copy only the first ndim-1 eigenvectors */
200 for (j = 0; j < ndim - 1; j++)
203 for (i = 0; i < 2 * ndim; i++)
212 for (i = 0; i < ndim; i++)
214 vh[j][i] = M_SQRT2 * om[i][index];
215 vk[j][i] = M_SQRT2 * om[i + ndim][index];
220 /* Calculate the last eigenvector as the outer-product of the first two.
221 * This insures that the conformation is not mirrored and
222 * prevents problems with completely flat reference structures.
224 cprod(vh[0], vh[1], vh[2]);
225 cprod(vk[0], vk[1], vk[2]);
229 /* Calculate the last eigenvector from the first one */
230 vh[1][XX] = -vh[0][YY];
231 vh[1][YY] = vh[0][XX];
232 vk[1][XX] = -vk[0][YY];
233 vk[1][YY] = vk[0][XX];
238 for (r = 0; r < ndim; r++)
240 for (c = 0; c < ndim; c++)
242 for (s = 0; s < ndim; s++)
244 R[r][c] += vk[s][r] * vh[s][c];
248 for (r = ndim; r < DIM; r++)
253 for (i = 0; i < 2 * ndim; i++)
262 void do_fit_ndim(int ndim, int natoms, real* w_rls, const rvec* xp, rvec* x)
268 /* Calculate the rotation matrix R */
269 calc_fit_R(ndim, natoms, w_rls, xp, x, R);
272 for (j = 0; j < natoms; j++)
274 for (m = 0; m < DIM; m++)
278 for (r = 0; r < DIM; r++)
281 for (c = 0; c < DIM; c++)
283 x[j][r] += R[r][c] * x_old[c];
289 void do_fit(int natoms, real* w_rls, const rvec* xp, rvec* x)
291 do_fit_ndim(3, natoms, w_rls, xp, x);
294 void reset_x_ndim(int ndim, int ncm, const int* ind_cm, int nreset, const int* ind_reset, rvec x[], const real mass[])
302 gmx_incons("More than 3 dimensions not supported.");
306 if (ind_cm != nullptr)
308 for (i = 0; i < ncm; i++)
312 for (m = 0; m < ndim; m++)
314 xcm[m] += mm * x[ai][m];
321 for (i = 0; i < ncm; i++)
324 for (m = 0; m < ndim; m++)
326 xcm[m] += mm * x[i][m];
331 for (m = 0; m < ndim; m++)
336 if (ind_reset != nullptr)
338 for (i = 0; i < nreset; i++)
340 rvec_dec(x[ind_reset[i]], xcm);
345 for (i = 0; i < nreset; i++)
352 void reset_x(int ncm, const int* ind_cm, int nreset, const int* ind_reset, rvec x[], const real mass[])
354 reset_x_ndim(3, ncm, ind_cm, nreset, ind_reset, x, mass);