Removed include of types/ifunc.h from typedefs.h
[gromacs.git] / src / gromacs / gmxpreprocess / topshake.cpp
blob587cce4dbe5f93c980717617a47133c4f00d24ab
1 /*
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, 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! */
38 #include "gmxpre.h"
40 #include "topshake.h"
42 #include <ctype.h>
44 #include <cmath>
46 #include "gromacs/gmxpreprocess/readir.h"
47 #include "gromacs/gmxpreprocess/topdirs.h"
48 #include "gromacs/gmxpreprocess/toppush.h"
49 #include "gromacs/gmxpreprocess/toputil.h"
50 #include "gromacs/legacyheaders/typedefs.h"
51 #include "gromacs/legacyheaders/types/ifunc.h"
52 #include "gromacs/math/units.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/smalloc.h"
56 static void copy_bond (t_params *pr, int to, int from)
57 /* copies an entry in a bond list to another position.
58 * does no allocing or freeing of memory
61 /*memcpy((char*) &(pr->param[to]),(char*) &(pr->param[from]),
62 (size_t)sizeof(pr->param[from]));*/
63 int i;
65 if (to != from)
67 range_check(to, 0, pr->nr);
68 range_check(from, 0, pr->nr);
69 for (i = 0; (i < MAXATOMLIST); i++)
71 pr->param[to].a[i] = pr->param[from].a[i];
73 for (i = 0; (i < MAXFORCEPARAM); i++)
75 pr->param[to].c[i] = pr->param[from].c[i];
77 for (i = 0; (i < MAXSLEN); i++)
79 pr->param[to].s[i] = pr->param[from].s[i];
84 static int count_hydrogens (char ***atomname, int nra, atom_id a[])
86 int i, nh;
88 if (!atomname)
90 gmx_fatal(FARGS, "Cannot call count_hydrogens with no atomname (%s %d)",
91 __FILE__, __LINE__);
94 nh = 0;
95 for (i = 0; (i < nra); i++)
97 if (toupper(**(atomname[a[i]])) == 'H')
99 nh++;
103 return nh;
106 void make_shake (t_params plist[], t_atoms *atoms, int nshake)
108 char ***info = atoms->atomname;
109 t_params *pr;
110 t_params *bonds;
111 t_param p, *bond, *ang;
112 real b_ij, b_jk;
113 int i, j, ftype, ftype_a;
114 gmx_bool bFound;
116 if (nshake != eshNONE)
118 switch (nshake)
120 case eshHBONDS:
121 printf("turning H bonds into constraints...\n");
122 break;
123 case eshALLBONDS:
124 printf("turning all bonds into constraints...\n");
125 break;
126 case eshHANGLES:
127 printf("turning all bonds and H angles into constraints...\n");
128 break;
129 case eshALLANGLES:
130 printf("turning all bonds and angles into constraints...\n");
131 break;
132 default:
133 gmx_fatal(FARGS, "Invalid option for make_shake (%d)", nshake);
136 if ((nshake == eshHANGLES) || (nshake == eshALLANGLES))
138 /* Add all the angles with hydrogens to the shake list
139 * and remove them from the bond list
141 for (ftype = 0; (ftype < F_NRE); ftype++)
143 if (interaction_function[ftype].flags & IF_BTYPE)
145 bonds = &(plist[ftype]);
147 for (ftype_a = 0; (bonds->nr > 0 && ftype_a < F_NRE); ftype_a++)
149 if (interaction_function[ftype_a].flags & IF_ATYPE)
151 pr = &(plist[ftype_a]);
153 for (i = 0; (i < pr->nr); )
155 int numhydrogens;
157 ang = &(pr->param[i]);
158 #ifdef DEBUG
159 printf("Angle: %d-%d-%d\n", ang->ai(), ang->aj(), ang->ak());
160 #endif
161 numhydrogens = count_hydrogens(info, 3, ang->a);
162 if ((nshake == eshALLANGLES) ||
163 (numhydrogens > 1) ||
164 (numhydrogens == 1 && toupper(**(info[ang->a[1]])) == 'O'))
166 /* Can only add hydrogen angle shake, if the two bonds
167 * are constrained.
168 * append this angle to the shake list
170 p.ai() = ang->ai();
171 p.aj() = ang->ak();
173 /* Calculate length of constraint */
174 bFound = FALSE;
175 b_ij = b_jk = 0.0;
176 for (j = 0; !bFound && (j < bonds->nr); j++)
178 bond = &(bonds->param[j]);
179 if (((bond->ai() == ang->ai()) &&
180 (bond->aj() == ang->aj())) ||
181 ((bond->ai() == ang->aj()) &&
182 (bond->aj() == ang->ai())))
184 b_ij = bond->c0();
186 if (((bond->ai() == ang->ak()) &&
187 (bond->aj() == ang->aj())) ||
188 ((bond->ai() == ang->aj()) &&
189 (bond->aj() == ang->ak())))
191 b_jk = bond->c0();
193 bFound = (b_ij != 0.0) && (b_jk != 0.0);
195 if (bFound)
197 /* apply law of cosines */
198 p.c0() = std::sqrt( b_ij*b_ij + b_jk*b_jk -
199 2.0*b_ij*b_jk*cos(DEG2RAD*ang->c0()) );
200 p.c1() = p.c0();
201 #ifdef DEBUG
202 printf("p: %d, q: %d, dist: %12.5e\n", p.ai(), p.aj(), p.c0());
203 #endif
204 add_param_to_list (&(plist[F_CONSTR]), &p);
205 /* move the last bond to this position */
206 copy_bond (pr, i, pr->nr-1);
207 /* should free memory here!! */
208 pr->nr--;
211 else
213 i++;
216 } /* if IF_ATYPE */
217 } /* for ftype_A */
218 } /* if IF_BTYPE */
219 } /* for ftype */
220 } /* if shake angles */
222 /* Add all the bonds with hydrogens to the shake list
223 * and remove them from the bond list
225 for (ftype = 0; (ftype < F_NRE); ftype++)
227 if (interaction_function[ftype].flags & IF_BTYPE)
229 pr = &(plist[ftype]);
230 j = 0;
231 for (i = 0; i < pr->nr; i++)
233 if ( (nshake != eshHBONDS) ||
234 (count_hydrogens (info, 2, pr->param[i].a) > 0) )
236 /* append this bond to the shake list */
237 p.ai() = pr->param[i].ai();
238 p.aj() = pr->param[i].aj();
239 p.c0() = pr->param[i].c0();
240 p.c1() = pr->param[i].c2();
241 add_param_to_list (&(plist[F_CONSTR]), &p);
243 else
245 copy_bond(pr, j++, i);
248 pr->nr = j;