3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35 /* This file is completely threadsafe - keep it that way! */
41 #include "hackblock.h"
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 };
50 static void free_t_bonded(t_rbonded *rb)
54 for (i=0; i<MAXATOMLIST; i++)
59 static void free_t_bondeds(t_rbondeds *rbs)
63 for(i=0; i<rbs->nb; i++)
64 free_t_bonded(&rbs->b[i]);
70 void free_t_restp(int nrtp, t_restp **rtp)
74 for(i=0; i<nrtp; i++) {
75 sfree((*rtp)[i].resname);
76 sfree((*rtp)[i].atom);
77 for(j=0; j<(*rtp)[i].natom; j++) {
78 sfree(*(*rtp)[i].atomname[j]);
79 sfree((*rtp)[i].atomname[j]);
81 sfree((*rtp)[i].atomname);
82 sfree((*rtp)[i].cgnr);
83 for(j=0; j<ebtsNR; j++)
84 free_t_bondeds(&(*rtp)[i].rb[j]);
89 void free_t_hack(int nh, t_hack **h)
104 void free_t_hackblock(int nhb, t_hackblock **hb)
108 for(i=0; i<nhb; i++) {
109 sfree((*hb)[i].name);
110 free_t_hack((*hb)[i].nhack, &(*hb)[i].hack);
111 for(j=0; j<ebtsNR; j++)
112 free_t_bondeds(&(*hb)[i].rb[j]);
117 void clear_t_hackblock(t_hackblock *hb)
125 for(i=0; i<ebtsNR; i++) {
131 void clear_t_hack(t_hack *hack)
145 hack->newx[i] = NOTSET;
148 #define safe_strdup(str) ((str != NULL) ? strdup(str) : NULL)
150 static void copy_t_rbonded(t_rbonded *s, t_rbonded *d)
154 for(i=0; i<MAXATOMLIST; i++)
155 d->a[i] = safe_strdup(s->a[i]);
156 d->s = safe_strdup(s->s);
159 static gmx_bool contains_char(t_rbonded *s,char c)
165 for(i=0; i<MAXATOMLIST; i++)
166 if (s->a[i] && s->a[i][0]==c)
172 gmx_bool merge_t_bondeds(t_rbondeds s[], t_rbondeds d[],gmx_bool bMin,gmx_bool bPlus)
175 gmx_bool bBondsRemoved;
177 bBondsRemoved = FALSE;
178 for(i=0; i < ebtsNR; i++) {
181 srenew(d[i].b, d[i].nb + s[i].nb);
182 for(j=0; j < s[i].nb; j++)
183 if (!(bMin && contains_char(&s[i].b[j],'-'))
184 && !(bPlus && contains_char(&s[i].b[j],'+'))) {
185 copy_t_rbonded(&s[i].b[j], &d[i].b[ d[i].nb ]);
187 } else if (i == ebtsBONDS) {
188 bBondsRemoved = TRUE;
193 return bBondsRemoved;
196 void copy_t_restp(t_restp *s, t_restp *d)
201 d->resname = safe_strdup(s->resname);
202 snew(d->atom, s->natom);
203 for(i=0; i < s->natom; i++)
204 d->atom[i] = s->atom[i];
205 snew(d->atomname, s->natom);
206 for(i=0; i < s->natom; i++) {
207 snew(d->atomname[i],1);
208 *d->atomname[i] = safe_strdup(*s->atomname[i]);
210 snew(d->cgnr, s->natom);
211 for(i=0; i < s->natom; i++)
212 d->cgnr[i] = s->cgnr[i];
213 for(i=0; i < ebtsNR; i++) {
214 d->rb[i].type = s->rb[i].type;
218 merge_t_bondeds(s->rb, d->rb,FALSE,FALSE);
221 void copy_t_hack(t_hack *s, t_hack *d)
226 d->oname = safe_strdup(s->oname);
227 d->nname = safe_strdup(s->nname);
230 *(d->atom) = *(s->atom);
234 d->a[i] = safe_strdup(s->a[i]);
235 copy_rvec(s->newx, d->newx);
238 void merge_hacks_lo(int ns, t_hack *s, int *nd, t_hack **d)
243 srenew(*d, *nd + ns);
244 for(i=0; i < ns; i++)
245 copy_t_hack(&s[i], &(*d)[*nd + i]);
250 void merge_hacks(t_hackblock *s, t_hackblock *d)
252 merge_hacks_lo(s->nhack, s->hack, &d->nhack, &d->hack);
255 void merge_t_hackblock(t_hackblock *s, t_hackblock *d)
258 merge_t_bondeds(s->rb, d->rb,FALSE,FALSE);
261 void copy_t_hackblock(t_hackblock *s, t_hackblock *d)
266 d->name = safe_strdup(s->name);
269 for(i=0; i<ebtsNR; i++) {
273 merge_t_hackblock(s, d);
278 void dump_hb(FILE *out, int nres, t_hackblock hb[])
282 #define SS(s) (s) ? (s) : "-"
283 #define SA(s) (s) ? "+" : ""
284 fprintf(out,"t_hackblock\n");
285 for(i=0; i<nres; i++) {
286 fprintf(out, "%3d %4s %2d %2d\n",
287 i, SS(hb[i].name), hb[i].nhack, hb[i].maxhack);
289 for(j=0; j<hb[i].nhack; j++) {
290 fprintf(out, "%d: %d %4s %4s %1s %2d %d %4s %4s %4s %4s\n",
292 SS(hb[i].hack[j].oname), SS(hb[i].hack[j].nname),
293 SA(hb[i].hack[j].atom), hb[i].hack[j].tp, hb[i].hack[j].cgnr,
294 SS(hb[i].hack[j].AI), SS(hb[i].hack[j].AJ),
295 SS(hb[i].hack[j].AK), SS(hb[i].hack[j].AL) );
297 for(j=0; j<ebtsNR; j++)
298 if (hb[i].rb[j].nb) {
299 fprintf(out, " %c %d:", btsNames[j][0], hb[i].rb[j].nb);
300 for(k=0; k<hb[i].rb[j].nb; k++) {
302 for(l=0; l<btsNiatoms[j]; l++)
303 fprintf(out, " %s",hb[i].rb[j].b[k].a[l]);
304 fprintf(out, " %s]",SS(hb[i].rb[j].b[k].s));
314 void init_t_protonate(t_protonate *protonate)
316 protonate->bInit=FALSE;