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, 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.
37 /* This file is completely threadsafe - keep it that way! */
55 #include "gmx_fatal.h"
57 static void copy_bond (t_params *pr, int to, int from)
58 /* copies an entry in a bond list to another position.
59 * does no allocing or freeing of memory
62 /*memcpy((char*) &(pr->param[to]),(char*) &(pr->param[from]),
63 (size_t)sizeof(pr->param[from]));*/
68 range_check(to, 0, pr->nr);
69 range_check(from, 0, pr->nr);
70 for (i = 0; (i < MAXATOMLIST); i++)
72 pr->param[to].a[i] = pr->param[from].a[i];
74 for (i = 0; (i < MAXFORCEPARAM); i++)
76 pr->param[to].c[i] = pr->param[from].c[i];
78 for (i = 0; (i < MAXSLEN); i++)
80 pr->param[to].s[i] = pr->param[from].s[i];
85 static int count_hydrogens (char ***atomname, int nra, atom_id a[])
91 gmx_fatal(FARGS, "Cannot call count_hydrogens with no atomname (%s %d)",
96 for (i = 0; (i < nra); i++)
98 if (toupper(**(atomname[a[i]])) == 'H')
107 void make_shake (t_params plist[], t_atoms *atoms, int nshake)
109 char ***info = atoms->atomname;
112 t_param p, *bond, *ang;
114 int nb, b, i, j, ftype, ftype_a;
117 if (nshake != eshNONE)
122 printf("turning H bonds into constraints...\n");
125 printf("turning all bonds into constraints...\n");
128 printf("turning all bonds and H angles into constraints...\n");
131 printf("turning all bonds and angles into constraints...\n");
134 gmx_fatal(FARGS, "Invalid option for make_shake (%d)", nshake);
137 if ((nshake == eshHANGLES) || (nshake == eshALLANGLES))
139 /* Add all the angles with hydrogens to the shake list
140 * and remove them from the bond list
142 for (ftype = 0; (ftype < F_NRE); ftype++)
144 if (interaction_function[ftype].flags & IF_BTYPE)
146 bonds = &(plist[ftype]);
148 for (ftype_a = 0; (bonds->nr > 0 && ftype_a < F_NRE); ftype_a++)
150 if (interaction_function[ftype_a].flags & IF_ATYPE)
152 pr = &(plist[ftype_a]);
154 for (i = 0; (i < pr->nr); )
158 ang = &(pr->param[i]);
160 printf("Angle: %d-%d-%d\n", ang->AI, ang->AJ, ang->AK);
162 numhydrogens = count_hydrogens(info, 3, ang->a);
163 if ((nshake == eshALLANGLES) ||
164 (numhydrogens > 1) ||
165 (numhydrogens == 1 && toupper(**(info[ang->a[1]])) == 'O'))
167 /* Can only add hydrogen angle shake, if the two bonds
169 * append this angle to the shake list
174 /* Calculate length of constraint */
177 for (j = 0; !bFound && (j < bonds->nr); j++)
179 bond = &(bonds->param[j]);
180 if (((bond->AI == ang->AI) &&
181 (bond->AJ == ang->AJ)) ||
182 ((bond->AI == ang->AJ) &&
183 (bond->AJ == ang->AI)))
187 if (((bond->AI == ang->AK) &&
188 (bond->AJ == ang->AJ)) ||
189 ((bond->AI == ang->AJ) &&
190 (bond->AJ == ang->AK)))
194 bFound = (b_ij != 0.0) && (b_jk != 0.0);
198 /* apply law of cosines */
199 p.C0 = sqrt( b_ij*b_ij + b_jk*b_jk -
200 2.0*b_ij*b_jk*cos(DEG2RAD*ang->C0) );
203 printf("p: %d, q: %d, dist: %12.5e\n", p.AI, p.AJ, p.C0);
205 add_param_to_list (&(plist[F_CONSTR]), &p);
206 /* move the last bond to this position */
207 copy_bond (pr, i, pr->nr-1);
208 /* should free memory here!! */
221 } /* if shake angles */
223 /* Add all the bonds with hydrogens to the shake list
224 * and remove them from the bond list
226 for (ftype = 0; (ftype < F_NRE); ftype++)
228 if (interaction_function[ftype].flags & IF_BTYPE)
230 pr = &(plist[ftype]);
232 for (i = 0; i < pr->nr; i++)
234 if ( (nshake != eshHBONDS) ||
235 (count_hydrogens (info, 2, pr->param[i].a) > 0) )
237 /* append this bond to the shake list */
238 p.AI = pr->param[i].AI;
239 p.AJ = pr->param[i].AJ;
240 p.C0 = pr->param[i].C0;
241 p.C1 = pr->param[i].C2;
242 add_param_to_list (&(plist[F_CONSTR]), &p);
246 copy_bond(pr, j++, i);