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.
39 #include "gromacs/math/utilities.h"
45 #include "gromacs/linearalgebra/nrjac.h"
46 #include "gmx_fatal.h"
47 #include "gromacs/utility/smalloc.h"
49 real calc_similar_ind(gmx_bool bRho, int nind, atom_id *index, real mass[],
53 real m, tm, xs, xd, rs, rd;
58 for (j = 0; j < nind; j++)
70 for (d = 0; d < DIM; d++)
72 xd = x[i][d] - xp[i][d];
76 xs = x[i][d] + xp[i][d];
91 real rmsdev_ind(int nind, atom_id index[], real mass[], rvec x[], rvec xp[])
93 return calc_similar_ind(FALSE, nind, index, mass, x, xp);
96 real rmsdev(int natoms, real mass[], rvec x[], rvec xp[])
98 return calc_similar_ind(FALSE, natoms, NULL, mass, x, xp);
101 real rhodev_ind(int nind, atom_id index[], real mass[], rvec x[], rvec xp[])
103 return calc_similar_ind(TRUE, nind, index, mass, x, xp);
106 real rhodev(int natoms, real mass[], rvec x[], rvec xp[])
108 return calc_similar_ind(TRUE, natoms, NULL, mass, x, xp);
111 void calc_fit_R(int ndim, int natoms, real *w_rls, rvec *xp, rvec *x, matrix R)
113 int c, r, n, j, m, i, irot, s;
114 double **omega, **om;
115 double d[2*DIM], xnr, xpc;
121 if (ndim != 3 && ndim != 2)
123 gmx_fatal(FARGS, "calc_fit_R called with ndim=%d instead of 3 or 2", ndim);
128 for (i = 0; i < 2*ndim; i++)
130 snew(omega[i], 2*ndim);
134 for (i = 0; i < 2*ndim; i++)
137 for (j = 0; j < 2*ndim; j++)
144 /*calculate the matrix U*/
146 for (n = 0; (n < natoms); n++)
148 if ((mn = w_rls[n]) != 0.0)
150 for (c = 0; (c < ndim); c++)
153 for (r = 0; (r < ndim); r++)
156 u[c][r] += mn*xnr*xpc;
163 /*omega is symmetric -> omega==omega' */
164 for (r = 0; r < 2*ndim; r++)
166 for (c = 0; c <= r; c++)
168 if (r >= ndim && c < ndim)
170 omega[r][c] = u[r-ndim][c];
171 omega[c][r] = u[r-ndim][c];
181 /*determine h and k*/
182 jacobi(omega, 2*ndim, d, om, &irot);
183 /*real **omega = input matrix a[0..n-1][0..n-1] must be symmetric
184 * int natoms = number of rows and columns
185 * real NULL = d[0]..d[n-1] are the eigenvalues of a[][]
186 * real **v = v[0..n-1][0..n-1] contains the vectors in columns
187 * int *irot = number of jacobi rotations
190 if (debug && irot == 0)
192 fprintf(debug, "IROT=0\n");
195 index = 0; /* For the compiler only */
197 /* Copy only the first ndim-1 eigenvectors */
198 for (j = 0; j < ndim-1; j++)
201 for (i = 0; i < 2*ndim; i++)
210 for (i = 0; i < ndim; i++)
212 vh[j][i] = M_SQRT2*om[i][index];
213 vk[j][i] = M_SQRT2*om[i+ndim][index];
218 /* Calculate the last eigenvector as the outer-product of the first two.
219 * This insures that the conformation is not mirrored and
220 * prevents problems with completely flat reference structures.
222 cprod(vh[0], vh[1], vh[2]);
223 cprod(vk[0], vk[1], vk[2]);
227 /* Calculate the last eigenvector from the first one */
228 vh[1][XX] = -vh[0][YY];
229 vh[1][YY] = vh[0][XX];
230 vk[1][XX] = -vk[0][YY];
231 vk[1][YY] = vk[0][XX];
236 for (r = 0; r < ndim; r++)
238 for (c = 0; c < ndim; c++)
240 for (s = 0; s < ndim; s++)
242 R[r][c] += vk[s][r]*vh[s][c];
246 for (r = ndim; r < DIM; r++)
251 for (i = 0; i < 2*ndim; i++)
260 void do_fit_ndim(int ndim, int natoms, real *w_rls, rvec *xp, rvec *x)
266 /* Calculate the rotation matrix R */
267 calc_fit_R(ndim, natoms, w_rls, xp, x, R);
270 for (j = 0; j < natoms; j++)
272 for (m = 0; m < DIM; m++)
276 for (r = 0; r < DIM; r++)
279 for (c = 0; c < DIM; c++)
281 x[j][r] += R[r][c]*x_old[c];
287 void do_fit(int natoms, real *w_rls, rvec *xp, rvec *x)
289 do_fit_ndim(3, natoms, w_rls, xp, x);
292 void reset_x_ndim(int ndim, int ncm, const atom_id *ind_cm,
293 int nreset, const atom_id *ind_reset,
294 rvec x[], const real mass[])
302 gmx_incons("More than 3 dimensions not supported.");
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 != NULL)
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 atom_id *ind_cm,
353 int nreset, const atom_id *ind_reset,
354 rvec x[], const real mass[])
356 reset_x_ndim(3, ncm, ind_cm, nreset, ind_reset, x, mass);