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, 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.
41 #include "gromacs/linearalgebra/nrjac.h"
42 #include "gromacs/math/utilities.h"
43 #include "gromacs/math/vec.h"
44 #include "gromacs/utility/fatalerror.h"
45 #include "gromacs/utility/smalloc.h"
47 real calc_similar_ind(gmx_bool bRho, int nind, atom_id *index, real mass[],
51 real m, tm, xs, xd, rs, rd;
56 for (j = 0; j < nind; j++)
68 for (d = 0; d < DIM; d++)
70 xd = x[i][d] - xp[i][d];
74 xs = x[i][d] + xp[i][d];
89 real rmsdev_ind(int nind, atom_id index[], real mass[], rvec x[], rvec xp[])
91 return calc_similar_ind(FALSE, nind, index, mass, x, xp);
94 real rmsdev(int natoms, real mass[], rvec x[], rvec xp[])
96 return calc_similar_ind(FALSE, natoms, NULL, mass, x, xp);
99 real rhodev_ind(int nind, atom_id index[], real mass[], rvec x[], rvec xp[])
101 return calc_similar_ind(TRUE, nind, index, mass, x, xp);
104 real rhodev(int natoms, real mass[], rvec x[], rvec xp[])
106 return calc_similar_ind(TRUE, natoms, NULL, mass, x, xp);
109 void calc_fit_R(int ndim, int natoms, real *w_rls, rvec *xp, rvec *x, matrix R)
111 int c, r, n, j, m, i, irot, s;
112 double **omega, **om;
113 double d[2*DIM], xnr, xpc;
119 if (ndim != 3 && ndim != 2)
121 gmx_fatal(FARGS, "calc_fit_R called with ndim=%d instead of 3 or 2", ndim);
126 for (i = 0; i < 2*ndim; i++)
128 snew(omega[i], 2*ndim);
132 for (i = 0; i < 2*ndim; i++)
135 for (j = 0; j < 2*ndim; j++)
142 /*calculate the matrix U*/
144 for (n = 0; (n < natoms); n++)
146 if ((mn = w_rls[n]) != 0.0)
148 for (c = 0; (c < ndim); c++)
151 for (r = 0; (r < ndim); r++)
154 u[c][r] += mn*xnr*xpc;
161 /*omega is symmetric -> omega==omega' */
162 for (r = 0; r < 2*ndim; r++)
164 for (c = 0; c <= r; c++)
166 if (r >= ndim && c < ndim)
168 omega[r][c] = u[r-ndim][c];
169 omega[c][r] = u[r-ndim][c];
179 /*determine h and k*/
180 jacobi(omega, 2*ndim, d, om, &irot);
181 /*real **omega = input matrix a[0..n-1][0..n-1] must be symmetric
182 * int natoms = number of rows and columns
183 * real NULL = d[0]..d[n-1] are the eigenvalues of a[][]
184 * real **v = v[0..n-1][0..n-1] contains the vectors in columns
185 * int *irot = number of jacobi rotations
188 if (debug && irot == 0)
190 fprintf(debug, "IROT=0\n");
193 index = 0; /* For the compiler only */
195 /* Copy only the first ndim-1 eigenvectors */
196 for (j = 0; j < ndim-1; j++)
199 for (i = 0; i < 2*ndim; i++)
208 for (i = 0; i < ndim; i++)
210 vh[j][i] = M_SQRT2*om[i][index];
211 vk[j][i] = M_SQRT2*om[i+ndim][index];
216 /* Calculate the last eigenvector as the outer-product of the first two.
217 * This insures that the conformation is not mirrored and
218 * prevents problems with completely flat reference structures.
220 cprod(vh[0], vh[1], vh[2]);
221 cprod(vk[0], vk[1], vk[2]);
225 /* Calculate the last eigenvector from the first one */
226 vh[1][XX] = -vh[0][YY];
227 vh[1][YY] = vh[0][XX];
228 vk[1][XX] = -vk[0][YY];
229 vk[1][YY] = vk[0][XX];
234 for (r = 0; r < ndim; r++)
236 for (c = 0; c < ndim; c++)
238 for (s = 0; s < ndim; s++)
240 R[r][c] += vk[s][r]*vh[s][c];
244 for (r = ndim; r < DIM; r++)
249 for (i = 0; i < 2*ndim; i++)
258 void do_fit_ndim(int ndim, int natoms, real *w_rls, rvec *xp, rvec *x)
264 /* Calculate the rotation matrix R */
265 calc_fit_R(ndim, natoms, w_rls, xp, x, R);
268 for (j = 0; j < natoms; j++)
270 for (m = 0; m < DIM; m++)
274 for (r = 0; r < DIM; r++)
277 for (c = 0; c < DIM; c++)
279 x[j][r] += R[r][c]*x_old[c];
285 void do_fit(int natoms, real *w_rls, rvec *xp, rvec *x)
287 do_fit_ndim(3, natoms, w_rls, xp, x);
290 void reset_x_ndim(int ndim, int ncm, const atom_id *ind_cm,
291 int nreset, const atom_id *ind_reset,
292 rvec x[], const real mass[])
300 gmx_incons("More than 3 dimensions not supported.");
306 for (i = 0; i < ncm; i++)
310 for (m = 0; m < ndim; m++)
312 xcm[m] += mm*x[ai][m];
319 for (i = 0; i < ncm; i++)
322 for (m = 0; m < ndim; m++)
324 xcm[m] += mm*x[i][m];
329 for (m = 0; m < ndim; m++)
334 if (ind_reset != NULL)
336 for (i = 0; i < nreset; i++)
338 rvec_dec(x[ind_reset[i]], xcm);
343 for (i = 0; i < nreset; i++)
350 void reset_x(int ncm, const atom_id *ind_cm,
351 int nreset, const atom_id *ind_reset,
352 rvec x[], const real mass[])
354 reset_x_ndim(3, ncm, ind_cm, nreset, ind_reset, x, mass);