Update copyright statements and change license to LGPL
[alexxy/gromacs.git] / src / kernel / 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  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 /* This file is completely threadsafe - keep it that way! */
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
42
43 #include <string.h>
44 #include "hackblock.h"
45 #include "smalloc.h"
46 #include "vec.h"
47 #include "string2.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     sfree(rb->a[i]);
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     free_t_bonded(&rbs->b[i]);
68   sfree(rbs->b);
69   rbs->b=NULL;
70   rbs->nb=0;
71 }
72
73 void free_t_restp(int nrtp, t_restp **rtp)
74 {
75   int i,j;
76   
77   for(i=0; i<nrtp; i++) {
78     sfree((*rtp)[i].resname);
79     sfree((*rtp)[i].atom);
80     for(j=0; j<(*rtp)[i].natom; j++) {
81       sfree(*(*rtp)[i].atomname[j]);
82       sfree((*rtp)[i].atomname[j]);
83     }
84     sfree((*rtp)[i].atomname);
85     sfree((*rtp)[i].cgnr);
86     for(j=0; j<ebtsNR; j++)
87       free_t_bondeds(&(*rtp)[i].rb[j]);
88   }
89   free(*rtp);
90 }
91
92 void free_t_hack(int nh, t_hack **h)
93 {
94   int i, j;
95   
96   for(i=0; i<nh; i++) {
97     sfree((*h)[i].oname);
98     sfree((*h)[i].nname);
99     sfree((*h)[i].atom);
100     for(j=0; j<4; j++)
101       sfree((*h)[i].a[j]);
102   }
103   sfree(*h);
104   *h=NULL;
105 }
106
107 void free_t_hackblock(int nhb, t_hackblock **hb)
108 {
109   int i, j;
110   
111   for(i=0; i<nhb; i++) {
112     sfree((*hb)[i].name);
113     free_t_hack((*hb)[i].nhack, &(*hb)[i].hack);
114     for(j=0; j<ebtsNR; j++)
115       free_t_bondeds(&(*hb)[i].rb[j]);
116   }
117   sfree(*hb);
118 }
119
120 void clear_t_hackblock(t_hackblock *hb)
121 {
122   int i;
123   
124   hb->name   = NULL;
125   hb->nhack  = 0;
126   hb->maxhack= 0;
127   hb->hack   = NULL;
128   for(i=0; i<ebtsNR; i++) {
129     hb->rb[i].nb=0;
130     hb->rb[i].b=NULL;
131   }
132 }
133
134 void clear_t_hack(t_hack *hack)
135 {
136   int i;
137   
138   hack->nr    = 0;
139   hack->oname = NULL;
140   hack->nname = NULL;
141   hack->atom  = NULL;
142   hack->cgnr  = NOTSET;
143   hack->tp    = 0;
144   hack->nctl  = 0;
145   for(i=0; i<4; i++)
146     hack->a[i]  = NULL;
147   for(i=0; i<DIM; i++)
148     hack->newx[i] = NOTSET;
149 }
150
151 #define safe_strdup(str) ((str != NULL) ? strdup(str) : NULL)
152
153 static void copy_t_rbonded(t_rbonded *s, t_rbonded *d)
154 {
155   int i;
156   
157   for(i=0; i<MAXATOMLIST; i++)
158     d->a[i] = safe_strdup(s->a[i]);
159   d->s = safe_strdup(s->s);
160 }
161
162 static gmx_bool contains_char(t_rbonded *s,char c)
163 {
164   int i;
165   gmx_bool bRet;
166   
167   bRet = FALSE;
168   for(i=0; i<MAXATOMLIST; i++)
169     if (s->a[i] && s->a[i][0]==c)
170       bRet = TRUE;
171   
172   return bRet;
173 }
174
175 gmx_bool merge_t_bondeds(t_rbondeds s[], t_rbondeds d[],gmx_bool bMin,gmx_bool bPlus)
176 {
177   int i, j;
178   gmx_bool bBondsRemoved;
179   
180   bBondsRemoved = FALSE;
181   for(i=0; i < ebtsNR; i++) {
182     if ( s[i].nb > 0 ) {
183       /* make space */
184       srenew(d[i].b, d[i].nb + s[i].nb);
185       for(j=0; j < s[i].nb; j++)
186         if (!(bMin && contains_char(&s[i].b[j],'-'))
187             && !(bPlus && contains_char(&s[i].b[j],'+'))) {
188           copy_t_rbonded(&s[i].b[j], &d[i].b[ d[i].nb ]);
189           d[i].nb ++;
190         } else if (i == ebtsBONDS) {
191           bBondsRemoved = TRUE;
192         }
193     }
194   }
195
196   return bBondsRemoved;
197 }
198
199 void copy_t_restp(t_restp *s, t_restp *d)
200 {
201   int i;
202   
203   *d = *s;
204   d->resname = safe_strdup(s->resname);
205   snew(d->atom, s->natom);
206   for(i=0; i < s->natom; i++)
207     d->atom[i] = s->atom[i];
208   snew(d->atomname, s->natom);
209   for(i=0; i < s->natom; i++) {
210     snew(d->atomname[i],1);
211     *d->atomname[i] = safe_strdup(*s->atomname[i]);
212   }
213   snew(d->cgnr, s->natom);
214   for(i=0; i < s->natom; i++)
215     d->cgnr[i] = s->cgnr[i];
216   for(i=0; i < ebtsNR; i++) {
217     d->rb[i].type = s->rb[i].type;
218     d->rb[i].nb = 0;
219     d->rb[i].b = NULL;
220   }
221   merge_t_bondeds(s->rb, d->rb,FALSE,FALSE);
222 }
223
224 void copy_t_hack(t_hack *s, t_hack *d)
225 {
226   int i;
227   
228   *d = *s;
229   d->oname = safe_strdup(s->oname);
230   d->nname = safe_strdup(s->nname);
231   if (s->atom) {
232     snew(d->atom, 1);
233     *(d->atom) = *(s->atom);
234   } else
235     d->atom = NULL;
236   for(i=0; i<4; i++)
237     d->a[i] = safe_strdup(s->a[i]);
238   copy_rvec(s->newx, d->newx);
239 }
240
241 void merge_hacks_lo(int ns, t_hack *s, int *nd, t_hack **d)
242 {
243   int i;
244   
245   if (ns) {
246     srenew(*d, *nd + ns);
247     for(i=0; i < ns; i++)
248       copy_t_hack(&s[i], &(*d)[*nd + i]);
249     (*nd) += ns;
250   }
251 }
252
253 void merge_hacks(t_hackblock *s, t_hackblock *d)
254 {
255   merge_hacks_lo(s->nhack, s->hack, &d->nhack, &d->hack);
256 }
257
258 void merge_t_hackblock(t_hackblock *s, t_hackblock *d)
259 {
260   merge_hacks(s, d);
261   merge_t_bondeds(s->rb, d->rb,FALSE,FALSE);
262 }
263
264 void copy_t_hackblock(t_hackblock *s, t_hackblock *d)
265 {
266   int i;
267   
268   *d       = *s;
269   d->name  = safe_strdup(s->name);
270   d->nhack = 0;
271   d->hack  = NULL;
272   for(i=0; i<ebtsNR; i++) {
273     d->rb[i].nb=0;
274     d->rb[i].b=NULL;
275   }
276   merge_t_hackblock(s, d);
277 }
278
279 #undef safe_strdup
280
281 void dump_hb(FILE *out, int nres, t_hackblock hb[])
282 {
283   int i,j,k,l;
284   
285 #define SS(s) (s) ? (s) : "-"
286 #define SA(s) (s) ? "+" : ""
287   fprintf(out,"t_hackblock\n");
288   for(i=0; i<nres; i++) {
289     fprintf(out, "%3d %4s %2d %2d\n",
290             i, SS(hb[i].name), hb[i].nhack, hb[i].maxhack);
291     if (hb[i].nhack)
292       for(j=0; j<hb[i].nhack; j++) {
293         fprintf(out, "%d: %d %4s %4s %1s %2d %d %4s %4s %4s %4s\n", 
294                 j, hb[i].hack[j].nr, 
295                 SS(hb[i].hack[j].oname), SS(hb[i].hack[j].nname),
296                 SA(hb[i].hack[j].atom), hb[i].hack[j].tp, hb[i].hack[j].cgnr,
297                 SS(hb[i].hack[j].AI), SS(hb[i].hack[j].AJ),
298                 SS(hb[i].hack[j].AK), SS(hb[i].hack[j].AL) );
299       }
300     for(j=0; j<ebtsNR; j++)
301       if (hb[i].rb[j].nb) {
302         fprintf(out, " %c %d:", btsNames[j][0], hb[i].rb[j].nb);
303         for(k=0; k<hb[i].rb[j].nb; k++) {
304           fprintf(out, " [");
305           for(l=0; l<btsNiatoms[j]; l++)
306             fprintf(out, " %s",hb[i].rb[j].b[k].a[l]);
307           fprintf(out, " %s]",SS(hb[i].rb[j].b[k].s));
308         }
309         fprintf(out,"\n");
310       }
311     fprintf(out,"\n");
312   }
313 #undef SS
314 #undef SA
315 }
316
317 void init_t_protonate(t_protonate *protonate)
318 {
319   protonate->bInit=FALSE;
320 }