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