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, 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! */
50 #include "gromacs/fileio/confio.h"
51 #include "gromacs/gmxpreprocess/gpp_nextnb.h"
52 #include "gromacs/gmxpreprocess/notset.h"
53 #include "gromacs/gmxpreprocess/pgutil.h"
54 #include "gromacs/gmxpreprocess/resall.h"
55 #include "gromacs/gmxpreprocess/topio.h"
56 #include "gromacs/gmxpreprocess/toputil.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/topology/ifunc.h"
59 #include "gromacs/utility/cstringutil.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/smalloc.h"
63 #define DIHEDRAL_WAS_SET_IN_RTP 0
64 static gmx_bool
was_dihedral_set_in_rtp(const t_param
*dih
)
66 return dih
->c
[MAXFORCEPARAM
-1] == DIHEDRAL_WAS_SET_IN_RTP
;
69 typedef gmx_bool (*peq
)(t_param
*p1
, t_param
*p2
);
71 static int acomp(const void *a1
, const void *a2
)
73 const t_param
*p1
, *p2
;
76 p1
= static_cast<const t_param
*>(a1
);
77 p2
= static_cast<const t_param
*>(a2
);
78 if ((ac
= (p1
->aj()-p2
->aj())) != 0)
82 else if ((ac
= (p1
->ai()-p2
->ai())) != 0)
88 return (p1
->ak()-p2
->ak());
92 static int pcomp(const void *a1
, const void *a2
)
94 const t_param
*p1
, *p2
;
97 p1
= static_cast<const t_param
*>(a1
);
98 p2
= static_cast<const t_param
*>(a2
);
99 if ((pc
= (p1
->ai()-p2
->ai())) != 0)
105 return (p1
->aj()-p2
->aj());
109 static int dcomp(const void *d1
, const void *d2
)
111 const t_param
*p1
, *p2
;
114 p1
= static_cast<const t_param
*>(d1
);
115 p2
= static_cast<const t_param
*>(d2
);
116 /* First sort by J & K (the two central) atoms */
117 if ((dc
= (p1
->aj()-p2
->aj())) != 0)
121 else if ((dc
= (p1
->ak()-p2
->ak())) != 0)
125 /* Then make sure to put rtp dihedrals before generated ones */
126 else if (was_dihedral_set_in_rtp(p1
) &&
127 !was_dihedral_set_in_rtp(p2
))
131 else if (!was_dihedral_set_in_rtp(p1
) &&
132 was_dihedral_set_in_rtp(p2
))
136 /* Finally, sort by I and J (two outer) atoms */
137 else if ((dc
= (p1
->ai()-p2
->ai())) != 0)
143 return (p1
->al()-p2
->al());
148 static gmx_bool
is_dihedral_on_same_bond(t_param
*p1
, t_param
*p2
)
150 if (((p1
->aj() == p2
->aj()) && (p1
->ak() == p2
->ak())) ||
151 ((p1
->aj() == p2
->ak()) && (p1
->ak() == p2
->aj())))
162 static gmx_bool
preq(t_param
*p1
, t_param
*p2
)
164 if ((p1
->ai() == p2
->ai()) && (p1
->aj() == p2
->aj()))
174 static void rm2par(t_param p
[], int *np
, peq eq
)
187 for (i
= 1; (i
< (*np
)); i
++)
189 if (!eq(&p
[i
], &p
[i
-1]))
194 /* Index now holds pointers to all the non-equal params,
195 * this only works when p is sorted of course
197 for (i
= 0; (i
< nind
); i
++)
199 for (j
= 0; (j
< MAXATOMLIST
); j
++)
201 p
[i
].a
[j
] = p
[index
[i
]].a
[j
];
203 for (j
= 0; (j
< MAXFORCEPARAM
); j
++)
205 p
[i
].c
[j
] = p
[index
[i
]].c
[j
];
207 if (p
[index
[i
]].a
[0] == p
[index
[i
]].a
[1])
212 "Something VERY strange is going on in rm2par (gen_ad.c)\n"
213 "a[0] %d a[1] %d a[2] %d a[3] %d\n",
214 p
[i
].a
[0], p
[i
].a
[1], p
[i
].a
[2], p
[i
].a
[3]);
218 else if (index
[i
] > i
)
220 /* Copy the string only if it comes from somewhere else
221 * otherwise we will end up copying a random (newly freed) pointer.
222 * Since the index is sorted we only have to test for index[i] > i.
224 strcpy(p
[i
].s
, p
[index
[i
]].s
);
232 static void cppar(t_param p
[], int np
, t_params plist
[], int ftype
)
234 int i
, j
, nral
, nrfp
;
243 for (i
= 0; (i
< np
); i
++)
245 for (j
= 0; (j
< nral
); j
++)
247 ps
->param
[ps
->nr
].a
[j
] = p
[i
].a
[j
];
249 for (j
= 0; (j
< nrfp
); j
++)
251 ps
->param
[ps
->nr
].c
[j
] = p
[i
].c
[j
];
253 for (j
= 0; (j
< MAXSLEN
); j
++)
255 ps
->param
[ps
->nr
].s
[j
] = p
[i
].s
[j
];
261 static void cpparam(t_param
*dest
, t_param
*src
)
265 for (j
= 0; (j
< MAXATOMLIST
); j
++)
267 dest
->a
[j
] = src
->a
[j
];
269 for (j
= 0; (j
< MAXFORCEPARAM
); j
++)
271 dest
->c
[j
] = src
->c
[j
];
273 for (j
= 0; (j
< MAXSLEN
); j
++)
275 dest
->s
[j
] = src
->s
[j
];
279 static void set_p(t_param
*p
, int ai
[4], real
*c
, char *s
)
283 for (j
= 0; (j
< 4); j
++)
287 for (j
= 0; (j
< MAXFORCEPARAM
); j
++)
302 static int int_comp(const void *a
, const void *b
)
304 return (*(int *)a
) - (*(int *)b
);
307 static int idcomp(const void *a
, const void *b
)
309 const t_param
*pa
, *pb
;
312 pa
= static_cast<const t_param
*>(a
);
313 pb
= static_cast<const t_param
*>(b
);
314 if ((d
= (pa
->a
[0]-pb
->a
[0])) != 0)
318 else if ((d
= (pa
->a
[3]-pb
->a
[3])) != 0)
322 else if ((d
= (pa
->a
[1]-pb
->a
[1])) != 0)
328 return (int) (pa
->a
[2]-pb
->a
[2]);
332 static void sort_id(int nr
, t_param ps
[])
336 /* First swap order of atoms around if necessary */
337 for (i
= 0; (i
< nr
); i
++)
339 if (ps
[i
].a
[3] < ps
[i
].a
[0])
341 tmp
= ps
[i
].a
[3]; ps
[i
].a
[3] = ps
[i
].a
[0]; ps
[i
].a
[0] = tmp
;
342 tmp
= ps
[i
].a
[2]; ps
[i
].a
[2] = ps
[i
].a
[1]; ps
[i
].a
[1] = tmp
;
348 qsort(ps
, nr
, (size_t)sizeof(ps
[0]), idcomp
);
352 static int n_hydro(int a
[], char ***atomname
)
357 for (i
= 0; (i
< 4); i
+= 3)
359 aname
= *atomname
[a
[i
]];
360 c0
= toupper(aname
[0]);
365 else if (((int)strlen(aname
) > 1) && (c0
>= '0') && (c0
<= '9'))
367 c1
= toupper(aname
[1]);
377 /* Clean up the dihedrals (both generated and read from the .rtp
379 static void clean_dih(t_param
*dih
, int *ndih
, t_param improper
[], int nimproper
,
380 t_atoms
*atoms
, gmx_bool bKeepAllGeneratedDihedrals
,
381 gmx_bool bRemoveDihedralIfWithImproper
)
386 /* Construct the list of the indices of the dihedrals
387 * (i.e. generated or read) that might be kept. */
388 snew(index
, *ndih
+1);
389 if (bKeepAllGeneratedDihedrals
)
391 fprintf(stderr
, "Keeping all generated dihedrals\n");
393 for (i
= 0; i
< nind
; i
++)
402 /* Check if generated dihedral i should be removed. The
403 * dihedrals have been sorted by dcomp() above, so all those
404 * on the same two central atoms are together, with those from
405 * the .rtp file preceding those that were automatically
406 * generated. We remove the latter if the former exist. */
407 for (i
= 0; i
< *ndih
; i
++)
409 /* Keep the dihedrals that were defined in the .rtp file,
410 * and the dihedrals that were generated and different
411 * from the last one (whether it was generated or not). */
412 if (was_dihedral_set_in_rtp(&dih
[i
]) ||
414 !is_dihedral_on_same_bond(&dih
[i
], &dih
[i
-1]))
423 for (i
= 0; i
< nind
; i
++)
425 gmx_bool bWasSetInRTP
= was_dihedral_set_in_rtp(&dih
[index
[i
]]);
426 gmx_bool bKeep
= TRUE
;
427 if (!bWasSetInRTP
&& bRemoveDihedralIfWithImproper
)
429 /* Remove the dihedral if there is an improper on the same
431 for (j
= 0; j
< nimproper
&& bKeep
; j
++)
433 bKeep
= !is_dihedral_on_same_bond(&dih
[index
[i
]], &improper
[j
]);
439 /* If we don't want all dihedrals, we want to select the
440 * ones with the fewest hydrogens. Note that any generated
441 * dihedrals on the same bond as an .rtp dihedral may have
442 * been already pruned above in the construction of
443 * index[]. However, their parameters are still present,
444 * and l is looping over this dihedral and all of its
445 * pruned siblings. */
446 int bestl
= index
[i
];
447 if (!bKeepAllGeneratedDihedrals
&& !bWasSetInRTP
)
449 /* Minimum number of hydrogens for i and l atoms */
453 is_dihedral_on_same_bond(&dih
[index
[i
]], &dih
[l
]));
456 int nh
= n_hydro(dih
[l
].a
, atoms
->atomname
);
470 cpparam(&dih
[k
], &dih
[bestl
]);
476 for (i
= k
; i
< *ndih
; i
++)
478 strcpy(dih
[i
].s
, "");
485 static int get_impropers(t_atoms
*atoms
, t_hackblock hb
[], t_param
**improper
,
486 gmx_bool bAllowMissing
)
488 t_rbondeds
*impropers
;
489 int nimproper
, i
, j
, k
, start
, ninc
, nalloc
;
495 snew(*improper
, nalloc
);
497 /* Add all the impropers from the residue database to the list. */
502 for (i
= 0; (i
< atoms
->nres
); i
++)
504 impropers
= &hb
[i
].rb
[ebtsIDIHS
];
505 for (j
= 0; (j
< impropers
->nb
); j
++)
508 for (k
= 0; (k
< 4) && !bStop
; k
++)
510 ai
[k
] = search_atom(impropers
->b
[j
].a
[k
], start
,
512 "improper", bAllowMissing
);
520 if (nimproper
== nalloc
)
523 srenew(*improper
, nalloc
);
526 set_p(&((*improper
)[nimproper
]), ai
, nullptr, impropers
->b
[j
].s
);
530 while ((start
< atoms
->nr
) && (atoms
->atom
[start
].resind
== i
))
540 static int nb_dist(t_nextnb
*nnb
, int ai
, int aj
)
552 nrexcl
= nnb
->nrexcl
[ai
];
553 for (nre
= 1; (nre
< nnb
->nrex
); nre
++)
556 for (nrx
= 0; (nrx
< nrexcl
[nre
]); nrx
++)
558 if ((aj
== a
[nrx
]) && (NRE
== -1))
567 static gmx_bool
is_hydro(t_atoms
*atoms
, int ai
)
569 return ((*(atoms
->atomname
[ai
]))[0] == 'H');
572 static void get_atomnames_min(int n
, char **anm
,
573 int resind
, t_atoms
*atoms
, int *a
)
577 /* Assume ascending residue numbering */
578 for (m
= 0; m
< n
; m
++)
580 if (atoms
->atom
[a
[m
]].resind
< resind
)
584 else if (atoms
->atom
[a
[m
]].resind
> resind
)
592 strcat(anm
[m
], *(atoms
->atomname
[a
[m
]]));
596 static void gen_excls(t_atoms
*atoms
, t_excls
*excls
, t_hackblock hb
[],
597 gmx_bool bAllowMissing
)
600 int a
, astart
, i1
, i2
, itmp
;
606 for (a
= 0; a
< atoms
->nr
; a
++)
608 r
= atoms
->atom
[a
].resind
;
609 if (a
== atoms
->nr
-1 || atoms
->atom
[a
+1].resind
!= r
)
611 hbexcl
= &hb
[r
].rb
[ebtsEXCLS
];
613 for (e
= 0; e
< hbexcl
->nb
; e
++)
615 anm
= hbexcl
->b
[e
].a
[0];
616 i1
= search_atom(anm
, astart
, atoms
,
617 "exclusion", bAllowMissing
);
618 anm
= hbexcl
->b
[e
].a
[1];
619 i2
= search_atom(anm
, astart
, atoms
,
620 "exclusion", bAllowMissing
);
621 if (i1
!= -1 && i2
!= -1)
629 srenew(excls
[i1
].e
, excls
[i1
].nr
+1);
630 excls
[i1
].e
[excls
[i1
].nr
] = i2
;
639 for (a
= 0; a
< atoms
->nr
; a
++)
643 qsort(excls
[a
].e
, excls
[a
].nr
, (size_t)sizeof(int), int_comp
);
648 static void remove_excl(t_excls
*excls
, int remove
)
652 for (i
= remove
+1; i
< excls
->nr
; i
++)
654 excls
->e
[i
-1] = excls
->e
[i
];
660 void clean_excls(t_nextnb
*nnb
, int nrexcl
, t_excls excls
[])
662 int i
, j
, j1
, k
, k1
, l
, l1
, e
;
667 /* extract all i-j-k-l neighbours from nnb struct */
668 for (i
= 0; (i
< nnb
->nr
); i
++)
670 /* For all particles */
673 for (j
= 0; (j
< nnb
->nrexcl
[i
][1]); j
++)
675 /* For all first neighbours */
676 j1
= nnb
->a
[i
][1][j
];
678 for (e
= 0; e
< excl
->nr
; e
++)
680 if (excl
->e
[e
] == j1
)
682 remove_excl(excl
, e
);
688 for (k
= 0; (k
< nnb
->nrexcl
[j1
][1]); k
++)
690 /* For all first neighbours of j1 */
691 k1
= nnb
->a
[j1
][1][k
];
693 for (e
= 0; e
< excl
->nr
; e
++)
695 if (excl
->e
[e
] == k1
)
697 remove_excl(excl
, e
);
703 for (l
= 0; (l
< nnb
->nrexcl
[k1
][1]); l
++)
705 /* For all first neighbours of k1 */
706 l1
= nnb
->a
[k1
][1][l
];
708 for (e
= 0; e
< excl
->nr
; e
++)
710 if (excl
->e
[e
] == l1
)
712 remove_excl(excl
, e
);
724 void generate_excls(t_nextnb
*nnb
, int nrexcl
, t_excls excls
[])
729 for (N
= 1; (N
< std::min(nrexcl
, nnb
->nrex
)); N
++)
731 /* extract all i-j-k-l neighbours from nnb struct */
732 for (i
= 0; (i
< nnb
->nr
); i
++)
734 /* For all particles */
737 excl
->nr
+= nnb
->nrexcl
[i
][N
];
738 srenew(excl
->e
, excl
->nr
);
739 for (j
= 0; (j
< nnb
->nrexcl
[i
][N
]); j
++)
741 /* For all first neighbours */
742 if (nnb
->a
[i
][N
][j
] != i
)
744 excl
->e
[n
++] = nnb
->a
[i
][N
][j
];
751 /* Generate pairs, angles and dihedrals from .rtp settings */
752 void gen_pad(t_nextnb
*nnb
, t_atoms
*atoms
, t_restp rtp
[],
753 t_params plist
[], t_excls excls
[], t_hackblock hb
[],
754 gmx_bool bAllowMissing
)
756 t_param
*ang
, *dih
, *pai
, *improper
;
757 t_rbondeds
*hbang
, *hbdih
;
760 int res
, minres
, maxres
;
761 int i
, j
, j1
, k
, k1
, l
, l1
, m
, n
, i1
, i2
;
762 int ninc
, maxang
, maxdih
, maxpai
;
763 int nang
, ndih
, npai
, nimproper
, nbd
;
765 gmx_bool bFound
, bExcl
;
767 /* These are the angles, dihedrals and pairs that we generate
768 * from the bonds. The ones that are already there from the rtp file
775 maxang
= maxdih
= maxpai
= ninc
;
781 for (i
= 0; i
< 4; i
++)
788 gen_excls(atoms
, excls
, hb
, bAllowMissing
);
789 /* mark all entries as not matched yet */
790 for (i
= 0; i
< atoms
->nres
; i
++)
792 for (j
= 0; j
< ebtsNR
; j
++)
794 for (k
= 0; k
< hb
[i
].rb
[j
].nb
; k
++)
796 hb
[i
].rb
[j
].b
[k
].match
= FALSE
;
802 /* Extract all i-j-k-l neighbours from nnb struct to generate all
803 * angles and dihedrals. */
804 for (i
= 0; (i
< nnb
->nr
); i
++)
806 /* For all particles */
807 for (j
= 0; (j
< nnb
->nrexcl
[i
][1]); j
++)
809 /* For all first neighbours */
810 j1
= nnb
->a
[i
][1][j
];
811 for (k
= 0; (k
< nnb
->nrexcl
[j1
][1]); k
++)
813 /* For all first neighbours of j1 */
814 k1
= nnb
->a
[j1
][1][k
];
817 /* Generate every angle only once */
828 ang
[nang
].c0() = NOTSET
;
829 ang
[nang
].c1() = NOTSET
;
830 set_p_string(&(ang
[nang
]), "");
833 minres
= atoms
->atom
[ang
[nang
].a
[0]].resind
;
835 for (m
= 1; m
< 3; m
++)
837 minres
= std::min(minres
, atoms
->atom
[ang
[nang
].a
[m
]].resind
);
838 maxres
= std::max(maxres
, atoms
->atom
[ang
[nang
].a
[m
]].resind
);
840 res
= 2*minres
-maxres
;
843 res
+= maxres
-minres
;
844 get_atomnames_min(3, anm
, res
, atoms
, ang
[nang
].a
);
845 hbang
= &hb
[res
].rb
[ebtsANGLES
];
846 for (l
= 0; (l
< hbang
->nb
); l
++)
848 if (strcmp(anm
[1], hbang
->b
[l
].aj()) == 0)
851 for (m
= 0; m
< 3; m
+= 2)
854 ((strcmp(anm
[m
], hbang
->b
[l
].ai()) == 0) &&
855 (strcmp(anm
[2-m
], hbang
->b
[l
].ak()) == 0)));
859 set_p_string(&(ang
[nang
]), hbang
->b
[l
].s
);
860 /* Mark that we found a match for this entry */
861 hbang
->b
[l
].match
= TRUE
;
866 while (res
< maxres
);
870 /* Generate every dihedral, 1-4 exclusion and 1-4 interaction
874 for (l
= 0; (l
< nnb
->nrexcl
[k1
][1]); l
++)
876 /* For all first neighbours of k1 */
877 l1
= nnb
->a
[k1
][1][l
];
878 if ((l1
!= i
) && (l1
!= j1
))
889 for (m
= 0; m
< MAXFORCEPARAM
; m
++)
891 dih
[ndih
].c
[m
] = NOTSET
;
893 set_p_string(&(dih
[ndih
]), "");
897 minres
= atoms
->atom
[dih
[ndih
].a
[0]].resind
;
899 for (m
= 1; m
< 4; m
++)
901 minres
= std::min(minres
, atoms
->atom
[dih
[ndih
].a
[m
]].resind
);
902 maxres
= std::max(maxres
, atoms
->atom
[dih
[ndih
].a
[m
]].resind
);
904 res
= 2*minres
-maxres
;
907 res
+= maxres
-minres
;
908 get_atomnames_min(4, anm
, res
, atoms
, dih
[ndih
].a
);
909 hbdih
= &hb
[res
].rb
[ebtsPDIHS
];
910 for (n
= 0; (n
< hbdih
->nb
); n
++)
913 for (m
= 0; m
< 2; m
++)
916 ((strcmp(anm
[3*m
], hbdih
->b
[n
].ai()) == 0) &&
917 (strcmp(anm
[1+m
], hbdih
->b
[n
].aj()) == 0) &&
918 (strcmp(anm
[2-m
], hbdih
->b
[n
].ak()) == 0) &&
919 (strcmp(anm
[3-3*m
], hbdih
->b
[n
].al()) == 0)));
923 set_p_string(&dih
[ndih
], hbdih
->b
[n
].s
);
924 /* Mark that we found a match for this entry */
925 hbdih
->b
[n
].match
= TRUE
;
927 /* Set the last parameter to be able to see
928 if the dihedral was in the rtp list.
930 dih
[ndih
].c
[MAXFORCEPARAM
-1] = DIHEDRAL_WAS_SET_IN_RTP
;
933 /* Set the next direct in case the rtp contains
934 multiple entries for this dihedral.
945 for (m
= 0; m
< MAXFORCEPARAM
; m
++)
947 dih
[ndih
].c
[m
] = NOTSET
;
952 while (res
< maxres
);
965 for (m
= 0; m
< MAXFORCEPARAM
; m
++)
967 dih
[ndih
].c
[m
] = NOTSET
;
969 set_p_string(&(dih
[ndih
]), "");
973 nbd
= nb_dist(nnb
, i
, l1
);
976 fprintf(debug
, "Distance (%d-%d) = %d\n", i
+1, l1
+1, nbd
);
980 i1
= std::min(i
, l1
);
981 i2
= std::max(i
, l1
);
983 for (m
= 0; m
< excls
[i1
].nr
; m
++)
985 bExcl
= bExcl
|| excls
[i1
].e
[m
] == i2
;
989 if (rtp
[0].bGenerateHH14Interactions
||
990 !(is_hydro(atoms
, i1
) && is_hydro(atoms
, i2
)))
999 pai
[npai
].c0() = NOTSET
;
1000 pai
[npai
].c1() = NOTSET
;
1001 set_p_string(&(pai
[npai
]), "");
1016 /* The above approach is great in that we double-check that e.g. an angle
1017 * really corresponds to three atoms connected by bonds, but this is not
1018 * generally true. Go through the angle and dihedral hackblocks to add
1019 * entries that we have not yet marked as matched when going through bonds.
1021 for (i
= 0; i
< atoms
->nres
; i
++)
1023 /* Add remaining angles from hackblock */
1024 hbang
= &hb
[i
].rb
[ebtsANGLES
];
1025 for (j
= 0; j
< hbang
->nb
; j
++)
1027 if (hbang
->b
[j
].match
== TRUE
)
1029 /* We already used this entry, continue to the next */
1032 /* Hm - entry not used, let's see if we can find all atoms */
1036 srenew(ang
, maxang
);
1039 for (k
= 0; k
< 3 && bFound
; k
++)
1041 p
= hbang
->b
[j
].a
[k
];
1048 else if (p
[0] == '+')
1053 ang
[nang
].a
[k
] = search_res_atom(p
, res
, atoms
, "angle", TRUE
);
1054 bFound
= (ang
[nang
].a
[k
] != -1);
1056 ang
[nang
].c0() = NOTSET
;
1057 ang
[nang
].c1() = NOTSET
;
1061 set_p_string(&(ang
[nang
]), hbang
->b
[j
].s
);
1062 hbang
->b
[j
].match
= TRUE
;
1063 /* Incrementing nang means we save this angle */
1068 /* Add remaining dihedrals from hackblock */
1069 hbdih
= &hb
[i
].rb
[ebtsPDIHS
];
1070 for (j
= 0; j
< hbdih
->nb
; j
++)
1072 if (hbdih
->b
[j
].match
== TRUE
)
1074 /* We already used this entry, continue to the next */
1077 /* Hm - entry not used, let's see if we can find all atoms */
1081 srenew(dih
, maxdih
);
1084 for (k
= 0; k
< 4 && bFound
; k
++)
1086 p
= hbdih
->b
[j
].a
[k
];
1093 else if (p
[0] == '+')
1098 dih
[ndih
].a
[k
] = search_res_atom(p
, res
, atoms
, "dihedral", TRUE
);
1099 bFound
= (dih
[ndih
].a
[k
] != -1);
1101 for (m
= 0; m
< MAXFORCEPARAM
; m
++)
1103 dih
[ndih
].c
[m
] = NOTSET
;
1108 set_p_string(&(dih
[ndih
]), hbdih
->b
[j
].s
);
1109 hbdih
->b
[j
].match
= TRUE
;
1110 /* Incrementing ndih means we save this dihedral */
1117 /* Sort angles with respect to j-i-k (middle atom first) */
1120 qsort(ang
, nang
, (size_t)sizeof(ang
[0]), acomp
);
1123 /* Sort dihedrals with respect to j-k-i-l (middle atoms first) */
1126 qsort(dih
, ndih
, (size_t)sizeof(dih
[0]), dcomp
);
1129 /* Sort the pairs */
1132 qsort(pai
, npai
, (size_t)sizeof(pai
[0]), pcomp
);
1136 /* Remove doubles, could occur in 6-rings, such as phenyls,
1137 maybe one does not want this when fudgeQQ < 1.
1139 fprintf(stderr
, "Before cleaning: %d pairs\n", npai
);
1140 rm2par(pai
, &npai
, preq
);
1143 /* Get the impropers from the database */
1144 nimproper
= get_impropers(atoms
, hb
, &improper
, bAllowMissing
);
1146 /* Sort the impropers */
1147 sort_id(nimproper
, improper
);
1151 fprintf(stderr
, "Before cleaning: %d dihedrals\n", ndih
);
1152 clean_dih(dih
, &ndih
, improper
, nimproper
, atoms
,
1153 rtp
[0].bKeepAllGeneratedDihedrals
,
1154 rtp
[0].bRemoveDihedralIfWithImproper
);
1157 /* Now we have unique lists of angles and dihedrals
1158 * Copy them into the destination struct
1160 cppar(ang
, nang
, plist
, F_ANGLES
);
1161 cppar(dih
, ndih
, plist
, F_PDIHS
);
1162 cppar(improper
, nimproper
, plist
, F_IDIHS
);
1163 cppar(pai
, npai
, plist
, F_LJ14
);
1165 /* Remove all exclusions which are within nrexcl */
1166 clean_excls(nnb
, rtp
[0].nrexcl
, excls
);