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,2015,2016,2017,2018, 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.
37 /* This file is completely threadsafe - keep it that way! */
49 #include "gromacs/fileio/confio.h"
50 #include "gromacs/gmxpreprocess/gpp_nextnb.h"
51 #include "gromacs/gmxpreprocess/notset.h"
52 #include "gromacs/gmxpreprocess/pgutil.h"
53 #include "gromacs/gmxpreprocess/resall.h"
54 #include "gromacs/gmxpreprocess/topio.h"
55 #include "gromacs/gmxpreprocess/toputil.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/topology/ifunc.h"
58 #include "gromacs/utility/cstringutil.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/smalloc.h"
62 #define DIHEDRAL_WAS_SET_IN_RTP 0
63 static bool was_dihedral_set_in_rtp(const t_param *dih)
65 return dih->c[MAXFORCEPARAM-1] == DIHEDRAL_WAS_SET_IN_RTP;
68 typedef bool (*peq)(t_param *p1, t_param *p2);
70 static int acomp(const void *a1, const void *a2)
72 const t_param *p1, *p2;
75 p1 = static_cast<const t_param *>(a1);
76 p2 = static_cast<const t_param *>(a2);
77 if ((ac = (p1->aj()-p2->aj())) != 0)
81 else if ((ac = (p1->ai()-p2->ai())) != 0)
87 return (p1->ak()-p2->ak());
91 static int pcomp(const void *a1, const void *a2)
93 const t_param *p1, *p2;
96 p1 = static_cast<const t_param *>(a1);
97 p2 = static_cast<const t_param *>(a2);
98 if ((pc = (p1->ai()-p2->ai())) != 0)
104 return (p1->aj()-p2->aj());
108 static int dcomp(const void *d1, const void *d2)
110 const t_param *p1, *p2;
113 p1 = static_cast<const t_param *>(d1);
114 p2 = static_cast<const t_param *>(d2);
115 /* First sort by J & K (the two central) atoms */
116 if ((dc = (p1->aj()-p2->aj())) != 0)
120 else if ((dc = (p1->ak()-p2->ak())) != 0)
124 /* Then make sure to put rtp dihedrals before generated ones */
125 else if (was_dihedral_set_in_rtp(p1) &&
126 !was_dihedral_set_in_rtp(p2))
130 else if (!was_dihedral_set_in_rtp(p1) &&
131 was_dihedral_set_in_rtp(p2))
135 /* Then sort by I and J (two outer) atoms */
136 else if ((dc = (p1->ai()-p2->ai())) != 0)
140 else if ((dc = (p1->al()-p2->al())) != 0)
146 // AMBER force fields with type 9 dihedrals can reach here, where we sort on
147 // the contents of the string that names the macro for the parameters.
148 return strcmp(p1->s, p2->s);
153 static bool is_dihedral_on_same_bond(t_param *p1, t_param *p2)
155 return ((p1->aj() == p2->aj()) && (p1->ak() == p2->ak())) ||
156 ((p1->aj() == p2->ak()) && (p1->ak() == p2->aj()));
160 static bool preq(t_param *p1, t_param *p2)
162 return (p1->ai() == p2->ai()) && (p1->aj() == p2->aj());
165 static void rm2par(t_param p[], int *np, peq eq)
178 for (i = 1; (i < (*np)); i++)
180 if (!eq(&p[i], &p[i-1]))
185 /* Index now holds pointers to all the non-equal params,
186 * this only works when p is sorted of course
188 for (i = 0; (i < nind); i++)
190 for (j = 0; (j < MAXATOMLIST); j++)
192 p[i].a[j] = p[index[i]].a[j];
194 for (j = 0; (j < MAXFORCEPARAM); j++)
196 p[i].c[j] = p[index[i]].c[j];
198 if (p[index[i]].a[0] == p[index[i]].a[1])
202 else if (index[i] > i)
204 /* Copy the string only if it comes from somewhere else
205 * otherwise we will end up copying a random (newly freed) pointer.
206 * Since the index is sorted we only have to test for index[i] > i.
208 strcpy(p[i].s, p[index[i]].s);
216 static void cppar(t_param p[], int np, t_params plist[], int ftype)
218 int i, j, nral, nrfp;
227 for (i = 0; (i < np); i++)
229 for (j = 0; (j < nral); j++)
231 ps->param[ps->nr].a[j] = p[i].a[j];
233 for (j = 0; (j < nrfp); j++)
235 ps->param[ps->nr].c[j] = p[i].c[j];
237 for (j = 0; (j < MAXSLEN); j++)
239 ps->param[ps->nr].s[j] = p[i].s[j];
245 static void cpparam(t_param *dest, t_param *src)
249 for (j = 0; (j < MAXATOMLIST); j++)
251 dest->a[j] = src->a[j];
253 for (j = 0; (j < MAXFORCEPARAM); j++)
255 dest->c[j] = src->c[j];
257 for (j = 0; (j < MAXSLEN); j++)
259 dest->s[j] = src->s[j];
263 static void set_p(t_param *p, const int ai[4], const real *c, char *s)
267 for (j = 0; (j < 4); j++)
271 for (j = 0; (j < MAXFORCEPARAM); j++)
286 static int idcomp(const void *a, const void *b)
288 const t_param *pa, *pb;
291 pa = static_cast<const t_param *>(a);
292 pb = static_cast<const t_param *>(b);
293 if ((d = (pa->a[0]-pb->a[0])) != 0)
297 else if ((d = (pa->a[3]-pb->a[3])) != 0)
301 else if ((d = (pa->a[1]-pb->a[1])) != 0)
307 return (pa->a[2]-pb->a[2]);
311 static void sort_id(int nr, t_param ps[])
315 /* First swap order of atoms around if necessary */
316 for (i = 0; (i < nr); i++)
318 if (ps[i].a[3] < ps[i].a[0])
320 tmp = ps[i].a[3]; ps[i].a[3] = ps[i].a[0]; ps[i].a[0] = tmp;
321 tmp = ps[i].a[2]; ps[i].a[2] = ps[i].a[1]; ps[i].a[1] = tmp;
327 qsort(ps, nr, static_cast<size_t>(sizeof(ps[0])), idcomp);
331 static int n_hydro(const int a[], char ***atomname)
336 for (i = 0; (i < 4); i += 3)
338 aname = *atomname[a[i]];
339 c0 = toupper(aname[0]);
344 else if ((static_cast<int>(strlen(aname)) > 1) && (c0 >= '0') && (c0 <= '9'))
346 c1 = toupper(aname[1]);
356 /* Clean up the dihedrals (both generated and read from the .rtp
358 static void clean_dih(t_param *dih, int *ndih, t_param improper[], int nimproper,
359 t_atoms *atoms, bool bKeepAllGeneratedDihedrals,
360 bool bRemoveDihedralIfWithImproper)
365 /* Construct the list of the indices of the dihedrals
366 * (i.e. generated or read) that might be kept. */
367 snew(index, *ndih+1);
368 if (bKeepAllGeneratedDihedrals)
370 fprintf(stderr, "Keeping all generated dihedrals\n");
372 for (i = 0; i < nind; i++)
381 /* Check if generated dihedral i should be removed. The
382 * dihedrals have been sorted by dcomp() above, so all those
383 * on the same two central atoms are together, with those from
384 * the .rtp file preceding those that were automatically
385 * generated. We remove the latter if the former exist. */
386 for (i = 0; i < *ndih; i++)
388 /* Keep the dihedrals that were defined in the .rtp file,
389 * and the dihedrals that were generated and different
390 * from the last one (whether it was generated or not). */
391 if (was_dihedral_set_in_rtp(&dih[i]) ||
393 !is_dihedral_on_same_bond(&dih[i], &dih[i-1]))
402 for (i = 0; i < nind; i++)
404 bool bWasSetInRTP = was_dihedral_set_in_rtp(&dih[index[i]]);
406 if (!bWasSetInRTP && bRemoveDihedralIfWithImproper)
408 /* Remove the dihedral if there is an improper on the same
410 for (j = 0; j < nimproper && bKeep; j++)
412 bKeep = !is_dihedral_on_same_bond(&dih[index[i]], &improper[j]);
418 /* If we don't want all dihedrals, we want to select the
419 * ones with the fewest hydrogens. Note that any generated
420 * dihedrals on the same bond as an .rtp dihedral may have
421 * been already pruned above in the construction of
422 * index[]. However, their parameters are still present,
423 * and l is looping over this dihedral and all of its
424 * pruned siblings. */
425 int bestl = index[i];
426 if (!bKeepAllGeneratedDihedrals && !bWasSetInRTP)
428 /* Minimum number of hydrogens for i and l atoms */
432 is_dihedral_on_same_bond(&dih[index[i]], &dih[l]));
435 int nh = n_hydro(dih[l].a, atoms->atomname);
449 cpparam(&dih[k], &dih[bestl]);
455 for (i = k; i < *ndih; i++)
457 strcpy(dih[i].s, "");
464 static int get_impropers(t_atoms *atoms, t_hackblock hb[], t_param **improper,
467 t_rbondeds *impropers;
468 int nimproper, i, j, k, start, ninc, nalloc;
474 snew(*improper, nalloc);
476 /* Add all the impropers from the residue database to the list. */
481 for (i = 0; (i < atoms->nres); i++)
483 impropers = &hb[i].rb[ebtsIDIHS];
484 for (j = 0; (j < impropers->nb); j++)
487 for (k = 0; (k < 4) && !bStop; k++)
489 ai[k] = search_atom(impropers->b[j].a[k], start,
491 "improper", bAllowMissing);
499 if (nimproper == nalloc)
502 srenew(*improper, nalloc);
505 set_p(&((*improper)[nimproper]), ai, nullptr, impropers->b[j].s);
509 while ((start < atoms->nr) && (atoms->atom[start].resind == i))
519 static int nb_dist(t_nextnb *nnb, int ai, int aj)
531 nrexcl = nnb->nrexcl[ai];
532 for (nre = 1; (nre < nnb->nrex); nre++)
535 for (nrx = 0; (nrx < nrexcl[nre]); nrx++)
537 if ((aj == a[nrx]) && (NRE == -1))
546 static bool is_hydro(t_atoms *atoms, int ai)
548 return ((*(atoms->atomname[ai]))[0] == 'H');
551 static void get_atomnames_min(int n, char **anm,
552 int resind, t_atoms *atoms, const int *a)
556 /* Assume ascending residue numbering */
557 for (m = 0; m < n; m++)
559 if (atoms->atom[a[m]].resind < resind)
563 else if (atoms->atom[a[m]].resind > resind)
571 strcat(anm[m], *(atoms->atomname[a[m]]));
575 static void gen_excls(t_atoms *atoms, t_excls *excls, t_hackblock hb[],
579 int a, astart, i1, i2, itmp;
585 for (a = 0; a < atoms->nr; a++)
587 r = atoms->atom[a].resind;
588 if (a == atoms->nr-1 || atoms->atom[a+1].resind != r)
590 hbexcl = &hb[r].rb[ebtsEXCLS];
592 for (e = 0; e < hbexcl->nb; e++)
594 anm = hbexcl->b[e].a[0];
595 i1 = search_atom(anm, astart, atoms,
596 "exclusion", bAllowMissing);
597 anm = hbexcl->b[e].a[1];
598 i2 = search_atom(anm, astart, atoms,
599 "exclusion", bAllowMissing);
600 if (i1 != -1 && i2 != -1)
608 srenew(excls[i1].e, excls[i1].nr+1);
609 excls[i1].e[excls[i1].nr] = i2;
618 for (a = 0; a < atoms->nr; a++)
622 std::sort(excls[a].e, excls[a].e+excls[a].nr);
627 static void remove_excl(t_excls *excls, int remove)
631 for (i = remove+1; i < excls->nr; i++)
633 excls->e[i-1] = excls->e[i];
639 void clean_excls(t_nextnb *nnb, int nrexcl, t_excls excls[])
641 int i, j, j1, k, k1, l, l1, e;
646 /* extract all i-j-k-l neighbours from nnb struct */
647 for (i = 0; (i < nnb->nr); i++)
649 /* For all particles */
652 for (j = 0; (j < nnb->nrexcl[i][1]); j++)
654 /* For all first neighbours */
655 j1 = nnb->a[i][1][j];
657 for (e = 0; e < excl->nr; e++)
659 if (excl->e[e] == j1)
661 remove_excl(excl, e);
667 for (k = 0; (k < nnb->nrexcl[j1][1]); k++)
669 /* For all first neighbours of j1 */
670 k1 = nnb->a[j1][1][k];
672 for (e = 0; e < excl->nr; e++)
674 if (excl->e[e] == k1)
676 remove_excl(excl, e);
682 for (l = 0; (l < nnb->nrexcl[k1][1]); l++)
684 /* For all first neighbours of k1 */
685 l1 = nnb->a[k1][1][l];
687 for (e = 0; e < excl->nr; e++)
689 if (excl->e[e] == l1)
691 remove_excl(excl, e);
703 void generate_excls(t_nextnb *nnb, int nrexcl, t_excls excls[])
708 for (N = 1; (N < std::min(nrexcl, nnb->nrex)); N++)
710 /* extract all i-j-k-l neighbours from nnb struct */
711 for (i = 0; (i < nnb->nr); i++)
713 /* For all particles */
716 excl->nr += nnb->nrexcl[i][N];
717 srenew(excl->e, excl->nr);
718 for (j = 0; (j < nnb->nrexcl[i][N]); j++)
720 /* For all first neighbours */
721 if (nnb->a[i][N][j] != i)
723 excl->e[n++] = nnb->a[i][N][j];
730 /* Generate pairs, angles and dihedrals from .rtp settings */
731 void gen_pad(t_nextnb *nnb, t_atoms *atoms, t_restp rtp[],
732 t_params plist[], t_excls excls[], t_hackblock hb[],
735 t_param *ang, *dih, *pai, *improper;
736 t_rbondeds *hbang, *hbdih;
739 int res, minres, maxres;
740 int i, j, j1, k, k1, l, l1, m, n, i1, i2;
741 int ninc, maxang, maxdih, maxpai;
742 int nang, ndih, npai, nimproper, nbd;
746 /* These are the angles, dihedrals and pairs that we generate
747 * from the bonds. The ones that are already there from the rtp file
754 maxang = maxdih = maxpai = ninc;
760 for (i = 0; i < 4; i++)
767 gen_excls(atoms, excls, hb, bAllowMissing);
768 /* mark all entries as not matched yet */
769 for (i = 0; i < atoms->nres; i++)
771 for (j = 0; j < ebtsNR; j++)
773 for (k = 0; k < hb[i].rb[j].nb; k++)
775 hb[i].rb[j].b[k].match = FALSE;
781 /* Extract all i-j-k-l neighbours from nnb struct to generate all
782 * angles and dihedrals. */
783 for (i = 0; (i < nnb->nr); i++)
785 /* For all particles */
786 for (j = 0; (j < nnb->nrexcl[i][1]); j++)
788 /* For all first neighbours */
789 j1 = nnb->a[i][1][j];
790 for (k = 0; (k < nnb->nrexcl[j1][1]); k++)
792 /* For all first neighbours of j1 */
793 k1 = nnb->a[j1][1][k];
796 /* Generate every angle only once */
807 ang[nang].c0() = NOTSET;
808 ang[nang].c1() = NOTSET;
809 set_p_string(&(ang[nang]), "");
812 minres = atoms->atom[ang[nang].a[0]].resind;
814 for (m = 1; m < 3; m++)
816 minres = std::min(minres, atoms->atom[ang[nang].a[m]].resind);
817 maxres = std::max(maxres, atoms->atom[ang[nang].a[m]].resind);
819 res = 2*minres-maxres;
822 res += maxres-minres;
823 get_atomnames_min(3, anm, res, atoms, ang[nang].a);
824 hbang = &hb[res].rb[ebtsANGLES];
825 for (l = 0; (l < hbang->nb); l++)
827 if (strcmp(anm[1], hbang->b[l].aj()) == 0)
830 for (m = 0; m < 3; m += 2)
833 ((strcmp(anm[m], hbang->b[l].ai()) == 0) &&
834 (strcmp(anm[2-m], hbang->b[l].ak()) == 0)));
838 set_p_string(&(ang[nang]), hbang->b[l].s);
839 /* Mark that we found a match for this entry */
840 hbang->b[l].match = TRUE;
845 while (res < maxres);
849 /* Generate every dihedral, 1-4 exclusion and 1-4 interaction
853 for (l = 0; (l < nnb->nrexcl[k1][1]); l++)
855 /* For all first neighbours of k1 */
856 l1 = nnb->a[k1][1][l];
857 if ((l1 != i) && (l1 != j1))
868 for (m = 0; m < MAXFORCEPARAM; m++)
870 dih[ndih].c[m] = NOTSET;
872 set_p_string(&(dih[ndih]), "");
876 minres = atoms->atom[dih[ndih].a[0]].resind;
878 for (m = 1; m < 4; m++)
880 minres = std::min(minres, atoms->atom[dih[ndih].a[m]].resind);
881 maxres = std::max(maxres, atoms->atom[dih[ndih].a[m]].resind);
883 res = 2*minres-maxres;
886 res += maxres-minres;
887 get_atomnames_min(4, anm, res, atoms, dih[ndih].a);
888 hbdih = &hb[res].rb[ebtsPDIHS];
889 for (n = 0; (n < hbdih->nb); n++)
892 for (m = 0; m < 2; m++)
895 ((strcmp(anm[3*m], hbdih->b[n].ai()) == 0) &&
896 (strcmp(anm[1+m], hbdih->b[n].aj()) == 0) &&
897 (strcmp(anm[2-m], hbdih->b[n].ak()) == 0) &&
898 (strcmp(anm[3-3*m], hbdih->b[n].al()) == 0)));
902 set_p_string(&dih[ndih], hbdih->b[n].s);
903 /* Mark that we found a match for this entry */
904 hbdih->b[n].match = TRUE;
906 /* Set the last parameter to be able to see
907 if the dihedral was in the rtp list.
909 dih[ndih].c[MAXFORCEPARAM-1] = DIHEDRAL_WAS_SET_IN_RTP;
912 /* Set the next direct in case the rtp contains
913 multiple entries for this dihedral.
924 for (m = 0; m < MAXFORCEPARAM; m++)
926 dih[ndih].c[m] = NOTSET;
931 while (res < maxres);
944 for (m = 0; m < MAXFORCEPARAM; m++)
946 dih[ndih].c[m] = NOTSET;
948 set_p_string(&(dih[ndih]), "");
952 nbd = nb_dist(nnb, i, l1);
955 i1 = std::min(i, l1);
956 i2 = std::max(i, l1);
958 for (m = 0; m < excls[i1].nr; m++)
960 bExcl = bExcl || excls[i1].e[m] == i2;
964 if (rtp[0].bGenerateHH14Interactions ||
965 !(is_hydro(atoms, i1) && is_hydro(atoms, i2)))
974 pai[npai].c0() = NOTSET;
975 pai[npai].c1() = NOTSET;
976 set_p_string(&(pai[npai]), "");
991 /* The above approach is great in that we double-check that e.g. an angle
992 * really corresponds to three atoms connected by bonds, but this is not
993 * generally true. Go through the angle and dihedral hackblocks to add
994 * entries that we have not yet marked as matched when going through bonds.
996 for (i = 0; i < atoms->nres; i++)
998 /* Add remaining angles from hackblock */
999 hbang = &hb[i].rb[ebtsANGLES];
1000 for (j = 0; j < hbang->nb; j++)
1002 if (hbang->b[j].match == TRUE)
1004 /* We already used this entry, continue to the next */
1007 /* Hm - entry not used, let's see if we can find all atoms */
1011 srenew(ang, maxang);
1014 for (k = 0; k < 3 && bFound; k++)
1016 p = hbang->b[j].a[k];
1023 else if (p[0] == '+')
1028 ang[nang].a[k] = search_res_atom(p, res, atoms, "angle", TRUE);
1029 bFound = (ang[nang].a[k] != -1);
1031 ang[nang].c0() = NOTSET;
1032 ang[nang].c1() = NOTSET;
1036 set_p_string(&(ang[nang]), hbang->b[j].s);
1037 hbang->b[j].match = TRUE;
1038 /* Incrementing nang means we save this angle */
1043 /* Add remaining dihedrals from hackblock */
1044 hbdih = &hb[i].rb[ebtsPDIHS];
1045 for (j = 0; j < hbdih->nb; j++)
1047 if (hbdih->b[j].match == TRUE)
1049 /* We already used this entry, continue to the next */
1052 /* Hm - entry not used, let's see if we can find all atoms */
1056 srenew(dih, maxdih);
1059 for (k = 0; k < 4 && bFound; k++)
1061 p = hbdih->b[j].a[k];
1068 else if (p[0] == '+')
1073 dih[ndih].a[k] = search_res_atom(p, res, atoms, "dihedral", TRUE);
1074 bFound = (dih[ndih].a[k] != -1);
1076 for (m = 0; m < MAXFORCEPARAM; m++)
1078 dih[ndih].c[m] = NOTSET;
1083 set_p_string(&(dih[ndih]), hbdih->b[j].s);
1084 hbdih->b[j].match = TRUE;
1085 /* Incrementing ndih means we save this dihedral */
1092 /* Sort angles with respect to j-i-k (middle atom first) */
1095 qsort(ang, nang, static_cast<size_t>(sizeof(ang[0])), acomp);
1098 /* Sort dihedrals with respect to j-k-i-l (middle atoms first) */
1101 qsort(dih, ndih, static_cast<size_t>(sizeof(dih[0])), dcomp);
1104 /* Sort the pairs */
1107 qsort(pai, npai, static_cast<size_t>(sizeof(pai[0])), pcomp);
1111 /* Remove doubles, could occur in 6-rings, such as phenyls,
1112 maybe one does not want this when fudgeQQ < 1.
1114 fprintf(stderr, "Before cleaning: %d pairs\n", npai);
1115 rm2par(pai, &npai, preq);
1118 /* Get the impropers from the database */
1119 nimproper = get_impropers(atoms, hb, &improper, bAllowMissing);
1121 /* Sort the impropers */
1122 sort_id(nimproper, improper);
1126 fprintf(stderr, "Before cleaning: %d dihedrals\n", ndih);
1127 clean_dih(dih, &ndih, improper, nimproper, atoms,
1128 rtp[0].bKeepAllGeneratedDihedrals,
1129 rtp[0].bRemoveDihedralIfWithImproper);
1132 /* Now we have unique lists of angles and dihedrals
1133 * Copy them into the destination struct
1135 cppar(ang, nang, plist, F_ANGLES);
1136 cppar(dih, ndih, plist, F_PDIHS);
1137 cppar(improper, nimproper, plist, F_IDIHS);
1138 cppar(pai, npai, plist, F_LJ14);
1140 /* Remove all exclusions which are within nrexcl */
1141 clean_excls(nnb, rtp[0].nrexcl, excls);