Replaced all occurences of str(n)casecmp with gmx_str(n)casecmp.
[gromacs/rigid-bodies.git] / src / gmxlib / do_fit.c
blob74ea0f22bcb2cd02f9c3fd4b810803955074d4b2
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * GROningen Mixture of Alchemy and Childrens' Stories
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include "maths.h"
41 #include "sysstuff.h"
42 #include "typedefs.h"
43 #include "nrjac.h"
44 #include "vec.h"
45 #include "txtdump.h"
46 #include "smalloc.h"
47 #include "do_fit.h"
49 #define EPS 1.0e-09
51 real calc_similar_ind(bool bRho,int nind,atom_id *index,real mass[],
52 rvec x[],rvec xp[])
54 int i, j, d;
55 real m, tm, xs, xd, rs, rd;
57 tm=0;
58 rs=0;
59 rd=0;
60 for(j=0; j<nind; j++) {
61 if (index)
62 i = index[j];
63 else
64 i = j;
65 m = mass[i];
66 tm += m;
67 for(d=0 ; d<DIM; d++) {
68 xd = x[i][d] - xp[i][d];
69 rd += m * sqr(xd);
70 if (bRho) {
71 xs = x[i][d] + xp[i][d];
72 rs += m * sqr(xs);
76 if (bRho)
77 return 2*sqrt(rd/rs);
78 else
79 return sqrt(rd/tm);
82 real rmsdev_ind(int nind,atom_id index[],real mass[],rvec x[],rvec xp[])
84 return calc_similar_ind(FALSE, nind, index, mass, x, xp);
87 real rmsdev(int natoms,real mass[],rvec x[],rvec xp[])
89 return calc_similar_ind(FALSE, natoms, NULL, mass, x, xp);
92 real rhodev_ind(int nind,atom_id index[],real mass[],rvec x[],rvec xp[])
94 return calc_similar_ind(TRUE, nind, index, mass, x, xp);
97 real rhodev(int natoms,real mass[],rvec x[],rvec xp[])
99 return calc_similar_ind(TRUE, natoms, NULL, mass, x, xp);
102 void calc_fit_R(int ndim,int natoms,real *w_rls,rvec *xp,rvec *x,matrix R)
104 int c,r,n,j,m,i,irot,s;
105 double **omega,**om;
106 double d[2*DIM],xnr,xpc;
107 matrix vh,vk,u;
108 real mn;
109 int index;
110 real max_d;
112 if (ndim != 3 && ndim != 2)
113 gmx_fatal(FARGS,"calc_fit_R called with ndim=%d instead of 3 or 2",ndim);
115 snew(omega,2*ndim);
116 snew(om,2*ndim);
117 for(i=0; i<2*ndim; i++) {
118 snew(omega[i],2*ndim);
119 snew(om[i],2*ndim);
122 for(i=0; i<2*ndim; i++) {
123 d[i]=0;
124 for(j=0; j<2*ndim; j++) {
125 omega[i][j]=0;
126 om[i][j]=0;
130 /*calculate the matrix U*/
131 clear_mat(u);
132 for(n=0;(n<natoms);n++)
133 if ((mn = w_rls[n]) != 0.0)
134 for(c=0; (c<ndim); c++) {
135 xpc=xp[n][c];
136 for(r=0; (r<ndim); r++) {
137 xnr=x[n][r];
138 u[c][r]+=mn*xnr*xpc;
142 /*construct omega*/
143 /*omega is symmetric -> omega==omega' */
144 for(r=0; r<2*ndim; r++)
145 for(c=0; c<=r; c++)
146 if (r>=ndim && c<ndim) {
147 omega[r][c]=u[r-ndim][c];
148 omega[c][r]=u[r-ndim][c];
149 } else {
150 omega[r][c]=0;
151 omega[c][r]=0;
154 /*determine h and k*/
155 jacobi(omega,2*ndim,d,om,&irot);
156 /*real **omega = input matrix a[0..n-1][0..n-1] must be symmetric
157 *int natoms = number of rows and columns
158 *real NULL = d[0]..d[n-1] are the eigenvalues of a[][]
159 *real **v = v[0..n-1][0..n-1] contains the vectors in columns
160 *int *irot = number of jacobi rotations
163 if (debug && irot==0) fprintf(debug,"IROT=0\n");
165 index=0; /* For the compiler only */
167 /* Copy only the first ndim-1 eigenvectors */
168 for(j=0; j<ndim-1; j++) {
169 max_d=-1000;
170 for(i=0; i<2*ndim; i++)
171 if (d[i]>max_d) {
172 max_d=d[i];
173 index=i;
175 d[index]=-10000;
176 for(i=0; i<ndim; i++) {
177 vh[j][i]=M_SQRT2*om[i][index];
178 vk[j][i]=M_SQRT2*om[i+ndim][index];
181 if (ndim == 3) {
182 /* Calculate the last eigenvector as the outer-product of the first two.
183 * This insures that the conformation is not mirrored and
184 * prevents problems with completely flat reference structures.
186 cprod(vh[0],vh[1],vh[2]);
187 cprod(vk[0],vk[1],vk[2]);
188 } else if (ndim == 2) {
189 /* Calculate the last eigenvector from the first one */
190 vh[1][XX] = -vh[0][YY];
191 vh[1][YY] = vh[0][XX];
192 vk[1][XX] = -vk[0][YY];
193 vk[1][YY] = vk[0][XX];
196 /* determine R */
197 clear_mat(R);
198 for(r=0; r<ndim; r++)
199 for(c=0; c<ndim; c++)
200 for(s=0; s<ndim; s++)
201 R[r][c] += vk[s][r]*vh[s][c];
202 for(r=ndim; r<DIM; r++)
203 R[r][r] = 1;
205 for(i=0; i<2*ndim; i++) {
206 sfree(omega[i]);
207 sfree(om[i]);
209 sfree(omega);
210 sfree(om);
213 void do_fit_ndim(int ndim,int natoms,real *w_rls,rvec *xp,rvec *x)
215 int i,j,m,r,c;
216 matrix R;
217 rvec x_old;
219 /* Calculate the rotation matrix R */
220 calc_fit_R(ndim,natoms,w_rls,xp,x,R);
222 /*rotate X*/
223 for(j=0; j<natoms; j++) {
224 for(m=0; m<DIM; m++)
225 x_old[m]=x[j][m];
226 for(r=0; r<DIM; r++) {
227 x[j][r]=0;
228 for(c=0; c<DIM; c++)
229 x[j][r]+=R[r][c]*x_old[c];
234 void do_fit(int natoms,real *w_rls,rvec *xp,rvec *x)
236 do_fit_ndim(3,natoms,w_rls,xp,x);
239 void reset_x_ndim(int ndim,int ncm,const atom_id *ind_cm,
240 int nreset,const atom_id *ind_reset,
241 rvec x[],const real mass[])
243 int i,m,ai;
244 rvec xcm;
245 real tm,mm;
247 tm = 0.0;
248 clear_rvec(xcm);
249 if (ind_cm != NULL)
251 for(i=0; i<ncm; i++)
253 ai = ind_cm[i];
254 mm = mass[ai];
255 for(m=0; m<ndim; m++)
257 xcm[m] += mm*x[ai][m];
259 tm += mm;
262 else
264 for(i=0; i<ncm; i++)
266 mm = mass[i];
267 for(m=0; m<ndim; m++)
269 xcm[m] += mm*x[i][m];
271 tm += mm;
274 for(m=0; m<ndim; m++)
276 xcm[m] /= tm;
279 if (ind_reset != NULL)
281 for(i=0; i<nreset; i++)
283 rvec_dec(x[ind_reset[i]],xcm);
286 else
288 for(i=0; i<nreset; i++)
290 rvec_dec(x[i],xcm);
295 void reset_x(int ncm,const atom_id *ind_cm,
296 int nreset,const atom_id *ind_reset,
297 rvec x[],const real mass[])
299 reset_x_ndim(3,ncm,ind_cm,nreset,ind_reset,x,mass);