3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * GROningen Mixture of Alchemy and Childrens' Stories
35 /* This file is completely threadsafe - keep it that way! */
48 static void m_op(matrix mat, rvec x)
53 for (m = 0; (m < DIM); m++)
55 xp[m] = mat[m][XX]*x[XX]+mat[m][YY]*x[YY]+mat[m][ZZ]*x[ZZ];
57 fprintf(stderr, "x %8.3f %8.3f %8.3f\n", x[XX], x[YY], x[ZZ]);
58 fprintf(stderr, "xp %8.3f %8.3f %8.3f\n", xp[XX], xp[YY], xp[ZZ]);
59 fprintf(stderr, "fac %8.3f %8.3f %8.3f\n", xp[XX]/x[XX], xp[YY]/x[YY],
66 static void ptrans(char *s, real **inten, real d[], real e[])
70 for (m = 1; (m < NDIM); m++)
76 fprintf(stderr, "%8s %8.3f %8.3f %8.3f, norm:%8.3f, d:%8.3f, e:%8.3f\n",
77 s, x, y, z, sqrt(n), d[m], e[m]);
79 fprintf(stderr, "\n");
82 void t_trans(matrix trans, real d[], real **ev)
86 for (j = 0; (j < DIM); j++)
92 fprintf(stderr, "d[%d]=%g\n", j, d[j+1]);
97 void principal_comp(int n, atom_id index[], t_atom atom[], rvec x[],
100 int i, j, ai, m, nrot;
102 double **inten, dd[NDIM], tvec[NDIM], **ev;
110 for (i = 0; (i < NDIM); i++)
112 snew(inten[i], NDIM);
120 for (i = 0; (i < NDIM); i++)
122 for (m = 0; (m < NDIM); m++)
127 for (i = 0; (i < n); i++)
134 inten[0][0] += mm*(sqr(ry)+sqr(rz));
135 inten[1][1] += mm*(sqr(rx)+sqr(rz));
136 inten[2][2] += mm*(sqr(rx)+sqr(ry));
137 inten[1][0] -= mm*(ry*rx);
138 inten[2][0] -= mm*(rx*rz);
139 inten[2][1] -= mm*(rz*ry);
141 inten[0][1] = inten[1][0];
142 inten[0][2] = inten[2][0];
143 inten[1][2] = inten[2][1];
145 ptrans("initial", inten, dd, e);
148 for (i = 0; (i < DIM); i++)
150 for (m = 0; (m < DIM); m++)
152 trans[i][m] = inten[i][m];
156 /* Call numerical recipe routines */
157 jacobi(inten, 3, dd, ev, &nrot);
159 ptrans("jacobi", ev, dd, e);
162 /* Sort eigenvalues in ascending order */
164 if (fabs(dd[i+1]) < fabs(dd[i])) { \
166 for (j = 0; (j < NDIM); j++) { tvec[j] = ev[j][i]; } \
168 for (j = 0; (j < NDIM); j++) { ev[j][i] = ev[j][i+1]; } \
170 for (j = 0; (j < NDIM); j++) { ev[j][i+1] = tvec[j]; } \
176 ptrans("swap", ev, dd, e);
177 t_trans(trans, dd, ev);
180 for (i = 0; (i < DIM); i++)
183 for (m = 0; (m < DIM); m++)
185 trans[i][m] = ev[m][i];
189 for (i = 0; (i < NDIM); i++)
198 void rotate_atoms(int gnx, atom_id *index, rvec x[], matrix trans)
203 for (i = 0; (i < gnx); i++)
205 ii = index ? index[i] : i;
209 x[ii][XX] = trans[XX][XX]*xt+trans[XX][YY]*yt+trans[XX][ZZ]*zt;
210 x[ii][YY] = trans[YY][XX]*xt+trans[YY][YY]*yt+trans[YY][ZZ]*zt;
211 x[ii][ZZ] = trans[ZZ][XX]*xt+trans[ZZ][YY]*yt+trans[ZZ][ZZ]*zt;
215 real calc_xcm(rvec x[], int gnx, atom_id *index, t_atom *atom, rvec xcm,
223 for (i = 0; (i < gnx); i++)
225 ii = index ? index[i] : i;
230 m0 = fabs(atom[ii].q);
242 for (m = 0; (m < DIM); m++)
244 xcm[m] += m0*x[ii][m];
247 for (m = 0; (m < DIM); m++)
255 real sub_xcm(rvec x[], int gnx, atom_id *index, t_atom atom[], rvec xcm,
261 tm = calc_xcm(x, gnx, index, atom, xcm, bQ);
262 for (i = 0; (i < gnx); i++)
264 ii = index ? index[i] : i;
265 rvec_dec(x[ii], xcm);
270 void add_xcm(rvec x[], int gnx, atom_id *index, rvec xcm)
274 for (i = 0; (i < gnx); i++)
276 ii = index ? index[i] : i;
277 rvec_inc(x[ii], xcm);
281 void orient_princ(t_atoms *atoms, int isize, atom_id *index,
282 int natoms, rvec x[], rvec *v, rvec d)
288 calc_xcm(x, isize, index, atoms->atom, xcm, FALSE);
289 for (i = 0; i < natoms; i++)
293 principal_comp(isize, index, atoms->atom, x, trans, prcomp);
296 copy_rvec(prcomp, d);
299 /* Check whether this trans matrix mirrors the molecule */
302 for (m = 0; (m < DIM); m++)
304 trans[ZZ][m] = -trans[ZZ][m];
307 rotate_atoms(natoms, NULL, x, trans);
310 rotate_atoms(natoms, NULL, v, trans);
313 for (i = 0; i < natoms; i++)