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