3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35 /* This file is completely threadsafe - keep it that way! */
52 #include "gmx_fatal.h"
54 static void copy_bond (t_params *pr, int to, int from)
55 /* copies an entry in a bond list to another position.
56 * does no allocing or freeing of memory
59 /*memcpy((char*) &(pr->param[to]),(char*) &(pr->param[from]),
60 (size_t)sizeof(pr->param[from]));*/
65 range_check(to, 0, pr->nr);
66 range_check(from, 0, pr->nr);
67 for (i = 0; (i < MAXATOMLIST); i++)
69 pr->param[to].a[i] = pr->param[from].a[i];
71 for (i = 0; (i < MAXFORCEPARAM); i++)
73 pr->param[to].c[i] = pr->param[from].c[i];
75 for (i = 0; (i < MAXSLEN); i++)
77 pr->param[to].s[i] = pr->param[from].s[i];
82 static int count_hydrogens (char ***atomname, int nra, atom_id a[])
88 gmx_fatal(FARGS, "Cannot call count_hydrogens with no atomname (%s %d)",
93 for (i = 0; (i < nra); i++)
95 if (toupper(**(atomname[a[i]])) == 'H')
104 void make_shake (t_params plist[], t_atoms *atoms, int nshake)
106 char ***info = atoms->atomname;
109 t_param p, *bond, *ang;
111 int nb, b, i, j, ftype, ftype_a;
114 if (nshake != eshNONE)
119 printf("turning H bonds into constraints...\n");
122 printf("turning all bonds into constraints...\n");
125 printf("turning all bonds and H angles into constraints...\n");
128 printf("turning all bonds and angles into constraints...\n");
131 gmx_fatal(FARGS, "Invalid option for make_shake (%d)", nshake);
134 if ((nshake == eshHANGLES) || (nshake == eshALLANGLES))
136 /* Add all the angles with hydrogens to the shake list
137 * and remove them from the bond list
139 for (ftype = 0; (ftype < F_NRE); ftype++)
141 if (interaction_function[ftype].flags & IF_BTYPE)
143 bonds = &(plist[ftype]);
145 for (ftype_a = 0; (bonds->nr > 0 && ftype_a < F_NRE); ftype_a++)
147 if (interaction_function[ftype_a].flags & IF_ATYPE)
149 pr = &(plist[ftype_a]);
151 for (i = 0; (i < pr->nr); )
155 ang = &(pr->param[i]);
157 printf("Angle: %d-%d-%d\n", ang->AI, ang->AJ, ang->AK);
159 numhydrogens = count_hydrogens(info, 3, ang->a);
160 if ((nshake == eshALLANGLES) ||
161 (numhydrogens > 1) ||
162 (numhydrogens == 1 && toupper(**(info[ang->a[1]])) == 'O'))
164 /* Can only add hydrogen angle shake, if the two bonds
166 * append this angle to the shake list
171 /* Calculate length of constraint */
174 for (j = 0; !bFound && (j < bonds->nr); j++)
176 bond = &(bonds->param[j]);
177 if (((bond->AI == ang->AI) &&
178 (bond->AJ == ang->AJ)) ||
179 ((bond->AI == ang->AJ) &&
180 (bond->AJ == ang->AI)))
184 if (((bond->AI == ang->AK) &&
185 (bond->AJ == ang->AJ)) ||
186 ((bond->AI == ang->AJ) &&
187 (bond->AJ == ang->AK)))
191 bFound = (b_ij != 0.0) && (b_jk != 0.0);
195 /* apply law of cosines */
196 p.C0 = sqrt( b_ij*b_ij + b_jk*b_jk -
197 2.0*b_ij*b_jk*cos(DEG2RAD*ang->C0) );
200 printf("p: %d, q: %d, dist: %12.5e\n", p.AI, p.AJ, p.C0);
202 add_param_to_list (&(plist[F_CONSTR]), &p);
203 /* move the last bond to this position */
204 copy_bond (pr, i, pr->nr-1);
205 /* should free memory here!! */
218 } /* if shake angles */
220 /* Add all the bonds with hydrogens to the shake list
221 * and remove them from the bond list
223 for (ftype = 0; (ftype < F_NRE); ftype++)
225 if (interaction_function[ftype].flags & IF_BTYPE)
227 pr = &(plist[ftype]);
229 for (i = 0; i < pr->nr; i++)
231 if ( (nshake != eshHBONDS) ||
232 (count_hydrogens (info, 2, pr->param[i].a) > 0) )
234 /* append this bond to the shake list */
235 p.AI = pr->param[i].AI;
236 p.AJ = pr->param[i].AJ;
237 p.C0 = pr->param[i].C0;
238 p.C1 = pr->param[i].C2;
239 add_param_to_list (&(plist[F_CONSTR]), &p);
243 copy_bond(pr, j++, i);