1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROningen Mixture of Alchemy and Childrens' Stories
51 real calc_similar_ind(gmx_bool bRho, int nind, atom_id *index, real mass[],
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];
78 xs = x[i][d] + xp[i][d];
93 real rmsdev_ind(int nind, atom_id 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, NULL, mass, x, xp);
103 real rhodev_ind(int nind, atom_id 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, NULL, mass, x, xp);
113 void calc_fit_R(int ndim, int natoms, real *w_rls, rvec *xp, rvec *x, matrix R)
115 int c, r, n, j, m, 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);
130 for (i = 0; i < 2*ndim; i++)
132 snew(omega[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, 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, rvec *xp, rvec *x)
291 do_fit_ndim(3, natoms, w_rls, xp, x);
294 void reset_x_ndim(int ndim, int ncm, const atom_id *ind_cm,
295 int nreset, const atom_id *ind_reset,
296 rvec x[], const real mass[])
304 gmx_incons("More than 3 dimensions not supported.");
310 for (i = 0; i < ncm; i++)
314 for (m = 0; m < ndim; m++)
316 xcm[m] += mm*x[ai][m];
323 for (i = 0; i < ncm; i++)
326 for (m = 0; m < ndim; m++)
328 xcm[m] += mm*x[i][m];
333 for (m = 0; m < ndim; m++)
338 if (ind_reset != NULL)
340 for (i = 0; i < nreset; i++)
342 rvec_dec(x[ind_reset[i]], xcm);
347 for (i = 0; i < nreset; i++)
354 void reset_x(int ncm, const atom_id *ind_cm,
355 int nreset, const atom_id *ind_reset,
356 rvec x[], const real mass[])
358 reset_x_ndim(3, ncm, ind_cm, nreset, ind_reset, x, mass);