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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * 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! */
51 static void m_op(matrix mat,rvec x)
56 for(m=0; (m<DIM); m++)
57 xp[m]=mat[m][XX]*x[XX]+mat[m][YY]*x[YY]+mat[m][ZZ]*x[ZZ];
58 fprintf(stderr,"x %8.3f %8.3f %8.3f\n",x[XX],x[YY],x[ZZ]);
59 fprintf(stderr,"xp %8.3f %8.3f %8.3f\n",xp[XX],xp[YY],xp[ZZ]);
60 fprintf(stderr,"fac %8.3f %8.3f %8.3f\n",xp[XX]/x[XX],xp[YY]/x[YY],
67 static void ptrans(char *s,real **inten,real d[],real e[])
71 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]);
82 void t_trans(matrix trans,real d[],real **ev)
86 for(j=0; (j<DIM); j++) {
91 fprintf(stderr,"d[%d]=%g\n",j,d[j+1]);
96 void principal_comp(int n,atom_id index[],t_atom atom[],rvec x[],
101 double **inten,dd[NDIM],tvec[NDIM],**ev;
109 for(i=0; (i<NDIM); i++) {
118 for(i=0; (i<NDIM); i++)
119 for(m=0; (m<NDIM); m++)
121 for(i=0; (i<n); i++) {
127 inten[0][0]+=mm*(sqr(ry)+sqr(rz));
128 inten[1][1]+=mm*(sqr(rx)+sqr(rz));
129 inten[2][2]+=mm*(sqr(rx)+sqr(ry));
130 inten[1][0]-=mm*(ry*rx);
131 inten[2][0]-=mm*(rx*rz);
132 inten[2][1]-=mm*(rz*ry);
134 inten[0][1]=inten[1][0];
135 inten[0][2]=inten[2][0];
136 inten[1][2]=inten[2][1];
138 ptrans("initial",inten,dd,e);
141 for(i=0; (i<DIM); i++) {
142 for(m=0; (m<DIM); m++)
143 trans[i][m]=inten[i][m];
146 /* Call numerical recipe routines */
147 jacobi(inten,3,dd,ev,&nrot);
149 ptrans("jacobi",ev,dd,e);
152 /* Sort eigenvalues in ascending order */
154 if (fabs(dd[i+1]) < fabs(dd[i])) { \
156 for(j=0; (j<NDIM); j++) tvec[j]=ev[j][i];\
158 for(j=0; (j<NDIM); j++) ev[j][i]=ev[j][i+1]; \
160 for(j=0; (j<NDIM); j++) ev[j][i+1]=tvec[j]; \
166 ptrans("swap",ev,dd,e);
167 t_trans(trans,dd,ev);
170 for(i=0; (i<DIM); i++) {
172 for(m=0; (m<DIM); m++)
173 trans[i][m]=ev[m][i];
176 for(i=0; (i<NDIM); i++) {
184 void rotate_atoms(int gnx,atom_id *index,rvec x[],matrix trans)
189 for(i=0; (i<gnx); i++) {
190 ii=index ? index[i] : i;
194 x[ii][XX]=trans[XX][XX]*xt+trans[XX][YY]*yt+trans[XX][ZZ]*zt;
195 x[ii][YY]=trans[YY][XX]*xt+trans[YY][YY]*yt+trans[YY][ZZ]*zt;
196 x[ii][ZZ]=trans[ZZ][XX]*xt+trans[ZZ][YY]*yt+trans[ZZ][ZZ]*zt;
200 real calc_xcm(rvec x[],int gnx,atom_id *index,t_atom *atom,rvec xcm,
208 for(i=0; (i<gnx); i++) {
209 ii=index ? index[i] : i;
219 for(m=0; (m<DIM); m++)
222 for(m=0; (m<DIM); m++)
228 real sub_xcm(rvec x[],int gnx,atom_id *index,t_atom atom[],rvec xcm,
234 tm=calc_xcm(x,gnx,index,atom,xcm,bQ);
235 for(i=0; (i<gnx); i++) {
236 ii=index ? index[i] : i;
242 void add_xcm(rvec x[],int gnx,atom_id *index,rvec xcm)
246 for(i=0; (i<gnx); i++) {
247 ii=index ? index[i] : i;
252 void orient_princ(t_atoms *atoms,int isize,atom_id *index,
253 int natoms, rvec x[], rvec *v, rvec d)
259 calc_xcm(x,isize,index,atoms->atom,xcm,FALSE);
260 for(i=0; i<natoms; i++)
262 principal_comp(isize,index,atoms->atom,x,trans,prcomp);
264 copy_rvec(prcomp, d);
266 /* Check whether this trans matrix mirrors the molecule */
267 if (det(trans) < 0) {
268 for(m=0; (m<DIM); m++)
269 trans[ZZ][m] = -trans[ZZ][m];
271 rotate_atoms(natoms,NULL,x,trans);
272 if (v) rotate_atoms(natoms,NULL,v,trans);
274 for(i=0; i<natoms; i++)