Fixed a bug in the pdb-writing code.
[gromacs.git] / src / tools / fitahx.c
blob885ec0468df713d98eb91610fcd4980edf8707a6
1 /*
2 * $Id$
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.1
11 * Copyright (c) 1991-2001, University of Groningen, The Netherlands
12 * This program is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU General Public License
14 * as published by the Free Software Foundation; either version 2
15 * of the License, or (at your option) any later version.
17 * If you want to redistribute modifications, please consider that
18 * scientific software is very special. Version control is crucial -
19 * bugs must be traceable. We will be happy to consider code for
20 * inclusion in the official distribution, but derived work must not
21 * be called official GROMACS. Details are found in the README & COPYING
22 * files - if they are missing, get the official version at www.gromacs.org.
24 * To help us fund GROMACS development, we humbly ask that you cite
25 * the papers on the package - you can find them in the top README file.
27 * For more info, check our website at http://www.gromacs.org
29 * And Hey:
30 * Gromacs Runs One Microsecond At Cannonball Speeds
32 static char *SRCID_fitahx_c = "$Id$";
33 #include "fitahx.h"
34 #include "vec.h"
35 #include "do_fit.h"
36 #include "smalloc.h"
38 static void my_calc_xcm(int nbb,atom_id bbind[],rvec x[],rvec xcm)
40 int i,m,ai;
42 clear_rvec(xcm);
43 for(i=0; (i<nbb); i++) {
44 ai=bbind[i];
45 rvec_inc(xcm,x[ai]);
47 for(m=0; (m<DIM); m++)
48 xcm[m]/=(nbb);
51 static void my_sub_xcm(int nbb,atom_id bbind[],rvec x[],rvec xcm)
53 int i,ai;
55 for(i=0; (i<nbb); i++) {
56 ai=bbind[i];
57 rvec_dec(x[ai],xcm);
61 real fit_ahx(int nres,t_bb bb[],int natoms,int nall,atom_id allindex[],
62 rvec x[],int nca,
63 atom_id caindex[],matrix box,bool bFit)
65 static rvec *xref=NULL;
66 static real *mass=NULL;
67 const real d=0.15; /* Rise per residue (nm) */
68 const real tw=1.745; /* Twist per residue (rad) */
69 const real rad=0.23; /* Radius of the helix (nm) */
70 real phi0,trms,rms;
71 rvec dx,xcm;
72 int ai,i,nmass;
74 if (nca < 3)
75 fatal_error(0,"Need at least 3 Calphas to fit to, (now %d)...\n",nca);
77 if (xref == NULL) {
78 snew(xref,natoms);
79 snew(mass,natoms);
81 phi0=0;
82 for(i=0; (i<nca); i++) {
83 ai=caindex[i];
84 xref[ai][XX]=rad*cos(phi0);
85 xref[ai][YY]=rad*sin(phi0);
86 xref[ai][ZZ]=d*i;
88 /* Set the mass to select that this Calpha contributes to fitting */
89 mass[ai]=10.0;
90 /*#define DEBUG*/
91 #ifdef DEBUG
92 fprintf(stderr,"%5d %8.3f %8.3f %8.3f %8.3f %8.3f %8.3f\n",ai,
93 x[ai][XX],x[ai][YY],x[ai][ZZ],
94 xref[ai][XX],xref[ai][YY],xref[ai][ZZ]);
95 #endif
96 phi0+=tw;
99 /* Center the referece around the origin */
100 my_calc_xcm(nca,caindex,xref,xcm);
101 my_sub_xcm(nca,caindex,xref,xcm);
103 if (bFit) {
104 /* Now center the to-be-fitted coords around the origin */
105 my_calc_xcm(nca,caindex,x,xcm);
106 my_sub_xcm(nall,allindex,x,xcm);
108 #ifdef DEBUG
109 dump_ahx(nres,bb,xref,box,0);
110 #endif
112 nmass=0;
113 for(i=0; (i<natoms); i++)
114 if (mass[i] > 0)
115 nmass++;
116 if (nmass != nca)
117 fatal_error(0,"nmass=%d, nca=%d\n",nmass,nca);
119 /* Now call the fitting routine */
120 if (bFit)
121 do_fit(natoms,mass,xref,x);
123 /* Reset masses and calc rms */
124 trms=0.0;
125 for(i=0; (i<nres); i++) {
126 ai=bb[i].CA;
128 if (mass[ai] > 0.0) {
129 rvec_sub(x[ai],xref[ai],dx);
130 rms=iprod(dx,dx);
131 bb[i].rmsa+=sqrt(rms);
132 bb[i].nrms++;
133 trms+=rms;
134 mass[ai]=0.0;
137 return sqrt(trms/nca);