2 * This file is part of the GROMACS molecular simulation package.
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 by the GROMACS development team.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
44 #include "gromacs/math/do_fit.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/smalloc.h"
49 static void my_calc_xcm(int nbb, const int bbind[], rvec x[], rvec xcm)
54 for (i = 0; (i < nbb); i++)
59 for (m = 0; (m < DIM); m++)
65 static void my_sub_xcm(int nbb, const int bbind[], rvec x[], rvec xcm)
69 for (i = 0; (i < nbb); i++)
76 real fit_ahx(int nres, t_bb bb[], int natoms, int nall, int allindex[], rvec x[], int nca, int caindex[], gmx_bool bFit)
78 static rvec* xref = nullptr;
79 static real* mass = nullptr;
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) */
89 gmx_fatal(FARGS, "Need at least 3 Calphas to fit to, (now %d)...\n", nca);
98 for (i = 0; (i < nca); i++)
101 xref[ai][XX] = rad * std::cos(phi0);
102 xref[ai][YY] = rad * std::sin(phi0);
103 xref[ai][ZZ] = d * i;
105 /* Set the mass to select that this Calpha contributes to fitting */
110 "%5d %8.3f %8.3f %8.3f %8.3f %8.3f %8.3f\n",
122 /* Center the referece around the origin */
123 my_calc_xcm(nca, caindex, xref, xcm);
124 my_sub_xcm(nca, caindex, xref, xcm);
128 /* Now center the to-be-fitted coords around the origin */
129 my_calc_xcm(nca, caindex, x, xcm);
130 my_sub_xcm(nall, allindex, x, xcm);
133 dump_ahx(nres, bb, xref, box, 0);
137 for (i = 0; (i < natoms); i++)
146 gmx_fatal(FARGS, "nmass=%d, nca=%d\n", nmass, nca);
149 /* Now call the fitting routine */
152 do_fit(natoms, mass, xref, x);
155 /* Reset masses and calc rms */
157 for (i = 0; (i < nres); i++)
163 rvec_sub(x[ai], xref[ai], dx);
165 bb[i].rmsa += std::sqrt(rms);
171 return std::sqrt(trms / nca);