Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / add_par.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
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.
14
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.
19  *
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.
26  *
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.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <string.h>
41 #include "typedefs.h"
42 #include "smalloc.h"
43 #include "grompp.h"
44 #include "macros.h"
45 #include "toputil.h"
46 #include "hackblock.h"
47 #include "string2.h"
48 #include "gmx_fatal.h"
49
50 static void clear_atom_list(int i0, atom_id a[])
51 {
52     int i;
53
54     for (i = i0; i < MAXATOMLIST; i++)
55     {
56         a[i] = -1;
57     }
58 }
59
60 static void clear_force_param(int i0, real c[])
61 {
62     int i;
63
64     for (i = i0; i < MAXFORCEPARAM; i++)
65     {
66         c[i] = NOTSET;
67     }
68 }
69
70 void add_param(t_params *ps, int ai, int aj, real *c, char *s)
71 {
72     int i;
73
74     if ((ai < 0) || (aj < 0))
75     {
76         gmx_fatal(FARGS, "Trying to add impossible atoms: ai=%d, aj=%d", ai, aj);
77     }
78     pr_alloc(1, ps);
79     ps->param[ps->nr].AI = ai;
80     ps->param[ps->nr].AJ = aj;
81     clear_atom_list(2, ps->param[ps->nr].a);
82     if (c == NULL)
83     {
84         clear_force_param(0, ps->param[ps->nr].c);
85     }
86     else
87     {
88         for (i = 0; (i < MAXFORCEPARAM); i++)
89         {
90             ps->param[ps->nr].c[i] = c[i];
91         }
92     }
93     set_p_string(&(ps->param[ps->nr]), s);
94     ps->nr++;
95 }
96
97 void add_imp_param(t_params *ps, int ai, int aj, int ak, int al, real c0, real c1,
98                    char *s)
99 {
100     pr_alloc(1, ps);
101     ps->param[ps->nr].AI = ai;
102     ps->param[ps->nr].AJ = aj;
103     ps->param[ps->nr].AK = ak;
104     ps->param[ps->nr].AL = al;
105     clear_atom_list  (4, ps->param[ps->nr].a);
106     ps->param[ps->nr].C0 = c0;
107     ps->param[ps->nr].C1 = c1;
108     clear_force_param(2, ps->param[ps->nr].c);
109     set_p_string(&(ps->param[ps->nr]), s);
110     ps->nr++;
111 }
112
113 void add_dih_param(t_params *ps, int ai, int aj, int ak, int al, real c0, real c1,
114                    real c2, char *s)
115 {
116     pr_alloc(1, ps);
117     ps->param[ps->nr].AI = ai;
118     ps->param[ps->nr].AJ = aj;
119     ps->param[ps->nr].AK = ak;
120     ps->param[ps->nr].AL = al;
121     clear_atom_list  (4, ps->param[ps->nr].a);
122     ps->param[ps->nr].C0 = c0;
123     ps->param[ps->nr].C1 = c1;
124     ps->param[ps->nr].C2 = c2;
125     clear_force_param(3, ps->param[ps->nr].c);
126     set_p_string(&(ps->param[ps->nr]), s);
127     ps->nr++;
128 }
129
130 void add_cmap_param(t_params *ps, int ai, int aj, int ak, int al, int am, char *s)
131 {
132     pr_alloc(1, ps);
133     ps->param[ps->nr].AI = ai;
134     ps->param[ps->nr].AJ = aj;
135     ps->param[ps->nr].AK = ak;
136     ps->param[ps->nr].AL = al;
137     ps->param[ps->nr].AM = am;
138     clear_atom_list(5, ps->param[ps->nr].a);
139     clear_force_param(0, ps->param[ps->nr].c);
140     set_p_string(&(ps->param[ps->nr]), s);
141     ps->nr++;
142 }
143
144 void add_vsite2_atoms(t_params *ps, int ai, int aj, int ak)
145 {
146     pr_alloc(1, ps);
147     ps->param[ps->nr].AI = ai;
148     ps->param[ps->nr].AJ = aj;
149     ps->param[ps->nr].AK = ak;
150     clear_atom_list  (3, ps->param[ps->nr].a);
151     clear_force_param(0, ps->param[ps->nr].c);
152     set_p_string(&(ps->param[ps->nr]), "");
153     ps->nr++;
154 }
155
156 void add_vsite2_param(t_params *ps, int ai, int aj, int ak, real c0)
157 {
158     pr_alloc(1, ps);
159     ps->param[ps->nr].AI = ai;
160     ps->param[ps->nr].AJ = aj;
161     ps->param[ps->nr].AK = ak;
162     clear_atom_list  (3, ps->param[ps->nr].a);
163     ps->param[ps->nr].C0 = c0;
164     clear_force_param(1, ps->param[ps->nr].c);
165     set_p_string(&(ps->param[ps->nr]), "");
166     ps->nr++;
167 }
168
169 void add_vsite3_param(t_params *ps, int ai, int aj, int ak, int al,
170                       real c0, real c1)
171 {
172     pr_alloc(1, ps);
173     ps->param[ps->nr].AI = ai;
174     ps->param[ps->nr].AJ = aj;
175     ps->param[ps->nr].AK = ak;
176     ps->param[ps->nr].AL = al;
177     clear_atom_list  (4, ps->param[ps->nr].a);
178     ps->param[ps->nr].C0 = c0;
179     ps->param[ps->nr].C1 = c1;
180     clear_force_param(2, ps->param[ps->nr].c);
181     set_p_string(&(ps->param[ps->nr]), "");
182     ps->nr++;
183 }
184
185 void add_vsite3_atoms(t_params *ps, int ai, int aj, int ak, int al, gmx_bool bSwapParity)
186 {
187     pr_alloc(1, ps);
188     ps->param[ps->nr].AI = ai;
189     ps->param[ps->nr].AJ = aj;
190     ps->param[ps->nr].AK = ak;
191     ps->param[ps->nr].AL = al;
192     clear_atom_list  (4, ps->param[ps->nr].a);
193     clear_force_param(0, ps->param[ps->nr].c);
194     if (bSwapParity)
195     {
196         ps->param[ps->nr].C1 = -1;
197     }
198     set_p_string(&(ps->param[ps->nr]), "");
199     ps->nr++;
200 }
201
202 void add_vsite4_atoms(t_params *ps, int ai, int aj, int ak, int al, int am)
203 {
204     pr_alloc(1, ps);
205     ps->param[ps->nr].AI = ai;
206     ps->param[ps->nr].AJ = aj;
207     ps->param[ps->nr].AK = ak;
208     ps->param[ps->nr].AL = al;
209     ps->param[ps->nr].AM = am;
210     clear_atom_list  (5, ps->param[ps->nr].a);
211     clear_force_param(0, ps->param[ps->nr].c);
212     set_p_string(&(ps->param[ps->nr]), "");
213     ps->nr++;
214 }
215
216 int search_jtype(t_restp *rtp, char *name, gmx_bool bNterm)
217 {
218     int   niter, iter, j, k, kmax, jmax, minstrlen;
219     char *rtpname, searchname[12];
220
221     strcpy(searchname, name);
222
223     /* Do a best match comparison */
224     /* for protein N-terminus, allow renaming of H1, H2 and H3 to H */
225     if (bNterm && (strlen(searchname) == 2) && (searchname[0] == 'H') &&
226         ( (searchname[1] == '1') || (searchname[1] == '2') ||
227           (searchname[1] == '3') ) )
228     {
229         niter = 2;
230     }
231     else
232     {
233         niter = 1;
234     }
235     kmax = 0;
236     jmax = -1;
237     for (iter = 0; (iter < niter && jmax == -1); iter++)
238     {
239         if (iter == 1)
240         {
241             /* Try without the hydrogen number in the N-terminus */
242             searchname[1] = '\0';
243         }
244         for (j = 0; (j < rtp->natom); j++)
245         {
246             rtpname = *(rtp->atomname[j]);
247             if (gmx_strcasecmp(searchname, rtpname) == 0)
248             {
249                 jmax = j;
250                 kmax = strlen(searchname);
251                 break;
252             }
253             if (iter == niter - 1)
254             {
255                 minstrlen = min(strlen(searchname), strlen(rtpname));
256                 for (k = 0; k < minstrlen; k++)
257                 {
258                     if (searchname[k] != rtpname[k])
259                     {
260                         break;
261                     }
262                 }
263                 if (k > kmax)
264                 {
265                     kmax = k;
266                     jmax = j;
267                 }
268             }
269         }
270     }
271     if (jmax == -1)
272     {
273         gmx_fatal(FARGS, "Atom %s not found in rtp database in residue %s",
274                   searchname, rtp->resname);
275     }
276     if (kmax != strlen(searchname))
277     {
278         gmx_fatal(FARGS, "Atom %s not found in rtp database in residue %s, "
279                   "it looks a bit like %s",
280                   searchname, rtp->resname, *(rtp->atomname[jmax]));
281     }
282     return jmax;
283 }