Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / hackblock.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 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <string.h>
43 #include "hackblock.h"
44 #include "gromacs/utility/smalloc.h"
45 #include "gromacs/math/vec.h"
46 #include "names.h"
47
48 /* these MUST correspond to the enum in hackblock.h */
49 const char *btsNames[ebtsNR] = { "bonds", "angles", "dihedrals", "impropers", "exclusions", "cmap" };
50 const int   btsNiatoms[ebtsNR] = { 2,       3,        4,           4,           2,             5 };
51
52 static void free_t_bonded(t_rbonded *rb)
53 {
54     int i;
55
56     for (i = 0; i < MAXATOMLIST; i++)
57     {
58         sfree(rb->a[i]);
59     }
60     sfree(rb->s);
61 }
62
63 static void free_t_bondeds(t_rbondeds *rbs)
64 {
65     int i;
66
67     for (i = 0; i < rbs->nb; i++)
68     {
69         free_t_bonded(&rbs->b[i]);
70     }
71     sfree(rbs->b);
72     rbs->b  = NULL;
73     rbs->nb = 0;
74 }
75
76 void free_t_restp(int nrtp, t_restp **rtp)
77 {
78     int i, j;
79
80     for (i = 0; i < nrtp; i++)
81     {
82         sfree((*rtp)[i].resname);
83         sfree((*rtp)[i].atom);
84         for (j = 0; j < (*rtp)[i].natom; j++)
85         {
86             sfree(*(*rtp)[i].atomname[j]);
87             sfree((*rtp)[i].atomname[j]);
88         }
89         sfree((*rtp)[i].atomname);
90         sfree((*rtp)[i].cgnr);
91         for (j = 0; j < ebtsNR; j++)
92         {
93             free_t_bondeds(&(*rtp)[i].rb[j]);
94         }
95     }
96     sfree(*rtp);
97 }
98
99 void free_t_hack(int nh, t_hack **h)
100 {
101     int i, j;
102
103     for (i = 0; i < nh; i++)
104     {
105         sfree((*h)[i].oname);
106         sfree((*h)[i].nname);
107         sfree((*h)[i].atom);
108         for (j = 0; j < 4; j++)
109         {
110             sfree((*h)[i].a[j]);
111         }
112     }
113     sfree(*h);
114     *h = NULL;
115 }
116
117 void free_t_hackblock(int nhb, t_hackblock **hb)
118 {
119     int i, j;
120
121     for (i = 0; i < nhb; i++)
122     {
123         sfree((*hb)[i].name);
124         free_t_hack((*hb)[i].nhack, &(*hb)[i].hack);
125         for (j = 0; j < ebtsNR; j++)
126         {
127             free_t_bondeds(&(*hb)[i].rb[j]);
128         }
129     }
130     sfree(*hb);
131 }
132
133 void clear_t_hackblock(t_hackblock *hb)
134 {
135     int i;
136
137     hb->name    = NULL;
138     hb->nhack   = 0;
139     hb->maxhack = 0;
140     hb->hack    = NULL;
141     for (i = 0; i < ebtsNR; i++)
142     {
143         hb->rb[i].nb = 0;
144         hb->rb[i].b  = NULL;
145     }
146 }
147
148 void clear_t_hack(t_hack *hack)
149 {
150     int i;
151
152     hack->nr    = 0;
153     hack->oname = NULL;
154     hack->nname = NULL;
155     hack->atom  = NULL;
156     hack->cgnr  = NOTSET;
157     hack->tp    = 0;
158     hack->nctl  = 0;
159     for (i = 0; i < 4; i++)
160     {
161         hack->a[i]  = NULL;
162     }
163     for (i = 0; i < DIM; i++)
164     {
165         hack->newx[i] = NOTSET;
166     }
167 }
168
169 #define safe_strdup(str) ((str != NULL) ? strdup(str) : NULL)
170
171 static void copy_t_rbonded(t_rbonded *s, t_rbonded *d)
172 {
173     int i;
174
175     for (i = 0; i < MAXATOMLIST; i++)
176     {
177         d->a[i] = safe_strdup(s->a[i]);
178     }
179     d->s     = safe_strdup(s->s);
180     d->match = s->match;
181 }
182
183 static gmx_bool contains_char(t_rbonded *s, char c)
184 {
185     int      i;
186     gmx_bool bRet;
187
188     bRet = FALSE;
189     for (i = 0; i < MAXATOMLIST; i++)
190     {
191         if (s->a[i] && s->a[i][0] == c)
192         {
193             bRet = TRUE;
194         }
195     }
196
197     return bRet;
198 }
199
200 gmx_bool rbonded_atoms_exist_in_list(t_rbonded *b, t_rbonded blist[], int nlist, int natoms)
201 {
202     int      i, k;
203     gmx_bool matchFound = FALSE;
204     gmx_bool atomsMatch;
205
206     for (i = 0; i < nlist && !matchFound; i++)
207     {
208         atomsMatch = TRUE;
209         for (k = 0; k < natoms && atomsMatch; k++)
210         {
211             atomsMatch = atomsMatch && !strcmp(b->a[k], blist[i].a[k]);
212         }
213         /* Try reverse if forward match did not work */
214         if (!atomsMatch)
215         {
216             atomsMatch = TRUE;
217             for (k = 0; k < natoms && atomsMatch; k++)
218             {
219                 atomsMatch = atomsMatch && !strcmp(b->a[k], blist[i].a[natoms-1-k]);
220             }
221         }
222         matchFound = atomsMatch;
223     }
224     return matchFound;
225 }
226
227 gmx_bool merge_t_bondeds(t_rbondeds s[], t_rbondeds d[], gmx_bool bMin, gmx_bool bPlus)
228 {
229     int      i, j;
230     gmx_bool bBondsRemoved;
231
232     bBondsRemoved = FALSE;
233     for (i = 0; i < ebtsNR; i++)
234     {
235         if (s[i].nb > 0)
236         {
237             /* make space */
238             srenew(d[i].b, d[i].nb + s[i].nb);
239             for (j = 0; j < s[i].nb; j++)
240             {
241                 /* Check if this bonded string already exists before adding.
242                  * We are merging from the main rtp to the hackblocks, so this
243                  * will mean the hackblocks overwrite the man rtp, as intended.
244                  */
245                 if (!rbonded_atoms_exist_in_list(&s[i].b[j], d[i].b, d[i].nb, btsNiatoms[i]))
246                 {
247                     if (!(bMin && contains_char(&s[i].b[j], '-'))
248                         && !(bPlus && contains_char(&s[i].b[j], '+')))
249                     {
250                         copy_t_rbonded(&s[i].b[j], &d[i].b[ d[i].nb ]);
251                         d[i].nb++;
252                     }
253                     else if (i == ebtsBONDS)
254                     {
255                         bBondsRemoved = TRUE;
256                     }
257                 }
258             }
259         }
260     }
261     return bBondsRemoved;
262 }
263
264 void copy_t_restp(t_restp *s, t_restp *d)
265 {
266     int i;
267
268     *d         = *s;
269     d->resname = safe_strdup(s->resname);
270     snew(d->atom, s->natom);
271     for (i = 0; i < s->natom; i++)
272     {
273         d->atom[i] = s->atom[i];
274     }
275     snew(d->atomname, s->natom);
276     for (i = 0; i < s->natom; i++)
277     {
278         snew(d->atomname[i], 1);
279         *d->atomname[i] = safe_strdup(*s->atomname[i]);
280     }
281     snew(d->cgnr, s->natom);
282     for (i = 0; i < s->natom; i++)
283     {
284         d->cgnr[i] = s->cgnr[i];
285     }
286     for (i = 0; i < ebtsNR; i++)
287     {
288         d->rb[i].type = s->rb[i].type;
289         d->rb[i].nb   = 0;
290         d->rb[i].b    = NULL;
291     }
292     merge_t_bondeds(s->rb, d->rb, FALSE, FALSE);
293 }
294
295 void copy_t_hack(t_hack *s, t_hack *d)
296 {
297     int i;
298
299     *d       = *s;
300     d->oname = safe_strdup(s->oname);
301     d->nname = safe_strdup(s->nname);
302     if (s->atom)
303     {
304         snew(d->atom, 1);
305         *(d->atom) = *(s->atom);
306     }
307     else
308     {
309         d->atom = NULL;
310     }
311     for (i = 0; i < 4; i++)
312     {
313         d->a[i] = safe_strdup(s->a[i]);
314     }
315     copy_rvec(s->newx, d->newx);
316 }
317
318 void merge_hacks_lo(int ns, t_hack *s, int *nd, t_hack **d)
319 {
320     int i;
321
322     if (ns)
323     {
324         srenew(*d, *nd + ns);
325         for (i = 0; i < ns; i++)
326         {
327             copy_t_hack(&s[i], &(*d)[*nd + i]);
328         }
329         (*nd) += ns;
330     }
331 }
332
333 void merge_hacks(t_hackblock *s, t_hackblock *d)
334 {
335     merge_hacks_lo(s->nhack, s->hack, &d->nhack, &d->hack);
336 }
337
338 void merge_t_hackblock(t_hackblock *s, t_hackblock *d)
339 {
340     merge_hacks(s, d);
341     merge_t_bondeds(s->rb, d->rb, FALSE, FALSE);
342 }
343
344 void copy_t_hackblock(t_hackblock *s, t_hackblock *d)
345 {
346     int i;
347
348     *d       = *s;
349     d->name  = safe_strdup(s->name);
350     d->nhack = 0;
351     d->hack  = NULL;
352     for (i = 0; i < ebtsNR; i++)
353     {
354         d->rb[i].nb = 0;
355         d->rb[i].b  = NULL;
356     }
357     merge_t_hackblock(s, d);
358 }
359
360 #undef safe_strdup
361
362 void dump_hb(FILE *out, int nres, t_hackblock hb[])
363 {
364     int i, j, k, l;
365
366 #define SS(s) (s) ? (s) : "-"
367 #define SA(s) (s) ? "+" : ""
368     fprintf(out, "t_hackblock\n");
369     for (i = 0; i < nres; i++)
370     {
371         fprintf(out, "%3d %4s %2d %2d\n",
372                 i, SS(hb[i].name), hb[i].nhack, hb[i].maxhack);
373         if (hb[i].nhack)
374         {
375             for (j = 0; j < hb[i].nhack; j++)
376             {
377                 fprintf(out, "%d: %d %4s %4s %1s %2d %d %4s %4s %4s %4s\n",
378                         j, hb[i].hack[j].nr,
379                         SS(hb[i].hack[j].oname), SS(hb[i].hack[j].nname),
380                         SA(hb[i].hack[j].atom), hb[i].hack[j].tp, hb[i].hack[j].cgnr,
381                         SS(hb[i].hack[j].AI), SS(hb[i].hack[j].AJ),
382                         SS(hb[i].hack[j].AK), SS(hb[i].hack[j].AL) );
383             }
384         }
385         for (j = 0; j < ebtsNR; j++)
386         {
387             if (hb[i].rb[j].nb)
388             {
389                 fprintf(out, " %c %d:", btsNames[j][0], hb[i].rb[j].nb);
390                 for (k = 0; k < hb[i].rb[j].nb; k++)
391                 {
392                     fprintf(out, " [");
393                     for (l = 0; l < btsNiatoms[j]; l++)
394                     {
395                         fprintf(out, " %s", hb[i].rb[j].b[k].a[l]);
396                     }
397                     fprintf(out, " %s]", SS(hb[i].rb[j].b[k].s));
398                 }
399                 fprintf(out, " Entry matched: %s\n", yesno_names[hb[i].rb[j].b[k].match]);
400             }
401         }
402         fprintf(out, "\n");
403     }
404 #undef SS
405 #undef SA
406 }
407
408 void init_t_protonate(t_protonate *protonate)
409 {
410     protonate->bInit = FALSE;
411 }