added Verlet scheme and NxN non-bonded functionality
[gromacs.git] / src / tools / fitahx.c
blob42b40d5660cd3a6c52ce289cf308c620f56f423e
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * Green Red Orange Magenta Azure Cyan Skyblue
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
39 #include "fitahx.h"
40 #include "vec.h"
41 #include "do_fit.h"
42 #include "smalloc.h"
44 static void my_calc_xcm(int nbb,atom_id bbind[],rvec x[],rvec xcm)
46 int i,m,ai;
48 clear_rvec(xcm);
49 for(i=0; (i<nbb); i++) {
50 ai=bbind[i];
51 rvec_inc(xcm,x[ai]);
53 for(m=0; (m<DIM); m++)
54 xcm[m]/=(nbb);
57 static void my_sub_xcm(int nbb,atom_id bbind[],rvec x[],rvec xcm)
59 int i,ai;
61 for(i=0; (i<nbb); i++) {
62 ai=bbind[i];
63 rvec_dec(x[ai],xcm);
67 real fit_ahx(int nres,t_bb bb[],int natoms,int nall,atom_id allindex[],
68 rvec x[],int nca,
69 atom_id caindex[],matrix box,gmx_bool bFit)
71 static rvec *xref=NULL;
72 static real *mass=NULL;
73 const real d=0.15; /* Rise per residue (nm) */
74 const real tw=1.745; /* Twist per residue (rad) */
75 const real rad=0.23; /* Radius of the helix (nm) */
76 real phi0,trms,rms;
77 rvec dx,xcm;
78 int ai,i,nmass;
80 if (nca < 3)
81 gmx_fatal(FARGS,"Need at least 3 Calphas to fit to, (now %d)...\n",nca);
83 if (xref == NULL) {
84 snew(xref,natoms);
85 snew(mass,natoms);
87 phi0=0;
88 for(i=0; (i<nca); i++) {
89 ai=caindex[i];
90 xref[ai][XX]=rad*cos(phi0);
91 xref[ai][YY]=rad*sin(phi0);
92 xref[ai][ZZ]=d*i;
94 /* Set the mass to select that this Calpha contributes to fitting */
95 mass[ai]=10.0;
96 /*#define DEBUG*/
97 #ifdef DEBUG
98 fprintf(stderr,"%5d %8.3f %8.3f %8.3f %8.3f %8.3f %8.3f\n",ai,
99 x[ai][XX],x[ai][YY],x[ai][ZZ],
100 xref[ai][XX],xref[ai][YY],xref[ai][ZZ]);
101 #endif
102 phi0+=tw;
105 /* Center the referece around the origin */
106 my_calc_xcm(nca,caindex,xref,xcm);
107 my_sub_xcm(nca,caindex,xref,xcm);
109 if (bFit) {
110 /* Now center the to-be-fitted coords around the origin */
111 my_calc_xcm(nca,caindex,x,xcm);
112 my_sub_xcm(nall,allindex,x,xcm);
114 #ifdef DEBUG
115 dump_ahx(nres,bb,xref,box,0);
116 #endif
118 nmass=0;
119 for(i=0; (i<natoms); i++)
120 if (mass[i] > 0)
121 nmass++;
122 if (nmass != nca)
123 gmx_fatal(FARGS,"nmass=%d, nca=%d\n",nmass,nca);
125 /* Now call the fitting routine */
126 if (bFit)
127 do_fit(natoms,mass,xref,x);
129 /* Reset masses and calc rms */
130 trms=0.0;
131 for(i=0; (i<nres); i++) {
132 ai=bb[i].CA;
134 if (mass[ai] > 0.0) {
135 rvec_sub(x[ai],xref[ai],dx);
136 rms=iprod(dx,dx);
137 bb[i].rmsa+=sqrt(rms);
138 bb[i].nrms++;
139 trms+=rms;
140 mass[ai]=0.0;
143 return sqrt(trms/nca);