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