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, 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.
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.
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.
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.
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.
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.
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"
48 static void my_calc_xcm(int nbb
, const int bbind
[], rvec x
[], rvec xcm
)
53 for (i
= 0; (i
< nbb
); i
++)
58 for (m
= 0; (m
< DIM
); m
++)
64 static void my_sub_xcm(int nbb
, const int bbind
[], rvec x
[], rvec xcm
)
68 for (i
= 0; (i
< nbb
); i
++)
75 real
fit_ahx(int nres
, t_bb bb
[], int natoms
, int nall
, int allindex
[],
77 int caindex
[], gmx_bool bFit
)
79 static rvec
*xref
= nullptr;
80 static real
*mass
= nullptr;
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) */
90 gmx_fatal(FARGS
, "Need at least 3 Calphas to fit to, (now %d)...\n", nca
);
99 for (i
= 0; (i
< nca
); i
++)
102 xref
[ai
][XX
] = rad
*std::cos(phi0
);
103 xref
[ai
][YY
] = rad
*std::sin(phi0
);
106 /* Set the mass to select that this Calpha contributes to fitting */
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
]);
117 /* Center the referece around the origin */
118 my_calc_xcm(nca
, caindex
, xref
, xcm
);
119 my_sub_xcm(nca
, caindex
, xref
, xcm
);
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
);
128 dump_ahx(nres
, bb
, xref
, box
, 0);
132 for (i
= 0; (i
< natoms
); i
++)
141 gmx_fatal(FARGS
, "nmass=%d, nca=%d\n", nmass
, nca
);
144 /* Now call the fitting routine */
147 do_fit(natoms
, mass
, xref
, x
);
150 /* Reset masses and calc rms */
152 for (i
= 0; (i
< nres
); i
++)
158 rvec_sub(x
[ai
], xref
[ai
], dx
);
160 bb
[i
].rmsa
+= std::sqrt(rms
);
166 return std::sqrt(trms
/nca
);