16700170b2a4e8e0189bd3ecd26532e141a212f1
[alexxy/gromacs.git] / src / gromacs / gmxana / fitahx.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "config.h"
40
41 #include "fitahx.h"
42 #include "gromacs/math/vec.h"
43
44 #include "gromacs/math/do_fit.h"
45 #include "gromacs/utility/fatalerror.h"
46 #include "gromacs/utility/smalloc.h"
47
48 static void my_calc_xcm(int nbb, atom_id bbind[], rvec x[], rvec xcm)
49 {
50     int    i, m, ai;
51
52     clear_rvec(xcm);
53     for (i = 0; (i < nbb); i++)
54     {
55         ai = bbind[i];
56         rvec_inc(xcm, x[ai]);
57     }
58     for (m = 0; (m < DIM); m++)
59     {
60         xcm[m] /= (nbb);
61     }
62 }
63
64 static void my_sub_xcm(int nbb, atom_id bbind[], rvec x[], rvec xcm)
65 {
66     int i, ai;
67
68     for (i = 0; (i < nbb); i++)
69     {
70         ai = bbind[i];
71         rvec_dec(x[ai], xcm);
72     }
73 }
74
75 real fit_ahx(int nres, t_bb bb[], int natoms, int nall, atom_id allindex[],
76              rvec x[], int nca,
77              atom_id caindex[], gmx_bool bFit)
78 {
79     static rvec *xref = NULL;
80     static real *mass = NULL;
81     const  real  d    = 0.15;  /* Rise per residue (nm)    */
82     const  real  tw   = 1.745; /* Twist per residue (rad)  */
83     const  real  rad  = 0.23;  /* Radius of the helix (nm) */
84     real         phi0, trms, rms;
85     rvec         dx, xcm;
86     int          ai, i, nmass;
87
88     if (nca < 3)
89     {
90         gmx_fatal(FARGS, "Need at least 3 Calphas to fit to, (now %d)...\n", nca);
91     }
92
93     if (xref == NULL)
94     {
95         snew(xref, natoms);
96         snew(mass, natoms);
97     }
98     phi0 = 0;
99     for (i = 0; (i < nca); i++)
100     {
101         ai           = caindex[i];
102         xref[ai][XX] = rad*cos(phi0);
103         xref[ai][YY] = rad*sin(phi0);
104         xref[ai][ZZ] = d*i;
105
106         /* Set the mass to select that this Calpha contributes to fitting */
107         mass[ai] = 10.0;
108 /*#define DEBUG*/
109 #ifdef DEBUG
110         fprintf(stderr, "%5d  %8.3f  %8.3f  %8.3f  %8.3f  %8.3f  %8.3f\n", ai,
111                 x[ai][XX], x[ai][YY], x[ai][ZZ],
112                 xref[ai][XX], xref[ai][YY], xref[ai][ZZ]);
113 #endif
114         phi0 += tw;
115     }
116
117     /* Center the referece around the origin */
118     my_calc_xcm(nca, caindex, xref, xcm);
119     my_sub_xcm(nca, caindex, xref, xcm);
120
121     if (bFit)
122     {
123         /* Now center the to-be-fitted coords around the origin */
124         my_calc_xcm(nca, caindex, x, xcm);
125         my_sub_xcm(nall, allindex, x, xcm);
126     }
127 #ifdef DEBUG
128     dump_ahx(nres, bb, xref, box, 0);
129 #endif
130
131     nmass = 0;
132     for (i = 0; (i < natoms); i++)
133     {
134         if (mass[i] > 0)
135         {
136             nmass++;
137         }
138     }
139     if (nmass != nca)
140     {
141         gmx_fatal(FARGS, "nmass=%d, nca=%d\n", nmass, nca);
142     }
143
144     /* Now call the fitting routine */
145     if (bFit)
146     {
147         do_fit(natoms, mass, xref, x);
148     }
149
150     /* Reset masses and calc rms */
151     trms = 0.0;
152     for (i = 0; (i < nres); i++)
153     {
154         ai = bb[i].CA;
155
156         if (mass[ai] > 0.0)
157         {
158             rvec_sub(x[ai], xref[ai], dx);
159             rms         = iprod(dx, dx);
160             bb[i].rmsa += sqrt(rms);
161             bb[i].nrms++;
162             trms    += rms;
163             mass[ai] = 0.0;
164         }
165     }
166     return sqrt(trms/nca);
167 }