2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2012,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.
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.
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.
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.
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.
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.
44 #include "gromacs/utility/smalloc.h"
50 #include "gmx_fatal.h"
51 #include "gpp_atomtype.h"
55 void set_p_string(t_param *p, const char *s)
59 if (strlen(s) < sizeof(p->s)-1)
61 strncpy(p->s, s, sizeof(p->s));
65 gmx_fatal(FARGS, "Increase MAXSLEN in include/grompp-impl.h to at least %d,"
66 " or shorten your definition of bonds like %s to at most %d",
67 strlen(s)+1, s, MAXSLEN-1);
76 void pr_alloc (int extra, t_params *pr)
80 /* get new space for arrays */
83 gmx_fatal(FARGS, "Trying to make array smaller.\n");
89 if ((pr->nr == 0) && (pr->param != NULL))
91 fprintf(stderr, "Warning: dangling pointer at %lx\n",
92 (unsigned long)pr->param);
95 if (pr->nr+extra > pr->maxnr)
97 pr->maxnr = max(1.2*pr->maxnr, pr->maxnr + extra);
98 srenew(pr->param, pr->maxnr);
99 for (i = pr->nr; (i < pr->maxnr); i++)
101 for (j = 0; (j < MAXATOMLIST); j++)
103 pr->param[i].a[j] = 0;
105 for (j = 0; (j < MAXFORCEPARAM); j++)
107 pr->param[i].c[j] = 0;
109 set_p_string(&(pr->param[i]), "");
114 void init_plist(t_params plist[])
118 for (i = 0; (i < F_NRE); i++)
122 plist[i].param = NULL;
126 plist[i].cmap = NULL;
127 plist[i].grid_spacing = 0;
130 plist[i].cmap_types = NULL;
134 void cp_param(t_param *dest, t_param *src)
138 for (j = 0; (j < MAXATOMLIST); j++)
140 dest->a[j] = src->a[j];
142 for (j = 0; (j < MAXFORCEPARAM); j++)
144 dest->c[j] = src->c[j];
146 strncpy(dest->s, src->s, sizeof(dest->s));
149 void add_param_to_list(t_params *list, t_param *b)
153 /* allocate one position extra */
156 /* fill the arrays */
157 for (j = 0; (j < MAXFORCEPARAM); j++)
159 list->param[list->nr].c[j] = b->c[j];
161 for (j = 0; (j < MAXATOMLIST); j++)
163 list->param[list->nr].a[j] = b->a[j];
165 memset(list->param[list->nr].s, 0, sizeof(list->param[list->nr].s));
171 void init_molinfo(t_molinfo *mol)
174 mol->excl_set = FALSE;
175 mol->bProcessed = FALSE;
176 init_plist(mol->plist);
177 init_block(&mol->cgs);
178 init_block(&mol->mols);
179 init_blocka(&mol->excls);
180 init_atom(&mol->atoms);
185 void done_bt (t_params *pl)
190 void done_mi(t_molinfo *mi)
194 done_atom (&(mi->atoms));
195 done_block(&(mi->cgs));
196 done_block(&(mi->mols));
197 for (i = 0; (i < F_NRE); i++)
199 done_bt(&(mi->plist[i]));
203 /* PRINTING STRUCTURES */
205 void print_bt(FILE *out, directive d, gpp_atomtype_t at,
206 int ftype, int fsubtype, t_params plist[],
209 /* This dihp is a DIRTY patch because the dih-types do not use
210 * all four atoms to determine the type.
212 const int dihp[2][2] = { { 1, 2 }, { 0, 3 } };
214 int i, j, f, nral, nrfp;
215 gmx_bool bDih = FALSE, bSwapParity;
217 bt = &(plist[ftype]);
245 case F_CROSS_BOND_ANGLES:
248 case F_CROSS_BOND_BONDS:
298 fprintf(out, "[ %s ]\n", dir2str(d));
302 fprintf (out, "%3s %4s", "ai", "aj");
303 for (j = 2; (j < nral); j++)
305 fprintf (out, " %3c%c", 'a', 'i'+j);
310 for (j = 0; (j < 2); j++)
312 fprintf (out, "%3c%c", 'a', 'i'+dihp[f][j]);
316 fprintf (out, " funct");
317 for (j = 0; (j < nrfp); j++)
319 fprintf (out, " %12c%1d", 'c', j);
323 /* print bondtypes */
324 for (i = 0; (i < bt->nr); i++)
326 bSwapParity = (bt->param[i].C0 == NOTSET) && (bt->param[i].C1 == -1);
329 for (j = 0; (j < nral); j++)
331 fprintf (out, "%5s ", get_atomtype_name(bt->param[i].a[j], at));
336 for (j = 0; (j < 2); j++)
338 fprintf (out, "%5s ", get_atomtype_name(bt->param[i].a[dihp[f][j]], at));
341 fprintf (out, "%5d ", bSwapParity ? -f-1 : f+1);
343 if (bt->param[i].s[0])
345 fprintf(out, " %s", bt->param[i].s);
349 for (j = 0; (j < nrfp && (bt->param[i].c[j] != NOTSET)); j++)
351 fprintf (out, "%13.6e ", bt->param[i].c[j]);
361 void print_blocka(FILE *out, const char *szName,
362 const char *szIndex, const char *szA,
367 fprintf (out, "; %s\n", szName);
368 fprintf (out, "; %4s %s\n", szIndex, szA);
369 for (i = 0; (i < block->nr); i++)
371 for (i = 0; (i < block->nr); i++)
373 fprintf (out, "%6d", i+1);
374 for (j = block->index[i]; (j < ((int)block->index[i+1])); j++)
376 fprintf (out, "%5d", block->a[j]+1);
384 void print_excl(FILE *out, int natoms, t_excls excls[])
391 for (i = 0; i < natoms && !have_excl; i++)
393 have_excl = (excls[i].nr > 0);
398 fprintf (out, "[ %s ]\n", dir2str(d_exclusions));
399 fprintf (out, "; %4s %s\n", "i", "excluded from i");
400 for (i = 0; i < natoms; i++)
404 fprintf (out, "%6d ", i+1);
405 for (j = 0; j < excls[i].nr; j++)
407 fprintf (out, " %5d", excls[i].e[j]+1);
417 static double get_residue_charge(const t_atoms *atoms, int at)
422 ri = atoms->atom[at].resind;
424 while (at < atoms->nr && atoms->atom[at].resind == ri)
426 q += atoms->atom[at].q;
433 void print_atoms(FILE *out, gpp_atomtype_t atype, t_atoms *at, int *cgnr,
434 gmx_bool bRTPresname)
442 as = dir2str(d_atoms);
443 fprintf(out, "[ %s ]\n", as);
444 fprintf(out, "; %4s %10s %6s %7s%6s %6s %10s %10s %6s %10s %10s\n",
445 "nr", "type", "resnr", "residue", "atom", "cgnr", "charge", "mass", "typeB", "chargeB", "massB");
451 fprintf(debug, "This molecule has %d atoms and %d residues\n",
457 /* if the information is present... */
458 for (i = 0; (i < at->nr); i++)
460 ri = at->atom[i].resind;
461 if ((i == 0 || ri != at->atom[i-1].resind) &&
462 at->resinfo[ri].rtp != NULL)
464 qres = get_residue_charge(at, i);
465 fprintf(out, "; residue %3d %-3s rtp %-4s q ",
467 *at->resinfo[ri].name,
468 *at->resinfo[ri].rtp);
469 if (fabs(qres) < 0.001)
471 fprintf(out, " %s", "0.0");
475 fprintf(out, "%+3.1f", qres);
479 tpA = at->atom[i].type;
480 if ((tpnmA = get_atomtype_name(tpA, atype)) == NULL)
482 gmx_fatal(FARGS, "tpA = %d, i= %d in print_atoms", tpA, i);
485 fprintf(out, "%6d %10s %6d%c %5s %6s %6d %10g %10g",
490 *(at->resinfo[at->atom[i].resind].rtp) :
491 *(at->resinfo[at->atom[i].resind].name),
492 *(at->atomname[i]), cgnr[i],
493 at->atom[i].q, at->atom[i].m);
494 if (PERTURBED(at->atom[i]))
496 tpB = at->atom[i].typeB;
497 if ((tpnmB = get_atomtype_name(tpB, atype)) == NULL)
499 gmx_fatal(FARGS, "tpB = %d, i= %d in print_atoms", tpB, i);
501 fprintf(out, " %6s %10g %10g",
502 tpnmB, at->atom[i].qB, at->atom[i].mB);
504 qtot += (double)at->atom[i].q;
505 if (fabs(qtot) < 4*GMX_REAL_EPS)
509 fprintf(out, " ; qtot %.4g\n", qtot);
516 void print_bondeds(FILE *out, int natoms, directive d,
517 int ftype, int fsubtype, t_params plist[])
520 gpp_atomtype_t atype;
525 atype = init_atomtype();
529 for (i = 0; (i < natoms); i++)
532 sprintf(buf, "%4d", (i+1));
533 add_atomtype(atype, &stab, a, buf, param, 0, 0, 0, 0, 0, 0, 0);
535 print_bt(out, d, atype, ftype, fsubtype, plist, TRUE);
540 done_atomtype(atype);