Update instructions in containers.rst
[gromacs.git] / src / gromacs / gmxana / princ.cpp
blob8f32afa7b342e5ecc3367c253a3b1417caeae03f
1 /*
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,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 /* This file is completely threadsafe - keep it that way! */
39 #include "gmxpre.h"
41 #include "princ.h"
43 #include <cmath>
45 #include "gromacs/linearalgebra/nrjac.h"
46 #include "gromacs/math/functions.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/topology/topology.h"
49 #include "gromacs/utility/smalloc.h"
51 #define NDIM 4
53 #ifdef DEBUG
54 static void m_op(matrix mat, rvec x)
56 rvec xp;
57 int m;
59 for (m = 0; (m < DIM); m++)
61 xp[m] = mat[m][XX] * x[XX] + mat[m][YY] * x[YY] + mat[m][ZZ] * x[ZZ];
63 fprintf(stderr, "x %8.3f %8.3f %8.3f\n", x[XX], x[YY], x[ZZ]);
64 fprintf(stderr, "xp %8.3f %8.3f %8.3f\n", xp[XX], xp[YY], xp[ZZ]);
65 fprintf(stderr, "fac %8.3f %8.3f %8.3f\n", xp[XX] / x[XX], xp[YY] / x[YY], xp[ZZ] / x[ZZ]);
68 static void ptrans(char* s, real** inten, real d[], real e[])
70 int m;
71 real n, x, y, z;
72 for (m = 1; (m < NDIM); m++)
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", s, x, y, z,
79 std::sqrt(n), d[m], e[m]);
81 fprintf(stderr, "\n");
84 void t_trans(matrix trans, real d[], real** ev)
86 rvec x;
87 int j;
88 for (j = 0; (j < DIM); j++)
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]);
97 #endif
99 void principal_comp(int n, const int index[], t_atom atom[], rvec x[], matrix trans, rvec d)
101 int i, j, ai, m, nrot;
102 real mm, rx, ry, rz;
103 double **inten, dd[NDIM], tvec[NDIM], **ev;
104 #ifdef DEBUG
105 real e[NDIM];
106 #endif
107 real temp;
109 snew(inten, NDIM);
110 snew(ev, NDIM);
111 for (i = 0; (i < NDIM); i++)
113 snew(inten[i], NDIM);
114 snew(ev[i], NDIM);
115 dd[i] = 0.0;
116 #ifdef DEBUG
117 e[i] = 0.0;
118 #endif
121 for (i = 0; (i < NDIM); i++)
123 for (m = 0; (m < NDIM); m++)
125 inten[i][m] = 0;
128 for (i = 0; (i < n); i++)
130 ai = index[i];
131 mm = atom[ai].m;
132 rx = x[ai][XX];
133 ry = x[ai][YY];
134 rz = x[ai][ZZ];
135 inten[0][0] += mm * (gmx::square(ry) + gmx::square(rz));
136 inten[1][1] += mm * (gmx::square(rx) + gmx::square(rz));
137 inten[2][2] += mm * (gmx::square(rx) + gmx::square(ry));
138 inten[1][0] -= mm * (ry * rx);
139 inten[2][0] -= mm * (rx * rz);
140 inten[2][1] -= mm * (rz * ry);
142 inten[0][1] = inten[1][0];
143 inten[0][2] = inten[2][0];
144 inten[1][2] = inten[2][1];
145 #ifdef DEBUG
146 ptrans("initial", inten, dd, e);
147 #endif
149 for (i = 0; (i < DIM); i++)
151 for (m = 0; (m < DIM); m++)
153 trans[i][m] = inten[i][m];
157 /* Call numerical recipe routines */
158 jacobi(inten, 3, dd, ev, &nrot);
159 #ifdef DEBUG
160 ptrans("jacobi", ev, dd, e);
161 #endif
163 /* Sort eigenvalues in ascending order */
164 #define SWAPPER(i) \
165 if (std::abs(dd[(i) + 1]) < std::abs(dd[i])) \
167 temp = dd[i]; \
168 for (j = 0; (j < NDIM); j++) \
170 tvec[j] = ev[j][i]; \
172 dd[i] = dd[(i) + 1]; \
173 for (j = 0; (j < NDIM); j++) \
175 ev[j][i] = ev[j][(i) + 1]; \
177 dd[(i) + 1] = temp; \
178 for (j = 0; (j < NDIM); j++) \
180 ev[j][(i) + 1] = tvec[j]; \
183 SWAPPER(0)
184 SWAPPER(1)
185 SWAPPER(0)
186 #ifdef DEBUG
187 ptrans("swap", ev, dd, e);
188 t_trans(trans, dd, ev);
189 #endif
191 for (i = 0; (i < DIM); i++)
193 d[i] = dd[i];
194 for (m = 0; (m < DIM); m++)
196 trans[i][m] = ev[m][i];
200 for (i = 0; (i < NDIM); i++)
202 sfree(inten[i]);
203 sfree(ev[i]);
205 sfree(inten);
206 sfree(ev);
209 void rotate_atoms(int gnx, const int* index, rvec x[], matrix trans)
211 real xt, yt, zt;
212 int i, ii;
214 for (i = 0; (i < gnx); i++)
216 ii = index ? index[i] : i;
217 xt = x[ii][XX];
218 yt = x[ii][YY];
219 zt = x[ii][ZZ];
220 x[ii][XX] = trans[XX][XX] * xt + trans[XX][YY] * yt + trans[XX][ZZ] * zt;
221 x[ii][YY] = trans[YY][XX] * xt + trans[YY][YY] * yt + trans[YY][ZZ] * zt;
222 x[ii][ZZ] = trans[ZZ][XX] * xt + trans[ZZ][YY] * yt + trans[ZZ][ZZ] * zt;
226 real calc_xcm(const rvec x[], int gnx, const int* index, const t_atom* atom, rvec xcm, gmx_bool bQ)
228 int i, ii, m;
229 real m0, tm;
231 clear_rvec(xcm);
232 tm = 0;
233 for (i = 0; (i < gnx); i++)
235 ii = index ? index[i] : i;
236 if (atom)
238 if (bQ)
240 m0 = std::abs(atom[ii].q);
242 else
244 m0 = atom[ii].m;
247 else
249 m0 = 1;
251 tm += m0;
252 for (m = 0; (m < DIM); m++)
254 xcm[m] += m0 * x[ii][m];
257 for (m = 0; (m < DIM); m++)
259 xcm[m] /= tm;
262 return tm;
265 real sub_xcm(rvec x[], int gnx, const int* index, const t_atom atom[], rvec xcm, gmx_bool bQ)
267 int i, ii;
268 real tm;
270 tm = calc_xcm(x, gnx, index, atom, xcm, bQ);
271 for (i = 0; (i < gnx); i++)
273 ii = index ? index[i] : i;
274 rvec_dec(x[ii], xcm);
276 return tm;
279 void add_xcm(rvec x[], int gnx, const int* index, rvec xcm)
281 int i, ii;
283 for (i = 0; (i < gnx); i++)
285 ii = index ? index[i] : i;
286 rvec_inc(x[ii], xcm);
290 void orient_princ(const t_atoms* atoms, int isize, const int* index, int natoms, rvec x[], rvec* v, rvec d)
292 int i, m;
293 rvec xcm, prcomp;
294 matrix trans;
296 calc_xcm(x, isize, index, atoms->atom, xcm, FALSE);
297 for (i = 0; i < natoms; i++)
299 rvec_dec(x[i], xcm);
301 principal_comp(isize, index, atoms->atom, x, trans, prcomp);
302 if (d)
304 copy_rvec(prcomp, d);
307 /* Check whether this trans matrix mirrors the molecule */
308 if (det(trans) < 0)
310 for (m = 0; (m < DIM); m++)
312 trans[ZZ][m] = -trans[ZZ][m];
315 rotate_atoms(natoms, nullptr, x, trans);
316 if (v)
318 rotate_atoms(natoms, nullptr, v, trans);
321 for (i = 0; i < natoms; i++)
323 rvec_inc(x[i], xcm);