Move smalloc.h to utility/
[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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include "fitahx.h"
42 #include "vec.h"
43 #include "gromacs/utility/smalloc.h"
44
45 #include "gromacs/math/do_fit.h"
46
47 static void my_calc_xcm(int nbb, atom_id bbind[], rvec x[], rvec xcm)
48 {
49     int    i, m, ai;
50
51     clear_rvec(xcm);
52     for (i = 0; (i < nbb); i++)
53     {
54         ai = bbind[i];
55         rvec_inc(xcm, x[ai]);
56     }
57     for (m = 0; (m < DIM); m++)
58     {
59         xcm[m] /= (nbb);
60     }
61 }
62
63 static void my_sub_xcm(int nbb, atom_id bbind[], rvec x[], rvec xcm)
64 {
65     int i, ai;
66
67     for (i = 0; (i < nbb); i++)
68     {
69         ai = bbind[i];
70         rvec_dec(x[ai], xcm);
71     }
72 }
73
74 real fit_ahx(int nres, t_bb bb[], int natoms, int nall, atom_id allindex[],
75              rvec x[], int nca,
76              atom_id caindex[], gmx_bool bFit)
77 {
78     static rvec *xref = NULL;
79     static real *mass = NULL;
80     const  real  d    = 0.15;  /* Rise per residue (nm)    */
81     const  real  tw   = 1.745; /* Twist per residue (rad)  */
82     const  real  rad  = 0.23;  /* Radius of the helix (nm) */
83     real         phi0, trms, rms;
84     rvec         dx, xcm;
85     int          ai, i, nmass;
86
87     if (nca < 3)
88     {
89         gmx_fatal(FARGS, "Need at least 3 Calphas to fit to, (now %d)...\n", nca);
90     }
91
92     if (xref == NULL)
93     {
94         snew(xref, natoms);
95         snew(mass, natoms);
96     }
97     phi0 = 0;
98     for (i = 0; (i < nca); i++)
99     {
100         ai           = caindex[i];
101         xref[ai][XX] = rad*cos(phi0);
102         xref[ai][YY] = rad*sin(phi0);
103         xref[ai][ZZ] = d*i;
104
105         /* Set the mass to select that this Calpha contributes to fitting */
106         mass[ai] = 10.0;
107 /*#define DEBUG*/
108 #ifdef DEBUG
109         fprintf(stderr, "%5d  %8.3f  %8.3f  %8.3f  %8.3f  %8.3f  %8.3f\n", ai,
110                 x[ai][XX], x[ai][YY], x[ai][ZZ],
111                 xref[ai][XX], xref[ai][YY], xref[ai][ZZ]);
112 #endif
113         phi0 += tw;
114     }
115
116     /* Center the referece around the origin */
117     my_calc_xcm(nca, caindex, xref, xcm);
118     my_sub_xcm(nca, caindex, xref, xcm);
119
120     if (bFit)
121     {
122         /* Now center the to-be-fitted coords around the origin */
123         my_calc_xcm(nca, caindex, x, xcm);
124         my_sub_xcm(nall, allindex, x, xcm);
125     }
126 #ifdef DEBUG
127     dump_ahx(nres, bb, xref, box, 0);
128 #endif
129
130     nmass = 0;
131     for (i = 0; (i < natoms); i++)
132     {
133         if (mass[i] > 0)
134         {
135             nmass++;
136         }
137     }
138     if (nmass != nca)
139     {
140         gmx_fatal(FARGS, "nmass=%d, nca=%d\n", nmass, nca);
141     }
142
143     /* Now call the fitting routine */
144     if (bFit)
145     {
146         do_fit(natoms, mass, xref, x);
147     }
148
149     /* Reset masses and calc rms */
150     trms = 0.0;
151     for (i = 0; (i < nres); i++)
152     {
153         ai = bb[i].CA;
154
155         if (mass[ai] > 0.0)
156         {
157             rvec_sub(x[ai], xref[ai], dx);
158             rms         = iprod(dx, dx);
159             bb[i].rmsa += sqrt(rms);
160             bb[i].nrms++;
161             trms    += rms;
162             mass[ai] = 0.0;
163         }
164     }
165     return sqrt(trms/nca);
166 }