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,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! */
44 #include "gromacs/gmxpreprocess/grompp-impl.h"
45 #include "gromacs/gmxpreprocess/hackblock.h"
46 #include "gromacs/gmxpreprocess/toputil.h"
47 #include "gromacs/legacyheaders/macros.h"
48 #include "gromacs/legacyheaders/typedefs.h"
49 #include "gromacs/utility/cstringutil.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/smalloc.h"
53 static void clear_atom_list(int i0, atom_id a[])
57 for (i = i0; i < MAXATOMLIST; i++)
63 static void clear_force_param(int i0, real c[])
67 for (i = i0; i < MAXFORCEPARAM; i++)
73 void add_param(t_params *ps, int ai, int aj, real *c, char *s)
77 if ((ai < 0) || (aj < 0))
79 gmx_fatal(FARGS, "Trying to add impossible atoms: ai=%d, aj=%d", ai, aj);
82 ps->param[ps->nr].AI = ai;
83 ps->param[ps->nr].AJ = aj;
84 clear_atom_list(2, ps->param[ps->nr].a);
87 clear_force_param(0, ps->param[ps->nr].c);
91 for (i = 0; (i < MAXFORCEPARAM); i++)
93 ps->param[ps->nr].c[i] = c[i];
96 set_p_string(&(ps->param[ps->nr]), s);
100 void add_imp_param(t_params *ps, int ai, int aj, int ak, int al, real c0, real c1,
104 ps->param[ps->nr].AI = ai;
105 ps->param[ps->nr].AJ = aj;
106 ps->param[ps->nr].AK = ak;
107 ps->param[ps->nr].AL = al;
108 clear_atom_list (4, ps->param[ps->nr].a);
109 ps->param[ps->nr].C0 = c0;
110 ps->param[ps->nr].C1 = c1;
111 clear_force_param(2, ps->param[ps->nr].c);
112 set_p_string(&(ps->param[ps->nr]), s);
116 void add_dih_param(t_params *ps, int ai, int aj, int ak, int al, real c0, real c1,
120 ps->param[ps->nr].AI = ai;
121 ps->param[ps->nr].AJ = aj;
122 ps->param[ps->nr].AK = ak;
123 ps->param[ps->nr].AL = al;
124 clear_atom_list (4, ps->param[ps->nr].a);
125 ps->param[ps->nr].C0 = c0;
126 ps->param[ps->nr].C1 = c1;
127 ps->param[ps->nr].C2 = c2;
128 clear_force_param(3, ps->param[ps->nr].c);
129 set_p_string(&(ps->param[ps->nr]), s);
133 void add_cmap_param(t_params *ps, int ai, int aj, int ak, int al, int am, char *s)
136 ps->param[ps->nr].AI = ai;
137 ps->param[ps->nr].AJ = aj;
138 ps->param[ps->nr].AK = ak;
139 ps->param[ps->nr].AL = al;
140 ps->param[ps->nr].AM = am;
141 clear_atom_list(5, ps->param[ps->nr].a);
142 clear_force_param(0, ps->param[ps->nr].c);
143 set_p_string(&(ps->param[ps->nr]), s);
147 void add_vsite2_atoms(t_params *ps, int ai, int aj, int ak)
150 ps->param[ps->nr].AI = ai;
151 ps->param[ps->nr].AJ = aj;
152 ps->param[ps->nr].AK = ak;
153 clear_atom_list (3, ps->param[ps->nr].a);
154 clear_force_param(0, ps->param[ps->nr].c);
155 set_p_string(&(ps->param[ps->nr]), "");
159 void add_vsite2_param(t_params *ps, int ai, int aj, int ak, real c0)
162 ps->param[ps->nr].AI = ai;
163 ps->param[ps->nr].AJ = aj;
164 ps->param[ps->nr].AK = ak;
165 clear_atom_list (3, ps->param[ps->nr].a);
166 ps->param[ps->nr].C0 = c0;
167 clear_force_param(1, ps->param[ps->nr].c);
168 set_p_string(&(ps->param[ps->nr]), "");
172 void add_vsite3_param(t_params *ps, int ai, int aj, int ak, int al,
176 ps->param[ps->nr].AI = ai;
177 ps->param[ps->nr].AJ = aj;
178 ps->param[ps->nr].AK = ak;
179 ps->param[ps->nr].AL = al;
180 clear_atom_list (4, ps->param[ps->nr].a);
181 ps->param[ps->nr].C0 = c0;
182 ps->param[ps->nr].C1 = c1;
183 clear_force_param(2, ps->param[ps->nr].c);
184 set_p_string(&(ps->param[ps->nr]), "");
188 void add_vsite3_atoms(t_params *ps, int ai, int aj, int ak, int al, gmx_bool bSwapParity)
191 ps->param[ps->nr].AI = ai;
192 ps->param[ps->nr].AJ = aj;
193 ps->param[ps->nr].AK = ak;
194 ps->param[ps->nr].AL = al;
195 clear_atom_list (4, ps->param[ps->nr].a);
196 clear_force_param(0, ps->param[ps->nr].c);
199 ps->param[ps->nr].C1 = -1;
201 set_p_string(&(ps->param[ps->nr]), "");
205 void add_vsite4_atoms(t_params *ps, int ai, int aj, int ak, int al, int am)
208 ps->param[ps->nr].AI = ai;
209 ps->param[ps->nr].AJ = aj;
210 ps->param[ps->nr].AK = ak;
211 ps->param[ps->nr].AL = al;
212 ps->param[ps->nr].AM = am;
213 clear_atom_list (5, ps->param[ps->nr].a);
214 clear_force_param(0, ps->param[ps->nr].c);
215 set_p_string(&(ps->param[ps->nr]), "");
219 int search_jtype(t_restp *rtp, char *name, gmx_bool bNterm)
221 int niter, iter, j, jmax;
222 size_t k, kmax, minstrlen;
223 char *rtpname, searchname[12];
225 strcpy(searchname, name);
227 /* Do a best match comparison */
228 /* for protein N-terminus, allow renaming of H1, H2 and H3 to H */
229 if (bNterm && (strlen(searchname) == 2) && (searchname[0] == 'H') &&
230 ( (searchname[1] == '1') || (searchname[1] == '2') ||
231 (searchname[1] == '3') ) )
241 for (iter = 0; (iter < niter && jmax == -1); iter++)
245 /* Try without the hydrogen number in the N-terminus */
246 searchname[1] = '\0';
248 for (j = 0; (j < rtp->natom); j++)
250 rtpname = *(rtp->atomname[j]);
251 if (gmx_strcasecmp(searchname, rtpname) == 0)
254 kmax = strlen(searchname);
257 if (iter == niter - 1)
259 minstrlen = std::min(strlen(searchname), strlen(rtpname));
260 for (k = 0; k < minstrlen; k++)
262 if (searchname[k] != rtpname[k])
277 gmx_fatal(FARGS, "Atom %s not found in rtp database in residue %s",
278 searchname, rtp->resname);
280 if (kmax != strlen(searchname))
282 gmx_fatal(FARGS, "Atom %s not found in rtp database in residue %s, "
283 "it looks a bit like %s",
284 searchname, rtp->resname, *(rtp->atomname[jmax]));