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! */
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"
54 static void m_op(matrix mat
, rvec x
)
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
[])
72 for (m
= 1; (m
< NDIM
); m
++)
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
)
88 for (j
= 0; (j
< DIM
); j
++)
94 fprintf(stderr
, "d[%d]=%g\n", j
, d
[j
+ 1]);
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
;
103 double **inten
, dd
[NDIM
], tvec
[NDIM
], **ev
;
111 for (i
= 0; (i
< NDIM
); i
++)
113 snew(inten
[i
], NDIM
);
121 for (i
= 0; (i
< NDIM
); i
++)
123 for (m
= 0; (m
< NDIM
); m
++)
128 for (i
= 0; (i
< n
); i
++)
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];
146 ptrans("initial", inten
, dd
, e
);
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
);
160 ptrans("jacobi", ev
, dd
, e
);
163 /* Sort eigenvalues in ascending order */
165 if (std::abs(dd[(i) + 1]) < std::abs(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]; \
187 ptrans("swap", ev
, dd
, e
);
188 t_trans(trans
, dd
, ev
);
191 for (i
= 0; (i
< DIM
); i
++)
194 for (m
= 0; (m
< DIM
); m
++)
196 trans
[i
][m
] = ev
[m
][i
];
200 for (i
= 0; (i
< NDIM
); i
++)
209 void rotate_atoms(int gnx
, const int* index
, rvec x
[], matrix trans
)
214 for (i
= 0; (i
< gnx
); i
++)
216 ii
= index
? index
[i
] : i
;
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
)
233 for (i
= 0; (i
< gnx
); i
++)
235 ii
= index
? index
[i
] : i
;
240 m0
= std::abs(atom
[ii
].q
);
252 for (m
= 0; (m
< DIM
); m
++)
254 xcm
[m
] += m0
* x
[ii
][m
];
257 for (m
= 0; (m
< DIM
); m
++)
265 real
sub_xcm(rvec x
[], int gnx
, const int* index
, const t_atom atom
[], rvec xcm
, gmx_bool bQ
)
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
);
279 void add_xcm(rvec x
[], int gnx
, const int* index
, rvec xcm
)
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
)
296 calc_xcm(x
, isize
, index
, atoms
->atom
, xcm
, FALSE
);
297 for (i
= 0; i
< natoms
; i
++)
301 principal_comp(isize
, index
, atoms
->atom
, x
, trans
, prcomp
);
304 copy_rvec(prcomp
, d
);
307 /* Check whether this trans matrix mirrors the molecule */
310 for (m
= 0; (m
< DIM
); m
++)
312 trans
[ZZ
][m
] = -trans
[ZZ
][m
];
315 rotate_atoms(natoms
, nullptr, x
, trans
);
318 rotate_atoms(natoms
, nullptr, v
, trans
);
321 for (i
= 0; i
< natoms
; i
++)