b202935ecba0c56ce0bb171fc260144c8d94c54d
[alexxy/gromacs.git] / src / gromacs / gmxana / princ.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include "config.h"
41
42 #include "gromacs/legacyheaders/typedefs.h"
43 #include "gromacs/math/vec.h"
44 #include "gromacs/utility/smalloc.h"
45 #include "gromacs/legacyheaders/txtdump.h"
46 #include "princ.h"
47
48 #include "gromacs/linearalgebra/nrjac.h"
49
50 static void m_op(matrix mat, rvec x)
51 {
52     rvec xp;
53     int  m;
54
55     for (m = 0; (m < DIM); m++)
56     {
57         xp[m] = mat[m][XX]*x[XX]+mat[m][YY]*x[YY]+mat[m][ZZ]*x[ZZ];
58     }
59     fprintf(stderr, "x    %8.3f  %8.3f  %8.3f\n", x[XX], x[YY], x[ZZ]);
60     fprintf(stderr, "xp   %8.3f  %8.3f  %8.3f\n", xp[XX], xp[YY], xp[ZZ]);
61     fprintf(stderr, "fac  %8.3f  %8.3f  %8.3f\n", xp[XX]/x[XX], xp[YY]/x[YY],
62             xp[ZZ]/x[ZZ]);
63 }
64
65 #define NDIM 4
66
67 #ifdef DEBUG
68 static void ptrans(char *s, real **inten, real d[], real e[])
69 {
70     int  m;
71     real n, x, y, z;
72     for (m = 1; (m < NDIM); m++)
73     {
74         x = inten[m][1];
75         y = inten[m][2];
76         z = inten[m][3];
77         n = x*x+y*y+z*z;
78         fprintf(stderr, "%8s %8.3f %8.3f %8.3f, norm:%8.3f, d:%8.3f, e:%8.3f\n",
79                 s, x, y, z, sqrt(n), d[m], e[m]);
80     }
81     fprintf(stderr, "\n");
82 }
83
84 void t_trans(matrix trans, real d[], real **ev)
85 {
86     rvec x;
87     int  j;
88     for (j = 0; (j < DIM); j++)
89     {
90         x[XX] = ev[1][j+1];
91         x[YY] = ev[2][j+1];
92         x[ZZ] = ev[3][j+1];
93         m_op(trans, x);
94         fprintf(stderr, "d[%d]=%g\n", j, d[j+1]);
95     }
96 }
97 #endif
98
99 void principal_comp(int n, atom_id index[], t_atom atom[], rvec x[],
100                     matrix trans, rvec d)
101 {
102     int      i, j, ai, m, nrot;
103     real     mm, rx, ry, rz;
104     double **inten, dd[NDIM], tvec[NDIM], **ev;
105 #ifdef DEBUG
106     real     e[NDIM];
107 #endif
108     real     temp;
109
110     snew(inten, NDIM);
111     snew(ev, NDIM);
112     for (i = 0; (i < NDIM); i++)
113     {
114         snew(inten[i], NDIM);
115         snew(ev[i], NDIM);
116         dd[i] = 0.0;
117 #ifdef DEBUG
118         e[i] = 0.0;
119 #endif
120     }
121
122     for (i = 0; (i < NDIM); i++)
123     {
124         for (m = 0; (m < NDIM); m++)
125         {
126             inten[i][m] = 0;
127         }
128     }
129     for (i = 0; (i < n); i++)
130     {
131         ai           = index[i];
132         mm           = atom[ai].m;
133         rx           = x[ai][XX];
134         ry           = x[ai][YY];
135         rz           = x[ai][ZZ];
136         inten[0][0] += mm*(sqr(ry)+sqr(rz));
137         inten[1][1] += mm*(sqr(rx)+sqr(rz));
138         inten[2][2] += mm*(sqr(rx)+sqr(ry));
139         inten[1][0] -= mm*(ry*rx);
140         inten[2][0] -= mm*(rx*rz);
141         inten[2][1] -= mm*(rz*ry);
142     }
143     inten[0][1] = inten[1][0];
144     inten[0][2] = inten[2][0];
145     inten[1][2] = inten[2][1];
146 #ifdef DEBUG
147     ptrans("initial", inten, dd, e);
148 #endif
149
150     for (i = 0; (i < DIM); i++)
151     {
152         for (m = 0; (m < DIM); m++)
153         {
154             trans[i][m] = inten[i][m];
155         }
156     }
157
158     /* Call numerical recipe routines */
159     jacobi(inten, 3, dd, ev, &nrot);
160 #ifdef DEBUG
161     ptrans("jacobi", ev, dd, e);
162 #endif
163
164     /* Sort eigenvalues in ascending order */
165 #define SWAPPER(i)          \
166     if (fabs(dd[i+1]) < fabs(dd[i])) {    \
167         temp = dd[i];         \
168         for (j = 0; (j < NDIM); j++) { tvec[j] = ev[j][i]; } \
169         dd[i] = dd[i+1];          \
170         for (j = 0; (j < NDIM); j++) { ev[j][i] = ev[j][i+1]; }        \
171         dd[i+1] = temp;           \
172         for (j = 0; (j < NDIM); j++) { ev[j][i+1] = tvec[j]; }         \
173     }
174     SWAPPER(0)
175     SWAPPER(1)
176     SWAPPER(0)
177 #ifdef DEBUG
178     ptrans("swap", ev, dd, e);
179     t_trans(trans, dd, ev);
180 #endif
181
182     for (i = 0; (i < DIM); i++)
183     {
184         d[i] = dd[i];
185         for (m = 0; (m < DIM); m++)
186         {
187             trans[i][m] = ev[m][i];
188         }
189     }
190
191     for (i = 0; (i < NDIM); i++)
192     {
193         sfree(inten[i]);
194         sfree(ev[i]);
195     }
196     sfree(inten);
197     sfree(ev);
198 }
199
200 void rotate_atoms(int gnx, atom_id *index, rvec x[], matrix trans)
201 {
202     real   xt, yt, zt;
203     int    i, ii;
204
205     for (i = 0; (i < gnx); i++)
206     {
207         ii        = index ? index[i] : i;
208         xt        = x[ii][XX];
209         yt        = x[ii][YY];
210         zt        = x[ii][ZZ];
211         x[ii][XX] = trans[XX][XX]*xt+trans[XX][YY]*yt+trans[XX][ZZ]*zt;
212         x[ii][YY] = trans[YY][XX]*xt+trans[YY][YY]*yt+trans[YY][ZZ]*zt;
213         x[ii][ZZ] = trans[ZZ][XX]*xt+trans[ZZ][YY]*yt+trans[ZZ][ZZ]*zt;
214     }
215 }
216
217 real calc_xcm(rvec x[], int gnx, atom_id *index, t_atom *atom, rvec xcm,
218               gmx_bool bQ)
219 {
220     int  i, ii, m;
221     real m0, tm;
222
223     clear_rvec(xcm);
224     tm = 0;
225     for (i = 0; (i < gnx); i++)
226     {
227         ii = index ? index[i] : i;
228         if (atom)
229         {
230             if (bQ)
231             {
232                 m0 = fabs(atom[ii].q);
233             }
234             else
235             {
236                 m0 = atom[ii].m;
237             }
238         }
239         else
240         {
241             m0 = 1;
242         }
243         tm += m0;
244         for (m = 0; (m < DIM); m++)
245         {
246             xcm[m] += m0*x[ii][m];
247         }
248     }
249     for (m = 0; (m < DIM); m++)
250     {
251         xcm[m] /= tm;
252     }
253
254     return tm;
255 }
256
257 real sub_xcm(rvec x[], int gnx, atom_id *index, t_atom atom[], rvec xcm,
258              gmx_bool bQ)
259 {
260     int  i, ii;
261     real tm;
262
263     tm = calc_xcm(x, gnx, index, atom, xcm, bQ);
264     for (i = 0; (i < gnx); i++)
265     {
266         ii = index ? index[i] : i;
267         rvec_dec(x[ii], xcm);
268     }
269     return tm;
270 }
271
272 void add_xcm(rvec x[], int gnx, atom_id *index, rvec xcm)
273 {
274     int  i, ii;
275
276     for (i = 0; (i < gnx); i++)
277     {
278         ii = index ? index[i] : i;
279         rvec_inc(x[ii], xcm);
280     }
281 }
282
283 void orient_princ(t_atoms *atoms, int isize, atom_id *index,
284                   int natoms, rvec x[], rvec *v, rvec d)
285 {
286     int     i, m;
287     rvec    xcm, prcomp;
288     matrix  trans;
289
290     calc_xcm(x, isize, index, atoms->atom, xcm, FALSE);
291     for (i = 0; i < natoms; i++)
292     {
293         rvec_dec(x[i], xcm);
294     }
295     principal_comp(isize, index, atoms->atom, x, trans, prcomp);
296     if (d)
297     {
298         copy_rvec(prcomp, d);
299     }
300
301     /* Check whether this trans matrix mirrors the molecule */
302     if (det(trans) < 0)
303     {
304         for (m = 0; (m < DIM); m++)
305         {
306             trans[ZZ][m] = -trans[ZZ][m];
307         }
308     }
309     rotate_atoms(natoms, NULL, x, trans);
310     if (v)
311     {
312         rotate_atoms(natoms, NULL, v, trans);
313     }
314
315     for (i = 0; i < natoms; i++)
316     {
317         rvec_inc(x[i], xcm);
318     }
319 }