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