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) 2011,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! */
42 #include "gromacs/gmxpreprocess/grompp-impl.h"
43 #include "gromacs/gmxpreprocess/hackblock.h"
44 #include "gromacs/gmxpreprocess/toputil.h"
45 #include "gromacs/legacyheaders/macros.h"
46 #include "gromacs/legacyheaders/typedefs.h"
47 #include "gromacs/utility/cstringutil.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/smalloc.h"
51 static void clear_atom_list(int i0, atom_id a[])
55 for (i = i0; i < MAXATOMLIST; i++)
61 static void clear_force_param(int i0, real c[])
65 for (i = i0; i < MAXFORCEPARAM; i++)
71 void add_param(t_params *ps, int ai, int aj, real *c, char *s)
75 if ((ai < 0) || (aj < 0))
77 gmx_fatal(FARGS, "Trying to add impossible atoms: ai=%d, aj=%d", ai, aj);
80 ps->param[ps->nr].AI = ai;
81 ps->param[ps->nr].AJ = aj;
82 clear_atom_list(2, ps->param[ps->nr].a);
85 clear_force_param(0, ps->param[ps->nr].c);
89 for (i = 0; (i < MAXFORCEPARAM); i++)
91 ps->param[ps->nr].c[i] = c[i];
94 set_p_string(&(ps->param[ps->nr]), s);
98 void add_imp_param(t_params *ps, int ai, int aj, int ak, int al, real c0, real c1,
102 ps->param[ps->nr].AI = ai;
103 ps->param[ps->nr].AJ = aj;
104 ps->param[ps->nr].AK = ak;
105 ps->param[ps->nr].AL = al;
106 clear_atom_list (4, ps->param[ps->nr].a);
107 ps->param[ps->nr].C0 = c0;
108 ps->param[ps->nr].C1 = c1;
109 clear_force_param(2, ps->param[ps->nr].c);
110 set_p_string(&(ps->param[ps->nr]), s);
114 void add_dih_param(t_params *ps, int ai, int aj, int ak, int al, real c0, real c1,
118 ps->param[ps->nr].AI = ai;
119 ps->param[ps->nr].AJ = aj;
120 ps->param[ps->nr].AK = ak;
121 ps->param[ps->nr].AL = al;
122 clear_atom_list (4, ps->param[ps->nr].a);
123 ps->param[ps->nr].C0 = c0;
124 ps->param[ps->nr].C1 = c1;
125 ps->param[ps->nr].C2 = c2;
126 clear_force_param(3, ps->param[ps->nr].c);
127 set_p_string(&(ps->param[ps->nr]), s);
131 void add_cmap_param(t_params *ps, int ai, int aj, int ak, int al, int am, char *s)
134 ps->param[ps->nr].AI = ai;
135 ps->param[ps->nr].AJ = aj;
136 ps->param[ps->nr].AK = ak;
137 ps->param[ps->nr].AL = al;
138 ps->param[ps->nr].AM = am;
139 clear_atom_list(5, ps->param[ps->nr].a);
140 clear_force_param(0, ps->param[ps->nr].c);
141 set_p_string(&(ps->param[ps->nr]), s);
145 void add_vsite2_atoms(t_params *ps, int ai, int aj, int ak)
148 ps->param[ps->nr].AI = ai;
149 ps->param[ps->nr].AJ = aj;
150 ps->param[ps->nr].AK = ak;
151 clear_atom_list (3, ps->param[ps->nr].a);
152 clear_force_param(0, ps->param[ps->nr].c);
153 set_p_string(&(ps->param[ps->nr]), "");
157 void add_vsite2_param(t_params *ps, int ai, int aj, int ak, real c0)
160 ps->param[ps->nr].AI = ai;
161 ps->param[ps->nr].AJ = aj;
162 ps->param[ps->nr].AK = ak;
163 clear_atom_list (3, ps->param[ps->nr].a);
164 ps->param[ps->nr].C0 = c0;
165 clear_force_param(1, ps->param[ps->nr].c);
166 set_p_string(&(ps->param[ps->nr]), "");
170 void add_vsite3_param(t_params *ps, int ai, int aj, int ak, int al,
174 ps->param[ps->nr].AI = ai;
175 ps->param[ps->nr].AJ = aj;
176 ps->param[ps->nr].AK = ak;
177 ps->param[ps->nr].AL = al;
178 clear_atom_list (4, ps->param[ps->nr].a);
179 ps->param[ps->nr].C0 = c0;
180 ps->param[ps->nr].C1 = c1;
181 clear_force_param(2, ps->param[ps->nr].c);
182 set_p_string(&(ps->param[ps->nr]), "");
186 void add_vsite3_atoms(t_params *ps, int ai, int aj, int ak, int al, gmx_bool bSwapParity)
189 ps->param[ps->nr].AI = ai;
190 ps->param[ps->nr].AJ = aj;
191 ps->param[ps->nr].AK = ak;
192 ps->param[ps->nr].AL = al;
193 clear_atom_list (4, ps->param[ps->nr].a);
194 clear_force_param(0, ps->param[ps->nr].c);
197 ps->param[ps->nr].C1 = -1;
199 set_p_string(&(ps->param[ps->nr]), "");
203 void add_vsite4_atoms(t_params *ps, int ai, int aj, int ak, int al, int am)
206 ps->param[ps->nr].AI = ai;
207 ps->param[ps->nr].AJ = aj;
208 ps->param[ps->nr].AK = ak;
209 ps->param[ps->nr].AL = al;
210 ps->param[ps->nr].AM = am;
211 clear_atom_list (5, ps->param[ps->nr].a);
212 clear_force_param(0, ps->param[ps->nr].c);
213 set_p_string(&(ps->param[ps->nr]), "");
217 int search_jtype(t_restp *rtp, char *name, gmx_bool bNterm)
219 int niter, iter, j, k, kmax, jmax, minstrlen;
220 char *rtpname, searchname[12];
222 strcpy(searchname, name);
224 /* Do a best match comparison */
225 /* for protein N-terminus, allow renaming of H1, H2 and H3 to H */
226 if (bNterm && (strlen(searchname) == 2) && (searchname[0] == 'H') &&
227 ( (searchname[1] == '1') || (searchname[1] == '2') ||
228 (searchname[1] == '3') ) )
238 for (iter = 0; (iter < niter && jmax == -1); iter++)
242 /* Try without the hydrogen number in the N-terminus */
243 searchname[1] = '\0';
245 for (j = 0; (j < rtp->natom); j++)
247 rtpname = *(rtp->atomname[j]);
248 if (gmx_strcasecmp(searchname, rtpname) == 0)
251 kmax = strlen(searchname);
254 if (iter == niter - 1)
256 minstrlen = min(strlen(searchname), strlen(rtpname));
257 for (k = 0; k < minstrlen; k++)
259 if (searchname[k] != rtpname[k])
274 gmx_fatal(FARGS, "Atom %s not found in rtp database in residue %s",
275 searchname, rtp->resname);
277 if (kmax != strlen(searchname))
279 gmx_fatal(FARGS, "Atom %s not found in rtp database in residue %s, "
280 "it looks a bit like %s",
281 searchname, rtp->resname, *(rtp->atomname[jmax]));