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