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