Clean up string2.h-related 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 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <string.h>
43 #include "hackblock.h"
44 #include "smalloc.h"
45 #include "vec.h"
46 #include "macros.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     free(*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 }
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 merge_t_bondeds(t_rbondeds s[], t_rbondeds d[], gmx_bool bMin, gmx_bool bPlus)
200 {
201     int      i, j;
202     gmx_bool bBondsRemoved;
203
204     bBondsRemoved = FALSE;
205     for (i = 0; i < ebtsNR; i++)
206     {
207         if (s[i].nb > 0)
208         {
209             /* make space */
210             srenew(d[i].b, d[i].nb + s[i].nb);
211             for (j = 0; j < s[i].nb; j++)
212             {
213                 if (!(bMin && contains_char(&s[i].b[j], '-'))
214                     && !(bPlus && contains_char(&s[i].b[j], '+')))
215                 {
216                     copy_t_rbonded(&s[i].b[j], &d[i].b[ d[i].nb ]);
217                     d[i].nb++;
218                 }
219                 else if (i == ebtsBONDS)
220                 {
221                     bBondsRemoved = TRUE;
222                 }
223             }
224         }
225     }
226
227     return bBondsRemoved;
228 }
229
230 void copy_t_restp(t_restp *s, t_restp *d)
231 {
232     int i;
233
234     *d         = *s;
235     d->resname = safe_strdup(s->resname);
236     snew(d->atom, s->natom);
237     for (i = 0; i < s->natom; i++)
238     {
239         d->atom[i] = s->atom[i];
240     }
241     snew(d->atomname, s->natom);
242     for (i = 0; i < s->natom; i++)
243     {
244         snew(d->atomname[i], 1);
245         *d->atomname[i] = safe_strdup(*s->atomname[i]);
246     }
247     snew(d->cgnr, s->natom);
248     for (i = 0; i < s->natom; i++)
249     {
250         d->cgnr[i] = s->cgnr[i];
251     }
252     for (i = 0; i < ebtsNR; i++)
253     {
254         d->rb[i].type = s->rb[i].type;
255         d->rb[i].nb   = 0;
256         d->rb[i].b    = NULL;
257     }
258     merge_t_bondeds(s->rb, d->rb, FALSE, FALSE);
259 }
260
261 void copy_t_hack(t_hack *s, t_hack *d)
262 {
263     int i;
264
265     *d       = *s;
266     d->oname = safe_strdup(s->oname);
267     d->nname = safe_strdup(s->nname);
268     if (s->atom)
269     {
270         snew(d->atom, 1);
271         *(d->atom) = *(s->atom);
272     }
273     else
274     {
275         d->atom = NULL;
276     }
277     for (i = 0; i < 4; i++)
278     {
279         d->a[i] = safe_strdup(s->a[i]);
280     }
281     copy_rvec(s->newx, d->newx);
282 }
283
284 void merge_hacks_lo(int ns, t_hack *s, int *nd, t_hack **d)
285 {
286     int i;
287
288     if (ns)
289     {
290         srenew(*d, *nd + ns);
291         for (i = 0; i < ns; i++)
292         {
293             copy_t_hack(&s[i], &(*d)[*nd + i]);
294         }
295         (*nd) += ns;
296     }
297 }
298
299 void merge_hacks(t_hackblock *s, t_hackblock *d)
300 {
301     merge_hacks_lo(s->nhack, s->hack, &d->nhack, &d->hack);
302 }
303
304 void merge_t_hackblock(t_hackblock *s, t_hackblock *d)
305 {
306     merge_hacks(s, d);
307     merge_t_bondeds(s->rb, d->rb, FALSE, FALSE);
308 }
309
310 void copy_t_hackblock(t_hackblock *s, t_hackblock *d)
311 {
312     int i;
313
314     *d       = *s;
315     d->name  = safe_strdup(s->name);
316     d->nhack = 0;
317     d->hack  = NULL;
318     for (i = 0; i < ebtsNR; i++)
319     {
320         d->rb[i].nb = 0;
321         d->rb[i].b  = NULL;
322     }
323     merge_t_hackblock(s, d);
324 }
325
326 #undef safe_strdup
327
328 void dump_hb(FILE *out, int nres, t_hackblock hb[])
329 {
330     int i, j, k, l;
331
332 #define SS(s) (s) ? (s) : "-"
333 #define SA(s) (s) ? "+" : ""
334     fprintf(out, "t_hackblock\n");
335     for (i = 0; i < nres; i++)
336     {
337         fprintf(out, "%3d %4s %2d %2d\n",
338                 i, SS(hb[i].name), hb[i].nhack, hb[i].maxhack);
339         if (hb[i].nhack)
340         {
341             for (j = 0; j < hb[i].nhack; j++)
342             {
343                 fprintf(out, "%d: %d %4s %4s %1s %2d %d %4s %4s %4s %4s\n",
344                         j, hb[i].hack[j].nr,
345                         SS(hb[i].hack[j].oname), SS(hb[i].hack[j].nname),
346                         SA(hb[i].hack[j].atom), hb[i].hack[j].tp, hb[i].hack[j].cgnr,
347                         SS(hb[i].hack[j].AI), SS(hb[i].hack[j].AJ),
348                         SS(hb[i].hack[j].AK), SS(hb[i].hack[j].AL) );
349             }
350         }
351         for (j = 0; j < ebtsNR; j++)
352         {
353             if (hb[i].rb[j].nb)
354             {
355                 fprintf(out, " %c %d:", btsNames[j][0], hb[i].rb[j].nb);
356                 for (k = 0; k < hb[i].rb[j].nb; k++)
357                 {
358                     fprintf(out, " [");
359                     for (l = 0; l < btsNiatoms[j]; l++)
360                     {
361                         fprintf(out, " %s", hb[i].rb[j].b[k].a[l]);
362                     }
363                     fprintf(out, " %s]", SS(hb[i].rb[j].b[k].s));
364                 }
365                 fprintf(out, "\n");
366             }
367         }
368         fprintf(out, "\n");
369     }
370 #undef SS
371 #undef SA
372 }
373
374 void init_t_protonate(t_protonate *protonate)
375 {
376     protonate->bInit = FALSE;
377 }