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) 2013,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.
46 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/utility/futil.h"
50 #include "gromacs/utility/fatalerror.h"
52 #include "gpp_nextnb.h"
59 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/fileio/pdbio.h"
63 #include "gromacs/fileio/filenm.h"
65 #include "gen_vsite.h"
68 #include "fflibutil.h"
71 #include "gromacs/fileio/strdb.h"
72 #include "gromacs/utility/exceptions.h"
73 #include "gromacs/utility/programcontext.h"
75 /* this must correspond to enum in pdb2top.h */
76 const char *hh[ehisNR] = { "HISD", "HISE", "HISH", "HIS1" };
78 static int missing_atoms(t_restp *rp, int resind, t_atoms *at, int i0, int i)
85 for (j = 0; j < rp->natom; j++)
87 name = *(rp->atomname[j]);
89 for (k = i0; k < i; k++)
91 bFound = (bFound || !gmx_strcasecmp(*(at->atomname[k]), name));
96 fprintf(stderr, "\nWARNING: "
97 "atom %s is missing in residue %s %d in the pdb file\n",
98 name, *(at->resinfo[resind].name), at->resinfo[resind].nr);
99 if (name[0] == 'H' || name[0] == 'h')
101 fprintf(stderr, " You might need to add atom %s to the hydrogen database of building block %s\n"
102 " in the file %s.hdb (see the manual)\n",
103 name, *(at->resinfo[resind].rtp), rp->filebase);
105 fprintf(stderr, "\n");
112 gmx_bool is_int(double x)
114 const double tol = 1e-4;
123 return (fabs(x-ix) < tol);
126 static void swap_strings(char **s, int i, int j)
136 choose_ff(const char *ffsel,
137 char *forcefield, int ff_maxlen,
138 char *ffdir, int ffdir_maxlen)
141 char **ffdirs, **ffs, **ffs_dir, *ptr;
142 int i, j, sel, cwdsel, nfound;
143 char buf[STRLEN], **desc;
147 nff = fflib_search_file_in_dirend(fflib_forcefield_itp(),
148 fflib_forcefield_dir_ext(),
153 gmx_fatal(FARGS, "No force fields found (files with name '%s' in subdirectories ending on '%s')",
154 fflib_forcefield_itp(), fflib_forcefield_dir_ext());
157 /* Replace with unix path separators */
158 if (DIR_SEPARATOR != '/')
160 for (i = 0; i < nff; i++)
162 while ( (ptr = strchr(ffdirs[i], DIR_SEPARATOR)) != NULL)
169 /* Store the force field names in ffs */
172 for (i = 0; i < nff; i++)
174 /* Remove the path from the ffdir name - use our unix standard here! */
175 ptr = strrchr(ffdirs[i], '/');
178 ffs[i] = strdup(ffdirs[i]);
179 ffs_dir[i] = low_gmxlibfn(ffdirs[i], FALSE, FALSE);
180 if (ffs_dir[i] == NULL)
182 gmx_fatal(FARGS, "Can no longer find file '%s'", ffdirs[i]);
187 ffs[i] = strdup(ptr+1);
188 ffs_dir[i] = strdup(ffdirs[i]);
190 ffs_dir[i][strlen(ffs_dir[i])-strlen(ffs[i])-1] = '\0';
191 /* Remove the extension from the ffdir name */
192 ffs[i][strlen(ffs[i])-strlen(fflib_forcefield_dir_ext())] = '\0';
200 for (i = 0; i < nff; i++)
202 if (strcmp(ffs[i], ffsel) == 0)
204 /* Matching ff name */
208 if (strncmp(ffs_dir[i], ".", 1) == 0)
225 "Force field '%s' occurs in %d places. pdb2gmx is using the one in the\n"
226 "current directory. Use interactive selection (not the -ff option) if\n"
227 "you would prefer a different one.\n", ffsel, nfound);
232 "Force field '%s' occurs in %d places, but not in the current directory.\n"
233 "Run without the -ff switch and select the force field interactively.", ffsel, nfound);
236 else if (nfound == 0)
238 gmx_fatal(FARGS, "Could not find force field '%s' in current directory, install tree or GMXDATA path.", ffsel);
244 for (i = 0; (i < nff); i++)
246 sprintf(buf, "%s%c%s%s%c%s",
247 ffs_dir[i], DIR_SEPARATOR,
248 ffs[i], fflib_forcefield_dir_ext(), DIR_SEPARATOR,
249 fflib_forcefield_doc());
252 /* We don't use fflib_open, because we don't want printf's */
253 fp = gmx_ffopen(buf, "r");
254 snew(desc[i], STRLEN);
255 get_a_line(fp, desc[i], STRLEN);
260 desc[i] = strdup(ffs[i]);
263 /* Order force fields from the same dir alphabetically
264 * and put deprecated force fields at the end.
266 for (i = 0; (i < nff); i++)
268 for (j = i+1; (j < nff); j++)
270 if (strcmp(ffs_dir[i], ffs_dir[j]) == 0 &&
271 ((desc[i][0] == '[' && desc[j][0] != '[') ||
272 ((desc[i][0] == '[' || desc[j][0] != '[') &&
273 gmx_strcasecmp(desc[i], desc[j]) > 0)))
275 swap_strings(ffdirs, i, j);
276 swap_strings(ffs, i, j);
277 swap_strings(desc, i, j);
282 printf("\nSelect the Force Field:\n");
283 for (i = 0; (i < nff); i++)
285 if (i == 0 || strcmp(ffs_dir[i-1], ffs_dir[i]) != 0)
287 if (strcmp(ffs_dir[i], ".") == 0)
289 printf("From current directory:\n");
293 printf("From '%s':\n", ffs_dir[i]);
296 printf("%2d: %s\n", i+1, desc[i]);
304 pret = fgets(buf, STRLEN, stdin);
308 sel = strtol(buf, NULL, 10);
312 while (pret == NULL || (sel < 0) || (sel >= nff));
314 /* Check for a current limitation of the fflib code.
315 * It will always read from the first ff directory in the list.
316 * This check assumes that the order of ffs matches the order
317 * in which fflib_open searches ff library files.
319 for (i = 0; i < sel; i++)
321 if (strcmp(ffs[i], ffs[sel]) == 0)
323 gmx_fatal(FARGS, "Can only select the first of multiple force field entries with directory name '%s%s' in the list. If you want to use the next entry, run pdb2gmx in a different directory or rename or move the force field directory present in the current working directory.",
324 ffs[sel], fflib_forcefield_dir_ext());
333 if (strlen(ffs[sel]) >= (size_t)ff_maxlen)
335 gmx_fatal(FARGS, "Length of force field name (%d) >= maxlen (%d)",
336 strlen(ffs[sel]), ff_maxlen);
338 strcpy(forcefield, ffs[sel]);
340 if (strlen(ffdirs[sel]) >= (size_t)ffdir_maxlen)
342 gmx_fatal(FARGS, "Length of force field dir (%d) >= maxlen (%d)",
343 strlen(ffdirs[sel]), ffdir_maxlen);
345 strcpy(ffdir, ffdirs[sel]);
347 for (i = 0; (i < nff); i++)
358 void choose_watermodel(const char *wmsel, const char *ffdir,
361 const char *fn_watermodels = "watermodels.dat";
362 char fn_list[STRLEN];
369 if (strcmp(wmsel, "none") == 0)
375 else if (strcmp(wmsel, "select") != 0)
377 *watermodel = strdup(wmsel);
382 sprintf(fn_list, "%s%c%s", ffdir, DIR_SEPARATOR, fn_watermodels);
384 if (!fflib_fexist(fn_list))
386 fprintf(stderr, "No file '%s' found, will not include a water model\n",
393 fp = fflib_open(fn_list);
394 printf("\nSelect the Water Model:\n");
397 while (get_a_line(fp, buf, STRLEN))
399 srenew(model, nwm+1);
400 snew(model[nwm], STRLEN);
401 sscanf(buf, "%s%n", model[nwm], &i);
405 fprintf(stderr, "%2d: %s\n", nwm+1, buf+i);
414 fprintf(stderr, "%2d: %s\n", nwm+1, "None");
419 pret = fgets(buf, STRLEN, stdin);
423 sel = strtol(buf, NULL, 10);
427 while (pret == NULL || sel < 0 || sel > nwm);
435 *watermodel = strdup(model[sel]);
438 for (i = 0; i < nwm; i++)
445 static int name2type(t_atoms *at, int **cgnr, gpp_atomtype_t atype,
446 t_restp restp[], gmx_residuetype_t rt)
448 int i, j, prevresind, resind, i0, prevcg, cg, curcg;
466 for (i = 0; (i < at->nr); i++)
469 if (at->atom[i].resind != resind)
472 resind = at->atom[i].resind;
473 bProt = gmx_residuetype_is_protein(rt, *(at->resinfo[resind].name));
474 bNterm = bProt && (resind == 0);
477 nmissat += missing_atoms(&restp[prevresind], prevresind, at, i0, i);
481 if (at->atom[i].m == 0)
485 fprintf(debug, "atom %d%s: curcg=%d, prevcg=%d, cg=%d\n",
486 i+1, *(at->atomname[i]), curcg, prevcg,
487 j == NOTSET ? NOTSET : restp[resind].cgnr[j]);
491 name = *(at->atomname[i]);
492 j = search_jtype(&restp[resind], name, bNterm);
493 at->atom[i].type = restp[resind].atom[j].type;
494 at->atom[i].q = restp[resind].atom[j].q;
495 at->atom[i].m = get_atomtype_massA(restp[resind].atom[j].type,
497 cg = restp[resind].cgnr[j];
498 /* A charge group number -1 signals a separate charge group
501 if ( (cg == -1) || (cg != prevcg) || (resind != prevresind) )
510 fprintf(debug, "atom %d%s: curcg=%d, qt=%g, is_int=%d\n",
511 i+1, *(at->atomname[i]), curcg, qt, is_int(qt));
522 at->atom[i].typeB = at->atom[i].type;
523 at->atom[i].qB = at->atom[i].q;
524 at->atom[i].mB = at->atom[i].m;
526 nmissat += missing_atoms(&restp[resind], resind, at, i0, i);
531 static void print_top_heavy_H(FILE *out, real mHmult)
535 fprintf(out, "; Using deuterium instead of hydrogen\n\n");
537 else if (mHmult == 4.0)
539 fprintf(out, "#define HEAVY_H\n\n");
541 else if (mHmult != 1.0)
543 fprintf(stderr, "WARNING: unsupported proton mass multiplier (%g) "
544 "in pdb2top\n", mHmult);
548 void print_top_comment(FILE *out,
549 const char *filename,
553 char ffdir_parent[STRLEN];
556 nice_header(out, filename);
557 fprintf(out, ";\tThis is a %s topology file\n;\n", bITP ? "include" : "standalone");
560 gmx::BinaryInformationSettings settings;
561 settings.generatedByHeader(true);
562 settings.linePrefix(";\t");
563 gmx::printBinaryInformation(out, gmx::getProgramContext(), settings);
565 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
567 if (strchr(ffdir, '/') == NULL)
569 fprintf(out, ";\tForce field was read from the standard Gromacs share directory.\n;\n\n");
571 else if (ffdir[0] == '.')
573 fprintf(out, ";\tForce field was read from current directory or a relative path - path added.\n;\n\n");
577 strncpy(ffdir_parent, ffdir, STRLEN-1);
578 ffdir_parent[STRLEN-1] = '\0'; /*make sure it is 0-terminated even for long string*/
579 p = strrchr(ffdir_parent, '/');
584 ";\tForce field data was read from:\n"
588 ";\tThis might be a non-standard force field location. When you use this topology, the\n"
589 ";\tforce field must either be present in the current directory, or the location\n"
590 ";\tspecified in the GMXLIB path variable or with the 'include' mdp file option.\n;\n\n",
595 void print_top_header(FILE *out, const char *filename,
596 gmx_bool bITP, const char *ffdir, real mHmult)
600 print_top_comment(out, filename, ffdir, bITP);
602 print_top_heavy_H(out, mHmult);
603 fprintf(out, "; Include forcefield parameters\n");
605 p = strrchr(ffdir, '/');
606 p = (ffdir[0] == '.' || p == NULL) ? ffdir : p+1;
608 fprintf(out, "#include \"%s/%s\"\n\n", p, fflib_forcefield_itp());
611 static void print_top_posre(FILE *out, const char *pr)
613 fprintf(out, "; Include Position restraint file\n");
614 fprintf(out, "#ifdef POSRES\n");
615 fprintf(out, "#include \"%s\"\n", pr);
616 fprintf(out, "#endif\n\n");
619 static void print_top_water(FILE *out, const char *ffdir, const char *water)
624 fprintf(out, "; Include water topology\n");
626 p = strrchr(ffdir, '/');
627 p = (ffdir[0] == '.' || p == NULL) ? ffdir : p+1;
628 fprintf(out, "#include \"%s/%s.itp\"\n", p, water);
631 fprintf(out, "#ifdef POSRES_WATER\n");
632 fprintf(out, "; Position restraint for each water oxygen\n");
633 fprintf(out, "[ position_restraints ]\n");
634 fprintf(out, ";%3s %5s %9s %10s %10s\n", "i", "funct", "fcx", "fcy", "fcz");
635 fprintf(out, "%4d %4d %10g %10g %10g\n", 1, 1, 1000.0, 1000.0, 1000.0);
636 fprintf(out, "#endif\n");
639 sprintf(buf, "%s/ions.itp", p);
641 if (fflib_fexist(buf))
643 fprintf(out, "; Include topology for ions\n");
644 fprintf(out, "#include \"%s\"\n", buf);
649 static void print_top_system(FILE *out, const char *title)
651 fprintf(out, "[ %s ]\n", dir2str(d_system));
652 fprintf(out, "; Name\n");
653 fprintf(out, "%s\n\n", title[0] ? title : "Protein");
656 void print_top_mols(FILE *out,
657 const char *title, const char *ffdir, const char *water,
658 int nincl, char **incls, int nmol, t_mols *mols)
665 fprintf(out, "; Include chain topologies\n");
666 for (i = 0; (i < nincl); i++)
668 incl = strrchr(incls[i], DIR_SEPARATOR);
675 /* Remove the path from the include name */
678 fprintf(out, "#include \"%s\"\n", incl);
685 print_top_water(out, ffdir, water);
687 print_top_system(out, title);
691 fprintf(out, "[ %s ]\n", dir2str(d_molecules));
692 fprintf(out, "; %-15s %5s\n", "Compound", "#mols");
693 for (i = 0; (i < nmol); i++)
695 fprintf(out, "%-15s %5d\n", mols[i].name, mols[i].nr);
700 void write_top(FILE *out, char *pr, char *molname,
701 t_atoms *at, gmx_bool bRTPresname,
702 int bts[], t_params plist[], t_excls excls[],
703 gpp_atomtype_t atype, int *cgnr, int nrexcl)
704 /* NOTE: nrexcl is not the size of *excl! */
706 if (at && atype && cgnr)
708 fprintf(out, "[ %s ]\n", dir2str(d_moleculetype));
709 fprintf(out, "; %-15s %5s\n", "Name", "nrexcl");
710 fprintf(out, "%-15s %5d\n\n", molname ? molname : "Protein", nrexcl);
712 print_atoms(out, atype, at, cgnr, bRTPresname);
713 print_bondeds(out, at->nr, d_bonds, F_BONDS, bts[ebtsBONDS], plist);
714 print_bondeds(out, at->nr, d_constraints, F_CONSTR, 0, plist);
715 print_bondeds(out, at->nr, d_constraints, F_CONSTRNC, 0, plist);
716 print_bondeds(out, at->nr, d_pairs, F_LJ14, 0, plist);
717 print_excl(out, at->nr, excls);
718 print_bondeds(out, at->nr, d_angles, F_ANGLES, bts[ebtsANGLES], plist);
719 print_bondeds(out, at->nr, d_dihedrals, F_PDIHS, bts[ebtsPDIHS], plist);
720 print_bondeds(out, at->nr, d_dihedrals, F_IDIHS, bts[ebtsIDIHS], plist);
721 print_bondeds(out, at->nr, d_cmap, F_CMAP, bts[ebtsCMAP], plist);
722 print_bondeds(out, at->nr, d_polarization, F_POLARIZATION, 0, plist);
723 print_bondeds(out, at->nr, d_thole_polarization, F_THOLE_POL, 0, plist);
724 print_bondeds(out, at->nr, d_vsites2, F_VSITE2, 0, plist);
725 print_bondeds(out, at->nr, d_vsites3, F_VSITE3, 0, plist);
726 print_bondeds(out, at->nr, d_vsites3, F_VSITE3FD, 0, plist);
727 print_bondeds(out, at->nr, d_vsites3, F_VSITE3FAD, 0, plist);
728 print_bondeds(out, at->nr, d_vsites3, F_VSITE3OUT, 0, plist);
729 print_bondeds(out, at->nr, d_vsites4, F_VSITE4FD, 0, plist);
730 print_bondeds(out, at->nr, d_vsites4, F_VSITE4FDN, 0, plist);
734 print_top_posre(out, pr);
739 static atom_id search_res_atom(const char *type, int resind,
741 const char *bondtype, gmx_bool bAllowMissing)
745 for (i = 0; (i < atoms->nr); i++)
747 if (atoms->atom[i].resind == resind)
749 return search_atom(type, i, atoms, bondtype, bAllowMissing);
756 static void do_ssbonds(t_params *ps, t_atoms *atoms,
757 int nssbonds, t_ssbond *ssbonds, gmx_bool bAllowMissing)
762 for (i = 0; (i < nssbonds); i++)
764 ri = ssbonds[i].res1;
765 rj = ssbonds[i].res2;
766 ai = search_res_atom(ssbonds[i].a1, ri, atoms,
767 "special bond", bAllowMissing);
768 aj = search_res_atom(ssbonds[i].a2, rj, atoms,
769 "special bond", bAllowMissing);
770 if ((ai == NO_ATID) || (aj == NO_ATID))
772 gmx_fatal(FARGS, "Trying to make impossible special bond (%s-%s)!",
773 ssbonds[i].a1, ssbonds[i].a2);
775 add_param(ps, ai, aj, NULL, NULL);
779 static void at2bonds(t_params *psb, t_hackblock *hb,
782 real long_bond_dist, real short_bond_dist)
786 real dist2, long_bond_dist2, short_bond_dist2;
789 long_bond_dist2 = sqr(long_bond_dist);
790 short_bond_dist2 = sqr(short_bond_dist);
801 fprintf(stderr, "Making bonds...\n");
803 for (resind = 0; (resind < atoms->nres) && (i < atoms->nr); resind++)
805 /* add bonds from list of bonded interactions */
806 for (j = 0; j < hb[resind].rb[ebtsBONDS].nb; j++)
808 /* Unfortunately we can not issue errors or warnings
809 * for missing atoms in bonds, as the hydrogens and terminal atoms
810 * have not been added yet.
812 ai = search_atom(hb[resind].rb[ebtsBONDS].b[j].a[0], i, atoms,
814 aj = search_atom(hb[resind].rb[ebtsBONDS].b[j].a[1], i, atoms,
816 if (ai != NO_ATID && aj != NO_ATID)
818 dist2 = distance2(x[ai], x[aj]);
819 if (dist2 > long_bond_dist2)
821 fprintf(stderr, "Warning: Long Bond (%d-%d = %g nm)\n",
822 ai+1, aj+1, sqrt(dist2));
824 else if (dist2 < short_bond_dist2)
826 fprintf(stderr, "Warning: Short Bond (%d-%d = %g nm)\n",
827 ai+1, aj+1, sqrt(dist2));
829 add_param(psb, ai, aj, NULL, hb[resind].rb[ebtsBONDS].b[j].s);
832 /* add bonds from list of hacks (each added atom gets a bond) */
833 while ( (i < atoms->nr) && (atoms->atom[i].resind == resind) )
835 for (j = 0; j < hb[resind].nhack; j++)
837 if ( ( hb[resind].hack[j].tp > 0 ||
838 hb[resind].hack[j].oname == NULL ) &&
839 strcmp(hb[resind].hack[j].a[0], *(atoms->atomname[i])) == 0)
841 switch (hb[resind].hack[j].tp)
843 case 9: /* COOH terminus */
844 add_param(psb, i, i+1, NULL, NULL); /* C-O */
845 add_param(psb, i, i+2, NULL, NULL); /* C-OA */
846 add_param(psb, i+2, i+3, NULL, NULL); /* OA-H */
849 for (k = 0; (k < hb[resind].hack[j].nr); k++)
851 add_param(psb, i, i+k+1, NULL, NULL);
858 /* we're now at the start of the next residue */
862 static int pcompar(const void *a, const void *b)
869 d = pa->a[0] - pb->a[0];
872 d = pa->a[1] - pb->a[1];
876 return strlen(pb->s) - strlen(pa->s);
884 static void clean_bonds(t_params *ps)
891 /* swap atomnumbers in bond if first larger than second: */
892 for (i = 0; (i < ps->nr); i++)
894 if (ps->param[i].a[1] < ps->param[i].a[0])
896 a = ps->param[i].a[0];
897 ps->param[i].a[0] = ps->param[i].a[1];
898 ps->param[i].a[1] = a;
903 qsort(ps->param, ps->nr, (size_t)sizeof(ps->param[0]), pcompar);
905 /* remove doubles, keep the first one always. */
907 for (i = 1; (i < ps->nr); i++)
909 if ((ps->param[i].a[0] != ps->param[j-1].a[0]) ||
910 (ps->param[i].a[1] != ps->param[j-1].a[1]) )
914 cp_param(&(ps->param[j]), &(ps->param[i]));
919 fprintf(stderr, "Number of bonds was %d, now %d\n", ps->nr, j);
924 fprintf(stderr, "No bonds\n");
928 void print_sums(t_atoms *atoms, gmx_bool bSystem)
936 where = " in system";
945 for (i = 0; (i < atoms->nr); i++)
947 m += atoms->atom[i].m;
948 qtot += atoms->atom[i].q;
951 fprintf(stderr, "Total mass%s %.3f a.m.u.\n", where, m);
952 fprintf(stderr, "Total charge%s %.3f e\n", where, qtot);
955 static void check_restp_type(const char *name, int t1, int t2)
959 gmx_fatal(FARGS, "Residues in one molecule have a different '%s' type: %d and %d", name, t1, t2);
963 static void check_restp_types(t_restp *r0, t_restp *r1)
967 check_restp_type("all dihedrals", r0->bKeepAllGeneratedDihedrals, r1->bKeepAllGeneratedDihedrals);
968 check_restp_type("nrexcl", r0->nrexcl, r1->nrexcl);
969 check_restp_type("HH14", r0->bGenerateHH14Interactions, r1->bGenerateHH14Interactions);
970 check_restp_type("remove dihedrals", r0->bRemoveDihedralIfWithImproper, r1->bRemoveDihedralIfWithImproper);
972 for (i = 0; i < ebtsNR; i++)
974 check_restp_type(btsNames[i], r0->rb[i].type, r1->rb[i].type);
978 void add_atom_to_restp(t_restp *restp, int at_start, const t_hack *hack)
982 const char *Hnum = "123456";
986 fprintf(debug,"adding atom(s) %s to atom %s in res %d%s in rtp\n",
988 * restp->atomname[at_start], resnr, restp->resname);
990 strcpy(buf, hack->nname);
991 buf[strlen(buf)+1] = '\0';
994 buf[strlen(buf)] = '-';
997 restp->natom += hack->nr;
998 srenew(restp->atom, restp->natom);
999 srenew(restp->atomname, restp->natom);
1000 srenew(restp->cgnr, restp->natom);
1001 /* shift the rest */
1002 for (k = restp->natom-1; k > at_start+hack->nr; k--)
1005 restp->atom [k - hack->nr];
1006 restp->atomname[k] =
1007 restp->atomname[k - hack->nr];
1009 restp->cgnr [k - hack->nr];
1012 for (k = 0; k < hack->nr; k++)
1014 /* set counter in atomname */
1017 buf[strlen(buf)-1] = Hnum[k];
1019 snew( restp->atomname[at_start+1+k], 1);
1020 restp->atom [at_start+1+k] = *hack->atom;
1021 *restp->atomname[at_start+1+k] = strdup(buf);
1022 if (hack->cgnr != NOTSET)
1024 restp->cgnr [at_start+1+k] = hack->cgnr;
1028 restp->cgnr [at_start+1+k] = restp->cgnr[at_start];
1033 void get_hackblocks_rtp(t_hackblock **hb, t_restp **restp,
1034 int nrtp, t_restp rtp[],
1035 int nres, t_resinfo *resinfo,
1037 t_hackblock **ntdb, t_hackblock **ctdb,
1048 /* first the termini */
1049 for (i = 0; i < nterpairs; i++)
1051 if (rn[i] >= 0 && ntdb[i] != NULL)
1053 copy_t_hackblock(ntdb[i], &(*hb)[rn[i]]);
1055 if (rc[i] >= 0 && ctdb[i] != NULL)
1057 merge_t_hackblock(ctdb[i], &(*hb)[rc[i]]);
1061 /* then the whole rtp */
1062 for (i = 0; i < nres; i++)
1064 /* Here we allow a mismatch of one character when looking for the rtp entry.
1065 * For such a mismatch there should be only one mismatching name.
1066 * This is mainly useful for small molecules such as ions.
1067 * Note that this will usually not work for protein, DNA and RNA,
1068 * since there the residue names should be listed in residuetypes.dat
1069 * and an error will have been generated earlier in the process.
1071 key = *resinfo[i].rtp;
1072 snew(resinfo[i].rtp, 1);
1073 *resinfo[i].rtp = search_rtp(key, nrtp, rtp);
1074 res = get_restp(*resinfo[i].rtp, nrtp, rtp);
1075 copy_t_restp(res, &(*restp)[i]);
1077 /* Check that we do not have different bonded types in one molecule */
1078 check_restp_types(&(*restp)[0], &(*restp)[i]);
1081 for (j = 0; j < nterpairs && tern == -1; j++)
1089 for (j = 0; j < nterpairs && terc == -1; j++)
1096 bRM = merge_t_bondeds(res->rb, (*hb)[i].rb, tern >= 0, terc >= 0);
1098 if (bRM && ((tern >= 0 && ntdb[tern] == NULL) ||
1099 (terc >= 0 && ctdb[terc] == NULL)))
1101 gmx_fatal(FARGS, "There is a dangling bond at at least one of the terminal ends and the force field does not provide terminal entries or files. Fix your terminal residues so that they match the residue database (.rtp) entries, or provide terminal database entries (.tdb).");
1103 if (bRM && ((tern >= 0 && ntdb[tern]->nhack == 0) ||
1104 (terc >= 0 && ctdb[terc]->nhack == 0)))
1106 gmx_fatal(FARGS, "There is a dangling bond at at least one of the terminal ends. Fix your coordinate file, add a new terminal database entry (.tdb), or select the proper existing terminal entry.");
1110 /* now perform t_hack's on t_restp's,
1111 i.e. add's and deletes from termini database will be
1112 added to/removed from residue topology
1113 we'll do this on one big dirty loop, so it won't make easy reading! */
1114 for (i = 0; i < nres; i++)
1116 for (j = 0; j < (*hb)[i].nhack; j++)
1118 if ( (*hb)[i].hack[j].nr)
1120 /* find atom in restp */
1121 for (l = 0; l < (*restp)[i].natom; l++)
1123 if ( ( (*hb)[i].hack[j].oname == NULL &&
1124 strcmp((*hb)[i].hack[j].a[0], *(*restp)[i].atomname[l]) == 0 ) ||
1125 ( (*hb)[i].hack[j].oname != NULL &&
1126 strcmp((*hb)[i].hack[j].oname, *(*restp)[i].atomname[l]) == 0 ) )
1131 if (l == (*restp)[i].natom)
1133 /* If we are doing an atom rename only, we don't need
1134 * to generate a fatal error if the old name is not found
1137 /* Deleting can happen also only on the input atoms,
1138 * not necessarily always on the rtp entry.
1140 if (!((*hb)[i].hack[j].oname != NULL &&
1141 (*hb)[i].hack[j].nname != NULL) &&
1142 !((*hb)[i].hack[j].oname != NULL &&
1143 (*hb)[i].hack[j].nname == NULL))
1146 "atom %s not found in buiding block %d%s "
1147 "while combining tdb and rtp",
1148 (*hb)[i].hack[j].oname != NULL ?
1149 (*hb)[i].hack[j].oname : (*hb)[i].hack[j].a[0],
1150 i+1, *resinfo[i].rtp);
1155 if ( (*hb)[i].hack[j].oname == NULL)
1158 add_atom_to_restp(&(*restp)[i], l, &(*hb)[i].hack[j]);
1163 if ( (*hb)[i].hack[j].nname == NULL)
1165 /* we're deleting */
1168 fprintf(debug, "deleting atom %s from res %d%s in rtp\n",
1169 *(*restp)[i].atomname[l],
1170 i+1, (*restp)[i].resname);
1172 /* shift the rest */
1173 (*restp)[i].natom--;
1174 for (k = l; k < (*restp)[i].natom; k++)
1176 (*restp)[i].atom [k] = (*restp)[i].atom [k+1];
1177 (*restp)[i].atomname[k] = (*restp)[i].atomname[k+1];
1178 (*restp)[i].cgnr [k] = (*restp)[i].cgnr [k+1];
1180 /* give back space */
1181 srenew((*restp)[i].atom, (*restp)[i].natom);
1182 srenew((*restp)[i].atomname, (*restp)[i].natom);
1183 srenew((*restp)[i].cgnr, (*restp)[i].natom);
1185 else /* nname != NULL */
1186 { /* we're replacing */
1189 fprintf(debug, "replacing atom %s by %s in res %d%s in rtp\n",
1190 *(*restp)[i].atomname[l], (*hb)[i].hack[j].nname,
1191 i+1, (*restp)[i].resname);
1193 snew( (*restp)[i].atomname[l], 1);
1194 (*restp)[i].atom[l] = *(*hb)[i].hack[j].atom;
1195 *(*restp)[i].atomname[l] = strdup((*hb)[i].hack[j].nname);
1196 if ( (*hb)[i].hack[j].cgnr != NOTSET)
1198 (*restp)[i].cgnr [l] = (*hb)[i].hack[j].cgnr;
1208 static gmx_bool atomname_cmp_nr(const char *anm, t_hack *hack, int *nr)
1215 return (gmx_strcasecmp(anm, hack->nname) == 0);
1219 if (isdigit(anm[strlen(anm)-1]))
1221 *nr = anm[strlen(anm)-1] - '0';
1227 if (*nr <= 0 || *nr > hack->nr)
1233 return (strlen(anm) == strlen(hack->nname) + 1 &&
1234 gmx_strncasecmp(anm, hack->nname, strlen(hack->nname)) == 0);
1239 static gmx_bool match_atomnames_with_rtp_atom(t_atoms *pdba, rvec *x, int atind,
1240 t_restp *rptr, t_hackblock *hbr,
1245 char *oldnm, *newnm;
1247 char *start_at, buf[STRLEN];
1249 gmx_bool bReplaceReplace, bFoundInAdd;
1252 oldnm = *pdba->atomname[atind];
1253 resnr = pdba->resinfo[pdba->atom[atind].resind].nr;
1256 for (j = 0; j < hbr->nhack; j++)
1258 if (hbr->hack[j].oname != NULL && hbr->hack[j].nname != NULL &&
1259 gmx_strcasecmp(oldnm, hbr->hack[j].oname) == 0)
1261 /* This is a replace entry. */
1262 /* Check if we are not replacing a replaced atom. */
1263 bReplaceReplace = FALSE;
1264 for (k = 0; k < hbr->nhack; k++)
1267 hbr->hack[k].oname != NULL && hbr->hack[k].nname != NULL &&
1268 gmx_strcasecmp(hbr->hack[k].nname, hbr->hack[j].oname) == 0)
1270 /* The replace in hack[j] replaces an atom that
1271 * was already replaced in hack[k], we do not want
1272 * second or higher level replaces at this stage.
1274 bReplaceReplace = TRUE;
1277 if (bReplaceReplace)
1279 /* Skip this replace. */
1283 /* This atom still has the old name, rename it */
1284 newnm = hbr->hack[j].nname;
1285 for (k = 0; k < rptr->natom; k++)
1287 if (gmx_strcasecmp(newnm, *rptr->atomname[k]) == 0)
1292 if (k == rptr->natom)
1294 /* The new name is not present in the rtp.
1295 * We need to apply the replace also to the rtp entry.
1298 /* We need to find the add hack that can add this atom
1299 * to find out after which atom it should be added.
1301 bFoundInAdd = FALSE;
1302 for (k = 0; k < hbr->nhack; k++)
1304 if (hbr->hack[k].oname == NULL &&
1305 hbr->hack[k].nname != NULL &&
1306 atomname_cmp_nr(newnm, &hbr->hack[k], &anmnr))
1310 start_at = hbr->hack[k].a[0];
1314 sprintf(buf, "%s%d", hbr->hack[k].nname, anmnr-1);
1317 for (start_nr = 0; start_nr < rptr->natom; start_nr++)
1319 if (gmx_strcasecmp(start_at, (*rptr->atomname[start_nr])) == 0)
1324 if (start_nr == rptr->natom)
1326 gmx_fatal(FARGS, "Could not find atom '%s' in residue building block '%s' to add atom '%s' to",
1327 start_at, rptr->resname, newnm);
1329 /* We can add the atom after atom start_nr */
1330 add_atom_to_restp(rptr, start_nr, &hbr->hack[j]);
1338 gmx_fatal(FARGS, "Could not find an 'add' entry for atom named '%s' corresponding to the 'replace' entry from atom name '%s' to '%s' for tdb or hdb database of residue type '%s'",
1340 hbr->hack[j].oname, hbr->hack[j].nname,
1347 printf("Renaming atom '%s' in residue '%s' %d to '%s'\n",
1348 oldnm, rptr->resname, resnr, newnm);
1350 /* Rename the atom in pdba */
1351 snew(pdba->atomname[atind], 1);
1352 *pdba->atomname[atind] = strdup(newnm);
1354 else if (hbr->hack[j].oname != NULL && hbr->hack[j].nname == NULL &&
1355 gmx_strcasecmp(oldnm, hbr->hack[j].oname) == 0)
1357 /* This is a delete entry, check if this atom is present
1358 * in the rtp entry of this residue.
1360 for (k = 0; k < rptr->natom; k++)
1362 if (gmx_strcasecmp(oldnm, *rptr->atomname[k]) == 0)
1367 if (k == rptr->natom)
1369 /* This atom is not present in the rtp entry,
1370 * delete is now from the input pdba.
1374 printf("Deleting atom '%s' in residue '%s' %d\n",
1375 oldnm, rptr->resname, resnr);
1377 /* We should free the atom name,
1378 * but it might be used multiple times in the symtab.
1379 * sfree(pdba->atomname[atind]);
1381 for (k = atind+1; k < pdba->nr; k++)
1383 pdba->atom[k-1] = pdba->atom[k];
1384 pdba->atomname[k-1] = pdba->atomname[k];
1385 copy_rvec(x[k], x[k-1]);
1396 void match_atomnames_with_rtp(t_restp restp[], t_hackblock hb[],
1397 t_atoms *pdba, rvec *x,
1404 for (i = 0; i < pdba->nr; i++)
1406 oldnm = *pdba->atomname[i];
1407 rptr = &restp[pdba->atom[i].resind];
1408 for (j = 0; (j < rptr->natom); j++)
1410 if (gmx_strcasecmp(oldnm, *(rptr->atomname[j])) == 0)
1415 if (j == rptr->natom)
1417 /* Not found yet, check if we have to rename this atom */
1418 if (match_atomnames_with_rtp_atom(pdba, x, i,
1419 rptr, &(hb[pdba->atom[i].resind]),
1422 /* We deleted this atom, decrease the atom counter by 1. */
1429 #define NUM_CMAP_ATOMS 5
1430 static void gen_cmap(t_params *psb, t_restp *restp, t_atoms *atoms, gmx_residuetype_t rt)
1432 int residx, i, j, k;
1434 t_resinfo *resinfo = atoms->resinfo;
1435 int nres = atoms->nres;
1437 atom_id cmap_atomid[NUM_CMAP_ATOMS];
1438 int cmap_chainnum = -1, this_residue_index;
1449 fprintf(stderr, "Making cmap torsions...");
1451 /* End loop at nres-1, since the very last residue does not have a +N atom, and
1452 * therefore we get a valgrind invalid 4 byte read error with atom am */
1453 for (residx = 0; residx < nres-1; residx++)
1455 /* Add CMAP terms from the list of CMAP interactions */
1456 for (j = 0; j < restp[residx].rb[ebtsCMAP].nb; j++)
1459 /* Loop over atoms in a candidate CMAP interaction and
1460 * check that they exist, are from the same chain and are
1461 * from residues labelled as protein. */
1462 for (k = 0; k < NUM_CMAP_ATOMS && bAddCMAP; k++)
1464 cmap_atomid[k] = search_atom(restp[residx].rb[ebtsCMAP].b[j].a[k],
1465 i, atoms, ptr, TRUE);
1466 bAddCMAP = bAddCMAP && (cmap_atomid[k] != NO_ATID);
1469 /* This break is necessary, because cmap_atomid[k]
1470 * == NO_ATID cannot be safely used as an index
1471 * into the atom array. */
1474 this_residue_index = atoms->atom[cmap_atomid[k]].resind;
1477 cmap_chainnum = resinfo[this_residue_index].chainnum;
1481 /* Does the residue for this atom have the same
1482 * chain number as the residues for previous
1484 bAddCMAP = bAddCMAP &&
1485 cmap_chainnum == resinfo[this_residue_index].chainnum;
1487 bAddCMAP = bAddCMAP && gmx_residuetype_is_protein(rt, *(resinfo[this_residue_index].name));
1492 add_cmap_param(psb, cmap_atomid[0], cmap_atomid[1], cmap_atomid[2], cmap_atomid[3], cmap_atomid[4], restp[residx].rb[ebtsCMAP].b[j].s);
1496 if (residx < nres-1)
1498 while (atoms->atom[i].resind < residx+1)
1505 /* Start the next residue */
1509 scrub_charge_groups(int *cgnr, int natoms)
1513 for (i = 0; i < natoms; i++)
1520 void pdb2top(FILE *top_file, char *posre_fn, char *molname,
1521 t_atoms *atoms, rvec **x, gpp_atomtype_t atype, t_symtab *tab,
1522 int nrtp, t_restp rtp[],
1523 t_restp *restp, t_hackblock *hb,
1524 gmx_bool bAllowMissing,
1525 gmx_bool bVsites, gmx_bool bVsiteAromatics,
1528 int nssbonds, t_ssbond *ssbonds,
1529 real long_bond_dist, real short_bond_dist,
1530 gmx_bool bDeuterate, gmx_bool bChargeGroups, gmx_bool bCmap,
1531 gmx_bool bRenumRes, gmx_bool bRTPresname)
1537 t_params plist[F_NRE];
1544 gmx_residuetype_t rt;
1547 gmx_residuetype_init(&rt);
1551 print_resall(debug, atoms->nres, restp, atype);
1552 dump_hb(debug, atoms->nres, hb);
1556 at2bonds(&(plist[F_BONDS]), hb,
1558 long_bond_dist, short_bond_dist);
1560 /* specbonds: disulphide bonds & heme-his */
1561 do_ssbonds(&(plist[F_BONDS]),
1562 atoms, nssbonds, ssbonds,
1565 nmissat = name2type(atoms, &cgnr, atype, restp, rt);
1570 fprintf(stderr, "There were %d missing atoms in molecule %s\n",
1575 gmx_fatal(FARGS, "There were %d missing atoms in molecule %s, if you want to use this incomplete topology anyhow, use the option -missing",
1580 /* Cleanup bonds (sort and rm doubles) */
1581 clean_bonds(&(plist[F_BONDS]));
1583 snew(vsite_type, atoms->nr);
1584 for (i = 0; i < atoms->nr; i++)
1586 vsite_type[i] = NOTSET;
1590 /* determine which atoms will be vsites and add dummy masses
1591 also renumber atom numbers in plist[0..F_NRE]! */
1592 do_vsites(nrtp, rtp, atype, atoms, tab, x, plist,
1593 &vsite_type, &cgnr, mHmult, bVsiteAromatics, ffdir);
1596 /* Make Angles and Dihedrals */
1597 fprintf(stderr, "Generating angles, dihedrals and pairs...\n");
1598 snew(excls, atoms->nr);
1599 init_nnb(&nnb, atoms->nr, 4);
1600 gen_nnb(&nnb, plist);
1601 print_nnb(&nnb, "NNB");
1602 gen_pad(&nnb, atoms, restp, plist, excls, hb, bAllowMissing);
1608 gen_cmap(&(plist[F_CMAP]), restp, atoms, rt);
1609 if (plist[F_CMAP].nr > 0)
1611 fprintf(stderr, "There are %4d cmap torsion pairs\n",
1616 /* set mass of all remaining hydrogen atoms */
1619 do_h_mass(&(plist[F_BONDS]), vsite_type, atoms, mHmult, bDeuterate);
1623 /* Cleanup bonds (sort and rm doubles) */
1624 /* clean_bonds(&(plist[F_BONDS]));*/
1627 "There are %4d dihedrals, %4d impropers, %4d angles\n"
1628 " %4d pairs, %4d bonds and %4d virtual sites\n",
1629 plist[F_PDIHS].nr, plist[F_IDIHS].nr, plist[F_ANGLES].nr,
1630 plist[F_LJ14].nr, plist[F_BONDS].nr,
1631 plist[F_VSITE2].nr +
1632 plist[F_VSITE3].nr +
1633 plist[F_VSITE3FD].nr +
1634 plist[F_VSITE3FAD].nr +
1635 plist[F_VSITE3OUT].nr +
1636 plist[F_VSITE4FD].nr +
1637 plist[F_VSITE4FDN].nr );
1639 print_sums(atoms, FALSE);
1641 if (FALSE == bChargeGroups)
1643 scrub_charge_groups(cgnr, atoms->nr);
1648 for (i = 0; i < atoms->nres; i++)
1650 atoms->resinfo[i].nr = i + 1;
1651 atoms->resinfo[i].ic = ' ';
1657 fprintf(stderr, "Writing topology\n");
1658 /* We can copy the bonded types from the first restp,
1659 * since the types have to be identical for all residues in one molecule.
1661 for (i = 0; i < ebtsNR; i++)
1663 bts[i] = restp[0].rb[i].type;
1665 write_top(top_file, posre_fn, molname,
1667 bts, plist, excls, atype, cgnr, restp[0].nrexcl);
1671 free_t_hackblock(atoms->nres, &hb);
1672 free_t_restp(atoms->nres, &restp);
1673 gmx_residuetype_destroy(rt);
1675 /* we should clean up hb and restp here, but that is a *L*O*T* of work! */
1677 for (i = 0; i < F_NRE; i++)
1679 sfree(plist[i].param);
1681 for (i = 0; i < atoms->nr; i++)