Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / gmxana / fitahx.cpp
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,2015,2017,2018,2019, 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 "fitahx.h"
40
41 #include <cmath>
42
43 #include "gromacs/math/do_fit.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/utility/fatalerror.h"
46 #include "gromacs/utility/smalloc.h"
47
48 static void my_calc_xcm(int nbb, const int 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, const int 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, int allindex[], rvec x[], int nca, int caindex[], gmx_bool bFit)
76 {
77     static rvec* xref = nullptr;
78     static real* mass = nullptr;
79     const real   d    = 0.15;  /* Rise per residue (nm)    */
80     const real   tw   = 1.745; /* Twist per residue (rad)  */
81     const real   rad  = 0.23;  /* Radius of the helix (nm) */
82     real         phi0, trms, rms;
83     rvec         dx, xcm;
84     int          ai, i, nmass;
85
86     if (nca < 3)
87     {
88         gmx_fatal(FARGS, "Need at least 3 Calphas to fit to, (now %d)...\n", nca);
89     }
90
91     if (xref == nullptr)
92     {
93         snew(xref, natoms);
94         snew(mass, natoms);
95     }
96     phi0 = 0;
97     for (i = 0; (i < nca); i++)
98     {
99         ai           = caindex[i];
100         xref[ai][XX] = rad * std::cos(phi0);
101         xref[ai][YY] = rad * std::sin(phi0);
102         xref[ai][ZZ] = d * i;
103
104         /* Set the mass to select that this Calpha contributes to fitting */
105         mass[ai] = 10.0;
106 /*#define DEBUG*/
107 #ifdef DEBUG
108         fprintf(stderr, "%5d  %8.3f  %8.3f  %8.3f  %8.3f  %8.3f  %8.3f\n", ai, x[ai][XX], x[ai][YY],
109                 x[ai][ZZ], xref[ai][XX], xref[ai][YY], xref[ai][ZZ]);
110 #endif
111         phi0 += tw;
112     }
113
114     /* Center the referece around the origin */
115     my_calc_xcm(nca, caindex, xref, xcm);
116     my_sub_xcm(nca, caindex, xref, xcm);
117
118     if (bFit)
119     {
120         /* Now center the to-be-fitted coords around the origin */
121         my_calc_xcm(nca, caindex, x, xcm);
122         my_sub_xcm(nall, allindex, x, xcm);
123     }
124 #ifdef DEBUG
125     dump_ahx(nres, bb, xref, box, 0);
126 #endif
127
128     nmass = 0;
129     for (i = 0; (i < natoms); i++)
130     {
131         if (mass[i] > 0)
132         {
133             nmass++;
134         }
135     }
136     if (nmass != nca)
137     {
138         gmx_fatal(FARGS, "nmass=%d, nca=%d\n", nmass, nca);
139     }
140
141     /* Now call the fitting routine */
142     if (bFit)
143     {
144         do_fit(natoms, mass, xref, x);
145     }
146
147     /* Reset masses and calc rms */
148     trms = 0.0;
149     for (i = 0; (i < nres); i++)
150     {
151         ai = bb[i].CA;
152
153         if (mass[ai] > 0.0)
154         {
155             rvec_sub(x[ai], xref[ai], dx);
156             rms = iprod(dx, dx);
157             bb[i].rmsa += std::sqrt(rms);
158             bb[i].nrms++;
159             trms += rms;
160             mass[ai] = 0.0;
161         }
162     }
163     return std::sqrt(trms / nca);
164 }