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, 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.
53 real calc_similar_ind(gmx_bool bRho,int nind,atom_id *index,real mass[],
57 real m, tm, xs, xd, rs, rd;
62 for(j=0; j<nind; j++) {
69 for(d=0 ; d<DIM; d++) {
70 xd = x[i][d] - xp[i][d];
73 xs = x[i][d] + xp[i][d];
84 real rmsdev_ind(int nind,atom_id index[],real mass[],rvec x[],rvec xp[])
86 return calc_similar_ind(FALSE, nind, index, mass, x, xp);
89 real rmsdev(int natoms,real mass[],rvec x[],rvec xp[])
91 return calc_similar_ind(FALSE, natoms, NULL, mass, x, xp);
94 real rhodev_ind(int nind,atom_id index[],real mass[],rvec x[],rvec xp[])
96 return calc_similar_ind(TRUE, nind, index, mass, x, xp);
99 real rhodev(int natoms,real mass[],rvec x[],rvec xp[])
101 return calc_similar_ind(TRUE, natoms, NULL, mass, x, xp);
104 void calc_fit_R(int ndim,int natoms,real *w_rls,rvec *xp,rvec *x,matrix R)
106 int c,r,n,j,m,i,irot,s;
108 double d[2*DIM],xnr,xpc;
114 if (ndim != 3 && ndim != 2)
115 gmx_fatal(FARGS,"calc_fit_R called with ndim=%d instead of 3 or 2",ndim);
119 for(i=0; i<2*ndim; i++) {
120 snew(omega[i],2*ndim);
124 for(i=0; i<2*ndim; i++) {
126 for(j=0; j<2*ndim; j++) {
132 /*calculate the matrix U*/
134 for(n=0;(n<natoms);n++)
135 if ((mn = w_rls[n]) != 0.0)
136 for(c=0; (c<ndim); c++) {
138 for(r=0; (r<ndim); r++) {
145 /*omega is symmetric -> omega==omega' */
146 for(r=0; r<2*ndim; r++)
148 if (r>=ndim && c<ndim) {
149 omega[r][c]=u[r-ndim][c];
150 omega[c][r]=u[r-ndim][c];
156 /*determine h and k*/
157 jacobi(omega,2*ndim,d,om,&irot);
158 /*real **omega = input matrix a[0..n-1][0..n-1] must be symmetric
159 *int natoms = number of rows and columns
160 *real NULL = d[0]..d[n-1] are the eigenvalues of a[][]
161 *real **v = v[0..n-1][0..n-1] contains the vectors in columns
162 *int *irot = number of jacobi rotations
165 if (debug && irot==0) fprintf(debug,"IROT=0\n");
167 index=0; /* For the compiler only */
169 /* Copy only the first ndim-1 eigenvectors */
170 for(j=0; j<ndim-1; j++) {
172 for(i=0; i<2*ndim; i++)
178 for(i=0; i<ndim; i++) {
179 vh[j][i]=M_SQRT2*om[i][index];
180 vk[j][i]=M_SQRT2*om[i+ndim][index];
184 /* Calculate the last eigenvector as the outer-product of the first two.
185 * This insures that the conformation is not mirrored and
186 * prevents problems with completely flat reference structures.
188 cprod(vh[0],vh[1],vh[2]);
189 cprod(vk[0],vk[1],vk[2]);
190 } else if (ndim == 2) {
191 /* Calculate the last eigenvector from the first one */
192 vh[1][XX] = -vh[0][YY];
193 vh[1][YY] = vh[0][XX];
194 vk[1][XX] = -vk[0][YY];
195 vk[1][YY] = vk[0][XX];
200 for(r=0; r<ndim; r++)
201 for(c=0; c<ndim; c++)
202 for(s=0; s<ndim; s++)
203 R[r][c] += vk[s][r]*vh[s][c];
204 for(r=ndim; r<DIM; r++)
207 for(i=0; i<2*ndim; i++) {
215 void do_fit_ndim(int ndim,int natoms,real *w_rls,rvec *xp,rvec *x)
221 /* Calculate the rotation matrix R */
222 calc_fit_R(ndim,natoms,w_rls,xp,x,R);
225 for(j=0; j<natoms; j++) {
228 for(r=0; r<DIM; r++) {
231 x[j][r]+=R[r][c]*x_old[c];
236 void do_fit(int natoms,real *w_rls,rvec *xp,rvec *x)
238 do_fit_ndim(3,natoms,w_rls,xp,x);
241 void reset_x_ndim(int ndim,int ncm,const atom_id *ind_cm,
242 int nreset,const atom_id *ind_reset,
243 rvec x[],const real mass[])
251 gmx_incons("More than 3 dimensions not supported.");
261 for(m=0; m<ndim; m++)
263 xcm[m] += mm*x[ai][m];
273 for(m=0; m<ndim; m++)
275 xcm[m] += mm*x[i][m];
280 for(m=0; m<ndim; m++)
285 if (ind_reset != NULL)
287 for(i=0; i<nreset; i++)
289 rvec_dec(x[ind_reset[i]],xcm);
294 for(i=0; i<nreset; i++)
301 void reset_x(int ncm,const atom_id *ind_cm,
302 int nreset,const atom_id *ind_reset,
303 rvec x[],const real mass[])
305 reset_x_ndim(3,ncm,ind_cm,nreset,ind_reset,x,mass);